Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Classes | Public Types | Public Member Functions | Protected Member Functions
itk::PSMEntropyModelFilter< TImage, TShapeMatrix > Class Template Reference

This the most basic of all PSM model optimization filters. This filter assembles all of the necessary components to perform an entropy-based correspondence optimization for constructing a point-based shape model of shape surfaces. More...

#include <itkPSMEntropyModelFilter.h>

+ Inheritance diagram for itk::PSMEntropyModelFilter< TImage, TShapeMatrix >:
+ Collaboration diagram for itk::PSMEntropyModelFilter< TImage, TShapeMatrix >:

Classes

struct  CuttingPlaneType
 

Public Types

typedef PSMEntropyModelFilter Self
 
typedef InPlaceImageFilter< TImage, TImage > Superclass
 
typedef SmartPointer< SelfPointer
 
typedef SmartPointer< const SelfConstPointer
 
typedef PSMParticleSystem< Dimension > ParticleSystemType
 
typedef TShapeMatrix ShapeMatrixType
 
typedef ShapeMatrixType::Pointer ShapeMatrixPointerType
 
typedef TImage ImageType
 
typedef ImageType::PointType PointType
 
typedef PSMGradientDescentOptimizer< typename ImageType::PixelType, Dimension > OptimizerType
 

Public Member Functions

 itkStaticConstMacro (Dimension, unsigned int, TImage::ImageDimension)
 
 itkNewMacro (Self)
 
 itkTypeMacro (PSMEntropyModelFilter, InPlaceImageFilter)
 
void SetInput (const std::string &s, itk::DataObject *o)
 
void SetInput (const TImage *image)
 
void SetInput (unsigned int index, const TImage *image)
 
void SetNumberOfScales (unsigned int n)
 
 itkGetMacro (NumberOfScales, unsigned int)
 
 itkSetMacro (CurrentScale, unsigned int)
 
 itkGetMacro (CurrentScale, unsigned int)
 
void SetInputCorrespondencePoints (unsigned int index, const std::vector< PointType > &corr)
 
 itkGetObjectMacro (ParticleSystem, ParticleSystemType)
 
 itkGetConstObjectMacro (ParticleSystem, ParticleSystemType)
 
PSMParticleEntropyFunction< typename ImageType::PixelType, Dimension > * GetParticleEntropyFunction ()
 
PSMShapeEntropyFunction< Dimension > * ShapeEntropyFunction ()
 
 itkGetObjectMacro (Optimizer, OptimizerType)
 
 itkGetConstObjectMacro (Optimizer, OptimizerType)
 
void SetMaximumNumberOfIterations (const std::vector< unsigned int > &n)
 
void SetMaximumNumberOfIterations (unsigned int n, unsigned int i=0)
 
unsigned int GetMaximumNumberOfIterations (unsigned int i=0) const
 
unsigned int GetNumberOfElapsedIterations () const
 
 itkGetMacro (Initialized, bool)
 
 itkGetMacro (Initializing, bool)
 
 itkSetMacro (Initializing, bool)
 
void SetShapeEntropyWeighting (double w)
 
double GetShapeEntropyWeighting () const
 
const std::vector< double > & GetShapePCAVariances () const
 
void SetRegularizationInitial (const std::vector< double > &v)
 
void SetRegularizationInitial (double t, unsigned int i=0)
 
double GetRegularizationInitial (unsigned int i=0) const
 
void SetRegularizationFinal (const std::vector< double > &v)
 
void SetRegularizationFinal (double t, unsigned int i=0)
 
double GetRegularizationFinal (unsigned int i=0) const
 
void SetRegularizationDecaySpan (const std::vector< double > &v)
 
void SetRegularizationDecaySpan (double t, unsigned int i=0)
 
double GetRegularizationDecaySpan (unsigned int i=0) const
 
void SetRegularizationConstant (double r)
 
void SetMinimumVariance (double r)
 
double GetRegularizationConstant () const
 
double GetMinimumVariance () const
 
void SetTolerance (const std::vector< double > &v)
 
void SetTolerance (double t, unsigned int i=0)
 
double GetTolerance (unsigned int i=0) const
 
 itkGetMacro (ShapeMatrix, ShapeMatrixPointerType)
 
const ShapeMatrixTypeGetShapeMatrix () const
 
bool CreateSingleCorrespondence ()
 
void SetDomainName (const std::string &s, unsigned int i)
 
unsigned int GetDomainIndexByName (const std::string &name)
 
void AddCuttingPlane (const vnl_vector_fixed< double, 3 > &x, const vnl_vector_fixed< double, 3 > &y, const vnl_vector_fixed< double, 3 > &z, unsigned int domain)
 
void AddCuttingPlane (const vnl_vector_fixed< double, 3 > &x, const vnl_vector_fixed< double, 3 > &y, const vnl_vector_fixed< double, 3 > &z, const std::string &s)
 

Protected Member Functions

void PrintSelf (std::ostream &os, Indent indent) const
 
void GenerateData ()
 
virtual void AllocateDataCaches ()
 
virtual void AllocateWorkingImages ()
 
virtual void AllocateDomainsAndNeighborhoods ()
 
void InitializeCorrespondences ()
 
 itkSetMacro (Initialized, bool)
 
void OptimizerIterateCallback (Object *, const EventObject &)
 
PSMImplicitSurfaceDomain< typename ImageType::PixelType, Dimension >::Pointer & GetDomainByName (const std::string &name)
 
PSMImplicitSurfaceDomain< typename ImageType::PixelType, Dimension >::Pointer & GetDomain (unsigned int i)
 

Detailed Description

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
class itk::PSMEntropyModelFilter< TImage, TShapeMatrix >

This the most basic of all PSM model optimization filters. This filter assembles all of the necessary components to perform an entropy-based correspondence optimization for constructing a point-based shape model of shape surfaces.

This class combines a PSMShapeEntropyFunction with a PSMParticleEntropyFunction within a PSMTwoCostFunction to produce an optimization of particle positions with respect to both their positional entropy on shape surfaces and their entropy within a higher-dimensional shape space (the entropy of the distributions of shapes). The relative weighting of these two terms may be adjusted using the SetShapeEntropyWeighting method, which increases or decreases the influence of the entropy of the shape model on the optimization. The default value for ShapeEntropyWeighting is 1.0.

NEED CITATION HERE

HOW DO YOU USE THIS FILTER? First you need to specify the input images. Input images are in the form of distance transforms that contain an implicit surface representation at the zero level set within the image. The zero level set is the surface on which the particles are constrained during the optimization. You will specify multiple input images, one for each shape in your shape cohort. To specify the inputs, use the method SetInput(i,image *), where i is the number of the input and image * is an itk::Image::Pointer to the distance transform image.

This filter can be run in "initialization mode", which causes all data structures to allocate and initialize, but does not run the actual optimization. To run in initialization mode, set the Initializing flag to true.

Author
Josh Cates

Definition at line 78 of file itkPSMEntropyModelFilter.h.

Member Typedef Documentation

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
typedef TImage itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::ImageType

Type of the input/output image.

Definition at line 105 of file itkPSMEntropyModelFilter.h.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
typedef PSMGradientDescentOptimizer<typename ImageType::PixelType, Dimension> itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::OptimizerType

Type of the optimizer

Definition at line 111 of file itkPSMEntropyModelFilter.h.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
typedef PSMParticleSystem<Dimension> itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::ParticleSystemType

Type of the particle system

Definition at line 92 of file itkPSMEntropyModelFilter.h.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
typedef ImageType::PointType itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::PointType

Expose the point type

Definition at line 108 of file itkPSMEntropyModelFilter.h.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
typedef PSMEntropyModelFilter itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::Self

Standard class typedefs.

Definition at line 83 of file itkPSMEntropyModelFilter.h.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
typedef TShapeMatrix itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::ShapeMatrixType

The type of the Shape Matrix, which defines optimization behavior.

Definition at line 95 of file itkPSMEntropyModelFilter.h.

Member Function Documentation

template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::AddCuttingPlane ( const vnl_vector_fixed< double, 3 > &  x,
const vnl_vector_fixed< double, 3 > &  y,
const vnl_vector_fixed< double, 3 > &  z,
unsigned int  domain 
)

Adds cutting plane information for given domain i. Currently, only one cutting plane per domain is supported.

Definition at line 136 of file itkPSMEntropyModelFilter.hxx.

