Chaste Release::3.1
ExtendedBidomainProblem< DIM > Class Template Reference

#include <ExtendedBidomainProblem.hpp>

Inheritance diagram for ExtendedBidomainProblem< DIM >:
Collaboration diagram for ExtendedBidomainProblem< DIM >:

List of all members.

Public Member Functions

 ExtendedBidomainProblem (AbstractCardiacCellFactory< DIM > *pCellFactory, AbstractCardiacCellFactory< DIM > *pSecondCellFactory, bool hasBath=false)
 ExtendedBidomainProblem ()
 ~ExtendedBidomainProblem ()
void SetFixedExtracellularPotentialNodes (std::vector< unsigned > nodes)
void SetIntracellularConductivitiesForSecondCell (c_vector< double, DIM > conductivities)
void SetExtendedBidomainParameters (double Am1, double Am2, double AmGap, double Cm1, double Cm2, double Ggap)
void SetGgapHeterogeneities (std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &rGgapHeterogeneityRegions, std::vector< double > &rGgapValues)
void SetExtracellularStimulusFactory (AbstractStimulusFactory< DIM > *pFactory)
void SetNodeForAverageOfPhiZeroed (unsigned node)
ExtendedBidomainTissue< DIM > * GetExtendedBidomainTissue ()
void WriteInfo (double time)
virtual void DefineWriterColumns (bool extending)
virtual void WriteOneStep (double time, Vec voltageVec)
void PreSolveChecks ()
bool GetHasBath ()
void SetHasBath (bool hasBath)

Protected Member Functions

Vec CreateInitialCondition ()
void AnalyseMeshForBath ()
void ProcessExtracellularStimulus ()
virtual AbstractCardiacTissue
< DIM > * 
CreateCardiacTissue ()
virtual
AbstractDynamicLinearPdeSolver
< DIM, DIM, 3 > * 
CreateSolver ()

Protected Attributes

AbstractCardiacCellFactory
< DIM, DIM > * 
mpSecondCellFactory
ExtendedBidomainTissue< DIM > * mpExtendedBidomainTissue
std::vector< unsignedmFixedExtracellularPotentialNodes
c_vector< double, DIM > mIntracellularConductivitiesSecondCell
unsigned mVoltageColumnId_Vm1
unsigned mVoltageColumnId_Vm2
unsigned mVoltageColumnId_Phie
std::vector< signed int > mVariablesIDs
bool mUserSpecifiedSecondCellConductivities
bool mUserHasSetBidomainValuesExplicitly
double mAmFirstCell
double mAmSecondCell
double mAmGap
double mCmFirstCell
double mCmSecondCell
double mGGap
std::vector< boost::shared_ptr
< AbstractChasteRegion< DIM > > > 
mGgapHeterogeneityRegions
std::vector< doublemGgapHeterogenousValues
AbstractStimulusFactory< DIM > * mpExtracellularStimulusFactory
int mRowForAverageOfPhiZeroed
unsigned mApplyAveragePhieZeroConstraintAfterSolving
bool mUserSuppliedExtracellularStimulus
bool mHasBath
AbstractExtendedBidomainSolver
< DIM, DIM > * 
mpSolver

Private Member Functions

template<class Archive >
void save (Archive &archive, const unsigned int version) const
template<class Archive >
void load (Archive &archive, const unsigned int version)

Friends

class boost::serialization::access
class TestArchivingExtendedBidomain

Detailed Description

template<unsigned DIM>
class ExtendedBidomainProblem< DIM >

Class which specifies and solves an extended bidomain problem.

See Buist ML, Poh YC. An Extended Bidomain Framework Incorporating Multiple Cell Types. Biophysical Journal, Volume 99, Issue 1, 13-18, 7 July 2010.

Briefly, the problems consits of 3 equations, 2 parabolic and one elliptic. The space is divided in 3 compartments:

  • First cell
  • Second cell
  • Extracellular space.

