ITK  4.1.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
itk::fem::SolverCrankNicolson< TDimension > Class Template Reference

#include <itkFEMSolverCrankNicolson.h>

+ Inheritance diagram for itk::fem::SolverCrankNicolson< TDimension >:
+ Collaboration diagram for itk::fem::SolverCrankNicolson< TDimension >:

List of all members.

Public Types

typedef SmartPointer< const SelfConstPointer
typedef Element::Float Float
typedef SmartPointer< SelfPointer
typedef SolverCrankNicolson Self
typedef Solver< TDimension > Superclass

Public Member Functions

void AddToDisplacements (Float optimum=1.0)
void AverageLastTwoDisplacements (Float t=0.5)
Float BrentsMethod (Float tol=0.01, unsigned int MaxIters=25)
Float EvaluateResidual (Float t=1.0)
void FindBracketingTriplet (Float *a, Float *b, Float *c)
Float GetCurrentMaxSolution ()
Float GetDeformationEnergy (Float t=1.0)
LinearSystemWrapperGetLS ()
virtual Float GetTimeStep (void) const
Float GoldenSection (Float tol=0.01, unsigned int MaxIters=25)
Float GSMax (Float a, Float b)
Float GSSign (Float a, Float b)
 itkGetConstMacro (Iterations, unsigned int)
 itkGetMacro (TotalSolutionIndex, unsigned int)
 itkGetMacro (SolutionTMinus1Index, unsigned int)
 itkNewMacro (Self)
 itkTypeMacro (SolverCrankNicolson, Solver< TDimension >)
void PrintDisplacements ()
void PrintForce ()
void PrintMinMaxOfSolution ()
void RecomputeForceVector (unsigned int index)
void ResetIterations (void)
void SetEnergyToMin (Float xmin)
virtual void SetTimeStep (Float dt)
void ZeroVector (int which=0)
 itkSetMacro (UseMassMatrix, bool)
 itkGetMacro (UseMassMatrix, bool)
 itkSetMacro (Alpha, Float)
 itkGetMacro (Alpha, Float)
 itkSetMacro (Rho, Float)
 itkGetMacro (Rho, Float)

Protected Member Functions

void AssembleFforTimeStep (int dim=0)
void AssembleKandM ()
void GenerateData ()
void InitializeForSolution ()
virtual void RunSolver (void)
 SolverCrankNicolson ()
 ~SolverCrankNicolson ()

Protected Attributes

Float m_Alpha
Float m_CurrentMaxSolution
unsigned int m_DifferenceMatrixIndex
unsigned int m_DiffMatrixBySolutionTMinus1Index
unsigned int m_ForceTIndex
unsigned int m_ForceTMinus1Index
unsigned int m_ForceTotalIndex
unsigned int m_Iterations
Float m_Rho
unsigned int m_SolutionTIndex
unsigned int m_SolutionTMinus1Index
unsigned int m_SolutionVectorTMinus1Index
unsigned int m_SumMatrixIndex
Float m_TimeStep
unsigned int m_TotalSolutionIndex
bool m_UseMassMatrix

Private Member Functions

void operator= (const Self &)
 SolverCrankNicolson (const Self &)

Detailed Description

template<unsigned int TDimension = 3>
class itk::fem::SolverCrankNicolson< TDimension >

FEM Solver for time dependent problems; uses Crank-Nicolson implicit discretization scheme.

This is the main class used for solving FEM time-dependent problems. It solves the following problem:

\[ ( M + \alpha*dt* K )*U_t=(M - (1.- \alpha)*dt* K)* U_{t-1} + dt*(\alpha*f_{n+1} + (1-\alpha)*f_n) \]

which is the Crank-Nicolson formulation of the static problem if $\alpha=0.5$. The static solution is gained if : $\rho = 0.0$; $\alpha = 1.0$; $dt = 1.0$; Practically, it is good to set rho to something small (for the itpack solver). The advantage of choosing $\alpha=0.5$ is that the solution is then stable for any choice of time step, dt. This class inherits and uses most of the Solver class functionality.

