OperatorSplittingMonodomainSolver.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 "OperatorSplittingMonodomainSolver.hpp"
00030 
00031 
00032 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00033 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00034 {
00035     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00036     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00037 
00039     // set up LHS matrix (and mass matrix)
00041     if(computeMatrix)
00042     {
00043         mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00044         mpMonodomainAssembler->AssembleMatrix();
00045 
00046         MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00047         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00048         mass_matrix_assembler.Assemble();
00049 
00050         this->mpLinearSystem->FinaliseLhsMatrix();
00051         PetscMatTools::Finalise(mMassMatrix);
00052     }
00053 
00054     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00055 
00057     // Set up z in b=Mz
00059     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00060     // dist stripe for the current Voltage
00061     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00062     // dist stripe for z (return value)
00063     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00064 
00065     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00066     double Cm  = HeartConfig::Instance()->GetCapacitance();
00067 
00068     for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00069          index!= dist_vec_matrix_based.End();
00070          ++index)
00071     {
00072         double V = distributed_current_solution[index];
00073         // in the main solver, the nodal ionic current and stimuli is computed and used.
00074         // However in operator splitting, this part of the solve is diffusion only, no reaction terms
00075         //double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00076         //           - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00077 
00078         dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse();
00079     }
00080     dist_vec_matrix_based.Restore();
00081 
00083     // b = Mz
00085     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00086 
00087     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00088     // the event will be begun again inside mpMonodomainAssembler->AssembleVector();
00089     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00090 
00092     // apply Neumann boundary conditions
00094     mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00095     mpNeumannSurfaceTermsAssembler->AssembleVector();
00096 
00097     // finalise
00098     this->mpLinearSystem->FinaliseRhsVector();
00099 }
00100 
00101 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00102 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec currentSolution)
00103 {
00104     double time = PdeSimulationTime::GetTime();
00105     double dt = PdeSimulationTime::GetPdeTimeStep();
00106     mpMonodomainTissue->SolveCellSystems(currentSolution, time, time+dt/2.0, true);
00107 }
00108 
00109 
00110 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00111 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::FollowingSolveLinearSystem(Vec currentSolution)
00112 {
00113     // solve cell models for second half timestep
00114     double time = PdeSimulationTime::GetTime();
00115     double dt = PdeSimulationTime::GetPdeTimeStep();
00116     mpMonodomainTissue->SolveCellSystems(currentSolution, time + dt/2, time+dt, true);
00117 }
00118 
00119 
00120 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00121 void OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00122 {
00123     if (this->mpLinearSystem != NULL)
00124     {
00125         return;
00126     }
00127 
00128     // call base class version...
00129     AbstractLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>::InitialiseForSolve(initialSolution);
00130 
00131     //..then do a bit extra
00132     if(HeartConfig::Instance()->GetUseAbsoluteTolerance())
00133     {
00134         this->mpLinearSystem->SetAbsoluteTolerance(HeartConfig::Instance()->GetAbsoluteTolerance());
00135     }
00136     else
00137     {
00138         NEVER_REACHED;
00139         // re-implement when needed
00140         //this->mpLinearSystem->SetRelativeTolerance(HeartConfig::Instance()->GetRelativeTolerance());
00141     }
00142 
00143     this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00144     this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00145     this->mpLinearSystem->SetMatrixIsSymmetric(true);
00146     this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
00147 
00148     // initialise matrix-based RHS vector and matrix, and use the linear
00149     // system rhs as a template
00150     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00151     VecDuplicate(r_template, &mVecForConstructingRhs);
00152     PetscInt ownership_range_lo;
00153     PetscInt ownership_range_hi;
00154     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00155     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00156     PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00157                          this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00158                          local_size, local_size);
00159 }
00160 
00161 
00162 
00163 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00164 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::OperatorSplittingMonodomainSolver(
00165             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00166             MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00167             BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00168             unsigned numQuadPoints)
00169     : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,1>(pMesh),
00170       mpBoundaryConditions(pBoundaryConditions),
00171       mpMonodomainTissue(pTissue),
00172       mNumQuadPoints(numQuadPoints)
00173 {
00174     assert(pTissue);
00175     assert(pBoundaryConditions);
00176     this->mMatrixIsConstant = true;
00177 
00178     mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00179     mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,1>(pMesh,pBoundaryConditions);
00180 
00181     // Tell tissue there's no need to replicate ionic caches
00182     pTissue->SetCacheReplication(false);
00183     mVecForConstructingRhs = NULL;
00184 }
00185 
00186 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00187 OperatorSplittingMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~OperatorSplittingMonodomainSolver()
00188 {
00189     delete mpMonodomainAssembler;
00190     delete mpNeumannSurfaceTermsAssembler;
00191 
00192     if(mVecForConstructingRhs)
00193     {
00194         VecDestroy(mVecForConstructingRhs);
00195         MatDestroy(mMassMatrix);
00196     }
00197 }
00198 
00199 
00200 
00202 // explicit instantiation
00204 
00205 template class OperatorSplittingMonodomainSolver<1,1>;
00206 template class OperatorSplittingMonodomainSolver<1,2>;
00207 template class OperatorSplittingMonodomainSolver<1,3>;
00208 template class OperatorSplittingMonodomainSolver<2,2>;
00209 template class OperatorSplittingMonodomainSolver<3,3>;
00210 
Generated on Thu Dec 22 13:00:06 2011 for Chaste by  doxygen 1.6.3