The three unkowns are the the intracellular potential of the two cells and the extracellular potential. This class allows the user to specify the parameters specific for the these simulations and also different intracellular conductivities for the two cells.

The solution vector used for calculations is arranged as a striped vector with this order:

  • Intracellular potential of the first cell
  • Intracellular potential of the second cell
  • Extracellular potential

However, the OUTPUT is rearranged as follows. The rearrangament occurs AFTER the solution is calculated:

  • Transmembrane potential of the first cell (Intracellular potential of the first cell - Extracellular potential)
  • Transmembrane potential of the second cell (Intracellular potential of the second cell - Extracellular potential)
  • Extracellular potential

Unlike a bidomain problem, a node-wise extracellular stimulus can be set up in extended bidomain problems in absence of a bath. This is done by setting a stimulus factory via the method SetExtracellularStimulusFactory. See documentation of AbstractStimulusFactory and ElectrodesStimulusFactory for more details. Note that compatibility conditions will be taken care of by this class by calling specific methods within the stimulus factory class.

Definition at line 91 of file ExtendedBidomainProblem.hpp.


Constructor & Destructor Documentation

template<unsigned DIM>
ExtendedBidomainProblem< DIM >::ExtendedBidomainProblem ( AbstractCardiacCellFactory< DIM > *  pCellFactory,
AbstractCardiacCellFactory< DIM > *  pSecondCellFactory,
bool  hasBath = false 
)

Constructor

Parameters:
pCellFactoryUser defined cell factory which shows how the pde should create cells.
pSecondCellFactoryUser defined cell factory which shows how the pde should create cells.
hasBathWhether the simulation has a bath (if this is true, all elements with attribute = 1 will be set to be bath elements (the rest should have attribute = 0)).

Definition at line 48 of file ExtendedBidomainProblem.cpp.

References ExtendedBidomainProblem< DIM >::mFixedExtracellularPotentialNodes.

template<unsigned DIM>
ExtendedBidomainProblem< DIM >::ExtendedBidomainProblem ( )

Archiving constructor

Definition at line 65 of file ExtendedBidomainProblem.cpp.

References ExtendedBidomainProblem< DIM >::mFixedExtracellularPotentialNodes.

template<unsigned DIM>
ExtendedBidomainProblem< DIM >::~ExtendedBidomainProblem ( )

Destructor

Definition at line 266 of file ExtendedBidomainProblem.cpp.


Member Function Documentation

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::AnalyseMeshForBath ( ) [protected]

Annotate bath nodes with the correct region code, if a bath is present. Will throw if mHasBath is set but no bath is present in the mesh.

template<unsigned DIM>
AbstractCardiacTissue< DIM > * ExtendedBidomainProblem< DIM >::CreateCardiacTissue ( ) [protected, virtual]
template<unsigned DIM>
Vec ExtendedBidomainProblem< DIM >::CreateInitialCondition ( ) [protected, virtual]

Create normal initial condition but overwrite V to zero for bath nodes, if there are any.

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 81 of file ExtendedBidomainProblem.cpp.

References DistributedVector::Begin(), DistributedVectorFactory::CreateDistributedVector(), DistributedVectorFactory::CreateVec(), DistributedVector::End(), and DistributedVector::Restore().

template<unsigned DIM>
AbstractDynamicLinearPdeSolver< DIM, DIM, 3 > * ExtendedBidomainProblem< DIM >::CreateSolver ( ) [protected, virtual]
template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::DefineWriterColumns ( bool  extending) [virtual]

Define what variables are written to the primary results file. Hardcoded variable names are "V", "V_2" and "Phi_e" for the three variables. If you request any extra variable (via HeartConfig), this method will also define the extra variables by calling a method in the parent class.

Parameters:
extendingwhether we are extending an existing results file

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 349 of file ExtendedBidomainProblem.cpp.

References AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DefineExtraVariablesWriterColumns().

template<unsigned DIM>
ExtendedBidomainTissue< DIM > * ExtendedBidomainProblem< DIM >::GetExtendedBidomainTissue ( )

