Skip to content
GitLab
Projects
Groups
Snippets
Help
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Open sidebar
OpenWalnut
OpenWalnut Core
Commits
09bfdc24
Commit
09bfdc24
authored
May 22, 2013
by
Sebastian Eichelbaum
Browse files
[MERGE]
parents
09735ed9
0a80d9bb
Changes
12
Hide whitespace changes
Inline
Side-by-side
Showing
12 changed files
with
359 additions
and
75 deletions
+359
-75
src/core/common/WHistogram2D.cpp
src/core/common/WHistogram2D.cpp
+1
-1
src/core/common/WHistogram2D.h
src/core/common/WHistogram2D.h
+0
-1
src/core/graphicsEngine/WGETextureUtils.cpp
src/core/graphicsEngine/WGETextureUtils.cpp
+142
-0
src/core/graphicsEngine/WGETextureUtils.h
src/core/graphicsEngine/WGETextureUtils.h
+16
-0
src/core/graphicsEngine/WTriangleMesh.cpp
src/core/graphicsEngine/WTriangleMesh.cpp
+63
-29
src/core/graphicsEngine/WTriangleMesh.h
src/core/graphicsEngine/WTriangleMesh.h
+23
-9
src/modules/openIGTLink/CMakeLists.txt
src/modules/openIGTLink/CMakeLists.txt
+0
-2
tools/cmake/FindEigen3.cmake
tools/cmake/FindEigen3.cmake
+81
-0
tools/cmake/Findeigen3.cmake
tools/cmake/Findeigen3.cmake
+0
-31
tools/cmake/Findniftilib.cmake
tools/cmake/Findniftilib.cmake
+1
-1
tools/cmake/OpenWalnut.cmake
tools/cmake/OpenWalnut.cmake
+1
-1
tools/release/packaging/debian/copyright
tools/release/packaging/debian/copyright
+31
-0
No files found.
src/core/common/WHistogram2D.cpp
View file @
09bfdc24
...
...
@@ -149,6 +149,6 @@ WGETexture2D::RPtr WHistogram2D::getTexture()
data
[
i
+
j
*
imageWidth
]
=
static_cast
<
float
>
(
m_bins
(
i
,
j
)
)
/
maxCount
;
}
}
;
return
WGETexture2D
::
RPtr
(
new
WGETexture2D
(
image
)
);
}
src/core/common/WHistogram2D.h
View file @
09bfdc24
...
...
@@ -42,7 +42,6 @@
class
WHistogram2D
:
public
WHistogramND
<
2
>
// NOLINT
{
public:
/**
* Convenience type for a shared_ptr on this type.
*/
...
...
src/core/graphicsEngine/WGETextureUtils.cpp
View file @
09bfdc24
...
...
@@ -22,6 +22,11 @@
//
//---------------------------------------------------------------------------
#include <algorithm>
#include <vector>
#include <boost/random.hpp>
#include "../common/exceptions/WPreconditionNotMet.h"
#include "WGETexture.h"
...
...
@@ -114,3 +119,140 @@ osg::ref_ptr< osg::Image > wge::genWhiteNoiseImage( size_t sizeX, size_t sizeY,
return
randImage
;
}
osg
::
ref_ptr
<
WGETexture
<
osg
::
Texture3D
>
>
wge
::
genTuringNoiseTexture
(
std
::
size_t
sizeX
,
std
::
size_t
sizeY
,
std
::
size_t
sizeZ
,
std
::
size_t
channels
)
{
WPrecond
(
channels
==
1
||
channels
==
3
||
channels
==
4
,
"Invalid number of channels. Valid are: 1, 3 and 4."
);
// some constants, maybe change to parameters
float
const
spotIrregularity
=
0.05
f
;
// 0.0 - 1.0
std
::
size_t
const
iterations
=
200
;
float
const
spotSize
=
0.5
;
float
const
spotFactor
=
(
0.02
f
+
0.58
f
*
(
1.0
f
-
spotSize
)
)
/
15.0
f
;
float
const
d1
=
0.125
f
;
float
const
d2
=
0.03125
f
;
float
const
speed
=
1.0
f
;
osg
::
ref_ptr
<
osg
::
Image
>
img
=
new
osg
::
Image
;
GLenum
type
=
GL_LUMINANCE
;
if
(
channels
==
3
)
{
type
=
GL_RGB
;
}
else
if
(
channels
==
4
)
{
type
=
GL_RGBA
;
}
std
::
vector
<
float
>
concentration1
(
sizeX
*
sizeY
*
sizeZ
,
4.0
);
std
::
vector
<
float
>
concentration2
(
sizeX
*
sizeY
*
sizeZ
,
4.0
);
std
::
vector
<
float
>
delta1
(
sizeX
*
sizeY
*
sizeZ
,
0.0
);
std
::
vector
<
float
>
delta2
(
sizeX
*
sizeY
*
sizeZ
,
0.0
);
std
::
vector
<
float
>
noise
(
sizeX
*
sizeY
*
sizeZ
);
boost
::
mt19937
generator
(
std
::
time
(
0
)
);
float
noiseRange
=
0.1
f
+
4.9
f
*
spotIrregularity
;
boost
::
uniform_real
<
float
>
dist
(
12.0
-
noiseRange
,
12.0
+
noiseRange
);
boost
::
variate_generator
<
boost
::
mt19937
&
,
boost
::
uniform_real
<
float
>
>
rand
(
generator
,
dist
);
// initialization step
for
(
std
::
size_t
i
=
0
;
i
<
sizeX
;
++
i
)
{
for
(
std
::
size_t
j
=
0
;
j
<
sizeY
;
++
j
)
{
for
(
std
::
size_t
k
=
0
;
k
<
sizeZ
;
++
k
)
{
std
::
size_t
idx
=
i
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
;
noise
[
idx
]
=
rand
();
}
}
}
// iteration
for
(
std
::
size_t
iter
=
0
;
iter
<
iterations
;
++
iter
)
{
std
::
cout
<<
"iterations: "
<<
iter
<<
std
::
endl
;
for
(
std
::
size_t
i
=
0
;
i
<
sizeX
;
++
i
)
{
std
::
size_t
iNext
=
(
i
+
1
)
%
sizeX
;
std
::
size_t
iPrev
=
(
i
+
sizeX
-
1
)
%
sizeX
;
for
(
std
::
size_t
j
=
0
;
j
<
sizeY
;
++
j
)
{
std
::
size_t
jNext
=
(
j
+
1
)
%
sizeY
;
std
::
size_t
jPrev
=
(
j
+
sizeY
-
1
)
%
sizeY
;
for
(
std
::
size_t
k
=
0
;
k
<
sizeZ
;
++
k
)
{
std
::
size_t
kNext
=
(
k
+
1
)
%
sizeZ
;
std
::
size_t
kPrev
=
(
k
+
sizeZ
-
1
)
%
sizeZ
;
std
::
size_t
idx
=
i
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
;
// estimate change in concentrations
// we use a 3d laplace filter here instead of the 2d filter as in Eichelbaum et al.
float
dc1
=
0.0
;
dc1
+=
concentration1
[
iPrev
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc1
+=
concentration1
[
iNext
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc1
+=
concentration1
[
i
+
jPrev
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc1
+=
concentration1
[
i
+
jNext
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc1
+=
concentration1
[
i
+
j
*
sizeX
+
kPrev
*
sizeX
*
sizeY
];
dc1
+=
concentration1
[
i
+
j
*
sizeX
+
kNext
*
sizeX
*
sizeY
];
dc1
-=
6.0
f
*
concentration1
[
idx
];
float
dc2
=
0.0
;
dc2
+=
concentration2
[
iPrev
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc2
+=
concentration2
[
iNext
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc2
+=
concentration2
[
i
+
jPrev
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc2
+=
concentration2
[
i
+
jNext
*
sizeX
+
k
*
sizeX
*
sizeY
];
dc2
+=
concentration2
[
i
+
j
*
sizeX
+
kPrev
*
sizeX
*
sizeY
];
dc2
+=
concentration2
[
i
+
j
*
sizeX
+
kNext
*
sizeX
*
sizeY
];
dc2
-=
6.0
f
*
concentration2
[
idx
];
// diffusion
delta1
[
idx
]
=
spotFactor
*
(
16.0
f
-
concentration1
[
idx
]
*
concentration2
[
idx
]
)
+
d1
*
dc1
;
delta2
[
idx
]
=
spotFactor
*
(
concentration1
[
idx
]
*
concentration2
[
idx
]
-
concentration2
[
idx
]
-
noise
[
idx
]
)
+
d2
*
dc2
;
}
}
}
for
(
std
::
size_t
i
=
0
;
i
<
sizeX
;
++
i
)
{
for
(
std
::
size_t
j
=
0
;
j
<
sizeY
;
++
j
)
{
for
(
std
::
size_t
k
=
0
;
k
<
sizeZ
;
++
k
)
{
std
::
size_t
idx
=
i
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
;
concentration1
[
idx
]
+=
speed
*
delta1
[
idx
];
concentration2
[
idx
]
+=
speed
*
delta2
[
idx
];
}
}
}
}
img
->
allocateImage
(
sizeX
,
sizeY
,
sizeZ
,
type
,
GL_UNSIGNED_BYTE
);
// find min and max
float
c1min
=
*
std
::
min_element
(
concentration1
.
begin
(),
concentration1
.
end
()
);
float
c1max
=
*
std
::
max_element
(
concentration1
.
begin
(),
concentration1
.
end
()
);
// copy to image
for
(
std
::
size_t
i
=
0
;
i
<
sizeX
;
++
i
)
{
for
(
std
::
size_t
j
=
0
;
j
<
sizeY
;
++
j
)
{
for
(
std
::
size_t
k
=
0
;
k
<
sizeZ
;
++
k
)
{
std
::
size_t
idx
=
i
+
j
*
sizeX
+
k
*
sizeX
*
sizeY
;
img
->
data
()[
idx
]
=
255.0
f
*
(
concentration1
[
idx
]
-
c1min
)
/
(
c1max
-
c1min
);
}
}
}
return
osg
::
ref_ptr
<
WGETexture
<
osg
::
Texture3D
>
>
(
new
WGETexture
<
osg
::
Texture3D
>
(
img
)
);
}
src/core/graphicsEngine/WGETextureUtils.h
View file @
09bfdc24
...
...
@@ -145,6 +145,22 @@ namespace wge
* \return the generated image.
*/
osg
::
ref_ptr
<
osg
::
Image
>
genWhiteNoiseImage
(
size_t
sizeX
,
size_t
sizeY
,
size_t
sizeZ
,
size_t
channels
=
1
);
/**
* Generates a 3D turing noise texture. For details, see:
*
* Eichelbaum, Sebastian, et al. "Fabric-like visualization of tensor field data on arbitrary surfaces
* in image space." New Developments in the Visualization and Processing of Tensor Fields. Springer Berlin Heidelberg, 2012. 71-92.
*
* \param sizeX The size of the textures in voxels, should be a power of 2.
* \param sizeY The size of the textures in voxels, should be a power of 2.
* \param sizeZ The size of the textures in voxels, should be a power of 2.
* \param channels The number of channels; either 1, 3 or 4.
*
* \return The texture.
*/
osg
::
ref_ptr
<
WGETexture
<
osg
::
Texture3D
>
>
genTuringNoiseTexture
(
std
::
size_t
sizeX
,
std
::
size_t
sizeY
,
std
::
size_t
sizeZ
,
std
::
size_t
channels
=
1
);
}
template
<
typename
T
>
...
...
src/core/graphicsEngine/WTriangleMesh.cpp
View file @
09bfdc24
...
...
@@ -31,10 +31,15 @@
#include <string>
#include <utility>
#include <vector>
#include <limits>
#include <osg/io_utils>
#include <Eigen/Sparse>
#include <Eigen/Eigenvalues>
#include "../common/WAssert.h"
#include "../common/WLogger.h"
#include "../common/math/WMath.h"
#include "../common/datastructures/WUnionFind.h"
#include "WTriangleMesh.h"
...
...
@@ -929,6 +934,8 @@ void WTriangleMesh::estimateCurvature()
n
+=
m_triangleNormals
->
operator
[]
(
triId
)
*
w
;
}
WAssert
(
n
.
length
()
>
0.0001
,
"Invalid normal!"
);
n
.
normalize
();
normals
[
vtxId
]
=
n
;
}
...
...
@@ -958,10 +965,9 @@ void WTriangleMesh::estimateCurvature()
}
}
WAssert
(
neighbors
.
size
()
>
2
,
""
);
WAssert
(
neighbors
.
size
()
>
2
,
"
Vertex has too few neighbors! Does this mesh have holes?
"
);
// choose tangent space coordinate system
double
maxCurvature
=
0.0
;
double
maxCurvature
=
-
std
::
numeric_limits
<
double
>::
infinity
();
std
::
vector
<
double
>
curvatures
;
osg
::
Vec3
maxCurvatureTangent
(
0.0
,
0.0
,
0.0
);
...
...
@@ -973,16 +979,13 @@ void WTriangleMesh::estimateCurvature()
osg
::
Vec3
const
&
neighbPos
=
m_verts
->
operator
[]
(
*
it
);
osg
::
Vec3
const
&
neighbNormal
=
normals
[
*
it
];
//
calc
projecte
d
tangent
// project
( neighbPos - p ) onto th
e tangent
plane
osg
::
Vec3
tangent
=
(
neighbPos
-
p
)
-
normal
*
(
(
neighbPos
-
p
)
*
normal
);
tangent
.
normalize
();
// approximate normal curvature in tangent direction
double
ncurv
=
-
1.0
*
(
(
neighbPos
-
p
)
*
(
neighbNormal
-
normal
)
)
/
(
(
neighbPos
-
p
)
*
(
neighbPos
-
p
)
);
// WAssert( 0.0 <= ncurv, "" );
// WAssert( ncurv <= 1.0, "" );
if
(
ncurv
>
maxCurvature
)
{
maxCurvature
=
ncurv
;
...
...
@@ -993,7 +996,9 @@ void WTriangleMesh::estimateCurvature()
curvatures
.
push_back
(
ncurv
);
}
// part 2: calculate coordinate system
WAssert
(
maxCurvatureTangent
.
length
()
>
0.0001
,
"Invalid tangent length!"
);
// part 2: choose a coordinate system in the tangent plane
osg
::
Vec3
const
&
e1
=
maxCurvatureTangent
;
osg
::
Vec3
e2
=
e1
^
normal
;
...
...
@@ -1009,9 +1014,10 @@ void WTriangleMesh::estimateCurvature()
}
}
// curvatures were all almost zero
if
(
!
significantCurvature
)
{
// curvatures were all almost zero
// the mesh is flat at this point, write values accordingly
m_mainNormalCurvature
->
operator
[]
(
vtxId
)
=
0.0
;
m_mainCurvaturePrincipalDirection
->
operator
[]
(
vtxId
)
=
e1
;
m_secondaryNormalCurvature
->
operator
[]
(
vtxId
)
=
0.0
;
...
...
@@ -1020,32 +1026,48 @@ void WTriangleMesh::estimateCurvature()
continue
;
}
// calculate coefficients
double
a11
=
0.0
;
double
a12
=
0.0
;
// a21 == a12
double
a22
=
0.0
;
double
a13
=
0.0
;
double
a23
=
0.0
;
// calculate coefficients of the ellipse a * cos²(theta) + b * cos(theta) * sin(theta) * c * sin²(theta)
// this is done by estimating a as the largest curvature amoung the estimated curvatures in all tangent
// direction belonging to points that share a triangle with the current one
double
const
&
a
=
maxCurvature
;
Eigen
::
Matrix
<
double
,
-
1
,
-
1
>
X
(
tangents
.
size
(),
2
);
Eigen
::
Matrix
<
double
,
-
1
,
-
1
>
y
(
tangents
.
size
(),
1
);
for
(
std
::
size_t
k
=
0
;
k
<
tangents
.
size
();
++
k
)
{
double
theta
=
calcAngleBetweenNormalizedVectors
(
tangents
[
k
],
e1
);
double
sintheta
=
sin
(
theta
);
double
costheta
=
cos
(
theta
);
X
(
k
,
0
)
=
cos
(
theta
)
*
sin
(
theta
);
X
(
k
,
1
)
=
sin
(
theta
)
*
sin
(
theta
);
y
(
k
,
0
)
=
curvatures
[
k
]
-
a
*
cos
(
theta
)
*
cos
(
theta
);
}
// use LU decomposition to calculate rank
Eigen
::
FullPivLU
<
Eigen
::
Matrix
<
double
,
-
1
,
-
1
>
>
lu
(
X
);
// we need a rank of at least 2 to calculate the coeffs
if
(
lu
.
rank
()
<
2
)
{
wlog
::
warn
(
"WTriangleMesh::estimateCurvature"
)
<<
"Rank too low, cannot estimate curvature for this vertex!"
;
a11
+=
costheta
*
costheta
*
sintheta
*
sintheta
;
a12
+=
costheta
*
sintheta
*
sintheta
*
sintheta
;
a22
+=
sintheta
*
sintheta
*
sintheta
*
sintheta
;
a13
+=
(
curvatures
[
k
]
-
a
*
costheta
*
costheta
)
*
costheta
*
sintheta
;
a23
+=
(
curvatures
[
k
]
-
a
*
costheta
*
costheta
)
*
sintheta
*
sintheta
;
m_mainNormalCurvature
->
operator
[]
(
vtxId
)
=
a
;
m_mainCurvaturePrincipalDirection
->
operator
[]
(
vtxId
)
=
e1
;
m_secondaryNormalCurvature
->
operator
[]
(
vtxId
)
=
0.0
;
m_secondaryCurvaturePrincipalDirection
->
operator
[]
(
vtxId
)
=
e2
;
continue
;
}
double
b
=
(
a13
*
a22
-
a23
*
a12
)
/
(
a11
*
a22
-
a12
*
a12
);
double
c
=
(
a11
*
a23
-
a12
*
a13
)
/
(
a11
*
a22
-
a12
*
a12
);
// do least squares
Eigen
::
Matrix
<
double
,
-
1
,
-
1
>
bb
=
(
X
.
transpose
()
*
X
).
inverse
()
*
X
.
transpose
()
*
y
;
// the missing coeffs b and c are now:
double
b
=
bb
(
0
,
0
);
double
c
=
bb
(
1
,
0
);
// this calculates the maximum and minimum normal curvatures
// (as eigenvalues)
double
Kg
=
a
*
c
-
0.25
*
b
*
b
;
double
H
=
0.5
*
(
a
+
c
);
double
s
=
H
*
H
-
Kg
;
...
...
@@ -1060,11 +1082,13 @@ void WTriangleMesh::estimateCurvature()
if
(
fabs
(
k1
-
k2
)
<
0.000001
)
{
// if the curvatures are equal, there is no single principal direction
m_mainNormalCurvature
->
operator
[]
(
vtxId
)
=
a
;
m_mainCurvaturePrincipalDirection
->
operator
[]
(
vtxId
)
=
e1
;
}
else
{
// if the curvatures differ, we can now find the direction of maximum curvature
double
temp
=
b
/
(
k2
-
k1
);
if
(
temp
>
1.0
)
...
...
@@ -1128,22 +1152,32 @@ double WTriangleMesh::calcAngleBetweenNormalizedVectors( osg::Vec3 const& v1, os
return
acos
(
temp
);
}
double
WTriangleMesh
::
getMainCurvature
(
std
::
size_t
vtxId
)
double
WTriangleMesh
::
getMainCurvature
(
std
::
size_t
vtxId
)
const
{
return
m_mainNormalCurvature
->
operator
[]
(
vtxId
);
}
double
WTriangleMesh
::
getSecondaryCurvature
(
std
::
size_t
vtxId
)
double
WTriangleMesh
::
getSecondaryCurvature
(
std
::
size_t
vtxId
)
const
{
return
m_secondaryNormalCurvature
->
operator
[]
(
vtxId
);
}
osg
::
Vec3
WTriangleMesh
::
getCurvatureMainPrincipalDirection
(
std
::
size_t
vtxId
)
boost
::
shared_ptr
<
std
::
vector
<
float
>
>
const
&
WTriangleMesh
::
getMainCurvatures
()
const
{
return
m_mainNormalCurvature
;
}
boost
::
shared_ptr
<
std
::
vector
<
float
>
>
const
&
WTriangleMesh
::
getSecondaryCurvatures
()
const
{
return
m_secondaryNormalCurvature
;
}
osg
::
Vec3
WTriangleMesh
::
getCurvatureMainPrincipalDirection
(
std
::
size_t
vtxId
)
const
{
return
m_mainCurvaturePrincipalDirection
->
operator
[]
(
vtxId
);
}
osg
::
Vec3
WTriangleMesh
::
getCurvatureSecondaryPrincipalDirection
(
std
::
size_t
vtxId
)
osg
::
Vec3
WTriangleMesh
::
getCurvatureSecondaryPrincipalDirection
(
std
::
size_t
vtxId
)
const
{
return
m_secondaryCurvaturePrincipalDirection
->
operator
[]
(
vtxId
);
}
...
...
src/core/graphicsEngine/WTriangleMesh.h
View file @
09bfdc24
...
...
@@ -428,7 +428,7 @@ public:
*
* \return The estimated main normal curvature.
*/
double
getMainCurvature
(
std
::
size_t
vtxId
);
double
getMainCurvature
(
std
::
size_t
vtxId
)
const
;
/**
* Retreive the secondary (minimum) curvature of a vertex.
...
...
@@ -437,7 +437,21 @@ public:
*
* \return The estimated secondary normal curvature.
*/
double
getSecondaryCurvature
(
std
::
size_t
vtxId
);
double
getSecondaryCurvature
(
std
::
size_t
vtxId
)
const
;
/**
* Get the vector of main curvature values.
*
* \return The curvature values for all the vertices.
*/
boost
::
shared_ptr
<
std
::
vector
<
float
>
>
const
&
getMainCurvatures
()
const
;
/**
* Get the vector of secondary curvature values.
*
* \return The curvature values for all the vertices.
*/
boost
::
shared_ptr
<
std
::
vector
<
float
>
>
const
&
getSecondaryCurvatures
()
const
;
/**
* Retreive the 3d principal direction of curvature of a vertex.
...
...
@@ -446,7 +460,7 @@ public:
*
* \return The first principal duirection.
*/
osg
::
Vec3
getCurvatureMainPrincipalDirection
(
std
::
size_t
vtxId
);
osg
::
Vec3
getCurvatureMainPrincipalDirection
(
std
::
size_t
vtxId
)
const
;
/**
* Retreive the 3d principal direction of curvature for the minimum normal curvature of a vertex.
...
...
@@ -455,7 +469,7 @@ public:
*
* \return The second principal duirection.
*/
osg
::
Vec3
getCurvatureSecondaryPrincipalDirection
(
std
::
size_t
vtxId
);
osg
::
Vec3
getCurvatureSecondaryPrincipalDirection
(
std
::
size_t
vtxId
)
const
;
/**
* Retreive the array of principal directions e.g. for use as texture coords.
...
...
@@ -471,6 +485,11 @@ public:
*/
osg
::
ref_ptr
<
osg
::
Vec3Array
>
getSecondaryPrincipalCurvatureDirectionArray
();
/**
* recalculates the vertex normals
*/
void
recalcVertNormals
();
/**
* Set this to true to force automatic normal calculation. Set it to false if you specify your own normals.
*
...
...
@@ -501,11 +520,6 @@ private:
*/
void
removeTriangle
(
size_t
index
);
/**
* recalculates the vertex normals
*/
void
recalcVertNormals
();
/**
* calculates a normal from the 3 points in space defining a triangle
*
...
...
src/modules/openIGTLink/CMakeLists.txt
View file @
09bfdc24
...
...
@@ -55,6 +55,4 @@ IF( OpenIGTLink_FOUND )
${
OpenIGTLink_LIBRARIES
}
# does your module need additional libs to compile?
""
# do you want to exclude files from stylechecking? (externals for example)
)
ELSE
(
OpenIGTLink_FOUND
)
MESSAGE
(
WARNING
"OpenIGTLink library wasn't found, don't compile openIGTLink module."
)
ENDIF
(
OpenIGTLink_FOUND
)
tools/cmake/FindEigen3.cmake
0 → 100644
View file @
09bfdc24
# - Try to find Eigen3 lib
#
# This module supports requiring a minimum version, e.g. you can do
# find_package(Eigen3 3.1.2)
# to require version 3.1.2 or newer of Eigen3.
#
# Once done this will define
#
# EIGEN3_FOUND - system has eigen lib with correct version
# EIGEN3_INCLUDE_DIR - the eigen include directory
# EIGEN3_VERSION - eigen version
# Copyright (c) 2006, 2007 Montel Laurent, <montel@kde.org>
# Copyright (c) 2008, 2009 Gael Guennebaud, <g.gael@free.fr>
# Copyright (c) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
# Redistribution and use is allowed according to the terms of the 2-clause BSD license.
if
(
NOT Eigen3_FIND_VERSION
)
if
(
NOT Eigen3_FIND_VERSION_MAJOR
)
set
(
Eigen3_FIND_VERSION_MAJOR 2
)
endif
(
NOT Eigen3_FIND_VERSION_MAJOR
)
if
(
NOT Eigen3_FIND_VERSION_MINOR
)
set
(
Eigen3_FIND_VERSION_MINOR 91
)
endif
(
NOT Eigen3_FIND_VERSION_MINOR
)
if
(
NOT Eigen3_FIND_VERSION_PATCH
)
set
(
Eigen3_FIND_VERSION_PATCH 0
)
endif
(
NOT Eigen3_FIND_VERSION_PATCH
)
set
(
Eigen3_FIND_VERSION
"
${
Eigen3_FIND_VERSION_MAJOR
}
.
${
Eigen3_FIND_VERSION_MINOR
}
.
${
Eigen3_FIND_VERSION_PATCH
}
"
)
endif
(
NOT Eigen3_FIND_VERSION
)
macro
(
_eigen3_check_version
)
file
(
READ
"
${
EIGEN3_INCLUDE_DIR
}
/Eigen/src/Core/util/Macros.h"
_eigen3_version_header
)
string
(
REGEX MATCH
"define[
\t
]+EIGEN_WORLD_VERSION[
\t
]+([0-9]+)"
_eigen3_world_version_match
"
${
_eigen3_version_header
}
"
)
set
(
EIGEN3_WORLD_VERSION
"
${
CMAKE_MATCH_1
}
"
)
string
(
REGEX MATCH
"define[
\t
]+EIGEN_MAJOR_VERSION[
\t
]+([0-9]+)"
_eigen3_major_version_match
"
${
_eigen3_version_header
}
"
)
set
(
EIGEN3_MAJOR_VERSION
"
${
CMAKE_MATCH_1
}
"
)
string
(
REGEX MATCH
"define[
\t
]+EIGEN_MINOR_VERSION[
\t
]+([0-9]+)"
_eigen3_minor_version_match
"
${
_eigen3_version_header
}
"
)
set
(
EIGEN3_MINOR_VERSION
"
${
CMAKE_MATCH_1
}
"
)
set
(
EIGEN3_VERSION
${
EIGEN3_WORLD_VERSION
}
.
${
EIGEN3_MAJOR_VERSION
}
.
${
EIGEN3_MINOR_VERSION
}
)
if
(
${
EIGEN3_VERSION
}
VERSION_LESS
${
Eigen3_FIND_VERSION
}
)
set
(
EIGEN3_VERSION_OK FALSE
)
else
(
${
EIGEN3_VERSION
}
VERSION_LESS
${
Eigen3_FIND_VERSION
}
)
set
(
EIGEN3_VERSION_OK TRUE
)
endif
(
${
EIGEN3_VERSION
}
VERSION_LESS
${
Eigen3_FIND_VERSION
}
)
if
(
NOT EIGEN3_VERSION_OK
)
message
(
STATUS
"Eigen3 version
${
EIGEN3_VERSION
}
found in
${
EIGEN3_INCLUDE_DIR
}
, "
"but at least version
${
Eigen3_FIND_VERSION
}
is required"
)
endif
(
NOT EIGEN3_VERSION_OK
)
endmacro
(
_eigen3_check_version
)
if
(
EIGEN3_INCLUDE_DIR
)
# in cache already
_eigen3_check_version
()
set
(
EIGEN3_FOUND
${
EIGEN3_VERSION_OK
}
)
else
(
EIGEN3_INCLUDE_DIR
)
find_path
(
EIGEN3_INCLUDE_DIR NAMES signature_of_eigen3_matrix_library
PATHS
${
CMAKE_INSTALL_PREFIX
}
/include
${
KDE4_INCLUDE_DIR
}
PATH_SUFFIXES eigen3 eigen
)
if
(
EIGEN3_INCLUDE_DIR
)
_eigen3_check_version
()
endif
(
EIGEN3_INCLUDE_DIR
)
include
(
FindPackageHandleStandardArgs
)
find_package_handle_standard_args
(
Eigen3 DEFAULT_MSG EIGEN3_INCLUDE_DIR EIGEN3_VERSION_OK
)
mark_as_advanced
(
EIGEN3_INCLUDE_DIR
)
endif
(
EIGEN3_INCLUDE_DIR
)
tools/cmake/Findeigen3.cmake
deleted
100644 → 0
View file @
09735ed9
# This script searches eigen3. See http://eigen.tuxfamily.org
#
# The following variables will be filled:
# * EIGEN3_FOUND - if Eigen (header) was found
# * EIGEN3_INCLUDE_DIR - the path of Eigne header if found
#
# You can set the environment variable "EIGEN3_INCLUDE_DIR" to help this script finding it if
# you placed it in some unusual location.
FIND_PATH
(
EIGEN3_INCLUDE_DIR Eigen/Core $ENV{EIGEN3_INCLUDE_DIR}
$ENV{HOME}/.local/include/eigen3
/usr/include/eigen3
/usr/local/include/eigen3
/opt/local/include/eigen3
/sw/include/eigen3
)
SET
(
EIGEN3_FOUND FALSE
)
IF
(
EIGEN3_INCLUDE_DIR
)
SET
(
EIGEN3_FOUND TRUE
)
ENDIF
()
IF
(
EIGEN3_FOUND
)
IF
(
NOT eigen3_FIND_QUIETLY
)
MESSAGE
(
STATUS
"Found eigen3:
${
EIGEN3_INCLUDE_DIR
}
"
)
ENDIF
()
ELSE
()
IF
(
eigen3_FIND_REQUIRED
)
MESSAGE
(
FATAL_ERROR
"Could not find eigen3. Check your distribution or download from http://eigen.tuxfamily.org"
)
ENDIF
()
ENDIF
()
tools/cmake/Findniftilib.cmake
View file @
09bfdc24
...
...
@@ -6,7 +6,7 @@
# * NIFTILIB_LIBRARY - the path to the library
#
FIND_PATH
(
NIFTILIB_INCLUDE_DIR nifti1.h
/usr/include/nifti /usr/local/include/nifti /opt/local/include/
nifti
)
FIND_PATH
(
NIFTILIB_INCLUDE_DIR nifti1.h
PATH_SUFFIXES
nifti
)
# This hack is inspired by FindBoost.cmake. It ensures that only shared objects are found. Linking a SO with a static lib is not possible
# in Linux. On other systems, this should be no problem.
...
...
tools/cmake/OpenWalnut.cmake
View file @
09bfdc24
...
...
@@ -410,7 +410,7 @@ MARK_AS_ADVANCED( FORCE OPENTHREADS_LIBRARY_DEBUG )
# Eigen3, at least 3.0.0
# See http://eigen.tuxfamily.org
FIND_PACKAGE
(
e
igen3 REQUIRED
)
FIND_PACKAGE
(
E
igen3 REQUIRED
)
IF
(
EIGEN3_FOUND
)
INCLUDE_DIRECTORIES
(
SYSTEM
${
EIGEN3_INCLUDE_DIR
}
)
...
...