]> Creatis software - clitk.git/blob - itk/clitkInvertVFFilter.txx
itk4 compatibility
[clitk.git] / itk / clitkInvertVFFilter.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://www.centreleonberard.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 #ifndef __clitkInvertVFFilter_txx
19 #define __clitkInvertVFFilter_txx
20
21 namespace
22 {
23
24 //=========================================================================================================================
25 //helper class 1 to allow a threaded execution: add contributions of input to output and update weights
26 //=========================================================================================================================
27 template<class InputImageType, class OutputImageType> class ITK_EXPORT HelperClass1 : public itk::ImageToImageFilter<InputImageType, OutputImageType>
28 {
29
30 public:
31   /** Standard class typedefs. */
32   typedef HelperClass1  Self;
33   typedef itk::ImageToImageFilter<InputImageType,OutputImageType> Superclass;
34   typedef itk::SmartPointer<Self>         Pointer;
35   typedef itk::SmartPointer<const Self>   ConstPointer;
36
37   /** Method for creation through the object factory. */
38   itkNewMacro(Self);
39
40   /** Run-time type information (and related methods) */
41   itkTypeMacro( HelperClass1, ImageToImageFilter );
42
43   /** Constants for the image dimensions */
44   itkStaticConstMacro(ImageDimension, unsigned int,InputImageType::ImageDimension);
45
46
47   //Typedefs
48   typedef typename OutputImageType::PixelType        PixelType;
49   typedef itk::Image<double, ImageDimension > WeightsImageType;
50   typedef itk::Image<itk::SimpleFastMutexLock, ImageDimension > MutexImageType;
51
52   //===================================================================================
53   //Set methods
54   void SetWeights(const typename WeightsImageType::Pointer input) {
55     m_Weights = input;
56     this->Modified();
57   }
58   void SetMutexImage(const typename MutexImageType::Pointer input) {
59     m_MutexImage=input;
60     this->Modified();
61     m_ThreadSafe=true;
62   }
63
64   //Get methods
65   typename  WeightsImageType::Pointer GetWeights() {
66     return m_Weights;
67   }
68
69   /** Typedef to describe the output image region type. */
70   typedef typename OutputImageType::RegionType OutputImageRegionType;
71
72 protected:
73   HelperClass1();
74   ~HelperClass1() {};
75
76   //the actual processing
77   void BeforeThreadedGenerateData();
78 #if ITK_VERSION_MAJOR >= 4
79   void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId );
80 #else
81   void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, int threadId );
82 #endif
83
84   //member data
85   typename  WeightsImageType::Pointer m_Weights;
86   typename MutexImageType::Pointer m_MutexImage;
87   bool m_ThreadSafe;
88
89 };
90
91
92
93 //=========================================================================================================================
94 //Member functions of the helper class 1
95 //=========================================================================================================================
96
97
98 //=========================================================================================================================
99 //Empty constructor
100 template<class InputImageType, class OutputImageType >
101 HelperClass1<InputImageType, OutputImageType>::HelperClass1()
102 {
103   m_ThreadSafe=false;
104 }
105
106 //=========================================================================================================================
107 //Before threaded data
108 template<class InputImageType, class OutputImageType >
109 void HelperClass1<InputImageType, OutputImageType>::BeforeThreadedGenerateData()
110 {
111   //std::cout << "HelperClass1::BeforeThreadedGenerateData - IN" << std::endl;
112   //Since we will add, put to zero!
113   this->GetOutput()->FillBuffer(itk::NumericTraits<double>::Zero);
114   this->GetWeights()->FillBuffer(itk::NumericTraits<double>::Zero);
115 }
116
117 //=========================================================================================================================
118 //update the output for the outputRegionForThread
119 template<class InputImageType, class OutputImageType>
120 #if ITK_VERSION_MAJOR >= 4
121 void HelperClass1<InputImageType, OutputImageType>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId )
122 #else
123 void HelperClass1<InputImageType, OutputImageType>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, int threadId )
124 #endif
125 {
126   //std::cout << "HelperClass1::ThreadedGenerateData - IN" << std::endl;
127   //Get pointer to the input
128   typename InputImageType::ConstPointer inputPtr = this->GetInput();
129
130   //Get pointer to the output
131   typename OutputImageType::Pointer outputPtr = this->GetOutput();
132   //typename OutputImageType::SizeType size=outputPtr->GetLargestPossibleRegion().GetSize();
133
134   //Iterator over input
135   typedef itk::ImageRegionConstIteratorWithIndex<InputImageType> InputImageIteratorType;
136
137   //define them over the outputRegionForThread
138   InputImageIteratorType inputIt(inputPtr, outputRegionForThread);
139
140   //Initialize
141   typename InputImageType::IndexType index;
142   itk::ContinuousIndex<double,ImageDimension> contIndex, inContIndex;
143   typename InputImageType::PointType ipoint;
144   typename OutputImageType::PointType opoint;
145   typedef typename OutputImageType::PixelType DisplacementType;
146   DisplacementType displacement;
147   inputIt.GoToBegin();
148
149   //define some temp variables
150   signed long baseIndex[ImageDimension];
151   double distance[ImageDimension];
152   unsigned int dim, counter, upper;
153   double totalOverlap,overlap;
154   typename OutputImageType::IndexType neighIndex;
155
156   //Find the number of neighbors
157   unsigned int neighbors =  1 << ImageDimension;
158
159   //==================================================================================================
160   //Loop over the output region and add the intensities from the input to the output and the weight to the weights
161   //==================================================================================================
162   while( !inputIt.IsAtEnd() ) {
163     // get the input image index
164     index = inputIt.GetIndex();
165     inputPtr->TransformIndexToPhysicalPoint( index,ipoint );
166
167     // get the required displacement
168     displacement = inputIt.Get();
169
170     // compute the required output image point
171     for(unsigned int j = 0; j < ImageDimension; j++ ) opoint[j] = ipoint[j] + (double)displacement[j];
172
173     // Update the output and the weights
174     if(outputPtr->TransformPhysicalPointToContinuousIndex(opoint, contIndex ) ) {
175       for(dim = 0; dim < ImageDimension; dim++) {
176         // The following  block is equivalent to the following line without
177         // having to call floor. (Only for positive inputs, we already now that is in the image)
178         // baseIndex[dim] = (long) vcl_floor(contIndex[dim] );
179
180         baseIndex[dim] = (long) contIndex[dim];
181         distance[dim] = contIndex[dim] - double( baseIndex[dim] );
182       }
183
184       //Add contribution for each neighbor
185       totalOverlap = itk::NumericTraits<double>::Zero;
186       for( counter = 0; counter < neighbors ; counter++ ) {
187         overlap = 1.0;          // fraction overlap
188         upper = counter;  // each bit indicates upper/lower neighbour
189
190         // get neighbor index and overlap fraction
191         for( dim = 0; dim < ImageDimension; dim++ ) {
192           if ( upper & 1 ) {
193             neighIndex[dim] = baseIndex[dim] + 1;
194             overlap *= distance[dim];
195           } else {
196             neighIndex[dim] = baseIndex[dim];
197             overlap *= 1.0 - distance[dim];
198           }
199           upper >>= 1;
200         }
201
202
203
204         //Set neighbor value only if overlap is not zero
205         if( (overlap>0.0)) // &&
206           //                    (static_cast<unsigned int>(neighIndex[0])<size[0]) &&
207           //                    (static_cast<unsigned int>(neighIndex[1])<size[1]) &&
208           //                    (static_cast<unsigned int>(neighIndex[2])<size[2]) &&
209           //                    (neighIndex[0]>=0) &&
210           //                    (neighIndex[1]>=0) &&
211           //                    (neighIndex[2]>=0) )
212         {
213           //what to store? the original displacement vector?
214           if (! m_ThreadSafe) {
215             //Set the pixel and weight at neighIndex
216             outputPtr->SetPixel(neighIndex, outputPtr->GetPixel(neighIndex) - (displacement*overlap));
217             m_Weights->SetPixel(neighIndex, m_Weights->GetPixel(neighIndex) + overlap);
218           }
219
220           else {
221             //Entering critilal section: shared memory
222             m_MutexImage->GetPixel(neighIndex).Lock();
223
224             //Set the pixel and weight at neighIndex
225             outputPtr->SetPixel(neighIndex, outputPtr->GetPixel(neighIndex) - (displacement*overlap));
226             m_Weights->SetPixel(neighIndex, m_Weights->GetPixel(neighIndex) + overlap);
227
228             //Unlock
229             m_MutexImage->GetPixel(neighIndex).Unlock();
230
231           }
232           //Add to total overlap
233           totalOverlap += overlap;
234         }
235
236         if( totalOverlap == 1.0 ) {
237           // finished
238           break;
239         }
240       }
241     }
242
243     ++inputIt;
244   }
245
246   //std::cout << "HelperClass1::ThreadedGenerateData - OUT" << std::endl;
247 }
248
249
250
251 //=========================================================================================================================
252 //helper class 2 to allow a threaded execution of normalisation by the weights
253 //=========================================================================================================================
254 template<class InputImageType, class OutputImageType> class HelperClass2 : public itk::ImageToImageFilter<InputImageType, OutputImageType>
255 {
256
257 public:
258   /** Standard class typedefs. */
259   typedef HelperClass2  Self;
260   typedef itk::ImageToImageFilter<InputImageType,OutputImageType> Superclass;
261   typedef itk::SmartPointer<Self>         Pointer;
262   typedef itk::SmartPointer<const Self>   ConstPointer;
263
264   /** Method for creation through the object factory. */
265   itkNewMacro(Self);
266
267   /** Run-time type information (and related methods) */
268   itkTypeMacro( HelperClass2, ImageToImageFilter );
269
270   /** Constants for the image dimensions */
271   itkStaticConstMacro(ImageDimension, unsigned int,InputImageType::ImageDimension);
272
273   //Typedefs
274   typedef typename OutputImageType::PixelType        PixelType;
275   typedef itk::Image<double,ImageDimension> WeightsImageType;
276
277   //Set methods
278   void SetWeights(const typename WeightsImageType::Pointer input) {
279     m_Weights = input;
280     this->Modified();
281   }
282   void SetEdgePaddingValue(PixelType value) {
283     m_EdgePaddingValue = value;
284     this->Modified();
285   }
286
287   /** Typedef to describe the output image region type. */
288   typedef typename OutputImageType::RegionType OutputImageRegionType;
289
290 protected:
291   HelperClass2();
292   ~HelperClass2() {};
293
294
295   //the actual processing
296 #if ITK_VERSION_MAJOR >= 4
297   void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId );
298 #else
299   void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, int threadId );
300 #endif
301
302   //member data
303   typename     WeightsImageType::Pointer m_Weights;
304   PixelType m_EdgePaddingValue;
305
306 } ;
307
308
309
310 //=========================================================================================================================
311 //Member functions of the helper class 2
312 //=========================================================================================================================
313
314
315 //=========================================================================================================================
316 //Empty constructor
317 template<class InputImageType, class OutputImageType > HelperClass2<InputImageType, OutputImageType>::HelperClass2()
318 {
319   m_EdgePaddingValue=itk::NumericTraits<PixelType>::Zero;
320 }
321
322
323 //=========================================================================================================================
324 //update the output for the outputRegionForThread
325 #if ITK_VERSION_MAJOR >= 4
326 template<class InputImageType, class OutputImageType > void HelperClass2<InputImageType, OutputImageType>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId )
327 #else
328 template<class InputImageType, class OutputImageType > void HelperClass2<InputImageType, OutputImageType>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, int threadId )
329 #endif
330 {
331   //std::cout << "HelperClass2::ThreadedGenerateData - IN" << std::endl;
332   
333   //Get pointer to the input
334   typename InputImageType::ConstPointer inputPtr = this->GetInput();
335
336   //Get pointer to the output
337   typename OutputImageType::Pointer outputPtr = this->GetOutput();
338
339   //Iterators over input, weigths  and output
340   typedef itk::ImageRegionConstIterator<InputImageType> InputImageIteratorType;
341   typedef itk::ImageRegionIterator<OutputImageType> OutputImageIteratorType;
342   typedef itk::ImageRegionIterator<WeightsImageType> WeightsImageIteratorType;
343
344   //define them over the outputRegionForThread
345   OutputImageIteratorType outputIt(outputPtr, outputRegionForThread);
346   InputImageIteratorType inputIt(inputPtr, outputRegionForThread);
347   WeightsImageIteratorType weightsIt(m_Weights, outputRegionForThread);
348
349
350   //==================================================================================================
351   //loop over the output and normalize the input, remove holes
352   PixelType neighValue;
353   double  zero = itk::NumericTraits<double>::Zero;
354   while (!outputIt.IsAtEnd()) {
355     //the weight is not zero
356     if (weightsIt.Get() != zero) {
357       //divide by the weight
358       outputIt.Set(static_cast<PixelType>(inputIt.Get()/weightsIt.Get()));
359     }
360
361     //copy the value of the  neighbour that was just processed
362     else {
363       if(!outputIt.IsAtBegin()) {
364         //go back
365         --outputIt;
366
367         //Neighbour cannot have zero weight because it should be filled already
368         neighValue=outputIt.Get();
369         ++outputIt;
370         outputIt.Set(neighValue);
371         //DD("hole filled");
372       } else {
373         //DD("is at begin, setting edgepadding value");
374         outputIt.Set(m_EdgePaddingValue);
375       }
376     }
377     ++weightsIt;
378     ++outputIt;
379     ++inputIt;
380
381   }//end while
382   
383   //std::cout << "HelperClass2::ThreadedGenerateData - OUT" << std::endl;
384   
385 }//end member
386
387
388 }//end nameless namespace
389
390
391
392 namespace clitk
393 {
394
395 //=========================================================================================================================
396 // The rest is the InvertVFFilter
397 //=========================================================================================================================
398
399 //=========================================================================================================================
400 //constructor
401 template <class InputImageType, class OutputImageType>
402 InvertVFFilter<InputImageType, OutputImageType>::InvertVFFilter()
403 {
404   m_EdgePaddingValue=itk::NumericTraits<PixelType>::Zero; //no other reasonable value?
405   m_ThreadSafe=false;
406   m_Verbose=false;
407 }
408
409 //=========================================================================================================================
410 //Update
411 template <class InputImageType, class OutputImageType> void InvertVFFilter<InputImageType, OutputImageType>::GenerateData()
412 {
413   //std::cout << "InvertVFFilter::GenerateData - IN" << std::endl;
414
415   //Get the properties of the input
416   typename InputImageType::ConstPointer inputPtr=this->GetInput();
417   typename WeightsImageType::RegionType region = inputPtr->GetLargestPossibleRegion();
418
419   //Allocate the weights
420   typename WeightsImageType::Pointer weights=WeightsImageType::New();
421   weights->SetOrigin(inputPtr->GetOrigin());
422   weights->SetRegions(region);
423   weights->Allocate();
424   weights->SetSpacing(inputPtr->GetSpacing());
425
426   //===========================================================================
427   //Inversion is divided in in two loops, for each we will call a threaded helper class
428   //1. add contribution of input to output and update weights
429   //2. normalize the output by the weight and remove holes
430   //===========================================================================
431
432
433   //===========================================================================
434   //1. add contribution of input to output and update weights
435
436   //Define an internal image type
437
438   typedef itk::Image<itk::Vector<double,ImageDimension>, ImageDimension > InternalImageType;
439
440   //Call threaded helper class 1
441   typedef HelperClass1<InputImageType, InternalImageType > HelperClass1Type;
442   typename HelperClass1Type::Pointer helper1=HelperClass1Type::New();
443
444   //Set input
445   if(m_NumberOfThreadsIsGiven)helper1->SetNumberOfThreads(m_NumberOfThreads);
446   helper1->SetInput(inputPtr);
447   helper1->SetWeights(weights);
448
449   //Threadsafe?
450   if(m_ThreadSafe) {
451     //Allocate the mutex image
452     typename MutexImageType::Pointer mutex=InvertVFFilter::MutexImageType::New();
453     mutex->SetRegions(region);
454     mutex->Allocate();
455     mutex->SetSpacing(inputPtr->GetSpacing());
456     helper1->SetMutexImage(mutex);
457     if (m_Verbose) std::cout <<"Inverting using a thread-safe algorithm" <<std::endl;
458   } else  if(m_Verbose)std::cout <<"Inverting using a thread-unsafe algorithm" <<std::endl;
459
460   //Execute helper class
461   helper1->Update();
462
463   //Get the output
464   typename InternalImageType::Pointer temp= helper1->GetOutput();
465   weights=helper1->GetWeights();
466
467
468   //===========================================================================
469   //2. Normalize the output by the weights and remove holes
470   //Call threaded helper class
471   typedef HelperClass2<InternalImageType, OutputImageType> HelperClass2Type;
472   typename HelperClass2Type::Pointer helper2=HelperClass2Type::New();
473
474   //Set temporary output as input
475   helper2->SetInput(temp);
476   helper2->SetWeights(weights);
477   helper2->SetEdgePaddingValue(m_EdgePaddingValue);
478
479   //Execute helper class
480   if (m_Verbose) std::cout << "Normalizing the output VF..."<<std::endl;
481   helper2->Update();
482
483   //Set the output
484   this->SetNthOutput(0, helper2->GetOutput());
485 }
486
487
488
489 }
490
491 #endif