Access the tissue object. Can only be called after Initialise()

Definition at line 312 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
bool ExtendedBidomainProblem< DIM >::GetHasBath ( ) [virtual]

Return whether this is a problem with bath or not

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 440 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
template<class Archive >
void ExtendedBidomainProblem< DIM >::load ( Archive &  archive,
const unsigned int  version 
) [inline, private]

Load the member variables from an archive.

Parameters:
archive
version

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 198 of file ExtendedBidomainProblem.hpp.

References DistributedVector::Begin(), PetscTools::Destroy(), DistributedVector::End(), ArchiveLocationInfo::GetArchiveDirectory(), ArchiveLocationInfo::GetArchiveRelativePath(), Hdf5DataReader::GetVariableOverNodes(), FileFinder::IsAbsolutePath(), ExtendedBidomainProblem< DIM >::mAmFirstCell, ExtendedBidomainProblem< DIM >::mAmGap, ExtendedBidomainProblem< DIM >::mAmSecondCell, ExtendedBidomainProblem< DIM >::mApplyAveragePhieZeroConstraintAfterSolving, ExtendedBidomainProblem< DIM >::mCmFirstCell, ExtendedBidomainProblem< DIM >::mCmSecondCell, ExtendedBidomainProblem< DIM >::mFixedExtracellularPotentialNodes, ExtendedBidomainProblem< DIM >::mGGap, ExtendedBidomainProblem< DIM >::mGgapHeterogeneityRegions, ExtendedBidomainProblem< DIM >::mGgapHeterogenousValues, ExtendedBidomainProblem< DIM >::mHasBath, ExtendedBidomainProblem< DIM >::mIntracellularConductivitiesSecondCell, ExtendedBidomainProblem< DIM >::mpExtendedBidomainTissue, AbstractCardiacProblem< DIM, DIM, 3 >::mpMesh, ExtendedBidomainProblem< DIM >::mRowForAverageOfPhiZeroed, AbstractCardiacProblem< DIM, DIM, 3 >::mSolution, ExtendedBidomainProblem< DIM >::mUserHasSetBidomainValuesExplicitly, ExtendedBidomainProblem< DIM >::mUserSpecifiedSecondCellConductivities, ExtendedBidomainProblem< DIM >::mUserSuppliedExtracellularStimulus, ExtendedBidomainProblem< DIM >::mVariablesIDs, and DistributedVector::Restore().

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::PreSolveChecks ( ) [virtual]

Performs a series of checks before solving. Checks that a suitable method of resolving the singularity is being used.

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 420 of file ExtendedBidomainProblem.cpp.

References EXCEPTION, HeartConfig::GetUseRelativeTolerance(), HeartConfig::Instance(), and AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::PreSolveChecks().

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::ProcessExtracellularStimulus ( ) [protected]

Small helper tidy-up method that incorporates everything related to the initialisation of the extracellular stimulus.

It checks whether there is any extracellular stimulus, If not, it creates a one (with default implementation to zero stimulus). In any case, it passes the mesh into the stimulus, calls the SetCompatibleExtracellularStimulus method and also checks if there are any grounded nodes. If so, it calls SetFixedExtracellularPotentialNodes.

Definition at line 108 of file ExtendedBidomainProblem.cpp.

References PetscTools::Barrier().

template<unsigned DIM>
template<class Archive >
void ExtendedBidomainProblem< DIM >::save ( Archive &  archive,
const unsigned int  version 
) const [inline, private]

Save the member variables to an archive.

Parameters:
archive
version

Todo:
#1369

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 105 of file ExtendedBidomainProblem.hpp.