Updated: The calls to to AssembleKandM (or AssembleK) and AssembleFforTimeStep (or AssembleF) are now handled internally by calling Update().

FIXME: 1) We should also account for the contribution to the force from essential BCs. Basically there are terms involving $ M * (\dot g_b) $ and $ K * g_b $ where $ g_b$ is the essential BC vector.

Definition at line 72 of file itkFEMSolverCrankNicolson.h.


Member Typedef Documentation

template<unsigned int TDimension = 3>
typedef SmartPointer<const Self> itk::fem::SolverCrankNicolson< TDimension >::ConstPointer

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 78 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
typedef Element::Float itk::fem::SolverCrankNicolson< TDimension >::Float

Some convenient typedefs.

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 86 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
typedef SmartPointer<Self> itk::fem::SolverCrankNicolson< TDimension >::Pointer

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 77 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
typedef SolverCrankNicolson itk::fem::SolverCrankNicolson< TDimension >::Self

Standard class typedefs.

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 75 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
typedef Solver<TDimension> itk::fem::SolverCrankNicolson< TDimension >::Superclass

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 76 of file itkFEMSolverCrankNicolson.h.


Constructor & Destructor Documentation

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::SolverCrankNicolson ( ) [protected]
template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::~SolverCrankNicolson ( ) [inline, protected]

Definition at line 207 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::SolverCrankNicolson ( const Self ) [private]

Member Function Documentation

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::AddToDisplacements ( Float  optimum = 1.0)

Add solution vector u to the corresponding nodal values, which are stored in node objects). This is standard post processing of the solution

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::AssembleFforTimeStep ( int  dim = 0) [protected]

Assemble the master force vector at a given time.

Parameters:
dimThis is a parameter that can be passed to the function and is normally used with isotropic elements to specify the dimension for which the master force vector should be assembled.
template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::AssembleKandM ( ) [protected]

Assemble the master stiffness and mass matrix. We actually assemble the right hand side and left hand side of the implicit scheme equation.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::AverageLastTwoDisplacements ( Float  t = 0.5)
template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::BrentsMethod ( Float  tol = 0.01,
unsigned int  MaxIters = 25 
)
template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::EvaluateResidual ( Float  t = 1.0)
template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::FindBracketingTriplet ( Float a,
Float b,
Float c 
)
template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::GenerateData ( ) [protected, virtual]

Method invoked by the pipeline in order to trigger the computation of the registration.

Reimplemented from itk::fem::Solver< TDimension >.

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::GetCurrentMaxSolution ( ) [inline]
template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::GetDeformationEnergy ( Float  t = 1.0)
template<unsigned int TDimension = 3>
LinearSystemWrapper* itk::fem::SolverCrankNicolson< TDimension >::GetLS ( ) [inline]

Definition at line 190 of file itkFEMSolverCrankNicolson.h.

References itk::fem::Solver< TDimension >::m_ls.

template<unsigned int TDimension = 3>
virtual Float itk::fem::SolverCrankNicolson< TDimension >::GetTimeStep ( void  ) const [inline, virtual]

Returns the time step used for dynamic problems.

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 142 of file itkFEMSolverCrankNicolson.h.

References itk::fem::SolverCrankNicolson< TDimension >::m_TimeStep.

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::GoldenSection ( Float  tol = 0.01,
unsigned int  MaxIters = 25 
)

Finds the optimum value between the last two solutions and sets the current solution to that value. Uses Evaluate Residual;

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::GSMax ( Float  a,
Float  b 
) [inline]

Definition at line 183 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::GSSign ( Float  a,
Float  b 
) [inline]

Definition at line 179 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::InitializeForSolution ( ) [protected]

helper initialization function before assembly but after generate GFN.

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkGetConstMacro ( Iterations  ,
unsigned  int 
)

Get the number of iterations run for the solver

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkGetMacro ( UseMassMatrix  ,
bool   
)

Get/Set the use of the Mass Matrix for the solution

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkGetMacro ( TotalSolutionIndex  ,
unsigned  int 
)

