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 "OperatorSplittingMonodomainSolver.hpp"
00030
00031
00032 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00033 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00034 {
00035 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00036 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00037
00038 if(!this->mpMonodomainAssembler)
00039 {
00040 this->mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00041 }
00042
00044
00046 if(computeMatrix)
00047 {
00048 this->mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00049 this->mpMonodomainAssembler->AssembleMatrix();
00050
00051 MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00052 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00053 mass_matrix_assembler.Assemble();
00054
00055 this->mpLinearSystem->AssembleFinalLhsMatrix();
00056 PetscMatTools::AssembleFinal(mMassMatrix);
00057 }
00058
00059 HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00060
00062
00064 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00065
00066 DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00067
00068 DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00069
00070 double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00071 double Cm = HeartConfig::Instance()->GetCapacitance();
00072
00073 for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00074 index!= dist_vec_matrix_based.End();
00075 ++index)
00076 {
00077 double V = distributed_current_solution[index];
00078
00079
00080
00081
00082
00083 dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse();
00084 }
00085 dist_vec_matrix_based.Restore();
00086
00088
00090 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00091
00092
00093
00094 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00095
00097
00099 this->mpMonodomainAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00100 this->mpMonodomainAssembler->SetApplyNeummanBoundaryConditionsToVector(this->mpBoundaryConditions);
00101 this->mpMonodomainAssembler->OnlyAssembleOnSurfaceElements();
00102
00103
00104 this->mpMonodomainAssembler->AssembleVector();
00105
00106
00107 this->mpLinearSystem->AssembleRhsVector();
00108 }
00109
00110 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00111 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution)
00112 {
00113 double time = PdeSimulationTime::GetTime();
00114 double dt = PdeSimulationTime::GetPdeTimeStep();
00115 mpMonodomainTissue->SolveCellSystems(currentSolution, time, time+dt/2.0, true);
00116 }
00117
00118
00119 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00120 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::FollowingSolveLinearSystem(Vec currentSolution)
00121 {
00122
00123 double time = PdeSimulationTime::GetTime();
00124 double dt = PdeSimulationTime::GetPdeTimeStep();
00125 mpMonodomainTissue->SolveCellSystems(currentSolution, time + dt/2, time+dt, true);
00126 }
00127
00128
00129 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00130 void OperatorSplittingMonodomainSolver<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 NEVER_REACHED;
00148
00149
00150 }
00151
00152 this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00153 this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00154 this->mpLinearSystem->SetMatrixIsSymmetric(true);
00155 this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
00156
00157
00158
00159 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00160 VecDuplicate(r_template, &mVecForConstructingRhs);
00161 PetscInt ownership_range_lo;
00162 PetscInt ownership_range_hi;
00163 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00164 PetscInt local_size = ownership_range_hi - ownership_range_lo;
00165 PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00166 this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00167 local_size, local_size);
00168 }
00169
00170
00171
00172 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00173 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::OperatorSplittingMonodomainSolver(
00174 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00175 MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00176 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00177 unsigned numQuadPoints)
00178 : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
00179 mpBoundaryConditions(pBoundaryConditions),
00180 mpMonodomainTissue(pTissue),
00181 mNumQuadPoints(numQuadPoints)
00182 {
00183 assert(pTissue);
00184 assert(pBoundaryConditions);
00185 this->mMatrixIsConstant = true;
00186 mpMonodomainAssembler = NULL;
00187
00188
00189 pTissue->SetCacheReplication(false);
00190 mVecForConstructingRhs = NULL;
00191 }
00192
00193 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00194 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~OperatorSplittingMonodomainSolver()
00195 {
00196 if(mpMonodomainAssembler)
00197 {
00198 delete mpMonodomainAssembler;
00199 }
00200
00201 if(mVecForConstructingRhs)
00202 {
00203 VecDestroy(mVecForConstructingRhs);
00204 MatDestroy(mMassMatrix);
00205 }
00206 }
00207
00208
00209
00211
00213
00214 template class OperatorSplittingMonodomainSolver<1,1>;
00215 template class OperatorSplittingMonodomainSolver<1,2>;
00216 template class OperatorSplittingMonodomainSolver<1,3>;
00217 template class OperatorSplittingMonodomainSolver<2,2>;
00218 template class OperatorSplittingMonodomainSolver<3,3>;
00219