References DistributedVector::Begin(), Hdf5DataWriter::DefineFixedDimension(), PetscTools::Destroy(), DistributedVector::End(), ArchiveLocationInfo::GetArchiveDirectory(), ArchiveLocationInfo::GetArchiveRelativePath(), HeartConfig::Instance(), ExtendedBidomainProblem< DIM >::mAmFirstCell, ExtendedBidomainProblem< DIM >::mAmGap, ExtendedBidomainProblem< DIM >::mAmSecondCell, ExtendedBidomainProblem< DIM >::mApplyAveragePhieZeroConstraintAfterSolving, ExtendedBidomainProblem< DIM >::mCmFirstCell, ExtendedBidomainProblem< DIM >::mCmSecondCell, ExtendedBidomainProblem< DIM >::mFixedExtracellularPotentialNodes, ExtendedBidomainProblem< DIM >::mGGap, ExtendedBidomainProblem< DIM >::mGgapHeterogeneityRegions, ExtendedBidomainProblem< DIM >::mGgapHeterogenousValues, ExtendedBidomainProblem< DIM >::mHasBath, ExtendedBidomainProblem< DIM >::mIntracellularConductivitiesSecondCell, ExtendedBidomainProblem< DIM >::mpExtendedBidomainTissue, AbstractCardiacProblem< DIM, DIM, 3 >::mpMesh, ExtendedBidomainProblem< DIM >::mRowForAverageOfPhiZeroed, AbstractCardiacProblem< DIM, DIM, 3 >::mSolution, ExtendedBidomainProblem< DIM >::mUserHasSetBidomainValuesExplicitly, ExtendedBidomainProblem< DIM >::mUserSpecifiedSecondCellConductivities, ExtendedBidomainProblem< DIM >::mUserSuppliedExtracellularStimulus, ExtendedBidomainProblem< DIM >::mVariablesIDs, and DistributedVector::Restore().

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetExtendedBidomainParameters ( double  Am1,
double  Am2,
double  AmGap,
double  Cm1,
double  Cm2,
double  Ggap 
)

Allow the user to specify different values for the bidomain parameters in the extened bidomain framework. This method switches the mUserHasSetBidomainValuesExplicitly to true and tells the problem to overwrite the corresponding parameters in HeartConfig. It needs to be called before Initialise() to have any effect.

Parameters:
Am1value of the surface-to-volume ratio for the first cell
Am2value of the surface-to-volume ratio for the second cell
AmGapvalue of the surface-to-volume ratio for the gap junction
Cm1value of the capacitance for the first cell
Cm2value of the capacitance for the second cell
Ggapvalue, in mS of the gap junction conductance

Definition at line 200 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetExtracellularStimulusFactory ( AbstractStimulusFactory< DIM > *  pFactory)

Sets a stimulus factory as the extracellular one.

Parameters:
pFactorythe factory to be set.

Definition at line 224 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetFixedExtracellularPotentialNodes ( std::vector< unsigned nodes)

Set the nodes at which phi_e (the extracellular potential) is fixed to zero. This does not necessarily have to be called. If it is not, phi_e is only defined up to a constant.

Parameters:
nodesthe nodes to be fixed.
Note:
currently, the value of phi_e at the fixed nodes cannot be set to be anything other than zero.

Definition at line 285 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetGgapHeterogeneities ( std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &  rGgapHeterogeneityRegions,
std::vector< double > &  rGgapValues 
)

Set the values of mCellHeterogeneityRegions and mGgapValues for the heterogeneities of Ggap. It just sets the member variables here that will later be passed on to the tissue object. It also checks that the two have the same size. Throws otherwise.

Parameters:
rGgapHeterogeneityRegionsa vector of heterogeneity regions for gap junctions
rGgapValuesa vector (of the same size as rGgapHeterogeneityRegions) with the respective values of Ggap for every region.

Definition at line 213 of file ExtendedBidomainProblem.cpp.

References EXCEPTION.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetHasBath ( bool  hasBath)

Set whether there is a bath or not. Useful for covering exceptions only (for example if bath problems are not supported and there is no method to analyse the mesh for a bath).

Parameters:
hasBathwhether we want the problem to have bath or not.