Get the index for the current solution

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkGetMacro ( SolutionTMinus1Index  ,
unsigned  int 
)

Get the index for the previous solution

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkGetMacro ( Alpha  ,
Float   
)

Set stability step for the solution. Initialized to 0.5

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkGetMacro ( Rho  ,
Float   
)

Set density constant.

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkNewMacro ( Self  )

Method for creation through the object factory.

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkSetMacro ( UseMassMatrix  ,
bool   
)

Get/Set the use of the Mass Matrix for the solution

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkSetMacro ( Alpha  ,
Float   
)

Set stability step for the solution. Initialized to 0.5

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkSetMacro ( Rho  ,
Float   
)

Set density constant.

template<unsigned int TDimension = 3>
itk::fem::SolverCrankNicolson< TDimension >::itkTypeMacro ( SolverCrankNicolson< TDimension >  ,
Solver< TDimension >   
)

Run-time type information (and related methods)

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::operator= ( const Self ) [private]

Make a DataObject of the correct type to be used as the specified output.

Reimplemented from itk::fem::Solver< TDimension >.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::PrintDisplacements ( )
template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::PrintForce ( )
template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::PrintMinMaxOfSolution ( )

Compute and print the minimum and maximum of the total solution and the last solution.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::RecomputeForceVector ( unsigned int  index)

compute the current state of the right hand side and store the current force for the next iteration.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::ResetIterations ( void  ) [inline]

Reset the number of iterations for the solver. This will prompt the Solver to Assemble the master stiffness and mass matrix again. This is only generated before the first iteration.

Definition at line 106 of file itkFEMSolverCrankNicolson.h.

References itk::fem::SolverCrankNicolson< TDimension >::m_Iterations.

template<unsigned int TDimension = 3>
virtual void itk::fem::SolverCrankNicolson< TDimension >::RunSolver ( void  ) [protected, virtual]

Solve for the displacement vector u at a given time. Update the total solution as well.

Reimplemented from itk::fem::Solver< TDimension >.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::SetEnergyToMin ( Float  xmin)
template<unsigned int TDimension = 3>
virtual void itk::fem::SolverCrankNicolson< TDimension >::SetTimeStep ( Float  dt) [inline, virtual]

Sets the time step used for dynamic problems.

Parameters:
dtNew time step.

Reimplemented from itk::fem::Solver< TDimension >.

Definition at line 152 of file itkFEMSolverCrankNicolson.h.

References itk::fem::SolverCrankNicolson< TDimension >::m_TimeStep.

template<unsigned int TDimension = 3>
void itk::fem::SolverCrankNicolson< TDimension >::ZeroVector ( int  which = 0)

Member Data Documentation

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::m_Alpha [protected]

Definition at line 240 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::m_CurrentMaxSolution [protected]
template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_DifferenceMatrixIndex [protected]

Definition at line 253 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_DiffMatrixBySolutionTMinus1Index [protected]

Definition at line 255 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_ForceTIndex [protected]

Definition at line 246 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_ForceTMinus1Index [protected]

Definition at line 248 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_ForceTotalIndex [protected]

Definition at line 247 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_Iterations [protected]
template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::m_Rho [protected]

Definition at line 239 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_SolutionTIndex [protected]

Definition at line 249 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_SolutionTMinus1Index [protected]

Definition at line 250 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_SolutionVectorTMinus1Index [protected]

Definition at line 251 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_SumMatrixIndex [protected]

Definition at line 254 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
Float itk::fem::SolverCrankNicolson< TDimension >::m_TimeStep [protected]
template<unsigned int TDimension = 3>
unsigned int itk::fem::SolverCrankNicolson< TDimension >::m_TotalSolutionIndex [protected]

Definition at line 252 of file itkFEMSolverCrankNicolson.h.

template<unsigned int TDimension = 3>
bool itk::fem::SolverCrankNicolson< TDimension >::m_UseMassMatrix [protected]

Definition at line 243 of file itkFEMSolverCrankNicolson.h.


The documentation for this class was generated from the following file: