#include <CellBasedPdeSolver.hpp>
Inherits SimpleLinearEllipticSolver< DIM, DIM >.

Public Member Functions | |
| CellBasedPdeSolver (TetrahedralMesh< DIM, DIM > *pMesh, AbstractLinearEllipticPde< DIM, DIM > *pPde, BoundaryConditionsContainer< DIM, DIM, 1 > *pBoundaryConditions) | |
| virtual | ~CellBasedPdeSolver () |
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) |
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 CellBasedPdeSolver.hpp.
| CellBasedPdeSolver< DIM >::CellBasedPdeSolver | ( | TetrahedralMesh< DIM, DIM > * | pMesh, | |
| AbstractLinearEllipticPde< DIM, DIM > * | pPde, | |||
| BoundaryConditionsContainer< DIM, DIM, 1 > * | pBoundaryConditions | |||
| ) | [inline] |
Constructor.
| pMesh | pointer to the mesh | |
| pPde | pointer to the PDE | |
| pBoundaryConditions | pointer to the boundary conditions |
Definition at line 41 of file CellBasedPdeSolver.cpp.
| CellBasedPdeSolver< DIM >::~CellBasedPdeSolver | ( | ) | [inline, virtual] |
Destructor.
Definition at line 49 of file CellBasedPdeSolver.cpp.
| c_matrix< double, 1 *(DIM+1), 1 *(DIM+1)> CellBasedPdeSolver< DIM >::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 | |||
| ) | [inline, protected, virtual] |
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 66 of file CellBasedPdeSolver.cpp.
References CellBasedPdeSolver< DIM >::mLinearInUCoeffInSourceTerm, and SimpleLinearEllipticSolver< DIM, DIM >::mpEllipticPde.
| c_vector< double, 1 *(DIM+1)> CellBasedPdeSolver< DIM >::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 | |||
| ) | [inline, protected, virtual] |
The SimpleLinearEllipticSolver version of this method is overloaded using the interpolated source term.
| rPhi | ||
| rGradPhi | ||
| rX | ||
| rU | ||
| rGradU | ||
| pElement |
Definition at line 54 of file CellBasedPdeSolver.cpp.
References CellBasedPdeSolver< DIM >::mConstantInUSourceTerm.
| void CellBasedPdeSolver< DIM >::IncrementInterpolatedQuantities | ( | double | phiI, | |
| const Node< DIM > * | pNode | |||
| ) | [inline, protected] |
Overridden IncrementInterpolatedQuantities() method.
| phiI | ||
| pNode |
Definition at line 96 of file CellBasedPdeSolver.cpp.
References CellBasedPdeSolver< DIM >::mConstantInUSourceTerm, CellBasedPdeSolver< DIM >::mLinearInUCoeffInSourceTerm, and SimpleLinearEllipticSolver< DIM, DIM >::mpEllipticPde.
| void CellBasedPdeSolver< DIM >::InitialiseForSolve | ( | Vec | initialSolution | ) | [inline, protected, virtual] |
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 SimpleLinearEllipticSolver< DIM, DIM >.
Definition at line 103 of file CellBasedPdeSolver.cpp.
References AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpLinearSystem, and LinearSystem::SetMatrixIsSymmetric().
| void CellBasedPdeSolver< DIM >::ResetInterpolatedQuantities | ( | void | ) | [inline, protected, virtual] |
Overridden ResetInterpolatedQuantities() method.
Reimplemented from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >.
Definition at line 89 of file CellBasedPdeSolver.cpp.
References CellBasedPdeSolver< DIM >::mConstantInUSourceTerm, and CellBasedPdeSolver< DIM >::mLinearInUCoeffInSourceTerm.
double CellBasedPdeSolver< DIM >::mConstantInUSourceTerm [private] |
The constant in u part of the source term, interpolated onto the current point.
Definition at line 54 of file CellBasedPdeSolver.hpp.
Referenced by CellBasedPdeSolver< DIM >::ComputeVectorTerm(), CellBasedPdeSolver< DIM >::IncrementInterpolatedQuantities(), and CellBasedPdeSolver< DIM >::ResetInterpolatedQuantities().
double CellBasedPdeSolver< DIM >::mLinearInUCoeffInSourceTerm [private] |
The linear in u part of the source term, interpolated onto the current point.
Definition at line 57 of file CellBasedPdeSolver.hpp.
Referenced by CellBasedPdeSolver< DIM >::ComputeMatrixTerm(), CellBasedPdeSolver< DIM >::IncrementInterpolatedQuantities(), and CellBasedPdeSolver< DIM >::ResetInterpolatedQuantities().
1.6.2