Chaste Commit::baa90ac2819b962188b7562f2326be23c47859a7
|
#include <SimpleLinearEllipticSolver.hpp>
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 Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL > | |
void | ComputeTransformedBasisFunctionDerivatives (const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rReturnValue) |
void | DoAssemble () |
virtual c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(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, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement) |
virtual c_vector< double, PROBLEM_DIM *(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, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement) |
virtual void | AssembleOnElement (Element< ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1) > &rAElem, c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> &rBElem) |
virtual bool | ElementAssemblyCriterion (Element< ELEMENT_DIM, SPACE_DIM > &rElement) |
Protected Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > | |
virtual double | GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown) |
virtual void | ResetInterpolatedQuantities () |
virtual void | IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode) |
virtual void | IncrementInterpolatedGradientQuantities (const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode) |
Additional Inherited Members | |
Protected Types inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL > | |
typedef LinearBasisFunction< ELEMENT_DIM > | BasisFunction |
Solver for solving AbstractLinearEllipticPdes.
Definition at line 49 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 | ||
) |
Constructor.
pMesh | pointer to the mesh |
pPde | pointer to the PDE |
pBoundaryConditions | pointer to the boundary conditions |
Definition at line 74 of file SimpleLinearEllipticSolver.cpp.
References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::mpEllipticPde.
|
protectedvirtual |
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 39 of file SimpleLinearEllipticSolver.cpp.
|
protectedvirtual |
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 62 of file SimpleLinearEllipticSolver.cpp.
|
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 >.
Definition at line 85 of file SimpleLinearEllipticSolver.cpp.
References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve().
Referenced by CellBasedEllipticPdeSolver< DIM >::InitialiseForSolve().
|
inlineprotectedvirtual |
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 110 of file SimpleLinearEllipticSolver.hpp.
References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpLinearSystem, and AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, 1, NORMAL >::SetupGivenLinearSystem().
|
protected |
The PDE to be solved.
Definition at line 56 of file SimpleLinearEllipticSolver.hpp.
Referenced by SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::SimpleLinearEllipticSolver().