ITK  5.0.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:
71 
77 
79  itkNewMacro( Self );
80 
83 
85  using TransformType = TTransform;
86  using TransformPointer = typename TransformType::Pointer;
87  using ParametersType = typename TransformType::ParametersType;
88  using ParametersValueType = typename TransformType::ParametersValueType;
89  using FixedParametersType = typename TransformType::FixedParametersType;
90  using FixedParametersValueType = typename TransformType::FixedParametersValueType;
91 
92  using TimeVaryingVelocityFieldControlPointLatticeType = typename TransformType::TimeVaryingVelocityFieldControlPointLatticeType;
93  using TimeVaryingVelocityFieldControlPointLatticePointer = typename TimeVaryingVelocityFieldControlPointLatticeType::Pointer;
96  using VectorType = typename TimeVaryingVelocityFieldControlPointLatticeType::PixelType;
98  using SpacingType = typename TimeVaryingVelocityFieldControlPointLatticeType::SpacingType;
103 
105  static constexpr unsigned int TotalDimension = TransformType::Dimension + 1;
106 
108  itkSetMacro( SplineOrder, SizeValueType );
109 
111  itkGetConstMacro( SplineOrder, SizeValueType );
112 
114  void SetRequiredTransformDomainMeshSize( const MeshSizeType & );
115 
117  itkGetConstReferenceMacro( RequiredTransformDomainMeshSize, MeshSizeType );
118 
120  void SetRequiredTransformDomainSize( const SizeType & );
121 
123  itkGetConstReferenceMacro( RequiredTransformDomainSize, SizeType );
124 
126  void SetRequiredTransformDomainSpacing( const SpacingType & );
127 
129  itkGetConstReferenceMacro( RequiredTransformDomainSpacing, SpacingType );
130 
132  void SetRequiredTransformDomainOrigin( const OriginType & );
133 
135  itkGetConstReferenceMacro( RequiredTransformDomainOrigin, OriginType );
136 
138  void SetRequiredTransformDomainDirection( const DirectionType & );
139 
141  itkGetConstReferenceMacro( RequiredTransformDomainDirection, DirectionType );
142 
145  {
146  OriginType requiredLatticeOrigin;
147  for( SizeValueType i = 0; i < TotalDimension; i++ )
148  {
149  requiredLatticeOrigin[i] = this->m_RequiredFixedParameters[TotalDimension + i];
150  }
151  return requiredLatticeOrigin;
152  }
154 
157  {
158  SpacingType requiredLatticeSpacing;
159  for( SizeValueType i = 0; i < TotalDimension; i++ )
160  {
161  FixedParametersValueType domainPhysicalDimensions = static_cast<FixedParametersValueType>( this->m_RequiredTransformDomainSize[i] - 1.0 ) *
162  this->m_RequiredTransformDomainSpacing[i];
163  requiredLatticeSpacing[i] = domainPhysicalDimensions / static_cast<FixedParametersValueType>( this->m_RequiredTransformDomainMeshSize[i] );
164  }
165  return requiredLatticeSpacing;
166  }
168 
171  {
172  SizeType requiredLatticeSize;
173  for( SizeValueType i = 0; i < TotalDimension; i++ )
174  {
175  requiredLatticeSize[i] = static_cast<SizeValueType>( this->m_RequiredFixedParameters[i] );
176  }
177  return requiredLatticeSize;
178  }
180 
183  {
184  return this->m_RequiredTransformDomainDirection;
185  }
186 
188  void AdaptTransformParameters() override;
189 
190  void SetRequiredFixedParameters( const FixedParametersType ) override;
191 
192 protected:
195 
196  void PrintSelf( std::ostream& os, Indent indent ) const override;
197 
198 private:
200  void UpdateRequiredFixedParameters();
201 
207 
209 
210 }; //class TimeVaryingBSplineVelocityFieldTransformParametersAdaptor
211 } // namespace itk
212 
213 #ifndef ITK_MANUAL_INSTANTIATION
214 #include "itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.hxx"
215 #endif
216 
217 #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:83
typename TimeVaryingVelocityFieldControlPointLatticeType::Pointer TimeVaryingVelocityFieldControlPointLatticePointer
typename TimeVaryingVelocityFieldControlPointLatticeType::DirectionType DirectionType
typename TimeVaryingVelocityFieldControlPointLatticeType::SizeValueType SizeValueType
typename TransformType::TimeVaryingVelocityFieldControlPointLatticeType TimeVaryingVelocityFieldControlPointLatticeType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
TimeVaryingBSplineVelocityFieldTransformParametersAdaptor is a helper class intended to definition...