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