Orfeo Toolbox  3.16
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | Protected Attributes
itk::fem::SolverCrankNicolson Class Reference

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

#include <itkFEMSolverCrankNicolson.h>

+ Inheritance diagram for itk::fem::SolverCrankNicolson:
+ Collaboration diagram for itk::fem::SolverCrankNicolson:

List of all members.

Public Types

typedef Element::ArrayType ElementArray
typedef Element::Float Float
typedef itk::Image
< Element::ConstPointer,
MaxGridDimensions
InterpolationGridType
typedef Load::ArrayType LoadArray
typedef Material::ArrayType MaterialArray
typedef Node::ArrayType NodeArray
typedef Element::VectorType VectorType

Public Member Functions

 SolverCrankNicolson ()
 ~SolverCrankNicolson ()
void AddToDisplacements (Float optimum=1.0)
void ApplyBC (int dim=0, unsigned int matrix=0)
virtual void AssembleElementMatrix (Element::Pointer e)
void AssembleF (int dim=0)
void AssembleFforTimeStep (int dim=0)
void AssembleK (void)
void AssembleKandM ()
virtual void AssembleLandmarkContribution (Element::Pointer e, float)
void AverageLastTwoDisplacements (Float t=0.5)
Float BrentsMethod (Float tol=0.01, unsigned int MaxIters=25)
virtual void Clear (void)
void DecomposeK (void)
Float EvaluateResidual (Float t=1.0)
virtual void FinalizeMatrixAfterAssembly (void)
void FindBracketingTriplet (Float *a, Float *b, Float *c)
void GenerateGFN (void)
Float GetCurrentMaxSolution ()
Float GetDeformationEnergy (Float t=1.0)
Float GetDeformationEnergy (unsigned int SolutionIndex=0)
const ElementGetElementAtPoint (const VectorType &pt) const
const InterpolationGridTypeGetInterpolationGrid (void) const
LinearSystemWrapper::Pointer GetLinearSystemWrapper ()
LinearSystemWrapperGetLS ()
unsigned int GetNumberOfDegreesOfFreedom (void)
Float GetSolution (unsigned int i, unsigned int which=0)
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)
void InitializeForSolution ()
void InitializeInterpolationGrid (const VectorType &size, const VectorType &bb1, const VectorType &bb2)
void InitializeInterpolationGrid (const VectorType &size)
virtual void InitializeLinearSystemWrapper (void)
virtual void InitializeMatrixForAssembly (unsigned int N)
void PrintDisplacements ()
void PrintForce ()
void PrintMinMaxOfSolution ()
void Read (std::istream &f)
void RecomputeForceVector (unsigned int index)
void SetAlpha (Float a=0.5)
void SetDeltatT (Float T)
void SetEnergyToMin (Float xmin)
void SetLinearSystemWrapper (LinearSystemWrapper::Pointer ls)
void SetRho (Float rho)
virtual void SetTimeStep (Float)
void Solve ()
void UpdateDisplacements (void)
void Write (std::ostream &f)
void ZeroVector (int which=0)

Public Attributes

unsigned int DifferenceMatrixIndex
unsigned int DiffMatrixBySolutionTMinus1Index
ElementArray el
unsigned int ForceTIndex
unsigned int ForceTMinus1Index
unsigned int ForceTotalIndex
LoadArray load
Float m_alpha
Float m_CurrentMaxSolution
Float m_deltaT
Float m_rho
MaterialArray mat
NodeArray node
unsigned int SolutionTIndex
unsigned int SolutionTMinus1Index
unsigned int SolutionVectorTMinus1Index
unsigned int SumMatrixIndex
unsigned int TotalSolutionIndex

Static Public Attributes

static const unsigned int MaxGridDimensions = 3

Protected Attributes

LinearSystemWrapper::Pointer m_ls
unsigned int NGFN
unsigned int NMFC

Detailed Description

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. One must call AssembleKandM instead of AssembleK and AssembleFforTimeStep instead of AssembleF. FIXMEs: 1) Members should be privatized, etc. 2) 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 64 of file itkFEMSolverCrankNicolson.h.


Member Typedef Documentation

Array that holds pointers to all elements. since we want to be able to manipulate the array we have to use special pointers

Definition at line 54 of file itkFEMSolver.h.

