ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.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 #ifndef itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor_h
19 #define itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor_h
20 
22 
23 namespace itk
24 {
65 template<typename TTransform>
67 : public TransformParametersAdaptor<TTransform>
68 {
69 public:
70 
76 
78  itkNewMacro( Self );
79 
82 
84  typedef TTransform TransformType;
85  typedef typename TransformType::Pointer TransformPointer;
86  typedef typename TransformType::ParametersType ParametersType;
87  typedef typename TransformType::ParametersValueType ParametersValueType;
88  typedef typename TransformType::FixedParametersType FixedParametersType;
89  typedef typename TransformType::FixedParametersValueType FixedParametersValueType;
90 
91  typedef typename TransformType::TimeVaryingVelocityFieldControlPointLatticeType TimeVaryingVelocityFieldControlPointLatticeType;
92  typedef typename TimeVaryingVelocityFieldControlPointLatticeType::Pointer TimeVaryingVelocityFieldControlPointLatticePointer;
93  typedef typename TimeVaryingVelocityFieldControlPointLatticeType::RegionType RegionType;
95  typedef typename TimeVaryingVelocityFieldControlPointLatticeType::PixelType VectorType;
97  typedef typename TimeVaryingVelocityFieldControlPointLatticeType::SpacingType SpacingType;
102 
104  itkStaticConstMacro( TotalDimension, unsigned int, TransformType::Dimension + 1 );
105 
107  itkSetMacro( SplineOrder, SizeValueType );
108 
110  itkGetConstMacro( SplineOrder, SizeValueType );
111 
113  void SetRequiredTransformDomainMeshSize( const MeshSizeType & );
114 
116  itkGetConstReferenceMacro( RequiredTransformDomainMeshSize, MeshSizeType );
117 
119  void SetRequiredTransformDomainSize( const SizeType & );
120 
122  itkGetConstReferenceMacro( RequiredTransformDomainSize, SizeType );
123 
125  void SetRequiredTransformDomainSpacing( const SpacingType & );
126 
128  itkGetConstReferenceMacro( RequiredTransformDomainSpacing, SpacingType );
129 
131  void SetRequiredTransformDomainOrigin( const OriginType & );
132 
134  itkGetConstReferenceMacro( RequiredTransformDomainOrigin, OriginType );
135 
137  void SetRequiredTransformDomainDirection( const DirectionType & );
138 
140  itkGetConstReferenceMacro( RequiredTransformDomainDirection, DirectionType );
141 
144  {
145  OriginType requiredLatticeOrigin;
146  for( SizeValueType i = 0; i < TotalDimension; i++ )
147  {
148  requiredLatticeOrigin[i] = this->m_RequiredFixedParameters[TotalDimension + i];
149  }
150  return requiredLatticeOrigin;
151  }
153 
156  {
157  SpacingType requiredLatticeSpacing;
158  for( SizeValueType i = 0; i < TotalDimension; i++ )
159  {
160  FixedParametersValueType domainPhysicalDimensions = static_cast<FixedParametersValueType>( this->m_RequiredTransformDomainSize[i] - 1.0 ) *
161  this->m_RequiredTransformDomainSpacing[i];
162  requiredLatticeSpacing[i] = domainPhysicalDimensions / static_cast<FixedParametersValueType>( this->m_RequiredTransformDomainMeshSize[i] );
163  }
164  return requiredLatticeSpacing;
165  }
167 
170  {
171  SizeType requiredLatticeSize;
172  for( SizeValueType i = 0; i < TotalDimension; i++ )
173  {
174  requiredLatticeSize[i] = static_cast<SizeValueType>( this->m_RequiredFixedParameters[i] );
175  }
176  return requiredLatticeSize;
177  }
179 
182  {
183  return this->m_RequiredTransformDomainDirection;
184  }
185 
187  virtual void AdaptTransformParameters() ITK_OVERRIDE;
188 
189  virtual void SetRequiredFixedParameters( const FixedParametersType ) ITK_OVERRIDE;
190 
191 protected:
193  ~TimeVaryingBSplineVelocityFieldTransformParametersAdaptor() ITK_OVERRIDE;
194 
195  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
196 
197 private:
198  ITK_DISALLOW_COPY_AND_ASSIGN(TimeVaryingBSplineVelocityFieldTransformParametersAdaptor);
199 
201  void UpdateRequiredFixedParameters();
202 
203  MeshSizeType m_RequiredTransformDomainMeshSize;
204  OriginType m_RequiredTransformDomainOrigin;
205  DirectionType m_RequiredTransformDomainDirection;
206  SpacingType m_RequiredTransformDomainSpacing;
207  SizeType m_RequiredTransformDomainSize;
208 
209  SizeValueType m_SplineOrder;
210 
211 }; //class TimeVaryingBSplineVelocityFieldTransformParametersAdaptor
212 } // namespace itk
213 
214 #ifndef ITK_MANUAL_INSTANTIATION
215 #include "itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.hxx"
216 #endif
217 
218 #endif /* itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor_h */
Light weight base class for most itk classes.
Base helper class intended for multi-resolution image registration.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
TransformType::TimeVaryingVelocityFieldControlPointLatticeType TimeVaryingVelocityFieldControlPointLatticeType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
TimeVaryingBSplineVelocityFieldTransformParametersAdaptor is a helper class intended to definition...