Definition at line 447 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetIntracellularConductivitiesForSecondCell ( c_vector< double, DIM >  conductivities)

Allow the user to specify different values for the conductivities of the second cell type

Parameters:
conductivities: the intracellular conductivities of the second cell that you wish to set

Definition at line 275 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::SetNodeForAverageOfPhiZeroed ( unsigned  node)

Set which row of the linear system should be used to enforce the condition that the average of phi_e is zero. If not called, this condition will not be used.

Parameters:
nodethe mesh node index giving the row at which to impose the constraint

Definition at line 298 of file ExtendedBidomainProblem.cpp.

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::WriteInfo ( double  time) [virtual]

Print out time and max/min of the intracellular potential of the two cells and phi_e values at current time.

Parameters:
timecurrent time.

Implements AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 319 of file ExtendedBidomainProblem.cpp.

References PetscTools::AmMaster().

template<unsigned DIM>
void ExtendedBidomainProblem< DIM >::WriteOneStep ( double  time,
Vec  voltageVec 
) [virtual]

Write one timestep of output data to the primary results file. The solution of the problem is in term of Phi_i_1, Phi_i_2 and Phi_e (i.e., intracellular potential of the two cells plus extracellular potential). This method works out the transmembrane potential of the two cells (V and V_2) calculated as V = Phi_i_1 - Phi_e and V_2 = Phi_i_2-Phi_e. It then writes to file V, V_2 and Phi_e.

It also writes any extra variable (defined by HeartConfig). Note that it does the job by calling the method in the parent class. This implies that the extra variable must be in the first cell (not the second) because the generic method to write extra variables only looks into the first cell factory. TODO: write a specific method for extended problems that looks in both cell factories

Parameters:
timethe current time
voltageVecthe solution vector to write

Implements AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 384 of file ExtendedBidomainProblem.cpp.

References DistributedVector::Begin(), PetscTools::Destroy(), DistributedVector::End(), DistributedVector::Restore(), and AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteExtraVariablesOneStep().


Friends And Related Function Documentation

template<unsigned DIM>
friend class boost::serialization::access [friend]

Needed for serialization.

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 94 of file ExtendedBidomainProblem.hpp.


Member Data Documentation

template<unsigned DIM>
double ExtendedBidomainProblem< DIM >::mAmFirstCell [protected]

Am for the first cell, set by the user and, if so, set into the PDE

Definition at line 307 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
double ExtendedBidomainProblem< DIM >::mAmGap [protected]

Am for the gap junctions, set by the user and, if so, set into the PDE

Definition at line 311 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
double ExtendedBidomainProblem< DIM >::mAmSecondCell [protected]

Am for the second cell, set by the user and, if so, set into the PDE

Definition at line 309 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

An alternative method for constraining the solution and resolving singularity. Singular system is solved and THEN the constraint is applied on the phi_e (after calculating.

This variable is checked within the specialization of "OnEndOfTimestep" method. False by default as initialised in the constructor.

Definition at line 346 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
double ExtendedBidomainProblem< DIM >::mCmFirstCell [protected]

Cm for the first cell, set by the user and, if so, set into the PDE

Definition at line 313 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
double ExtendedBidomainProblem< DIM >::mCmSecondCell [protected]

Cm for the second cell, set by the user and, if so, set into the PDE

Definition at line 315 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
std::vector<unsigned> ExtendedBidomainProblem< DIM >::mFixedExtracellularPotentialNodes [protected]

Nodes at which the extracellular voltage is fixed to zero (replicated)

Definition at line 286 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::ExtendedBidomainProblem(), ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
double ExtendedBidomainProblem< DIM >::mGGap [protected]

Conductance, in mS of the gap junction conductance, set by the user and, if so, set into the PDE (otherwise its default value is 0

Definition at line 317 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
std::vector<boost::shared_ptr<AbstractChasteRegion<DIM> > > ExtendedBidomainProblem< DIM >::mGgapHeterogeneityRegions [protected]

Vector of Ggap heterogeneity regions. Set by the user, it is passed on to the tissue object. If empty (i.e., the user did not specify any heterogeneities) then the empty vector will still be passed to the tissue object, which, seeing the empty vector, will apply mGgap everywhere

Definition at line 324 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
std::vector<double> ExtendedBidomainProblem< DIM >::mGgapHeterogenousValues [protected]

Corresponding vector of values of Ggap for every heterogeneous region in mGgapHeterogeneityRegions

Definition at line 327 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
bool ExtendedBidomainProblem< DIM >::mHasBath [protected]

Whether the mesh has a bath, ie whether this is a bath simulation

Definition at line 352 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
c_vector<double, DIM> ExtendedBidomainProblem< DIM >::mIntracellularConductivitiesSecondCell [protected]

stores the values of the conductivities for the second cell.

Definition at line 289 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
ExtendedBidomainTissue<DIM>* ExtendedBidomainProblem< DIM >::mpExtendedBidomainTissue [protected]

The bidomain tissue object for the extended problem

Definition at line 283 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
AbstractStimulusFactory<DIM>* ExtendedBidomainProblem< DIM >::mpExtracellularStimulusFactory [protected]

Factory to generate the stimulus

Definition at line 331 of file ExtendedBidomainProblem.hpp.

template<unsigned DIM>
AbstractCardiacCellFactory<DIM,DIM>* ExtendedBidomainProblem< DIM >::mpSecondCellFactory [protected]

The cell factory creates the cells for each node

Definition at line 280 of file ExtendedBidomainProblem.hpp.

template<unsigned DIM>
AbstractExtendedBidomainSolver<DIM,DIM>* ExtendedBidomainProblem< DIM >::mpSolver [protected]

We need to save the assembler that is being used to be able to switch off the electrodes (by adding default boundary conditions to the assembler)

Reimplemented from AbstractCardiacProblem< DIM, DIM, 3 >.

Definition at line 382 of file ExtendedBidomainProblem.hpp.

template<unsigned DIM>
int ExtendedBidomainProblem< DIM >::mRowForAverageOfPhiZeroed [protected]

Another method of resolving the singularity in the bidomain equations. Specifies a row of the matrix at which to impose an extra condition.

Definition at line 337 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
bool ExtendedBidomainProblem< DIM >::mUserHasSetBidomainValuesExplicitly [protected]

flag used to check whwther the user wanted specific, out-of-heartconfig values

Definition at line 304 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
bool ExtendedBidomainProblem< DIM >::mUserSpecifiedSecondCellConductivities [protected]

Flag for checking that the user specified conductivities for the second cell. Get method available for this variable.

Definition at line 301 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
bool ExtendedBidomainProblem< DIM >::mUserSuppliedExtracellularStimulus [protected]

keeps track of whether the user set an extracellular stimulus. True if the user did.

Definition at line 349 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
std::vector<signed int> ExtendedBidomainProblem< DIM >::mVariablesIDs [protected]

Used by the writer, stores the variable output names

Definition at line 298 of file ExtendedBidomainProblem.hpp.

Referenced by ExtendedBidomainProblem< DIM >::load(), and ExtendedBidomainProblem< DIM >::save().

template<unsigned DIM>
unsigned ExtendedBidomainProblem< DIM >::mVoltageColumnId_Phie [protected]

Used by the writer, extracellular potential

Definition at line 296 of file ExtendedBidomainProblem.hpp.

template<unsigned DIM>
unsigned ExtendedBidomainProblem< DIM >::mVoltageColumnId_Vm1 [protected]

Used by the writer, transmembrane potential of the first cell

Definition at line 292 of file ExtendedBidomainProblem.hpp.

template<unsigned DIM>
unsigned ExtendedBidomainProblem< DIM >::mVoltageColumnId_Vm2 [protected]

Used by the writer, transmembrane potential of the second cell

Definition at line 294 of file ExtendedBidomainProblem.hpp.


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