ExtendedBidomainSolver.cpp
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
00037 #include "ExtendedBidomainSolver.hpp"
00038 #include "ExtendedBidomainMassMatrixAssembler.hpp"
00039 #include "ExtendedBidomainAssembler.hpp"
00040
00041 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00042 void ExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00043 {
00044 if (this->mpLinearSystem != NULL)
00045 {
00046 return;
00047 }
00048 AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00049
00050
00051
00052 Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00053 VecDuplicate(r_template, &mVecForConstructingRhs);
00054 PetscInt ownership_range_lo;
00055 PetscInt ownership_range_hi;
00056 VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00057 PetscInt local_size = ownership_range_hi - ownership_range_lo;
00058 PetscTools::SetupMat(mMassMatrix, 3*this->mpMesh->GetNumNodes(), 3*this->mpMesh->GetNumNodes(),
00059 3*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00060 local_size, local_size);
00061 }
00062
00063
00064 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00065 void ExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(
00066 Vec currentSolution,
00067 bool computeMatrix)
00068 {
00069 assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00070 assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00071 assert(currentSolution != NULL);
00072
00073
00075
00077 if(computeMatrix)
00078 {
00079 mpExtendedBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00080 mpExtendedBidomainAssembler->AssembleMatrix();
00081
00082
00083
00084 assert(SPACE_DIM==ELEMENT_DIM);
00085 ExtendedBidomainMassMatrixAssembler<SPACE_DIM> mass_matrix_assembler(this->mpMesh);
00086 mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00087 mass_matrix_assembler.Assemble();
00088
00089 this->mpLinearSystem->SwitchWriteModeLhsMatrix();
00090 PetscMatTools::Finalise(mMassMatrix);
00091 }
00092
00093
00094 HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00095
00097
00099 DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00100
00101
00102 double Am1 = this->mpExtendedBidomainTissue->GetAmFirstCell();
00103 double Am2 = this->mpExtendedBidomainTissue->GetAmSecondCell();
00104 double AmGap = this->mpExtendedBidomainTissue->GetAmGap();
00105 double Cm1 = this->mpExtendedBidomainTissue->GetCmFirstCell();
00106 double Cm2 = this->mpExtendedBidomainTissue->GetCmSecondCell();
00107
00108
00109 DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00110 DistributedVector::Stripe distributed_current_solution_v_first_cell(distributed_current_solution, 0);
00111 DistributedVector::Stripe distributed_current_solution_v_second_cell(distributed_current_solution, 1);
00112 DistributedVector::Stripe distributed_current_solution_phi_e(distributed_current_solution, 2);
00113
00114
00115 DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00116 DistributedVector::Stripe dist_vec_matrix_based_v_first_cell(dist_vec_matrix_based, 0);
00117 DistributedVector::Stripe dist_vec_matrix_based_v_second_cell(dist_vec_matrix_based, 1);
00118 DistributedVector::Stripe dist_vec_matrix_based_phi_e(dist_vec_matrix_based, 2);
00119
00120 for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00121 index!= dist_vec_matrix_based.End();
00122 ++index)
00123 {
00124 double V_first_cell = distributed_current_solution_v_first_cell[index];
00125 double V_second_Cell = distributed_current_solution_v_second_cell[index];
00126
00127 double i_ionic_first_cell = this->mpExtendedBidomainTissue->rGetIionicCacheReplicated()[index.Global];
00128 double i_ionic_second_cell = this->mpExtendedBidomainTissue->rGetIionicCacheReplicatedSecondCell()[index.Global];
00129 double intracellular_stimulus_first_cell = this->mpExtendedBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00130 double intracellular_stimulus_second_cell = this->mpExtendedBidomainTissue->rGetIntracellularStimulusCacheReplicatedSecondCell()[index.Global];
00131 double extracellular_stimulus = this->mpExtendedBidomainTissue->rGetExtracellularStimulusCacheReplicated()[index.Global];
00132 double g_gap = this->mpExtendedBidomainTissue->rGetGgapCacheReplicated()[index.Global];
00133 double delta_t = PdeSimulationTime::GetPdeTimeStep();
00134 dist_vec_matrix_based_v_first_cell[index] = Am1*Cm1*V_first_cell/delta_t - Am1*i_ionic_first_cell + AmGap*g_gap*(V_second_Cell - V_first_cell) - intracellular_stimulus_first_cell;
00135 dist_vec_matrix_based_v_second_cell[index] = Am2*Cm2*V_second_Cell/delta_t - Am2*i_ionic_second_cell + AmGap*g_gap*(V_first_cell - V_second_Cell) - intracellular_stimulus_second_cell;
00136
00137 if (this->mpExtendedBidomainTissue->HasTheUserSuppliedExtracellularStimulus() )
00138 {
00139 assert((fabs(intracellular_stimulus_first_cell) < 1e-12)
00140 && (fabs(intracellular_stimulus_second_cell) < 1e-12));
00141
00146 dist_vec_matrix_based_phi_e[index] = -extracellular_stimulus;
00147 }
00148 else
00149 {
00150 dist_vec_matrix_based_phi_e[index] = 0.0;
00151 }
00152 }
00153
00154
00155 dist_vec_matrix_based.Restore();
00156
00158
00160
00161 MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00162
00163
00164
00165 HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00166
00168
00170 mpExtendedBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions);
00171 mpExtendedBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false);
00172 mpExtendedBidomainNeumannSurfaceTermAssembler->AssembleVector();
00173
00174 this->mpLinearSystem->FinaliseRhsVector();
00175
00176 this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
00177
00178 if(computeMatrix)
00179 {
00180 this->mpLinearSystem->FinaliseLhsMatrix();
00181 }
00182 this->mpLinearSystem->FinaliseRhsVector();
00183 }
00184
00185
00186 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00187 ExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::ExtendedBidomainSolver(
00188 bool bathSimulation,
00189 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00190 ExtendedBidomainTissue<SPACE_DIM>* pTissue,
00191 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,3>* pBoundaryConditions)
00192 : AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
00193 {
00194
00195 pTissue->SetCacheReplication(false);
00196 mVecForConstructingRhs = NULL;
00197
00198
00199 if(this->mBathSimulation)
00200 {
00201
00202 EXCEPTION("Bath simulations are not yet supported for extended bidomain problems");
00203 }
00204 else
00205 {
00206 mpExtendedBidomainAssembler = new ExtendedBidomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpExtendedBidomainTissue);
00207 }
00208
00209 mpExtendedBidomainNeumannSurfaceTermAssembler = new ExtendedBidomainNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM>(pMesh,pBoundaryConditions);
00210
00211 }
00212
00213 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00214 ExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::~ExtendedBidomainSolver()
00215 {
00216 delete mpExtendedBidomainAssembler;
00217 delete mpExtendedBidomainNeumannSurfaceTermAssembler;
00218
00219 if(mVecForConstructingRhs)
00220 {
00221 PetscTools::Destroy(mVecForConstructingRhs);
00222 PetscTools::Destroy(mMassMatrix);
00223 }
00224 }
00225
00227
00229
00230 template class ExtendedBidomainSolver<1,1>;
00231 template class ExtendedBidomainSolver<2,2>;
00232 template class ExtendedBidomainSolver<3,3>;
00233