Local float type

Definition at line 48 of file itkFEMSolver.h.

Type used to store interpolation grid

Definition at line 95 of file itkFEMSolver.h.

Array that holds special pointers to all external loads

Definition at line 66 of file itkFEMSolver.h.

Array that holds pointers to the materials

Definition at line 72 of file itkFEMSolver.h.

Array that holds special pointers to the nodes

Definition at line 60 of file itkFEMSolver.h.

VectorType from the Element base class

Definition at line 78 of file itkFEMSolver.h.


Constructor & Destructor Documentation

itk::fem::SolverCrankNicolson::SolverCrankNicolson ( )
inline

Default constructor which sets the indices for the matrix and vector storage. Time step and other parameters are also initialized.

Definition at line 142 of file itkFEMSolverCrankNicolson.h.

References DifferenceMatrixIndex, DiffMatrixBySolutionTMinus1Index, ForceTIndex, ForceTMinus1Index, ForceTotalIndex, m_alpha, m_CurrentMaxSolution, m_deltaT, m_rho, SolutionTIndex, SolutionTMinus1Index, SolutionVectorTMinus1Index, SumMatrixIndex, and TotalSolutionIndex.

itk::fem::SolverCrankNicolson::~SolverCrankNicolson ( )
inline

Definition at line 162 of file itkFEMSolverCrankNicolson.h.


Member Function Documentation

void itk::fem::SolverCrankNicolson::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

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

Copy the resulting displacements from solution vector back to node objects.

Definition at line 613 of file itkFEMSolverCrankNicolson.cxx.

References itk::fem::LinearSystemWrapper::AddSolutionValue(), itk::fem::LinearSystemWrapper::AddVectorValue(), ForceTIndex, ForceTMinus1Index, ForceTotalIndex, itk::fem::Solution::GetSolutionValue(), itk::fem::LinearSystemWrapper::GetVectorValue(), m_CurrentMaxSolution, itk::fem::Solver::m_ls, itk::fem::Solver::NGFN, itk::fem::LinearSystemWrapper::SetSolutionValue(), itk::fem::LinearSystemWrapper::SetVectorValue(), SolutionTIndex, SolutionTMinus1Index, SolutionVectorTMinus1Index, and TotalSolutionIndex.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::IterativeSolve().

void itk::fem::Solver::ApplyBC ( int  dim = 0,
unsigned int  matrix = 0 
)
inherited

Apply the boundary conditions to the system.

Note:
This function must be called after AssembleK().
Parameters:
matrixIndex of a matrix, to which the BCs should be applied (master stiffness matrix). Normally this is zero, but in derived classes many matrices may be used and this index must be specified.
dimThis is a parameter that can be passed to the function and is normally used with isotropic elements to specify the dimension in which the DOF is fixed.

Apply the boundary conditions to the system.

   Store a temporary pointer to load object for later,
   so that we don't have to access it via the iterator
   Apply boundary conditions in form of MFC loads.

   We add the multi freedom constraints contribution to the master
   stiffness matrix using the lagrange multipliers. Basically we only
   change the last couple of rows and columns in K.
   Apply essential boundary conditions

Definition at line 643 of file itkFEMSolver.cxx.

References itk::fem::LinearSystemWrapper::AddVectorValue(), itk::fem::LinearSystemWrapper::DestroyVector(), itk::fem::LinearSystemWrapper::GetColumnsOfNonZeroMatrixElementsInRow(), itk::fem::LinearSystemWrapper::GetMatrixValue(), itk::fem::LinearSystemWrapper::InitializeVector(), itk::fem::LinearSystemWrapper::IsVectorInitialized(), itk::fem::Solver::load, itk::fem::Solver::m_ls, itk::fem::Solver::NGFN, and itk::fem::LinearSystemWrapper::SetMatrixValue().

Referenced by AssembleKandM(), and itk::fem::Solver::FinalizeMatrixAfterAssembly().

void itk::fem::Solver::AssembleElementMatrix ( Element::Pointer  e)
virtualinherited

Copy the element stiffness matrix into the correct position in the master stiffess matrix. Since more complex Solver classes may need to assemble many matrices and may also do some funky stuff to them, this function is virtual and can be overriden in a derived solver class.

