MatrixBasedMonodomainSolver.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 "MatrixBasedMonodomainSolver.hpp"
00030 #include "MassMatrixAssembler.hpp"
00031 #include "PetscMatTools.hpp"
00032 
00033 
00034 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00035 void MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(Vec currentSolution, bool computeMatrix)
00036 {
00037     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00038     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00039 
00040     if(!this->mpMonodomainAssembler)
00041     {
00042         this->mpMonodomainAssembler = new MonodomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00043     }        
00044 
00045 
00047     // set up LHS matrix (and mass matrix)
00049     if(computeMatrix)
00050     {
00051         this->mpMonodomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00052         this->mpMonodomainAssembler->AssembleMatrix();
00053 
00054         MassMatrixAssembler<ELEMENT_DIM,SPACE_DIM> mass_matrix_assembler(this->mpMesh, HeartConfig::Instance()->GetUseMassLumping());
00055         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00056         mass_matrix_assembler.Assemble();
00057 
00058         this->mpLinearSystem->AssembleFinalLhsMatrix();
00059         PetscMatTools::AssembleFinal(mMassMatrix);
00060         
00061         if (HeartConfig::Instance()->GetUseMassLumpingForPrecond() && !HeartConfig::Instance()->GetUseMassLumping())
00062         {
00063             this->mpLinearSystem->SetPrecondMatrixIsDifferentFromLhs();
00064             
00065             MonodomainAssembler<ELEMENT_DIM,SPACE_DIM> lumped_mass_assembler(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);            
00066             lumped_mass_assembler.SetMatrixToAssemble(this->mpLinearSystem->rGetPrecondMatrix());
00067 
00068             HeartConfig::Instance()->SetUseMassLumping(true);
00069             lumped_mass_assembler.AssembleMatrix();
00070             HeartConfig::Instance()->SetUseMassLumping(false);
00071             
00072             this->mpLinearSystem->AssembleFinalPrecondMatrix();
00073         }
00074         
00075     }
00076 
00077     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00078 
00080     // Set up z in b=Mz
00082     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00083     // dist stripe for the current Voltage
00084     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00085     // dist stripe for z (return value)
00086     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00087 
00088     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00089     double Cm  = HeartConfig::Instance()->GetCapacitance();
00090 
00091     for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00092          index!= dist_vec_matrix_based.End();
00093          ++index)
00094     {
00095         double V = distributed_current_solution[index];
00096         double F = - Am*this->mpMonodomainTissue->rGetIionicCacheReplicated()[index.Global]
00097                    - this->mpMonodomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00098 
00099         dist_vec_matrix_based[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00100     }
00101     dist_vec_matrix_based.Restore();
00102     
00104     // b = Mz
00106     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00107 
00108     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00109     // the event will be begun again inside this->mpMonodomainAssembler->AssembleVector();
00110     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00111 
00113     // apply Neumann boundary conditions
00115     this->mpMonodomainAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00116     this->mpMonodomainAssembler->SetApplyNeummanBoundaryConditionsToVector(this->mpBoundaryConditions);
00117     this->mpMonodomainAssembler->OnlyAssembleOnSurfaceElements();
00118     // note: don't need this for neumann bcs, would introduce parallel replication overhead
00119     //this->mpMonodomainAssembler->SetCurrentSolution(currentSolution);
00120     this->mpMonodomainAssembler->AssembleVector();
00121   
00123     // apply correction term
00125     if(mpMonodomainCorrectionTermAssembler)
00126     {
00127         mpMonodomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00128         // don't need to set current solution
00129         mpMonodomainCorrectionTermAssembler->AssembleVector();
00130     }
00131   
00132     // finalise 
00133     this->mpLinearSystem->AssembleRhsVector();
00134 }
00135 
00136 
00137 
00138 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00139 void MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00140 {
00141     if (this->mpLinearSystem != NULL)
00142     {
00143         return;
00144     }
00145     AbstractMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00146 
00147     // initialise matrix-based RHS vector and matrix, and use the linear
00148     // system rhs as a template
00149     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00150     VecDuplicate(r_template, &mVecForConstructingRhs);
00151     PetscInt ownership_range_lo;
00152     PetscInt ownership_range_hi;
00153     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00154     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00155     PetscTools::SetupMat(mMassMatrix, this->mpMesh->GetNumNodes(), this->mpMesh->GetNumNodes(),
00156                          this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00157                          local_size, local_size);
00158 }
00159 
00160 
00161 
00162 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00163 MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::MatrixBasedMonodomainSolver(
00164             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00165             MonodomainTissue<ELEMENT_DIM,SPACE_DIM>* pTissue,
00166             BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,1>* pBoundaryConditions,
00167             unsigned numQuadPoints)
00168     : AbstractMonodomainSolver<ELEMENT_DIM,SPACE_DIM>(pMesh,pTissue,pBoundaryConditions,numQuadPoints)
00169 {
00170     // Tell tissue there's no need to replicate ionic caches
00171     pTissue->SetCacheReplication(false);
00172     mVecForConstructingRhs = NULL;
00173     
00174     if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00175     {
00176         mpMonodomainCorrectionTermAssembler 
00177             = new MonodomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpMonodomainTissue,this->mNumQuadPoints);
00178         //We are going to need those caches after all
00179         pTissue->SetCacheReplication(true);
00180     }
00181     else
00182     {
00183         mpMonodomainCorrectionTermAssembler = NULL;
00184     }
00185 }
00186 
00187 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00188 MatrixBasedMonodomainSolver<ELEMENT_DIM,SPACE_DIM>::~MatrixBasedMonodomainSolver()
00189 {
00190     if(mVecForConstructingRhs)
00191     {
00192         VecDestroy(mVecForConstructingRhs);
00193         MatDestroy(mMassMatrix);
00194     }
00195     
00196     if(mpMonodomainCorrectionTermAssembler)
00197     {
00198         delete mpMonodomainCorrectionTermAssembler;
00199     }
00200 }
00201 
00202 
00204 // explicit instantiation
00206 
00207 template class MatrixBasedMonodomainSolver<1,1>;
00208 template class MatrixBasedMonodomainSolver<1,2>;
00209 template class MatrixBasedMonodomainSolver<1,3>;
00210 template class MatrixBasedMonodomainSolver<2,2>;
00211 template class MatrixBasedMonodomainSolver<3,3>;
00212 

Generated on Mon Apr 18 11:35:28 2011 for Chaste by  doxygen 1.5.5