]> Creatis software - cpPlugins.git/blob - plugins/ImageMeshFilters/RasterMeshFilter.cxx
20a6769299597c2443dbbfb8ff02269c8d998b72
[cpPlugins.git] / plugins / ImageMeshFilters / RasterMeshFilter.cxx
1 #include <ImageMeshFilters/RasterMeshFilter.h>
2 #include <cpPlugins/DataObjects/BoundingBox.h>
3 #include <cpPlugins/DataObjects/Mesh.h>
4
5 #include <itkTriangleMeshToBinaryImageFilter.h>
6 #include <cpExtensions/Algorithms/RasterContourFilter.h>
7
8 // -------------------------------------------------------------------------
9 cpPluginsImageMeshFilters::RasterMeshFilter::
10 RasterMeshFilter( )
11   : Superclass( )
12 {
13   typedef cpPlugins::BaseObjects::DataObject _TDataObject;
14   typedef cpPlugins::DataObjects::Image      _TImage;
15   typedef cpPlugins::DataObjects::Mesh       _TMesh;
16
17   this->_ConfigureInput< _TMesh >( "Input", true, false );
18   this->_ConfigureInput< _TDataObject >( "Template", false, false );
19   this->_ConfigureOutput< _TImage >( "Output" );
20
21   this->m_Parameters.ConfigureAsUint( "InsideValue" );
22   this->m_Parameters.ConfigureAsUint( "OutsideValue" );
23   this->m_Parameters.ConfigureAsUint( "MinimumSize" );
24
25   this->m_Parameters.SetUint( "InsideValue", 1 );
26   this->m_Parameters.SetUint( "OutsideValue", 0 );
27   this->m_Parameters.SetUint( "MinimumSize", 100 );
28 }
29
30 // -------------------------------------------------------------------------
31 cpPluginsImageMeshFilters::RasterMeshFilter::
32 ~RasterMeshFilter( )
33 {
34 }
35
36 // -------------------------------------------------------------------------
37 void cpPluginsImageMeshFilters::RasterMeshFilter::
38 _GenerateData( )
39 {
40   typedef itk::Mesh< float, 3 >  _3F;
41   typedef itk::Mesh< double, 3 > _3D;
42
43   bool is_2d = false;
44   auto pd = this->GetInputData< vtkPolyData >( "Input" );
45   if( pd != NULL )
46   {
47     double bounds[ 6 ];
48     pd->GetBounds( bounds );
49     is_2d = ( bounds[ 4 ] == bounds[ 5 ] );
50
51   } // fi
52   if( !is_2d )
53   {
54     auto _3f = this->GetInputData< _3F >( "Input" );
55     auto _3d = this->GetInputData< _3D >( "Input" );
56     if     ( _3f != NULL ) this->_GD0_3D( _3f );
57     else if( _3d != NULL ) this->_GD0_3D( _3d );
58     else this->_Error( "No valid input mesh." );
59   }
60   else
61     this->_GD0_2D( pd );
62 }
63
64 // -------------------------------------------------------------------------
65 template< class _TMesh >
66 void cpPluginsImageMeshFilters::RasterMeshFilter::
67 _GD0_2D( _TMesh* mesh )
68 {
69   typedef unsigned char _TPixel;
70   typedef itk::ImageBase< 2 > _TImageBase;
71   typedef itk::Image< _TPixel, 2 > _TImage;
72   typedef cpExtensions::Algorithms::RasterContourFilter< _TImage > _TFilter;
73
74   auto in_im = this->GetInputData< _TImageBase >( "Template" );
75   if( in_im == NULL )
76     this->_Error( "A template is needed for 2D raster (this is temporal)." );
77   _TPixel inside = _TPixel( this->m_Parameters.GetUint( "InsideValue" ) );
78   _TPixel outside = _TPixel( this->m_Parameters.GetUint( "OutsideValue" ) );
79
80   auto filter = this->_CreateITK< _TFilter >( );
81   filter->ClearPoints( );
82   double pnt[ 3 ];
83   for( long i = 0; i < mesh->GetNumberOfPoints( ); ++i )
84   {
85     mesh->GetPoint( i, pnt );
86     filter->AddPoint( pnt[ 0 ], pnt[ 1 ] );
87
88   } // rof
89   filter->SetTemplate( in_im );
90   filter->SetInsideValue( inside );
91   filter->SetOutsideValue( outside );
92   filter->Update( );
93   this->GetOutput( "Output" )->SetITK( filter->GetOutput( ) );
94 }
95
96 // -------------------------------------------------------------------------
97 template< class _TMesh >
98 void cpPluginsImageMeshFilters::RasterMeshFilter::
99 _GD0_3D( _TMesh* mesh )
100 {
101   typedef unsigned char _TPixel;
102   typedef cpPlugins::DataObjects::BoundingBox _TBB;
103   typedef itk::ImageBase< _TMesh::PointDimension > _TImageBase;
104   typedef itk::Image< _TPixel, _TMesh::PointDimension > _TImage;
105   typedef itk::TriangleMeshToBinaryImageFilter< _TMesh, _TImage > _TFilter;
106   typedef typename _TImage::PointType _TPoint;
107
108   static const unsigned int PAD = 10;
109
110   _TFilter* filter = this->_CreateITK< _TFilter >( );
111
112   // Get configuration values
113   typename _TImage::SpacingType spac;
114   typename _TImage::SizeType size;
115   typename _TImage::PointType origin;
116   typename _TImage::DirectionType direction;
117   typename _TImage::IndexType index;
118
119   auto in_bb = this->GetInput< _TBB >( "Template" );
120   auto in_im = this->GetInputData< _TImageBase >( "Template" );
121   if( in_im != NULL )
122   {
123     spac = in_im->GetSpacing( );
124     size = in_im->GetLargestPossibleRegion( ).GetSize( );
125     origin = in_im->GetOrigin( );
126     direction = in_im->GetDirection( );
127     index = in_im->GetLargestPossibleRegion( ).GetIndex( );
128   }
129   else
130   {
131     _TPoint minBB, maxBB;
132     if( in_bb != NULL )
133     {
134       minBB = in_bb->GetMinimum< _TPoint >( );
135       maxBB = in_bb->GetMaximum< _TPoint >( );
136     }
137     else
138     {
139       auto bb = mesh->GetBoundingBox( );
140       minBB = bb->GetMinimum( );
141       maxBB = bb->GetMaximum( );
142
143     } // fi
144
145     double lx = double( maxBB[ 0 ] - minBB[ 0 ] );
146     double ly = double( maxBB[ 1 ] - minBB[ 1 ] );
147     double lz = double( maxBB[ 2 ] - minBB[ 2 ] );
148     double lm = ( lx < ly )? lx: ly;
149     lm = ( lm < lz )? lm: lz;
150
151     // Compute spacing
152     spac.Fill( lm / double( this->m_Parameters.GetUint( "MinimumSize" ) ) );
153
154     // Compute size
155     size[ 0 ] = ( unsigned long )( std::ceil( lx / spac[ 0 ] ) );
156     size[ 1 ] = ( unsigned long )( std::ceil( ly / spac[ 1 ] ) );
157     size[ 2 ] = ( unsigned long )( std::ceil( lz / spac[ 2 ] ) );
158
159     // ... add some padding pixels
160     size[ 0 ] += PAD;
161     size[ 1 ] += PAD;
162     size[ 2 ] += PAD;
163
164     // Set origin
165     origin = minBB;
166     origin[ 0 ] -= double( PAD >> 1 ) * spac[ 0 ];
167     origin[ 1 ] -= double( PAD >> 1 ) * spac[ 1 ];
168     origin[ 2 ] -= double( PAD >> 1 ) * spac[ 2 ];
169
170     // Remaining values
171     direction.SetIdentity( );
172     index.Fill( 0 );
173
174   } // fi
175
176   // Execute
177   filter->SetInput( mesh );
178   filter->SetSpacing( spac );
179   filter->SetSize( size );
180   filter->SetOrigin( origin );
181   filter->SetDirection( direction );
182   filter->SetIndex( index );
183   filter->SetInsideValue( this->m_Parameters.GetUint( "InsideValue" ) );
184   filter->SetOutsideValue( this->m_Parameters.GetUint( "OutsideValue" ) );
185   filter->Update( );
186   this->GetOutput( "Output" )->SetITK( filter->GetOutput( ) );
187 }
188
189 // eof - $RCSfile$