Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
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 >:

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, bool clearDirectory=false)
 
void SetSamplingTimeStep (double samplingTimeStep)
 
void SolveAndWriteResultsToFile ()
 
void WriteVtkResultsToFile (Vec solution, unsigned numTimeStepsElapsed)
 
AbstractOdeSystemForCoupledPdeSystemGetOdeSystemAtNode (unsigned index)
 
- Public Member Functions inherited from AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >
 AbstractAssemblerSolverHybrid (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBoundaryConditions)
 
virtual ~AbstractAssemblerSolverHybrid ()
 
void SetupGivenLinearSystem (Vec currentSolution, bool computeMatrix, LinearSystem *pLinearSystem)
 
- Public Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
 AbstractFeVolumeIntegralAssembler (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
virtual ~AbstractFeVolumeIntegralAssembler ()
 
- Public Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
 AbstractFeAssemblerCommon ()
 
void SetCurrentSolution (Vec currentSolution)
 
virtual ~AbstractFeAssemblerCommon ()
 
- Public Member Functions inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
 AbstractFeAssemblerInterface ()
 
void SetMatrixToAssemble (Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
 
void SetVectorToAssemble (Vec &rVecToAssemble, bool zeroVectorBeforeAssembly)
 
void Assemble ()
 
void AssembleMatrix ()
 
void AssembleVector ()
 
virtual ~AbstractFeAssemblerInterface ()
 
- Public Member Functions inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
 AbstractDynamicLinearPdeSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
void SetTimes (double tStart, double tEnd)
 
void SetTimeStep (double dt)
 
void SetInitialCondition (Vec initialCondition)
 
virtual Vec Solve ()
 
void SetMatrixIsNotAssembled ()
 
void SetTimeAdaptivityController (AbstractTimeAdaptivityController *pTimeAdaptivityController)
 
void SetOutputToVtk (bool output)
 
void SetOutputToParallelVtk (bool output)
 
void SetOutputToTxt (bool output)
 
void SetOutputDirectoryAndPrefix (std::string outputDirectory, std::string prefix)
 
void SetPrintingTimestepMultiple (unsigned multiple)
 
- Public Member Functions inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
 AbstractLinearPdeSolver (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
virtual ~AbstractLinearPdeSolver ()
 
virtual void FinaliseLinearSystem (Vec currentSolution)
 
virtual void FollowingSolveLinearSystem (Vec currentSolution)
 
LinearSystemGetLinearSystem ()
 

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< AbstractIvpOdeSolvermpOdeSolver
 
double mSamplingTimeStep
 
bool mOdeSystemsPresent
 
out_stream mpVtkMetaFile
 
bool mClearOutputDirectory
 

Additional Inherited Members

- Protected Types inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
typedef LinearBasisFunction< ELEMENT_DIM > BasisFunction
 
- 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 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 IncrementInterpolatedGradientQuantities (const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode)
 
- Protected Member Functions inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
void InitialiseHdf5Writer ()
 
void WriteOneStep (double time, Vec solution)
 
- Protected Attributes inherited from AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > mNaturalNeumannSurfaceTermAssembler
 
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * mpBoundaryConditions
 
- Protected Attributes inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, INTERPOLATION_LEVEL >
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
 
GaussianQuadratureRule< ELEMENT_DIM > * mpQuadRule
 
- Protected Attributes inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >
ReplicatableVector mCurrentSolutionOrGuessReplicated
 
- Protected Attributes inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
Vec mVectorToAssemble
 
Mat mMatrixToAssemble
 
bool mAssembleMatrix
 
bool mAssembleVector
 
bool mZeroMatrixBeforeAssembly
 
bool mZeroVectorBeforeAssembly
 
PetscInt mOwnershipRangeLo
 
PetscInt mOwnershipRangeHi
 
- Protected Attributes inherited from AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
double mTstart
 
double mTend
 
bool mTimesSet
 
Vec mInitialCondition
 
bool mMatrixIsAssembled
 
bool mMatrixIsConstant
 
double mIdealTimeStep
 
double mLastWorkingTimeStep
 
AbstractTimeAdaptivityControllermpTimeAdaptivityController
 
bool mOutputToVtk
 
bool mOutputToParallelVtk
 
bool mOutputToTxt
 
std::string mOutputDirectory
 
std::string mFilenamePrefix
 
unsigned mPrintingTimestepMultiple
 
Hdf5DataWritermpHdf5Writer
 
std::vector< int > mVariableColumnIds
 
- Protected Attributes inherited from AbstractLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
LinearSystemmpLinearSystem
 
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
 

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

◆ LinearParabolicPdeSystemWithCoupledOdeSystemSolver()

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>() 
)

◆ ~LinearParabolicPdeSystemWithCoupledOdeSystemSolver()

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

Member Function Documentation

◆ ComputeMatrixTerm()

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 
)
privatevirtual
Returns
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 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().

◆ ComputeVectorTerm()

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 
)
privatevirtual
Returns
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 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().

◆ GetOdeSystemAtNode()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractOdeSystemForCoupledPdeSystem * LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetOdeSystemAtNode ( unsigned  index)

Get a pointer to the ODE system defined at a given node.

Parameters
indexthe global index of a node in the mpMesh
Returns
mOdeSystemsAtNodes[index]

Definition at line 631 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

◆ IncrementInterpolatedQuantities()

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 
)
privatevirtual

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 325 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

References Node< SPACE_DIM >::GetIndex().

◆ InitialiseForSolve()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseForSolve ( Vec  initialSolution = NULL)
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.

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 339 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::CalculateMaximumContainingElementsPerProcess().

◆ PrepareForSetupLinearSystem()

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 427 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

References PdeSimulationTime::GetNextTime(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), PdeSimulationTime::GetPdeTimeStep(), and PdeSimulationTime::GetTime().

◆ ResetInterpolatedQuantities()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ResetInterpolatedQuantities ( void  )
privatevirtual

◆ SetOutputDirectory()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetOutputDirectory ( std::string  outputDirectory,
bool  clearDirectory = false 
)

Set mOutputDirectory.

Parameters
outputDirectorythe output directory to use
clearDirectorywhether 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.

◆ SetSamplingTimeStep()

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 465 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

◆ SetupLinearSystem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetupLinearSystem ( Vec  currentSolution,
bool  computeMatrix 
)
privatevirtual

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 365 of file LinearParabolicPdeSystemWithCoupledOdeSystemSolver.hpp.

◆ SolveAndWriteResultsToFile()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void LinearParabolicPdeSystemWithCoupledOdeSystemSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SolveAndWriteResultsToFile ( )

◆ WriteVtkResultsToFile() [1/2]

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.

◆ WriteVtkResultsToFile() [2/2]

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 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().

Member Data Documentation

◆ mClearOutputDirectory

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

◆ mInterpolatedOdeStateVariables

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.

◆ mOdeSystemsAtNodes

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

◆ mOdeSystemsPresent

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().

◆ mpMesh

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

◆ mpOdeSolver

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

◆ mpPdeSystem

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.

◆ mpVtkMetaFile

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.

◆ mSamplingTimeStep

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: