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 "MatrixBasedMonodomainSolver.hpp"
00030 #include "MassMatrixAssembler.hpp"
00031 #include "PetscMatTools.hpp"
00032
00033
00034 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00035 void MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00036 {
00037 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00038 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00039
00040 if(!this->mpMonodomainAssembler)
00041 {
00042 this->mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00043 }
00044
00045
00047
00049 if(computeMatrix)
00050 {
00051 this->mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00052 this->mpMonodomainAssembler->AssembleMatrix();
00053
00054 MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00055 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00056 mass_matrix_assembler.Assemble();
00057
00058 this->mpLinearSystem->AssembleFinalLhsMatrix();
00059 PetscMatTools::AssembleFinal(mMassMatrix);
00060
00061 if (HeartConfig::Instance()->GetUseMassLumpingForPrecond() && !HeartConfig::Instance()->GetUseMassLumping())
00062 {
00063 this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
00064
00065 MonodomainAssembler<ELEMENT_DIM,SPACE_DIM> lumped_mass_assembler(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00066 lumped_mass_assembler.SetMatrixToAssemble(this->mpLinearSystem->rGetPrecondMatrix());
00067
00068 HeartConfig::Instance()->SetUseMassLumping(true);
00069 lumped_mass_assembler.AssembleMatrix();
00070 HeartConfig::Instance()->SetUseMassLumping(false);
00071
00072 this->mpLinearSystem->AssembleFinalPrecondMatrix();
00073 }
00074
00075 }
00076
00077 HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00078
00080
00082 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00083
00084 DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00085
00086 DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00087
00088 double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00089 double Cm = HeartConfig::Instance()->GetCapacitance();
00090
00091 for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00092 index!= dist_vec_matrix_based.End();
00093 ++index)
00094 {
00095 double V = distributed_current_solution[index];
00096 double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00097 - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00098
00099 dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00100 }
00101 dist_vec_matrix_based.Restore();
00102
00104
00106 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00107
00108
00109
00110 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00111
00113
00115 this->mpMonodomainAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00116 this->mpMonodomainAssembler->SetApplyNeummanBoundaryConditionsToVector(this->mpBoundaryConditions);
00117 this->mpMonodomainAssembler->OnlyAssembleOnSurfaceElements();
00118
00119
00120 this->mpMonodomainAssembler->AssembleVector();
00121
00123
00125 if(mpMonodomainCorrectionTermAssembler)
00126 {
00127 mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00128
00129 mpMonodomainCorrectionTermAssembler->AssembleVector();
00130 }
00131
00132
00133 this->mpLinearSystem->AssembleRhsVector();
00134 }
00135
00136
00137
00138 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00139 void MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00140 {
00141 if (this->mpLinearSystem != NULL)
00142 {
00143 return;
00144 }
00145 AbstractMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00146
00147
00148
00149 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00150 VecDuplicate(r_template, &mVecForConstructingRhs);
00151 PetscInt ownership_range_lo;
00152 PetscInt ownership_range_hi;
00153 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00154 PetscInt local_size = ownership_range_hi - ownership_range_lo;
00155 PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00156 this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00157 local_size, local_size);
00158 }
00159
00160
00161
00162 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00163 MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::MatrixBasedMonodomainSolver(
00164 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00165 MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00166 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00167 unsigned numQuadPoints)
00168 : AbstractMonodomainSolver<ELEMENT_DIM,SPACE_DIM>(pMesh,pTissue,pBoundaryConditions,numQuadPoints)
00169 {
00170
00171 pTissue->SetCacheReplication(false);
00172 mVecForConstructingRhs = NULL;
00173
00174 if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00175 {
00176 mpMonodomainCorrectionTermAssembler
00177 = new MonodomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00178
00179 pTissue->SetCacheReplication(true);
00180 }
00181 else
00182 {
00183 mpMonodomainCorrectionTermAssembler = NULL;
00184 }
00185 }
00186
00187 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00188 MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~MatrixBasedMonodomainSolver()
00189 {
00190 if(mVecForConstructingRhs)
00191 {
00192 VecDestroy(mVecForConstructingRhs);
00193 MatDestroy(mMassMatrix);
00194 }
00195
00196 if(mpMonodomainCorrectionTermAssembler)
00197 {
00198 delete mpMonodomainCorrectionTermAssembler;
00199 }
00200 }
00201
00202
00204
00206
00207 template class MatrixBasedMonodomainSolver<1,1>;
00208 template class MatrixBasedMonodomainSolver<1,2>;
00209 template class MatrixBasedMonodomainSolver<1,3>;
00210 template class MatrixBasedMonodomainSolver<2,2>;
00211 template class MatrixBasedMonodomainSolver<3,3>;
00212