Chaste Release::3.1
SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM > Class Template Reference

#include <SimpleLinearEllipticSolver.hpp>

Inheritance diagram for SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >:
Collaboration diagram for SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >:

List of all members.

Public Member Functions

 SimpleLinearEllipticSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, AbstractLinearEllipticPde< ELEMENT_DIM, SPACE_DIM > *pPde, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions, unsigned numQuadPoints=2)
void InitialiseForSolve (Vec initialSolution=NULL)

Protected Member Functions

virtual c_matrix< double,
1 *(ELEMENT_DIM+1),
1 *(ELEMENT_DIM+1)> 
ComputeMatrixTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, 1 > &rU, c_matrix< double, 1, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
virtual c_vector< double,
1 *(ELEMENT_DIM+1)> 
ComputeVectorTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, 1 > &rU, c_matrix< double, 1, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement)
void SetupLinearSystem (Vec currentSolution, bool computeMatrix)

Protected Attributes

AbstractLinearEllipticPde
< ELEMENT_DIM, SPACE_DIM > * 
mpEllipticPde

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
class SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >

SimpleLinearEllipticSolver.

Solver for solving AbstractLinearEllipticPdes.

Definition at line 49 of file SimpleLinearEllipticSolver.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SimpleLinearEllipticSolver ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
AbstractLinearEllipticPde< ELEMENT_DIM, SPACE_DIM > *  pPde,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *  pBoundaryConditions,
unsigned  numQuadPoints = 2 
)

Constructor.

Parameters:
pMeshpointer to the mesh
pPdepointer to the PDE
pBoundaryConditionspointer to the boundary conditions
numQuadPointsnumber of quadrature points (defaults to 2)

Definition at line 75 of file SimpleLinearEllipticSolver.cpp.

References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_matrix< double, 1 *(ELEMENT_DIM+1), 1 *(ELEMENT_DIM+1)> SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm ( c_vector< double, ELEMENT_DIM+1 > &  rPhi,
c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, 1 > &  rU,
c_matrix< double, 1, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
) [protected, virtual]

The term to be added to the element stiffness matrix - see AbstractFeVolumeIntegralAssembler

grad_phi[row] . ( pde_diffusion_term * grad_phi[col])

Parameters:
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer to the element

Definition at line 39 of file SimpleLinearEllipticSolver.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
c_vector< double, 1 *(ELEMENT_DIM+1)> SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm ( c_vector< double, ELEMENT_DIM+1 > &  rPhi,
c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, 1 > &  rU,
c_matrix< double, 1, SPACE_DIM > &  rGradU,
Element< ELEMENT_DIM, SPACE_DIM > *  pElement 
) [protected, virtual]

The term arising from boundary conditions to be added to the element stiffness vector - see AbstractFeVolumeIntegralAssembler

Parameters:
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer to the element

Definition at line 62 of file SimpleLinearEllipticSolver.cpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve ( Vec  initialSolution = NULL) [virtual]

Overloaded InitaliseForSolve() which just calls the base class but also sets the matrix as symmetric and sets Conjugate Gradients as the solver

Parameters:
initialSolutioninitialSolution (used in base class version of this method)

Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Reimplemented in CellBasedPdeSolver< DIM >.

Definition at line 87 of file SimpleLinearEllipticSolver.cpp.

References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve().

Referenced by CellBasedPdeSolver< DIM >::InitialiseForSolve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
void SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem ( Vec  currentSolution,
bool  computeMatrix 
) [inline, protected, virtual]

Delegate to AbstractAssemblerSolverHybrid::SetupGivenLinearSystem.

Parameters:
currentSolutionThe current solution which can be used in setting up the linear system if needed (NULL if there isn't a current solution)
computeMatrixWhether to compute the LHS matrix of the linear system (mainly for dynamic solves).

Implements AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 112 of file SimpleLinearEllipticSolver.hpp.


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
AbstractLinearEllipticPde<ELEMENT_DIM,SPACE_DIM>* SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde [protected]

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