Here we finaly update the corresponding element in the master stiffness matrix. We first check if element in Ke is zero, to prevent zeros from being allocated in sparse matrix.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 347 of file itkFEMSolver.cxx.

References itk::fem::LinearSystemWrapper::AddMatrixValue(), itk::fem::Element::GetDegreeOfFreedom(), itk::fem::Element::GetNumberOfDegreesOfFreedom(), itk::fem::Element::GetStiffnessMatrix(), itk::fem::Solver::m_ls, and itk::fem::Solver::NGFN.

Referenced by itk::fem::Solver::AssembleK().

void itk::fem::Solver::AssembleF ( int  dim = 0)
inherited

Assemble the master force vector.

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.

Assemble the master force vector

 Convert the external loads to the nodal loads and
 add them to the master force vector F.
   Store a temporary pointer to load object for later,
   so that we don't have to access it via the iterator
   Pass the vector to the solution to the Load object.
   Here we only handle Nodal loads
       If using the extra dim parameter, we can apply the force to
       different isotropic dimension.

       FIXME: We assume that the implementation of force vector
       inside the LoadNode class is correct for given number of
       dimensions
   Element loads...
       If array of element pointers is not empty,
       we apply the load to all elements in that array
       If the list of element pointers in load object is empty,
       we apply the load to all elements in a system.
   Handle boundary conditions in form of MFC loads are handled next.
   Handle essential boundary conditions.
   If we got here, we were unable to handle that class of Load object.
   We do nothing...
 Adjust the master force vector for essential boundary
 conditions as required.

Definition at line 389 of file itkFEMSolver.cxx.

References itk::fem::LinearSystemWrapper::AddVectorValue(), itk::fem::Solver::el, itk::fem::Element::GetDegreeOfFreedom(), itk::fem::Element::GetLoadVector(), itk::fem::Element::GetNumberOfDegreesOfFreedom(), itk::fem::LinearSystemWrapper::GetVectorValue(), itk::fem::LinearSystemWrapper::InitializeVector(), itk::fem::LinearSystemWrapper::IsVectorInitialized(), itk::fem::Solver::load, itk::fem::Solver::m_ls, itk::fem::Solver::NGFN, itk::fem::Solver::NMFC, and itk::fem::LinearSystemWrapper::SetVectorValue().

