]> Creatis software - clitk.git/blob - itk/clitkSliceBySliceRelativePositionFilter.txx
add autocrop and options
[clitk.git] / itk / clitkSliceBySliceRelativePositionFilter.txx
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to: 
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://oncora1.lyon.fnclcc.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17   ======================================================================-====*/
18
19 // clitk
20 #include "clitkSegmentationUtils.h"
21 #include "clitkExtractSliceFilter.h"
22 #include "clitkResampleImageWithOptionsFilter.h"
23
24 // itk
25 #include <itkJoinSeriesImageFilter.h>
26
27 //--------------------------------------------------------------------
28 template <class ImageType>
29 clitk::SliceBySliceRelativePositionFilter<ImageType>::
30 SliceBySliceRelativePositionFilter():
31   clitk::FilterBase(),
32   itk::ImageToImageFilter<ImageType, ImageType>()
33 {
34   this->SetNumberOfRequiredInputs(2);
35   SetDirection(2);
36   SetObjectBackgroundValue(0);  
37   SetFuzzyThreshold(0.6);
38   SetOrientationTypeString("Left");
39   SetIntermediateSpacing(10);
40   ResampleBeforeRelativePositionFilterOff();
41   UniqueConnectedComponentBySliceOff();
42   NotFlagOff();
43 }
44 //--------------------------------------------------------------------
45
46
47 //--------------------------------------------------------------------
48 template <class ImageType>
49 void 
50 clitk::SliceBySliceRelativePositionFilter<ImageType>::
51 SetInput(const ImageType * image) 
52 {
53   // Process object is not const-correct so the const casting is required.
54   this->SetNthInput(0, const_cast<ImageType *>(image));
55 }
56 //--------------------------------------------------------------------
57   
58
59 //--------------------------------------------------------------------
60 template <class ImageType>
61 void 
62 clitk::SliceBySliceRelativePositionFilter<ImageType>::
63 SetInputObject(const ImageType * image) 
64 {
65   // Process object is not const-correct so the const casting is required.
66   this->SetNthInput(1, const_cast<ImageType *>(image));
67 }
68 //--------------------------------------------------------------------
69   
70
71 //--------------------------------------------------------------------
72 template <class ImageType>
73 void 
74 clitk::SliceBySliceRelativePositionFilter<ImageType>::
75 GenerateInputRequestedRegion() 
76 {
77   // Call default
78   itk::ImageToImageFilter<ImageType, ImageType>::GenerateInputRequestedRegion();
79   // Get input pointers and set requested region to common region
80   ImagePointer input1 = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(0));
81   ImagePointer input2 = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(1));
82   input1->SetRequestedRegion(input1->GetLargestPossibleRegion());
83   input2->SetRequestedRegion(input2->GetLargestPossibleRegion());
84 }
85 //--------------------------------------------------------------------
86
87
88 //--------------------------------------------------------------------
89 template <class ImageType>
90 void 
91 clitk::SliceBySliceRelativePositionFilter<ImageType>::
92 GenerateOutputInformation() 
93 {
94   // Get input pointer
95   input = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(0));
96   object = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(1));
97
98   //--------------------------------------------------------------------
99   // Resample object to the same spacing than input
100   if (!clitk::HaveSameSpacing<ImageType, ImageType>(object, input)) {
101     StartNewStep("Resample object to the same spacing than input");
102     m_working_object = clitk::ResampleImageSpacing<ImageType>(object, input->GetSpacing());
103     StopCurrentStep<ImageType>(m_working_object);
104   }
105   else {
106     m_working_object = object;
107   }
108   
109   //--------------------------------------------------------------------
110   // Pad object to the same size than input
111   if (!clitk::HaveSameSizeAndSpacing<ImageType, ImageType>(m_working_object, input)) {
112     StartNewStep("Pad object to the same size than input");
113     m_working_object = clitk::ResizeImageLike<ImageType>(m_working_object, 
114                                                           input, 
115                                                           GetObjectBackgroundValue());
116     StopCurrentStep<ImageType>(m_working_object);
117   }
118   else {
119   }
120
121   /*
122     - extract vector of slices in input, in object
123     - slice by slice rel position
124     - joint result
125     - post process
126   */
127
128
129   //--------------------------------------------------------------------
130   // Extract input slices
131   StartNewStep("Extract input slices");
132   typedef clitk::ExtractSliceFilter<ImageType> ExtractSliceFilterType;
133   typename ExtractSliceFilterType::Pointer extractSliceFilter = ExtractSliceFilterType::New();
134   extractSliceFilter->SetInput(input);
135   extractSliceFilter->SetDirection(GetDirection());
136   extractSliceFilter->Update();
137   typedef typename ExtractSliceFilterType::SliceType SliceType;
138   std::vector<typename SliceType::Pointer> mInputSlices;
139   extractSliceFilter->GetOutputSlices(mInputSlices);
140   StopCurrentStep<SliceType>(mInputSlices[0]);
141   
142   //--------------------------------------------------------------------
143   // Extract object slices
144   StartNewStep("Extract object slices");
145   extractSliceFilter = ExtractSliceFilterType::New();
146   extractSliceFilter->SetInput(m_working_object);//object);
147   extractSliceFilter->SetDirection(GetDirection());
148   extractSliceFilter->Update();
149   std::vector<typename SliceType::Pointer> mObjectSlices;
150   extractSliceFilter->GetOutputSlices(mObjectSlices);
151   StopCurrentStep<SliceType>(mObjectSlices[0]);
152
153   //--------------------------------------------------------------------
154   // Perform slice by slice relative position
155   StartNewStep("Perform slice by slice relative position");
156   for(unsigned int i=0; i<mInputSlices.size(); i++) {
157     // Select main CC in each object slice (required ?)
158     mObjectSlices[i] = Labelize<SliceType>(mObjectSlices[i], 0, true, 1);
159     mObjectSlices[i] = KeepLabels<SliceType>(mObjectSlices[i], 0, 1, 1, 1, true);
160
161     // Relative position
162     typedef clitk::AddRelativePositionConstraintToLabelImageFilter<SliceType> RelPosFilterType;
163     typename RelPosFilterType::Pointer relPosFilter = RelPosFilterType::New();
164     relPosFilter->VerboseStepOff();
165     relPosFilter->WriteStepOff();
166     relPosFilter->SetCurrentStepBaseId(this->GetCurrentStepId());
167     relPosFilter->SetInput(mInputSlices[i]); 
168     relPosFilter->SetInputObject(mObjectSlices[i]); 
169     relPosFilter->SetNotFlag(GetNotFlag());
170     relPosFilter->SetOrientationTypeString(this->GetOrientationTypeString());
171     relPosFilter->SetIntermediateSpacing(this->GetIntermediateSpacing());
172     relPosFilter->SetResampleBeforeRelativePositionFilter(this->GetResampleBeforeRelativePositionFilter());
173     relPosFilter->SetFuzzyThreshold(this->GetFuzzyThreshold());
174     relPosFilter->AutoCropFlagOff(); // important ! because we join the slices after this loop
175     relPosFilter->Update();
176     mInputSlices[i] = relPosFilter->GetOutput();
177
178     // Select main CC if needed
179     if (GetUniqueConnectedComponentBySlice()) {
180       mInputSlices[i] = Labelize<SliceType>(mInputSlices[i], 0, true, 1);
181       mInputSlices[i] = KeepLabels<SliceType>(mInputSlices[i], 0, 1, 1, 1, true);
182     }
183
184   }
185
186   typedef itk::JoinSeriesImageFilter<SliceType, ImageType> JoinSeriesFilterType;
187   typename JoinSeriesFilterType::Pointer joinFilter = JoinSeriesFilterType::New();
188   joinFilter->SetOrigin(input->GetOrigin()[GetDirection()]);
189   joinFilter->SetSpacing(input->GetSpacing()[GetDirection()]);
190   for(unsigned int i=0; i<mInputSlices.size(); i++) {
191     joinFilter->PushBackInput(mInputSlices[i]);
192   }
193   joinFilter->Update();
194   m_working_input = joinFilter->GetOutput();
195   StopCurrentStep<ImageType>(m_working_input);
196
197   //--------------------------------------------------------------------
198   // Step 7: autocrop
199   if (GetAutoCropFlag()) {
200     StartNewStep("Final AutoCrop");
201     typedef clitk::AutoCropFilter<ImageType> CropFilterType;
202     typename CropFilterType::Pointer cropFilter = CropFilterType::New();
203     cropFilter->SetInput(m_working_input);
204     cropFilter->ReleaseDataFlagOff();
205     cropFilter->Update();   
206     m_working_input = cropFilter->GetOutput();
207     StopCurrentStep<ImageType>(m_working_input);    
208   }
209
210   // Update output info
211   this->GetOutput(0)->SetRegions(m_working_input->GetLargestPossibleRegion());  
212 }
213 //--------------------------------------------------------------------
214
215
216 //--------------------------------------------------------------------
217 template <class ImageType>
218 void 
219 clitk::SliceBySliceRelativePositionFilter<ImageType>::
220 GenerateData() 
221 {
222   // Get input pointer
223   //--------------------------------------------------------------------
224   //--------------------------------------------------------------------  
225   // Final Step -> set output
226   //this->SetNthOutput(0, m_working_input);
227   this->GraftOutput(m_working_input);
228   return;
229 }
230 //--------------------------------------------------------------------
231