1 /*=========================================================================
2 Program: vv http://www.creatis.insa-lyon.fr/rio/vv
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
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.
13 It is distributed under dual licence
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_tx
20 #include "clitkBSplineDeformableTransform.h"
23 #include "itkContinuousIndex.h"
24 #include "itkImageRegionIterator.h"
25 #include "itkImageRegionConstIterator.h"
26 #include "itkImageRegionConstIteratorWithIndex.h"
27 #include "itkIdentityTransform.h"
32 // Constructor with default arguments
33 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
34 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
35 ::BSplineDeformableTransform():Superclass(OutputDimension,0)
39 //JV default spline order
40 for ( i=0;i<InputDimension; i++)
43 //JV default LUTSamplingfactor
44 for ( i=0;i<InputDimension; i++)
45 m_LUTSamplingFactors[i]=20;
47 for (i=0;i<InputDimension; i++)
48 m_SupportSize[i] = m_SplineOrders[i]+1;
50 //Instantiate a Vector Bspline Interpolator
51 m_VectorInterpolator =VectorInterpolatorType::New();
52 m_VectorInterpolator->SetLUTSamplingFactors(m_LUTSamplingFactors);
53 m_VectorInterpolator->SetSplineOrders(m_SplineOrders);
55 // Set Bulk transform to NULL (no bulk performed)
56 m_BulkTransform = NULL;
61 // Default grid size is zero
62 typename RegionType::SizeType size;
63 typename RegionType::IndexType index;
66 m_GridRegion.SetSize( size );
67 m_GridRegion.SetIndex( index );
69 m_GridOrigin.Fill( 0.0 ); // default origin is all zeros
70 m_GridSpacing.Fill( 1.0 ); // default spacing is all ones
71 m_GridDirection.SetIdentity(); // default spacing is all ones
73 m_InternalParametersBuffer = ParametersType(0);
74 // Make sure the parameters pointer is not NULL after construction.
75 m_InputParametersPointer = &m_InternalParametersBuffer;
77 // Initialize coeffient images
78 m_WrappedImage = CoefficientImageType::New();
79 m_WrappedImage->SetRegions( m_GridRegion );
80 m_WrappedImage->SetOrigin( m_GridOrigin.GetDataPointer() );
81 m_WrappedImage->SetSpacing( m_GridSpacing.GetDataPointer() );
82 m_WrappedImage->SetDirection( m_GridDirection );
83 m_CoefficientImage = NULL;
85 // Variables for computing interpolation
86 for (i=0; i <InputDimension;i++)
88 m_Offset[i] = m_SplineOrders[i] / 2;
89 m_SplineOrderOdd[i] = m_SplineOrders[i] % 2;
91 m_ValidRegion = m_GridRegion;
93 // Initialize jacobian images
94 for (i=0; i <OutputDimension;i++)
96 m_JacobianImage[i] = JacobianImageType::New();
97 m_JacobianImage[i]->SetRegions( m_GridRegion );
98 m_JacobianImage[i]->SetOrigin( m_GridOrigin.GetDataPointer() );
99 m_JacobianImage[i]->SetSpacing( m_GridSpacing.GetDataPointer() );
100 m_JacobianImage[i]->SetDirection( m_GridDirection );
103 /** Fixed Parameters store the following information:
108 //JV we add the splineOrders, LUTsamplingfactor, m_Mask and m_BulkTransform
114 The size of these is equal to the NInputDimensions
115 *********************************************************/
116 this->m_FixedParameters.SetSize ( NInputDimensions * (NInputDimensions + 5)+2 );
117 this->m_FixedParameters.Fill ( 0.0 );
118 for ( i=0; i<NInputDimensions; i++)
120 this->m_FixedParameters[2*NInputDimensions+i] = m_GridSpacing[i];
122 for (unsigned int di=0; di<NInputDimensions; di++)
124 for (unsigned int dj=0; dj<NInputDimensions; dj++)
126 this->m_FixedParameters[3*NInputDimensions+(di*NInputDimensions+dj)] = m_GridDirection[di][dj];
130 //JV add splineOrders
131 for ( i=0; i<NInputDimensions; i++)
133 this->m_FixedParameters[ ( (3+NInputDimensions)*NInputDimensions)+i] = (this->GetSplineOrders())[i];
136 //JV add LUTsamplingFactors
137 for ( i=0; i<NInputDimensions; i++)
139 this->m_FixedParameters[( (4+NInputDimensions)*NInputDimensions)+i ] = this->GetLUTSamplingFactors()[i];
142 // JV add the mask pointer
143 this->m_FixedParameters[( (5+NInputDimensions)*NInputDimensions)]=(double)((size_t)m_Mask);
145 // JV add the bulkTransform pointer
146 this->m_FixedParameters[( (5+NInputDimensions)*NInputDimensions) +1]=(double)((size_t)m_BulkTransform);
149 // Calculate the PointToIndex matrices
151 for( unsigned int i=0; i<OutputDimension; i++)
153 scale[i][i] = m_GridSpacing[i];
156 m_IndexToPoint = m_GridDirection * scale;
157 m_PointToIndex = m_IndexToPoint.GetInverse();
160 m_LastJacobianIndex = m_ValidRegion.GetIndex();
162 //JV initialize some variables for jacobian calculation
163 m_SupportRegion.SetSize(m_SupportSize);
164 m_SupportIndex.Fill(0);
165 m_SupportRegion.SetIndex(m_SupportIndex);
166 for ( i = 0; i < OutputDimension; i++ )
167 m_ZeroVector[i]=itk::NumericTraits<JacobianValueType>::Zero;
172 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
173 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
174 ::~BSplineDeformableTransform()
180 // JV set Spline Order
181 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
183 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
184 ::SetSplineOrder(const unsigned int & splineOrder)
186 SizeType splineOrders;
187 for (unsigned int i=0;i<InputDimension;i++)splineOrders[i]=splineOrder;
189 this->SetSplineOrders(splineOrders);
193 // JV set Spline Orders
194 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
196 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
197 ::SetSplineOrders(const SizeType & splineOrders)
199 if (m_SplineOrders!=splineOrders)
201 m_SplineOrders=splineOrders;
203 //update the interpolation function
204 m_VectorInterpolator->SetSplineOrders(m_SplineOrders);
206 //update the variables for computing interpolation
207 for (unsigned int i=0; i <InputDimension;i++)
209 m_SupportSize[i] = m_SplineOrders[i]+1;
210 m_Offset[i] = m_SplineOrders[i] / 2;
211 m_SplineOrderOdd[i] = m_SplineOrders[i] % 2;
218 // JV set sampling factor
219 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
221 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
222 ::SetLUTSamplingFactor( const int & samplingFactor)
224 SizeType samplingFactors;
225 for (unsigned int i=0; i<NInputDimensions; i++)
226 samplingFactors[i]=samplingFactor;
228 //update the interpolation function
229 this->SetLUTSamplingFactors(samplingFactors);
233 // JV set sampling factors
234 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
236 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
237 ::SetLUTSamplingFactors( const SizeType & samplingFactors)
239 if(m_LUTSamplingFactors!=samplingFactors)
241 for (unsigned int i=0; i<NInputDimensions; i++)
242 m_LUTSamplingFactors[i]=samplingFactors[i];
244 //update the interpolation function
245 m_VectorInterpolator->SetLUTSamplingFactors(m_LUTSamplingFactors);
252 // Get the number of parameters
253 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
255 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
256 ::GetNumberOfParameters(void) const
258 // The number of parameters equal OutputDimension * number of
259 // of pixels in the grid region.
260 return ( static_cast<unsigned int>( OutputDimension ) *
261 static_cast<unsigned int>( m_GridRegion.GetNumberOfPixels() ) );
266 // Get the number of parameters per dimension
267 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
269 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
270 ::GetNumberOfParametersPerDimension(void) const
272 // The number of parameters per dimension equal number of
273 // of pixels in the grid region.
274 return ( static_cast<unsigned int>( m_GridRegion.GetNumberOfPixels() ) );
279 // Set the grid region
280 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
282 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
283 ::SetGridRegion( const RegionType & region )
285 if ( m_GridRegion != region )
287 m_GridRegion = region;
289 // set regions for each coefficient and jacobian image
290 m_WrappedImage->SetRegions( m_GridRegion );
291 for (unsigned int j=0; j <OutputDimension;j++)
292 m_JacobianImage[j]->SetRegions( m_GridRegion );
294 // Set the valid region
295 // If the grid spans the interval [start, last].
296 // The valid interval for evaluation is [start+offset, last-offset]
297 // when spline order is even.
298 // The valid interval for evaluation is [start+offset, last-offset)
299 // when spline order is odd.
300 // Where offset = vcl_floor(spline / 2 ).
301 // Note that the last pixel is not included in the valid region
302 // with odd spline orders.
303 typename RegionType::SizeType size = m_GridRegion.GetSize();
304 typename RegionType::IndexType index = m_GridRegion.GetIndex();
305 for ( unsigned int j = 0; j < NInputDimensions; j++ )
308 static_cast< typename RegionType::IndexValueType >( m_Offset[j] );
310 static_cast< typename RegionType::SizeValueType> ( 2 * m_Offset[j] );
311 m_ValidRegionLast[j] = index[j] +
312 static_cast< typename RegionType::IndexValueType >( size[j] ) - 1;
314 m_ValidRegion.SetSize( size );
315 m_ValidRegion.SetIndex( index );
317 // If we are using the default parameters, update their size and set to identity.
318 // Input parameters point to internal buffer => using default parameters.
319 if (m_InputParametersPointer == &m_InternalParametersBuffer)
321 // Check if we need to resize the default parameter buffer.
322 if ( m_InternalParametersBuffer.GetSize() != this->GetNumberOfParameters() )
324 m_InternalParametersBuffer.SetSize( this->GetNumberOfParameters() );
325 // Fill with zeros for identity.
326 m_InternalParametersBuffer.Fill( 0 );
335 // Set the grid spacing
336 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
338 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
339 ::SetGridSpacing( const SpacingType & spacing )
341 if ( m_GridSpacing != spacing )
343 m_GridSpacing = spacing;
345 // Set spacing for each coefficient and jacobian image
346 m_WrappedImage->SetSpacing( m_GridSpacing.GetDataPointer() );
347 for (unsigned int j=0; j <OutputDimension; j++) m_JacobianImage[j]->SetSpacing( m_GridSpacing.GetDataPointer() );
351 for( unsigned int i=0; i<OutputDimension; i++)
353 scale[i][i] = m_GridSpacing[i];
356 m_IndexToPoint = m_GridDirection * scale;
357 m_PointToIndex = m_IndexToPoint.GetInverse();
364 // Set the grid direction
365 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
367 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
368 ::SetGridDirection( const DirectionType & direction )
370 if ( m_GridDirection != direction )
372 m_GridDirection = direction;
374 // Set direction for each coefficient and jacobian image
375 m_WrappedImage->SetDirection( m_GridDirection );
376 for (unsigned int j=0; j <OutputDimension; j++) m_JacobianImage[j]->SetDirection( m_GridDirection );
380 for( unsigned int i=0; i<OutputDimension; i++)
382 scale[i][i] = m_GridSpacing[i];
385 m_IndexToPoint = m_GridDirection * scale;
386 m_PointToIndex = m_IndexToPoint.GetInverse();
394 // Set the grid origin
395 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
397 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
398 ::SetGridOrigin( const OriginType& origin )
400 if ( m_GridOrigin != origin )
402 m_GridOrigin = origin;
404 // Set origin for each coefficient and jacobianimage
405 m_WrappedImage->SetOrigin( m_GridOrigin.GetDataPointer() );
406 for (unsigned int j=0; j <OutputDimension; j++) m_JacobianImage[j]->SetOrigin( m_GridOrigin.GetDataPointer() );
414 // Set the parameters
415 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
417 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
420 if( m_InputParametersPointer )
422 ParametersType * parameters =
423 const_cast<ParametersType *>( m_InputParametersPointer );
424 parameters->Fill( 0.0 );
429 itkExceptionMacro( << "Input parameters for the spline haven't been set ! "
430 << "Set them using the SetParameters or SetCoefficientImage method first." );
434 // Set the parameters
435 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
437 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
438 ::SetParameters( const ParametersType & parameters )
441 // Check if the number of parameters match the
442 // Expected number of parameters
443 if ( parameters.Size() != this->GetNumberOfParameters() )
445 itkExceptionMacro(<<"Mismatched between parameters size "
447 << " and region size "
448 << m_GridRegion.GetNumberOfPixels() );
451 // Clean up buffered parameters
452 m_InternalParametersBuffer = ParametersType( 0 );
454 // Keep a reference to the input parameters
455 m_InputParametersPointer = ¶meters;
457 // Wrap flat array as images of coefficients
458 this->WrapAsImages();
460 //Set input to vector interpolator
461 m_VectorInterpolator->SetInputImage(this->GetCoefficientImage());
463 // Modified is always called since we just have a pointer to the
464 // parameters and cannot know if the parameters have changed.
469 // Set the Fixed Parameters
470 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
472 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
473 ::SetFixedParameters( const ParametersType & parameters )
476 // JV number should be exact, no defaults for spacing
477 if ( parameters.Size() != NInputDimensions * (5 + NInputDimensions)+2 )
479 itkExceptionMacro(<< "Mismatched between parameters size "
481 << " and number of fixed parameters "
482 << NInputDimensions * (5 + NInputDimensions)+2 );
484 /*********************************************************
485 Fixed Parameters store the following information:
490 // JV we add the splineOrders, LUTsamplingfactor, mask pointer and bulktransform pointer
496 The size of these is equal to the NInputDimensions
497 *********************************************************/
499 /** Set the Grid Parameters */
501 for (unsigned int i=0; i<NInputDimensions; i++)
503 gridSize[i] = static_cast<int> (parameters[i]);
505 RegionType bsplineRegion;
506 bsplineRegion.SetSize( gridSize );
508 /** Set the Origin Parameters */
510 for (unsigned int i=0; i<NInputDimensions; i++)
512 origin[i] = parameters[NInputDimensions+i];
515 /** Set the Spacing Parameters */
517 for (unsigned int i=0; i<NInputDimensions; i++)
519 spacing[i] = parameters[2*NInputDimensions+i];
522 /** Set the Direction Parameters */
523 DirectionType direction;
524 for (unsigned int di=0; di<NInputDimensions; di++)
526 for (unsigned int dj=0; dj<NInputDimensions; dj++)
528 direction[di][dj] = parameters[3*NInputDimensions+(di*NInputDimensions+dj)];
532 //JV add the spline orders
533 SizeType splineOrders;
534 for (unsigned int i=0; i<NInputDimensions; i++)
536 splineOrders[i]=parameters[(3+NInputDimensions)*NInputDimensions+i];
539 //JV add the LUT sampling factor
540 SizeType samplingFactors;
541 for (unsigned int i=0; i<NInputDimensions; i++)
543 samplingFactors[i]=parameters[(4+NInputDimensions)*NInputDimensions+i];
546 //JV add the MaskPointer
547 m_Mask=((MaskPointer)( (size_t)parameters[(5+NInputDimensions)*NInputDimensions]));
549 //JV add the MaskPointer
550 m_BulkTransform=((BulkTransformPointer)( (size_t)parameters[(5+NInputDimensions)*NInputDimensions+1]));
553 this->SetGridSpacing( spacing );
554 this->SetGridDirection( direction );
555 this->SetGridOrigin( origin );
556 this->SetGridRegion( bsplineRegion );
558 //JV add the LUT parameters
559 this->SetSplineOrders( splineOrders );
560 this->SetLUTSamplingFactors( samplingFactors );
565 // Wrap flat parameters as images
566 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
568 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
571 //JV Wrap parameter array in vectorial image, changed parameter order: A1x A1y A1z, A2x ....
572 PixelType * dataPointer =reinterpret_cast<PixelType *>( const_cast<double *>(m_InputParametersPointer->data_block() )) ;
573 unsigned int numberOfPixels = m_GridRegion.GetNumberOfPixels();
575 m_WrappedImage->GetPixelContainer()->SetImportPointer( dataPointer,numberOfPixels);//InputDimension
576 m_CoefficientImage = m_WrappedImage;
578 //JV Wrap jacobian into OutputDimension X Vectorial images
579 this->m_Jacobian.set_size( OutputDimension, this->GetNumberOfParameters() );
581 // Use memset to set the memory
582 JacobianPixelType * jacobianDataPointer = reinterpret_cast<JacobianPixelType *>(this->m_Jacobian.data_block());
583 memset(jacobianDataPointer, 0, OutputDimension*numberOfPixels*sizeof(JacobianPixelType));
584 m_LastJacobianIndex = m_ValidRegion.GetIndex();
585 m_NeedResetJacobian = false;
587 for (unsigned int j=0; j<OutputDimension; j++)
589 m_JacobianImage[j]->GetPixelContainer()->
590 SetImportPointer( jacobianDataPointer, numberOfPixels );
591 jacobianDataPointer += numberOfPixels;
596 // Set the parameters by value
597 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
599 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
600 ::SetParametersByValue( const ParametersType & parameters )
603 // Check if the number of parameters match the
604 // Expected number of parameters
605 if ( parameters.Size() != this->GetNumberOfParameters() )
607 itkExceptionMacro(<<"Mismatched between parameters size "
609 << " and region size "
610 << m_GridRegion.GetNumberOfPixels() );
614 m_InternalParametersBuffer = parameters;
615 m_InputParametersPointer = &m_InternalParametersBuffer;
617 // Wrap flat array as images of coefficients
618 this->WrapAsImages();
620 // Modified is always called since we just have a pointer to the
621 // Parameters and cannot know if the parameters have changed.
626 // Get the parameters
627 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
629 typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
631 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
632 ::GetParameters( void ) const
634 /** NOTE: For efficiency, this class does not keep a copy of the parameters -
635 * it just keeps pointer to input parameters.
637 if (NULL == m_InputParametersPointer)
639 itkExceptionMacro( <<"Cannot GetParameters() because m_InputParametersPointer is NULL. Perhaps SetCoefficientImage() has been called causing the NULL pointer." );
642 return (*m_InputParametersPointer);
646 // Get the parameters
647 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
649 typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
651 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
652 ::GetFixedParameters( void ) const
654 RegionType resRegion = this->GetGridRegion( );
656 for (unsigned int i=0; i<NInputDimensions; i++)
658 this->m_FixedParameters[i] = (resRegion.GetSize())[i];
660 for (unsigned int i=0; i<NInputDimensions; i++)
662 this->m_FixedParameters[NInputDimensions+i] = (this->GetGridOrigin())[i];
664 for (unsigned int i=0; i<NInputDimensions; i++)
666 this->m_FixedParameters[2*NInputDimensions+i] = (this->GetGridSpacing())[i];
668 for (unsigned int di=0; di<NInputDimensions; di++)
670 for (unsigned int dj=0; dj<NInputDimensions; dj++)
672 this->m_FixedParameters[3*NInputDimensions+(di*NInputDimensions+dj)] = (this->GetGridDirection())[di][dj];
676 //JV add splineOrders
677 for (unsigned int i=0; i<NInputDimensions; i++)
679 this->m_FixedParameters[(3+NInputDimensions)*NInputDimensions+i] = (this->GetSplineOrders())[i];
682 //JV add LUTsamplingFactor
683 for (unsigned int i=0; i<NInputDimensions; i++)
685 this->m_FixedParameters[(4+NInputDimensions)*NInputDimensions+i] = (this->GetLUTSamplingFactors())[i];
689 this->m_FixedParameters[(5+NInputDimensions)*NInputDimensions]=(double)((size_t) m_Mask);
691 //JV add the bulktransform pointer
692 this->m_FixedParameters[(5+NInputDimensions)*NInputDimensions+1]=(double)((size_t) m_BulkTransform);
694 return (this->m_FixedParameters);
698 // Set the B-Spline coefficients using input images
699 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
701 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
702 ::SetCoefficientImage( CoefficientImagePointer image )
705 this->SetGridRegion( image->GetBufferedRegion() );
706 this->SetGridSpacing( image->GetSpacing() );
707 this->SetGridDirection( image->GetDirection() );
708 this->SetGridOrigin( image->GetOrigin() );
709 m_CoefficientImage = image;
711 // Update the interpolator
712 m_VectorInterpolator->SetInputImage(this->GetCoefficientImage());
714 // Clean up buffered parameters
715 m_InternalParametersBuffer = ParametersType( 0 );
716 m_InputParametersPointer = NULL;
722 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
724 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
725 ::PrintSelf(std::ostream &os, itk::Indent indent) const
728 this->Superclass::PrintSelf(os, indent);
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;
737 os << indent << "CoefficientImage: [ ";
738 os << m_CoefficientImage.GetPointer() << " ]" << std::endl;
740 os << indent << "WrappedImage: [ ";
741 os << m_WrappedImage.GetPointer() << " ]" << std::endl;
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;
750 if ( m_BulkTransform )
752 os << indent << "BulkTransformType: "
753 << m_BulkTransform->GetNameOfClass() << std::endl;
755 os << indent << "VectorBSplineInterpolator: ";
756 os << m_VectorInterpolator.GetPointer() << std::endl;
757 os << indent << "Mask: ";
758 os << m_Mask<< std::endl;
763 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
765 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
766 ::InsideValidRegion( const ContinuousIndexType& index ) const
770 if ( !m_ValidRegion.IsInside( index ) )
774 //JV verify for each input dimension
777 typedef typename ContinuousIndexType::ValueType ValueType;
778 for( unsigned int j = 0; j < InputDimension; j++ )
780 if (m_SplineOrderOdd[j])
782 if ( index[j] >= static_cast<ValueType>( m_ValidRegionLast[j] ) )
796 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
797 typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
799 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
800 ::TransformPoint(const InputPointType &inputPoint) const
803 InputPointType transformedPoint;
804 OutputPointType outputPoint;
807 if ( m_BulkTransform )
809 transformedPoint = m_BulkTransform->TransformPoint( inputPoint );
813 transformedPoint = inputPoint;
816 // Deformable transform
817 if ( m_CoefficientImage )
819 // Check if inside mask
820 if(m_Mask && !(m_Mask->IsInside(inputPoint) ) )
822 // Outside: no (deformable) displacement
823 return transformedPoint;
826 // Check if inside valid region
828 ContinuousIndexType index;
829 this->TransformPointToContinuousIndex( inputPoint, index );
830 inside = this->InsideValidRegion( index );
833 // Outside: no (deformable) displacement
834 return transformedPoint;
837 // Call the vector interpolator
838 itk::Vector<TCoordRep,OutputDimension> displacement=m_VectorInterpolator->EvaluateAtContinuousIndex(index);
840 // Return the results
841 outputPoint = transformedPoint+displacement;
847 itkWarningMacro( << "B-spline coefficients have not been set" );
848 outputPoint = transformedPoint;
856 //JV Deformably transform a point
857 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
858 typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
860 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
861 ::DeformablyTransformPoint(const InputPointType &inputPoint) const
863 OutputPointType outputPoint;
864 if ( m_CoefficientImage )
867 // Check if inside mask
868 if(m_Mask && !(m_Mask->IsInside(inputPoint) ) )
870 // Outside: no (deformable) displacement
876 ContinuousIndexType index;
877 this->TransformPointToContinuousIndex( inputPoint, index );
878 inside = this->InsideValidRegion( index );
882 //outside: no (deformable) displacement
884 //return outputPoint;
887 // Call the vector interpolator
888 itk::Vector<TCoordRep,OutputDimension> displacement=m_VectorInterpolator->EvaluateAtContinuousIndex(index);
890 // Return the results
891 outputPoint = inputPoint+displacement;
894 // No coefficients available
897 itkWarningMacro( << "B-spline coefficients have not been set" );
898 outputPoint = inputPoint;
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>
909 typename BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
911 BSplineDeformableTransform<TCoordRep, NInputDimensions, NOutputDimensions>
912 ::GetJacobian( const InputPointType & point ) const
914 // Can only compute Jacobian if parameters are set via
915 // SetParameters or SetParametersByValue
916 // if( m_InputParametersPointer == NULL )
918 // itkExceptionMacro( <<"Cannot compute Jacobian: parameters not set" );
921 if (m_NeedResetJacobian)
924 //========================================================
925 // For each dimension, copy the weight to the support region
926 //========================================================
928 // Check if inside mask
929 if(m_Mask && !(m_Mask->IsInside(point) ) )
931 // Outside: no (deformable) displacement
932 return this->m_Jacobian;
936 this->TransformPointToContinuousIndex( point, m_Index );
938 // NOTE: if the support region does not lie totally within the grid
939 // we assume zero displacement and return the input point
940 if ( !this->InsideValidRegion( m_Index ) )
942 return this->m_Jacobian;
945 //Compute interpolation weights
946 const WeightsDataType *weights=NULL;
947 m_VectorInterpolator->EvaluateWeightsAtContinuousIndex( m_Index, &weights, m_LastJacobianIndex);
948 m_SupportRegion.SetIndex( m_LastJacobianIndex );
950 //Reset the iterators
952 for ( j = 0; j < OutputDimension; j++ )
953 m_Iterator[j] = IteratorType( m_JacobianImage[j], m_SupportRegion);
955 // For each dimension, copy the weight to the support region
956 while ( ! (m_Iterator[0]).IsAtEnd() )
958 //copy weight to jacobian image
959 for ( j = 0; j < OutputDimension; j++ )
961 m_ZeroVector[j]=*weights;
962 (m_Iterator[j]).Set( m_ZeroVector);
963 m_ZeroVector[j]=itk::NumericTraits<JacobianValueType>::Zero;
966 // go to next coefficient in the support region
969 m_NeedResetJacobian = true;
972 return this->m_Jacobian;
977 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
979 BSplineDeformableTransform<TCoordRep, NInputDimensions,NOutputDimensions>
980 ::TransformPointToContinuousIndex( const InputPointType & point, ContinuousIndexType & index ) const
984 itk::Vector<double, OutputDimension> tvector;
986 for ( j = 0; j < OutputDimension; j++ )
988 tvector[j] = point[j] - this->m_GridOrigin[j];
991 itk::Vector<double, OutputDimension> cvector;
993 cvector = m_PointToIndex * tvector;
995 for ( j = 0; j < OutputDimension; j++ )
997 index[j] = static_cast< typename ContinuousIndexType::CoordRepType >( cvector[j] );
1001 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
1003 BSplineDeformableTransform<TCoordRep, NInputDimensions,NOutputDimensions>
1004 ::SetJacobianImageData(JacobianPixelType * jacobianDataPointer, unsigned dim)
1006 unsigned int numberOfPixels = m_GridRegion.GetNumberOfPixels();
1007 m_JacobianImage[dim]->GetPixelContainer()->SetImportPointer(jacobianDataPointer, numberOfPixels);
1008 return numberOfPixels;
1011 template<class TCoordRep, unsigned int NInputDimensions, unsigned int NOutputDimensions>
1013 BSplineDeformableTransform<TCoordRep, NInputDimensions,NOutputDimensions>
1014 ::ResetJacobian() const
1016 //========================================================
1017 // Zero all components of jacobian
1018 //========================================================
1019 // JV not thread safe (m_LastJacobianIndex), instantiate N transforms
1020 // NOTE: for efficiency, we only need to zero out the coefficients
1021 // that got fill last time this method was called.
1025 //Define the region for each jacobian image
1026 m_SupportRegion.SetIndex(m_LastJacobianIndex);
1028 //Initialize the iterators
1029 for (j = 0; j < OutputDimension; j++)
1030 m_Iterator[j] = IteratorType(m_JacobianImage[j], m_SupportRegion);
1032 //Set the previously-set to zero
1033 while (!(m_Iterator[0]).IsAtEnd())
1035 for (j = 0; j < OutputDimension; j++)
1037 m_Iterator[j].Set(m_ZeroVector);
1041 m_NeedResetJacobian = false;