Chaste  Release::3.4
BidomainSolver.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 "BidomainSolver.hpp"
38 #include "BidomainAssembler.hpp"
39 #include "BidomainWithBathAssembler.hpp"
40 #include "PetscMatTools.hpp"
41 
42 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
44 {
45  if (this->mpLinearSystem != NULL)
46  {
47  return;
48  }
50 
51  // initialise matrix-based RHS vector and matrix, and use the linear
52  // system rhs as a template
53  Vec& r_template = this->mpLinearSystem->rGetRhsVector();
54  VecDuplicate(r_template, &mVecForConstructingRhs);
55  PetscInt ownership_range_lo;
56  PetscInt ownership_range_hi;
57  VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
58  PetscInt local_size = ownership_range_hi - ownership_range_lo;
59  PetscTools::SetupMat(mMassMatrix, 2*this->mpMesh->GetNumNodes(), 2*this->mpMesh->GetNumNodes(),
60  2*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
61  local_size, local_size);
62 }
63 
64 
65 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
67  Vec currentSolution,
68  bool computeMatrix)
69 {
70  assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
71  assert(this->mpLinearSystem->rGetRhsVector() != NULL);
72  assert(currentSolution != NULL);
73 
74 
76  // set up LHS matrix (and mass matrix)
78  if (computeMatrix)
79  {
80  mpBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
81  mpBidomainAssembler->AssembleMatrix();
82 
83  // the BidomainMassMatrixAssembler deals with the mass matrix
84  // for both bath and nonbath problems
85  assert(SPACE_DIM==ELEMENT_DIM);
86  BidomainMassMatrixAssembler<SPACE_DIM> mass_matrix_assembler(this->mpMesh);
87  mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
88  mass_matrix_assembler.Assemble();
89 
90  this->mpLinearSystem->SwitchWriteModeLhsMatrix();
91  PetscMatTools::Finalise(mMassMatrix);
92  }
93 
94 
95  HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
96 
98  // Set up z in b=Mz
100  DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
101 
102  // dist stripe for the current Voltage
103  DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
104  DistributedVector::Stripe distributed_current_solution_vm(distributed_current_solution, 0);
105 
106  // dist stripe for z
107  DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
108  DistributedVector::Stripe dist_vec_matrix_based_vm(dist_vec_matrix_based, 0);
109  DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
110 
112  double Cm = HeartConfig::Instance()->GetCapacitance();
113 
114  if(!(this->mBathSimulation))
115  {
116  for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
117  index!= dist_vec_matrix_based.End();
118  ++index)
119  {
120  double V = distributed_current_solution_vm[index];
121  double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
122  - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
123 
124  dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
125  dist_vec_matrix_based_phie[index] = 0.0;
126  }
127  }
128  else
129  {
130  for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
131  index!= dist_vec_matrix_based.End();
132  ++index)
133  {
134 
135  if ( !HeartRegionCode::IsRegionBath( this->mpMesh->GetNode(index.Global)->GetRegion() ))
136  {
137  double V = distributed_current_solution_vm[index];
138  double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
139  - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
140 
141  dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
142  }
143  else
144  {
145  dist_vec_matrix_based_vm[index] = 0.0;
146  }
147 
148  dist_vec_matrix_based_phie[index] = 0.0;
149 
150  }
151  }
152 
153  dist_vec_matrix_based.Restore();
154 
156  // b = Mz
158  MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
159 
160  // assembling RHS is not finished yet, as Neumann bcs are added below, but
161  // the event will be begun again inside mpBidomainAssembler->AssembleVector();
162  HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
163 
164 
166  // apply Neumann boundary conditions
168  mpBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions); // as the BCC can change
169  mpBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
170  mpBidomainNeumannSurfaceTermAssembler->AssembleVector();
171 
172 
174  // apply correction term
176  if(mpBidomainCorrectionTermAssembler)
177  {
178  mpBidomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
179  // don't need to set current solution
180  mpBidomainCorrectionTermAssembler->AssembleVector();
181  }
182 
183  this->mpLinearSystem->FinaliseRhsVector();
184 
185  this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
186 
187  if(this->mBathSimulation)
188  {
189  this->mpLinearSystem->FinaliseLhsMatrix();
190  this->FinaliseForBath(computeMatrix,true);
191  }
192 
193  if(computeMatrix)
194  {
195  this->mpLinearSystem->FinaliseLhsMatrix();
196  }
197  this->mpLinearSystem->FinaliseRhsVector();
198 }
199 
200 
201 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
203  bool bathSimulation,
205  BidomainTissue<SPACE_DIM>* pTissue,
207  : AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
208 {
209  // Tell tissue there's no need to replicate ionic caches
210  pTissue->SetCacheReplication(false);
211  mVecForConstructingRhs = NULL;
212 
213  // create assembler
214  if(bathSimulation)
215  {
217  }
218  else
219  {
221  }
222 
223 
225 
227  {
230  //We are going to need those caches after all
231  pTissue->SetCacheReplication(true);
232  }
233  else
234  {
236  }
237 }
238 
239 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
241 {
242  delete mpBidomainAssembler;
243  delete mpBidomainNeumannSurfaceTermAssembler;
244 
245  if(mVecForConstructingRhs)
246  {
247  PetscTools::Destroy(mVecForConstructingRhs);
248  PetscTools::Destroy(mMassMatrix);
249  }
250 
251  if(mpBidomainCorrectionTermAssembler)
252  {
253  delete mpBidomainCorrectionTermAssembler;
254  }
255 }
256 
258 // explicit instantiation
260 
261 template class BidomainSolver<1,1>;
262 template class BidomainSolver<2,2>;
263 template class BidomainSolver<3,3>;
264 
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
static bool IsRegionBath(HeartRegionType regionId)
BidomainAssembler< ELEMENT_DIM, SPACE_DIM > * mpBidomainAssembler
void InitialiseForSolve(Vec initialSolution)
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
BidomainNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM > * mpBidomainNeumannSurfaceTermAssembler
BidomainTissue< SPACE_DIM > * mpBidomainTissue
BidomainCorrectionTermAssembler< ELEMENT_DIM, SPACE_DIM > * mpBidomainCorrectionTermAssembler
bool GetUseStateVariableInterpolation() const
void SetupLinearSystem(Vec currentSolution, bool computeMatrix)
BidomainSolver(bool bathSimulation, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BidomainTissue< SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 2 > *pBoundaryConditions)
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 double GetPdeTimeStepInverse()
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
void SetCacheReplication(bool doCacheReplication)
static void Finalise(Mat matrix)
double GetCapacitance() const
double GetSurfaceAreaToVolumeRatio() const
void InitialiseForSolve(Vec initialSolution)
void SetMatrixToAssemble(Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
static HeartConfig * Instance()