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 #include "itkImage.h"
39 #include "itkImageRegion.h"
43 
44 namespace itk
45 {
46 
47 // Forward declarations for friendship
48 template< class TScalarType, unsigned int NDimensions, unsigned int VSplineOrder >
50 
128 template<
129 class TScalarType = double, // Data type for scalars
130 unsigned int NDimensions = 3, // Number of dimensions
131 unsigned int VSplineOrder = 3 >
132 // Spline order
134  public AdvancedBSplineDeformableTransformBase< TScalarType, NDimensions >
135 {
136 public:
137 
141  TScalarType, NDimensions > Superclass;
143  typedef SmartPointer< const Self > ConstPointer;
144 
146  itkNewMacro( Self );
147 
150 
152  itkStaticConstMacro( SpaceDimension, unsigned int, NDimensions );
153 
155  itkStaticConstMacro( SplineOrder, unsigned int, VSplineOrder );
156 
170  typedef typename Superclass::InputCovariantVectorType
174 
175  typedef typename Superclass
178  typedef typename Superclass
181  typedef typename Superclass
186 
191 
194 
196  typedef typename Superclass::SizeType SizeType;
201 
203  virtual void SetGridRegion( const RegionType & region );
204 
206  OutputPointType TransformPoint( const InputPointType & point ) const;
207 
210  itkGetStaticConstMacro( SpaceDimension ),
211  itkGetStaticConstMacro( SplineOrder ) > WeightsFunctionType;
216  ScalarType,
217  itkGetStaticConstMacro( SpaceDimension ),
218  itkGetStaticConstMacro( SplineOrder ) > DerivativeWeightsFunctionType;
221  ScalarType,
222  itkGetStaticConstMacro( SpaceDimension ),
223  itkGetStaticConstMacro( SplineOrder ) > SODerivativeWeightsFunctionType;
225 
228 
236  virtual void TransformPoint(
237  const InputPointType & inputPoint,
238  OutputPointType & outputPoint,
239  WeightsType & weights,
240  ParameterIndexArrayType & indices,
241  bool & inside ) const;
242 
244  unsigned long GetNumberOfWeights( void ) const
245  {
246  return this->m_WeightsFunction->GetNumberOfWeights();
247  }
248 
249 
250  unsigned int GetNumberOfAffectedWeights( void ) const;
251 
252  virtual NumberOfParametersType GetNumberOfNonZeroJacobianIndices( void ) const;
253 
255  virtual void GetJacobian(
256  const InputPointType & ipp,
257  JacobianType & j,
258  NonZeroJacobianIndicesType & nzji ) const;
259 
264  const InputPointType & ipp,
265  const MovingImageGradientType & movingImageGradient,
266  DerivativeType & imageJacobian,
267  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const;
268 
270  virtual void GetSpatialJacobian(
271  const InputPointType & ipp,
272  SpatialJacobianType & sj ) const;
273 
275  virtual void GetSpatialHessian(
276  const InputPointType & ipp,
277  SpatialHessianType & sh ) const;
278 
280  virtual void GetJacobianOfSpatialJacobian(
281  const InputPointType & ipp,
282  JacobianOfSpatialJacobianType & jsj,
283  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const;
284 
288  virtual void GetJacobianOfSpatialJacobian(
289  const InputPointType & ipp,
290  SpatialJacobianType & sj,
291  JacobianOfSpatialJacobianType & jsj,
292  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const;
293 
295  virtual void GetJacobianOfSpatialHessian(
296  const InputPointType & ipp,
297  JacobianOfSpatialHessianType & jsh,
298  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const;
299 
303  virtual void GetJacobianOfSpatialHessian(
304  const InputPointType & ipp,
305  SpatialHessianType & sh,
306  JacobianOfSpatialHessianType & jsh,
307  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const;
308 
309 protected:
310 
312  virtual void PrintSelf( std::ostream & os, Indent indent ) const;
313 
316 
318  // Why??
319  itkSetObjectMacro( WeightsFunction, WeightsFunctionType );
320  itkGetObjectMacro( WeightsFunction, WeightsFunctionType );
321 
323  void WrapAsImages( void );
324 
325  virtual void ComputeNonZeroJacobianIndices(
326  NonZeroJacobianIndicesType & nonZeroJacobianIndices,
327  const RegionType & supportRegion ) const;
328 
331 
336  WeightsFunctionPointer m_WeightsFunction;
337  std::vector< DerivativeWeightsFunctionPointer > m_DerivativeWeightsFunctions;
338  std::vector< std::vector< SODerivativeWeightsFunctionPointer > > m_SODerivativeWeightsFunctions;
339 
340 private:
341 
342  AdvancedBSplineDeformableTransform( const Self & ); // purposely not implemented
343  void operator=( const Self & ); // purposely not implemented
344 
345  friend class MultiBSplineDeformableTransformWithNormal< ScalarType,
346  itkGetStaticConstMacro( SpaceDimension ),
347  itkGetStaticConstMacro( SplineOrder ) >;
348 
349 };
350 
351 } // namespace itk
352 
353 #ifndef ITK_MANUAL_INSTANTIATION
354 #include "itkAdvancedBSplineDeformableTransform.hxx"
355 #endif
356 
357 #endif /* __itkAdvancedBSplineDeformableTransform_h */
Matrix< ScalarType, OutputSpaceDimension, InputSpaceDimension > SpatialJacobianType
std::vector< DerivativeWeightsFunctionPointer > m_DerivativeWeightsFunctions
Deformable transform using a B-spline representation.
itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions)
virtual void SetGridRegion(const RegionType &region)
Superclass::InputCovariantVectorType InputCovariantVectorType
virtual void GetJacobian(const InputPointType &ipp, JacobianType &j, NonZeroJacobianIndicesType &nzji) const
BSplineInterpolationSecondOrderDerivativeWeightFunction< ScalarType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SplineOrder) > SODerivativeWeightsFunctionType
virtual void GetJacobianOfSpatialJacobian(const InputPointType &ipp, JacobianOfSpatialJacobianType &jsj, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const
virtual void ComputeNonZeroJacobianIndices(NonZeroJacobianIndicesType &nonZeroJacobianIndices, const RegionType &supportRegion) const
Superclass::JacobianOfSpatialHessianType JacobianOfSpatialHessianType
ImageRegion< itkGetStaticConstMacro(SpaceDimension) > RegionType
Superclass::OutputCovariantVectorType OutputCovariantVectorType
BSplineInterpolationWeightFunction2< ScalarType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SplineOrder) > WeightsFunctionType
BSplineInterpolationDerivativeWeightFunction< ScalarType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SplineOrder) > DerivativeWeightsFunctionType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
AdvancedBSplineDeformableTransformBase< TScalarType, NDimensions > Superclass
virtual void PrintSelf(std::ostream &os, Indent indent) const
Image< JacobianPixelType, itkGetStaticConstMacro(SpaceDimension) > JacobianImageType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
virtual void EvaluateJacobianWithImageGradientProduct(const InputPointType &ipp, const MovingImageGradientType &movingImageGradient, DerivativeType &imageJacobian, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const
Superclass::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
std::vector< std::vector< SODerivativeWeightsFunctionPointer > > m_SODerivativeWeightsFunctions
WeightsFunctionType::ContinuousIndexType ContinuousIndexType
unsigned int GetNumberOfAffectedWeights(void) const
OutputPointType TransformPoint(const InputPointType &point) const
This transform is a composition of B-spline transformations, allowing sliding motion between differen...
Image< PixelType, itkGetStaticConstMacro(SpaceDimension) > ImageType
Superclass::JacobianOfSpatialJacobianType JacobianOfSpatialJacobianType
virtual void GetSpatialJacobian(const InputPointType &ipp, SpatialJacobianType &sj) const
virtual void GetSpatialHessian(const InputPointType &ipp, SpatialHessianType &sh) const
DerivativeWeightsFunctionType::Pointer DerivativeWeightsFunctionPointer
Superclass::MovingImageGradientValueType MovingImageGradientValueType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
SODerivativeWeightsFunctionType::Pointer SODerivativeWeightsFunctionPointer
Base class for deformable transform using a B-spline representation.
virtual void GetJacobianOfSpatialHessian(const InputPointType &ipp, JacobianOfSpatialHessianType &jsh, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const
FixedArray< Matrix< ScalarType, InputSpaceDimension, InputSpaceDimension >, OutputSpaceDimension > SpatialHessianType
virtual NumberOfParametersType GetNumberOfNonZeroJacobianIndices(void) const


Generated on 04-09-2015 for elastix by doxygen 1.8.9.1 elastix logo