]> Creatis software - clitk.git/blob - registration/clitkBSplineDeformableTransform.txx
Debug RTStruct conversion with empty struc
[clitk.git] / registration / clitkBSplineDeformableTransform.txx
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://www.centreleonberard.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_txx
19 #define __clitkBSplineDeformableTransform_txx
20 #include "clitkBSplineDeformableTransform.h"
21
22 //itk
23 #include "itkContinuousIndex.h"
24 #include "itkImageRegionIterator.h"
25 #include "itkImageRegionConstIterator.h"
26 #include "itkImageRegionConstIteratorWithIndex.h"
27 #include "itkIdentityTransform.h"
28
29 namespace clitk
30 {
31
32   // Constructor with default arguments
33   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
34   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>::BSplineDeformableTransform():Superclass(0)
35   {
36     unsigned int i;
37     
38     //JV default spline order
39     for ( i=0;i<InputDimension; i++)
40       m_SplineOrders[i]=3;
41
42     //JV default LUTSamplingfactor
43     for ( i=0;i<InputDimension; i++)
44       m_LUTSamplingFactors[i]=20;
45
46     for (i=0;i<InputDimension; i++)
47       m_SupportSize[i] = m_SplineOrders[i]+1;
48
49     //Instantiate a Vector Bspline Interpolator
50     m_VectorInterpolator =VectorInterpolatorType::New();
51     m_VectorInterpolator->SetLUTSamplingFactors(m_LUTSamplingFactors);
52     m_VectorInterpolator->SetSplineOrders(m_SplineOrders);
53
54     // Set Bulk transform to NULL (no bulk performed)
55     m_BulkTransform = NULL;
56
57     // Mask
58     m_Mask=NULL;
59
60     // Default grid size is zero
61     typename RegionType::SizeType size;
62     typename RegionType::IndexType index;
63     size.Fill( 0 );
64     index.Fill( 0 );
65     m_GridRegion.SetSize( size );
66     m_GridRegion.SetIndex( index );
67
68     m_GridOrigin.Fill( 0.0 );  // default origin is all zeros
69     m_GridSpacing.Fill( 1.0 ); // default spacing is all ones
70     m_GridDirection.SetIdentity(); // default spacing is all ones
71
72     m_InternalParametersBuffer = ParametersType(0);
73     // Make sure the parameters pointer is not NULL after construction.
74     m_InputParametersPointer = &m_InternalParametersBuffer;
75
76     // Initialize coeffient images
77     m_WrappedImage = CoefficientImageType::New();
78     m_WrappedImage->SetRegions( m_GridRegion );
79     m_WrappedImage->SetOrigin( m_GridOrigin.GetDataPointer() );
80     m_WrappedImage->SetSpacing( m_GridSpacing.GetDataPointer() );
81     m_WrappedImage->SetDirection( m_GridDirection );
82     m_CoefficientImage = NULL;
83   
84     // Variables for computing interpolation
85     for (i=0; i <InputDimension;i++)
86       {
87         m_Offset[i] = m_SplineOrders[i] / 2;
88         m_SplineOrderOdd[i] = m_SplineOrders[i] % 2;
89       }
90     m_ValidRegion = m_GridRegion;
91         
92     // Initialize jacobian images 
93     for (i=0; i <OutputDimension;i++)
94       {
95         m_JacobianImage[i] = JacobianImageType::New();
96         m_JacobianImage[i]->SetRegions( m_GridRegion );
97         m_JacobianImage[i]->SetOrigin( m_GridOrigin.GetDataPointer() );
98         m_JacobianImage[i]->SetSpacing( m_GridSpacing.GetDataPointer() );
99         m_JacobianImage[i]->SetDirection( m_GridDirection );
100       }
101     
102     /** Fixed Parameters store the following information:
103      *     Grid Size
104      *     Grid Origin
105      *     Grid Spacing
106      *     Grid Direction */
107     //JV we add the splineOrders, LUTsamplingfactor, m_Mask  and m_BulkTransform
108     /*
109       Spline orders
110       Sampling factors
111       m_Mask
112       m_BulkTransform
113       The size of these is equal to the  NInputDimensions
114     *********************************************************/
115     this->m_FixedParameters.SetSize ( NInputDimensions * (NInputDimensions + 5)+2 );
116     this->m_FixedParameters.Fill ( 0.0 );
117     for ( i=0; i<NInputDimensions; i++)
118       {
119         this->m_FixedParameters[2*NInputDimensions+i] = m_GridSpacing[i];
120       }
121     for (unsigned int di=0; di<NInputDimensions; di++)
122       {
123         for (unsigned int dj=0; dj<NInputDimensions; dj++)
124           {
125             this->m_FixedParameters[3*NInputDimensions+(di*NInputDimensions+dj)] = m_GridDirection[di][dj];
126           }
127       }
128
129     //JV add splineOrders
130     for ( i=0; i<NInputDimensions; i++)
131       {
132         this->m_FixedParameters[ ( (3+NInputDimensions)*NInputDimensions)+i] = (this->GetSplineOrders())[i];
133       }
134     
135     //JV add LUTsamplingFactors
136     for ( i=0; i<NInputDimensions; i++)
137       {
138         this->m_FixedParameters[( (4+NInputDimensions)*NInputDimensions)+i ] = this->GetLUTSamplingFactors()[i];
139       }
140
141     // JV add the mask pointer
142     this->m_FixedParameters[( (5+NInputDimensions)*NInputDimensions)]=(double)((size_t)m_Mask);
143
144     // JV add the bulkTransform pointer
145     this->m_FixedParameters[( (5+NInputDimensions)*NInputDimensions) +1]=(double)((size_t)m_BulkTransform);
146
147
148     // Calculate the PointToIndex matrices
149     DirectionType scale;
150     for( unsigned int i=0; i<OutputDimension; i++)
151       {
152         scale[i][i] = m_GridSpacing[i];
153       }
154
155     m_IndexToPoint = m_GridDirection * scale;
156     m_PointToIndex = m_IndexToPoint.GetInverse();
157
158     // Jacobian
159     m_LastJacobianIndex = m_ValidRegion.GetIndex();
160
161     //JV initialize some variables for jacobian calculation
162     m_SupportRegion.SetSize(m_SupportSize);
163     m_SupportIndex.Fill(0);
164     m_SupportRegion.SetIndex(m_SupportIndex);
165     for (  i = 0; i < OutputDimension; i++ ) 
166       m_ZeroVector[i]=itk::NumericTraits<JacobianValueType>::Zero;
167   }
168     
169
170   // Destructor
171   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
172   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
173   ::~BSplineDeformableTransform()
174   {
175
176   }
177
178
179   // JV set Spline Order
180   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
181   void
182   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
183   ::SetSplineOrder(const unsigned int & splineOrder)
184   {
185     SizeType splineOrders;
186     for (unsigned int i=0;i<InputDimension;i++)splineOrders[i]=splineOrder;
187
188     this->SetSplineOrders(splineOrders);
189   }
190    
191
192   // JV set Spline Orders
193   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
194   void
195   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
196   ::SetSplineOrders(const SizeType & splineOrders)
197   {
198     if (m_SplineOrders!=splineOrders)
199       {
200         m_SplineOrders=splineOrders;
201         
202         //update the interpolation function
203         m_VectorInterpolator->SetSplineOrders(m_SplineOrders);
204         
205         //update the variables for computing interpolation
206         for (unsigned int i=0; i <InputDimension;i++)
207         {
208           m_SupportSize[i] = m_SplineOrders[i]+1;
209           m_Offset[i] = m_SplineOrders[i] / 2;
210           m_SplineOrderOdd[i] = m_SplineOrders[i] % 2;
211         }
212         this->Modified();
213       }
214   }
215
216
217   // JV set sampling factor
218   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
219   void
220   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
221   ::SetLUTSamplingFactor( const int & samplingFactor)
222   {
223     SizeType samplingFactors;
224     for (unsigned int i=0; i<NInputDimensions; i++)
225       samplingFactors[i]=samplingFactor;
226     
227     //update the interpolation function
228     this->SetLUTSamplingFactors(samplingFactors);
229   }
230
231
232   // JV set sampling factors
233   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
234   void
235   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
236   ::SetLUTSamplingFactors( const SizeType & samplingFactors)
237   {
238     if(m_LUTSamplingFactors!=samplingFactors)
239       {
240         for (unsigned int i=0; i<NInputDimensions; i++)
241           m_LUTSamplingFactors[i]=samplingFactors[i];
242         
243         //update the interpolation function
244         m_VectorInterpolator->SetLUTSamplingFactors(m_LUTSamplingFactors);
245         
246         this->Modified();
247       }
248   }
249
250
251   // Get the number of parameters
252   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
253   typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>::NumberOfParametersType
254   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
255   ::GetNumberOfParameters(void) const
256   {
257     // The number of parameters equal OutputDimension * number of
258     // of pixels in the grid region.
259     return ( static_cast<unsigned int>( OutputDimension ) *
260              static_cast<unsigned int>( m_GridRegion.GetNumberOfPixels() ) );
261
262   }
263
264
265   // Get the number of parameters per dimension
266   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
267   unsigned int
268   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
269   ::GetNumberOfParametersPerDimension(void) const
270   {
271     // The number of parameters per dimension equal number of
272     // of pixels in the grid region.
273     return ( static_cast<unsigned int>( m_GridRegion.GetNumberOfPixels() ) );
274
275   }
276
277
278   // Set the grid region
279   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
280   void
281   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
282   ::SetGridRegion( const RegionType & region )
283   {
284     if ( m_GridRegion != region )
285       {
286         m_GridRegion = region;
287
288         // set regions for each coefficient and jacobian image
289         m_WrappedImage->SetRegions( m_GridRegion );
290         for (unsigned int j=0; j <OutputDimension;j++)
291           m_JacobianImage[j]->SetRegions( m_GridRegion );
292
293         // Set the valid region
294         // If the grid spans the interval [start, last].
295         // The valid interval for evaluation is [start+offset, last-offset]
296         // when spline order is even.
297         // The valid interval for evaluation is [start+offset, last-offset)
298         // when spline order is odd.
299         // Where offset = std::floor(spline / 2 ).
300         // Note that the last pixel is not included in the valid region
301         // with odd spline orders.
302         typename RegionType::SizeType size = m_GridRegion.GetSize();
303         typename RegionType::IndexType index = m_GridRegion.GetIndex();
304         for ( unsigned int j = 0; j < NInputDimensions; j++ )
305           {
306             index[j] += 
307               static_cast< typename RegionType::IndexValueType >( m_Offset[j] );
308             size[j] -= 
309               static_cast< typename RegionType::SizeValueType> ( 2 * m_Offset[j] );
310             m_ValidRegionLast[j] = index[j] +
311               static_cast< typename RegionType::IndexValueType >( size[j] ) - 1;
312           }
313         m_ValidRegion.SetSize( size );
314         m_ValidRegion.SetIndex( index );
315
316         // If we are using the default parameters, update their size and set to identity.
317         // Input parameters point to internal buffer => using default parameters.
318         if (m_InputParametersPointer == &m_InternalParametersBuffer)
319           {
320             // Check if we need to resize the default parameter buffer.
321             if ( m_InternalParametersBuffer.GetSize() != this->GetNumberOfParameters() )
322               {
323                 m_InternalParametersBuffer.SetSize( this->GetNumberOfParameters() );
324                 // Fill with zeros for identity.
325                 m_InternalParametersBuffer.Fill( 0 );
326               }
327           }
328
329         this->Modified();
330       }
331   }
332
333
334   // Set the grid spacing
335   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
336   void
337   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
338   ::SetGridSpacing( const SpacingType & spacing )
339   {
340     if ( m_GridSpacing != spacing )
341       {
342         m_GridSpacing = spacing;
343
344         // Set spacing for each coefficient and jacobian image
345         m_WrappedImage->SetSpacing( m_GridSpacing.GetDataPointer() );
346         for (unsigned int j=0; j <OutputDimension; j++) m_JacobianImage[j]->SetSpacing( m_GridSpacing.GetDataPointer() );
347         
348         // Set scale
349         DirectionType scale;
350         for( unsigned int i=0; i<OutputDimension; i++)
351           {
352             scale[i][i] = m_GridSpacing[i];
353           }
354
355         m_IndexToPoint = m_GridDirection * scale;
356         m_PointToIndex = m_IndexToPoint.GetInverse();
357
358         this->Modified();
359       }
360
361   }
362
363   // Set the grid direction
364   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
365   void
366   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
367   ::SetGridDirection( const DirectionType & direction )
368   {
369     if ( m_GridDirection != direction )
370       {
371         m_GridDirection = direction;
372
373         // Set direction for each coefficient and jacobian image
374         m_WrappedImage->SetDirection( m_GridDirection );
375         for (unsigned int j=0; j <OutputDimension; j++) m_JacobianImage[j]->SetDirection( m_GridDirection );
376         
377         // Set scale
378         DirectionType scale;
379         for( unsigned int i=0; i<OutputDimension; i++)
380           {
381             scale[i][i] = m_GridSpacing[i];
382           }
383
384         m_IndexToPoint = m_GridDirection * scale;
385         m_PointToIndex = m_IndexToPoint.GetInverse();
386
387         this->Modified();
388       }
389
390   }
391
392
393   // Set the grid origin
394   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
395   void
396   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
397   ::SetGridOrigin( const OriginType& origin )
398   {
399     if ( m_GridOrigin != origin )
400       {
401         m_GridOrigin = origin;
402
403         // Set origin for each coefficient and jacobianimage
404         m_WrappedImage->SetOrigin( m_GridOrigin.GetDataPointer() );
405         for (unsigned int j=0; j <OutputDimension; j++) m_JacobianImage[j]->SetOrigin( m_GridOrigin.GetDataPointer() );
406         
407         this->Modified();
408       }
409
410   }
411
412
413   // Set the parameters
414   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
415   void
416   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
417   ::SetIdentity()
418   {
419     if( m_InputParametersPointer )
420       {
421         ParametersType * parameters =
422           const_cast<ParametersType *>( m_InputParametersPointer );
423         parameters->Fill( 0.0 );
424         this->Modified();
425       }
426     else 
427       {
428         itkExceptionMacro( << "Input parameters for the spline haven't been set ! "
429                            << "Set them using the SetParameters or SetCoefficientImage method first." );
430       }
431   }
432
433   // Set the parameters
434   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
435   void
436   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
437   ::SetParameters( const ParametersType & parameters )
438   {
439
440     // Check if the number of parameters match the
441     // Expected number of parameters
442     if ( parameters.Size() != this->GetNumberOfParameters() )
443       {
444         itkExceptionMacro(<<"Mismatched between parameters size "
445                           << parameters.size() 
446                           << " and region size "
447                           << m_GridRegion.GetNumberOfPixels() );
448       }
449
450     // Clean up buffered parameters
451     m_InternalParametersBuffer = ParametersType( 0 );
452
453     // Keep a reference to the input parameters
454     m_InputParametersPointer = &parameters;
455
456     // Wrap flat array as images of coefficients
457     this->WrapAsImages();
458
459     //Set input to vector interpolator
460     m_VectorInterpolator->SetInputImage(this->GetCoefficientImage());
461
462     // Modified is always called since we just have a pointer to the
463     // parameters and cannot know if the parameters have changed.
464     this->Modified();
465   }
466
467
468   // Set the Fixed Parameters
469   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
470   void
471   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
472   ::SetFixedParameters( const ParametersType & parameters )
473   {
474  
475     // JV number should be exact, no defaults for spacing
476     if ( parameters.Size() != NInputDimensions * (5 + NInputDimensions)+2 )
477       {
478         itkExceptionMacro(<< "Mismatched between parameters size "
479                           << parameters.size() 
480                           << " and number of fixed parameters "
481                           << NInputDimensions * (5 + NInputDimensions)+2 );
482       }
483     /********************************************************* 
484     Fixed Parameters store the following information:
485         Grid Size
486         Grid Origin
487         Grid Spacing
488         Grid Direction */
489     // JV we add the splineOrders, LUTsamplingfactor, mask pointer and bulktransform pointer
490     /*
491       Spline orders
492       Sampling factors
493       m_Mask
494       m_BulkTransform
495     The size of these is equal to the  NInputDimensions
496     *********************************************************/
497   
498     /** Set the Grid Parameters */
499     SizeType   gridSize;
500     for (unsigned int i=0; i<NInputDimensions; i++)
501       {
502         gridSize[i] = static_cast<int> (parameters[i]);
503       }
504     RegionType bsplineRegion;
505     bsplineRegion.SetSize( gridSize );
506   
507     /** Set the Origin Parameters */
508     OriginType origin;
509     for (unsigned int i=0; i<NInputDimensions; i++)
510       {
511         origin[i] = parameters[NInputDimensions+i];
512       }
513   
514     /** Set the Spacing Parameters */
515     SpacingType spacing;
516     for (unsigned int i=0; i<NInputDimensions; i++)
517       {
518         spacing[i] = parameters[2*NInputDimensions+i];
519       }
520
521     /** Set the Direction Parameters */
522     DirectionType direction;
523     for (unsigned int di=0; di<NInputDimensions; di++)
524       {
525         for (unsigned int dj=0; dj<NInputDimensions; dj++)
526           {
527             direction[di][dj] = parameters[3*NInputDimensions+(di*NInputDimensions+dj)];
528           }
529       }
530
531     //JV add the spline orders
532     SizeType splineOrders;
533     for (unsigned int i=0; i<NInputDimensions; i++)
534       {
535         splineOrders[i]=parameters[(3+NInputDimensions)*NInputDimensions+i];
536       }
537
538     //JV add the LUT sampling factor
539     SizeType  samplingFactors;
540     for (unsigned int i=0; i<NInputDimensions; i++)
541       {
542         samplingFactors[i]=parameters[(4+NInputDimensions)*NInputDimensions+i];
543       }
544
545     //JV add the MaskPointer
546     m_Mask=((MaskPointer)( (size_t)parameters[(5+NInputDimensions)*NInputDimensions]));
547
548     //JV add the MaskPointer
549     m_BulkTransform=((BulkTransformPointer)( (size_t)parameters[(5+NInputDimensions)*NInputDimensions+1]));
550
551     //Set the members
552     this->SetGridSpacing( spacing );
553     this->SetGridDirection( direction );
554     this->SetGridOrigin( origin );
555     this->SetGridRegion( bsplineRegion );
556
557     //JV add the LUT parameters
558     this->SetSplineOrders( splineOrders );
559     this->SetLUTSamplingFactors( samplingFactors );
560       
561   }
562
563
564   // Wrap flat parameters as images
565   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
566   void
567   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
568   ::WrapAsImages()
569   {
570     //JV Wrap parameter array in vectorial image, changed parameter order: A1x A1y A1z, A2x ....
571     PixelType * dataPointer =reinterpret_cast<PixelType *>( const_cast<double *>(m_InputParametersPointer->data_block() )) ;
572     unsigned int numberOfPixels = m_GridRegion.GetNumberOfPixels();
573     
574     m_WrappedImage->GetPixelContainer()->SetImportPointer( dataPointer,numberOfPixels);//InputDimension
575     m_CoefficientImage = m_WrappedImage;
576     
577     //JV Wrap jacobian into OutputDimension X Vectorial images
578     this->m_SharedDataBSplineJacobian.set_size( OutputDimension, this->GetNumberOfParameters() );
579
580     // Use memset to set the memory
581     JacobianPixelType * jacobianDataPointer = reinterpret_cast<JacobianPixelType *>(this->m_SharedDataBSplineJacobian.data_block());
582
583     memset(jacobianDataPointer, 0,  OutputDimension*numberOfPixels*sizeof(JacobianPixelType));
584     m_LastJacobianIndex = m_ValidRegion.GetIndex();
585     m_NeedResetJacobian = false;
586
587     for (unsigned int j=0; j<OutputDimension; j++)
588       {
589         m_JacobianImage[j]->GetPixelContainer()->
590           SetImportPointer( jacobianDataPointer, numberOfPixels );
591         jacobianDataPointer += numberOfPixels;
592      }
593   }
594
595
596   // Set the parameters by value
597   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
598   void
599   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
600   ::SetParametersByValue( const ParametersType & parameters )
601   {
602
603     // Check if the number of parameters match the
604     // Expected number of parameters
605     if ( parameters.Size() != this->GetNumberOfParameters() )
606       {
607         itkExceptionMacro(<<"Mismatched between parameters size "
608                           << parameters.size() 
609                           << " and region size "
610                           << m_GridRegion.GetNumberOfPixels() );
611       }
612
613     // Copy it
614     m_InternalParametersBuffer = parameters;
615     m_InputParametersPointer = &m_InternalParametersBuffer;
616
617     // Wrap flat array as images of coefficients
618     this->WrapAsImages();
619
620     // Modified is always called since we just have a pointer to the
621     // Parameters and cannot know if the parameters have changed.
622     this->Modified();
623
624   }
625
626   // Get the parameters
627   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
628   const 
629   typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
630   ::ParametersType &
631   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
632   ::GetParameters( void ) const
633   {
634     /** NOTE: For efficiency, this class does not keep a copy of the parameters - 
635      * it just keeps pointer to input parameters. 
636      */
637     if (NULL == m_InputParametersPointer)
638       {
639         itkExceptionMacro( <<"Cannot GetParameters() because m_InputParametersPointer is NULL. Perhaps SetCoefficientImage() has been called causing the NULL pointer." );
640       }
641
642     return (*m_InputParametersPointer);
643   }
644
645
646   // Get the parameters
647   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
648   const 
649   typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
650   ::ParametersType &
651   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
652   ::GetFixedParameters( void ) const
653   {
654     RegionType resRegion = this->GetGridRegion(  );
655   
656     for (unsigned int i=0; i<NInputDimensions; i++)
657       {
658         this->m_FixedParameters[i] = (resRegion.GetSize())[i];
659       }
660     for (unsigned int i=0; i<NInputDimensions; i++)
661       {
662         this->m_FixedParameters[NInputDimensions+i] = (this->GetGridOrigin())[i];
663       } 
664     for (unsigned int i=0; i<NInputDimensions; i++)
665       {
666         this->m_FixedParameters[2*NInputDimensions+i] =  (this->GetGridSpacing())[i];
667       }
668     for (unsigned int di=0; di<NInputDimensions; di++)
669       {
670         for (unsigned int dj=0; dj<NInputDimensions; dj++)
671           {
672             this->m_FixedParameters[3*NInputDimensions+(di*NInputDimensions+dj)] = (this->GetGridDirection())[di][dj];
673           }
674       }
675
676     //JV add splineOrders
677     for (unsigned int i=0; i<NInputDimensions; i++)
678       {
679         this->m_FixedParameters[(3+NInputDimensions)*NInputDimensions+i] = (this->GetSplineOrders())[i];
680       }
681     
682     //JV add LUTsamplingFactor
683     for (unsigned int i=0; i<NInputDimensions; i++)
684       {
685         this->m_FixedParameters[(4+NInputDimensions)*NInputDimensions+i] = (this->GetLUTSamplingFactors())[i];
686       }
687     
688     //JV add the mask
689     this->m_FixedParameters[(5+NInputDimensions)*NInputDimensions]=(double)((size_t) m_Mask);
690
691     //JV add the bulktransform pointer
692     this->m_FixedParameters[(5+NInputDimensions)*NInputDimensions+1]=(double)((size_t) m_BulkTransform);
693     
694     return (this->m_FixedParameters);
695   }
696
697   
698   // Set the B-Spline coefficients using input images
699   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
700   void 
701   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
702   ::SetCoefficientImage( CoefficientImagePointer image )
703   {
704     
705     this->SetGridRegion( image->GetBufferedRegion() );
706     this->SetGridSpacing( image->GetSpacing() );
707     this->SetGridDirection( image->GetDirection() );
708     this->SetGridOrigin( image->GetOrigin() );
709     m_CoefficientImage = image;
710
711     // Update the interpolator    
712     m_VectorInterpolator->SetInputImage(this->GetCoefficientImage());
713
714     // Clean up buffered parameters
715     m_InternalParametersBuffer = ParametersType( 0 );
716     m_InputParametersPointer  = NULL;
717       
718   }  
719
720
721   // Print self
722   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
723   void
724   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
725   ::PrintSelf(std::ostream &os, itk::Indent indent) const
726   {
727
728     this->Superclass::PrintSelf(os, indent);
729
730     os << indent << "GridRegion: " << m_GridRegion << std::endl;
731     os << indent << "GridOrigin: " << m_GridOrigin << std::endl;
732     os << indent << "GridSpacing: " << m_GridSpacing << std::endl;
733     os << indent << "GridDirection: " << m_GridDirection << std::endl;
734     os << indent << "IndexToPoint: " << m_IndexToPoint << std::endl;
735     os << indent << "PointToIndex: " << m_PointToIndex << std::endl;
736
737     os << indent << "CoefficientImage: [ ";
738     os << m_CoefficientImage.GetPointer() << " ]" << std::endl;
739
740     os << indent << "WrappedImage: [ ";
741     os << m_WrappedImage.GetPointer() << " ]" << std::endl;
742  
743     os << indent << "InputParametersPointer: " 
744        << m_InputParametersPointer << std::endl;
745     os << indent << "ValidRegion: " << m_ValidRegion << std::endl;
746     os << indent << "LastJacobianIndex: " << m_LastJacobianIndex << std::endl;
747     os << indent << "BulkTransform: ";
748     os << m_BulkTransform << std::endl;
749
750     if ( m_BulkTransform )
751       {
752         os << indent << "BulkTransformType: " 
753            << m_BulkTransform->GetNameOfClass() << std::endl;
754       }
755     os << indent << "VectorBSplineInterpolator: ";
756     os << m_VectorInterpolator.GetPointer() << std::endl;
757     os << indent << "Mask: ";
758     os << m_Mask<< std::endl;
759   }
760
761
762   // Verify 
763   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
764   bool 
765   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
766   ::InsideValidRegion( const ContinuousIndexType& index ) const
767   {
768     bool inside = true;
769
770     if ( !m_ValidRegion.IsInside( index ) )
771       {
772         inside = false;
773       }
774     //JV verify for each input dimension
775     if ( inside)
776       {
777         typedef typename ContinuousIndexType::ValueType ValueType;
778         for( unsigned int j = 0; j < InputDimension; j++ )
779           {
780             if (m_SplineOrderOdd[j])
781               {
782                 if ( index[j] >= static_cast<ValueType>( m_ValidRegionLast[j] ) )
783                   { 
784                     inside = false;
785                     break;
786                   }
787               }
788           }
789       }
790     
791     return inside;
792   }
793
794
795   // Transform a point
796   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
797   typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
798   ::OutputPointType
799   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
800   ::TransformPoint(const InputPointType &inputPoint) const 
801   {
802   
803    InputPointType transformedPoint;
804    OutputPointType outputPoint;
805
806    // BulkTransform
807    if ( m_BulkTransform )
808      {
809        transformedPoint = m_BulkTransform->TransformPoint( inputPoint );
810      }
811    else
812       {
813         transformedPoint = inputPoint;
814       }
815    
816     // Deformable transform
817     if ( m_CoefficientImage )
818       {
819         // Check if inside mask
820         if(m_Mask &&  !(m_Mask->IsInside(inputPoint) ) )
821           {
822             // Outside: no (deformable) displacement
823             return transformedPoint;
824           }             
825         
826         // Check if inside valid region
827         bool inside = true;
828         ContinuousIndexType index;
829         this->TransformPointToContinuousIndex( inputPoint, index );
830         inside = this->InsideValidRegion( index );
831         if ( !inside )
832           {
833             // Outside: no (deformable) displacement
834             return transformedPoint;
835           }
836         
837         // Call the vector interpolator
838         itk::Vector<TCoordRep,OutputDimension> displacement=m_VectorInterpolator->EvaluateAtContinuousIndex(index);
839
840         // Return the results
841         outputPoint = transformedPoint+displacement;
842
843       }
844
845     else
846       {
847         itkWarningMacro( << "B-spline coefficients have not been set" );
848         outputPoint = transformedPoint;
849       }
850     
851     return outputPoint;
852   }
853
854
855
856   //JV Deformably transform a point
857   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
858   typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
859   ::OutputPointType
860   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
861   ::DeformablyTransformPoint(const InputPointType &inputPoint) const 
862   {
863     OutputPointType outputPoint;
864     if ( m_CoefficientImage )
865       {
866
867         // Check if inside mask
868         if(m_Mask && !(m_Mask->IsInside(inputPoint) ) )
869           {
870             // Outside: no (deformable) displacement
871             return inputPoint;
872           }     
873         
874         // Check if inside
875         bool inside = true;
876         ContinuousIndexType index;
877         this->TransformPointToContinuousIndex( inputPoint, index );
878         inside = this->InsideValidRegion( index );
879
880         if ( !inside )
881           {
882             //outside: no (deformable) displacement
883              return inputPoint;
884              //return outputPoint;
885           }
886
887         // Call the vector interpolator
888         itk::Vector<TCoordRep,OutputDimension> displacement=m_VectorInterpolator->EvaluateAtContinuousIndex(index);
889         
890         // Return the results
891         outputPoint = inputPoint+displacement;
892       }
893
894     // No coefficients available
895     else
896       {
897         itkWarningMacro( << "B-spline coefficients have not been set" );
898         outputPoint = inputPoint;
899       }
900    
901     return outputPoint;
902   }
903   
904
905   // JV weights are identical as for transformpoint, could be done simultaneously in metric!!!!
906   // Compute the Jacobian in one position 
907   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
908   void
909   BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
910   ::ComputeJacobianWithRespectToParameters(const InputPointType & point, JacobianType & jacobian) const
911   {
912     if (m_NeedResetJacobian)
913       ResetJacobian();
914
915     //========================================================
916     // For each dimension, copy the weight to the support region
917     //========================================================
918
919     // Check if inside mask
920     if (m_Mask &&  !(m_Mask->IsInside(point) ) )
921       {
922         // Outside: no (deformable) displacement
923         jacobian = m_SharedDataBSplineJacobian;
924         return;
925       }
926
927     //Get index
928     this->TransformPointToContinuousIndex( point, m_Index );
929
930     // NOTE: if the support region does not lie totally within the grid
931     // we assume zero displacement and return the input point
932     if ( !this->InsideValidRegion( m_Index ) )
933       {
934         jacobian = m_SharedDataBSplineJacobian;
935         return;
936       }
937
938     //Compute interpolation weights
939     const WeightsDataType *weights = NULL;
940     m_VectorInterpolator->EvaluateWeightsAtContinuousIndex( m_Index, &weights, m_LastJacobianIndex);
941     m_SupportRegion.SetIndex( m_LastJacobianIndex );
942
943     //Reset the iterators
944     unsigned int j = 0;
945     for ( j = 0; j < OutputDimension; j++ )
946       m_Iterator[j] = IteratorType( m_JacobianImage[j], m_SupportRegion);
947
948     // For each dimension, copy the weight to the support region
949     while ( ! (m_Iterator[0]).IsAtEnd() )
950       {
951         //copy weight to jacobian image
952         for ( j = 0; j < OutputDimension; j++ )
953           {
954             m_ZeroVector[j]=*weights;
955             (m_Iterator[j]).Set( m_ZeroVector);
956             m_ZeroVector[j]=itk::NumericTraits<JacobianValueType>::Zero;
957             ++(m_Iterator[j]);
958           }
959         // go to next coefficient in the support region
960         weights++;
961       }
962     m_NeedResetJacobian = true;
963
964     // Return the result
965     jacobian = m_SharedDataBSplineJacobian;
966
967   }
968
969
970
971   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
972   void
973   BSplineDeformableTransform<TCoordRep, NInputDimensions,NOutputDimensions>
974   ::TransformPointToContinuousIndex( const InputPointType & point, ContinuousIndexType & index ) const
975   {
976     unsigned int j;
977
978     itk::Vector<double, OutputDimension> tvector;
979
980     for ( j = 0; j < OutputDimension; j++ )
981       {
982         tvector[j] = point[j] - this->m_GridOrigin[j];
983       }
984
985     itk::Vector<double, OutputDimension> cvector;
986
987     cvector = m_PointToIndex * tvector;
988
989     for ( j = 0; j < OutputDimension; j++ )
990       {
991         index[j] = static_cast< typename ContinuousIndexType::CoordRepType >( cvector[j] );
992       }
993   }
994
995   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
996   unsigned
997   BSplineDeformableTransform<TCoordRep, NInputDimensions,NOutputDimensions>
998   ::SetJacobianImageData(JacobianPixelType * jacobianDataPointer, unsigned dim)
999   {
1000     unsigned int numberOfPixels = m_GridRegion.GetNumberOfPixels();
1001     m_JacobianImage[dim]->GetPixelContainer()->SetImportPointer(jacobianDataPointer, numberOfPixels);
1002     return numberOfPixels;
1003   }
1004
1005   template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
1006   void
1007   BSplineDeformableTransform<TCoordRep, NInputDimensions,NOutputDimensions>
1008   ::ResetJacobian() const
1009   {
1010     //========================================================
1011     // Zero all components of jacobian
1012     //========================================================
1013     // JV not thread safe (m_LastJacobianIndex), instantiate N transforms
1014     // NOTE: for efficiency, we only need to zero out the coefficients
1015     // that got fill last time this method was called.
1016
1017     unsigned int j = 0;
1018
1019     //Define the  region for each jacobian image
1020     m_SupportRegion.SetIndex(m_LastJacobianIndex);
1021
1022     //Initialize the iterators
1023     for (j = 0; j < OutputDimension; j++)
1024       m_Iterator[j] = IteratorType(m_JacobianImage[j], m_SupportRegion);
1025
1026     //Set the previously-set to zero
1027     while (!(m_Iterator[0]).IsAtEnd())
1028       {
1029         for (j = 0; j < OutputDimension; j++)
1030           {
1031             m_Iterator[j].Set(m_ZeroVector);
1032             ++(m_Iterator[j]);
1033           }
1034       }
1035     m_NeedResetJacobian = false;
1036   }
1037 } // namespace
1038
1039 #endif