#include <CardiacElectroMechanicsProblem.hpp>
Public Member Functions | |
CardiacElectroMechanicsProblem (ContractionModel contractionModel, TetrahedralMesh< DIM, DIM > *pElectricsMesh, QuadraticMesh< DIM > *pMechanicsMesh, std::vector< unsigned > fixedMechanicsNodes, AbstractCardiacCellFactory< DIM > *pCellFactory, double endTime, double electricsPdeTimeStep, double mechanicsSolveTimestep, double contractionModelOdeTimeStep, std::string outputDirectory) | |
virtual | ~CardiacElectroMechanicsProblem () |
void | Initialise () |
void | Solve () |
double | Max (std::vector< double > &vec) |
void | SetNoElectricsOutput () |
void | SetWatchedPosition (c_vector< double, DIM > watchedLocation) |
void | SetVariableFibreSheetDirectionsFile (std::string orthoFile, bool definedPerQuadPoint) |
std::vector< c_vector< double, DIM > > & | rGetDeformedPosition () |
void | UseMechanoElectricFeedback () |
c_matrix< double, DIM, DIM > & | rGetModifiedConductivityTensor (unsigned elementIndex, const c_matrix< double, DIM, DIM > &rOriginalConductivity) |
virtual void | PrepareForSolve () |
virtual void | OnEndOfTimeStep (unsigned counter) |
AbstractNonlinearElasticitySolver < DIM > * | GetCardiacMechanicsSolver () |
void | SetMaterialLaw (AbstractMaterialLaw< DIM > *pMaterialLaw) |
Protected Member Functions | |
void | DetermineWatchedNodes () |
void | WriteWatchedLocationData (double time, Vec voltage) |
Protected Attributes | |
ContractionModel | mContractionModel |
MonodomainProblem< DIM > * | mpMonodomainProblem |
AbstractCardiacMechanicsSolverInterface < DIM > * | mpCardiacMechSolver |
AbstractNonlinearElasticitySolver < DIM > * | mpMechanicsSolver |
double | mEndTime |
double | mElectricsTimeStep |
double | mMechanicsTimeStep |
unsigned | mNumElecTimestepsPerMechTimestep |
double | mContractionModelOdeTimeStep |
TetrahedralMesh< DIM, DIM > * | mpElectricsMesh |
QuadraticMesh< DIM > * | mpMechanicsMesh |
FineCoarseMeshPair< DIM > * | mpMeshPair |
std::string | mOutputDirectory |
std::string | mDeformationOutputDirectory |
bool | mWriteOutput |
bool | mNoElectricsOutput |
bool | mIsWatchedLocation |
c_vector< double, DIM > | mWatchedLocation |
unsigned | mWatchedElectricsNodeIndex |
unsigned | mWatchedMechanicsNodeIndex |
out_stream | mpWatchedLocationFile |
std::vector< unsigned > | mFixedNodes |
std::string | mFibreSheetDirectionsFile |
bool | mFibreSheetDirectionsDefinedPerQuadraturePoint |
std::vector< double > | mStretchesForEachMechanicsElement |
std::vector< c_matrix< double, DIM, DIM > > | mDeformationGradientsForEachMechanicsElement |
std::pair< unsigned, c_matrix < double, DIM, DIM > > | mLastModifiedConductivity |
bool | mConductivityAffectedByDeformationMef |
bool | mCellModelsAffectedByDeformationMef |
AbstractMaterialLaw< DIM > * | mpMaterialLaw |
bool | mAllocatedMaterialLawMemory |
Static Protected Attributes | |
static const int | WRITE_EVERY_NTH_TIME = 1 |
Friends | |
class | TestCardiacElectroMechanicsProblem |
For solving full electro-mechanical problems.
Solves a monodomain problem (diffusion plus cell models) on a (fine) electrics mesh, and a mechanics problem (finite elasticity plus contraction model) on a coarse mesh. An implicit scheme (Jon Whiteley's algorithm) be be used.
For solving problems on regular grids use CardiacElectroMechProbRegularGeom
The implicit algorithm:
Store the position in the electrics mesh of each quad point in the mechanics mesh For every time: Solve the monodomain problem (ie integrate ODEs, solve PDE) Get intracellular [Ca] at each electrics node and interpolate on each mechanics quad point Set [Ca] on each contraction model (one for each point) Solve static finite elasticity problem implicity
Definition at line 73 of file CardiacElectroMechanicsProblem.hpp.
CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem | ( | ContractionModel | contractionModel, | |
TetrahedralMesh< DIM, DIM > * | pElectricsMesh, | |||
QuadraticMesh< DIM > * | pMechanicsMesh, | |||
std::vector< unsigned > | fixedMechanicsNodes, | |||
AbstractCardiacCellFactory< DIM > * | pCellFactory, | |||
double | endTime, | |||
double | electricsPdeTimeStep, | |||
double | mechanicsSolveTimestep, | |||
double | contractionModelOdeTimeStep, | |||
std::string | outputDirectory = "" | |||
) | [inline] |
Constructor.
contractionModel | contraction model (see the enum "ContractionModel" for the options). | |
pElectricsMesh | Mesh on which to solve electrics (Monodomain) | |
pMechanicsMesh | Mesh (2nd order) on which to solve mechanics | |
fixedMechanicsNodes | Indices of those nodes which a pinned in space | |
pCellFactory | factory to use to create cells | |
endTime | the end time to use | |
electricsPdeTimeStep | timestep used in solving for the electrical activity | |
mechanicsSolveTimestep | how often the mechanics is solved for (should be a multiple of electrics PDE timestep) | |
contractionModelOdeTimeStep | Step size for contraction model (of active tension in cardiac cells) being used. | |
outputDirectory | the output directory |
Definition at line 226 of file CardiacElectroMechanicsProblem.cpp.
References GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), GenericEventHandler< 13, HeartEventHandler >::Disable(), EXCEPTION, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), LogFile::Instance(), HeartConfig::Instance(), CardiacElectroMechanicsProblem< DIM >::mAllocatedMaterialLawMemory, CardiacElectroMechanicsProblem< DIM >::mCellModelsAffectedByDeformationMef, CardiacElectroMechanicsProblem< DIM >::mConductivityAffectedByDeformationMef, CardiacElectroMechanicsProblem< DIM >::mContractionModel, CardiacElectroMechanicsProblem< DIM >::mContractionModelOdeTimeStep, CardiacElectroMechanicsProblem< DIM >::mDeformationOutputDirectory, CardiacElectroMechanicsProblem< DIM >::mElectricsTimeStep, CardiacElectroMechanicsProblem< DIM >::mEndTime, CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsFile, CardiacElectroMechanicsProblem< DIM >::mFixedNodes, CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mLastModifiedConductivity, CardiacElectroMechanicsProblem< DIM >::mMechanicsTimeStep, CardiacElectroMechanicsProblem< DIM >::mNoElectricsOutput, CardiacElectroMechanicsProblem< DIM >::mNumElecTimestepsPerMechTimestep, CardiacElectroMechanicsProblem< DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM >::mpElectricsMesh, CardiacElectroMechanicsProblem< DIM >::mpMaterialLaw, CardiacElectroMechanicsProblem< DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver, CardiacElectroMechanicsProblem< DIM >::mpMonodomainProblem, CardiacElectroMechanicsProblem< DIM >::mWatchedElectricsNodeIndex, CardiacElectroMechanicsProblem< DIM >::mWatchedMechanicsNodeIndex, CardiacElectroMechanicsProblem< DIM >::mWriteOutput, GenericEventHandler< 7, MechanicsEventHandler >::Reset(), LogFile::Set(), HeartConfig::SetOutputDirectory(), HeartConfig::SetOutputFilenamePrefix(), HeartConfig::SetSimulationDuration(), UNSIGNED_UNSET, and LogFile::WriteHeader().
CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem | ( | ) | [inline, virtual] |
Delete allocated memory and close the watched location file
NOTE if SetWatchedLocation but not Initialise has been called, mpWatchedLocationFile will be uninitialised and using it will cause a seg fault. Hence the mpMechanicsMesh!=NULL it is true if Initialise has been called.
NOTE if SetWatchedLocation but not Initialise has been called, mpWatchedLocationFile will be uninitialised and using it will cause a seg fault. Hence the mpMechanicsMesh!=NULL it is true if Initialise has been called.
Definition at line 343 of file CardiacElectroMechanicsProblem.cpp.
References LogFile::Close(), CardiacElectroMechanicsProblem< DIM >::mAllocatedMaterialLawMemory, CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM >::mpMaterialLaw, CardiacElectroMechanicsProblem< DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM >::mpMeshPair, CardiacElectroMechanicsProblem< DIM >::mpMonodomainProblem, and CardiacElectroMechanicsProblem< DIM >::mpWatchedLocationFile.
void CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes | ( | ) | [inline, protected] |
Determine which node is closest to the watched location
Definition at line 55 of file CardiacElectroMechanicsProblem.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM >::mpElectricsMesh, CardiacElectroMechanicsProblem< DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM >::mpWatchedLocationFile, CardiacElectroMechanicsProblem< DIM >::mWatchedElectricsNodeIndex, CardiacElectroMechanicsProblem< DIM >::mWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mWatchedMechanicsNodeIndex, NEVER_REACHED, OutputFileHandler::OpenOutputFile(), and UNSIGNED_UNSET.
Referenced by CardiacElectroMechanicsProblem< DIM >::Initialise().
void CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData | ( | double | time, | |
Vec | voltage | |||
) | [inline, protected] |
Write info (x, y, [z], V) for the watched node.
time | Time-step now, to write out | |
voltage | Vm vector (this is Monodomain) |
Definition at line 139 of file CardiacElectroMechanicsProblem.cpp.
References CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver, CardiacElectroMechanicsProblem< DIM >::mpWatchedLocationFile, CardiacElectroMechanicsProblem< DIM >::mWatchedElectricsNodeIndex, and CardiacElectroMechanicsProblem< DIM >::mWatchedMechanicsNodeIndex.
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
void CardiacElectroMechanicsProblem< DIM >::Initialise | ( | ) | [inline] |
Initialise the class. Initialises the MonodomainProblem and sets up the electrics mesh to mechanics mesh data.
Definition at line 370 of file CardiacElectroMechanicsProblem.cpp.
References CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), EXCEPTION, MonodomainProblem< ELEMENT_DIM, SPACE_DIM >::GetMonodomainTissue(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Initialise(), HeartConfig::Instance(), CardiacElectroMechanicsProblem< DIM >::mAllocatedMaterialLawMemory, CardiacElectroMechanicsProblem< DIM >::mCellModelsAffectedByDeformationMef, CardiacElectroMechanicsProblem< DIM >::mConductivityAffectedByDeformationMef, CardiacElectroMechanicsProblem< DIM >::mContractionModel, CardiacElectroMechanicsProblem< DIM >::mDeformationGradientsForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM >::mDeformationOutputDirectory, CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsDefinedPerQuadraturePoint, CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsFile, CardiacElectroMechanicsProblem< DIM >::mFixedNodes, CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM >::mpElectricsMesh, CardiacElectroMechanicsProblem< DIM >::mpMaterialLaw, CardiacElectroMechanicsProblem< DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver, CardiacElectroMechanicsProblem< DIM >::mpMeshPair, CardiacElectroMechanicsProblem< DIM >::mpMonodomainProblem, CardiacElectroMechanicsProblem< DIM >::mStretchesForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM >::mWriteOutput, AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM >::SetConductivityModifier(), HeartConfig::SetIntracellularConductivities(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetMesh(), and AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
void CardiacElectroMechanicsProblem< DIM >::Solve | ( | ) | [inline] |
Solve the electromechanics problem
Definition at line 475 of file CardiacElectroMechanicsProblem.cpp.
References Hdf5DataWriter::AdvanceAlongUnlimitedDimension(), TimeStepper::AdvanceOneTimeStep(), PetscTools::Barrier(), GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), Hdf5DataWriter::Close(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CreateInitialCondition(), MonodomainProblem< ELEMENT_DIM, SPACE_DIM >::CreateSolver(), GenericEventHandler< 7, MechanicsEventHandler >::EndEvent(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetDistributedVectorFactory(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement(), DistributedVectorFactory::GetHigh(), DistributedVectorFactory::GetLow(), TimeStepper::GetNextTime(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), HeartConfig::GetOutputDirectory(), TimeStepper::GetTime(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetTissue(), CardiacElectroMechanicsProblem< DIM >::Initialise(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseWriter(), LogFile::Instance(), HeartConfig::Instance(), TimeStepper::IsTimeAtEnd(), CardiacElectroMechanicsProblem< DIM >::Max(), CardiacElectroMechanicsProblem< DIM >::mCellModelsAffectedByDeformationMef, CardiacElectroMechanicsProblem< DIM >::mConductivityAffectedByDeformationMef, CardiacElectroMechanicsProblem< DIM >::mContractionModelOdeTimeStep, CardiacElectroMechanicsProblem< DIM >::mDeformationGradientsForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM >::mElectricsTimeStep, CardiacElectroMechanicsProblem< DIM >::mEndTime, CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM >::mMechanicsTimeStep, CardiacElectroMechanicsProblem< DIM >::mNoElectricsOutput, CardiacElectroMechanicsProblem< DIM >::mNumElecTimestepsPerMechTimestep, CardiacElectroMechanicsProblem< DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM >::mpElectricsMesh, CardiacElectroMechanicsProblem< DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver, CardiacElectroMechanicsProblem< DIM >::mpMeshPair, CardiacElectroMechanicsProblem< DIM >::mpMonodomainProblem, AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpWriter, CardiacElectroMechanicsProblem< DIM >::mStretchesForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM >::mWriteOutput, CardiacElectroMechanicsProblem< DIM >::OnEndOfTimeStep(), CardiacElectroMechanicsProblem< DIM >::PrepareForSolve(), HeartConfig::Reset(), CardiacElectroMechanicsProblem< DIM >::rGetDeformedPosition(), CmguiMeshWriter< DIM, DIM >::SetAdditionalFieldNames(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetBoundaryConditionsContainer(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetInitialCondition(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetMatrixIsNotAssembled(), HeartConfig::SetOutputDirectory(), HeartConfig::SetPrintingTimeStep(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetTimes(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetTimeStep(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve(), CardiacElectroMechanicsProblem< DIM >::WRITE_EVERY_NTH_TIME, CmguiDeformedSolutionsWriter< DIM >::WriteCmguiScript(), CmguiDeformedSolutionsWriter< DIM >::WriteDeformationPositions(), LogFile::WriteElapsedTime(), CmguiDeformedSolutionsWriter< DIM >::WriteInitialMesh(), MonodomainProblem< ELEMENT_DIM, SPACE_DIM >::WriteOneStep(), and CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData().
double CardiacElectroMechanicsProblem< DIM >::Max | ( | std::vector< double > & | vec | ) | [inline] |
Short helper function - the max of a std::vector
vec | a vector of doubles |
Definition at line 791 of file CardiacElectroMechanicsProblem.cpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
void CardiacElectroMechanicsProblem< DIM >::SetNoElectricsOutput | ( | ) | [inline] |
Call to not write out voltages
Definition at line 802 of file CardiacElectroMechanicsProblem.cpp.
References CardiacElectroMechanicsProblem< DIM >::mNoElectricsOutput.
void CardiacElectroMechanicsProblem< DIM >::SetWatchedPosition | ( | c_vector< double, DIM > | watchedLocation | ) | [inline] |
Set a location to be watched - for which lots of output is given. Should correspond to nodes in both meshes.
The watched file will have rows that look like: time x_pos y_pos [z_pos] voltage Ca_i_conc.
NOTE: for the Calcium - assumes LUO_RUDY IS USED
watchedLocation | location (x,y,z) in space. Watched node is the closest to this point. |
Definition at line 808 of file CardiacElectroMechanicsProblem.cpp.
References CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation, and CardiacElectroMechanicsProblem< DIM >::mWatchedLocation.
void CardiacElectroMechanicsProblem< DIM >::SetVariableFibreSheetDirectionsFile | ( | std::string | orthoFile, | |
bool | definedPerQuadPoint | |||
) | [inline] |
Set a variable fibre-sheet-normal direction (matrices), from file. If the second parameter is false, there should be one fibre-sheet definition for each element; otherwise there should be one fibre-sheet definition for each *quadrature point* in the mesh. In the first case, the file should be a .ortho file (ie each line has the fibre dir, sheet dir, normal dir for that element), in the second it should have .orthoquad as the format.
orthoFile | the file containing the fibre/sheet directions | |
definedPerQuadPoint | whether the fibre-sheet definitions are for each quadrature point in the mesh (if not, one for each element is assumed). |
Definition at line 815 of file CardiacElectroMechanicsProblem.cpp.
References CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsDefinedPerQuadraturePoint, and CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsFile.
std::vector< c_vector< double, DIM > > & CardiacElectroMechanicsProblem< DIM >::rGetDeformedPosition | ( | ) | [inline] |
Definition at line 822 of file CardiacElectroMechanicsProblem.cpp.
References CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver.
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
void CardiacElectroMechanicsProblem< DIM >::UseMechanoElectricFeedback | ( | ) | [inline] |
By default (at the moment), the deformation does not affect the electrophysiology in any way. Call this to allow it to, then (i) The stretch will be passed back to the cell models for use stretch-activated channels etc (ii) The deformation to alter the conductivity
EMTODO - check these Two things to note: (i) this can't be called if fibre-sheet directions have been defined from file for each quadrature point (as opposed to each mechanics element) - this is because if the stretch is to be passed back to the electric mesh nodes, the fibre direction has to be defined at those nodes (ii) currently the set-up stage (computing mechanics mesh elements and weights for electrics mesh nodes) is inefficiently implemented - setup will be very slow for big meshes
Definition at line 284 of file CardiacElectroMechanicsProblem.hpp.
References CardiacElectroMechanicsProblem< DIM >::mCellModelsAffectedByDeformationMef, and CardiacElectroMechanicsProblem< DIM >::mConductivityAffectedByDeformationMef.
c_matrix< double, DIM, DIM > & CardiacElectroMechanicsProblem< DIM >::rGetModifiedConductivityTensor | ( | unsigned | elementIndex, | |
const c_matrix< double, DIM, DIM > & | rOriginalConductivity | |||
) | [inline] |
The implementation of the pure method defined in the base class AbstractConductivityModifier. The tissue class will call this method.
elementIndex | Index of current element | |
rOriginalConductivity | Reference to the original (for example, undeformed) conductivity tensor |
Definition at line 164 of file CardiacElectroMechanicsProblem.cpp.
References Inverse(), CardiacElectroMechanicsProblem< DIM >::mDeformationGradientsForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM >::mLastModifiedConductivity, and CardiacElectroMechanicsProblem< DIM >::mpMeshPair.
virtual void CardiacElectroMechanicsProblem< DIM >::PrepareForSolve | ( | ) | [inline, virtual] |
Called in Solve() before the time loop
Definition at line 305 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
virtual void CardiacElectroMechanicsProblem< DIM >::OnEndOfTimeStep | ( | unsigned | counter | ) | [inline, virtual] |
Called in Solve() at the end of every time step
counter | time step |
Definition at line 311 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
AbstractNonlinearElasticitySolver<DIM>* CardiacElectroMechanicsProblem< DIM >::GetCardiacMechanicsSolver | ( | ) | [inline] |
Get the mechanics solver. Needs to be called after Initialise().
///
Definition at line 318 of file CardiacElectroMechanicsProblem.hpp.
References CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver.
void CardiacElectroMechanicsProblem< DIM >::SetMaterialLaw | ( | AbstractMaterialLaw< DIM > * | pMaterialLaw | ) | [inline] |
Set the material law. If this isn't called the default material law will be used. Only call before Initialise() and before Solve()
pMaterialLaw | the material law |
Definition at line 328 of file CardiacElectroMechanicsProblem.hpp.
References CardiacElectroMechanicsProblem< DIM >::mAllocatedMaterialLawMemory, and CardiacElectroMechanicsProblem< DIM >::mpMaterialLaw.
ContractionModel CardiacElectroMechanicsProblem< DIM >::mContractionModel [protected] |
Contraction model (from enumeration)
Definition at line 82 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Initialise().
MonodomainProblem<DIM>* CardiacElectroMechanicsProblem< DIM >::mpMonodomainProblem [protected] |
The cardiac problem class
Definition at line 84 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
AbstractCardiacMechanicsSolverInterface<DIM>* CardiacElectroMechanicsProblem< DIM >::mpCardiacMechSolver [protected] |
Definition at line 88 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
AbstractNonlinearElasticitySolver<DIM>* CardiacElectroMechanicsProblem< DIM >::mpMechanicsSolver [protected] |
The mechanics solver - a pointer to the part that sees the solid mechanics solver
Definition at line 90 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::GetCardiacMechanicsSolver(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::rGetDeformedPosition(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData().
double CardiacElectroMechanicsProblem< DIM >::mEndTime [protected] |
End time. The start time is assumed to be 0.0
Definition at line 93 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Solve().
double CardiacElectroMechanicsProblem< DIM >::mElectricsTimeStep [protected] |
The electrics timestep.
Definition at line 95 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Solve().
double CardiacElectroMechanicsProblem< DIM >::mMechanicsTimeStep [protected] |
The mechanics timestep. Needs to be a multiple of the electrics timestep
Definition at line 97 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Solve().
unsigned CardiacElectroMechanicsProblem< DIM >::mNumElecTimestepsPerMechTimestep [protected] |
The number of electrics timesteps per mechanics timestep
Definition at line 99 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Solve().
double CardiacElectroMechanicsProblem< DIM >::mContractionModelOdeTimeStep [protected] |
Timestep to use when solving contraction models
Definition at line 101 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Solve().
TetrahedralMesh<DIM,DIM>* CardiacElectroMechanicsProblem< DIM >::mpElectricsMesh [protected] |
The mesh for the electrics
Definition at line 104 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom(), CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), CardiacElectroMechanicsProblem< DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM >::Solve().
QuadraticMesh<DIM>* CardiacElectroMechanicsProblem< DIM >::mpMechanicsMesh [protected] |
The mesh for the mechanics
Definition at line 106 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom(), CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
FineCoarseMeshPair<DIM>* CardiacElectroMechanicsProblem< DIM >::mpMeshPair [protected] |
Class wrapping both meshes, useful for transferring information
Definition at line 109 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::rGetModifiedConductivityTensor(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
std::string CardiacElectroMechanicsProblem< DIM >::mOutputDirectory [protected] |
Output directory, relative to TEST_OUTPUT
Definition at line 112 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), CardiacElectroMechanicsProblem< DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM >::Solve().
std::string CardiacElectroMechanicsProblem< DIM >::mDeformationOutputDirectory [protected] |
Deformation output-sub-directory
Definition at line 114 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::Initialise().
bool CardiacElectroMechanicsProblem< DIM >::mWriteOutput [protected] |
Whether to write any output
Definition at line 116 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM >::Solve().
bool CardiacElectroMechanicsProblem< DIM >::mNoElectricsOutput [protected] |
Whether to not write out voltages
Definition at line 118 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::SetNoElectricsOutput(), and CardiacElectroMechanicsProblem< DIM >::Solve().
const int CardiacElectroMechanicsProblem< DIM >::WRITE_EVERY_NTH_TIME = 1 [static, protected] |
when to write output
Definition at line 121 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Solve().
bool CardiacElectroMechanicsProblem< DIM >::mIsWatchedLocation [protected] |
Whether any location has been set to be watched (lots of output for that location
Definition at line 124 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::SetWatchedPosition(), CardiacElectroMechanicsProblem< DIM >::Solve(), CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
c_vector<double,DIM> CardiacElectroMechanicsProblem< DIM >::mWatchedLocation [protected] |
The watched location if there is one
Definition at line 126 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), and CardiacElectroMechanicsProblem< DIM >::SetWatchedPosition().
unsigned CardiacElectroMechanicsProblem< DIM >::mWatchedElectricsNodeIndex [protected] |
The node in the electrics mesh corresponding to the watched location
Definition at line 128 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), and CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData().
unsigned CardiacElectroMechanicsProblem< DIM >::mWatchedMechanicsNodeIndex [protected] |
The node in the mechanics mesh corresponding to the watched location
Definition at line 130 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), and CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData().
out_stream CardiacElectroMechanicsProblem< DIM >::mpWatchedLocationFile [protected] |
File where watched location info is written
Definition at line 132 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::DetermineWatchedNodes(), CardiacElectroMechanicsProblem< DIM >::WriteWatchedLocationData(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
std::vector<unsigned> CardiacElectroMechanicsProblem< DIM >::mFixedNodes [protected] |
Nodes for which the deformation is fixed to zero
Definition at line 135 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom(), and CardiacElectroMechanicsProblem< DIM >::Initialise().
std::string CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsFile [protected] |
.ortho/.orthoquad file from which to read element-wise/quadpoint-wise fibre-sheet-normal-directions
Definition at line 137 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM >::SetVariableFibreSheetDirectionsFile().
bool CardiacElectroMechanicsProblem< DIM >::mFibreSheetDirectionsDefinedPerQuadraturePoint [protected] |
Whether the mFibreSheetDirectionsFile file gives the fibre info at each element, or each quadrature point
Definition at line 139 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM >::SetVariableFibreSheetDirectionsFile().
std::vector<double> CardiacElectroMechanicsProblem< DIM >::mStretchesForEachMechanicsElement [protected] |
A vector of stretches (in the fibre direction), one for each element in the mechanics mesh
Definition at line 142 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM >::Solve().
std::vector<c_matrix<double,DIM,DIM> > CardiacElectroMechanicsProblem< DIM >::mDeformationGradientsForEachMechanicsElement [protected] |
A vector of deformation gradients (each entry a matrix), one for each element in the mechanics mesh
Definition at line 144 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::rGetModifiedConductivityTensor(), and CardiacElectroMechanicsProblem< DIM >::Solve().
std::pair<unsigned, c_matrix<double,DIM,DIM> > CardiacElectroMechanicsProblem< DIM >::mLastModifiedConductivity [protected] |
A pair of (element_index, deformed_conductivity) for the last element on which the deformed conductivity sigma_def = F^{-1} sigma_undef F^{-T} has been computed. Used in rGetModifiedConductivityTensor().
Definition at line 149 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM >::rGetModifiedConductivityTensor().
bool CardiacElectroMechanicsProblem< DIM >::mConductivityAffectedByDeformationMef [protected] |
Whether the deformation is always to alter the conductivities (ie one part of MEF)
Definition at line 152 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::UseMechanoElectricFeedback().
bool CardiacElectroMechanicsProblem< DIM >::mCellModelsAffectedByDeformationMef [protected] |
Whether the deformation is always to affect the cell models (for example, for use in stretch-activated channels) (ie one part of MEF)
Definition at line 157 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM >::UseMechanoElectricFeedback().
AbstractMaterialLaw<DIM>* CardiacElectroMechanicsProblem< DIM >::mpMaterialLaw [protected] |
The material law to be used. Defaults to NashHunterPoleZero if an incompressible problem is being solved
Definition at line 160 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::SetMaterialLaw(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().
bool CardiacElectroMechanicsProblem< DIM >::mAllocatedMaterialLawMemory [protected] |
Whether a material law was passed in or the default used
Definition at line 163 of file CardiacElectroMechanicsProblem.hpp.
Referenced by CardiacElectroMechanicsProblem< DIM >::CardiacElectroMechanicsProblem(), CardiacElectroMechanicsProblem< DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM >::SetMaterialLaw(), and CardiacElectroMechanicsProblem< DIM >::~CardiacElectroMechanicsProblem().