Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

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:
Inheritance graph
[legend]
Collaboration diagram for itk::fem::SolverCrankNicolson:
Collaboration graph
[legend]

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

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)
virtual void Clear (void)
void DecomposeK (void)
virtual void FinalizeMatrixAfterAssembly (void)
void FindBracketingTriplet (Float *a, Float *b, Float *c)
void GenerateGFN (void)
Float GetCurrentMaxSolution ()
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
void InitializeForSolution ()
void InitializeInterpolationGrid (const VectorType &size, const VectorType &bb1, const VectorType &bb2)
virtual void InitializeLinearSystemWrapper (void)
virtual void InitializeMatrixForAssembly (unsigned int N)
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 ()
 SolverCrankNicolson ()
void UpdateDisplacements (void)
void Write (std::ostream &f)
 ~SolverCrankNicolson ()



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)
Float GetDeformationEnergy (Float t=1.0)
Float GoldenSection (Float tol=0.01, unsigned int MaxIters=25)
Float GSMax (Float a, Float b)
Float GSSign (Float a, Float b)
void PrintDisplacements ()
void PrintForce ()
void ZeroVector (int which=0)



void InitializeInterpolationGrid (const VectorType &size)

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 147 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 167 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

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:
matrix Index 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.
dim This 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.

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

virtual void itk::fem::Solver::AssembleElementMatrix ( Element::Pointer  e  )  [virtual, inherited]

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.

Reimplemented in itk::fem::SolverHyperbolic.

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

Assemble the master force vector.

Parameters:
dim This 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.
void itk::fem::SolverCrankNicolson::AssembleFforTimeStep ( int  dim = 0  ) 

Assemble the master force vector at a given time.

Parameters:
dim This 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.
void itk::fem::Solver::AssembleK ( void   )  [inherited]

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

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.

virtual void itk::fem::Solver::AssembleLandmarkContribution ( Element::Pointer  e,
float   
) [virtual, inherited]

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.

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

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

Float itk::fem::SolverCrankNicolson::BrentsMethod ( Float  tol = 0.01,
unsigned int  MaxIters = 25 
)

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

virtual void itk::fem::Solver::Clear ( void   )  [virtual, inherited]

Cleans all data members, and initializes the solver to initial state.

void itk::fem::Solver::DecomposeK ( void   )  [inherited]

Decompose matrix using svd, qr, whatever ...

