Shapeworks Studio  2.1
Shape analysis software suite
itkPSMMeanCurvatureAttribute.h
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 #ifndef __itkPSMMeanCurvatureAttribute_h
19 #define __itkPSMMeanCurvatureAttribute_h
20 
21 #include "itkDataObject.h"
22 #include "itkWeakPointer.h"
23 #include "itkPSMContainer.h"
24 #include "itkPSMContainerArrayAttribute.h"
25 #include "itkPSMImageDomainWithCurvature.h"
26 #include "itkPSMParticleSystem.h"
27 
28 namespace itk
29 {
33 template <class TNumericType, unsigned int VDimension>
34 class ITK_EXPORT PSMMeanCurvatureAttribute
35  : public PSMContainerArrayAttribute<TNumericType,VDimension>
36 {
37 public:
39  typedef TNumericType NumericType;
42  typedef SmartPointer<Self> Pointer;
43  typedef SmartPointer<const Self> ConstPointer;
44  typedef WeakPointer<const Self> ConstWeakPointer;
45 
49  typedef vnl_vector_fixed<TNumericType, VDimension> VnlVectorType;
50 
52  itkNewMacro(Self);
53 
56 
57  virtual void PositionAddEventCallback(Object *o, const EventObject &e)
58  {
59  Superclass::PositionAddEventCallback(o, e);
60  const ParticlePositionAddEvent &event = dynamic_cast<const ParticlePositionAddEvent &>(e);
61  const ParticleSystemType *ps = dynamic_cast<const ParticleSystemType *>(o);
62  this->ComputeMeanCurvature(ps, event.GetPositionIndex(), event.GetDomainIndex());
63  }
64 
65  virtual void PositionSetEventCallback(Object *o, const EventObject &e)
66  {
67  const ParticlePositionSetEvent &event = dynamic_cast<const ParticlePositionSetEvent &>(e);
68  const ParticleSystemType *ps= dynamic_cast<const ParticleSystemType *>(o);
69  this->ComputeMeanCurvature(ps, event.GetPositionIndex(), event.GetDomainIndex());
70  }
71 
72  virtual void DomainAddEventCallback(Object *o, const EventObject &e)
73  {
74  Superclass::DomainAddEventCallback(o, e);
75  m_MeanCurvatureList.push_back(0.0);
76  m_CurvatureStandardDeviationList.push_back(0.0);
77  const ParticleDomainAddEvent &event = dynamic_cast<const ParticleDomainAddEvent &>(e);
78  const ParticleSystemType *ps= dynamic_cast<const ParticleSystemType *>(o);
79  this->ComputeCurvatureStatistics(ps, event.GetDomainIndex());
80  }
81 
83  inline void ComputeMeanCurvature(const ParticleSystemType *system,
84  unsigned int idx, unsigned int dom)
85  {
86  // First we need a pointer to the domain. Assume we have
87  // Curvature information.
90  system->GetDomain(dom) );
91  // Get the position and index.
92  PointType pos = system->GetPosition(idx, dom);
93  this->operator[](dom)->operator[](idx) = domain->GetCurvature(pos);
94  }
95 
98  virtual void ComputeCurvatureStatistics(const ParticleSystemType *, unsigned int d);
99 
100  double GetMeanCurvature(int d)
101  { return m_MeanCurvatureList[d]; }
102  double GetCurvatureStandardDeviation(int d)
103  { return m_CurvatureStandardDeviationList[d];}
104 
105 protected:
107  {
108  this->m_DefinedCallbacks.PositionSetEvent = true;
109  this->m_DefinedCallbacks.DomainAddEvent = true;
110  }
111  virtual ~PSMMeanCurvatureAttribute() {};
112 
113  void PrintSelf(std::ostream& os, Indent indent) const
114  { Superclass::PrintSelf(os,indent); }
115 
116 private:
117  PSMMeanCurvatureAttribute(const Self&); //purposely not implemented
118  void operator=(const Self&); //purposely not implemented
119 
120  std::vector<double> m_MeanCurvatureList;
121  std::vector<double> m_CurvatureStandardDeviationList;
122 
123 };
124 
125 } // end namespace
126 
127 #ifndef ITK_MANUAL_INSTANTIATION
128 #include "itkPSMMeanCurvatureAttribute.hxx"
129 #endif
130 
131 #endif
PointType & GetPosition(unsigned long int k, unsigned int d=0)
virtual void DomainAddEventCallback(Object *o, const EventObject &e)
A facade class that manages interactions with a particle system.
DomainType * GetDomain(unsigned int i)
Base class for PSMParticleSystem attribute classes.
Definition: Shape.h:14
PSMParticleSystem< VDimension > ParticleSystemType