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

itk::fem::Solver Class Reference

Main FEM solver. More...

#include <itkFEMSolver.h>

Inheritance diagram for itk::fem::Solver:

Inheritance graph
[legend]
Collaboration diagram for itk::fem::Solver:

Collaboration graph
[legend]
List of all members.

Public Types

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

Public Member Functions

 itkStaticConstMacro (MaxGridDimensions, unsigned int, 3)
void InitializeInterpolationGrid (const VectorType &size, const VectorType &bb1, const VectorType &bb2)
const InterpolationGridTypeGetInterpolationGrid (void) const
const ElementGetElementAtPoint (const VectorType &pt) const
void Read (std::istream &f)
void Write (std::ostream &f)
virtual void Clear (void)
void GenerateGFN (void)
void AssembleK (void)
virtual void InitializeMatrixForAssembly (unsigned int N)
virtual void AssembleElementMatrix (Element::Pointer e)
virtual void AssembleLandmarkContribution (const Element::Pointer e, float)
void ApplyBC (int dim=0, unsigned int matrix=0)
void AssembleF (int dim=0)
void DecomposeK (void)
virtual void Solve (void)
void UpdateDisplacements (void)
Float GetSolution (unsigned int i, unsigned int which=0)
unsigned int GetNumberOfDegreesOfFreedom (void)
Float GetDeformationEnergy (unsigned int SolutionIndex=0)
 Solver ()
virtual ~Solver ()
void SetLinearSystemWrapper (LinearSystemWrapper::Pointer ls)
LinearSystemWrapper::Pointer GetLinearSystemWrapper ()
virtual void InitializeLinearSystemWrapper (void)
virtual Float GetTimeStep (void) const
virtual void SetTimeStep (Float)
void InitializeInterpolationGrid (const VectorType &size)
virtual void FinalizeMatrixAfterAssembly (void)

Public Attributes

ElementArray el
NodeArray node
LoadArray load
MaterialArray mat

Protected Attributes

unsigned int NGFN
unsigned int NMFC
LinearSystemWrapper::Pointer m_ls

Detailed Description

Main FEM solver.

This is the main class used for solving the FEM problems. It also stores all the objects that define the specific FEM problem. Normally there is one Solver object for each FEM problem.

Definition at line 41 of file itkFEMSolver.h.


Member Typedef Documentation

typedef Element::ArrayType itk::fem::Solver::ElementArray
 

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.

typedef Element::Float itk::fem::Solver::Float
 

Local float type Definition at line 48 of file itkFEMSolver.h.

typedef itk::Image<Element::ConstPointer,MaxGridDimensions> itk::fem::Solver::InterpolationGridType
 

Type used to store interpolation grid Definition at line 95 of file itkFEMSolver.h.

typedef Load::ArrayType itk::fem::Solver::LoadArray
 

Array that holds special pointers to all external loads Definition at line 66 of file itkFEMSolver.h.

typedef Material::ArrayType itk::fem::Solver::MaterialArray
 

Array that holds pointers to the materials Definition at line 72 of file itkFEMSolver.h.

typedef Node::ArrayType itk::fem::Solver::NodeArray
 

Array that holds special pointers to the nodes Definition at line 60 of file itkFEMSolver.h.

typedef Element::VectorType itk::fem::Solver::VectorType
 

VectorType from the Element base class Definition at line 78 of file itkFEMSolver.h.

Referenced by InitializeInterpolationGrid().


Constructor & Destructor Documentation

itk::fem::Solver::Solver  ) 
 

Default constructor sets Solver to use VNL linear system .

See also:
Solver::SetLinearSystemWrapper

virtual itk::fem::Solver::~Solver  )  [inline, virtual]
 

Virtual destructor Definition at line 276 of file itkFEMSolver.h.


Member Function Documentation

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

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 FinalizeMatrixAfterAssembly().

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

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  ) 
 

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::Solver::AssembleK void   ) 
 

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

virtual void itk::fem::Solver::AssembleLandmarkContribution const Element::Pointer  e,
float 
[virtual]
 

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.

virtual void itk::fem::Solver::Clear void   )  [virtual]
 

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

void itk::fem::Solver::DecomposeK void   ) 
 

Decompose matrix using svd, qr, whatever ...

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

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 ApplyBC().

void itk::fem::Solver::GenerateGFN void   ) 
 

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::Solver::GetDeformationEnergy unsigned int  SolutionIndex = 0  ) 
 

Get the total deformation energy using the chosen solution

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

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]
 

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.

LinearSystemWrapper::Pointer itk::fem::Solver::GetLinearSystemWrapper  )  [inline]
 

Gets the LinearSystemWrapper object.

See also:
SetLinearSystemWrapper
Definition at line 299 of file itkFEMSolver.h.

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

Definition at line 258 of file itkFEMSolver.h.

References NGFN.

Float itk::fem::Solver::GetSolution unsigned int  i,
unsigned int  which = 0
[inline]
 

Definition at line 253 of file itkFEMSolver.h.

References m_ls.

virtual Float itk::fem::Solver::GetTimeStep void   )  const [inline, virtual]
 

Returns the time step used for dynamic problems.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 310 of file itkFEMSolver.h.

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

Same as InitializeInterpolationGrid(size, {0,0...}, size); Definition at line 116 of file itkFEMSolver.h.

References InitializeInterpolationGrid(), and VectorType.

void itk::fem::Solver::InitializeInterpolationGrid const VectorType size,
const VectorType bb1,
const VectorType bb2
 

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 InitializeInterpolationGrid().

virtual void itk::fem::Solver::InitializeLinearSystemWrapper void   )  [virtual]
 

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]
 

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.

itk::fem::Solver::itkStaticConstMacro MaxGridDimensions  ,
unsigned  int,
 

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};

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

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

void itk::fem::Solver::SetLinearSystemWrapper LinearSystemWrapper::Pointer  ls  ) 
 

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.

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

Sets the time step used for dynamic problems.

Parameters:
dt New time step.

Reimplemented in itk::fem::SolverHyperbolic.

Definition at line 317 of file itkFEMSolver.h.

References m_ls, NGFN, NMFC, and itk::fem::LinearSystemWrapper::Pointer.

virtual void itk::fem::Solver::Solve void   )  [virtual]
 

Solve for the displacement vector u. May be overriden in derived classes.

Reimplemented in itk::fem::SolverHyperbolic.

void itk::fem::Solver::UpdateDisplacements void   ) 
 

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  ) 
 

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


Member Data Documentation

ElementArray itk::fem::Solver::el
 

Definition at line 55 of file itkFEMSolver.h.

LoadArray itk::fem::Solver::load
 

Definition at line 67 of file itkFEMSolver.h.

LinearSystemWrapper::Pointer itk::fem::Solver::m_ls [protected]
 

Pointer to LinearSystemWrapper object. Definition at line 333 of file itkFEMSolver.h.

Referenced by GetSolution(), and SetTimeStep().

MaterialArray itk::fem::Solver::mat
 

Definition at line 73 of file itkFEMSolver.h.

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

Number of global degrees of freedom in a system Definition at line 324 of file itkFEMSolver.h.

Referenced by GetNumberOfDegreesOfFreedom(), and SetTimeStep().

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

Number of multi freedom constraints in a system. This member is set in a AssembleK function. Definition at line 330 of file itkFEMSolver.h.

Referenced by SetTimeStep().

NodeArray itk::fem::Solver::node
 

Definition at line 61 of file itkFEMSolver.h.


The documentation for this class was generated from the following file:
Generated at Sun Apr 1 03:23:42 2007 for ITK by doxygen 1.3.8 written by Dimitri van Heesch, © 1997-2000