]> Creatis software - clitk.git/blob - registration/clitkBSplineDeformableTransform.h
*** empty log message ***
[clitk.git] / registration / clitkBSplineDeformableTransform.h
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to: 
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://oncora1.lyon.fnclcc.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17 ======================================================================-====*/
18 #ifndef __clitkBSplineDeformableTransform_h
19 #define __clitkBSplineDeformableTransform_h
20 #include "clitkVectorBSplineResampleImageFunctionWithLUT.h"
21
22 //itk include
23 #include "itkTransform.h"
24 #include "itkImage.h"
25 #include "itkImageRegion.h"
26 #include "itkSpatialObject.h"
27
28 namespace clitk
29 {
30
31   template <
32     class TCoordRep = double,               // Data type for scalars, coordinate representation,vectors
33     unsigned int NInputDimensions = 3,        // Number of input dimensions
34     unsigned int NOutputDimensions = 3 >      // Number of output dimensions
35   class ITK_EXPORT BSplineDeformableTransform : 
36     public itk::Transform< TCoordRep, NInputDimensions, NOutputDimensions >
37   {
38   public:
39
40     //====================================================================
41     // Typedefs
42     //====================================================================
43     typedef BSplineDeformableTransform                         Self;
44     typedef itk::Transform< TCoordRep, NInputDimensions, NOutputDimensions > Superclass;
45     typedef itk::SmartPointer<Self>                                 Pointer;
46     typedef itk::SmartPointer<const Self>                           ConstPointer;
47       
48     /** New macro for creation of through the object factory.*/
49     itkNewMacro( Self );
50
51     /** Run-time type information (and related methods). */
52     itkTypeMacro( BSplineDeformableTransform, Transform );
53
54     /** Dimension of the domain space. */
55     itkStaticConstMacro(OutputDimension, unsigned int, NOutputDimensions);
56
57     /** Dimension of the input model. */
58     itkStaticConstMacro(InputDimension, unsigned int, NInputDimensions);
59
60     /** Standard scalar type for this class. */
61     typedef typename Superclass::ScalarType ScalarType;
62
63     /** Standard parameters container. */
64     typedef typename Superclass::ParametersType ParametersType;
65
66     /** Standard Jacobian container. */
67     typedef typename Superclass::JacobianType JacobianType;
68
69     /** Standard vector type for this class. */
70     typedef itk::Vector<TCoordRep,
71                         InputDimension> InputVectorType;
72     typedef itk::Vector<TCoordRep,
73                         OutputDimension> OutputVectorType;
74   
75     /** Standard covariant vector type for this class. */
76     typedef itk::CovariantVector<TCoordRep,
77                                  InputDimension> InputCovariantVectorType;
78     typedef itk::CovariantVector<TCoordRep,
79                                  OutputDimension> OutputCovariantVectorType;
80
81     /** Standard vnl_vector type for this class. */
82     typedef vnl_vector_fixed<TCoordRep,
83                              InputDimension> InputVnlVectorType;
84     typedef vnl_vector_fixed<TCoordRep,
85                              OutputDimension> OutputVnlVectorType;
86
87     /** Standard coordinate point type for this class. */
88     typedef itk::Point<TCoordRep,
89                        InputDimension> InputPointType;
90     typedef itk::Point<TCoordRep,
91                        OutputDimension> OutputPointType;
92
93     //JV Parameters as images with OutputDimension number of components per Pixel
94     typedef typename ParametersType::ValueType                          ParametersValueType;
95     typedef typename itk::Vector<ParametersValueType, OutputDimension>  PixelType;
96     typedef itk::Image<PixelType, InputDimension>                       CoefficientImageType;
97     typedef typename CoefficientImageType::Pointer                      CoefficientImagePointer;
98
99   
100     /** Typedefs for specifying the extend to the grid. */
101     typedef itk::ImageRegion<InputDimension>                RegionType;
102     typedef typename RegionType::IndexType                  IndexType;
103     typedef typename RegionType::SizeType                   SizeType;
104     typedef typename CoefficientImageType::SpacingType      SpacingType;
105     typedef typename CoefficientImageType::DirectionType    DirectionType;
106     typedef typename CoefficientImageType::PointType        OriginType;
107     typedef itk::ContinuousIndex<TCoordRep, InputDimension> ContinuousIndexType;
108   
109     //JV added for the BLUT interpolator
110     //typedef itk::Vector<TCoordRep, InputDimension> OutputSpacingType;
111
112     //JV m_VectorInterpolator
113     typedef VectorBSplineResampleImageFunctionWithLUT
114     <CoefficientImageType, TCoordRep> VectorInterpolatorType;
115     typedef  typename  VectorInterpolatorType::CoefficientDataType CoefficientDataType;
116     typedef  typename  VectorInterpolatorType::CoefficientDataType WeightsDataType;
117
118     /** Typedef of the bulk transform. */
119     typedef itk::Transform<ScalarType, InputDimension, OutputDimension> BulkTransformType;
120     typedef BulkTransformType*                  BulkTransformPointer;
121
122     /** Typedef of the mask */
123     typedef itk::SpatialObject<  InputDimension >   MaskType;
124     typedef MaskType*    MaskPointer;
125
126     //====================================================================
127     // Set et Gets
128     //====================================================================
129     //JV  added for the BLUT interpolator
130     void SetSplineOrder(const unsigned int &  splineOrder);
131     void SetSplineOrders(const SizeType &  splineOrders);
132     itkGetMacro( SplineOrders, SizeType );
133     itkGetConstMacro( SplineOrders, SizeType );
134     void SetLUTSamplingFactor(const int &  samplingFactor);
135     void SetLUTSamplingFactors(const SizeType &  samplingFactors);
136     itkGetMacro( LUTSamplingFactors, SizeType );
137     itkGetConstMacro( LUTSamplingFactors,SizeType );
138     //void SetOutputSpacing(const OutputSpacingType & outputSpacing);
139     //itkGetMacro( OutputSpacing, OutputSpacingType );
140     //itkGetConstMacro( OutputSpacing, OutputSpacingType );
141
142     void SetParameters(const ParametersType & parameters);
143   
144     void SetFixedParameters(const ParametersType & parameters);
145
146     void SetParametersByValue(const ParametersType & parameters);
147
148     void SetIdentity();
149
150     /** Get the Transformation Parameters. */
151     virtual const ParametersType& GetParameters(void) const;
152
153     /** Get the Transformation Fixed Parameters. */
154     virtual const ParametersType& GetFixedParameters(void) const;
155   
156     // The coefficientImage
157     virtual CoefficientImagePointer GetCoefficientImage()
158     { return m_CoefficientImage; }
159     virtual const CoefficientImagePointer  GetCoefficientImage() const
160     { return m_CoefficientImage; }
161     virtual void SetCoefficientImage(CoefficientImagePointer image);  
162
163     /** This method specifies the region over which the grid resides. */
164     virtual void SetGridRegion( const RegionType& region );
165     itkGetMacro( GridRegion, RegionType );
166     itkGetConstMacro( GridRegion, RegionType );
167
168     /** This method specifies the grid spacing or resolution. */
169     virtual void SetGridSpacing( const SpacingType& spacing );
170     itkGetMacro( GridSpacing, SpacingType );
171     itkGetConstMacro( GridSpacing, SpacingType );
172
173     /** This method specifies the grid directions . */
174     virtual void SetGridDirection( const DirectionType & spacing );
175     itkGetMacro( GridDirection, DirectionType );
176     itkGetConstMacro( GridDirection, DirectionType );
177
178     /** This method specifies the grid origin. */
179     virtual void SetGridOrigin( const OriginType& origin );
180     itkGetMacro( GridOrigin, OriginType );
181     itkGetConstMacro( GridOrigin, OriginType );
182     
183     // Set the bulk transform, real pointer
184     // itkSetConstObjectMacro( BulkTransform, BulkTransformType );
185     // itkGetConstObjectMacro( BulkTransform, BulkTransformType );
186     void SetBulkTransform(BulkTransformPointer b){m_BulkTransform=b;}
187     BulkTransformPointer GetBulkTransform(void) {return m_BulkTransform;}
188
189     //Set mask, inside transform applies, outside zero, real pointer
190     void SetMask(MaskPointer m){m_Mask=m;}
191     MaskPointer GetMask(void){return m_Mask;}
192     // itkSetConstObjectMacro( Mask, MaskType );
193     // itkGetConstObjectMacro( Mask, MaskType );
194
195     /** Transform points by a BSpline deformable transformation. */
196     OutputPointType  TransformPoint(const InputPointType  &point ) const;
197   
198     // JV added for just the deformable part, without bulk
199     OutputPointType  DeformablyTransformPoint(const InputPointType  &point ) const;
200   
201     /** Parameter index array type. */
202     typedef itk::Array<unsigned long> ParameterIndexArrayType;
203
204     /** Transform points by a BSpline deformable transformation. 
205      * On return, weights contains the interpolation weights used to compute the 
206      * deformation and indices of the x (zeroth) dimension coefficient parameters
207      * in the support region used to compute the deformation.
208      * Parameter indices for the i-th dimension can be obtained by adding
209      * ( i * this->GetNumberOfParametersPerDimension() ) to the indices array.
210      */
211
212     // JV not implemented
213     //  virtual void TransformPoint( const InputPointType & inputPoint,
214     //                           OutputPointType & outputPoint,
215     //                           WeightsType & weights,
216     //                           ParameterIndexArrayType & indices, 
217     //                           bool & inside ) const;
218     //     virtual void DeformablyTransformPoint( const InputPointType & inputPoint,
219     //                                     OutputPointType & outputPoint,
220     //                                     WeightsType & weights,
221     //                                     ParameterIndexArrayType & indices, 
222     //                                     bool & inside ) const;
223     //     virtual void GetJacobian( const InputPointType & inputPoint,
224     //                                WeightsType & weights,
225     //                                ParameterIndexArrayType & indices
226     //                                ) const;
227    
228     /** Method to transform a vector - 
229      *  not applicable for this type of transform. */
230     virtual OutputVectorType TransformVector(const InputVectorType &) const
231     { 
232       itkExceptionMacro(<< "Method not applicable for deformable transform." );
233       return OutputVectorType(); 
234     }
235
236     /** Method to transform a vnl_vector - 
237      *  not applicable for this type of transform */
238     virtual OutputVnlVectorType TransformVector(const InputVnlVectorType &) const
239     { 
240       itkExceptionMacro(<< "Method not applicable for deformable transform. ");
241       return OutputVnlVectorType(); 
242     }
243
244     /** Method to transform a CovariantVector - 
245      *  not applicable for this type of transform */
246     virtual OutputCovariantVectorType TransformCovariantVector(
247                                                                const InputCovariantVectorType &) const
248     { 
249       itkExceptionMacro(<< "Method not applicable for deformable transfrom. ");
250       return OutputCovariantVectorType(); 
251     } 
252     
253     /** Compute the Jacobian Matrix of the transformation at one point */
254     virtual const JacobianType& GetJacobian(const InputPointType  &point ) const;
255
256     /** Return the number of parameters that completely define the Transfom */
257     virtual unsigned int GetNumberOfParameters(void) const;
258
259     /** Return the number of parameters per dimension */
260     unsigned int GetNumberOfParametersPerDimension(void) const;
261
262     /** Return the region of the grid wholly within the support region */
263     itkGetConstReferenceMacro( ValidRegion, RegionType );
264
265     /** Indicates that this transform is linear. That is, given two
266      * points P and Q, and scalar coefficients a and b, then
267      *
268      *           T( a*P + b*Q ) = a * T(P) + b * T(Q)
269      */
270     virtual bool IsLinear() const { return false; }
271
272    //unsigned int GetNumberOfAffectedWeights() const;
273
274   protected:
275     /** Print contents of an BSplineDeformableTransform. */
276     void PrintSelf(std::ostream &os, itk::Indent indent) const;
277
278
279     BSplineDeformableTransform();
280     virtual ~BSplineDeformableTransform();
281
282     /** Wrap flat array into images of coefficients. */
283     void WrapAsImages();
284
285     /** Convert an input point to a continuous index inside the BSpline grid */
286     void TransformPointToContinuousIndex( const InputPointType & point, ContinuousIndexType & index ) const;
287
288   private:
289     BSplineDeformableTransform(const Self&); //purposely not implemented
290     void operator=(const Self&); //purposely not implemented
291
292     /** The bulk transform. */
293     BulkTransformPointer  m_BulkTransform;
294     MaskPointer m_Mask;
295
296     // JV added for BLUT interpolator
297     SizeType m_SplineOrders;
298     SizeType m_LUTSamplingFactors;
299
300     /** Variables defining the coefficient grid extend. */
301     RegionType    m_GridRegion;
302     SpacingType   m_GridSpacing;
303     DirectionType m_GridDirection;
304     OriginType    m_GridOrigin;
305
306     DirectionType m_PointToIndex;
307     DirectionType m_IndexToPoint;
308   
309     RegionType    m_ValidRegion;
310
311     /** Variables defining the interpolation support region. */
312     SizeType m_Offset;
313     itk::FixedArray<bool,InputDimension>    m_SplineOrderOdd;
314     SizeType      m_SupportSize;
315     IndexType     m_ValidRegionLast;
316   
317     /** Array holding images wrapped from the flat parameters. */
318     CoefficientImagePointer   m_WrappedImage;
319
320     /** Vector image representing the B-spline coefficients 
321      *  in each dimension. */
322     CoefficientImagePointer   m_CoefficientImage;
323
324     /** Jacobian as OutputDimension number of images. */
325     typedef typename JacobianType::ValueType JacobianValueType;
326     typedef typename itk::Vector<JacobianValueType,OutputDimension> JacobianPixelType;
327     typedef itk::Image<JacobianPixelType, OutputDimension> JacobianImageType;
328     typename JacobianImageType::Pointer m_JacobianImage[OutputDimension];
329     typedef itk::ImageRegionIterator<JacobianImageType> IteratorType;
330
331     /** Keep track of last support region used in computing the Jacobian
332      * for fast resetting of Jacobian to zero.
333      */
334     //JV for J calculation
335     mutable RegionType    m_SupportRegion;
336     mutable IndexType     m_SupportIndex;
337     mutable IndexType m_LastJacobianIndex;
338     mutable IteratorType m_Iterator[OutputDimension];
339     mutable JacobianPixelType m_ZeroVector;
340     mutable ContinuousIndexType m_Index;
341
342     /** Keep a pointer to the input parameters. */
343     const ParametersType *  m_InputParametersPointer;
344
345     /** Internal parameters buffer. */
346     ParametersType          m_InternalParametersBuffer;
347
348     //JV  the BLUT interpolator
349     typename VectorInterpolatorType::Pointer m_VectorInterpolator;
350   
351     /** Check if a continuous index is inside the valid region. */
352     bool InsideValidRegion( const ContinuousIndexType& index ) const;
353
354
355   }; //class BSplineDeformableTransform
356
357
358 }  // namespace itk
359
360 #if ITK_TEMPLATE_TXX
361 # include "clitkBSplineDeformableTransform.txx"
362 #endif
363
364
365 #endif // __clitkBSplineDeformableTransform_h