]> Creatis software - clitk.git/blob - itk/clitkSegmentationUtils.txx
add findExtrema, jointSlices
[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
35 //--------------------------------------------------------------------
36 template<class ImageType>
37 void clitk::ComputeBBFromImageRegion(typename ImageType::Pointer image, 
38                                      typename ImageType::RegionType region,
39                                      typename itk::BoundingBox<unsigned long, 
40                                                                ImageType::ImageDimension>::Pointer bb) {
41   typedef typename ImageType::IndexType IndexType;
42   IndexType firstIndex;
43   IndexType lastIndex;
44   for(unsigned int i=0; i<image->GetImageDimension(); i++) {
45     firstIndex[i] = region.GetIndex()[i];
46     lastIndex[i] = firstIndex[i]+region.GetSize()[i];
47   }
48
49   typedef itk::BoundingBox<unsigned long, 
50                            ImageType::ImageDimension> BBType;
51   typedef typename BBType::PointType PointType;
52   PointType lastPoint;
53   PointType firstPoint;
54   image->TransformIndexToPhysicalPoint(firstIndex, firstPoint);
55   image->TransformIndexToPhysicalPoint(lastIndex, lastPoint);
56
57   bb->SetMaximum(lastPoint);
58   bb->SetMinimum(firstPoint);
59 }
60 //--------------------------------------------------------------------
61
62
63 //--------------------------------------------------------------------
64 template<int Dimension>
65 void clitk::ComputeBBIntersection(typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbo, 
66                                   typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbi1, 
67                                   typename itk::BoundingBox<unsigned long, Dimension>::Pointer bbi2) {
68
69   typedef itk::BoundingBox<unsigned long, Dimension> BBType;
70   typedef typename BBType::PointType PointType;
71   PointType lastPoint;
72   PointType firstPoint;
73
74   for(unsigned int i=0; i<Dimension; i++) {
75     firstPoint[i] = std::max(bbi1->GetMinimum()[i], 
76                              bbi2->GetMinimum()[i]);
77     lastPoint[i] = std::min(bbi1->GetMaximum()[i], 
78                             bbi2->GetMaximum()[i]);
79   }
80
81   bbo->SetMaximum(lastPoint);
82   bbo->SetMinimum(firstPoint);
83 }
84 //--------------------------------------------------------------------
85
86
87 //--------------------------------------------------------------------
88 template<class ImageType>
89 void clitk::ComputeRegionFromBB(typename ImageType::Pointer image, 
90                                 const typename itk::BoundingBox<unsigned long, 
91                                 ImageType::ImageDimension>::Pointer bb, 
92                                 typename ImageType::RegionType & region) {
93   // Types
94   typedef typename ImageType::IndexType  IndexType;
95   typedef typename ImageType::PointType  PointType;
96   typedef typename ImageType::RegionType RegionType;
97   typedef typename ImageType::SizeType   SizeType;
98
99   // Region starting point
100   IndexType regionStart;
101   PointType start = bb->GetMinimum();
102   image->TransformPhysicalPointToIndex(start, regionStart);
103     
104   // Region size
105   SizeType regionSize;
106   PointType maxs = bb->GetMaximum();
107   PointType mins = bb->GetMinimum();
108   for(unsigned int i=0; i<ImageType::ImageDimension; i++) {
109     // DD(maxs[i]);
110     // DD(mins[i]);
111     // DD((maxs[i] - mins[i])/image->GetSpacing()[i]);
112     regionSize[i] = lrint((maxs[i] - mins[i])/image->GetSpacing()[i]);
113     // DD(regionSize[i]);
114   }
115    
116   // Create region
117   region.SetIndex(regionStart);
118   region.SetSize(regionSize);
119 }
120 //--------------------------------------------------------------------
121
122 //--------------------------------------------------------------------
123 template<class ImageType, class TMaskImageType>
124 typename ImageType::Pointer
125 clitk::SetBackground(const ImageType * input, 
126                      const TMaskImageType * mask, 
127                      typename TMaskImageType::PixelType maskBG,
128                      typename ImageType::PixelType outValue, 
129                      bool inPlace) {
130   typedef clitk::SetBackgroundImageFilter<ImageType, TMaskImageType, ImageType> 
131     SetBackgroundImageFilterType;
132   typename SetBackgroundImageFilterType::Pointer setBackgroundFilter 
133     = SetBackgroundImageFilterType::New();
134   //  if (inPlace) setBackgroundFilter->ReleaseDataFlagOn(); // No seg fault
135   setBackgroundFilter->SetInPlace(inPlace); // This is important to keep memory low
136   setBackgroundFilter->SetInput(input);
137   setBackgroundFilter->SetInput2(mask);
138   setBackgroundFilter->SetMaskValue(maskBG);
139   setBackgroundFilter->SetOutsideValue(outValue);
140   setBackgroundFilter->Update();
141   return setBackgroundFilter->GetOutput();
142 }
143 //--------------------------------------------------------------------
144
145
146 //--------------------------------------------------------------------
147 template<class ImageType>
148 int clitk::GetNumberOfConnectedComponentLabels(typename ImageType::Pointer input, 
149                                                typename ImageType::PixelType BG, 
150                                                bool isFullyConnected) {
151   // Connected Component label 
152   typedef itk::ConnectedComponentImageFilter<ImageType, ImageType> ConnectFilterType;
153   typename ConnectFilterType::Pointer connectFilter = ConnectFilterType::New();
154   connectFilter->SetInput(input);
155   connectFilter->SetBackgroundValue(BG);
156   connectFilter->SetFullyConnected(isFullyConnected);
157   connectFilter->Update();
158   
159   // Return result
160   return connectFilter->GetObjectCount();
161 }
162 //--------------------------------------------------------------------
163
164 //--------------------------------------------------------------------
165 /*
166   Warning : in this cas, we consider outputType like inputType, not
167   InternalImageType. Be sure it fits.
168  */
169 template<class ImageType>
170 typename ImageType::Pointer
171 clitk::Labelize(const ImageType * input, 
172                 typename ImageType::PixelType BG, 
173                 bool isFullyConnected, 
174                 int minimalComponentSize) {
175   // InternalImageType for storing large number of component
176   typedef itk::Image<int, ImageType::ImageDimension> InternalImageType;
177   
178   // Connected Component label 
179   typedef itk::ConnectedComponentImageFilter<ImageType, InternalImageType> ConnectFilterType;
180   typename ConnectFilterType::Pointer connectFilter = ConnectFilterType::New();
181   //  connectFilter->ReleaseDataFlagOn(); 
182   connectFilter->SetInput(input);
183   connectFilter->SetBackgroundValue(BG);
184   connectFilter->SetFullyConnected(isFullyConnected);
185   
186   // Sort by size and remove too small area.
187   typedef itk::RelabelComponentImageFilter<InternalImageType, ImageType> RelabelFilterType;
188   typename RelabelFilterType::Pointer relabelFilter = RelabelFilterType::New();
189   //  relabelFilter->ReleaseDataFlagOn(); // if yes, fail when ExplosionControlledThresholdConnectedImageFilter ???
190   relabelFilter->SetInput(connectFilter->GetOutput());
191   relabelFilter->SetMinimumObjectSize(minimalComponentSize);
192   relabelFilter->Update();
193
194   // Return result
195   typename ImageType::Pointer output = relabelFilter->GetOutput();
196   return output;
197 }
198 //--------------------------------------------------------------------
199
200
201 //--------------------------------------------------------------------
202 template<class ImageType>
203 typename ImageType::Pointer
204 clitk::RemoveLabels(typename ImageType::Pointer input, 
205                     typename ImageType::PixelType BG,
206                     std::vector<typename ImageType::PixelType> & labelsToRemove) {
207   typename ImageType::Pointer working_image = input;
208   for (unsigned int i=0; i <labelsToRemove.size(); i++) {
209     typedef clitk::SetBackgroundImageFilter<ImageType, ImageType> SetBackgroundImageFilterType;
210     typename SetBackgroundImageFilterType::Pointer setBackgroundFilter = SetBackgroundImageFilterType::New();
211     setBackgroundFilter->SetInput(input);
212     setBackgroundFilter->SetInput2(input);
213     setBackgroundFilter->SetMaskValue(labelsToRemove[i]);
214     setBackgroundFilter->SetOutsideValue(BG);
215     setBackgroundFilter->Update();
216     working_image = setBackgroundFilter->GetOutput();
217   }
218   return working_image;
219 }
220 //--------------------------------------------------------------------
221
222
223 //--------------------------------------------------------------------
224 template<class ImageType>
225 typename ImageType::Pointer
226 clitk::KeepLabels(const ImageType * input, 
227                   typename ImageType::PixelType BG, 
228                   typename ImageType::PixelType FG, 
229                   typename ImageType::PixelType firstKeep, 
230                   typename ImageType::PixelType lastKeep, 
231                   bool useLastKeep) {
232   typedef itk::BinaryThresholdImageFilter<ImageType, ImageType> BinarizeFilterType; 
233   typename BinarizeFilterType::Pointer binarizeFilter = BinarizeFilterType::New();
234   binarizeFilter->SetInput(input);
235   binarizeFilter->SetLowerThreshold(firstKeep);
236   if (useLastKeep) binarizeFilter->SetUpperThreshold(lastKeep);
237   binarizeFilter->SetInsideValue(FG);
238   binarizeFilter->SetOutsideValue(BG);
239   binarizeFilter->Update();
240   return binarizeFilter->GetOutput();
241 }
242 //--------------------------------------------------------------------
243
244
245 //--------------------------------------------------------------------
246 template<class ImageType>
247 typename ImageType::Pointer
248 clitk::LabelizeAndSelectLabels(typename ImageType::Pointer input,
249                                typename ImageType::PixelType BG, 
250                                typename ImageType::PixelType FG, 
251                                bool isFullyConnected,
252                                int minimalComponentSize,
253                                LabelizeParameters<typename ImageType::PixelType> * param)
254 {
255   typename ImageType::Pointer working_image;
256   working_image = Labelize<ImageType>(input, BG, isFullyConnected, minimalComponentSize);
257   working_image = RemoveLabels<ImageType>(working_image, BG, param->GetLabelsToRemove());
258   working_image = KeepLabels<ImageType>(working_image, 
259                                         BG, FG, 
260                                         param->GetFirstKeep(), 
261                                         param->GetLastKeep(), 
262                                         param->GetUseLastKeep());
263   return working_image;
264 }
265 //--------------------------------------------------------------------
266
267
268 //--------------------------------------------------------------------
269 template<class ImageType>
270 typename ImageType::Pointer
271 clitk::ResizeImageLike(typename ImageType::Pointer input,
272                        typename ImageType::Pointer like, 
273                        typename ImageType::PixelType backgroundValue) 
274 {
275   typedef clitk::CropLikeImageFilter<ImageType> CropFilterType;
276   typename CropFilterType::Pointer cropFilter = CropFilterType::New();
277   cropFilter->SetInput(input);
278   cropFilter->SetCropLikeImage(like);
279   cropFilter->SetBackgroundValue(backgroundValue);
280   cropFilter->Update();
281   return cropFilter->GetOutput();  
282 }
283 //--------------------------------------------------------------------
284
285
286 //--------------------------------------------------------------------
287 template<class MaskImageType>
288 typename MaskImageType::Pointer
289 clitk::SliceBySliceRelativePosition(const MaskImageType * input,
290                                     const MaskImageType * object,
291                                     int direction, 
292                                     double threshold, 
293                                     std::string orientation, 
294                                     bool uniqueConnectedComponent, 
295                                     double spacing, 
296                                     bool inverseflag) 
297 {
298   typedef clitk::SliceBySliceRelativePositionFilter<MaskImageType> SliceRelPosFilterType;
299   typename SliceRelPosFilterType::Pointer sliceRelPosFilter = SliceRelPosFilterType::New();
300   sliceRelPosFilter->VerboseStepFlagOff();
301   sliceRelPosFilter->WriteStepFlagOff();
302   sliceRelPosFilter->SetInput(input);
303   sliceRelPosFilter->SetInputObject(object);
304   sliceRelPosFilter->SetDirection(direction);
305   sliceRelPosFilter->SetFuzzyThreshold(threshold);
306   sliceRelPosFilter->AddOrientationTypeString(orientation);
307   sliceRelPosFilter->SetResampleBeforeRelativePositionFilter((spacing != -1));
308   sliceRelPosFilter->SetIntermediateSpacing(spacing);
309   sliceRelPosFilter->SetUniqueConnectedComponentBySlice(uniqueConnectedComponent);
310   sliceRelPosFilter->SetInverseOrientationFlag(inverseflag);
311   //  sliceRelPosFilter->SetAutoCropFlag(true); ??
312   sliceRelPosFilter->Update();
313   return sliceRelPosFilter->GetOutput();
314 }
315 //--------------------------------------------------------------------
316
317 //--------------------------------------------------------------------
318 template<class ImageType>
319 bool
320 clitk::FindExtremaPointInAGivenDirection(const ImageType * input, 
321                                          typename ImageType::PixelType bg, 
322                                          int direction, bool opposite, 
323                                          typename ImageType::PointType & point)
324 {
325   typename ImageType::PointType dummy;
326   return clitk::FindExtremaPointInAGivenDirection(input, bg, direction, 
327                                                   opposite, dummy, 0, point);
328 }
329 //--------------------------------------------------------------------
330
331 //--------------------------------------------------------------------
332 template<class ImageType>
333 bool
334 clitk::FindExtremaPointInAGivenDirection(const ImageType * input, 
335                                          typename ImageType::PixelType bg, 
336                                          int direction, bool opposite, 
337                                          typename ImageType::PointType refpoint,
338                                          double distanceMax, 
339                                          typename ImageType::PointType & point)
340 {
341   /*
342     loop over input pixels, store the index in the fg that is max
343     according to the given direction. 
344   */
345   typedef itk::ImageRegionConstIteratorWithIndex<ImageType> IteratorType;
346   IteratorType iter(input, input->GetLargestPossibleRegion());
347   iter.GoToBegin();
348   typename ImageType::IndexType max = input->GetLargestPossibleRegion().GetIndex();
349   if (opposite) max = max+input->GetLargestPossibleRegion().GetSize();
350   bool found=false;
351   while (!iter.IsAtEnd()) {
352     if (iter.Get() != bg) {
353       bool test = iter.GetIndex()[direction] >  max[direction];
354       if (opposite) test = !test;
355       if (test) {
356         typename ImageType::PointType p;
357         input->TransformIndexToPhysicalPoint(iter.GetIndex(), p);
358         if ((distanceMax==0) || (p.EuclideanDistanceTo(refpoint) < distanceMax)) {
359           max = iter.GetIndex();
360           found = true;
361         }
362       }
363     }
364     ++iter;
365   }
366   if (!found) return false;
367   input->TransformIndexToPhysicalPoint(max, point);
368   return true;
369 }
370 //--------------------------------------------------------------------
371
372
373 //--------------------------------------------------------------------
374 template<class ImageType>
375 typename ImageType::Pointer
376 clitk::CropImageAbove(typename ImageType::Pointer image, 
377                       int dim, double min, 
378                       bool autoCrop,
379                       typename ImageType::PixelType BG) 
380 {
381   return clitk::CropImageAlongOneAxis<ImageType>(image, dim, 
382                                                  image->GetOrigin()[dim], 
383                                                  min,
384                                                  autoCrop, BG);
385 }
386 //--------------------------------------------------------------------
387
388
389 //--------------------------------------------------------------------
390 template<class ImageType>
391 typename ImageType::Pointer
392 clitk::CropImageBelow(typename ImageType::Pointer image, 
393                       int dim, double max, 
394                       bool autoCrop,
395                       typename ImageType::PixelType BG) 
396 {
397   typename ImageType::PointType p;
398   image->TransformIndexToPhysicalPoint(image->GetLargestPossibleRegion().GetIndex()+
399                                        image->GetLargestPossibleRegion().GetSize(), p);
400   return clitk::CropImageAlongOneAxis<ImageType>(image, dim, max, p[dim], autoCrop, BG);
401 }
402 //--------------------------------------------------------------------
403
404
405 //--------------------------------------------------------------------
406 template<class ImageType>
407 typename ImageType::Pointer
408 clitk::CropImageAlongOneAxis(typename ImageType::Pointer image, 
409                              int dim, double min, double max, 
410                              bool autoCrop,
411                              typename ImageType::PixelType BG) 
412 {
413   // Compute region size
414   typename ImageType::RegionType region;
415   typename ImageType::SizeType size = image->GetLargestPossibleRegion().GetSize();
416   typename ImageType::PointType p = image->GetOrigin();
417   p[dim] = min;
418   typename ImageType::IndexType start;
419   image->TransformPhysicalPointToIndex(p, start);
420   p[dim] = max;
421   typename ImageType::IndexType end;
422   image->TransformPhysicalPointToIndex(p, end);
423   size[dim] = fabs(end[dim]-start[dim]);
424   region.SetIndex(start);
425   region.SetSize(size);
426   
427   // Perform Crop
428   typedef itk::RegionOfInterestImageFilter<ImageType, ImageType> CropFilterType;
429   typename CropFilterType::Pointer cropFilter = CropFilterType::New();
430   cropFilter->SetInput(image);
431   cropFilter->SetRegionOfInterest(region);
432   cropFilter->Update();
433   typename ImageType::Pointer result = cropFilter->GetOutput();
434   
435   // Auto Crop
436   if (autoCrop) {
437     result = clitk::AutoCrop<ImageType>(result, BG);
438   }
439   return result;
440 }
441 //--------------------------------------------------------------------
442
443
444 //--------------------------------------------------------------------
445 template<class ImageType>
446 void
447 clitk::ComputeCentroids(typename ImageType::Pointer image, 
448                         typename ImageType::PixelType BG, 
449                         std::vector<typename ImageType::PointType> & centroids) 
450 {
451   typedef long LabelType;
452   static const unsigned int Dim = ImageType::ImageDimension;
453   typedef itk::ShapeLabelObject< LabelType, Dim > LabelObjectType;
454   typedef itk::LabelMap< LabelObjectType > LabelMapType;
455   typedef itk::LabelImageToLabelMapFilter<ImageType, LabelMapType> ImageToMapFilterType;
456   typename ImageToMapFilterType::Pointer imageToLabelFilter = ImageToMapFilterType::New(); 
457   typedef itk::ShapeLabelMapFilter<LabelMapType, ImageType> ShapeFilterType; 
458   typename ShapeFilterType::Pointer statFilter = ShapeFilterType::New();
459   imageToLabelFilter->SetBackgroundValue(BG);
460   imageToLabelFilter->SetInput(image);
461   statFilter->SetInput(imageToLabelFilter->GetOutput());
462   statFilter->Update();
463   typename LabelMapType::Pointer labelMap = statFilter->GetOutput();
464
465   centroids.clear();
466   typename ImageType::PointType dummy;
467   centroids.push_back(dummy); // label 0 -> no centroid, use dummy point
468   for(uint i=1; i<labelMap->GetNumberOfLabelObjects()+1; i++) {
469     centroids.push_back(labelMap->GetLabelObject(i)->GetCentroid());
470   } 
471 }
472 //--------------------------------------------------------------------
473
474
475 //--------------------------------------------------------------------
476 template<class ImageType>
477 void
478 clitk::ExtractSlices(typename ImageType::Pointer image, 
479                      int direction, 
480                      std::vector<typename itk::Image<typename ImageType::PixelType, 
481                      ImageType::ImageDimension-1>::Pointer > & slices) 
482 {
483   typedef clitk::ExtractSliceFilter<ImageType> ExtractSliceFilterType;
484   typedef typename ExtractSliceFilterType::SliceType SliceType;
485   typename ExtractSliceFilterType::Pointer 
486     extractSliceFilter = ExtractSliceFilterType::New();
487   extractSliceFilter->SetInput(image);
488   extractSliceFilter->SetDirection(direction);
489   extractSliceFilter->Update();
490   extractSliceFilter->GetOutputSlices(slices);
491 }
492 //--------------------------------------------------------------------
493
494
495 //--------------------------------------------------------------------
496 template<class ImageType>
497 typename ImageType::Pointer
498 clitk::JoinSlices(std::vector<typename itk::Image<typename ImageType::PixelType, 
499                   ImageType::ImageDimension-1>::Pointer > & slices, 
500                   typename ImageType::Pointer input, 
501                   int direction) {
502   typedef typename itk::Image<typename ImageType::PixelType, ImageType::ImageDimension-1> SliceType;
503   typedef itk::JoinSeriesImageFilter<SliceType, ImageType> JoinSeriesFilterType;
504   typename JoinSeriesFilterType::Pointer joinFilter = JoinSeriesFilterType::New();
505   joinFilter->SetOrigin(input->GetOrigin()[direction]);
506   joinFilter->SetSpacing(input->GetSpacing()[direction]);
507   for(unsigned int i=0; i<slices.size(); i++) {
508     joinFilter->PushBackInput(slices[i]);
509   }
510   joinFilter->Update();
511   return joinFilter->GetOutput();
512 }
513 //--------------------------------------------------------------------
514
515
516 //--------------------------------------------------------------------
517 template<class ImageType>
518 void
519 clitk::PointsUtils<ImageType>::Convert2DTo3D(const PointType2D & p, 
520                                              ImagePointer image, 
521                                              const int slice, 
522                                              PointType3D & p3D)  
523 {
524   p3D[0] = p[0]; 
525   p3D[1] = p[1];
526   p3D[2] = (image->GetLargestPossibleRegion().GetIndex()[2]+slice)*image->GetSpacing()[2] 
527     + image->GetOrigin()[2];
528 }
529 //--------------------------------------------------------------------
530
531
532 //--------------------------------------------------------------------
533 template<class ImageType>
534 void 
535 clitk::PointsUtils<ImageType>::Convert2DTo3DList(const MapPoint2DType & map, 
536                                                  ImagePointer image, 
537                                                  VectorPoint3DType & list)
538 {
539   typename MapPoint2DType::const_iterator iter = map.begin();
540   while (iter != map.end()) {
541     PointType3D p;
542     Convert2DTo3D(iter->second, image, iter->first, p);
543     list.push_back(p);
544     ++iter;
545   }
546 }
547 //--------------------------------------------------------------------
548
549 //--------------------------------------------------------------------
550 template<class ImageType>
551 void 
552 clitk::WriteListOfLandmarks(std::vector<typename ImageType::PointType> points, 
553                             std::string filename)
554 {
555   std::ofstream os; 
556   openFileForWriting(os, filename); 
557   os << "LANDMARKS1" << std::endl;  
558   for(uint i=0; i<points.size(); i++) {
559     const typename ImageType::PointType & p = points[i];
560     // Write it in the file
561     os << i << " " << p[0] << " " << p[1] << " " << p[2] << " 0 0 " << std::endl;
562   }
563   os.close();
564 }
565 //--------------------------------------------------------------------
566
567
568 //--------------------------------------------------------------------
569 template<class ImageType>
570 typename ImageType::Pointer 
571 clitk::Dilate(typename ImageType::Pointer image, 
572               double radiusInMM,               
573               typename ImageType::PixelType BG,
574               typename ImageType::PixelType FG,  
575               bool extendSupport)
576 {
577   typename ImageType::SizeType r;
578   for(uint i=0; i<ImageType::ImageDimension; i++) 
579     r[i] = (uint)lrint(radiusInMM/image->GetSpacing()[i]);
580   return clitk::Dilate<ImageType>(image, r, BG, FG, extendSupport);
581 }
582 //--------------------------------------------------------------------
583
584
585 //--------------------------------------------------------------------
586 template<class ImageType>
587 typename ImageType::Pointer 
588 clitk::Dilate(typename ImageType::Pointer image, 
589               typename ImageType::PointType radiusInMM, 
590               typename ImageType::PixelType BG, 
591               typename ImageType::PixelType FG, 
592               bool extendSupport)
593 {
594   typename ImageType::SizeType r;
595   for(uint i=0; i<ImageType::ImageDimension; i++) 
596     r[i] = (uint)lrint(radiusInMM[i]/image->GetSpacing()[i]);
597   return clitk::Dilate<ImageType>(image, r, BG, FG, extendSupport);
598 }
599 //--------------------------------------------------------------------
600
601
602 //--------------------------------------------------------------------
603 template<class ImageType>
604 typename ImageType::Pointer 
605 clitk::Dilate(typename ImageType::Pointer image, 
606               typename ImageType::SizeType radius, 
607               typename ImageType::PixelType BG, 
608               typename ImageType::PixelType FG, 
609               bool extendSupport)
610 {
611   // Create kernel for dilatation
612   typedef itk::BinaryBallStructuringElement<typename ImageType::PixelType, 
613                                             ImageType::ImageDimension> KernelType;
614   KernelType structuringElement;
615   structuringElement.SetRadius(radius);
616   structuringElement.CreateStructuringElement();
617
618   if (extendSupport) {
619     typedef itk::ConstantPadImageFilter<ImageType, ImageType> PadFilterType;
620     typename PadFilterType::Pointer padFilter = PadFilterType::New();
621     padFilter->SetInput(image);
622     typename ImageType::SizeType lower;
623     typename ImageType::SizeType upper;
624     for(uint i=0; i<3; i++) {
625       lower[i] = upper[i] = 2*(radius[i]+1);
626     }
627     padFilter->SetPadLowerBound(lower);
628     padFilter->SetPadUpperBound(upper);
629     padFilter->Update();
630     image = padFilter->GetOutput();
631   }
632
633   // Dilate  filter
634   typedef itk::BinaryDilateImageFilter<ImageType, ImageType , KernelType> DilateFilterType;
635   typename DilateFilterType::Pointer dilateFilter = DilateFilterType::New();
636   dilateFilter->SetBackgroundValue(BG);
637   dilateFilter->SetForegroundValue(FG);
638   dilateFilter->SetBoundaryToForeground(false);
639   dilateFilter->SetKernel(structuringElement);
640   dilateFilter->SetInput(image);
641   dilateFilter->Update();
642   return image = dilateFilter->GetOutput();
643 }
644 //--------------------------------------------------------------------
645
646
647 //--------------------------------------------------------------------
648 template<class ValueType, class VectorType>
649 void clitk::ConvertOption(std::string optionName, uint given, 
650                           ValueType * values, VectorType & p, 
651                           uint dim, bool required) 
652 {
653   if (required && (given == 0)) {
654     clitkExceptionMacro("The option --" << optionName << " must be set and have 1 or " 
655                         << dim << " values.");
656   }
657   if (given == 1) {
658     for(uint i=0; i<dim; i++) p[i] = values[0];
659     return;
660   }
661   if (given == dim) {
662     for(uint i=0; i<dim; i++) p[i] = values[i];
663     return;
664   }
665   if (given == 0) return;
666   clitkExceptionMacro("The option --" << optionName << " must have 1 or " 
667                       << dim << " values.");
668 }
669 //--------------------------------------------------------------------
670
671