Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
|
#include <LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp>
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 |
bool | mClearOutputDirectory |
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 371 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, INTERPOLATION_LEVEL >::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. If an ODE system is present, the pointers to the ODE system objects are deleted here.
Definition at line 415 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
privatevirtual |
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 251 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References AbstractLinearParabolicPdeSystemForCoupledOdeSystem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeDiffusionTerm(), AbstractLinearParabolicPdeSystemForCoupledOdeSystem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeDuDtCoefficientFunction(), and PdeSimulationTime::GetPdeTimeStepInverse().
|
privatevirtual |
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 283 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References AbstractLinearParabolicPdeSystemForCoupledOdeSystem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeDuDtCoefficientFunction(), AbstractLinearParabolicPdeSystemForCoupledOdeSystem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeSourceTerm(), and PdeSimulationTime::GetPdeTimeStepInverse().
AbstractOdeSystemForCoupledPdeSystem * LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetOdeSystemAtNode | ( | unsigned | index | ) |
Get a pointer to the ODE system defined at a given node.
index | the global index of a node in the mpMesh |
Definition at line 631 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
privatevirtual |
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 325 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References Node< SPACE_DIM >::GetIndex().
|
privatevirtual |
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 339 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMaximumContainingElementsPerProcess().
|
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 427 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References PdeSimulationTime::GetNextTime(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), PdeSimulationTime::GetPdeTimeStep(), and PdeSimulationTime::GetTime().
|
privatevirtual |
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 313 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetOutputDirectory | ( | std::string | outputDirectory, |
bool | clearDirectory = false |
||
) |
Set mOutputDirectory.
outputDirectory | the output directory to use |
clearDirectory | whether to clear outputDirectory or not. Note that the actual clearing happens when you call SolveAndWriteResultsToFile(). False by default. |
Definition at line 458 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 465 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
privatevirtual |
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 365 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 472 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References TimeStepper::AdvanceOneTimeStep(), PetscTools::Destroy(), DOUBLE_UNSET, EXCEPTION, TimeStepper::GetNextTime(), TimeStepper::GetTime(), TimeStepper::GetTotalTimeStepsTaken(), TimeStepper::IsTimeAtEnd(), and OutputFileHandler::OpenOutputFile().
|
private |
Write the current results to mpVtkMetaFile.
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 553 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
References VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddPointData(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), and VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
|
private |
Whether the output directory should be cleared before solve or not. False by default. Can be changed when setting the output directory
Definition at line 99 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
The values of the ODE system state variables, interpolated at a quadrature point.
Definition at line 77 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
Vector of pointers to ODE systems, defined at nodes.
Definition at line 74 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
Whether ODE systems are present (if not, then the system comprises coupled PDEs only).
Definition at line 90 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
Pointer to the mesh.
Definition at line 68 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
The ODE solver.
Definition at line 80 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
The PDE system to be solved.
Definition at line 71 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
private |
Meta results file for VTK.
Definition at line 93 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.
|
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.