140 {
141  if (m_CuttingPlanes.size() < domain+1)
142  {
143  m_CuttingPlanes.resize(domain+1);
144  }
145  CuttingPlaneType p;
146  p.x = x;
147  p.y = y;
148  p.z = z;
149  p.valid = true;
150  m_CuttingPlanes[domain] = p;
151 }
template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::AllocateDataCaches ( )
protectedvirtual

This method sets up temporary data structures needed for the optimization. It is called during the first run of GenerateData.

Definition at line 80 of file itkPSMEntropyModelFilter.hxx.

81 {
82  // Set up the various data caches that the optimization functions will use.
83  // Sigma cache is caching a parameter for the local particle entropy
84  // computation that is updated for each particle.
85  m_SigmaCache = PSMContainerArrayAttribute<double, Dimension>::New();
86  m_ParticleSystem->RegisterAttribute(m_SigmaCache);
87  m_ParticleEntropyFunction->SetSpatialSigmaCache(m_SigmaCache);
88 
89  // m_MeanCurvatureCache = PSMMeanCurvatureAttribute<typename ImageType::PixelType, Dimension>::New();
90  // m_ParticleSystem->RegisterAttribute(m_MeanCurvatureCache);
91 }
template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::AllocateDomainsAndNeighborhoods ( )
protectedvirtual

This method allocates data structures necessary to keep track of image domains and image neighborhoods. It is called during the first run of GenerateData.

Definition at line 107 of file itkPSMEntropyModelFilter.hxx.

108 {
109  // Allocate all the necessary domains and neighborhoods. This must be done
110  // *after* registering the attributes to the particle system since some of
111  // them respond to AddDomain.
112  for (unsigned int i = 0; i < this->GetNumberOfInputs(); i++)
113  {
114  m_DomainList.push_back( PSMImplicitSurfaceDomain<typename
115  ImageType::PixelType, Dimension>::New() );
116  m_NeighborhoodList.push_back( PSMSurfaceNeighborhood<ImageType>::New() );
117  m_DomainList[i]->SetSigma(m_WorkingImages[i]->GetSpacing()[0] * 2.0);
118  m_DomainList[i]->SetImage(m_WorkingImages[i]);
119 
120  if (m_CuttingPlanes.size() > i)
121  {
122  m_DomainList[i]->SetCuttingPlane(m_CuttingPlanes[i].x,
123  m_CuttingPlanes[i].y,
124  m_CuttingPlanes[i].z);
125  }
126 
127  m_ParticleSystem->AddDomain(m_DomainList[i]);
128  m_ParticleSystem->SetNeighborhood(i, m_NeighborhoodList[i]);
129  }
130 }
template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::AllocateWorkingImages ( )
protectedvirtual

This method allocates working images needed during the optimization. It is called during the first run of GenerateData.

Definition at line 95 of file itkPSMEntropyModelFilter.hxx.

96 {
97  m_WorkingImages.resize(this->GetNumberOfInputs());
98  for (unsigned int i = 0; i < this->GetNumberOfInputs(); i++)
99  {
100  m_WorkingImages[i] = const_cast<TImage *>(this->GetInput(i));
101  }
102 }
template<class TImage , class TShapeMatrix >
bool itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::CreateSingleCorrespondence ( )

Creates a single correspondence point to each domain that is placed at the surface position closest to the upper-left-hand corner of each of the input images.

Definition at line 347 of file itkPSMEntropyModelFilter.hxx.

348 {
349  bool ok = true;
350  for (unsigned int i = 0; i < m_ParticleSystem->GetNumberOfDomains(); i++)
351  {
352  typename itk::ZeroCrossingImageFilter<ImageType, ImageType>::Pointer zc =
353  itk::ZeroCrossingImageFilter<ImageType, ImageType>::New();
354 
355  zc->SetInput( dynamic_cast<PSMImplicitSurfaceDomain<typename ImageType::PixelType, Dimension> *>(
356  m_ParticleSystem->GetDomain(i))->GetImage());
357  zc->Update();
358  itk::ImageRegionConstIteratorWithIndex<ImageType> it(zc->GetOutput(),
359  zc->GetOutput()->GetRequestedRegion());
360 
361  bool done = false;
362  for (it.GoToReverseBegin(); !it.IsAtReverseEnd() && done == false; --it)
363  {
364  if (it.Get() == 1.0)
365  {
366  PointType pos;
367  dynamic_cast<PSMImplicitSurfaceDomain<typename ImageType::PixelType, Dimension> *>(
368  m_ParticleSystem->GetDomain(i))->GetImage()->TransformIndexToPhysicalPoint(it.GetIndex(), pos);
369  try
370  {
371  m_ParticleSystem->AddPosition(pos, i);
372  done = true;
373  }
374  catch(itk::ExceptionObject &)
375  {
376  done = false;
377  ok = false;
378  }
379  }
380  }
381  }
382  return ok;
383 }
template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GenerateData ( )
protected

