Chaste Release::3.1
LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > Class Template Reference

#include <LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp>

Inheritance diagram for LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:
Collaboration diagram for LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:

List of all members.

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< doublemInterpolatedOdeStateVariables
boost::shared_ptr
< AbstractIvpOdeSolver
mpOdeSolver
double mSamplingTimeStep
bool mOdeSystemsPresent
out_stream mpVtkMetaFile

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
class LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >

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.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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>() 
)
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~LinearParabolicPdeSystemWithCoupledOdeSystemSolver ( )

Destructor.

Definition at line 396 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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.

Parameters:
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer 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().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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.

Parameters:
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer 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().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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

Parameters:
phiI
pNodepointer 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().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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.

Parameters:
initialSolutionInitial 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.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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.

Parameters:
currentPdeSolutionthe 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().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ResetInterpolatedQuantities ( void  ) [private, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetOutputDirectory ( std::string  outputDirectory)

Set mOutputDirectory.

Parameters:
outputDirectorythe output directory to use

Definition at line 438 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetSamplingTimeStep ( double  samplingTimeStep)

Set mSamplingTimeStep.

Parameters:
samplingTimeStepthe sampling timestep to use

Definition at line 444 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
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.

Parameters:
currentSolutionThe current solution which can be used in setting up the linear system if needed (NULL if there isn't a current solution)
computeMatrixWhether 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.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SolveAndWriteResultsToFile ( )
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteVtkResultsToFile ( Vec  solution,
unsigned  numTimeStepsElapsed 
)

Write the solution to VTK. Called by SolveAndWriteResultsToFile().

Parameters:
solutionthe solution of the coupled PDE/ODE system
numTimeStepsElapsedthe number of timesteps that have elapsed

Definition at line 532 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

References VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddPointData().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteVtkResultsToFile ( ) [private]

Write the current results to mpVtkMetaFile.


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
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.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
std::vector<AbstractOdeSystemForCoupledPdeSystem*> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mOdeSystemsAtNodes [private]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
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.

Referenced by LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::LinearParabolicPdeSystemWithCoupledOdeSystemSolver().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh [private]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
boost::shared_ptr<AbstractIvpOdeSolver> LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpOdeSolver [private]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
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.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
out_stream LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpVtkMetaFile [private]

Meta results file for VTK.

Definition at line 93 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM = ELEMENT_DIM, unsigned PROBLEM_DIM = 1>
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.


The documentation for this class was generated from the following file: