00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #include "MonodomainSolver.hpp"
00037 #include "MassMatrixAssembler.hpp"
00038 #include "PetscMatTools.hpp"
00039
00040
00041 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00042 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00043 {
00044 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00045 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00046
00047
00049
00051 if(computeMatrix)
00052 {
00053 mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00054 mpMonodomainAssembler->AssembleMatrix();
00055
00056 MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00057 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00058 mass_matrix_assembler.Assemble();
00059
00060 this->mpLinearSystem->FinaliseLhsMatrix();
00061 PetscMatTools::Finalise(mMassMatrix);
00062
00063 if (HeartConfig::Instance()->GetUseMassLumpingForPrecond() && !HeartConfig::Instance()->GetUseMassLumping())
00064 {
00065 this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
00066
00067 MonodomainAssembler<ELEMENT_DIM,SPACE_DIM> lumped_mass_assembler(this->mpMesh,this->mpMonodomainTissue);
00068 lumped_mass_assembler.SetMatrixToAssemble(this->mpLinearSystem->rGetPrecondMatrix());
00069
00070 HeartConfig::Instance()->SetUseMassLumping(true);
00071 lumped_mass_assembler.AssembleMatrix();
00072 HeartConfig::Instance()->SetUseMassLumping(false);
00073
00074 this->mpLinearSystem->FinalisePrecondMatrix();
00075 }
00076
00077 }
00078
00079 HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00080
00082
00084 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00085
00086 DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00087
00088 DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00089
00090 double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00091 double Cm = HeartConfig::Instance()->GetCapacitance();
00092
00093 for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00094 index!= dist_vec_matrix_based.End();
00095 ++index)
00096 {
00097 double V = distributed_current_solution[index];
00098 double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00099 - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00100
00101 dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00102 }
00103 dist_vec_matrix_based.Restore();
00104
00106
00108 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00109
00110
00111
00112 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00113
00115
00117 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00118 mpNeumannSurfaceTermsAssembler->AssembleVector();
00119
00121
00123 if(mpMonodomainCorrectionTermAssembler)
00124 {
00125 mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00126
00127 mpMonodomainCorrectionTermAssembler->AssembleVector();
00128 }
00129
00130
00131 this->mpLinearSystem->FinaliseRhsVector();
00132 }
00133
00134
00135
00136 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00137 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00138 {
00139 if (this->mpLinearSystem != NULL)
00140 {
00141 return;
00142 }
00143
00144
00145 AbstractLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>::InitialiseForSolve(initialSolution);
00146
00147
00148 if(HeartConfig::Instance()->GetUseAbsoluteTolerance())
00149 {
00150 this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance());
00151 }
00152 else
00153 {
00154 this->mpLinearSystem->SetRelativeTolerance(HeartConfig::Instance()->GetRelativeTolerance());
00155 }
00156
00157 this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00158 this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00159 this->mpLinearSystem->SetMatrixIsSymmetric(true);
00160 this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
00161
00162
00163
00164 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00165 VecDuplicate(r_template, &mVecForConstructingRhs);
00166 PetscInt ownership_range_lo;
00167 PetscInt ownership_range_hi;
00168 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00169 PetscInt local_size = ownership_range_hi - ownership_range_lo;
00170 PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00171 this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00172 local_size, local_size);
00173 }
00174
00175 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00176 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution)
00177 {
00178
00179 mpMonodomainTissue->SolveCellSystems(currentSolution, PdeSimulationTime::GetTime(), PdeSimulationTime::GetNextTime());
00180 }
00181
00182 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00183 MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::MonodomainSolver(
00184 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00185 MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00186 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions)
00187 : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
00188 mpMonodomainTissue(pTissue),
00189 mpBoundaryConditions(pBoundaryConditions)
00190 {
00191 assert(pTissue);
00192 assert(pBoundaryConditions);
00193 this->mMatrixIsConstant = true;
00194
00195 mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue);
00196 mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>(pMesh,pBoundaryConditions);
00197
00198
00199
00200 pTissue->SetCacheReplication(false);
00201 mVecForConstructingRhs = NULL;
00202
00203 if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00204 {
00205 mpMonodomainCorrectionTermAssembler
00206 = new MonodomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue);
00207
00208 pTissue->SetCacheReplication(true);
00209 }
00210 else
00211 {
00212 mpMonodomainCorrectionTermAssembler = NULL;
00213 }
00214 }
00215
00216 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00217 MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~MonodomainSolver()
00218 {
00219 delete mpMonodomainAssembler;
00220 delete mpNeumannSurfaceTermsAssembler;
00221
00222 if(mVecForConstructingRhs)
00223 {
00224 PetscTools::Destroy(mVecForConstructingRhs);
00225 PetscTools::Destroy(mMassMatrix);
00226 }
00227
00228 if(mpMonodomainCorrectionTermAssembler)
00229 {
00230 delete mpMonodomainCorrectionTermAssembler;
00231 }
00232 }
00233
00234
00236
00238
00239 template class MonodomainSolver<1,1>;
00240 template class MonodomainSolver<1,2>;
00241 template class MonodomainSolver<1,3>;
00242 template class MonodomainSolver<2,2>;
00243 template class MonodomainSolver<3,3>;
00244