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 "OperatorSplittingMonodomainSolver.hpp"
00037
00038
00039 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00040 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00041 {
00042 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00043 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00044
00046
00048 if(computeMatrix)
00049 {
00050 mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00051 mpMonodomainAssembler->AssembleMatrix();
00052
00053 MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00054 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00055 mass_matrix_assembler.Assemble();
00056
00057 this->mpLinearSystem->FinaliseLhsMatrix();
00058 PetscMatTools::Finalise(mMassMatrix);
00059 }
00060
00061 HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00062
00064
00066 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00067
00068 DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00069
00070 DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00071
00072 double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00073 double Cm = HeartConfig::Instance()->GetCapacitance();
00074
00075 for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00076 index!= dist_vec_matrix_based.End();
00077 ++index)
00078 {
00079 double V = distributed_current_solution[index];
00080
00081
00082
00083
00084
00085 dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse();
00086 }
00087 dist_vec_matrix_based.Restore();
00088
00090
00092 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00093
00094
00095
00096 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00097
00099
00101 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00102 mpNeumannSurfaceTermsAssembler->AssembleVector();
00103
00104
00105 this->mpLinearSystem->FinaliseRhsVector();
00106 }
00107
00108 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00109 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution)
00110 {
00111 double time = PdeSimulationTime::GetTime();
00112 double dt = PdeSimulationTime::GetPdeTimeStep();
00113 mpMonodomainTissue->SolveCellSystems(currentSolution, time, time+dt/2.0, true);
00114 }
00115
00116
00117 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00118 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::FollowingSolveLinearSystem(Vec currentSolution)
00119 {
00120
00121 double time = PdeSimulationTime::GetTime();
00122 double dt = PdeSimulationTime::GetPdeTimeStep();
00123 mpMonodomainTissue->SolveCellSystems(currentSolution, time + dt/2, PdeSimulationTime::GetNextTime(), true);
00124 }
00125
00126
00127 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00128 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00129 {
00130 if (this->mpLinearSystem != NULL)
00131 {
00132 return;
00133 }
00134
00135
00136 AbstractLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>::InitialiseForSolve(initialSolution);
00137
00138
00139 if(HeartConfig::Instance()->GetUseAbsoluteTolerance())
00140 {
00141 this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance());
00142 }
00143 else
00144 {
00145 NEVER_REACHED;
00146
00147
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
00169
00170 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00171 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::OperatorSplittingMonodomainSolver(
00172 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00173 MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00174 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions)
00175 : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
00176 mpBoundaryConditions(pBoundaryConditions),
00177 mpMonodomainTissue(pTissue)
00178 {
00179 assert(pTissue);
00180 assert(pBoundaryConditions);
00181 this->mMatrixIsConstant = true;
00182
00183 mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue);
00184 mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>(pMesh,pBoundaryConditions);
00185
00186
00187 pTissue->SetCacheReplication(false);
00188 mVecForConstructingRhs = NULL;
00189 }
00190
00191 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00192 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~OperatorSplittingMonodomainSolver()
00193 {
00194 delete mpMonodomainAssembler;
00195 delete mpNeumannSurfaceTermsAssembler;
00196
00197 if(mVecForConstructingRhs)
00198 {
00199 PetscTools::Destroy(mVecForConstructingRhs);
00200 PetscTools::Destroy(mMassMatrix);
00201 }
00202 }
00203
00204
00205
00207
00209
00210 template class OperatorSplittingMonodomainSolver<1,1>;
00211 template class OperatorSplittingMonodomainSolver<1,2>;
00212 template class OperatorSplittingMonodomainSolver<1,3>;
00213 template class OperatorSplittingMonodomainSolver<2,2>;
00214 template class OperatorSplittingMonodomainSolver<3,3>;
00215