Chaste Release::3.1
MonodomainSolver.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include "MonodomainSolver.hpp"
00037 #include "MassMatrixAssembler.hpp"
00038 #include "PetscMatTools.hpp"
00039 
00040 
00041 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00042 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00043 {
00044     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00045     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00046 
00047 
00049     // set up LHS matrix (and mass matrix)
00051     if(computeMatrix)
00052     {
00053         mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00054         mpMonodomainAssembler->AssembleMatrix();
00055 
00056         MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00057         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00058         mass_matrix_assembler.Assemble();
00059 
00060         this->mpLinearSystem->FinaliseLhsMatrix();
00061         PetscMatTools::Finalise(mMassMatrix);
00062 
00063         if (HeartConfig::Instance()->GetUseMassLumpingForPrecond() && !HeartConfig::Instance()->GetUseMassLumping())
00064         {
00065             this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
00066 
00067             MonodomainAssembler<ELEMENT_DIM,SPACE_DIM> lumped_mass_assembler(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00068             lumped_mass_assembler.SetMatrixToAssemble(this->mpLinearSystem->rGetPrecondMatrix());
00069 
00070             HeartConfig::Instance()->SetUseMassLumping(true);
00071             lumped_mass_assembler.AssembleMatrix();
00072             HeartConfig::Instance()->SetUseMassLumping(false);
00073 
00074             this->mpLinearSystem->FinalisePrecondMatrix();
00075         }
00076 
00077     }
00078 
00079     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00080 
00082     // Set up z in b=Mz
00084     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00085     // dist stripe for the current Voltage
00086     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00087     // dist stripe for z (return value)
00088     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00089 
00090     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00091     double Cm  = HeartConfig::Instance()->GetCapacitance();
00092 
00093     for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00094          index!= dist_vec_matrix_based.End();
00095          ++index)
00096     {
00097         double V = distributed_current_solution[index];
00098         double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00099                    - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00100 
00101         dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00102     }
00103     dist_vec_matrix_based.Restore();
00104 
00106     // b = Mz
00108     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00109 
00110     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00111     // the event will be begun again inside mpMonodomainAssembler->AssembleVector();
00112     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00113 
00115     // apply Neumann boundary conditions
00117     mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00118     mpNeumannSurfaceTermsAssembler->AssembleVector();
00119 
00121     // apply correction term
00123     if(mpMonodomainCorrectionTermAssembler)
00124     {
00125         mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00126         // don't need to set current solution
00127         mpMonodomainCorrectionTermAssembler->AssembleVector();
00128     }
00129 
00130     // finalise
00131     this->mpLinearSystem->FinaliseRhsVector();
00132 }
00133 
00134 
00135 
00136 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00137 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00138 {
00139     if (this->mpLinearSystem != NULL)
00140     {
00141         return;
00142     }
00143 
00144     // call base class version...
00145     AbstractLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>::InitialiseForSolve(initialSolution);
00146 
00147     //..then do a bit extra
00148     if(HeartConfig::Instance()->GetUseAbsoluteTolerance())
00149     {
00150         this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance());
00151     }
00152     else
00153     {
00154         this->mpLinearSystem->SetRelativeTolerance(HeartConfig::Instance()->GetRelativeTolerance());
00155     }
00156 
00157     this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00158     this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00159     this->mpLinearSystem->SetMatrixIsSymmetric(true);
00160     this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
00161 
00162     // initialise matrix-based RHS vector and matrix, and use the linear
00163     // system rhs as a template
00164     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00165     VecDuplicate(r_template, &mVecForConstructingRhs);
00166     PetscInt ownership_range_lo;
00167     PetscInt ownership_range_hi;
00168     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00169     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00170     PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00171                          this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00172                          local_size, local_size);
00173 }
00174 
00175 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00176 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution)
00177 {
00178     // solve cell models
00179     double time = PdeSimulationTime::GetTime();
00180     double dt = PdeSimulationTime::GetPdeTimeStep();
00181     mpMonodomainTissue->SolveCellSystems(currentSolution, time, time+dt);
00182 }
00183 
00184 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00185 MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::MonodomainSolver(
00186             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00187             MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00188             BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00189             unsigned numQuadPoints)
00190     : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
00191       mpMonodomainTissue(pTissue),
00192       mNumQuadPoints(numQuadPoints),
00193       mpBoundaryConditions(pBoundaryConditions)
00194 {
00195     assert(pTissue);
00196     assert(pBoundaryConditions);
00197     this->mMatrixIsConstant = true;
00198 
00199     mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00200     mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>(pMesh,pBoundaryConditions);
00201 
00202 
00203     // Tell tissue there's no need to replicate ionic caches
00204     pTissue->SetCacheReplication(false);
00205     mVecForConstructingRhs = NULL;
00206 
00207     if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00208     {
00209         mpMonodomainCorrectionTermAssembler
00210             = new MonodomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00211         //We are going to need those caches after all
00212         pTissue->SetCacheReplication(true);
00213     }
00214     else
00215     {
00216         mpMonodomainCorrectionTermAssembler = NULL;
00217     }
00218 }
00219 
00220 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00221 MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~MonodomainSolver()
00222 {
00223     delete mpMonodomainAssembler;
00224     delete mpNeumannSurfaceTermsAssembler;
00225 
00226     if(mVecForConstructingRhs)
00227     {
00228         PetscTools::Destroy(mVecForConstructingRhs);
00229         PetscTools::Destroy(mMassMatrix);
00230     }
00231 
00232     if(mpMonodomainCorrectionTermAssembler)
00233     {
00234         delete mpMonodomainCorrectionTermAssembler;
00235     }
00236 }
00237 
00238 
00240 // explicit instantiation
00242 
00243 template class MonodomainSolver<1,1>;
00244 template class MonodomainSolver<1,2>;
00245 template class MonodomainSolver<1,3>;
00246 template class MonodomainSolver<2,2>;
00247 template class MonodomainSolver<3,3>;
00248