Skip to content

Commit ff7fb80

Browse files
authored
Merge pull request #46 from DVigneault/UseDirectionMatrix
BUG: Use Direction When Computing Offsets
2 parents b4f2ae5 + e6bb029 commit ff7fb80

File tree

3 files changed

+128
-4
lines changed

3 files changed

+128
-4
lines changed

include/itkCuberilleImageToMeshFilter.hxx

+4-4
Original file line numberDiff line numberDiff line change
@@ -265,10 +265,10 @@ CuberilleImageToMeshFilter<TInputImage,TOutputMesh,TInterpolator>
265265
{
266266
PointType vertex;
267267
image->TransformIndexToPhysicalPoint( index, vertex );
268-
SpacingType spacing = image->GetSpacing();
269-
vertex[0] -= ( spacing[0] / 2.0 );
270-
vertex[1] -= ( spacing[1] / 2.0 );
271-
vertex[2] -= ( spacing[2] / 2.0 );
268+
const auto spacing = image->GetSpacing();
269+
const auto direction = image->GetDirection();
270+
const auto offset = direction * spacing * 0.5;
271+
vertex -= offset;
272272
if( m_ProjectVerticesToIsoSurface )
273273
{
274274
ProjectVertexToIsoSurface( vertex );

test/CMakeLists.txt

+4
Original file line numberDiff line numberDiff line change
@@ -4,6 +4,7 @@ set(DATA_DIR ${CMAKE_CURRENT_SOURCE_DIR}/Input)
44

55
set(CuberilleTests
66
CuberilleTest01.cxx
7+
CuberilleTest02.cxx
78
)
89

910
CreateTestDriver(Cuberille "${Cuberille-Test_LIBRARIES}" "${CuberilleTests}")
@@ -347,3 +348,6 @@ add_test(
347348
# 0.85 # Step length relaxation factor
348349
# 25 # Maximum number of steps
349350
#)
351+
352+
itk_add_test(NAME CuberilleTestDirectionMatrix
353+
COMMAND ${itk-module}TestDriver CuberilleTest02 )

test/CuberilleTest02.cxx

+120
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,120 @@
1+
/*=========================================================================
2+
*
3+
* Copyright Insight Software Consortium
4+
*
5+
* Licensed under the Apache License, Version 2.0 (the "License");
6+
* you may not use this file except in compliance with the License.
7+
* You may obtain a copy of the License at
8+
*
9+
* http://www.apache.org/licenses/LICENSE-2.0.txt
10+
*
11+
* Unless required by applicable law or agreed to in writing, software
12+
* distributed under the License is distributed on an "AS IS" BASIS,
13+
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14+
* See the License for the specific language governing permissions and
15+
* limitations under the License.
16+
*
17+
*=========================================================================*/
18+
19+
// STD
20+
#include <array>
21+
22+
// ITK
23+
#include "itkTestingMacros.h"
24+
#include "itkImage.h"
25+
#include "itkMesh.h"
26+
#include "itkCuberilleImageToMeshFilter.h"
27+
#include "itkNearestNeighborInterpolateImageFunction.h"
28+
29+
const unsigned int Dimension = 3;
30+
using TPixel = unsigned char;
31+
using TImage = itk::Image<TPixel, Dimension>;
32+
using TMesh = itk::Mesh<double, Dimension>;
33+
using TImageToMesh = itk::CuberilleImageToMeshFilter< TImage, TMesh >;
34+
using TInterp = itk::NearestNeighborInterpolateImageFunction< TImage >;
35+
36+
TImage::Pointer
37+
CuberilleTestCreateImage(const std::array<bool, 3> &flips) {
38+
const auto image = TImage::New();
39+
40+
TImage::IndexType start;
41+
start.Fill( 0 );
42+
43+
TImage::SizeType size;
44+
size.Fill( 3 );
45+
46+
TImage::RegionType region;
47+
region.SetSize(size);
48+
region.SetIndex(start);
49+
50+
image->SetRegions(region);
51+
image->Allocate();
52+
image->FillBuffer( 0 );
53+
54+
TImage::IndexType index;
55+
index.Fill( 1 );
56+
image->SetPixel( index, 1 );
57+
58+
TImage::DirectionType direction;
59+
direction.SetIdentity();
60+
61+
for (size_t i = 0; i < flips.size(); ++i) {
62+
if (flips[i]) {
63+
direction(i, i) *= -1;
64+
}
65+
}
66+
67+
image->SetDirection( direction );
68+
69+
return image;
70+
}
71+
72+
void CuberilleTestHelper(TImage::Pointer input) {
73+
74+
const auto image_to_mesh = TImageToMesh::New();
75+
image_to_mesh->SetInput( input );
76+
image_to_mesh->ProjectVerticesToIsoSurfaceOff();
77+
image_to_mesh->Update();
78+
79+
const auto mesh = TMesh::New();
80+
mesh->Graft( image_to_mesh->GetOutput() );
81+
mesh->DisconnectPipeline();
82+
83+
TMesh::PointType center;
84+
center.Fill( 0.0 );
85+
86+
for (auto it = mesh->GetPoints()->Begin();
87+
it != mesh->GetPoints()->End();
88+
++it) {
89+
center[0] += it.Value()[0];
90+
center[1] += it.Value()[1];
91+
center[2] += it.Value()[2];
92+
}
93+
94+
center[0] /= 8.0;
95+
center[1] /= 8.0;
96+
center[2] /= 8.0;
97+
98+
const auto interp = TInterp::New();
99+
interp->SetInputImage( input );
100+
101+
if (1 != interp->Evaluate( center )) {
102+
throw 0;
103+
}
104+
105+
}
106+
107+
int CuberilleTest02 (int itkNotUsed(argc), char * itkNotUsed(argv) [] ) {
108+
109+
CuberilleTestHelper( CuberilleTestCreateImage({{false, false, false}}) );
110+
CuberilleTestHelper( CuberilleTestCreateImage({{false, false, true}}) );
111+
CuberilleTestHelper( CuberilleTestCreateImage({{false, true, false}}) );
112+
CuberilleTestHelper( CuberilleTestCreateImage({{false, true, true}}) );
113+
CuberilleTestHelper( CuberilleTestCreateImage({{true, false, false}}) );
114+
CuberilleTestHelper( CuberilleTestCreateImage({{true, false, true}}) );
115+
CuberilleTestHelper( CuberilleTestCreateImage({{true, true, false}}) );
116+
CuberilleTestHelper( CuberilleTestCreateImage({{true, true, true}}) );
117+
118+
return EXIT_SUCCESS;
119+
120+
}

0 commit comments

Comments
 (0)