go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkAdvancedBSplineDeformableTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright UMC Utrecht and contributors
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 
20  Program: Insight Segmentation & Registration Toolkit
21  Module: $RCSfile: itkAdvancedBSplineDeformableTransform.h,v $
22  Language: C++
23  Date: $Date: 2008-04-11 16:28:11 $
24  Version: $Revision: 1.38 $
25 
26  Copyright (c) Insight Software Consortium. All rights reserved.
27  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
28 
29  This software is distributed WITHOUT ANY WARRANTY; without even
30  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
31  PURPOSE. See the above copyright notices for more information.
32 
33 =========================================================================*/
34 #ifndef __itkAdvancedBSplineDeformableTransform_h
35 #define __itkAdvancedBSplineDeformableTransform_h
36 
38 
39 #include "itkImage.h"
40 #include "itkImageRegion.h"
44 
45 namespace itk
46 {
47 
48 // Forward declarations for friendship
49 template< class TScalarType, unsigned int NDimensions, unsigned int VSplineOrder >
50 class MultiBSplineDeformableTransformWithNormal;
51 
129 template<
130 class TScalarType = double, // Data type for scalars
131 unsigned int NDimensions = 3, // Number of dimensions
132 unsigned int VSplineOrder = 3 >
133 // Spline order
135  public AdvancedBSplineDeformableTransformBase< TScalarType, NDimensions >
136 {
137 public:
138 
142  TScalarType, NDimensions > Superclass;
143  typedef SmartPointer< Self > Pointer;
144  typedef SmartPointer< const Self > ConstPointer;
145 
147  itkNewMacro( Self );
148 
151 
153  itkStaticConstMacro( SpaceDimension, unsigned int, NDimensions );
154 
156  itkStaticConstMacro( SplineOrder, unsigned int, VSplineOrder );
157 
171  typedef typename Superclass::InputCovariantVectorType
175 
176  typedef typename Superclass
179  typedef typename Superclass
182  typedef typename Superclass
187 
189  typedef typename Superclass::PixelType PixelType;
192 
195 
197  typedef typename Superclass::SizeType SizeType;
202 
204  void SetGridRegion( const RegionType & region ) override;
205 
207  OutputPointType TransformPoint( const InputPointType & point ) const override;
208 
211  itkGetStaticConstMacro( SpaceDimension ),
212  itkGetStaticConstMacro( SplineOrder ) > WeightsFunctionType;
217  ScalarType,
218  itkGetStaticConstMacro( SpaceDimension ),
219  itkGetStaticConstMacro( SplineOrder ) > DerivativeWeightsFunctionType;
222  ScalarType,
223  itkGetStaticConstMacro( SpaceDimension ),
224  itkGetStaticConstMacro( SplineOrder ) > SODerivativeWeightsFunctionType;
226 
229 
237  virtual void TransformPoint(
238  const InputPointType & inputPoint,
239  OutputPointType & outputPoint,
240  WeightsType & weights,
241  ParameterIndexArrayType & indices,
242  bool & inside ) const;
243 
245  unsigned long GetNumberOfWeights( void ) const
246  {
247  return this->m_WeightsFunction->GetNumberOfWeights();
248  }
249 
250 
251  unsigned int GetNumberOfAffectedWeights( void ) const override;
252 
253  NumberOfParametersType GetNumberOfNonZeroJacobianIndices( void ) const override;
254 
257  const InputPointType & ipp,
258  JacobianType & j,
259  NonZeroJacobianIndicesType & nzji ) const override;
260 
265  const InputPointType & ipp,
266  const MovingImageGradientType & movingImageGradient,
267  DerivativeType & imageJacobian,
268  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const override;
269 
272  const InputPointType & ipp,
273  SpatialJacobianType & sj ) const override;
274 
277  const InputPointType & ipp,
278  SpatialHessianType & sh ) const override;
279 
282  const InputPointType & ipp,
284  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const override;
285 
290  const InputPointType & ipp,
291  SpatialJacobianType & sj,
293  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const override;
294 
297  const InputPointType & ipp,
299  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const override;
300 
305  const InputPointType & ipp,
306  SpatialHessianType & sh,
308  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const override;
309 
310 protected:
311 
313  void PrintSelf( std::ostream & os, Indent indent ) const override;
314 
317 
319  // Why??
320  itkSetObjectMacro( WeightsFunction, WeightsFunctionType );
322 
324  void WrapAsImages( void );
325 
327  NonZeroJacobianIndicesType & nonZeroJacobianIndices,
328  const RegionType & supportRegion ) const override;
329 
332 
338  std::vector< DerivativeWeightsFunctionPointer > m_DerivativeWeightsFunctions;
339  std::vector< std::vector< SODerivativeWeightsFunctionPointer > > m_SODerivativeWeightsFunctions;
340 
341 private:
342 
343  AdvancedBSplineDeformableTransform( const Self & ); // purposely not implemented
344  void operator=( const Self & ); // purposely not implemented
345 
347  itkGetStaticConstMacro( SpaceDimension ),
348  itkGetStaticConstMacro( SplineOrder ) >;
349 
350 };
351 
352 } // namespace itk
353 
354 #ifndef ITK_MANUAL_INSTANTIATION
355 #include "itkAdvancedBSplineDeformableTransform.hxx"
356 #endif
357 
358 #endif /* __itkAdvancedBSplineDeformableTransform_h */
Base class for deformable transform using a B-spline representation.
Superclass ::JacobianOfSpatialJacobianType JacobianOfSpatialJacobianType
Superclass ::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
Superclass::MovingImageGradientValueType MovingImageGradientValueType
Superclass ::JacobianOfSpatialHessianType JacobianOfSpatialHessianType
ImageRegion< itkGetStaticConstMacro(SpaceDimension) > RegionType
Deformable transform using a B-spline representation.
AdvancedBSplineDeformableTransformBase< TScalarType, NDimensions > Superclass
void GetSpatialHessian(const InputPointType &ipp, SpatialHessianType &sh) const override
void GetSpatialJacobian(const InputPointType &ipp, SpatialJacobianType &sj) const override
Superclass ::JacobianOfSpatialHessianType JacobianOfSpatialHessianType
BSplineInterpolationDerivativeWeightFunction< ScalarType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SplineOrder) > DerivativeWeightsFunctionType
virtual void TransformPoint(const InputPointType &inputPoint, OutputPointType &outputPoint, WeightsType &weights, ParameterIndexArrayType &indices, bool &inside) const
void GetJacobianOfSpatialHessian(const InputPointType &ipp, JacobianOfSpatialHessianType &jsh, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void SetGridRegion(const RegionType &region) override
Superclass::MovingImageGradientValueType MovingImageGradientValueType
DerivativeWeightsFunctionType::Pointer DerivativeWeightsFunctionPointer
std::vector< std::vector< SODerivativeWeightsFunctionPointer > > m_SODerivativeWeightsFunctions
itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType)
Superclass::InputCovariantVectorType InputCovariantVectorType
WeightsFunctionType::ContinuousIndexType ContinuousIndexType
void GetJacobianOfSpatialJacobian(const InputPointType &ipp, JacobianOfSpatialJacobianType &jsj, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
BSplineInterpolationSecondOrderDerivativeWeightFunction< ScalarType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SplineOrder) > SODerivativeWeightsFunctionType
void EvaluateJacobianWithImageGradientProduct(const InputPointType &ipp, const MovingImageGradientType &movingImageGradient, DerivativeType &imageJacobian, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void ComputeNonZeroJacobianIndices(NonZeroJacobianIndicesType &nonZeroJacobianIndices, const RegionType &supportRegion) const override
Superclass ::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
void GetJacobianOfSpatialJacobian(const InputPointType &ipp, SpatialJacobianType &sj, JacobianOfSpatialJacobianType &jsj, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
std::vector< DerivativeWeightsFunctionPointer > m_DerivativeWeightsFunctions
SODerivativeWeightsFunctionType::Pointer SODerivativeWeightsFunctionPointer
Superclass::OutputCovariantVectorType OutputCovariantVectorType
Superclass ::JacobianOfSpatialJacobianType JacobianOfSpatialJacobianType
NumberOfParametersType GetNumberOfNonZeroJacobianIndices(void) const override
void GetJacobianOfSpatialHessian(const InputPointType &ipp, SpatialHessianType &sh, JacobianOfSpatialHessianType &jsh, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
BSplineInterpolationWeightFunction2< ScalarType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SplineOrder) > WeightsFunctionType
void PrintSelf(std::ostream &os, Indent indent) const override
OutputPointType TransformPoint(const InputPointType &point) const override
void GetJacobian(const InputPointType &ipp, JacobianType &j, NonZeroJacobianIndicesType &nzji) const override
itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions)
itkStaticConstMacro(SplineOrder, unsigned int, VSplineOrder)
unsigned int GetNumberOfAffectedWeights(void) const override
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Returns the weights over the support region used for B-spline interpolation/reconstruction.
This transform is a composition of B-spline transformations, allowing sliding motion between differen...


Generated on OURCE_DATE_EPOCH for elastix by doxygen 1.9.1 elastix logo