Chaste  Release::3.4
MonodomainSolver.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 #include "MonodomainSolver.hpp"
37 #include "MassMatrixAssembler.hpp"
38 #include "PetscMatTools.hpp"
39 
40 
41 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
42 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
43 {
44  assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
45  assert(this->mpLinearSystem->rGetRhsVector() != NULL);
46 
47 
49  // set up LHS matrix (and mass matrix)
51  if(computeMatrix)
52  {
53  mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
54  mpMonodomainAssembler->AssembleMatrix();
55 
56  MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
57  mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
58  mass_matrix_assembler.Assemble();
59 
60  this->mpLinearSystem->FinaliseLhsMatrix();
61  PetscMatTools::Finalise(mMassMatrix);
62 
63  if (HeartConfig::Instance()->GetUseMassLumpingForPrecond() && !HeartConfig::Instance()->GetUseMassLumping())
64  {
65  this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
66 
67  MonodomainAssembler<ELEMENT_DIM,SPACE_DIM> lumped_mass_assembler(this->mpMesh,this->mpMonodomainTissue);
68  lumped_mass_assembler.SetMatrixToAssemble(this->mpLinearSystem->rGetPrecondMatrix());
69 
71  lumped_mass_assembler.AssembleMatrix();
73 
74  this->mpLinearSystem->FinalisePrecondMatrix();
75  }
76 
77  }
78 
79  HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
80 
82  // Set up z in b=Mz
84  DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
85  // dist stripe for the current Voltage
86  DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
87  // dist stripe for z (return value)
88  DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
89 
91  double Cm = HeartConfig::Instance()->GetCapacitance();
92 
93  for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
94  index!= dist_vec_matrix_based.End();
95  ++index)
96  {
97  double V = distributed_current_solution[index];
98  double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
99  - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
100 
101  dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
102  }
103  dist_vec_matrix_based.Restore();
104 
106  // b = Mz
108  MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
109 
110  // assembling RHS is not finished yet, as Neumann bcs are added below, but
111  // the event will be begun again inside mpMonodomainAssembler->AssembleVector();
112  HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
113 
115  // apply Neumann boundary conditions
117  mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
118  mpNeumannSurfaceTermsAssembler->AssembleVector();
119 
121  // apply correction term
123  if(mpMonodomainCorrectionTermAssembler)
124  {
125  mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
126  // don't need to set current solution
127  mpMonodomainCorrectionTermAssembler->AssembleVector();
128  }
129 
130  // finalise
131  this->mpLinearSystem->FinaliseRhsVector();
132 }
133 
134 
135 
136 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
138 {
139  if (this->mpLinearSystem != NULL)
140  {
141  return;
142  }
143 
144  // call base class version...
146 
147  //..then do a bit extra
148  if(HeartConfig::Instance()->GetUseAbsoluteTolerance())
149  {
150  this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance());
151  }
152  else
153  {
154  this->mpLinearSystem->SetRelativeTolerance(HeartConfig::Instance()->GetRelativeTolerance());
155  }
156 
157  this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
158  this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
159  this->mpLinearSystem->SetMatrixIsSymmetric(true);
160  this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
161 
162  // initialise matrix-based RHS vector and matrix, and use the linear
163  // system rhs as a template
164  Vec& r_template = this->mpLinearSystem->rGetRhsVector();
165  VecDuplicate(r_template, &mVecForConstructingRhs);
166  PetscInt ownership_range_lo;
167  PetscInt ownership_range_hi;
168  VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
169  PetscInt local_size = ownership_range_hi - ownership_range_lo;
170  PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
171  this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
172  local_size, local_size);
173 }
174 
175 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
177 {
178  // solve cell models
179  mpMonodomainTissue->SolveCellSystems(currentSolution, PdeSimulationTime::GetTime(), PdeSimulationTime::GetNextTime());
180 }
181 
182 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
187  : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
188  mpMonodomainTissue(pTissue),
189  mpBoundaryConditions(pBoundaryConditions)
190 {
191  assert(pTissue);
192  assert(pBoundaryConditions);
193  this->mMatrixIsConstant = true;
194 
197 
198 
199  // Tell tissue there's no need to replicate ionic caches
200  pTissue->SetCacheReplication(false);
201  mVecForConstructingRhs = NULL;
202 
203  if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
204  {
207  //We are going to need those caches after all
208  pTissue->SetCacheReplication(true);
209  }
210  else
211  {
213  }
214 }
215 
216 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
218 {
219  delete mpMonodomainAssembler;
220  delete mpNeumannSurfaceTermsAssembler;
221 
222  if(mVecForConstructingRhs)
223  {
224  PetscTools::Destroy(mVecForConstructingRhs);
225  PetscTools::Destroy(mMassMatrix);
226  }
227 
228  if(mpMonodomainCorrectionTermAssembler)
229  {
230  delete mpMonodomainCorrectionTermAssembler;
231  }
232 }
233 
234 
236 // explicit instantiation
238 
239 template class MonodomainSolver<1,1>;
240 template class MonodomainSolver<1,2>;
241 template class MonodomainSolver<1,3>;
242 template class MonodomainSolver<2,2>;
243 template class MonodomainSolver<3,3>;
244 
MonodomainCorrectionTermAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainCorrectionTermAssembler
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
NaturalNeumannSurfaceTermAssembler< ELEMENT_DIM, SPACE_DIM, 1 > * mpNeumannSurfaceTermsAssembler
MonodomainAssembler< ELEMENT_DIM, SPACE_DIM > * mpMonodomainAssembler
virtual void InitialiseForSolve(Vec initialSolution=NULL)
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
void SetupLinearSystem(Vec currentSolution, bool computeMatrix)
void SetUseMassLumping(bool useMassLumping=true)
static double GetNextTime()
virtual void InitialiseForSolve(Vec initialSolution)
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
MonodomainSolver(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, MonodomainTissue< ELEMENT_DIM, SPACE_DIM > *pTissue, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, 1 > *pBoundaryConditions)
static double GetPdeTimeStepInverse()
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
void PrepareForSetupLinearSystem(Vec currentSolution)
void SetCacheReplication(bool doCacheReplication)
static void Finalise(Mat matrix)
double GetCapacitance() const
double GetSurfaceAreaToVolumeRatio() const
void SetMatrixToAssemble(Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
static HeartConfig * Instance()
MonodomainTissue< ELEMENT_DIM, SPACE_DIM > * mpMonodomainTissue
static double GetTime()