Referenced by AssembleFforTimeStep(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::IterativeSolve().

void itk::fem::SolverCrankNicolson::AssembleFforTimeStep ( int  dim = 0)

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.

Assemble the master force vector

Definition at line 203 of file itkFEMSolverCrankNicolson.cxx.

References itk::fem::Solver::AssembleF(), DifferenceMatrixIndex, DiffMatrixBySolutionTMinus1Index, ForceTIndex, itk::fem::Solver::load, itk::fem::Solver::m_ls, itk::fem::LinearSystemWrapper::MultiplyMatrixVector(), itk::fem::Solver::NGFN, RecomputeForceVector(), itk::fem::LinearSystemWrapper::SetSolutionValue(), itk::fem::LinearSystemWrapper::SetVectorValue(), SolutionTMinus1Index, SolutionVectorTMinus1Index, and TotalSolutionIndex.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::IterativeSolve().

void itk::fem::Solver::AssembleK ( void  )
inherited

Assemble the master stiffness matrix (also apply the MFCs to K)

 Before we can start the assembly procedure, we need to know,
 how many boundary conditions if form of MFCs are there in a system.
 Now we can assemble the master stiffness matrix from
 element stiffness matrices.

 Since we're using the Lagrange multiplier method to apply the MFC,
 each constraint adds a new global DOF.
 Step over all elements
 Step over all the loads again to add the landmark contributions
 to the appropriate place in the stiffness matrix

Definition at line 241 of file itkFEMSolver.cxx.

References itk::fem::Solver::AssembleElementMatrix(), itk::fem::Solver::AssembleLandmarkContribution(), itk::fem::Solver::el, itk::fem::Solver::FinalizeMatrixAfterAssembly(), itk::fem::Solver::InitializeMatrixForAssembly(), itk::fem::Solver::load, itk::fem::Solver::NGFN, and itk::fem::Solver::NMFC.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MultiResSolve(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::RunRegistration().

void itk::fem::SolverCrankNicolson::AssembleKandM ( )

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

Assemble the master stiffness matrix (also apply the MFCs to K)

 Now we can assemble the master stiffness matrix
 from element stiffness matrices
 Step over all elements
 Step over all the loads to add the landmark contributions to the
 appropriate place in the stiffness matrix

Definition at line 57 of file itkFEMSolverCrankNicolson.cxx.

References itk::fem::LinearSystemWrapper::AddMatrixValue(), itk::fem::Solver::ApplyBC(), DifferenceMatrixIndex, itk::fem::Solver::el, itk::fem::Element::GetDegreeOfFreedom(), itk::fem::Element::GetLandmarkContributionMatrix(), itk::fem::Element::GetNumberOfDegreesOfFreedom(), InitializeForSolution(), itk::fem::Solver::load, m_alpha, m_deltaT, itk::fem::Solver::m_ls, m_rho, itk::fem::Solver::NGFN, itk::fem::Solver::NMFC, and SumMatrixIndex.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MultiResSolve(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::RunRegistration().

void itk::fem::Solver::AssembleLandmarkContribution ( Element::Pointer  e,
float  eta 
)
virtualinherited

Add the contribution of the landmark-containing elements to the correct position in the master stiffess matrix. Since more complex Solver classes may need to assemble many matrices and may also do some funky stuff to them, this function is virtual and can be overriden in a derived solver class.

Here we finaly update the corresponding element in the master stiffness matrix. We first check if element in Le is zero, to prevent zeros from being allocated in sparse matrix.

Definition at line 311 of file itkFEMSolver.cxx.

References itk::fem::LinearSystemWrapper::AddMatrixValue(), itk::fem::Element::GetDegreeOfFreedom(), itk::fem::Element::GetLandmarkContributionMatrix(), itk::fem::Element::GetNumberOfDegreesOfFreedom(), itk::fem::Solver::m_ls, and itk::fem::Solver::NGFN.

Referenced by itk::fem::Solver::AssembleK().

void itk::fem::SolverCrankNicolson::AverageLastTwoDisplacements ( Float  t = 0.5)

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

Definition at line 697 of file itkFEMSolverCrankNicolson.cxx.

References itk::fem::Solution::GetSolutionValue(), itk::fem::Solver::m_ls, itk::fem::Solver::NGFN, itk::fem::LinearSystemWrapper::SetSolutionValue(), itk::fem::LinearSystemWrapper::SetVectorValue(), SolutionTIndex, SolutionTMinus1Index, and SolutionVectorTMinus1Index.

Element::Float itk::fem::SolverCrankNicolson::BrentsMethod ( Float  tol = 0.01,
unsigned int  MaxIters = 25 
)
void itk::fem::Solver::Clear ( void  )
virtualinherited
void itk::fem::Solver::DecomposeK ( void  )
inherited

Decompose matrix using svd, qr, whatever ...

Decompose matrix using svd, qr, whatever ... if needed

Definition at line 586 of file itkFEMSolver.cxx.

Element::Float itk::fem::SolverCrankNicolson::EvaluateResidual ( Float  t = 1.0)
virtual void itk::fem::Solver::FinalizeMatrixAfterAssembly ( void  )
inlinevirtualinherited

This function is called after the assebly has been completed. In this class it is only used to apply the BCs. You may however use it to perform other stuff in derived solver classes.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 188 of file itkFEMSolver.h.

References itk::fem::Solver::ApplyBC().

Referenced by itk::fem::Solver::AssembleK().

void itk::fem::SolverCrankNicolson::FindBracketingTriplet ( Float a,
Float b,
Float c 
)

Definition at line 264 of file itkFEMSolverCrankNicolson.cxx.

References EvaluateResidual(), GSMax(), and GSSign().

Referenced by BrentsMethod(), and GoldenSection().

void itk::fem::Solver::GenerateGFN ( void  )
inherited

System solver functions. Call all six functions below (in listed order) to solve system. Assign a global freedom numbers to each DOF in a system. This must be done before any other solve function can be called.

Assign a global freedom number to each DOF in a system.

Assign new ID to every DOF in a system

Definition at line 186 of file itkFEMSolver.cxx.

References itk::fem::Solver::el, itk::fem::Element::InvalidDegreeOfFreedomID, itk::fem::Solver::NGFN, and itk::fem::Solver::node.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::CreateMesh().

Float itk::fem::SolverCrankNicolson::GetCurrentMaxSolution ( )
inline
Element::Float itk::fem::SolverCrankNicolson::GetDeformationEnergy ( Float  t = 1.0)
Solver::Float itk::fem::Solver::GetDeformationEnergy ( unsigned int  SolutionIndex = 0)
inherited

Get the total deformation energy using the chosen solution

Definition at line 620 of file itkFEMSolver.cxx.

References itk::fem::Solver::el, itk::fem::Solution::GetSolutionValue(), and itk::fem::Solver::m_ls.

const Element * itk::fem::Solver::GetElementAtPoint ( const VectorType pt) const
inherited

Returns the pointer to the element which contains global point pt.

Parameters:
ptPoint in global coordinate system.
Note:
Interpolation grid must be initializes before you can call this function.

Definition at line 883 of file itkFEMSolver.cxx.

References itk::fem::Solver::m_InterpolationGrid, and itk::fem::Solver::MaxGridDimensions.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::InterpolateVectorField().

const InterpolationGridType* itk::fem::Solver::GetInterpolationGrid ( void  ) const
inlineinherited

Returns pointer to interpolation grid, which is an itk::Image of pointers to Element objects. Normally you would use physical coordinates to get specific points (pointers to elements) from the image. You can then use the Elemenet::InterpolateSolution member function on the returned element to obtain the solution at this point.

Note:
Physical coordinates in an image correspond to the global coordinate system in which the mesh (nodes) are.

Definition at line 131 of file itkFEMSolver.h.

References itk::SmartPointer< TObjectType >::GetPointer(), and itk::fem::Solver::m_InterpolationGrid.

LinearSystemWrapper::Pointer itk::fem::Solver::GetLinearSystemWrapper ( )
inlineinherited
LinearSystemWrapper* itk::fem::SolverCrankNicolson::GetLS ( )
inline
unsigned int itk::fem::Solver::GetNumberOfDegreesOfFreedom ( void  )
inlineinherited
Float itk::fem::Solver::GetSolution ( unsigned int  i,
unsigned int  which = 0 
)
inlineinherited
virtual Float itk::fem::Solver::GetTimeStep ( void  ) const
inlinevirtualinherited

Returns the time step used for dynamic problems.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 308 of file itkFEMSolver.h.

Element::Float itk::fem::SolverCrankNicolson::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;

Definition at line 445 of file itkFEMSolverCrankNicolson.cxx.

References EvaluateResidual(), FindBracketingTriplet(), and SetEnergyToMin().

Float itk::fem::SolverCrankNicolson::GSMax ( Float  a,
Float  b 
)
inline
Float itk::fem::SolverCrankNicolson::GSSign ( Float  a,
Float  b 
)
inline
void itk::fem::SolverCrankNicolson::InitializeForSolution ( )
void itk::fem::Solver::InitializeInterpolationGrid ( const VectorType size,
const VectorType bb1,
const VectorType bb2 
)
inherited

Initialize the interpolation grid. The interpolation grid is used to find elements that containg specific points in a mesh. The interpolation grid stores pointers to elements for each point on a grid thereby providing a fast way (lookup table) to perform interpolation of results.

Note:
Interpolation grid must be reinitialized each time a mesh changes.
Parameters:
sizeVector that represents number of points on a grid in each dimension.
bb1Lower limit of a bounding box of a grid.
bb2Upper limit of a bounding box of a grid.
See also:
GetInterpolationGrid

Initialize the interpolation grid

Definition at line 758 of file itkFEMSolver.cxx.

References itk::fem::Solver::el, itk::ImageConstIterator< TImage >::GetIndex(), itk::ImageRegionConstIterator< TImage >::GoToBegin(), itk::ImageConstIterator< TImage >::IsAtEnd(), itk::fem::Solver::m_InterpolationGrid, itk::fem::Solver::MaxGridDimensions, itk::Image< TPixel, VImageDimension >::New(), and itk::ImageRegionIterator< TImage >::Set().

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::CreateMesh(), and itk::fem::Solver::InitializeInterpolationGrid().

void itk::fem::Solver::InitializeInterpolationGrid ( const VectorType size)
inlineinherited

Same as InitializeInterpolationGrid(size, {0,0...}, size);

Definition at line 116 of file itkFEMSolver.h.

References itk::fem::Solver::InitializeInterpolationGrid().

void itk::fem::Solver::InitializeLinearSystemWrapper ( void  )
virtualinherited

Performs any initialization needed for LinearSystemWrapper object i.e. sets the maximum number of matrices and vectors.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 72 of file itkFEMSolver.cxx.

References itk::fem::Solver::m_ls, itk::fem::LinearSystemWrapper::SetNumberOfMatrices(), itk::fem::LinearSystemWrapper::SetNumberOfSolutions(), and itk::fem::LinearSystemWrapper::SetNumberOfVectors().

Referenced by itk::fem::Solver::SetLinearSystemWrapper().

void itk::fem::Solver::InitializeMatrixForAssembly ( unsigned int  N)
virtualinherited

This function is called before assembling the matrices. You can override it in a derived class to account for special needs.

Parameters:
NSize of the matrix.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 303 of file itkFEMSolver.cxx.

References itk::fem::LinearSystemWrapper::InitializeMatrix(), itk::fem::Solver::m_ls, and itk::fem::LinearSystemWrapper::SetSystemOrder().

Referenced by itk::fem::Solver::AssembleK().

void itk::fem::SolverCrankNicolson::PrintDisplacements ( )
void itk::fem::SolverCrankNicolson::PrintForce ( )
void itk::fem::SolverCrankNicolson::PrintMinMaxOfSolution ( )

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

Compute maximum and minimum solution values.

Copy the resulting displacements from solution vector back to node objects.

Definition at line 674 of file itkFEMSolverCrankNicolson.cxx.

References itk::fem::Solution::GetSolutionValue(), itk::fem::Solver::m_ls, itk::fem::Solver::NGFN, SolutionTIndex, and TotalSolutionIndex.

void itk::fem::Solver::Read ( std::istream &  f)
inherited

Reads the whole system (nodes, materials and elements) from input stream

   If CreateFromStream returned 0, we're ok. That was the signal
   for the end of stream. Just continue reading... and consequently
   exit the do loop.
   Find out what kind of object did we read from stream
   and store it in the appropriate array
   If we got here, something strange was in the file...

Definition at line 84 of file itkFEMSolver.cxx.

References itk::fem::FEMLightObject::CreateFromStream(), itk::fem::Solver::el, itk::fem::Solver::load, itk::fem::Solver::mat, and itk::fem::Solver::node.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ApplyLoads(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::CreateMesh().

void itk::fem::SolverCrankNicolson::RecomputeForceVector ( unsigned int  index)

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

Definition at line 242 of file itkFEMSolverCrankNicolson.cxx.

References DiffMatrixBySolutionTMinus1Index, ForceTIndex, ForceTMinus1Index, itk::fem::LinearSystemWrapper::GetVectorValue(), m_alpha, m_deltaT, itk::fem::Solver::m_ls, and itk::fem::LinearSystemWrapper::SetVectorValue().

Referenced by AssembleFforTimeStep().

void itk::fem::SolverCrankNicolson::SetAlpha ( Float  a = 0.5)
inline
void itk::fem::SolverCrankNicolson::SetDeltatT ( Float  T)
inline
void itk::fem::SolverCrankNicolson::SetEnergyToMin ( Float  xmin)
void itk::fem::Solver::SetLinearSystemWrapper ( LinearSystemWrapper::Pointer  ls)
inherited

Sets the LinearSystemWrapper object that will be used when solving the master equation. If this function is not called, a default VNL linear system representation will be used (class LinearSystemWrapperVNL).

Parameters:
lsPointer to an object of class which is derived from LinearSystemWrapper.
Note:
Once the LinearSystemWrapper object is changed, it is used until the member function SetLinearSystemWrapper is called again. Since LinearSystemWrapper object was created outside the Solver class, it should also be destroyed outside. Solver class will not destroy it when the Solver object is destroyed.

Change the LinearSystemWrapper object used to solve system of equations.

Definition at line 64 of file itkFEMSolver.cxx.

References itk::fem::Solver::InitializeLinearSystemWrapper(), and itk::fem::Solver::m_ls.

Referenced by itk::fem::Solver::Clear(), itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MultiResSolve(), itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::RunRegistration(), and itk::fem::Solver::Solver().

void itk::fem::SolverCrankNicolson::SetRho ( Float  rho)
inline
virtual void itk::fem::Solver::SetTimeStep ( Float  )
inlinevirtualinherited

Sets the time step used for dynamic problems.

Parameters:
dtNew time step.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 315 of file itkFEMSolver.h.

void itk::fem::SolverCrankNicolson::Solve ( void  )
virtual

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

Solve for the displacement vector u

Reimplemented from itk::fem::Solver.

Definition at line 254 of file itkFEMSolverCrankNicolson.cxx.

References itk::fem::LinearSystemWrapper::InitializeSolution(), itk::fem::Solver::m_ls, SolutionTIndex, and itk::fem::LinearSystemWrapper::Solve().

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::IterativeSolve().

void itk::fem::Solver::UpdateDisplacements ( void  )
inherited

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

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

Definition at line 616 of file itkFEMSolver.cxx.

void itk::fem::Solver::Write ( std::ostream &  f)
inherited

Writes everything (nodes, materials and elements) to output stream

Definition at line 155 of file itkFEMSolver.cxx.

References itk::fem::Solver::el, itk::fem::Solver::load, itk::fem::Solver::mat, and itk::fem::Solver::node.

void itk::fem::SolverCrankNicolson::ZeroVector ( int  which = 0)

Member Data Documentation

unsigned int itk::fem::SolverCrankNicolson::DifferenceMatrixIndex
unsigned int itk::fem::SolverCrankNicolson::DiffMatrixBySolutionTMinus1Index
ElementArray itk::fem::Solver::el
inherited
unsigned int itk::fem::SolverCrankNicolson::ForceTIndex
unsigned int itk::fem::SolverCrankNicolson::ForceTMinus1Index
unsigned int itk::fem::SolverCrankNicolson::ForceTotalIndex
LoadArray itk::fem::Solver::load
inherited
Float itk::fem::SolverCrankNicolson::m_alpha
Float itk::fem::SolverCrankNicolson::m_CurrentMaxSolution
Float itk::fem::SolverCrankNicolson::m_deltaT
LinearSystemWrapper::Pointer itk::fem::Solver::m_ls
protectedinherited
Float itk::fem::SolverCrankNicolson::m_rho

Definition at line 165 of file itkFEMSolverCrankNicolson.h.

Referenced by AssembleKandM(), SetRho(), and SolverCrankNicolson().

MaterialArray itk::fem::Solver::mat
inherited
const unsigned int itk::fem::Solver::MaxGridDimensions = 3
staticinherited

Since the itk::Image is templated over the number of dimensions, we have to know this at compile time. Solver class, however, can handle elements in any number of dimensions. In order to be able to use the Image, we choose the maximum number of space dimension that this function will be able to handle. Any unused dimensions are filled with zero.

For example: If a 2D node coordinates are {1.0,3.0} then the corresponding phisycal point in an image is {1.0,3.0,0.0};

Definition at line 90 of file itkFEMSolver.h.

Referenced by itk::fem::Solver::GetElementAtPoint(), and itk::fem::Solver::InitializeInterpolationGrid().

unsigned int itk::fem::Solver::NGFN
protectedinherited
unsigned int itk::fem::Solver::NMFC
protectedinherited

Number of multi freedom constraints in a system. This member is set in a AssembleK function.

Definition at line 328 of file itkFEMSolver.h.

Referenced by itk::fem::Solver::AssembleF(), itk::fem::Solver::AssembleK(), AssembleKandM(), itk::fem::Solver::Clear(), and InitializeForSolution().

NodeArray itk::fem::Solver::node
inherited
unsigned int itk::fem::SolverCrankNicolson::SolutionTIndex
unsigned int itk::fem::SolverCrankNicolson::SolutionTMinus1Index
unsigned int itk::fem::SolverCrankNicolson::SolutionVectorTMinus1Index
unsigned int itk::fem::SolverCrankNicolson::SumMatrixIndex
unsigned int itk::fem::SolverCrankNicolson::TotalSolutionIndex

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

Generated at Sun Feb 3 2013 03:01:27 for Orfeo Toolbox with doxygen 1.8.1.1