Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
|
#include <CellBasedEllipticPdeSolver.hpp>
Protected Member Functions | |
virtual c_vector< double, 1 *(DIM+1)> | ComputeVectorTerm (c_vector< double, DIM+1 > &rPhi, c_matrix< double, DIM, DIM+1 > &rGradPhi, ChastePoint< DIM > &rX, c_vector< double, 1 > &rU, c_matrix< double, 1, DIM > &rGradU, Element< DIM, DIM > *pElement) |
virtual c_matrix< double, 1 *(DIM+1), 1 *(DIM+1)> | ComputeMatrixTerm (c_vector< double, DIM+1 > &rPhi, c_matrix< double, DIM, DIM+1 > &rGradPhi, ChastePoint< DIM > &rX, c_vector< double, 1 > &rU, c_matrix< double, 1, DIM > &rGradU, Element< DIM, DIM > *pElement) |
void | ResetInterpolatedQuantities () |
void | IncrementInterpolatedQuantities (double phiI, const Node< DIM > *pNode) |
void | InitialiseForSolve (Vec initialSolution) |
Protected Member Functions inherited from SimpleLinearEllipticSolver< DIM, DIM > | |
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 | 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) |
Private Attributes | |
double | mConstantInUSourceTerm |
double | mLinearInUCoeffInSourceTerm |
A purpose-made elliptic solver that interpolates the source terms from node onto Gauss points, as for a cell-based simulation with PDEs the source will only be known at the cells (nodes), not the Gauss points.
Definition at line 48 of file CellBasedEllipticPdeSolver.hpp.
CellBasedEllipticPdeSolver< DIM >::CellBasedEllipticPdeSolver | ( | TetrahedralMesh< DIM, DIM > * | pMesh, |
AbstractLinearEllipticPde< DIM, DIM > * | pPde, | ||
BoundaryConditionsContainer< DIM, DIM, 1 > * | pBoundaryConditions | ||
) |
Constructor.
pMesh | pointer to the mesh |
pPde | pointer to the PDE |
pBoundaryConditions | pointer to the boundary conditions |
Definition at line 39 of file CellBasedEllipticPdeSolver.cpp.
|
virtual |
Destructor.
Definition at line 47 of file CellBasedEllipticPdeSolver.cpp.
|
protectedvirtual |
The SimpleLinearEllipticSolver version of this method is overloaded using the interpolated source term.
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 64 of file CellBasedEllipticPdeSolver.cpp.
|
protectedvirtual |
The SimpleLinearEllipticSolver version of this method is overloaded using the interpolated source term.
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 52 of file CellBasedEllipticPdeSolver.cpp.
|
protected |
Overridden IncrementInterpolatedQuantities() method.
phiI | |
pNode |
Definition at line 94 of file CellBasedEllipticPdeSolver.cpp.
|
protectedvirtual |
Create the linear system object if it hasn't been already. Can use an initial solution as PETSc template, or base it on the mesh size.
initialSolution | an initial guess |
Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 101 of file CellBasedEllipticPdeSolver.cpp.
References SimpleLinearEllipticSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().
|
protectedvirtual |
Overridden ResetInterpolatedQuantities() method.
Reimplemented from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >.
Definition at line 87 of file CellBasedEllipticPdeSolver.cpp.
|
private |
The constant in u part of the source term, interpolated onto the current point.
Definition at line 53 of file CellBasedEllipticPdeSolver.hpp.
|
private |
The linear in u part of the source term, interpolated onto the current point.
Definition at line 56 of file CellBasedEllipticPdeSolver.hpp.