#include <SimpleLinearEllipticSolver.hpp>
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) |
virtual c_vector< double, ELEMENT_DIM > | ComputeVectorSurfaceTerm (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > &rSurfaceElement, c_vector< double, ELEMENT_DIM > &rPhi, ChastePoint< SPACE_DIM > &rX) |
void | SetupLinearSystem (Vec currentSolution, bool computeMatrix) |
Protected Attributes | |
AbstractLinearEllipticPde < ELEMENT_DIM, SPACE_DIM > * | mpEllipticPde |
Solver for solving AbstractLinearEllipticPdes.
Definition at line 45 of file SimpleLinearEllipticSolver.hpp.
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 | |||
) | [inline] |
Constructor.
pMesh | pointer to the mesh | |
pPde | pointer to the PDE | |
pBoundaryConditions | pointer to the boundary conditions | |
numQuadPoints | number of quadrature points (defaults to 2) |
Definition at line 84 of file SimpleLinearEllipticSolver.cpp.
References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde.
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 | |||
) | [inline, protected, virtual] |
The term to be added to the element stiffness matrix:
grad_phi[row] . ( pde_diffusion_term * grad_phi[col])
rPhi | The basis functions, rPhi(i) = phi_i, i=1..numBases | |
rGradPhi | Basis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i) | |
rX | The point in space | |
rU | The unknown as a vector, u(i) = u_i | |
rGradU | The gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j) | |
pElement | Pointer to the element |
Definition at line 33 of file SimpleLinearEllipticSolver.cpp.
References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde.
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 | |||
) | [inline, protected, virtual] |
The term arising from boundary conditions to be added to the element stiffness vector.
rPhi | The basis functions, rPhi(i) = phi_i, i=1..numBases | |
rGradPhi | Basis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i) | |
rX | The point in space | |
rU | The unknown as a vector, u(i) = u_i | |
rGradU | The gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j) | |
pElement | Pointer to the element |
Definition at line 57 of file SimpleLinearEllipticSolver.cpp.
References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde.
c_vector< double, ELEMENT_DIM > SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm | ( | const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > & | rSurfaceElement, | |
c_vector< double, ELEMENT_DIM > & | rPhi, | |||
ChastePoint< SPACE_DIM > & | rX | |||
) | [inline, protected, virtual] |
This method returns the vector to be added to element stiffness vector for a given gauss point in BoundaryElement. The arguments are the bases, x and current solution computed at the Gauss point. The returned vector will be multiplied by the gauss weight and jacobian determinent and added to the element stiffness matrix (see AssembleOnElement()).
rSurfaceElement | the element which is being considered. | |
rPhi | The basis functions, rPhi(i) = phi_i, i=1..numBases | |
rX | The point in space |
Reimplemented from AbstractFeObjectAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >.
Definition at line 71 of file SimpleLinearEllipticSolver.cpp.
void SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem | ( | Vec | currentSolution, | |
bool | computeMatrix | |||
) | [inline, protected, virtual] |
Delegate to AbstractAssemblerSolverHybrid::SetupGivenLinearSystem.
currentSolution | The current solution which can be used in setting up the linear system if needed (NULL if there isn't a current solution) | |
computeMatrix | Whether to compute the LHS matrix of the linear system (mainly for dynamic solves). |
Implements AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 116 of file SimpleLinearEllipticSolver.hpp.
void SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve | ( | Vec | initialSolution = NULL |
) | [inline, virtual] |
Overloaded InitaliseForSolve() which just calls the base class but also sets the matrix as symmetric and sets Conjugate Gradients as the solver
initialSolution | initialSolution (used in base class version of this method) |
Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Reimplemented in CellBasedSimulationWithPdesAssembler< DIM >.
Definition at line 97 of file SimpleLinearEllipticSolver.cpp.
References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve(), AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpLinearSystem, LinearSystem::SetKspType(), and LinearSystem::SetMatrixIsSymmetric().
Referenced by CellBasedSimulationWithPdesAssembler< DIM >::InitialiseForSolve().
AbstractLinearEllipticPde<ELEMENT_DIM,SPACE_DIM>* SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde [protected] |
The PDE to be solved.
Definition at line 52 of file SimpleLinearEllipticSolver.hpp.
Referenced by SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeMatrixTerm(), SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::ComputeVectorTerm(), and SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SimpleLinearEllipticSolver().