Chaste Release::3.1
|
#include <LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp>
Public Member Functions | |
LinearParabolicPdeSystemWithCoupledOdeSystemSolver (TetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, AbstractLinearParabolicPdeSystemForCoupledOdeSystem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pPdeSystem, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBoundaryConditions, std::vector< AbstractOdeSystemForCoupledPdeSystem * > odeSystemsAtNodes=std::vector< AbstractOdeSystemForCoupledPdeSystem * >(), boost::shared_ptr< AbstractIvpOdeSolver > pOdeSolver=boost::shared_ptr< AbstractIvpOdeSolver >()) | |
~LinearParabolicPdeSystemWithCoupledOdeSystemSolver () | |
void | PrepareForSetupLinearSystem (Vec currentPdeSolution) |
void | SetOutputDirectory (std::string outputDirectory) |
void | SetSamplingTimeStep (double samplingTimeStep) |
void | SolveAndWriteResultsToFile () |
void | WriteVtkResultsToFile (Vec solution, unsigned numTimeStepsElapsed) |
Private Member Functions | |
void | WriteVtkResultsToFile () |
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) |
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) |
void | ResetInterpolatedQuantities () |
void | IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode) |
void | InitialiseForSolve (Vec initialSolution=NULL) |
void | SetupLinearSystem (Vec currentSolution, bool computeMatrix) |
Private Attributes | |
AbstractTetrahedralMesh < ELEMENT_DIM, SPACE_DIM > * | mpMesh |
AbstractLinearParabolicPdeSystemForCoupledOdeSystem < ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | mpPdeSystem |
std::vector < AbstractOdeSystemForCoupledPdeSystem * > | mOdeSystemsAtNodes |
std::vector< double > | mInterpolatedOdeStateVariables |
boost::shared_ptr < AbstractIvpOdeSolver > | mpOdeSolver |
double | mSamplingTimeStep |
bool | mOdeSystemsPresent |
out_stream | mpVtkMetaFile |
A class for solving systems of parabolic PDEs and ODEs, which may be coupled via their source terms:
d/dt (u_i) = div (D(x) grad (u_i)) + f_i (x, u_1, ..., u_p, v_1, ..., v_q), i=1,...,p, d/dt (v_j) = g_j(x, u_1, ..., u_p, v_1, ..., v_q), j=1,...,q.
The solver class is templated over spatial dimension and PDE problem dimension (p).
Definition at line 61 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::LinearParabolicPdeSystemWithCoupledOdeSystemSolver | ( | TetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | pMesh, |
AbstractLinearParabolicPdeSystemForCoupledOdeSystem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | pPdeSystem, | ||
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | pBoundaryConditions, | ||
std::vector< AbstractOdeSystemForCoupledPdeSystem * > | odeSystemsAtNodes = std::vector<AbstractOdeSystemForCoupledPdeSystem*>() , |
||
boost::shared_ptr< AbstractIvpOdeSolver > | pOdeSolver = boost::shared_ptr<AbstractIvpOdeSolver>() |
||
) |
Constructor.
pMesh | pointer to the mesh |
pPdeSystem | pointer to the PDE system |
pBoundaryConditions | pointer to the boundary conditions. |
odeSystemsAtNodes | optional vector of pointers to ODE systems, defined at nodes |
pOdeSolver | optional pointer to an ODE solver (defaults to NULL) |
Definition at line 353 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mOdeSystemsAtNodes, LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mOdeSystemsPresent, AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NORMAL >::mpBoundaryConditions, LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, and LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpOdeSolver.
LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~LinearParabolicPdeSystemWithCoupledOdeSystemSolver | ( | ) |
Destructor.
Definition at line 396 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1)> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_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, PROBLEM_DIM > & | rU, | ||
c_matrix< double, PROBLEM_DIM, SPACE_DIM > & | rGradU, | ||
Element< ELEMENT_DIM, SPACE_DIM > * | pElement | ||
) | [private, virtual] |
The term to be added to the element stiffness matrix.
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 |
Reimplemented from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >.
Definition at line 234 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References PdeSimulationTime::GetPdeTimeStepInverse().
c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_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, PROBLEM_DIM > & | rU, | ||
c_matrix< double, PROBLEM_DIM, SPACE_DIM > & | rGradU, | ||
Element< ELEMENT_DIM, SPACE_DIM > * | pElement | ||
) | [private, virtual] |
The term 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 |
Reimplemented from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >.
Definition at line 266 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References PdeSimulationTime::GetPdeTimeStepInverse().
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::IncrementInterpolatedQuantities | ( | double | phiI, |
const Node< SPACE_DIM > * | pNode | ||
) | [private, virtual] |
Update the member variable mInterpolatedOdeStateVariables by computing the interpolated value of each ODE state variable at each Gauss point
phiI | |
pNode | pointer to a Node |
Reimplemented from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >.
Definition at line 307 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References Node< SPACE_DIM >::GetIndex().
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve | ( | Vec | initialSolution = NULL | ) | [private, virtual] |
Initialise method: sets up the linear system (using the mesh to determine the number of unknowns per row to preallocate) if it is not already set up. Can use an initial solution as PETSc template, or base it on the mesh size.
initialSolution | Initial solution (defaults to NULL) for PETSc to use as a template. |
Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 321 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::PrepareForSetupLinearSystem | ( | Vec | currentPdeSolution | ) | [virtual] |
Overridden PrepareForSetupLinearSystem() method. Pass the current solution to the PDE system to the ODE system and solve it over the next timestep.
currentPdeSolution | the solution to the PDE system at the current time |
Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 408 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References PdeSimulationTime::GetPdeTimeStep(), and PdeSimulationTime::GetTime().
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ResetInterpolatedQuantities | ( | void | ) | [private, virtual] |
Reset the member variable mInterpolatedOdeStateVariables.
Reimplemented from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >.
Definition at line 295 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetOutputDirectory | ( | std::string | outputDirectory | ) |
Set mOutputDirectory.
outputDirectory | the output directory to use |
Definition at line 438 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetSamplingTimeStep | ( | double | samplingTimeStep | ) |
Set mSamplingTimeStep.
samplingTimeStep | the sampling timestep to use |
Definition at line 444 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetupLinearSystem | ( | Vec | currentSolution, |
bool | computeMatrix | ||
) | [private, virtual] |
Completely set up the linear system that has to be solved each timestep.
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 347 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SolveAndWriteResultsToFile | ( | ) |
Solve the coupled PDE/ODE system over the pre-specified time interval, and record results using mSamplingTimeStep.
Definition at line 451 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References TimeStepper::AdvanceOneTimeStep(), PetscTools::Destroy(), DOUBLE_UNSET, EXCEPTION, TimeStepper::GetNextTime(), TimeStepper::GetTime(), TimeStepper::GetTotalTimeStepsTaken(), TimeStepper::IsTimeAtEnd(), and OutputFileHandler::OpenOutputFile().
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteVtkResultsToFile | ( | Vec | solution, |
unsigned | numTimeStepsElapsed | ||
) |
Write the solution to VTK. Called by SolveAndWriteResultsToFile().
solution | the solution of the coupled PDE/ODE system |
numTimeStepsElapsed | the number of timesteps that have elapsed |
Definition at line 532 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddPointData().
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteVtkResultsToFile | ( | ) | [private] |
Write the current results to mpVtkMetaFile.
std::vector<double> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mInterpolatedOdeStateVariables [private] |
The values of the ODE system state variables, interpolated at a quadrature point.
Definition at line 77 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
std::vector<AbstractOdeSystemForCoupledPdeSystem*> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mOdeSystemsAtNodes [private] |
Vector of pointers to ODE systems, defined at nodes.
Definition at line 74 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
bool LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mOdeSystemsPresent [private] |
Whether ODE systems are present (if not, then the system comprises coupled PDEs only).
Definition at line 90 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh [private] |
Pointer to the mesh.
Reimplemented from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.
Definition at line 68 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
boost::shared_ptr<AbstractIvpOdeSolver> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpOdeSolver [private] |
The ODE solver.
Definition at line 80 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
AbstractLinearParabolicPdeSystemForCoupledOdeSystem<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpPdeSystem [private] |
The PDE system to be solved.
Definition at line 71 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
out_stream LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpVtkMetaFile [private] |
Meta results file for VTK.
Definition at line 93 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
double LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mSamplingTimeStep [private] |
A sampling timestep for writing results to file. Set to PdeSimulationTime::GetPdeTimeStep() in the constructor; may be overwritten using the SetSamplingTimeStep() method.
Definition at line 87 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.