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
00032
00033 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00034 void MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00035 {
00036 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00037 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00038
00039 if(!this->mpMonodomainAssembler)
00040 {
00041 this->mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mDt,this->mNumQuadPoints);
00042 }
00043
00045
00047 if(computeMatrix)
00048 {
00049 this->mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00050 this->mpMonodomainAssembler->AssembleMatrix();
00051
00052 MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh);
00053 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00054 mass_matrix_assembler.Assemble();
00055
00056 this->mpLinearSystem->AssembleFinalLhsMatrix();
00057 MatAssemblyBegin(mMassMatrix, MAT_FINAL_ASSEMBLY);
00058 MatAssemblyEnd(mMassMatrix, MAT_FINAL_ASSEMBLY);
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 double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00081 - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00082
00083 dist_vec_matrix_based[index] = Am*Cm*V*this->mDtInverse + F;
00084 }
00085 dist_vec_matrix_based.Restore();
00087
00089 this->mpLinearSystem->ZeroRhsVector();
00090
00091 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00092
00093
00094
00095 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00096
00098
00100 this->mpMonodomainAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00101 this->mpMonodomainAssembler->SetApplyNeummanBoundaryConditionsToVector(this->mpBoundaryConditions);
00102 this->mpMonodomainAssembler->OnlyAssembleOnSurfaceElements();
00103
00104
00105 this->mpMonodomainAssembler->AssembleVector();
00106
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126 this->mpLinearSystem->AssembleRhsVector();
00127 }
00128
00129
00130
00131 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00132 void MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00133 {
00134 if (this->mpLinearSystem != NULL)
00135 {
00136 return;
00137 }
00138 AbstractMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00139
00140
00141
00142 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00143 VecDuplicate(r_template, &mVecForConstructingRhs);
00144 PetscInt ownership_range_lo;
00145 PetscInt ownership_range_hi;
00146 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00147 PetscInt local_size = ownership_range_hi - ownership_range_lo;
00148 PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00149 this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00150 local_size, local_size);
00151 }
00152
00153
00154
00155 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00156 MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::MatrixBasedMonodomainSolver(
00157 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00158 MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00159 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00160 unsigned numQuadPoints)
00161 : AbstractMonodomainSolver<ELEMENT_DIM,SPACE_DIM>(pMesh,pTissue,pBoundaryConditions,numQuadPoints)
00162 {
00163
00164 pTissue->SetCacheReplication(false);
00165 mVecForConstructingRhs = NULL;
00166
00168
00169
00170 }
00171
00172 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00173 MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~MatrixBasedMonodomainSolver()
00174 {
00175 if(mVecForConstructingRhs)
00176 {
00177 VecDestroy(mVecForConstructingRhs);
00178 MatDestroy(mMassMatrix);
00179 }
00180
00182
00183
00184
00185
00186
00187
00188
00189
00190 }
00191
00192
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00205
00207
00208 template class MatrixBasedMonodomainSolver<1,1>;
00209 template class MatrixBasedMonodomainSolver<1,2>;
00210 template class MatrixBasedMonodomainSolver<1,3>;
00211 template class MatrixBasedMonodomainSolver<2,2>;
00212 template class MatrixBasedMonodomainSolver<3,3>;
00213