Chaste  Release::3.4
ExtendedBidomainSolver.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 
37 #include "ExtendedBidomainSolver.hpp"
38 #include "ExtendedBidomainMassMatrixAssembler.hpp"
39 #include "ExtendedBidomainAssembler.hpp"
40 
41 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
43 {
44  if (this->mpLinearSystem != NULL)
45  {
46  return;
47  }
49 
50  // initialise matrix-based RHS vector and matrix, and use the linear
51  // system rhs as a template
52  Vec& r_template = this->mpLinearSystem->rGetRhsVector();
53  VecDuplicate(r_template, &mVecForConstructingRhs);
54  PetscInt ownership_range_lo;
55  PetscInt ownership_range_hi;
56  VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
57  PetscInt local_size = ownership_range_hi - ownership_range_lo;
58  PetscTools::SetupMat(mMassMatrix, 3*this->mpMesh->GetNumNodes(), 3*this->mpMesh->GetNumNodes(),
59  3*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
60  local_size, local_size);
61 }
62 
63 
64 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
66  Vec currentSolution,
67  bool computeMatrix)
68 {
69  assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
70  assert(this->mpLinearSystem->rGetRhsVector() != NULL);
71  assert(currentSolution != NULL);
72 
73 
75  // set up LHS matrix (and mass matrix)
77  if(computeMatrix)
78  {
79  mpExtendedBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
80  mpExtendedBidomainAssembler->AssembleMatrix();
81 
82  // the ExtendedBidomainMassMatrixAssembler deals with the mass matrix
83  // for both bath and nonbath problems
84  assert(SPACE_DIM==ELEMENT_DIM);
85  ExtendedBidomainMassMatrixAssembler<SPACE_DIM> mass_matrix_assembler(this->mpMesh);
86  mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
87  mass_matrix_assembler.Assemble();
88 
89  this->mpLinearSystem->SwitchWriteModeLhsMatrix();
90  PetscMatTools::Finalise(mMassMatrix);
91  }
92 
93 
94  HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
95 
97  // Set up z in b=Mz
99  DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
100 
101  // get bidomain parameters
102  double Am1 = this->mpExtendedBidomainTissue->GetAmFirstCell();
103  double Am2 = this->mpExtendedBidomainTissue->GetAmSecondCell();
104  double AmGap = this->mpExtendedBidomainTissue->GetAmGap();
105  double Cm1 = this->mpExtendedBidomainTissue->GetCmFirstCell();
106  double Cm2 = this->mpExtendedBidomainTissue->GetCmSecondCell();
107 
108  // dist stripe for the current Voltage
109  DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
110  DistributedVector::Stripe distributed_current_solution_v_first_cell(distributed_current_solution, 0);
111  DistributedVector::Stripe distributed_current_solution_v_second_cell(distributed_current_solution, 1);
112  DistributedVector::Stripe distributed_current_solution_phi_e(distributed_current_solution, 2);
113 
114  // dist stripe for z
115  DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
116  DistributedVector::Stripe dist_vec_matrix_based_v_first_cell(dist_vec_matrix_based, 0);
117  DistributedVector::Stripe dist_vec_matrix_based_v_second_cell(dist_vec_matrix_based, 1);
118  DistributedVector::Stripe dist_vec_matrix_based_phi_e(dist_vec_matrix_based, 2);
119 
120  for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
121  index!= dist_vec_matrix_based.End();
122  ++index)
123  {
124  double V_first_cell = distributed_current_solution_v_first_cell[index];
125  double V_second_Cell = distributed_current_solution_v_second_cell[index];
126 
127  double i_ionic_first_cell = this->mpExtendedBidomainTissue->rGetIionicCacheReplicated()[index.Global];
128  double i_ionic_second_cell = this->mpExtendedBidomainTissue->rGetIionicCacheReplicatedSecondCell()[index.Global];
129  double intracellular_stimulus_first_cell = this->mpExtendedBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
130  double intracellular_stimulus_second_cell = this->mpExtendedBidomainTissue->rGetIntracellularStimulusCacheReplicatedSecondCell()[index.Global];
131  double extracellular_stimulus = this->mpExtendedBidomainTissue->rGetExtracellularStimulusCacheReplicated()[index.Global];
132  double g_gap = this->mpExtendedBidomainTissue->rGetGgapCacheReplicated()[index.Global];
133  double delta_t = PdeSimulationTime::GetPdeTimeStep();
134  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;
135  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;
136 
137  if (this->mpExtendedBidomainTissue->HasTheUserSuppliedExtracellularStimulus() )
138  {
139  assert((fabs(intracellular_stimulus_first_cell) < 1e-12)
140  && (fabs(intracellular_stimulus_second_cell) < 1e-12));
141 
146  dist_vec_matrix_based_phi_e[index] = -extracellular_stimulus;
147  }
148  else
149  {
150  dist_vec_matrix_based_phi_e[index] = 0.0;
151  }
152  }
153 
154 
155  dist_vec_matrix_based.Restore();
156 
158  // b = Mz
160 
161  MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
162 
163  // assembling RHS is not finished yet, as Neumann bcs are added below, but
164  // the event will be begun again inside mpExtendedBidomainAssembler->AssembleVector();
165  HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
166 
168  // apply Neumann boundary conditions
170  mpExtendedBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions); // as the BCC can change
171  mpExtendedBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
172  mpExtendedBidomainNeumannSurfaceTermAssembler->AssembleVector();
173 
174  this->mpLinearSystem->FinaliseRhsVector();
175 
176  this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
177 
178  if(computeMatrix)
179  {
180  this->mpLinearSystem->FinaliseLhsMatrix();
181  }
182  this->mpLinearSystem->FinaliseRhsVector();
183 }
184 
185 
186 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
188  bool bathSimulation,
192  : AbstractExtendedBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
193 {
194  // Tell Tissue there's no need to replicate ionic caches
195  pTissue->SetCacheReplication(false);
196  mVecForConstructingRhs = NULL;
197 
198  // create assembler
199  if(this->mBathSimulation)
200  {
201  //this->mpExtendedExtendedBidomainAssembler = new ExtendedExtendedBidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpExtendedExtendedBidomainTissue,this->mDt);
202  EXCEPTION("Bath simulations are not yet supported for extended bidomain problems");
203  }
204  else
205  {
207  }
208 
210 
211 }
212 
213 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
215 {
216  delete mpExtendedBidomainAssembler;
217  delete mpExtendedBidomainNeumannSurfaceTermAssembler;
218 
219  if(mVecForConstructingRhs)
220  {
221  PetscTools::Destroy(mVecForConstructingRhs);
222  PetscTools::Destroy(mMassMatrix);
223  }
224 }
225 
227 // explicit instantiation
229 
230 template class ExtendedBidomainSolver<1,1>;
231 template class ExtendedBidomainSolver<2,2>;
232 template class ExtendedBidomainSolver<3,3>;
233 
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
ExtendedBidomainNeumannSurfaceTermAssembler< ELEM_DIM, SPACE_DIM > * mpExtendedBidomainNeumannSurfaceTermAssembler
ExtendedBidomainSolver(bool bathSimulation, AbstractTetrahedralMesh< ELEM_DIM, SPACE_DIM > *pMesh, ExtendedBidomainTissue< SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, 3 > *pBoundaryConditions)
ExtendedBidomainAssembler< ELEM_DIM, SPACE_DIM > * mpExtendedBidomainAssembler
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
#define EXCEPTION(message)
Definition: Exception.hpp:143
void SetupLinearSystem(Vec currentSolution, bool computeMatrix)
static void SetupMat(Mat &rMat, int numRows, int numColumns, unsigned rowPreallocation, int numLocalRows=PETSC_DECIDE, int numLocalColumns=PETSC_DECIDE, bool ignoreOffProcEntries=true, bool newAllocationError=true)
Definition: PetscTools.cpp:268
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
static double GetPdeTimeStep()
void SetCacheReplication(bool doCacheReplication)
static void Finalise(Mat matrix)
void SetMatrixToAssemble(Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
void InitialiseForSolve(Vec initialSolution)