36 #include "MonodomainSolver.hpp"
37 #include "MassMatrixAssembler.hpp"
38 #include "PetscMatTools.hpp"
41 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
44 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
45 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
53 mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
54 mpMonodomainAssembler->AssembleMatrix();
60 this->mpLinearSystem->FinaliseLhsMatrix();
65 this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
74 this->mpLinearSystem->FinalisePrecondMatrix();
94 index!= dist_vec_matrix_based.
End();
97 double V = distributed_current_solution[index];
98 double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
99 - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
103 dist_vec_matrix_based.
Restore();
108 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
117 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(),
false);
118 mpNeumannSurfaceTermsAssembler->AssembleVector();
123 if(mpMonodomainCorrectionTermAssembler)
125 mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(),
false);
127 mpMonodomainCorrectionTermAssembler->AssembleVector();
131 this->mpLinearSystem->FinaliseRhsVector();
136 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
139 if (this->mpLinearSystem != NULL)
159 this->mpLinearSystem->SetMatrixIsSymmetric(
true);
164 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
165 VecDuplicate(r_template, &mVecForConstructingRhs);
168 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
169 PetscInt local_size = ownership_range_hi - ownership_range_lo;
171 this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
172 local_size, local_size);
175 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
182 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
188 mpMonodomainTissue(pTissue),
189 mpBoundaryConditions(pBoundaryConditions)
192 assert(pBoundaryConditions);
216 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
219 delete mpMonodomainAssembler;
220 delete mpNeumannSurfaceTermsAssembler;
222 if(mVecForConstructingRhs)
228 if(mpMonodomainCorrectionTermAssembler)
230 delete mpMonodomainCorrectionTermAssembler;
MonodomainCorrectionTermAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainCorrectionTermAssembler
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, 1 > * mpNeumannSurfaceTermsAssembler
MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainAssembler
virtual void InitialiseForSolve(Vec initialSolution=NULL)
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
static void BeginEvent(unsigned event)
void SetupLinearSystem(Vec currentSolution, bool computeMatrix)
Vec mVecForConstructingRhs
void SetUseMassLumping(bool useMassLumping=true)
static double GetNextTime()
virtual void InitialiseForSolve(Vec initialSolution)
MonodomainSolver(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions)
static double GetPdeTimeStepInverse()
void PrepareForSetupLinearSystem(Vec currentSolution)
void SetCacheReplication(bool doCacheReplication)
double GetCapacitance() const
double GetSurfaceAreaToVolumeRatio() const
static void EndEvent(unsigned event)
void SetMatrixToAssemble(Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
static HeartConfig * Instance()
virtual ~MonodomainSolver()
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * mpMonodomainTissue