]> Creatis software - clitk.git/blob - itk/clitkSegmentationUtils.txx
replaces round by itk:Math:Round (because round not found in Msvc2008)
[clitk.git] / itk / clitkSegmentationUtils.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
19 // clitk
20 #include "clitkSetBackgroundImageFilter.h"
21 #include "clitkSliceBySliceRelativePositionFilter.h"
22 #include "clitkCropLikeImageFilter.h"
23 #include "clitkMemoryUsage.h"
24
25 // itk
26 #include <itkConnectedComponentImageFilter.h>
27 #include <itkRelabelComponentImageFilter.h>
28 #include <itkBinaryThresholdImageFilter.h>
29 #include <itkPasteImageFilter.h>
30 #include <itkStatisticsLabelMapFilter.h>
31 #include <itkBinaryBallStructuringElement.h>
32 #include <itkBinaryDilateImageFilter.h>
33 #include <itkConstantPadImageFilter.h>
34 #include <itkImageSliceIteratorWithIndex.h>
35 #include <itkBinaryMorphologicalOpeningImageFilter.h>
36 #include <itkImageDuplicator.h>
37
38 namespace clitk {
39
40   //--------------------------------------------------------------------
41   template<class ImageType>
42   void ComputeBBFromImageRegion(const ImageType * image, 
43                                 typename ImageType::RegionType region,
44                                 typename itk::BoundingBox<unsigned long, 
45                                 ImageType::ImageDimension>::Pointer bb) {
46     typedef typename ImageType::IndexType IndexType;
47     IndexType firstIndex;
48     IndexType lastIndex;
49     for(unsigned int i=0; i<image->GetImageDimension(); i++) {
50       firstIndex[i] = region.GetIndex()[i];
51       lastIndex[i] = firstIndex[i]+region.GetSize()[i];
52     }
53
54     typedef itk::BoundingBox<unsigned long, 
55                              ImageType::ImageDimension> BBType;
56     typedef typename BBType::PointType PointType;
57     PointType lastPoint;
58     PointType firstPoint;
59     image->TransformIndexToPhysicalPoint(firstIndex, firstPoint);
60     image->TransformIndexToPhysicalPoint(lastIndex, lastPoint);
61
62     bb->SetMaximum(lastPoint);
63     bb->SetMinimum(firstPoint);
64   }
65   //--------------------------------------------------------------------
66
67
68   //--------------------------------------------------------------------
69   template<int Dimension>
70   void ComputeBBIntersection(typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbo, 
71                              typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbi1, 
72                              typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbi2) {
73
74     typedef itk::BoundingBox<unsigned long, Dimension> BBType;
75     typedef typename BBType::PointType PointType;
76     PointType lastPoint;
77     PointType firstPoint;
78
79     for(unsigned int i=0; i<Dimension; i++) {
80       firstPoint[i] = std::max(bbi1->GetMinimum()[i], 
81                                bbi2->GetMinimum()[i]);
82       lastPoint[i] = std::min(bbi1->GetMaximum()[i], 
83                               bbi2->GetMaximum()[i]);
84     }
85
86     bbo->SetMaximum(lastPoint);
87     bbo->SetMinimum(firstPoint);
88   }
89   //--------------------------------------------------------------------
90
91
92   //--------------------------------------------------------------------
93   template<class ImageType>
94   void ComputeRegionFromBB(const ImageType * image, 
95                            const typename itk::BoundingBox<unsigned long, 
96                                                            ImageType::ImageDimension>::Pointer bb, 
97                            typename ImageType::RegionType & region) {
98     // Types
99     typedef typename ImageType::IndexType  IndexType;
100     typedef typename ImageType::PointType  PointType;
101     typedef typename ImageType::RegionType RegionType;
102     typedef typename ImageType::SizeType   SizeType;
103
104     // Region starting point
105     IndexType regionStart;
106     PointType start = bb->GetMinimum();
107     image->TransformPhysicalPointToIndex(start, regionStart);
108     
109     // Region size
110     SizeType regionSize;
111     PointType maxs = bb->GetMaximum();
112     PointType mins = bb->GetMinimum();
113     for(unsigned int i=0; i<ImageType::ImageDimension; i++) {
114       regionSize[i] = lrint((maxs[i] - mins[i])/image->GetSpacing()[i]);
115     }
116    
117     // Create region
118     region.SetIndex(regionStart);
119     region.SetSize(regionSize);
120   }
121   //--------------------------------------------------------------------
122
123   //--------------------------------------------------------------------
124   template<class ImageType, class TMaskImageType>
125   typename ImageType::Pointer
126   SetBackground(const ImageType * input, 
127                 const TMaskImageType * mask, 
128                 typename TMaskImageType::PixelType maskBG,
129                 typename ImageType::PixelType outValue, 
130                 bool inPlace) {
131     typedef SetBackgroundImageFilter<ImageType, TMaskImageType, ImageType> 
132       SetBackgroundImageFilterType;
133     typename SetBackgroundImageFilterType::Pointer setBackgroundFilter 
134       = SetBackgroundImageFilterType::New();
135     //  if (inPlace) setBackgroundFilter->ReleaseDataFlagOn(); // No seg fault
136     setBackgroundFilter->SetInPlace(inPlace); // This is important to keep memory low
137     setBackgroundFilter->SetInput(input);
138     setBackgroundFilter->SetInput2(mask);
139     setBackgroundFilter->SetMaskValue(maskBG);
140     setBackgroundFilter->SetOutsideValue(outValue);
141     setBackgroundFilter->Update();
142     return setBackgroundFilter->GetOutput();
143   }
144   //--------------------------------------------------------------------
145
146
147   //--------------------------------------------------------------------
148   template<class ImageType>
149   int GetNumberOfConnectedComponentLabels(const ImageType * input, 
150                                           typename ImageType::PixelType BG, 
151                                           bool isFullyConnected) {
152     // Connected Component label 
153     typedef itk::ConnectedComponentImageFilter<ImageType, ImageType> ConnectFilterType;
154     typename ConnectFilterType::Pointer connectFilter = ConnectFilterType::New();
155     connectFilter->SetInput(input);
156     connectFilter->SetBackgroundValue(BG);
157     connectFilter->SetFullyConnected(isFullyConnected);
158     connectFilter->Update();
159   
160     // Return result
161     return connectFilter->GetObjectCount();
162   }
163   //--------------------------------------------------------------------
164
165   //--------------------------------------------------------------------
166   /*
167     Warning : in this cas, we consider outputType like inputType, not
168     InternalImageType. Be sure it fits.
169   */
170   template<class ImageType>
171   typename ImageType::Pointer
172   Labelize(const ImageType * input, 
173            typename ImageType::PixelType BG, 
174            bool isFullyConnected, 
175            int minimalComponentSize) {
176     // InternalImageType for storing large number of component
177     typedef itk::Image<int, ImageType::ImageDimension> InternalImageType;
178   
179     // Connected Component label 
180     typedef itk::ConnectedComponentImageFilter<ImageType, InternalImageType> ConnectFilterType;
181     typename ConnectFilterType::Pointer connectFilter = ConnectFilterType::New();
182     //  connectFilter->ReleaseDataFlagOn(); 
183     connectFilter->SetInput(input);
184     connectFilter->SetBackgroundValue(BG);
185     connectFilter->SetFullyConnected(isFullyConnected);
186   
187     // Sort by size and remove too small area.
188     typedef itk::RelabelComponentImageFilter<InternalImageType, ImageType> RelabelFilterType;
189     typename RelabelFilterType::Pointer relabelFilter = RelabelFilterType::New();
190     //  relabelFilter->ReleaseDataFlagOn(); // if yes, fail when ExplosionControlledThresholdConnectedImageFilter ???
191     relabelFilter->SetInput(connectFilter->GetOutput());
192     relabelFilter->SetMinimumObjectSize(minimalComponentSize);
193     relabelFilter->Update();
194
195     // Return result
196     typename ImageType::Pointer output = relabelFilter->GetOutput();
197     return output;
198   }
199   //--------------------------------------------------------------------
200
201
202   //--------------------------------------------------------------------
203   /*
204     Warning : in this cas, we consider outputType like inputType, not
205     InternalImageType. Be sure it fits.
206   */
207   template<class ImageType>
208   typename ImageType::Pointer
209   LabelizeAndCountNumberOfObjects(const ImageType * input, 
210                                   typename ImageType::PixelType BG, 
211                                   bool isFullyConnected, 
212                                   int minimalComponentSize, 
213                                   int & nb) {
214     // InternalImageType for storing large number of component
215     typedef itk::Image<int, ImageType::ImageDimension> InternalImageType;
216   
217     // Connected Component label 
218     typedef itk::ConnectedComponentImageFilter<ImageType, InternalImageType> ConnectFilterType;
219     typename ConnectFilterType::Pointer connectFilter = ConnectFilterType::New();
220     //  connectFilter->ReleaseDataFlagOn(); 
221     connectFilter->SetInput(input);
222     connectFilter->SetBackgroundValue(BG);
223     connectFilter->SetFullyConnected(isFullyConnected);
224   
225     // Sort by size and remove too small area.
226     typedef itk::RelabelComponentImageFilter<InternalImageType, ImageType> RelabelFilterType;
227     typename RelabelFilterType::Pointer relabelFilter = RelabelFilterType::New();
228     //  relabelFilter->ReleaseDataFlagOn(); // if yes, fail when ExplosionControlledThresholdConnectedImageFilter ???
229     relabelFilter->SetInput(connectFilter->GetOutput());
230     relabelFilter->SetMinimumObjectSize(minimalComponentSize);
231     relabelFilter->Update();
232
233     nb = relabelFilter->GetNumberOfObjects();
234     // DD(relabelFilter->GetOriginalNumberOfObjects());
235     // DD(relabelFilter->GetSizeOfObjectsInPhysicalUnits()[0]);
236
237     // Return result
238     typename ImageType::Pointer output = relabelFilter->GetOutput();
239     return output;
240   }
241   //--------------------------------------------------------------------
242
243
244
245   //--------------------------------------------------------------------
246   template<class ImageType>
247   typename ImageType::Pointer
248   RemoveLabels(const ImageType * input, 
249                typename ImageType::PixelType BG,
250                std::vector<typename ImageType::PixelType> & labelsToRemove) {
251     assert(labelsToRemove.size() != 0);
252     typename ImageType::Pointer working_image;// = input;
253     for (unsigned int i=0; i <labelsToRemove.size(); i++) {
254       typedef SetBackgroundImageFilter<ImageType, ImageType> SetBackgroundImageFilterType;
255       typename SetBackgroundImageFilterType::Pointer setBackgroundFilter = SetBackgroundImageFilterType::New();
256       setBackgroundFilter->SetInput(input);
257       setBackgroundFilter->SetInput2(input);
258       setBackgroundFilter->SetMaskValue(labelsToRemove[i]);
259       setBackgroundFilter->SetOutsideValue(BG);
260       setBackgroundFilter->Update();
261       working_image = setBackgroundFilter->GetOutput();
262     }
263     return working_image;
264   }
265   //--------------------------------------------------------------------
266
267
268   //--------------------------------------------------------------------
269   template<class ImageType>
270   typename ImageType::Pointer
271   KeepLabels(const ImageType * input, 
272              typename ImageType::PixelType BG, 
273              typename ImageType::PixelType FG, 
274              typename ImageType::PixelType firstKeep, 
275              typename ImageType::PixelType lastKeep, 
276              bool useLastKeep) {
277     typedef itk::BinaryThresholdImageFilter<ImageType, ImageType> BinarizeFilterType; 
278     typename BinarizeFilterType::Pointer binarizeFilter = BinarizeFilterType::New();
279     binarizeFilter->SetInput(input);
280     binarizeFilter->SetLowerThreshold(firstKeep);
281     if (useLastKeep) binarizeFilter->SetUpperThreshold(lastKeep);
282     binarizeFilter->SetInsideValue(FG);
283     binarizeFilter->SetOutsideValue(BG);
284     binarizeFilter->Update();
285     return binarizeFilter->GetOutput();
286   }
287   //--------------------------------------------------------------------
288
289
290   //--------------------------------------------------------------------
291   template<class ImageType>
292   typename ImageType::Pointer
293   LabelizeAndSelectLabels(const ImageType * input,
294                           typename ImageType::PixelType BG, 
295                           typename ImageType::PixelType FG, 
296                           bool isFullyConnected,
297                           int minimalComponentSize,
298                           LabelizeParameters<typename ImageType::PixelType> * param)
299   {
300     typename ImageType::Pointer working_image;
301     working_image = Labelize<ImageType>(input, BG, isFullyConnected, minimalComponentSize);
302     if (param->GetLabelsToRemove().size() != 0)
303       working_image = RemoveLabels<ImageType>(working_image, BG, param->GetLabelsToRemove());
304     working_image = KeepLabels<ImageType>(working_image, 
305                                           BG, FG, 
306                                           param->GetFirstKeep(), 
307                                           param->GetLastKeep(), 
308                                           param->GetUseLastKeep());
309     return working_image;
310   }
311   //--------------------------------------------------------------------
312
313
314   //--------------------------------------------------------------------
315   template<class ImageType>
316   typename ImageType::Pointer
317   ResizeImageLike(const ImageType * input,                       
318                   const itk::ImageBase<ImageType::ImageDimension> * like, 
319                   typename ImageType::PixelType backgroundValue) 
320   {
321     typedef CropLikeImageFilter<ImageType> CropFilterType;
322     typename CropFilterType::Pointer cropFilter = CropFilterType::New();
323     cropFilter->SetInput(input);
324     cropFilter->SetCropLikeImage(like);
325     cropFilter->SetBackgroundValue(backgroundValue);
326     cropFilter->Update();
327     return cropFilter->GetOutput();  
328   }
329   //--------------------------------------------------------------------
330
331
332   //--------------------------------------------------------------------
333   template<class MaskImageType>
334   typename MaskImageType::Pointer
335   SliceBySliceRelativePosition(const MaskImageType * input,
336                                const MaskImageType * object,
337                                int direction, 
338                                double threshold, 
339                                std::string orientation, 
340                                bool uniqueConnectedComponent, 
341                                double spacing, 
342                                bool autocropFlag, 
343                                bool singleObjectCCL) 
344   {
345     typedef SliceBySliceRelativePositionFilter<MaskImageType> SliceRelPosFilterType;
346     typename SliceRelPosFilterType::Pointer sliceRelPosFilter = SliceRelPosFilterType::New();
347     sliceRelPosFilter->VerboseStepFlagOff();
348     sliceRelPosFilter->WriteStepFlagOff();
349     sliceRelPosFilter->SetInput(input);
350     sliceRelPosFilter->SetInputObject(object);
351     sliceRelPosFilter->SetDirection(direction);
352     sliceRelPosFilter->SetFuzzyThreshold(threshold);
353     sliceRelPosFilter->AddOrientationTypeString(orientation);
354     sliceRelPosFilter->SetIntermediateSpacingFlag((spacing != -1));
355     sliceRelPosFilter->SetIntermediateSpacing(spacing);
356     sliceRelPosFilter->SetUniqueConnectedComponentBySlice(uniqueConnectedComponent);
357     sliceRelPosFilter->CCLSelectionFlagOff();
358     sliceRelPosFilter->SetUseASingleObjectConnectedComponentBySliceFlag(singleObjectCCL);
359     //    sliceRelPosFilter->SetInverseOrientationFlag(inverseflag); 
360     sliceRelPosFilter->SetAutoCropFlag(autocropFlag); 
361     sliceRelPosFilter->IgnoreEmptySliceObjectFlagOn();
362     sliceRelPosFilter->Update();
363     return sliceRelPosFilter->GetOutput();
364   }
365   //--------------------------------------------------------------------
366
367   //--------------------------------------------------------------------
368   template<class ImageType>
369   bool
370   FindExtremaPointInAGivenDirection(const ImageType * input, 
371                                     typename ImageType::PixelType bg, 
372                                     int direction, bool opposite, 
373                                     typename ImageType::PointType & point)
374   {
375     typename ImageType::PointType dummy;
376     return FindExtremaPointInAGivenDirection(input, bg, direction, 
377                                              opposite, dummy, 0, point);
378   }
379   //--------------------------------------------------------------------
380
381   //--------------------------------------------------------------------
382   template<class ImageType>
383   bool
384   FindExtremaPointInAGivenDirection(const ImageType * input, 
385                                     typename ImageType::PixelType bg, 
386                                     int direction, bool opposite, 
387                                     typename ImageType::PointType refpoint,
388                                     double distanceMax, 
389                                     typename ImageType::PointType & point)
390   {
391     /*
392       loop over input pixels, store the index in the fg that is max
393       according to the given direction. 
394     */    
395     typedef itk::ImageRegionConstIteratorWithIndex<ImageType> IteratorType;
396     IteratorType iter(input, input->GetLargestPossibleRegion());
397     iter.GoToBegin();
398     typename ImageType::IndexType max = input->GetLargestPossibleRegion().GetIndex();
399     if (opposite) max = max+input->GetLargestPossibleRegion().GetSize();
400     bool found=false;
401     while (!iter.IsAtEnd()) {
402       if (iter.Get() != bg) {
403         bool test = iter.GetIndex()[direction] >  max[direction];
404         if (opposite) test = !test;
405         if (test) {
406           typename ImageType::PointType p;
407           input->TransformIndexToPhysicalPoint(iter.GetIndex(), p);
408           if ((distanceMax==0) || (p.EuclideanDistanceTo(refpoint) < distanceMax)) {
409             max = iter.GetIndex();
410             found = true;
411           }
412         }
413       }
414       ++iter;
415     }
416     if (!found) return false;
417     input->TransformIndexToPhysicalPoint(max, point);
418     return true;
419   }
420   //--------------------------------------------------------------------
421
422
423   //--------------------------------------------------------------------
424   template<class ImageType>
425   typename ImageType::Pointer
426   CropImageRemoveGreaterThan(const ImageType * image, 
427                  int dim, double min, bool autoCrop,
428                  typename ImageType::PixelType BG) 
429   {
430     return CropImageAlongOneAxis<ImageType>(image, dim, 
431                                             image->GetOrigin()[dim], 
432                                             min,
433                                             autoCrop, BG);
434   }
435   //--------------------------------------------------------------------
436
437
438   //--------------------------------------------------------------------
439   template<class ImageType>
440   typename ImageType::Pointer
441   CropImageRemoveLowerThan(const ImageType * image, 
442                  int dim, double max, bool autoCrop,
443                  typename ImageType::PixelType BG) 
444   {
445     typename ImageType::PointType p;
446     image->TransformIndexToPhysicalPoint(image->GetLargestPossibleRegion().GetIndex()+
447                                          image->GetLargestPossibleRegion().GetSize(), p);
448     return CropImageAlongOneAxis<ImageType>(image, dim, max, p[dim], autoCrop, BG);
449   }
450   //--------------------------------------------------------------------
451
452
453   //--------------------------------------------------------------------
454   template<class ImageType>
455   typename ImageType::Pointer
456   CropImageAlongOneAxis(const ImageType * image, 
457                         int dim, double min, double max, 
458                         bool autoCrop, typename ImageType::PixelType BG) 
459   {
460     // Compute region size
461     typename ImageType::RegionType region;
462     typename ImageType::SizeType size = image->GetLargestPossibleRegion().GetSize();
463     typename ImageType::PointType p = image->GetOrigin();
464     p[dim] = min;
465     typename ImageType::IndexType start;
466     image->TransformPhysicalPointToIndex(p, start);
467     p[dim] = max;
468     typename ImageType::IndexType end;
469     image->TransformPhysicalPointToIndex(p, end);
470     size[dim] = abs(end[dim]-start[dim]);
471     region.SetIndex(start);
472     region.SetSize(size);
473   
474     // Perform Crop
475     typedef itk::RegionOfInterestImageFilter<ImageType, ImageType> CropFilterType;
476     typename CropFilterType::Pointer cropFilter = CropFilterType::New();
477     cropFilter->SetInput(image);
478     cropFilter->SetRegionOfInterest(region);
479     cropFilter->Update();
480     typename ImageType::Pointer result = cropFilter->GetOutput();
481   
482     // Auto Crop
483     if (autoCrop) {
484       result = AutoCrop<ImageType>(result, BG);
485     }
486     return result;
487   }
488   //--------------------------------------------------------------------
489
490
491   //--------------------------------------------------------------------
492   template<class ImageType>
493   void
494   ComputeCentroids(const ImageType * image, 
495                    typename ImageType::PixelType BG, 
496                    std::vector<typename ImageType::PointType> & centroids) 
497   {
498     typedef long LabelType;
499     static const unsigned int Dim = ImageType::ImageDimension;
500     typedef itk::ShapeLabelObject< LabelType, Dim > LabelObjectType;
501     typedef itk::LabelMap< LabelObjectType > LabelMapType;
502     typedef itk::LabelImageToLabelMapFilter<ImageType, LabelMapType> ImageToMapFilterType;
503     typename ImageToMapFilterType::Pointer imageToLabelFilter = ImageToMapFilterType::New(); 
504     typedef itk::ShapeLabelMapFilter<LabelMapType, ImageType> ShapeFilterType; 
505     typename ShapeFilterType::Pointer statFilter = ShapeFilterType::New();
506     imageToLabelFilter->SetBackgroundValue(BG);
507     imageToLabelFilter->SetInput(image);
508     statFilter->SetInput(imageToLabelFilter->GetOutput());
509     statFilter->Update();
510     typename LabelMapType::Pointer labelMap = statFilter->GetOutput();
511
512     centroids.clear();
513     typename ImageType::PointType dummy;
514     centroids.push_back(dummy); // label 0 -> no centroid, use dummy point for BG 
515     //DS FIXME (not useful ! to change ..)
516     for(uint i=0; i<labelMap->GetNumberOfLabelObjects(); i++) {
517       int label = labelMap->GetLabels()[i];
518       centroids.push_back(labelMap->GetLabelObject(label)->GetCentroid());
519     } 
520   }
521   //--------------------------------------------------------------------
522
523
524   //--------------------------------------------------------------------
525   template<class ImageType, class LabelType>
526   typename itk::LabelMap< itk::ShapeLabelObject<LabelType, ImageType::ImageDimension> >::Pointer
527   ComputeLabelMap(const ImageType * image, 
528                   typename ImageType::PixelType BG, 
529                   bool computePerimeterFlag) 
530   {
531     static const unsigned int Dim = ImageType::ImageDimension;
532     typedef itk::ShapeLabelObject< LabelType, Dim > LabelObjectType;
533     typedef itk::LabelMap< LabelObjectType > LabelMapType;
534     typedef itk::LabelImageToLabelMapFilter<ImageType, LabelMapType> ImageToMapFilterType;
535     typename ImageToMapFilterType::Pointer imageToLabelFilter = ImageToMapFilterType::New(); 
536     typedef itk::ShapeLabelMapFilter<LabelMapType, ImageType> ShapeFilterType; 
537     typename ShapeFilterType::Pointer statFilter = ShapeFilterType::New();
538     imageToLabelFilter->SetBackgroundValue(BG);
539     imageToLabelFilter->SetInput(image);
540     statFilter->SetInput(imageToLabelFilter->GetOutput());
541     statFilter->SetComputePerimeter(computePerimeterFlag);
542     statFilter->Update();
543     return statFilter->GetOutput();
544   }
545   //--------------------------------------------------------------------
546
547
548   //--------------------------------------------------------------------
549   template<class ImageType>
550   void
551   ComputeCentroids2(const ImageType * image, 
552                    typename ImageType::PixelType BG, 
553                    std::vector<typename ImageType::PointType> & centroids) 
554   {
555     typedef long LabelType;
556     static const unsigned int Dim = ImageType::ImageDimension;
557     typedef itk::ShapeLabelObject< LabelType, Dim > LabelObjectType;
558     typedef itk::LabelMap< LabelObjectType > LabelMapType;
559     typedef itk::LabelImageToLabelMapFilter<ImageType, LabelMapType> ImageToMapFilterType;
560     typename ImageToMapFilterType::Pointer imageToLabelFilter = ImageToMapFilterType::New(); 
561     typedef itk::ShapeLabelMapFilter<LabelMapType, ImageType> ShapeFilterType; 
562     typename ShapeFilterType::Pointer statFilter = ShapeFilterType::New();
563     imageToLabelFilter->SetBackgroundValue(BG);
564     imageToLabelFilter->SetInput(image);
565     statFilter->SetInput(imageToLabelFilter->GetOutput());
566     statFilter->Update();
567     typename LabelMapType::Pointer labelMap = statFilter->GetOutput();
568
569     centroids.clear();
570     typename ImageType::PointType dummy;
571     centroids.push_back(dummy); // label 0 -> no centroid, use dummy point
572     for(uint i=1; i<labelMap->GetNumberOfLabelObjects()+1; i++) {
573       centroids.push_back(labelMap->GetLabelObject(i)->GetCentroid());
574     } 
575     
576     for(uint i=1; i<labelMap->GetNumberOfLabelObjects()+1; i++) {
577       DD(labelMap->GetLabelObject(i)->GetBinaryPrincipalAxes());
578       DD(labelMap->GetLabelObject(i)->GetBinaryFlatness());
579       DD(labelMap->GetLabelObject(i)->GetRoundness ());      
580
581       // search for the point on the boundary alog PA
582
583     }
584
585   }
586   //--------------------------------------------------------------------
587
588
589   //--------------------------------------------------------------------
590   template<class ImageType>
591   void
592   PointsUtils<ImageType>::Convert2DTo3D(const PointType2D & p2D, 
593                                         const ImageType * image, 
594                                         const int slice, 
595                                         PointType3D & p3D)  
596   {
597     IndexType3D index3D;
598     index3D[0] = index3D[1] = 0;
599     index3D[2] = image->GetLargestPossibleRegion().GetIndex()[2]+slice;
600     image->TransformIndexToPhysicalPoint(index3D, p3D);
601     p3D[0] = p2D[0]; 
602     p3D[1] = p2D[1];
603     //  p3D[2] = p[2];//(image->GetLargestPossibleRegion().GetIndex()[2]+slice)*image->GetSpacing()[2] 
604     //    + image->GetOrigin()[2];
605   }
606   //--------------------------------------------------------------------
607
608
609   //--------------------------------------------------------------------
610   template<class ImageType>
611   void 
612   PointsUtils<ImageType>::Convert2DMapTo3DList(const MapPoint2DType & map, 
613                                             const ImageType * image, 
614                                             VectorPoint3DType & list)
615   {
616     typename MapPoint2DType::const_iterator iter = map.begin();
617     while (iter != map.end()) {
618       PointType3D p;
619       Convert2DTo3D(iter->second, image, iter->first, p);
620       list.push_back(p);
621       ++iter;
622     }
623   }
624   //--------------------------------------------------------------------
625
626
627   //--------------------------------------------------------------------
628   template<class ImageType>
629   void 
630   PointsUtils<ImageType>::Convert2DListTo3DList(const VectorPoint2DType & p2D, 
631                                                 int slice,
632                                                 const ImageType * image, 
633                                                 VectorPoint3DType & list) 
634   {
635     for(uint i=0; i<p2D.size(); i++) {
636       PointType3D p;
637       Convert2DTo3D(p2D[i], image, slice, p);
638       list.push_back(p);
639     }
640   }
641   //--------------------------------------------------------------------
642
643
644   //--------------------------------------------------------------------
645   template<class ImageType>
646   void 
647   WriteListOfLandmarks(std::vector<typename ImageType::PointType> points, 
648                        std::string filename)
649   {
650     std::ofstream os; 
651     openFileForWriting(os, filename); 
652     os << "LANDMARKS1" << std::endl;  
653     for(uint i=0; i<points.size(); i++) {
654       const typename ImageType::PointType & p = points[i];
655       // Write it in the file
656       os << i << " " << p[0] << " " << p[1] << " " << p[2] << " 0 0 " << std::endl;
657     }
658     os.close();
659   }
660   //--------------------------------------------------------------------
661
662
663   //--------------------------------------------------------------------
664   template<class ImageType>
665   typename ImageType::Pointer 
666   Dilate(const ImageType * image, double radiusInMM,               
667          typename ImageType::PixelType BG,
668          typename ImageType::PixelType FG,  
669          bool extendSupport)
670   {
671     typename ImageType::SizeType r;
672     for(uint i=0; i<ImageType::ImageDimension; i++) 
673       r[i] = (uint)lrint(radiusInMM/image->GetSpacing()[i]);
674     return Dilate<ImageType>(image, r, BG, FG, extendSupport);
675   }
676   //--------------------------------------------------------------------
677
678
679   //--------------------------------------------------------------------
680   template<class ImageType>
681   typename ImageType::Pointer 
682   Dilate(const ImageType * image, typename ImageType::PointType radiusInMM, 
683          typename ImageType::PixelType BG, 
684          typename ImageType::PixelType FG, 
685          bool extendSupport)
686   {
687     typename ImageType::SizeType r;
688     for(uint i=0; i<ImageType::ImageDimension; i++) 
689       r[i] = (uint)lrint(radiusInMM[i]/image->GetSpacing()[i]);
690     return Dilate<ImageType>(image, r, BG, FG, extendSupport);
691   }
692   //--------------------------------------------------------------------
693
694
695   //--------------------------------------------------------------------
696   template<class ImageType>
697   typename ImageType::Pointer 
698   Dilate(const ImageType * image, typename ImageType::SizeType radius, 
699          typename ImageType::PixelType BG, 
700          typename ImageType::PixelType FG, 
701          bool extendSupport)
702   {
703     // Create kernel for dilatation
704     typedef itk::BinaryBallStructuringElement<typename ImageType::PixelType, 
705                                               ImageType::ImageDimension> KernelType;
706     KernelType structuringElement;
707     structuringElement.SetRadius(radius);
708     structuringElement.CreateStructuringElement();
709
710     typename ImageType::Pointer output;
711     if (extendSupport) {
712       typedef itk::ConstantPadImageFilter<ImageType, ImageType> PadFilterType;
713       typename PadFilterType::Pointer padFilter = PadFilterType::New();
714       padFilter->SetInput(image);
715       typename ImageType::SizeType lower;
716       typename ImageType::SizeType upper;
717       for(uint i=0; i<3; i++) {
718         lower[i] = upper[i] = 2*(radius[i]+1);
719       }
720       padFilter->SetPadLowerBound(lower);
721       padFilter->SetPadUpperBound(upper);
722       padFilter->Update();
723       output = padFilter->GetOutput();
724     }
725
726     // Dilate  filter
727     typedef itk::BinaryDilateImageFilter<ImageType, ImageType , KernelType> DilateFilterType;
728     typename DilateFilterType::Pointer dilateFilter = DilateFilterType::New();
729     dilateFilter->SetBackgroundValue(BG);
730     dilateFilter->SetForegroundValue(FG);
731     dilateFilter->SetBoundaryToForeground(false);
732     dilateFilter->SetKernel(structuringElement);
733     dilateFilter->SetInput(output);
734     dilateFilter->Update();
735     return dilateFilter->GetOutput();
736   }
737   //--------------------------------------------------------------------
738
739
740   //--------------------------------------------------------------------
741   template<class ImageType>
742   typename ImageType::Pointer 
743   Opening(const ImageType * image, typename ImageType::SizeType radius,
744          typename ImageType::PixelType BG,
745          typename ImageType::PixelType FG)
746   {
747     // Kernel 
748     typedef itk::BinaryBallStructuringElement<typename ImageType::PixelType, 
749                                               ImageType::ImageDimension> KernelType;    
750     KernelType structuringElement;
751     structuringElement.SetRadius(radius);
752     structuringElement.CreateStructuringElement();
753     
754     // Filter
755     typedef itk::BinaryMorphologicalOpeningImageFilter<ImageType, ImageType , KernelType> OpeningFilterType;
756     typename OpeningFilterType::Pointer open = OpeningFilterType::New();
757     open->SetInput(image);
758     open->SetBackgroundValue(BG);
759     open->SetForegroundValue(FG);
760     open->SetKernel(structuringElement);
761     open->Update();
762     return open->GetOutput();
763   }
764   //--------------------------------------------------------------------
765
766
767
768   //--------------------------------------------------------------------
769   template<class ValueType, class VectorType>
770   void ConvertOption(std::string optionName, uint given, 
771                      ValueType * values, VectorType & p, 
772                      uint dim, bool required) 
773   {
774     if (required && (given == 0)) {
775       clitkExceptionMacro("The option --" << optionName << " must be set and have 1 or " 
776                           << dim << " values.");
777     }
778     if (given == 1) {
779       for(uint i=0; i<dim; i++) p[i] = values[0];
780       return;
781     }
782     if (given == dim) {
783       for(uint i=0; i<dim; i++) p[i] = values[i];
784       return;
785     }
786     if (given == 0) return;
787     clitkExceptionMacro("The option --" << optionName << " must have 1 or " 
788                         << dim << " values.");
789   }
790   //--------------------------------------------------------------------
791
792
793   //--------------------------------------------------------------------
794   /*
795     http://www.gamedev.net/community/forums/topic.asp?topic_id=542870
796     Assuming the points are (Ax,Ay) (Bx,By) and (Cx,Cy), you need to compute:
797     (Bx - Ax) * (Cy - Ay) - (By - Ay) * (Cx - Ax)
798     This will equal zero if the point C is on the line formed by
799     points A and B, and will have a different sign depending on the
800     side. Which side this is depends on the orientation of your (x,y)
801     coordinates, but you can plug test values for A,B and C into this
802     formula to determine whether negative values are to the left or to
803     the right.
804     => to accelerate, start with formula, when change sign -> stop and fill
805
806     offsetToKeep = is used to determine which side of the line we
807     keep. The point along the mainDirection but 'offsetToKeep' mm away
808     is kept.
809   
810   */
811   template<class ImageType>
812   void 
813   SliceBySliceSetBackgroundFromLineSeparation(ImageType * input, 
814                                               std::vector<typename ImageType::PointType> & lA, 
815                                               std::vector<typename ImageType::PointType> & lB, 
816                                               typename ImageType::PixelType BG, 
817                                               int mainDirection, 
818                                               double offsetToKeep)
819   {
820     typedef itk::ImageSliceIteratorWithIndex<ImageType> SliceIteratorType;
821     SliceIteratorType siter = SliceIteratorType(input, 
822                                                 input->GetLargestPossibleRegion());
823     siter.SetFirstDirection(0);
824     siter.SetSecondDirection(1);
825     siter.GoToBegin();
826     uint i=0;
827     typename ImageType::PointType A;
828     typename ImageType::PointType B;
829     typename ImageType::PointType C;
830     assert(lA.size() == lB.size());
831     while ((i<lA.size()) && (!siter.IsAtEnd())) {
832       // Check that the current slice correspond to the current point
833       input->TransformIndexToPhysicalPoint(siter.GetIndex(), C);
834       // DD(C);
835       // DD(i);
836       // DD(lA[i]);
837       if ((fabs(C[2] - lA[i][2]))>0.01) { // is !equal with a tolerance of 0.01 mm
838       }
839       else {
840         // Define A,B,C points
841         A = lA[i];
842         B = lB[i];
843         C = A;
844         // DD(A);
845         // DD(B);
846         // DD(C);
847       
848         // Check that the line is not a point (A=B)
849         bool p = (A[0] == B[0]) && (A[1] == B[1]);
850       
851         if (!p) {
852           C[mainDirection] += offsetToKeep; // I know I must keep this point
853           double s = (B[0] - A[0]) * (C[1] - A[1]) - (B[1] - A[1]) * (C[0] - A[0]);
854           bool isPositive = s<0;
855           while (!siter.IsAtEndOfSlice()) {
856             while (!siter.IsAtEndOfLine()) {
857               // Very slow, I know ... but image should be very small
858               input->TransformIndexToPhysicalPoint(siter.GetIndex(), C);
859               double s = (B[0] - A[0]) * (C[1] - A[1]) - (B[1] - A[1]) * (C[0] - A[0]);
860               if (s == 0) siter.Set(BG); // on the line, we decide to remove
861               if (isPositive) {
862                 if (s > 0) siter.Set(BG);
863               }
864               else {
865                 if (s < 0) siter.Set(BG); 
866               }
867               ++siter;
868             }
869             siter.NextLine();
870           } // end loop slice
871         }      
872
873         ++i;
874       } // End of current slice
875       siter.NextSlice();
876     }
877   }                                                   
878   //--------------------------------------------------------------------
879
880
881   //--------------------------------------------------------------------
882   template<class ImageType>
883   void 
884   AndNot(ImageType * input, 
885          const ImageType * object, 
886          typename ImageType::PixelType BG)
887   {
888     typename ImageType::Pointer o;
889     bool resized=false;
890     if (!clitk::HaveSameSizeAndSpacing<ImageType, ImageType>(input, object)) {
891       o = clitk::ResizeImageLike<ImageType>(object, input, BG);
892       resized = true;
893     }
894
895     typedef clitk::BooleanOperatorLabelImageFilter<ImageType> BoolFilterType;
896     typename BoolFilterType::Pointer boolFilter = BoolFilterType::New(); 
897     boolFilter->InPlaceOn();
898     boolFilter->SetInput1(input);
899     if (resized) boolFilter->SetInput2(o);  
900     else boolFilter->SetInput2(object);
901     boolFilter->SetBackgroundValue1(BG);
902     boolFilter->SetBackgroundValue2(BG);
903     boolFilter->SetOperationType(BoolFilterType::AndNot);
904     boolFilter->Update();
905   }
906   //--------------------------------------------------------------------
907
908
909   //--------------------------------------------------------------------
910   template<class ImageType>
911   void 
912   And(ImageType * input, 
913       const ImageType * object, 
914       typename ImageType::PixelType BG)
915   {
916     typename ImageType::Pointer o;
917     bool resized=false;
918     if (!clitk::HaveSameSizeAndSpacing<ImageType, ImageType>(input, object)) {
919       o = clitk::ResizeImageLike<ImageType>(object, input, BG);
920       resized = true;
921     }
922
923     typedef clitk::BooleanOperatorLabelImageFilter<ImageType> BoolFilterType;
924     typename BoolFilterType::Pointer boolFilter = BoolFilterType::New(); 
925     boolFilter->InPlaceOn();
926     boolFilter->SetInput1(input);
927     if (resized) boolFilter->SetInput2(o);  
928     else boolFilter->SetInput2(object);
929     boolFilter->SetBackgroundValue1(BG);
930     boolFilter->SetBackgroundValue2(BG);
931     boolFilter->SetOperationType(BoolFilterType::And);
932     boolFilter->Update();
933   }
934   //--------------------------------------------------------------------
935
936
937   //--------------------------------------------------------------------
938   template<class ImageType>
939   typename ImageType::Pointer
940   Binarize(const ImageType * input, 
941            typename ImageType::PixelType lower, 
942            typename ImageType::PixelType upper, 
943            typename ImageType::PixelType BG,
944            typename ImageType::PixelType FG) 
945   {
946     typedef itk::BinaryThresholdImageFilter<ImageType, ImageType> BinaryThresholdFilterType;
947     typename BinaryThresholdFilterType::Pointer binarizeFilter = BinaryThresholdFilterType::New();
948     binarizeFilter->SetInput(input);
949     binarizeFilter->InPlaceOff();
950     binarizeFilter->SetLowerThreshold(lower);
951     binarizeFilter->SetUpperThreshold(upper);
952     binarizeFilter->SetInsideValue(FG);
953     binarizeFilter->SetOutsideValue(BG);
954     binarizeFilter->Update();
955     return binarizeFilter->GetOutput();
956   }
957   //--------------------------------------------------------------------
958
959
960   //--------------------------------------------------------------------
961   template<class ImageType>
962   void
963   GetMinMaxPointPosition(const ImageType * input, 
964                          typename ImageType::PointType & min,
965                          typename ImageType::PointType & max) 
966   {
967     typename ImageType::IndexType index = input->GetLargestPossibleRegion().GetIndex();
968     input->TransformIndexToPhysicalPoint(index, min);
969     index = index+input->GetLargestPossibleRegion().GetSize();
970     input->TransformIndexToPhysicalPoint(index, max);
971   }
972   //--------------------------------------------------------------------
973
974
975   //--------------------------------------------------------------------
976   template<class ImageType>
977   typename ImageType::PointType
978   FindExtremaPointInAGivenLine(const ImageType * input, 
979                                int dimension, 
980                                bool inverse, 
981                                typename ImageType::PointType p, 
982                                typename ImageType::PixelType BG, 
983                                double distanceMax) 
984   {
985     // Which direction ?  Increasing or decreasing.
986     int d=1;
987     if (inverse) d=-1;
988   
989     // Transform to pixel index
990     typename ImageType::IndexType index;
991     input->TransformPhysicalPointToIndex(p, index);
992
993     // Loop while inside the mask;
994     while (input->GetPixel(index) != BG) {
995       index[dimension] += d;
996     }
997
998     // Transform back to Physical Units
999     typename ImageType::PointType result;
1000     input->TransformIndexToPhysicalPoint(index, result);
1001
1002     // Check that is is not too far away
1003     double distance = p.EuclideanDistanceTo(result);
1004     if (distance > distanceMax) {
1005       result = p; // Get back to initial value
1006     }
1007
1008     return result;
1009   }
1010   //--------------------------------------------------------------------
1011
1012
1013   //--------------------------------------------------------------------
1014   template<class PointType>
1015   bool
1016   IsOnTheSameLineSide(PointType C, PointType A, PointType B, PointType like) 
1017   {
1018     // Look at the position of point 'like' according to the AB line
1019     double s = (B[0] - A[0]) * (like[1] - A[1]) - (B[1] - A[1]) * (like[0] - A[0]);
1020     bool negative = s<0;
1021   
1022     // Look the C position
1023     s = (B[0] - A[0]) * (C[1] - A[1]) - (B[1] - A[1]) * (C[0] - A[0]);
1024
1025     if (negative && (s<=0)) return true;
1026     if (!negative && (s>=0)) return true;
1027     return false;
1028   }
1029   //--------------------------------------------------------------------
1030
1031
1032   //--------------------------------------------------------------------
1033   /* Consider an input object, for each slice, find the extrema
1034      position according to a given direction and build a line segment
1035      passing throught this point in a given direction.  Output is a
1036      vector of line (from point A to B), for each slice;
1037    */
1038   template<class ImageType>
1039   void 
1040   SliceBySliceBuildLineSegmentAccordingToExtremaPosition(const ImageType * input, 
1041                                                          typename ImageType::PixelType BG, 
1042                                                          int sliceDimension, 
1043                                                          int extremaDirection, 
1044                                                          bool extremaOppositeFlag, 
1045                                                          int lineDirection,
1046                                                          double margin,
1047                                                          std::vector<typename ImageType::PointType> & A, 
1048                                                          std::vector<typename ImageType::PointType> & B)
1049   {
1050     // Type of a slice
1051     typedef typename itk::Image<typename ImageType::PixelType, ImageType::ImageDimension-1> SliceType;
1052     
1053     // Build the list of slices
1054     std::vector<typename SliceType::Pointer> slices;
1055     clitk::ExtractSlices<ImageType>(input, sliceDimension, slices);
1056
1057     // Build the list of 2D points
1058     std::map<int, typename SliceType::PointType> position2D;
1059     for(uint i=0; i<slices.size(); i++) {
1060       typename SliceType::PointType p;
1061       bool found = 
1062         clitk::FindExtremaPointInAGivenDirection<SliceType>(slices[i], BG, 
1063                                                             extremaDirection, extremaOppositeFlag, p);
1064       if (found) {
1065         position2D[i] = p;
1066       }    
1067     }
1068     
1069     // Convert 2D points in slice into 3D points
1070     clitk::PointsUtils<ImageType>::Convert2DMapTo3DList(position2D, input, A);
1071     
1072     // Create additional point just right to the previous ones, on the
1073     // given lineDirection, in order to create a horizontal/vertical line.
1074     for(uint i=0; i<A.size(); i++) {
1075       typename ImageType::PointType p = A[i];
1076       p[lineDirection] += 10;
1077       B.push_back(p);
1078       // Margins ?
1079       A[i][1] += margin;
1080       B[i][1] += margin;
1081     }
1082
1083   }
1084   //--------------------------------------------------------------------
1085
1086
1087   //--------------------------------------------------------------------
1088   template<class ImageType>
1089   typename ImageType::Pointer
1090   SliceBySliceKeepMainCCL(const ImageType * input, 
1091                           typename ImageType::PixelType BG,
1092                           typename ImageType::PixelType FG)  {
1093     
1094     // Extract slices
1095     const int d = ImageType::ImageDimension-1;
1096     typedef typename itk::Image<typename ImageType::PixelType, d> SliceType;
1097     std::vector<typename SliceType::Pointer> slices;
1098     clitk::ExtractSlices<ImageType>(input, d, slices);
1099     
1100     // Labelize and keep the main one
1101     std::vector<typename SliceType::Pointer> o;
1102     for(uint i=0; i<slices.size(); i++) {
1103       o.push_back(clitk::Labelize<SliceType>(slices[i], BG, false, 1));
1104       o[i] = clitk::KeepLabels<SliceType>(o[i], BG, FG, 1, 1, true);
1105     }
1106     
1107     // Join slices
1108     typename ImageType::Pointer output;
1109     output = clitk::JoinSlices<ImageType>(o, input, d);
1110     return output;
1111   }
1112   //--------------------------------------------------------------------
1113
1114
1115   //--------------------------------------------------------------------
1116   template<class ImageType>
1117   typename ImageType::Pointer
1118   Clone(const ImageType * input) {
1119     typedef itk::ImageDuplicator<ImageType> DuplicatorType;
1120     typename DuplicatorType::Pointer duplicator = DuplicatorType::New();
1121     duplicator->SetInputImage(input);
1122     duplicator->Update();
1123     return duplicator->GetOutput();
1124   }
1125   //--------------------------------------------------------------------
1126
1127
1128   //--------------------------------------------------------------------
1129   /* Consider an input object, start at A, for each slice (dim1): 
1130      - compute the intersection between the AB line and the current slice
1131      - remove what is at lower or greater according to dim2 of this point
1132      - stop at B
1133   */
1134   template<class ImageType>
1135   typename ImageType::Pointer
1136   SliceBySliceSetBackgroundFromSingleLine(const ImageType * input, 
1137                                           typename ImageType::PixelType BG, 
1138                                           typename ImageType::PointType & A, 
1139                                           typename ImageType::PointType & B, 
1140                                           int dim1, int dim2, bool removeLowerPartFlag)
1141     
1142   {
1143     // Extract slices
1144     typedef typename itk::Image<typename ImageType::PixelType, ImageType::ImageDimension-1> SliceType;
1145     typedef typename SliceType::Pointer SlicePointer;
1146     std::vector<SlicePointer> slices;
1147     clitk::ExtractSlices<ImageType>(input, dim1, slices);
1148
1149     // Start at slice that contains A, and stop at B
1150     typename ImageType::IndexType Ap;
1151     typename ImageType::IndexType Bp;
1152     input->TransformPhysicalPointToIndex(A, Ap);
1153     input->TransformPhysicalPointToIndex(B, Bp);
1154     
1155     // Determine slice largest region
1156     typename SliceType::RegionType region = slices[0]->GetLargestPossibleRegion();
1157     typename SliceType::SizeType size = region.GetSize();
1158     typename SliceType::IndexType index = region.GetIndex();
1159
1160     // Line slope
1161     double a = (Bp[dim2]-Ap[dim2])/(Bp[dim1]-Ap[dim1]);
1162     double b = Ap[dim2];
1163
1164     // Loop from slice A to slice B
1165     for(uint i=0; i<(Bp[dim1]-Ap[dim1]); i++) {
1166       // Compute intersection between line AB and current slice for the dim2
1167       double p = a*i+b;
1168       // Change region (lower than dim2)
1169       if (removeLowerPartFlag) {
1170         size[dim2] = p-Ap[dim2];
1171       }
1172       else {
1173         size[dim2] = slices[0]->GetLargestPossibleRegion().GetSize()[dim2]-p;
1174         index[dim2] = p;
1175       }
1176       region.SetSize(size);
1177       region.SetIndex(index);
1178       // Fill region with BG (simple region iterator)
1179       FillRegionWithValue<SliceType>(slices[i+Ap[dim1]], BG, region);
1180       /*
1181       typedef itk::ImageRegionIterator<SliceType> IteratorType;
1182       IteratorType iter(slices[i+Ap[dim1]], region);
1183       iter.GoToBegin();
1184       while (!iter.IsAtEnd()) {
1185         iter.Set(BG);
1186         ++iter;
1187       }
1188       */
1189       // Loop
1190     }
1191     
1192     // Merge slices
1193     typename ImageType::Pointer output;
1194     output = clitk::JoinSlices<ImageType>(slices, input, dim1);
1195     return output;
1196   }
1197   //--------------------------------------------------------------------
1198
1199   //--------------------------------------------------------------------
1200   template<class ImageType>
1201   void
1202   FillRegionWithValue(ImageType * input, typename ImageType::PixelType value, typename ImageType::RegionType & region)
1203   {
1204     typedef itk::ImageRegionIterator<ImageType> IteratorType;
1205     IteratorType iter(input, region);
1206     iter.GoToBegin();
1207     while (!iter.IsAtEnd()) {
1208       iter.Set(value);
1209       ++iter;
1210     }    
1211   }
1212   //--------------------------------------------------------------------
1213
1214 } // end of namespace
1215