]> Creatis software - clitk.git/blob - itk/clitkSegmentationUtils.txx
Add LabelizeAndCountNumberOfObjects
[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://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 "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
36 //--------------------------------------------------------------------
37 template<class ImageType>
38 void clitk::ComputeBBFromImageRegion(typename ImageType::Pointer image, 
39                                      typename ImageType::RegionType region,
40                                      typename itk::BoundingBox<unsigned long, 
41                                                                ImageType::ImageDimension>::Pointer bb) {
42   typedef typename ImageType::IndexType IndexType;
43   IndexType firstIndex;
44   IndexType lastIndex;
45   for(unsigned int i=0; i<image->GetImageDimension(); i++) {
46     firstIndex[i] = region.GetIndex()[i];
47     lastIndex[i] = firstIndex[i]+region.GetSize()[i];
48   }
49
50   typedef itk::BoundingBox<unsigned long, 
51                            ImageType::ImageDimension> BBType;
52   typedef typename BBType::PointType PointType;
53   PointType lastPoint;
54   PointType firstPoint;
55   image->TransformIndexToPhysicalPoint(firstIndex, firstPoint);
56   image->TransformIndexToPhysicalPoint(lastIndex, lastPoint);
57
58   bb->SetMaximum(lastPoint);
59   bb->SetMinimum(firstPoint);
60 }
61 //--------------------------------------------------------------------
62
63
64 //--------------------------------------------------------------------
65 template<int Dimension>
66 void clitk::ComputeBBIntersection(typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbo, 
67                                   typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbi1, 
68                                   typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbi2) {
69
70   typedef itk::BoundingBox<unsigned long, Dimension> BBType;
71   typedef typename BBType::PointType PointType;
72   PointType lastPoint;
73   PointType firstPoint;
74
75   for(unsigned int i=0; i<Dimension; i++) {
76     firstPoint[i] = std::max(bbi1->GetMinimum()[i], 
77                              bbi2->GetMinimum()[i]);
78     lastPoint[i] = std::min(bbi1->GetMaximum()[i], 
79                             bbi2->GetMaximum()[i]);
80   }
81
82   bbo->SetMaximum(lastPoint);
83   bbo->SetMinimum(firstPoint);
84 }
85 //--------------------------------------------------------------------
86
87
88 //--------------------------------------------------------------------
89 template<class ImageType>
90 void clitk::ComputeRegionFromBB(typename ImageType::Pointer image, 
91                                 const typename itk::BoundingBox<unsigned long, 
92                                 ImageType::ImageDimension>::Pointer bb, 
93                                 typename ImageType::RegionType & region) {
94   // Types
95   typedef typename ImageType::IndexType  IndexType;
96   typedef typename ImageType::PointType  PointType;
97   typedef typename ImageType::RegionType RegionType;
98   typedef typename ImageType::SizeType   SizeType;
99
100   // Region starting point
101   IndexType regionStart;
102   PointType start = bb->GetMinimum();
103   image->TransformPhysicalPointToIndex(start, regionStart);
104     
105   // Region size
106   SizeType regionSize;
107   PointType maxs = bb->GetMaximum();
108   PointType mins = bb->GetMinimum();
109   for(unsigned int i=0; i<ImageType::ImageDimension; i++) {
110     // DD(maxs[i]);
111     // DD(mins[i]);
112     // DD((maxs[i] - mins[i])/image->GetSpacing()[i]);
113     regionSize[i] = lrint((maxs[i] - mins[i])/image->GetSpacing()[i]);
114     // DD(regionSize[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 clitk::SetBackground(const ImageType * input, 
127                      const TMaskImageType * mask, 
128                      typename TMaskImageType::PixelType maskBG,
129                      typename ImageType::PixelType outValue, 
130                      bool inPlace) {
131   typedef clitk::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 clitk::GetNumberOfConnectedComponentLabels(typename ImageType::Pointer 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 clitk::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   // DD(relabelFilter->GetNumberOfObjects());
196   // DD(relabelFilter->GetOriginalNumberOfObjects());
197   // DD(relabelFilter->GetSizeOfObjectsInPhysicalUnits()[0]);
198
199   // Return result
200   typename ImageType::Pointer output = relabelFilter->GetOutput();
201   return output;
202 }
203 //--------------------------------------------------------------------
204
205
206 //--------------------------------------------------------------------
207 /*
208   Warning : in this cas, we consider outputType like inputType, not
209   InternalImageType. Be sure it fits.
210  */
211 template<class ImageType>
212 typename ImageType::Pointer
213 clitk::LabelizeAndCountNumberOfObjects(const ImageType * input, 
214                                        typename ImageType::PixelType BG, 
215                                        bool isFullyConnected, 
216                                        int minimalComponentSize, 
217                                        int & nb) {
218   // InternalImageType for storing large number of component
219   typedef itk::Image<int, ImageType::ImageDimension> InternalImageType;
220   
221   // Connected Component label 
222   typedef itk::ConnectedComponentImageFilter<ImageType, InternalImageType> ConnectFilterType;
223   typename ConnectFilterType::Pointer connectFilter = ConnectFilterType::New();
224   //  connectFilter->ReleaseDataFlagOn(); 
225   connectFilter->SetInput(input);
226   connectFilter->SetBackgroundValue(BG);
227   connectFilter->SetFullyConnected(isFullyConnected);
228   
229   // Sort by size and remove too small area.
230   typedef itk::RelabelComponentImageFilter<InternalImageType, ImageType> RelabelFilterType;
231   typename RelabelFilterType::Pointer relabelFilter = RelabelFilterType::New();
232   //  relabelFilter->ReleaseDataFlagOn(); // if yes, fail when ExplosionControlledThresholdConnectedImageFilter ???
233   relabelFilter->SetInput(connectFilter->GetOutput());
234   relabelFilter->SetMinimumObjectSize(minimalComponentSize);
235   relabelFilter->Update();
236
237   nb = relabelFilter->GetNumberOfObjects();
238   // DD(relabelFilter->GetNumberOfObjects());
239   // DD(relabelFilter->GetOriginalNumberOfObjects());
240   // DD(relabelFilter->GetSizeOfObjectsInPhysicalUnits()[0]);
241
242   // Return result
243   typename ImageType::Pointer output = relabelFilter->GetOutput();
244   return output;
245 }
246 //--------------------------------------------------------------------
247
248
249
250 //--------------------------------------------------------------------
251 template<class ImageType>
252 typename ImageType::Pointer
253 clitk::RemoveLabels(typename ImageType::Pointer input, 
254                     typename ImageType::PixelType BG,
255                     std::vector<typename ImageType::PixelType> & labelsToRemove) {
256   typename ImageType::Pointer working_image = input;
257   for (unsigned int i=0; i <labelsToRemove.size(); i++) {
258     typedef clitk::SetBackgroundImageFilter<ImageType, ImageType> SetBackgroundImageFilterType;
259     typename SetBackgroundImageFilterType::Pointer setBackgroundFilter = SetBackgroundImageFilterType::New();
260     setBackgroundFilter->SetInput(input);
261     setBackgroundFilter->SetInput2(input);
262     setBackgroundFilter->SetMaskValue(labelsToRemove[i]);
263     setBackgroundFilter->SetOutsideValue(BG);
264     setBackgroundFilter->Update();
265     working_image = setBackgroundFilter->GetOutput();
266   }
267   return working_image;
268 }
269 //--------------------------------------------------------------------
270
271
272 //--------------------------------------------------------------------
273 template<class ImageType>
274 typename ImageType::Pointer
275 clitk::KeepLabels(const ImageType * input, 
276                   typename ImageType::PixelType BG, 
277                   typename ImageType::PixelType FG, 
278                   typename ImageType::PixelType firstKeep, 
279                   typename ImageType::PixelType lastKeep, 
280                   bool useLastKeep) {
281   typedef itk::BinaryThresholdImageFilter<ImageType, ImageType> BinarizeFilterType; 
282   typename BinarizeFilterType::Pointer binarizeFilter = BinarizeFilterType::New();
283   binarizeFilter->SetInput(input);
284   binarizeFilter->SetLowerThreshold(firstKeep);
285   if (useLastKeep) binarizeFilter->SetUpperThreshold(lastKeep);
286   binarizeFilter->SetInsideValue(FG);
287   binarizeFilter->SetOutsideValue(BG);
288   binarizeFilter->Update();
289   return binarizeFilter->GetOutput();
290 }
291 //--------------------------------------------------------------------
292
293
294 //--------------------------------------------------------------------
295 template<class ImageType>
296 typename ImageType::Pointer
297 clitk::LabelizeAndSelectLabels(typename ImageType::Pointer input,
298                                typename ImageType::PixelType BG, 
299                                typename ImageType::PixelType FG, 
300                                bool isFullyConnected,
301                                int minimalComponentSize,
302                                LabelizeParameters<typename ImageType::PixelType> * param)
303 {
304   typename ImageType::Pointer working_image;
305   working_image = Labelize<ImageType>(input, BG, isFullyConnected, minimalComponentSize);
306   working_image = RemoveLabels<ImageType>(working_image, BG, param->GetLabelsToRemove());
307   working_image = KeepLabels<ImageType>(working_image, 
308                                         BG, FG, 
309                                         param->GetFirstKeep(), 
310                                         param->GetLastKeep(), 
311                                         param->GetUseLastKeep());
312   return working_image;
313 }
314 //--------------------------------------------------------------------
315
316
317 //--------------------------------------------------------------------
318 template<class ImageType>
319 typename ImageType::Pointer
320 clitk::ResizeImageLike(typename ImageType::Pointer input,
321                        typename ImageType::Pointer like, 
322                        typename ImageType::PixelType backgroundValue) 
323 {
324   typedef clitk::CropLikeImageFilter<ImageType> CropFilterType;
325   typename CropFilterType::Pointer cropFilter = CropFilterType::New();
326   cropFilter->SetInput(input);
327   cropFilter->SetCropLikeImage(like);
328   cropFilter->SetBackgroundValue(backgroundValue);
329   cropFilter->Update();
330   return cropFilter->GetOutput();  
331 }
332 //--------------------------------------------------------------------
333
334
335 //--------------------------------------------------------------------
336 template<class MaskImageType>
337 typename MaskImageType::Pointer
338 clitk::SliceBySliceRelativePosition(const MaskImageType * input,
339                                     const MaskImageType * object,
340                                     int direction, 
341                                     double threshold, 
342                                     std::string orientation, 
343                                     bool uniqueConnectedComponent, 
344                                     double spacing, 
345                                     bool inverseflag) 
346 {
347   typedef clitk::SliceBySliceRelativePositionFilter<MaskImageType> SliceRelPosFilterType;
348   typename SliceRelPosFilterType::Pointer sliceRelPosFilter = SliceRelPosFilterType::New();
349   sliceRelPosFilter->VerboseStepFlagOff();
350   sliceRelPosFilter->WriteStepFlagOff();
351   sliceRelPosFilter->SetInput(input);
352   sliceRelPosFilter->SetInputObject(object);
353   sliceRelPosFilter->SetDirection(direction);
354   sliceRelPosFilter->SetFuzzyThreshold(threshold);
355   sliceRelPosFilter->AddOrientationTypeString(orientation);
356   sliceRelPosFilter->SetResampleBeforeRelativePositionFilter((spacing != -1));
357   sliceRelPosFilter->SetIntermediateSpacing(spacing);
358   sliceRelPosFilter->SetUniqueConnectedComponentBySlice(uniqueConnectedComponent);
359   sliceRelPosFilter->SetInverseOrientationFlag(inverseflag);
360   //  sliceRelPosFilter->SetAutoCropFlag(true); ??
361   sliceRelPosFilter->Update();
362   return sliceRelPosFilter->GetOutput();
363 }
364 //--------------------------------------------------------------------
365
366 //--------------------------------------------------------------------
367 template<class ImageType>
368 bool
369 clitk::FindExtremaPointInAGivenDirection(const ImageType * input, 
370                                          typename ImageType::PixelType bg, 
371                                          int direction, bool opposite, 
372                                          typename ImageType::PointType & point)
373 {
374   typename ImageType::PointType dummy;
375   return clitk::FindExtremaPointInAGivenDirection(input, bg, direction, 
376                                                   opposite, dummy, 0, point);
377 }
378 //--------------------------------------------------------------------
379
380 //--------------------------------------------------------------------
381 template<class ImageType>
382 bool
383 clitk::FindExtremaPointInAGivenDirection(const ImageType * input, 
384                                          typename ImageType::PixelType bg, 
385                                          int direction, bool opposite, 
386                                          typename ImageType::PointType refpoint,
387                                          double distanceMax, 
388                                          typename ImageType::PointType & point)
389 {
390   /*
391     loop over input pixels, store the index in the fg that is max
392     according to the given direction. 
393   */
394   typedef itk::ImageRegionConstIteratorWithIndex<ImageType> IteratorType;
395   IteratorType iter(input, input->GetLargestPossibleRegion());
396   iter.GoToBegin();
397   typename ImageType::IndexType max = input->GetLargestPossibleRegion().GetIndex();
398   if (opposite) max = max+input->GetLargestPossibleRegion().GetSize();
399   bool found=false;
400   while (!iter.IsAtEnd()) {
401     if (iter.Get() != bg) {
402       bool test = iter.GetIndex()[direction] >  max[direction];
403       if (opposite) test = !test;
404       if (test) {
405         typename ImageType::PointType p;
406         input->TransformIndexToPhysicalPoint(iter.GetIndex(), p);
407         if ((distanceMax==0) || (p.EuclideanDistanceTo(refpoint) < distanceMax)) {
408           max = iter.GetIndex();
409           found = true;
410         }
411       }
412     }
413     ++iter;
414   }
415   if (!found) return false;
416   input->TransformIndexToPhysicalPoint(max, point);
417   return true;
418 }
419 //--------------------------------------------------------------------
420
421
422 //--------------------------------------------------------------------
423 template<class ImageType>
424 typename ImageType::Pointer
425 clitk::CropImageAbove(typename ImageType::Pointer image, 
426                       int dim, double min, 
427                       bool autoCrop,
428                       typename ImageType::PixelType BG) 
429 {
430   return clitk::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 clitk::CropImageBelow(typename ImageType::Pointer image, 
442                       int dim, double max, 
443                       bool autoCrop,
444                       typename ImageType::PixelType BG) 
445 {
446   typename ImageType::PointType p;
447   image->TransformIndexToPhysicalPoint(image->GetLargestPossibleRegion().GetIndex()+
448                                        image->GetLargestPossibleRegion().GetSize(), p);
449   return clitk::CropImageAlongOneAxis<ImageType>(image, dim, max, p[dim], autoCrop, BG);
450 }
451 //--------------------------------------------------------------------
452
453
454 //--------------------------------------------------------------------
455 template<class ImageType>
456 typename ImageType::Pointer
457 clitk::CropImageAlongOneAxis(typename ImageType::Pointer image, 
458                              int dim, double min, double max, 
459                              bool autoCrop,
460                              typename ImageType::PixelType BG) 
461 {
462   // Compute region size
463   typename ImageType::RegionType region;
464   typename ImageType::SizeType size = image->GetLargestPossibleRegion().GetSize();
465   typename ImageType::PointType p = image->GetOrigin();
466   p[dim] = min;
467   typename ImageType::IndexType start;
468   image->TransformPhysicalPointToIndex(p, start);
469   p[dim] = max;
470   typename ImageType::IndexType end;
471   image->TransformPhysicalPointToIndex(p, end);
472   size[dim] = fabs(end[dim]-start[dim]);
473   region.SetIndex(start);
474   region.SetSize(size);
475   
476   // Perform Crop
477   typedef itk::RegionOfInterestImageFilter<ImageType, ImageType> CropFilterType;
478   typename CropFilterType::Pointer cropFilter = CropFilterType::New();
479   cropFilter->SetInput(image);
480   cropFilter->SetRegionOfInterest(region);
481   cropFilter->Update();
482   typename ImageType::Pointer result = cropFilter->GetOutput();
483   
484   // Auto Crop
485   if (autoCrop) {
486     result = clitk::AutoCrop<ImageType>(result, BG);
487   }
488   return result;
489 }
490 //--------------------------------------------------------------------
491
492
493 //--------------------------------------------------------------------
494 template<class ImageType>
495 void
496 clitk::ComputeCentroids(typename ImageType::Pointer image, 
497                         typename ImageType::PixelType BG, 
498                         std::vector<typename ImageType::PointType> & centroids) 
499 {
500   typedef long LabelType;
501   static const unsigned int Dim = ImageType::ImageDimension;
502   typedef itk::ShapeLabelObject< LabelType, Dim > LabelObjectType;
503   typedef itk::LabelMap< LabelObjectType > LabelMapType;
504   typedef itk::LabelImageToLabelMapFilter<ImageType, LabelMapType> ImageToMapFilterType;
505   typename ImageToMapFilterType::Pointer imageToLabelFilter = ImageToMapFilterType::New(); 
506   typedef itk::ShapeLabelMapFilter<LabelMapType, ImageType> ShapeFilterType; 
507   typename ShapeFilterType::Pointer statFilter = ShapeFilterType::New();
508   imageToLabelFilter->SetBackgroundValue(BG);
509   imageToLabelFilter->SetInput(image);
510   statFilter->SetInput(imageToLabelFilter->GetOutput());
511   statFilter->Update();
512   typename LabelMapType::Pointer labelMap = statFilter->GetOutput();
513
514   centroids.clear();
515   typename ImageType::PointType dummy;
516   centroids.push_back(dummy); // label 0 -> no centroid, use dummy point
517   for(uint i=1; i<labelMap->GetNumberOfLabelObjects()+1; i++) {
518     centroids.push_back(labelMap->GetLabelObject(i)->GetCentroid());
519   } 
520 }
521 //--------------------------------------------------------------------
522
523
524 //--------------------------------------------------------------------
525 template<class ImageType>
526 void
527 clitk::ExtractSlices(typename ImageType::Pointer image, 
528                      int direction, 
529                      std::vector<typename itk::Image<typename ImageType::PixelType, 
530                      ImageType::ImageDimension-1>::Pointer > & slices) 
531 {
532   typedef clitk::ExtractSliceFilter<ImageType> ExtractSliceFilterType;
533   typedef typename ExtractSliceFilterType::SliceType SliceType;
534   typename ExtractSliceFilterType::Pointer 
535     extractSliceFilter = ExtractSliceFilterType::New();
536   extractSliceFilter->SetInput(image);
537   extractSliceFilter->SetDirection(direction);
538   extractSliceFilter->Update();
539   extractSliceFilter->GetOutputSlices(slices);
540 }
541 //--------------------------------------------------------------------
542
543
544 //--------------------------------------------------------------------
545 template<class ImageType>
546 typename ImageType::Pointer
547 clitk::JoinSlices(std::vector<typename itk::Image<typename ImageType::PixelType, 
548                   ImageType::ImageDimension-1>::Pointer > & slices, 
549                   typename ImageType::Pointer input, 
550                   int direction) {
551   typedef typename itk::Image<typename ImageType::PixelType, ImageType::ImageDimension-1> SliceType;
552   typedef itk::JoinSeriesImageFilter<SliceType, ImageType> JoinSeriesFilterType;
553   typename JoinSeriesFilterType::Pointer joinFilter = JoinSeriesFilterType::New();
554   joinFilter->SetOrigin(input->GetOrigin()[direction]);
555   joinFilter->SetSpacing(input->GetSpacing()[direction]);
556   for(unsigned int i=0; i<slices.size(); i++) {
557     joinFilter->PushBackInput(slices[i]);
558   }
559   joinFilter->Update();
560   return joinFilter->GetOutput();
561 }
562 //--------------------------------------------------------------------
563
564
565 //--------------------------------------------------------------------
566 template<class ImageType>
567 void
568 clitk::PointsUtils<ImageType>::Convert2DTo3D(const PointType2D & p, 
569                                              ImagePointer image, 
570                                              const int slice, 
571                                              PointType3D & p3D)  
572 {
573   p3D[0] = p[0]; 
574   p3D[1] = p[1];
575   p3D[2] = (image->GetLargestPossibleRegion().GetIndex()[2]+slice)*image->GetSpacing()[2] 
576     + image->GetOrigin()[2];
577 }
578 //--------------------------------------------------------------------
579
580
581 //--------------------------------------------------------------------
582 template<class ImageType>
583 void 
584 clitk::PointsUtils<ImageType>::Convert2DTo3DList(const MapPoint2DType & map, 
585                                                  ImagePointer image, 
586                                                  VectorPoint3DType & list)
587 {
588   typename MapPoint2DType::const_iterator iter = map.begin();
589   while (iter != map.end()) {
590     PointType3D p;
591     Convert2DTo3D(iter->second, image, iter->first, p);
592     list.push_back(p);
593     ++iter;
594   }
595 }
596 //--------------------------------------------------------------------
597
598 //--------------------------------------------------------------------
599 template<class ImageType>
600 void 
601 clitk::WriteListOfLandmarks(std::vector<typename ImageType::PointType> points, 
602                             std::string filename)
603 {
604   std::ofstream os; 
605   openFileForWriting(os, filename); 
606   os << "LANDMARKS1" << std::endl;  
607   for(uint i=0; i<points.size(); i++) {
608     const typename ImageType::PointType & p = points[i];
609     // Write it in the file
610     os << i << " " << p[0] << " " << p[1] << " " << p[2] << " 0 0 " << std::endl;
611   }
612   os.close();
613 }
614 //--------------------------------------------------------------------
615
616
617 //--------------------------------------------------------------------
618 template<class ImageType>
619 typename ImageType::Pointer 
620 clitk::Dilate(typename ImageType::Pointer image, 
621               double radiusInMM,               
622               typename ImageType::PixelType BG,
623               typename ImageType::PixelType FG,  
624               bool extendSupport)
625 {
626   typename ImageType::SizeType r;
627   for(uint i=0; i<ImageType::ImageDimension; i++) 
628     r[i] = (uint)lrint(radiusInMM/image->GetSpacing()[i]);
629   return clitk::Dilate<ImageType>(image, r, BG, FG, extendSupport);
630 }
631 //--------------------------------------------------------------------
632
633
634 //--------------------------------------------------------------------
635 template<class ImageType>
636 typename ImageType::Pointer 
637 clitk::Dilate(typename ImageType::Pointer image, 
638               typename ImageType::PointType radiusInMM, 
639               typename ImageType::PixelType BG, 
640               typename ImageType::PixelType FG, 
641               bool extendSupport)
642 {
643   typename ImageType::SizeType r;
644   for(uint i=0; i<ImageType::ImageDimension; i++) 
645     r[i] = (uint)lrint(radiusInMM[i]/image->GetSpacing()[i]);
646   return clitk::Dilate<ImageType>(image, r, BG, FG, extendSupport);
647 }
648 //--------------------------------------------------------------------
649
650
651 //--------------------------------------------------------------------
652 template<class ImageType>
653 typename ImageType::Pointer 
654 clitk::Dilate(typename ImageType::Pointer image, 
655               typename ImageType::SizeType radius, 
656               typename ImageType::PixelType BG, 
657               typename ImageType::PixelType FG, 
658               bool extendSupport)
659 {
660   // Create kernel for dilatation
661   typedef itk::BinaryBallStructuringElement<typename ImageType::PixelType, 
662                                             ImageType::ImageDimension> KernelType;
663   KernelType structuringElement;
664   structuringElement.SetRadius(radius);
665   structuringElement.CreateStructuringElement();
666
667   if (extendSupport) {
668     typedef itk::ConstantPadImageFilter<ImageType, ImageType> PadFilterType;
669     typename PadFilterType::Pointer padFilter = PadFilterType::New();
670     padFilter->SetInput(image);
671     typename ImageType::SizeType lower;
672     typename ImageType::SizeType upper;
673     for(uint i=0; i<3; i++) {
674       lower[i] = upper[i] = 2*(radius[i]+1);
675     }
676     padFilter->SetPadLowerBound(lower);
677     padFilter->SetPadUpperBound(upper);
678     padFilter->Update();
679     image = padFilter->GetOutput();
680   }
681
682   // Dilate  filter
683   typedef itk::BinaryDilateImageFilter<ImageType, ImageType , KernelType> DilateFilterType;
684   typename DilateFilterType::Pointer dilateFilter = DilateFilterType::New();
685   dilateFilter->SetBackgroundValue(BG);
686   dilateFilter->SetForegroundValue(FG);
687   dilateFilter->SetBoundaryToForeground(false);
688   dilateFilter->SetKernel(structuringElement);
689   dilateFilter->SetInput(image);
690   dilateFilter->Update();
691   return image = dilateFilter->GetOutput();
692 }
693 //--------------------------------------------------------------------
694
695
696 //--------------------------------------------------------------------
697 template<class ValueType, class VectorType>
698 void clitk::ConvertOption(std::string optionName, uint given, 
699                           ValueType * values, VectorType & p, 
700                           uint dim, bool required) 
701 {
702   if (required && (given == 0)) {
703     clitkExceptionMacro("The option --" << optionName << " must be set and have 1 or " 
704                         << dim << " values.");
705   }
706   if (given == 1) {
707     for(uint i=0; i<dim; i++) p[i] = values[0];
708     return;
709   }
710   if (given == dim) {
711     for(uint i=0; i<dim; i++) p[i] = values[i];
712     return;
713   }
714   if (given == 0) return;
715   clitkExceptionMacro("The option --" << optionName << " must have 1 or " 
716                       << dim << " values.");
717 }
718 //--------------------------------------------------------------------
719
720
721 //--------------------------------------------------------------------
722 /*
723   http://www.gamedev.net/community/forums/topic.asp?topic_id=542870
724   Assuming the points are (Ax,Ay) (Bx,By) and (Cx,Cy), you need to compute:
725   (Bx - Ax) * (Cy - Ay) - (By - Ay) * (Cx - Ax)
726   This will equal zero if the point C is on the line formed by
727   points A and B, and will have a different sign depending on the
728   side. Which side this is depends on the orientation of your (x,y)
729   coordinates, but you can plug test values for A,B and C into this
730   formula to determine whether negative values are to the left or to
731   the right.
732   => to accelerate, start with formula, when change sign -> stop and fill
733 */
734 template<class ImageType>
735 void 
736 clitk::SliceBySliceSetBackgroundFromLineSeparation(typename ImageType::Pointer input, 
737                                                    std::vector<typename ImageType::PointType> & lA, 
738                                                    std::vector<typename ImageType::PointType> & lB, 
739                                                    typename ImageType::PixelType BG, 
740                                                    int mainDirection, 
741                                                    double offsetToKeep)
742 {
743   
744   typedef itk::ImageSliceIteratorWithIndex<ImageType> SliceIteratorType;
745   SliceIteratorType siter = SliceIteratorType(input, 
746                                               input->GetLargestPossibleRegion());
747   siter.SetFirstDirection(0);
748   siter.SetSecondDirection(1);
749   siter.GoToBegin();
750   int i=0;
751   typename ImageType::PointType A;
752   typename ImageType::PointType B;
753   typename ImageType::PointType C;
754   while (!siter.IsAtEnd()) {
755     // Check that the current slice correspond to the current point
756     input->TransformIndexToPhysicalPoint(siter.GetIndex(), C);
757     if (C[2] != lA[i][2]) {
758       // DD(C);
759       // DD(lA[i]);
760     }
761     else {
762       // Define A,B,C points
763       A = lA[i];
764       B = lB[i];
765       C = A;
766       C[mainDirection] += offsetToKeep; // I know I must keep this point
767       double s = (B[0] - A[0]) * (C[1] - A[1]) - (B[1] - A[1]) * (C[0] - A[0]);
768       bool isPositive = s<0;
769       while (!siter.IsAtEndOfSlice()) {
770         while (!siter.IsAtEndOfLine()) {
771           // Very slow, I know ... but image should be very small
772           input->TransformIndexToPhysicalPoint(siter.GetIndex(), C);
773           double s = (B[0] - A[0]) * (C[1] - A[1]) - (B[1] - A[1]) * (C[0] - A[0]);
774           if (s == 0) siter.Set(BG); // on the line, we decide to remove
775           if (isPositive) {
776             if (s > 0) siter.Set(BG);
777           }
778           else {
779             if (s < 0) siter.Set(BG); 
780           }
781           ++siter;
782         }
783         siter.NextLine();
784       }
785       ++i;
786     }
787     siter.NextSlice();
788   }
789 }                                                   
790 //--------------------------------------------------------------------