This method does the work of the optimization.

Definition at line 220 of file itkPSMEntropyModelFilter.hxx.

221 {
222  // First, undo the default InPlaceImageFilter functionality that will
223  // release our inputs.
224  this->SetInPlace(false);
225 
226  // If this is the first call to GenerateData, then we need to allocate some structures
227  if (this->GetInitialized() == false)
228  {
229  this->AllocateWorkingImages();
230  this->AllocateDataCaches();
231  this->GetOptimizer()->SetCostFunction(m_CostFunction);
232 
234 
235  // Point the optimizer to the particle system.
236  this->GetOptimizer()->SetParticleSystem(this->GetParticleSystem());
237  // this->ReadTransforms();
239  this->InitializeOptimizationFunctions();
240 
241  this->SetInitialized(true);
242  }
243 
244  // This filter can be run in "initialization mode", which causes all
245  // data structures to allocate and initialize, but does not run the
246  // actual optimization.
247  if (this->GetInitializing() == true) return;
248 
249  // Multiscale optimization. This can be stopped and restarted by
250  // calling GenerateData again by maintaining m_CurrentScale.
251  for (unsigned int scale = m_CurrentScale; scale < m_NumberOfScales; scale++)
252  {
253  m_CurrentScale = scale;
254 
255  // Set up the optimization parameters for this scale, unless
256  // this is a restart from a previous call to GenerateData. This
257  // section deals with convergence criteria.
258  if (this->GetAbortGenerateData() == false)
259  {
260  // This section deals with convergence.
261  m_Optimizer->SetTolerance(m_Tolerances[scale]);
262  m_Optimizer->SetMaximumNumberOfIterations(m_MaximumNumberOfIterations[scale]);
263 
264  // Set up the exponentially-decreasing regularization constant. If
265  // the decay span is greater than 1 iteration, then we will set up
266  // the annealing approach. Otherwise, the optimizer will simply use
267  // its constant annealing parameter.
268  if (m_RegularizationDecaySpan[scale] >= 1.0f)
269  {
270  m_ShapeEntropyFunction->SetMinimumVarianceDecay(m_RegularizationInitial[scale],
271  m_RegularizationFinal[scale],
272  m_RegularizationDecaySpan[scale]);
273  }
274  else
275  {
276  m_ShapeEntropyFunction->SetMinimumVariance(m_RegularizationInitial[scale]);
277  m_ShapeEntropyFunction->SetHoldMinimumVariance(true);
278  }
279  }
280 
281  this->SetAbortGenerateData(false); // reset the abort flag
282 
283  // Finally, run the optimization. Abort flag is checked in the
284  // iteration callback registered with the optimizer. See
285  // this->OptimizerIterateCallback
286  this->GetOptimizer()->StartOptimization();
287 
288 
289  // If this is multiscale, split particles for the next iteration
290  // -- but not if this is the last iteration.
291  if ((m_NumberOfScales > 1) && (scale != m_NumberOfScales-1) )
292  {
293  m_ParticleSystem->SplitAllParticles(this->GetInput()->GetSpacing()[0] * 1.0);
294  }
295 
296  }
297 }
void SetMinimumVarianceDecay(double initial_value, double final_value, double time_period)
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
PSMImplicitSurfaceDomain< typename ImageType::PixelType, Dimension>::Pointer& itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GetDomain ( unsigned int  i)
inlineprotected

Returns a pointer to the ith domain in the domain list.

Definition at line 613 of file itkPSMEntropyModelFilter.h.

