17 #ifndef __itkKernelTransform_txx
18 #define __itkKernelTransform_txx
28 template <
class TScalarType,
unsigned int NDimensions>
40 this->
m_I.set_identity();
52 template <
class TScalarType,
unsigned int NDimensions>
62 template <
class TScalarType,
unsigned int NDimensions>
67 itkDebugMacro(
"setting SourceLandmarks to " << landmarks );
68 if (this->m_SourceLandmarks != landmarks)
70 this->m_SourceLandmarks = landmarks;
71 this->UpdateParameters();
80 template <
class TScalarType,
unsigned int NDimensions>
85 itkDebugMacro(
"setting TargetLandmarks to " << landmarks );
86 if (this->m_TargetLandmarks != landmarks)
88 this->m_TargetLandmarks = landmarks;
89 this->UpdateParameters();
99 #if !defined(ITK_LEGACY_REMOVE)
100 template <
class TScalarType,
unsigned int NDimensions>
105 itkLegacyReplaceBodyMacro( itkKernelTransform::ComputeG_vector, 3.6,
106 itkKernelTransform::ComputeG_vector_gmatrix );
114 template <
class TScalarType,
unsigned int NDimensions>
117 ComputeG(
const InputVectorType &, GMatrixType & itkNotUsed( gmatrix ) )
const
119 itkExceptionMacro( <<
"ComputeG(vector,gmatrix) must be reimplemented"
120 <<
" in subclasses of KernelTransform." );
126 template <
class TScalarType,
unsigned int NDimensions>
127 const typename KernelTransform<TScalarType, NDimensions>::GMatrixType &
131 m_GMatrix.fill( NumericTraits< TScalarType >::Zero );
132 m_GMatrix.fill_diagonal( m_Stiffness );
142 template <
class TScalarType,
unsigned int NDimensions>
149 unsigned long numberOfLandmarks = this->m_SourceLandmarks
150 ->GetNumberOfPoints();
152 PointsIterator sp = this->m_SourceLandmarks->GetPoints()->Begin();
156 for(
unsigned int lnd=0; lnd < numberOfLandmarks; lnd++ )
158 this->ComputeG( thisPoint - sp->Value(), Gmatrix );
159 for(
unsigned int dim=0; dim < NDimensions; dim++ )
161 for(
unsigned int odim=0; odim < NDimensions; odim++ )
163 result[ odim ] += Gmatrix(dim, odim ) * m_DMatrix(dim,lnd);
175 template <
class TScalarType,
unsigned int NDimensions>
179 unsigned long numberOfLandmarks = this->m_SourceLandmarks
180 ->GetNumberOfPoints();
182 PointsIterator sp = this->m_SourceLandmarks->GetPoints()->Begin();
183 PointsIterator tp = this->m_TargetLandmarks->GetPoints()->Begin();
186 this->m_Displacements->Reserve( numberOfLandmarks );
191 vt->
Value() = tp->Value() - sp->Value();
201 template <
class TScalarType,
unsigned int NDimensions>
206 typedef vnl_svd<TScalarType> SVDSolverType;
210 SVDSolverType svd( this->m_LMatrix, 1e-8 );
211 this->m_WMatrix = svd.solve( this->m_YMatrix );
220 template <
class TScalarType,
unsigned int NDimensions>
224 unsigned long numberOfLandmarks = this->m_SourceLandmarks
225 ->GetNumberOfPoints();
227 NDimensions*(NDimensions+1), 0);
232 this->m_LMatrix.set_size(
233 NDimensions*(numberOfLandmarks+NDimensions+1),
234 NDimensions*(numberOfLandmarks+NDimensions+1) );
236 this->m_LMatrix.fill( 0.0 );
238 this->m_LMatrix.update( this->m_KMatrix, 0, 0 );
239 this->m_LMatrix.update( this->m_PMatrix, 0, this->m_KMatrix.columns() );
240 this->m_LMatrix.update( this->m_PMatrix.transpose(),
241 this->m_KMatrix.rows(), 0);
242 this->m_LMatrix.update( O2, this->m_KMatrix.rows(),
243 this->m_KMatrix.columns());
251 template <
class TScalarType,
unsigned int NDimensions>
255 unsigned long numberOfLandmarks = this->m_SourceLandmarks
256 ->GetNumberOfPoints();
261 this->m_KMatrix.set_size( NDimensions * numberOfLandmarks,
262 NDimensions * numberOfLandmarks );
264 this->m_KMatrix.fill( 0.0 );
266 PointsIterator p1 = this->m_SourceLandmarks->GetPoints()->Begin();
278 G = this->ComputeReflexiveG(p1);
279 this->m_KMatrix.update(G, i*NDimensions, i*NDimensions);
287 this->ComputeG(s, G);
289 this->m_KMatrix.update(G, i*NDimensions, j*NDimensions);
290 this->m_KMatrix.update(G, j*NDimensions, i*NDimensions);
302 template <
class TScalarType,
unsigned int NDimensions>
306 unsigned long numberOfLandmarks = this->m_SourceLandmarks
307 ->GetNumberOfPoints();
313 this->m_PMatrix.set_size( NDimensions*numberOfLandmarks,
314 NDimensions*(NDimensions+1) );
315 this->m_PMatrix.fill( 0.0 );
316 for (
unsigned int i = 0; i < numberOfLandmarks; i++)
318 this->m_SourceLandmarks->GetPoint(i, &p);
319 for (
unsigned int j = 0; j < NDimensions; j++)
322 this->m_PMatrix.update(temp, i*NDimensions, j*NDimensions);
324 this->m_PMatrix.update(I, i*NDimensions, NDimensions*NDimensions);
331 template <
class TScalarType,
unsigned int NDimensions>
335 unsigned long numberOfLandmarks = this->m_SourceLandmarks
336 ->GetNumberOfPoints();
339 this->m_Displacements->Begin();
341 this->m_YMatrix.set_size( NDimensions*(numberOfLandmarks+NDimensions+1), 1);
343 this->m_YMatrix.fill( 0.0 );
345 for (
unsigned int i = 0; i < numberOfLandmarks; i++)
347 for (
unsigned int j = 0; j < NDimensions; j++)
349 this->m_YMatrix.put(i*NDimensions+j, 0, displacement.
Value()[j]);
354 for (
unsigned int i = 0; i < NDimensions*(NDimensions+1); i++)
356 this->m_YMatrix.put(numberOfLandmarks*NDimensions+i, 0, 0);
364 template <
class TScalarType,
unsigned int NDimensions>
369 unsigned long numberOfLandmarks = this->m_SourceLandmarks
370 ->GetNumberOfPoints();
373 this->m_DMatrix.set_size(NDimensions,numberOfLandmarks);
375 for(
unsigned int lnd=0; lnd < numberOfLandmarks; lnd++ )
377 for(
unsigned int dim=0; dim < NDimensions; dim++ )
379 this->m_DMatrix(dim,lnd) = this->m_WMatrix(ci++,0);
384 for(
unsigned int j=0; j < NDimensions; j++ )
386 for(
unsigned int i=0; i < NDimensions; i++ )
388 this->m_AMatrix(i,j) = this->m_WMatrix(ci++,0);
393 for(
unsigned int k=0; k < NDimensions; k++ )
395 this->m_BVector(k) = this->m_WMatrix(ci++,0);
406 template <
class TScalarType,
unsigned int NDimensions>
416 result.
Fill( NumericTraits< ValueType >::Zero );
418 this->ComputeDeformationContribution( thisPoint, result );
421 for(
unsigned int j=0; j < NDimensions; j++ )
423 for(
unsigned int i=0; i < NDimensions; i++ )
425 result[i] += this->m_AMatrix(i,j) * thisPoint[j];
432 for(
unsigned int k=0; k < NDimensions; k++ )
434 result[k] += this->m_BVector(k) + thisPoint[k];
442 template <
class TScalarType,
unsigned int NDimensions>
449 this->m_Jacobian.
Fill( 0.0 );
454 itkExceptionMacro(<<
"GetJacobian must be implemented in subclasses"
455 <<
" of KernelTransform.");
457 return this->m_Jacobian;
468 template <
class TScalarType,
unsigned int NDimensions>
473 typename PointsContainer::Pointer landmarks = PointsContainer::New();
474 const unsigned int numberOfLandmarks = parameters.Size() / NDimensions;
475 landmarks->Reserve( numberOfLandmarks );
482 unsigned int pcounter = 0;
485 for(
unsigned int dim=0; dim<NDimensions; dim++)
487 landMark[ dim ] = parameters[ pcounter ];
490 itr.Value() = landMark;
494 this->m_SourceLandmarks->SetPoints( landmarks );
507 template <
class TScalarType,
unsigned int NDimensions>
512 typename PointsContainer::Pointer landmarks = PointsContainer::New();
513 const unsigned int numberOfLandmarks = parameters.Size() / NDimensions;
515 landmarks->Reserve( numberOfLandmarks );
522 unsigned int pcounter = 0;
525 for(
unsigned int dim=0; dim<NDimensions; dim++)
527 landMark[ dim ] = parameters[ pcounter ];
530 itr.Value() = landMark;
534 this->m_TargetLandmarks->SetPoints( landmarks );
540 template <
class TScalarType,
unsigned int NDimensions>
546 ->GetNumberOfPoints()
549 PointsIterator itr = this->m_SourceLandmarks->GetPoints()->Begin();
552 unsigned int pcounter = 0;
556 for(
unsigned int dim=0; dim<NDimensions; dim++)
558 this->m_Parameters[ pcounter ] = landmark[ dim ];
567 template <
class TScalarType,
unsigned int NDimensions>
572 this->UpdateParameters();
573 return this->m_Parameters;
581 template <
class TScalarType,
unsigned int NDimensions>
587 ->GetNumberOfPoints()
590 PointsIterator itr = this->m_TargetLandmarks->GetPoints()->Begin();
593 unsigned int pcounter = 0;
597 for(
unsigned int dim=0; dim<NDimensions; dim++)
599 this->m_FixedParameters[ pcounter ] = landmark[ dim ];
605 return this->m_FixedParameters;
609 template <
class TScalarType,
unsigned int NDimensions>
614 Superclass::PrintSelf(os,indent);
615 if( this->m_SourceLandmarks )
617 os << indent <<
"SourceLandmarks: " << std::endl;
620 if( this->m_TargetLandmarks )
622 os << indent <<
"TargetLandmarks: " << std::endl;
625 if( this->m_Displacements )
627 os << indent <<
"Displacements: " << std::endl;
630 os << indent <<
"Stiffness: " << this->m_Stiffness << std::endl;