Float itk::fem::SolverCrankNicolson::EvaluateResidual ( Float  t = 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

virtual void itk::fem::Solver::FinalizeMatrixAfterAssembly ( void   )  [inline, virtual, inherited]

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 189 of file itkFEMSolver.h.

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

void itk::fem::SolverCrankNicolson::FindBracketingTriplet ( Float a,
Float b,
Float c 
)
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.

Float itk::fem::SolverCrankNicolson::GetCurrentMaxSolution (  )  [inline]

Definition at line 137 of file itkFEMSolverCrankNicolson.h.

References m_CurrentMaxSolution.

Float itk::fem::Solver::GetDeformationEnergy ( unsigned int  SolutionIndex = 0  )  [inherited]

Get the total deformation energy using the chosen solution

Float itk::fem::SolverCrankNicolson::GetDeformationEnergy ( Float  t = 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

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

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

Parameters:
pt Point in global coordinate system.
Note:
Interpolation grid must be initializes before you can call this function.
const InterpolationGridType* itk::fem::Solver::GetInterpolationGrid ( void   )  const [inline, inherited]

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 132 of file itkFEMSolver.h.

References itk::SmartPointer< TObjectType >::GetPointer().

LinearSystemWrapper::Pointer itk::fem::Solver::GetLinearSystemWrapper (  )  [inline, inherited]

Gets the LinearSystemWrapper object.

See also:
SetLinearSystemWrapper

Definition at line 298 of file itkFEMSolver.h.

References itk::fem::Solver::m_ls.

LinearSystemWrapper* itk::fem::SolverCrankNicolson::GetLS (  )  [inline]

Definition at line 135 of file itkFEMSolverCrankNicolson.h.

References itk::fem::Solver::m_ls.

unsigned int itk::fem::Solver::GetNumberOfDegreesOfFreedom ( void   )  [inline, inherited]

Definition at line 257 of file itkFEMSolver.h.

References itk::fem::Solver::NGFN.

Float itk::fem::Solver::GetSolution ( unsigned int  i,
unsigned int  which = 0 
) [inline, inherited]
virtual Float itk::fem::Solver::GetTimeStep ( void   )  const [inline, virtual, inherited]

Returns the time step used for dynamic problems.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 309 of file itkFEMSolver.h.

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;

Float itk::fem::SolverCrankNicolson::GSMax ( Float  a,
Float  b 
) [inline]

add 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 131 of file itkFEMSolverCrankNicolson.h.

Float itk::fem::SolverCrankNicolson::GSSign ( Float  a,
Float  b 
) [inline]

add 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 130 of file itkFEMSolverCrankNicolson.h.

void itk::fem::SolverCrankNicolson::InitializeForSolution (  ) 

helper initialization function before assembly but after generate GFN.

void itk::fem::Solver::InitializeInterpolationGrid ( const VectorType size  )  [inline, inherited]

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

Definition at line 116 of file itkFEMSolver.h.

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

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:
size Vector that represents number of points on a grid in each dimension.
bb1 Lower limit of a bounding box of a grid.
bb2 Upper limit of a bounding box of a grid.
See also:
GetInterpolationGrid

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

virtual void itk::fem::Solver::InitializeLinearSystemWrapper ( void   )  [virtual, inherited]

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

Reimplemented in itk::fem::SolverHyperbolic.

virtual void itk::fem::Solver::InitializeMatrixForAssembly ( unsigned int  N  )  [virtual, inherited]

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

Parameters:
N Size of the matrix.

Reimplemented in itk::fem::SolverHyperbolic.

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

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

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

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

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

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

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

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

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.

void itk::fem::SolverCrankNicolson::SetAlpha ( Float  a = 0.5  )  [inline]

Set stability step for the solution.

Definition at line 106 of file itkFEMSolverCrankNicolson.h.

References m_alpha.

void itk::fem::SolverCrankNicolson::SetDeltatT ( Float  T  )  [inline]

Set time step for the solution. Should be 1/2.

Definition at line 109 of file itkFEMSolverCrankNicolson.h.

References m_deltaT.

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:
ls Pointer 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.
void itk::fem::SolverCrankNicolson::SetRho ( Float  rho  )  [inline]

Set density constant.

Definition at line 112 of file itkFEMSolverCrankNicolson.h.

References m_rho.

virtual void itk::fem::Solver::SetTimeStep ( Float   )  [inline, virtual, inherited]

Sets the time step used for dynamic problems.

Parameters:
dt New time step.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 316 of file itkFEMSolver.h.

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

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

Reimplemented from itk::fem::Solver.

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

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

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

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

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


Member Data Documentation

Definition at line 181 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 183 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 55 of file itkFEMSolver.h.

Definition at line 174 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 176 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 175 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 67 of file itkFEMSolver.h.

Definition at line 171 of file itkFEMSolverCrankNicolson.h.

Referenced by SetAlpha(), and SolverCrankNicolson().

Definition at line 172 of file itkFEMSolverCrankNicolson.h.

Referenced by GetCurrentMaxSolution(), and SolverCrankNicolson().

Definition at line 169 of file itkFEMSolverCrankNicolson.h.

Referenced by SetDeltatT(), and SolverCrankNicolson().

Definition at line 170 of file itkFEMSolverCrankNicolson.h.

Referenced by SetRho(), and SolverCrankNicolson().

Definition at line 73 of file itkFEMSolver.h.

const unsigned int itk::fem::Solver::MaxGridDimensions = 3 [static, inherited]

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.

unsigned int itk::fem::Solver::NGFN [protected, inherited]

Number of global degrees of freedom in a system

Definition at line 323 of file itkFEMSolver.h.

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

unsigned int itk::fem::Solver::NMFC [protected, inherited]

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

Definition at line 329 of file itkFEMSolver.h.

Definition at line 61 of file itkFEMSolver.h.

Definition at line 177 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 178 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 179 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 182 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().

Definition at line 180 of file itkFEMSolverCrankNicolson.h.

Referenced by SolverCrankNicolson().


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

Generated at Sat Apr 17 01:57:59 2010 for ITK by doxygen 1.6.1 written by Dimitri van Heesch, © 1997-2000