ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkFEMRegistrationFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 
19 #ifndef itkFEMRegistrationFilter_h
20 #define itkFEMRegistrationFilter_h
21 
24 #include "itkFEMObject.h"
27 #include "itkFEMImageMetricLoad.h"
29 
30 #include "itkImage.h"
31 #include "itkVector.h"
36 #include "itkWarpImageFilter.h"
37 #include "itkImageToImageMetric.h"
40 #include "itkFixedArray.h"
41 
43 #include "itkFEMLoadLandmark.h"
44 
45 #include "vnl/vnl_vector.h"
46 #include "itkMath.h"
47 #include "vnl/vnl_vector_fixed.h"
48 
49 #include <iostream>
50 #include <string>
51 
52 namespace itk
53 {
54 namespace fem
55 {
56 
117 template <typename TMovingImage, typename TFixedImage, typename TFemObjectType>
118 class ITK_TEMPLATE_EXPORT FEMRegistrationFilter : public ImageToImageFilter<TMovingImage, TFixedImage>
119 {
120 public:
125 
127  itkNewMacro(Self);
128 
131 
132  typedef TMovingImage MovingImageType;
133  typedef TFixedImage FixedImageType;
134  typedef TFemObjectType FEMObjectType;
135  typedef typename FixedImageType::PixelType PixelType;
138 
140  itkStaticConstMacro(ImageDimension, unsigned int,
141  FixedImageType::ImageDimension);
142 
146 
147  enum Sign { positive = 1, negative = -1 };
148  typedef double Float;
150 
151  typedef std::vector<typename LoadLandmark::Pointer> LandmarkArrayType;
155 
160 
162 
164  typedef double CoordRepType;
167 
169 
173 
175 
178 
180 
186 
189 
191  void SetMovingImage(MovingImageType* R);
192 
198  {
199  return m_MovingImage;
200  }
201 
204  {
205  return m_OriginalMovingImage;
206  }
207 
209  void SetFixedImage(FixedImageType* T);
210 
212  {
213  return m_FixedImage;
214  }
215 
221  void SetInputFEMObject(FEMObjectType* F, unsigned int level = 0);
222 
223  FEMObjectType * GetInputFEMObject(unsigned int level = 0);
224 
226  void RunRegistration();
227 
229  void IterativeSolve(SolverType *S);
230 
232  void MultiResSolve();
233 
235  void WarpImage(const MovingImageType * R);
236 
240  {
241  return m_WarpedImage;
242  }
243 
245  void ComputeJacobian();
246 
249  {
250  return m_FloatImage;
251  }
252 
255  {
256  return m_Field;
257  }
258 
261  {
262  m_FieldSize = F->GetLargestPossibleRegion().GetSize();
263  m_Field = F;
264  }
266 
268  void AddLandmark(PointType source, PointType target);
269 
270  void InsertLandmark(unsigned int i, PointType source, PointType target);
271 
272  void DeleteLandmark(unsigned int i);
273 
274  void ClearLandmarks();
275 
276  void GetLandmark(unsigned int i, PointType& source, PointType& target);
277 
285  void EnforceDiffeomorphism(float thresh, SolverType *S, bool onlywriteimages);
286 
293  void SetMeshPixelsPerElementAtEachResolution(unsigned int i, unsigned int which = 0)
294  {
295  m_MeshPixelsPerElementAtEachResolution[which] = i;
296  }
298 
302  void SetNumberOfIntegrationPoints(unsigned int i, unsigned int which = 0)
303  {
304  m_NumberOfIntegrationPoints[which] = i;
305  }
306 
312  void SetWidthOfMetricRegion(unsigned int i, unsigned int which = 0)
313  {
314  m_MetricWidth[which] = i;
315  }
316 
317  unsigned int GetWidthOfMetricRegion(unsigned int which = 0)
318  {
319  return m_MetricWidth[which];
320  }
321 
326  void SetMaximumIterations(unsigned int i, unsigned int which)
327  {
328  m_Maxiters[which] = i;
329  }
330 
336  itkSetMacro(TimeStep, Float);
337  itkGetMacro(TimeStep, Float);
339 
344  itkSetMacro(Alpha, Float);
345  itkGetMacro(Alpha, Float);
347 
351  itkSetMacro(UseLandmarks, bool);
352  itkGetMacro(UseLandmarks, bool);
353  itkBooleanMacro(UseLandmarks);
354 #if !defined(ITK_LEGACY_REMOVE)
355 
357  itkLegacyMacro(void SetUseLandmarksOff())
358  {
359  SetUseLandmarks(false);
360  }
361 
363  itkLegacyMacro(void SetUseLandmarksOn())
364  {
365  SetUseLandmarks(true);
366  }
367 #endif
368 
372  itkSetMacro(UseMassMatrix, bool);
373  itkGetMacro(UseMassMatrix, bool);
374  itkBooleanMacro(UseMassMatrix);
376 
380  itkSetMacro(EnergyReductionFactor, Float);
381  itkGetMacro(EnergyReductionFactor, Float);
383 
385  void SetElasticity(Float i, unsigned int which = 0)
386  {
387  m_E[which] = i;
388  }
389 
391  Float GetElasticity(unsigned int which = 0)
392  {
393  return m_E[which];
394  }
395 
397  void SetRho(Float r, unsigned int which = 0)
398  {
399  m_Rho[which] = r;
400  }
401 
403  void SetGamma(Float r, unsigned int which = 0)
404  {
405  m_Gamma[which] = r;
406  }
407 
410  {
411  m_DescentDirection = positive;
412  }
413 
416  {
417  m_DescentDirection = negative;
418  }
419 
424  itkSetMacro(DoLineSearchOnImageEnergy, unsigned int);
425  itkGetMacro(DoLineSearchOnImageEnergy, unsigned int);
427 
428 
433  itkSetMacro(UseNormalizedGradient, bool);
434  itkGetMacro(UseNormalizedGradient, bool);
435  itkBooleanMacro(UseNormalizedGradient);
436 #if !defined(ITK_LEGACY_REMOVE)
437 
439  itkLegacyMacro(void SetUseNormalizedGradientOff())
440  {
441  SetUseNormalizedGradient(false);
442  }
443 
445  itkLegacyMacro(void SetUseNormalizedGradientOn())
446  {
447  SetUseNormalizedGradient(true);
448  }
449 #endif
450 
454  itkSetMacro(EmployRegridding, unsigned int);
455  itkGetMacro(EmployRegridding, unsigned int);
457 
462  itkSetMacro(LineSearchMaximumIterations, unsigned int);
463  itkGetMacro(LineSearchMaximumIterations, unsigned int);
465 
470  {
471  return m_FullImageSize;
472  }
473 
478  itkGetModifiableObjectMacro(Metric, MetricBaseType);
479  itkSetObjectMacro(Metric, MetricBaseType);
481 
490  void ChooseMetric( unsigned int whichmetric );
491 
495  unsigned int GetMetricType()
496  {
497  return m_WhichMetric;
498  }
499 
502  {
503  m_Element = e;
504  }
505 
508  {
509  m_Material = m;
510  }
511 
513  void PrintVectorField(unsigned int modnum = 1000);
514 
518  void SetMaxLevel(unsigned int level);
519  itkGetMacro(MaxLevel, unsigned int);
521 
526  itkSetMacro(CreateMeshFromImage, bool);
527  itkGetMacro(CreateMeshFromImage, bool);
528  itkBooleanMacro(CreateMeshFromImage);
529 #if !defined(ITK_LEGACY_REMOVE)
530 
532  itkLegacyMacro(void SetCreateMeshFromImageOn())
533  {
534  SetCreateMeshFromImage(true);
535  }
536 
538  itkLegacyMacro(void SetCreateMeshFromImageOff())
539  {
540  SetCreateMeshFromImage(false);
541  }
542 #endif
543 
545  itkSetObjectMacro( Interpolator, InterpolatorType );
546 
548  itkGetModifiableObjectMacro( Interpolator, InterpolatorType );
549 
553  itkSetMacro(StandardDeviations, StandardDeviationsType);
554  virtual void SetStandardDeviations(double value);
556 
559  itkGetConstReferenceMacro(StandardDeviations, StandardDeviationsType);
560 
563  itkSetMacro(MaximumKernelWidth, unsigned int);
564  itkGetConstMacro(MaximumKernelWidth, unsigned int);
566 
569  itkSetMacro(MaximumError, double);
570  itkGetConstMacro(MaximumError, double);
572 
573 protected:
574  FEMRegistrationFilter();
575  ~FEMRegistrationFilter() ITK_OVERRIDE;
576 
577  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
578 
580  void CreateMesh(unsigned int ElementsPerSide, SolverType *solver);
581 
583  void ApplyLoads(ImageSizeType Isz, double* spacing = ITK_NULLPTR);
584 
586  void ApplyImageLoads(MovingImageType* i1, FixedImageType* i2);
587 
588  // FIXME - Not implemented
591  void CreateLinearSystemSolver();
592 
593  // FIXME - Not implemented
595  Float EvaluateEnergy();
596 
601  void InterpolateVectorField(SolverType *S);
602 
603  // FIXME - Not implemented
606  FloatImageType * GetMetricImage(FieldType* F);
607 
609  FieldPointer ExpandVectorField(ExpandFactorsType* expandFactors, FieldType* f);
610 
612  void SampleVectorFieldAtNodes(SolverType *S);
613 
615  Float EvaluateResidual(SolverType *mySolver, Float t);
616 
617  // FIXME - Replace with BSD Code
618  /* Finds a triplet that brackets the energy minimum. From Numerical Recipes.*/
619  // void FindBracketingTriplet(SolverType& mySolver,Float* a, Float* b, Float* c);
620  void FindBracketingTriplet(SolverType *mySolver, Float* a, Float* b, Float* c);
621 
626  Float GoldenSection(SolverType *mySolver, Float tol = 0.01, unsigned int MaxIters = 25);
627 
629  itkSetObjectMacro( Load, ImageMetricLoadType );
630  itkGetModifiableObjectMacro(Load, ImageMetricLoadType );
632 
634  void SmoothDisplacementField();
635 
636 private:
637 
638  void InitializeField();
639 
640  ITK_DISALLOW_COPY_AND_ASSIGN(FEMRegistrationFilter);
641 
642  unsigned int m_DoLineSearchOnImageEnergy;
643  unsigned int m_LineSearchMaximumIterations;
644 
646  vnl_vector<unsigned int> m_NumberOfIntegrationPoints; // resolution of integration
647  vnl_vector<unsigned int> m_MetricWidth; // number of iterations at each level
648  vnl_vector<unsigned int> m_Maxiters;
649  unsigned int m_TotalIterations; // total number of iterations that were run
650  unsigned int m_MaxLevel;
651  unsigned int m_FileCount; // keeps track of number of files written
652  unsigned int m_CurrentLevel; // current resolution level
653 
654  typename FixedImageType::SizeType m_CurrentLevelImageSize;
655 
656  unsigned int m_WhichMetric;
657 
660  vnl_vector<unsigned int> m_MeshPixelsPerElementAtEachResolution;
661 
662  Float m_TimeStep;
663  vnl_vector<Float> m_E;
664  vnl_vector<Float> m_Rho;
665  vnl_vector<Float> m_Gamma;
666  Float m_Energy; // current value of energy
667  Float m_MinE; // minimum recorded energy
668  Float m_MinJacobian; // minimum recorded energy
669  Float m_Alpha;
670 
671  bool m_UseLandmarks;
672  bool m_UseMassMatrix;
673  bool m_UseNormalizedGradient;
674  bool m_CreateMeshFromImage;
675  unsigned int m_EmployRegridding;
676  Sign m_DescentDirection;
677  Float m_EnergyReductionFactor;
678  ImageSizeType m_FullImageSize;
679  ImageSizeType m_ImageOrigin;
680 
683  ImageSizeType m_ImageScaling;
684  ImageSizeType m_CurrentImageScaling;
685  typename FieldType::RegionType m_FieldRegion;
686  typename FieldType::SizeType m_FieldSize;
687  typename FieldType::Pointer m_Field;
688 
689  // Only use TotalField if re-gridding is employed.
690  typename FieldType::Pointer m_TotalField;
691  typename ImageMetricLoadType::Pointer m_Load; // Defines the load to use
692 
693  // Define the warper
694  typename WarperType::Pointer m_Warper;
695 
696  // Declare a new image to hold the warped reference
697  typename FixedImageType::Pointer m_WarpedImage;
698  typename FloatImageType::Pointer m_FloatImage;
699  typename FixedImageType::RegionType m_Wregion;
700  typename FixedImageType::IndexType m_Windex;
701 
702  // Declare images for target and reference
703  typename MovingImageType::Pointer m_MovingImage;
704  typename MovingImageType::Pointer m_OriginalMovingImage;
705  typename FixedImageType::Pointer m_FixedImage;
706 
707  // Element and metric pointers
708  typename Element::Pointer m_Element;
709  typename MaterialType::Pointer m_Material;
711  typename FEMObjectType::Pointer m_FEMObject;
712 
713  LandmarkArrayType m_LandmarkArray;
714  InterpolatorPointer m_Interpolator;
715 
716  double m_MaximumError;
717 
718  unsigned int m_MaximumKernelWidth;
719 
720  StandardDeviationsType m_StandardDeviations;
721 
722 };
723 
724 }
725 } // end namespace itk::fem
726 
727 #ifndef ITK_MANUAL_INSTANTIATION
728 #include "itkFEMRegistrationFilter.hxx"
729 #endif
730 
731 #endif
void SetElasticity(Float i, unsigned int which=0)
InterpolationGridType::SizeType InterpolationGridSizeType
Image< float, itkGetStaticConstMacro(ImageDimension)> FloatImageType
Light weight base class for most itk classes.
SmartPointer< const Self > ConstPointer
itk::Image< Element::ConstPointer, ImageDimension > InterpolationGridType
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
void SetMeshPixelsPerElementAtEachResolution(unsigned int i, unsigned int which=0)
itk::Image< VectorType, itkGetStaticConstMacro(ImageDimension)> FieldType
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:259
itk::VectorIndexSelectionCastImageFilter< FieldType, FloatImageType > IndexSelectCasterType
FEM Solver for time dependent problems; uses Crank-Nicolson implicit discretization scheme...
InterpolationGridType::PointType InterpolationGridPointType
itk::ImageRegionIteratorWithIndex< FloatImageType > FloatImageIterator
void SetNumberOfIntegrationPoints(unsigned int i, unsigned int which=0)
Generate a rectilinar mesh from an image. The result is stored in a FEMObject.
void SetWidthOfMetricRegion(unsigned int i, unsigned int which=0)
FixedArray< double, ImageDimension > StandardDeviationsType
Float GetElasticity(unsigned int which=0)
itk::fem::ImageToRectilinearFEMObjectFilter< TMovingImage > ImageToMeshType
itk::Vector< float, itkGetStaticConstMacro(ImageDimension)> VectorType
itk::VectorExpandImageFilter< FieldType, FieldType > ExpanderType
void SetRho(Float r, unsigned int which=0)
Expand the size of a vector image by an integer factor in each dimension.
Linear elasticity material class.
ExpanderType::ExpandFactorsType ExpandFactorsType
General image pair load that uses the itkFiniteDifferenceFunctions.
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
void SetMaterial(MaterialType::Pointer m)
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
InterpolatorType::Pointer InterpolatorPointer
ImageToImageFilter< TMovingImage, TFixedImage > Superclass
FiniteDifferenceFunctionLoad< MovingImageType, FixedImageType > ImageMetricLoadType
const SizeType & GetSize() const
Warps an image using an input displacement field.
FEM Image registration filter. The image registration problem is modelled here with the finite elemen...
VectorInterpolateImageFunction< FieldType, CoordRepType > InterpolatorType
itk::ImageRegionIteratorWithIndex< FixedImageType > ImageIterator
SolverCrankNicolson< ImageDimension > SolverType
unsigned int GetWidthOfMetricRegion(unsigned int which=0)
Array for FEMP objects.
Definition: itkFEMPArray.h:40
Abstract base element class.
itk::WarpImageFilter< MovingImageType, FixedImageType, FieldType > WarperType
VectorLinearInterpolateImageFunction< FieldType, CoordRepType > DefaultInterpolatorType
Base class for filters that take an image as input and produce an image as output.
itk::ImageRegionIteratorWithIndex< FieldType > FieldIterator
LinearSystemWrapperItpack LinearSystemSolverType
Base class for all vector image interpolaters.
Extracts the selected index of the vector that is the input pixel type.
LinearSystemWrapper class that uses Itpack numeric library functions to define and solve a sparse lin...
MetricBaseType::Pointer MetricBaseTypePointer
static ITK_CONSTEXPR_VAR double e
The base of the natural logarithm or Euler&#39;s number
Definition: itkMath.h:56
std::vector< typename LoadLandmark::Pointer > LandmarkArrayType
void SetGamma(Float r, unsigned int which=0)
void SetMaximumIterations(unsigned int i, unsigned int which)
Templated n-dimensional image class.
Definition: itkImage.h:75
PDEDeformableRegistrationFunction< FixedImageType, MovingImageType, FieldType > MetricBaseType
Linearly interpolate a vector image at specified positions.