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