614  {
615  if (i >= m_DomainList.size())
616  { itkExceptionMacro("Requested domain not available"); }
617  return m_DomainList[i];
618  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
PSMImplicitSurfaceDomain< typename ImageType::PixelType, Dimension>::Pointer& itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GetDomainByName ( const std::string &  name)
inlineprotected

Returns a pointer to the domain that is referenced by the given name

Definition at line 605 of file itkPSMEntropyModelFilter.h.

606  {
607  return this->GetDomain(this->GetDomainIndexByName(name));
608  }
unsigned int GetDomainIndexByName(const std::string &name)
PSMImplicitSurfaceDomain< typename ImageType::PixelType, Dimension >::Pointer & GetDomain(unsigned int i)
template<class TImage , class TShapeMatrix >
unsigned int itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GetDomainIndexByName ( const std::string &  name)

Returns the index associated with the given name. Throws an exception if the given name is not found in the domain name list.

Definition at line 168 of file itkPSMEntropyModelFilter.hxx.

169 {
170  // This is slow, but is not anticipated to be used very extensively.
171  // Better to have simple code than overengineer something not
172  // performance-critical.
173 
174  for (unsigned int i = 0; i < m_DomainListNames.size(); i++)
175  {
176  if (m_DomainListNames[i] == name)
177  {
178  return i;
179  }
180  }
181 
182  itkExceptionMacro("The domain name " + name + " was not found. ");
183  return 0; // algorithm will never reach this point
184 }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
unsigned int itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GetNumberOfElapsedIterations ( ) const
inline

Return the number of elapsed iterations of the solver. Note that this number is reset to zero with each call to GenerateData.

Definition at line 274 of file itkPSMEntropyModelFilter.h.

275  {
276  return m_Optimizer->GetNumberOfIterations();
277  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
PSMParticleEntropyFunction<typename ImageType::PixelType, Dimension>* itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GetParticleEntropyFunction ( )
inline

Returns a pointer to the ParticleEntropyFunction. This is the cost function that will maximize positional entropy of particle distributions on shape surfaces. This is the cost function term that ensures good shape representations.

Definition at line 220 of file itkPSMEntropyModelFilter.h.

221  {
222  return m_ParticleEntropyFunction;
223  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
const std::vector<double>& itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::GetShapePCAVariances ( ) const
inline

Returns a vector of the variances associated with the PCA decomposition of the shape space. There will be only number_of_samples - 1 total valid PCA modes. This method simply calls the method of the same name in PSMShapeEntropyFunction. It is included here for convenience.

Definition at line 342 of file itkPSMEntropyModelFilter.h.

343  {
344  return m_ShapeEntropyFunction->GetShapePCAVariances();
345  }
const std::vector< double > & GetShapePCAVariances() const
template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::InitializeCorrespondences ( )
protected

This method copies the correspondence point lists provided by SetInputCorrespondencePoints into the initialized particle system. It is called during the first run of GenerateData.

Definition at line 188 of file itkPSMEntropyModelFilter.hxx.

189 {
190  // If the user has specified correspondence points on the input,
191  // then make sure the number of lists of correspondence points
192  // matches the number of inputs (distance transforms). Everything
193  // should match this->GetNumberOfInputs().
194  if (m_InputCorrespondencePoints.size() > 0)
195  {
196  if (this->GetNumberOfInputs() != m_InputCorrespondencePoints.size() )
197  {
198  itkExceptionMacro("The number of inputs does not match the number of correspondence point lists.");
199  }
200  else
201  {
202  for (unsigned int i = 0; i < m_InputCorrespondencePoints.size(); i++)
203  {
204  this->GetParticleSystem()->AddPositionList(m_InputCorrespondencePoints[i],i);
205  }
206  }
207  }
208  else // No input correspondences are specified, so add a point to each surface.
209  {
211  }
212 
213  // Push position information out to all observers (necessary to
214  // correctly fill out the shape matrix)
215  this->GetParticleSystem()->SynchronizePositions();
216 }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkGetMacro ( Initialized  ,
bool   
)

The Initialized flag indicates whether GenerateData has previously been called and all the necessary temporary data structures have been allocated.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkGetMacro ( Initializing  ,
bool   
)

Get/Set the Initializing flag. When set to true, this flag will prevent GenerateData from starting the optimization. When m_Initializing is true, GenerateData will allocated all data structures and set necessary parameters, but return before running the optimization.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkGetMacro ( ShapeMatrix  ,
ShapeMatrixPointerType   
)

Returns shape matrix. Each column of the shape matrix is a different point-based representation of a shape in the cohort. There are several types of shape matrices, which may be defined differently for subclasses of the PSMEntropyModelFilter class.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkGetObjectMacro ( ParticleSystem  ,
ParticleSystemType   
)

Returns the particle system used in the shape model computation.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkGetObjectMacro ( Optimizer  ,
OptimizerType   
)

Return a pointer to the optimizer object.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkNewMacro ( Self  )

Method for creation through the object factory.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkSetMacro ( Initialized  ,
bool   
)
protected

The Initialized flag indicates whether GenerateData has previously been called and all the necessary temporary data structures have been allocated.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkStaticConstMacro ( Dimension  ,
unsigned  int,
TImage::ImageDimension   
)

Dimensionality of the domain of the particle system.

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::itkTypeMacro ( PSMEntropyModelFilter< TImage, TShapeMatrix >  ,
InPlaceImageFilter   
)

Run-time type information (and related methods).

template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::OptimizerIterateCallback ( Object *  ,
const EventObject &   
)
protected

A callback method that will be called after each iteration of the Solver.

Definition at line 328 of file itkPSMEntropyModelFilter.hxx.

329 {
330  // Update any optimization values based on potential changes from
331  // the user during the last iteration here. This allows for changes in a
332  // thread-safe manner.
333 
334  // Now invoke any iteration events registered with this filter.
335  this->InvokeEvent(itk::IterationEvent());
336 
337  // Request to terminate this filter.
338  if (this->GetAbortGenerateData() == true)
339  {
340  m_Optimizer->StopOptimization();
341  }
342 }
template<class TImage , class TShapeMatrix >
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetDomainName ( const std::string &  s,
unsigned int  i 
)

Sets the name of the ith domain. This functionality is only for convenience, but required if you want to lookup a domain index by its name.

Definition at line 156 of file itkPSMEntropyModelFilter.hxx.

157 {
158  if (m_DomainListNames.size() < (i+1))
159  {
160  m_DomainListNames.resize(i+1);
161  }
162  m_DomainListNames[i] = s;
163 }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetInput ( const std::string &  s,
itk::DataObject *  o 
)
inline

Override parent classes to expand input list on new inputs. This filter requires multiple input images. Each input to this filter represents a shape sample. Note that there may be hundreds of inputs to this filter!

Definition at line 132 of file itkPSMEntropyModelFilter.h.

133  {
134  Superclass::SetInput(s,o);
135  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetInputCorrespondencePoints ( unsigned int  index,
const std::vector< PointType > &  corr 
)
inline

This is a placeholder for a more formalized implementation of a point-based shape model. Set the Nth list of correspondences for the input model. These are the initial correspondence point positions for the Nth input surface (see SetInput).

Definition at line 202 of file itkPSMEntropyModelFilter.h.

203  {
204  if (m_InputCorrespondencePoints.size() < (index+1))
205  {
206  m_InputCorrespondencePoints.resize(index+1);
207  }
208  m_InputCorrespondencePoints[index] = corr;
209  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetMaximumNumberOfIterations ( const std::vector< unsigned int > &  n)
inline

Set / Get the maximum number of iterations that the solver is allowed to execute. If set to 0, this parameter is ignored and the solver will iterate until any other convergence criteria are met. A maximum number of iterations must be specified for each scale (second, optional parameter if NumberOfScales > 1).

Definition at line 243 of file itkPSMEntropyModelFilter.h.

244  {
245  m_MaximumNumberOfIterations = n;
246  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetNumberOfScales ( unsigned int  n)
inline

The number of scales for the optimization. Starting with the input correspondence points for Scale 0, the algorithm will run a full optimization with the parameters given for that scale (RegularizationInitial, RegularizationFinal, etc.). After optimization at each scale, the particle system size will be doubled by splitting each particle and optimization will proceed at the next scale. If no correspondence points are supplied on the inputs, then the system at Scale 0 will be initialized with a single correspondence point, placed at the surface position closest to the upper-left-hand corner of each of the input images.

If NumberOfScales = 1, then the filter will simply optimize the positions of the input correspondence points. Note that you must supply parameters separately for each scale of the optimization.

Definition at line 166 of file itkPSMEntropyModelFilter.h.

167  {
168  m_NumberOfScales = n;
169  if (m_Tolerances.size() < n)
170  {
171  m_Tolerances.resize(n);
172  }
173  if (m_MaximumNumberOfIterations.size() < n)
174  {
175  m_MaximumNumberOfIterations.resize(n);
176  }
177  if (m_RegularizationInitial.size() < n)
178  {
179  m_RegularizationInitial.resize(n);
180  }
181  if (m_RegularizationFinal.size() < n)
182  {
183  m_RegularizationFinal.resize(n);
184  }
185  if (m_RegularizationDecaySpan.size() < n)
186  {
187  m_RegularizationDecaySpan.resize(n);
188  }
189 
190  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetRegularizationConstant ( double  r)
inline

Directly set the regularization. This bypasses all other methods, but will only be valid until the current scale is optimized.

Definition at line 461 of file itkPSMEntropyModelFilter.h.

462  {
463  m_ShapeEntropyFunction->SetMinimumVariance(r);
464  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetRegularizationInitial ( const std::vector< double > &  v)
inline

A regularization constant is added to the diagonal elements of the covariance matrix of the particle positions in shape space during the optimization. It is called "minimum variance" in the PSMShapeEntropyFunction class because it enforces a minimum variance in each of the eigen-modes of the shape space during the optimization. Higher MinimumVariance produces larger gradients, thus allowing shapes to move larger distances in shape space during optimization of positions with this function. Higher MinimumVariance values minimize the effects of smaller modes of variation in the covariance matrix. Conversely, smaller MinimumVariance values increase the significance of the smaller modes of variation on the gradient calculation. The parameters RegularizationInitial, RegularizationFinal, and RegularizationDecaySpan define an "annealing" type optimization, where the minimum variance starts at RegularizationInitial and exponentially decreases over RegularizationDecaySpan to the RegularizationFinal value. If the RegularizationDecaySpan is set to 0, then only RegularizationConstant is used as the MinimumVariance parameter and is held fixed during the optimization.

Definition at line 367 of file itkPSMEntropyModelFilter.h.

368  {
369  if (v.size() != m_RegularizationInitial.size())
370  {
371  itkExceptionMacro("Number of variables does not match number of scales. Call SetNumberOfScales first.");
372  }
373  m_RegularizationInitial = v;
374  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetShapeEntropyWeighting ( double  w)
inline

This method sets/gets the relative weighting of terms in the two-term cost function used for the optimization. This parameter adjusts the weighting of the shape representation term with respect to the strength of the entropy minimization in shape space. Higher values will favor a more compact representation in shape space and produce less regular samplings of shape surfaces. Default value is an equal weighting between the two terms (1.0).

Definition at line 299 of file itkPSMEntropyModelFilter.h.

300  {
301  m_CostFunction->SetRelativeGradientScaling(w);
302  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
void itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::SetTolerance ( const std::vector< double > &  v)
inline

Set the tolerance for the optimization, below which the optimization will be considered to have converged. Tolerance is specified as a distance in the implicit surface images. The solver will converge when no pixel has moved more than the specified tolerance in a given iteration. An optional argument specifies the scale to which you want this tolerance applied (see SetNumberOfScales). If NumberOfScales > 1, then you MUST supply a value for each scale.

Definition at line 486 of file itkPSMEntropyModelFilter.h.

487  {
488  if (v.size() != m_Tolerances.size())
489  {
490  itkExceptionMacro("Number of variables does not match number of scales. Call SetNumberOfScales first.");
491  }
492  m_Tolerances = v;
493  }
template<class TImage, class TShapeMatrix = PSMShapeMatrixAttribute<double, TImage::ImageDimension>>
PSMShapeEntropyFunction<Dimension>* itk::PSMEntropyModelFilter< TImage, TShapeMatrix >::ShapeEntropyFunction ( )
inline

Returns a pointer to the ShapeEntropyFunction. This is the cost function that will minimize the entropy of the shape vectors in high-dimensional shape space. This is the cost function term that establishes correspondence across shapes.

Definition at line 229 of file itkPSMEntropyModelFilter.h.

230  {
231  return m_ShapeEntropyFunction;
232  }

The documentation for this class was generated from the following files: