MonodomainSolver.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include "MonodomainSolver.hpp"
00030 #include "MassMatrixAssembler.hpp"
00031 #include "PetscMatTools.hpp"
00032 
00033 
00034 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00035 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00036 {
00037     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00038     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00039 
00040 
00042     // set up LHS matrix (and mass matrix)
00044     if(computeMatrix)
00045     {
00046         mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00047         mpMonodomainAssembler->AssembleMatrix();
00048 
00049         MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00050         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00051         mass_matrix_assembler.Assemble();
00052 
00053         this->mpLinearSystem->FinaliseLhsMatrix();
00054         PetscMatTools::Finalise(mMassMatrix);
00055 
00056         if (HeartConfig::Instance()->GetUseMassLumpingForPrecond() && !HeartConfig::Instance()->GetUseMassLumping())
00057         {
00058             this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
00059 
00060             MonodomainAssembler<ELEMENT_DIM,SPACE_DIM> lumped_mass_assembler(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00061             lumped_mass_assembler.SetMatrixToAssemble(this->mpLinearSystem->rGetPrecondMatrix());
00062 
00063             HeartConfig::Instance()->SetUseMassLumping(true);
00064             lumped_mass_assembler.AssembleMatrix();
00065             HeartConfig::Instance()->SetUseMassLumping(false);
00066 
00067             this->mpLinearSystem->FinalisePrecondMatrix();
00068         }
00069 
00070     }
00071 
00072     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00073 
00075     // Set up z in b=Mz
00077     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00078     // dist stripe for the current Voltage
00079     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00080     // dist stripe for z (return value)
00081     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00082 
00083     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00084     double Cm  = HeartConfig::Instance()->GetCapacitance();
00085 
00086     for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00087          index!= dist_vec_matrix_based.End();
00088          ++index)
00089     {
00090         double V = distributed_current_solution[index];
00091         double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00092                    - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00093 
00094         dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00095     }
00096     dist_vec_matrix_based.Restore();
00097 
00099     // b = Mz
00101     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00102 
00103     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00104     // the event will be begun again inside mpMonodomainAssembler->AssembleVector();
00105     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00106 
00108     // apply Neumann boundary conditions
00110     mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00111     mpNeumannSurfaceTermsAssembler->AssembleVector();
00112 
00114     // apply correction term
00116     if(mpMonodomainCorrectionTermAssembler)
00117     {
00118         mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00119         // don't need to set current solution
00120         mpMonodomainCorrectionTermAssembler->AssembleVector();
00121     }
00122 
00123     // finalise
00124     this->mpLinearSystem->FinaliseRhsVector();
00125 }
00126 
00127 
00128 
00129 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00130 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00131 {
00132     if (this->mpLinearSystem != NULL)
00133     {
00134         return;
00135     }
00136 
00137     // call base class version...
00138     AbstractLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>::InitialiseForSolve(initialSolution);
00139 
00140     //..then do a bit extra
00141     if(HeartConfig::Instance()->GetUseAbsoluteTolerance())
00142     {
00143         this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance());
00144     }
00145     else
00146     {
00147         this->mpLinearSystem->SetRelativeTolerance(HeartConfig::Instance()->GetRelativeTolerance());
00148     }
00149 
00150     this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00151     this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00152     this->mpLinearSystem->SetMatrixIsSymmetric(true);
00153     this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
00154 
00155     // initialise matrix-based RHS vector and matrix, and use the linear
00156     // system rhs as a template
00157     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00158     VecDuplicate(r_template, &mVecForConstructingRhs);
00159     PetscInt ownership_range_lo;
00160     PetscInt ownership_range_hi;
00161     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00162     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00163     PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00164                          this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00165                          local_size, local_size);
00166 }
00167 
00168 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00169 void MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution)
00170 {
00171     // solve cell models
00172     double time = PdeSimulationTime::GetTime();
00173     double dt = PdeSimulationTime::GetPdeTimeStep();
00174     mpMonodomainTissue->SolveCellSystems(currentSolution, time, time+dt);
00175 }
00176 
00177 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00178 MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::MonodomainSolver(
00179             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00180             MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00181             BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00182             unsigned numQuadPoints)
00183     : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
00184       mpMonodomainTissue(pTissue),
00185       mNumQuadPoints(numQuadPoints),
00186       mpBoundaryConditions(pBoundaryConditions)
00187 {
00188     assert(pTissue);
00189     assert(pBoundaryConditions);
00190     this->mMatrixIsConstant = true;
00191 
00192     mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00193     mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>(pMesh,pBoundaryConditions);
00194 
00195 
00196     // Tell tissue there's no need to replicate ionic caches
00197     pTissue->SetCacheReplication(false);
00198     mVecForConstructingRhs = NULL;
00199 
00200     if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00201     {
00202         mpMonodomainCorrectionTermAssembler
00203             = new MonodomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00204         //We are going to need those caches after all
00205         pTissue->SetCacheReplication(true);
00206     }
00207     else
00208     {
00209         mpMonodomainCorrectionTermAssembler = NULL;
00210     }
00211 }
00212 
00213 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00214 MonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~MonodomainSolver()
00215 {
00216     delete mpMonodomainAssembler;
00217     delete mpNeumannSurfaceTermsAssembler;
00218 
00219     if(mVecForConstructingRhs)
00220     {
00221         VecDestroy(mVecForConstructingRhs);
00222         MatDestroy(mMassMatrix);
00223     }
00224 
00225     if(mpMonodomainCorrectionTermAssembler)
00226     {
00227         delete mpMonodomainCorrectionTermAssembler;
00228     }
00229 }
00230 
00231 
00233 // explicit instantiation
00235 
00236 template class MonodomainSolver<1,1>;
00237 template class MonodomainSolver<1,2>;
00238 template class MonodomainSolver<1,3>;
00239 template class MonodomainSolver<2,2>;
00240 template class MonodomainSolver<3,3>;
00241 
Generated on Thu Dec 22 13:00:06 2011 for Chaste by  doxygen 1.6.3