37 #include "ExtendedBidomainSolver.hpp" 38 #include "ExtendedBidomainMassMatrixAssembler.hpp" 39 #include "ExtendedBidomainAssembler.hpp" 41 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
44 if (this->mpLinearSystem != NULL)
52 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
53 VecDuplicate(r_template, &mVecForConstructingRhs);
56 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
57 PetscInt local_size = ownership_range_hi - ownership_range_lo;
59 3*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
60 local_size, local_size);
63 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
68 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
69 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
70 assert(currentSolution != NULL);
77 mpExtendedBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
78 mpExtendedBidomainAssembler->AssembleMatrix();
82 assert(SPACE_DIM==ELEMENT_DIM);
87 this->mpLinearSystem->SwitchWriteModeLhsMatrix();
100 double Am1 = this->mpExtendedBidomainTissue->GetAmFirstCell();
101 double Am2 = this->mpExtendedBidomainTissue->GetAmSecondCell();
102 double AmGap = this->mpExtendedBidomainTissue->GetAmGap();
103 double Cm1 = this->mpExtendedBidomainTissue->GetCmFirstCell();
104 double Cm2 = this->mpExtendedBidomainTissue->GetCmSecondCell();
119 index!= dist_vec_matrix_based.
End();
122 double V_first_cell = distributed_current_solution_v_first_cell[index];
123 double V_second_Cell = distributed_current_solution_v_second_cell[index];
125 double i_ionic_first_cell = this->mpExtendedBidomainTissue->rGetIionicCacheReplicated()[index.Global];
126 double i_ionic_second_cell = this->mpExtendedBidomainTissue->rGetIionicCacheReplicatedSecondCell()[index.Global];
127 double intracellular_stimulus_first_cell = this->mpExtendedBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
128 double intracellular_stimulus_second_cell = this->mpExtendedBidomainTissue->rGetIntracellularStimulusCacheReplicatedSecondCell()[index.Global];
129 double extracellular_stimulus = this->mpExtendedBidomainTissue->rGetExtracellularStimulusCacheReplicated()[index.Global];
130 double g_gap = this->mpExtendedBidomainTissue->rGetGgapCacheReplicated()[index.Global];
132 dist_vec_matrix_based_v_first_cell[index] = Am1*Cm1*V_first_cell/delta_t - Am1*i_ionic_first_cell + AmGap*g_gap*(V_second_Cell - V_first_cell) - intracellular_stimulus_first_cell;
133 dist_vec_matrix_based_v_second_cell[index] = Am2*Cm2*V_second_Cell/delta_t - Am2*i_ionic_second_cell + AmGap*g_gap*(V_first_cell - V_second_Cell) - intracellular_stimulus_second_cell;
135 if (this->mpExtendedBidomainTissue->HasTheUserSuppliedExtracellularStimulus() )
137 assert((fabs(intracellular_stimulus_first_cell) < 1e-12)
138 && (fabs(intracellular_stimulus_second_cell) < 1e-12));
144 dist_vec_matrix_based_phi_e[index] = -extracellular_stimulus;
148 dist_vec_matrix_based_phi_e[index] = 0.0;
153 dist_vec_matrix_based.
Restore();
159 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
168 mpExtendedBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions);
169 mpExtendedBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(),
false);
170 mpExtendedBidomainNeumannSurfaceTermAssembler->AssembleVector();
172 this->mpLinearSystem->FinaliseRhsVector();
174 this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
178 this->mpLinearSystem->FinaliseLhsMatrix();
180 this->mpLinearSystem->FinaliseRhsVector();
183 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
199 EXCEPTION(
"Bath simulations are not yet supported for extended bidomain problems");
210 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
ExtendedBidomainNeumannSurfaceTermAssembler< ELEM_DIM, SPACE_DIM > * mpExtendedBidomainNeumannSurfaceTermAssembler
ExtendedBidomainSolver(bool bathSimulation, AbstractTetrahedralMesh< ELEM_DIM, SPACE_DIM > *pMesh, ExtendedBidomainTissue< SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, 3 > *pBoundaryConditions)
ExtendedBidomainAssembler< ELEM_DIM, SPACE_DIM > * mpExtendedBidomainAssembler
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
#define EXCEPTION(message)
static void BeginEvent(unsigned event)
ExtendedBidomainTissue< SPACE_DIM > * mpExtendedBidomainTissue
void SetupLinearSystem(Vec currentSolution, bool computeMatrix)
Vec mVecForConstructingRhs
static double GetPdeTimeStep()
~ExtendedBidomainSolver()
void InitialiseForSolve(Vec initialSolution)
void SetCacheReplication(bool doCacheReplication)
static void EndEvent(unsigned event)
void SetMatrixToAssemble(Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
void InitialiseForSolve(Vec initialSolution)