BidomainSolver.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 
00030 #include "BidomainSolver.hpp"
00031 #include "BidomainAssembler.hpp"
00032 #include "BidomainWithBathAssembler.hpp"
00033 #include "PetscMatTools.hpp"
00034 
00035 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00036 void BidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00037 {
00038     if (this->mpLinearSystem != NULL)
00039     {
00040         return;
00041     }
00042     AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00043 
00044     // initialise matrix-based RHS vector and matrix, and use the linear
00045     // system rhs as a template
00046     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00047     VecDuplicate(r_template, &mVecForConstructingRhs);
00048     PetscInt ownership_range_lo;
00049     PetscInt ownership_range_hi;
00050     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00051     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00052     PetscTools::SetupMat(mMassMatrix, 2*this->mpMesh->GetNumNodes(), 2*this->mpMesh->GetNumNodes(),
00053                          2*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00054                          local_size, local_size);
00055 }
00056 
00057 
00058 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00059 void BidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(
00060         Vec currentSolution,
00061         bool computeMatrix)
00062 {
00063     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00064     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00065     assert(currentSolution != NULL);
00066 
00067 
00069     // set up LHS matrix (and mass matrix)
00071     if (computeMatrix)
00072     {
00073         mpBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00074         mpBidomainAssembler->AssembleMatrix();
00075 
00076         // the BidomainMassMatrixAssembler deals with the mass matrix
00077         // for both bath and nonbath problems
00078         assert(SPACE_DIM==ELEMENT_DIM);
00079         BidomainMassMatrixAssembler<SPACE_DIM> mass_matrix_assembler(this->mpMesh);
00080         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00081         mass_matrix_assembler.Assemble();
00082 
00083         this->mpLinearSystem->SwitchWriteModeLhsMatrix();
00084         PetscMatTools::Finalise(mMassMatrix);
00085     }
00086 
00087 
00088     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00089 
00091     // Set up z in b=Mz
00093     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00094 
00095     // dist stripe for the current Voltage
00096     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00097     DistributedVector::Stripe distributed_current_solution_vm(distributed_current_solution, 0);
00098 
00099     // dist stripe for z
00100     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00101     DistributedVector::Stripe dist_vec_matrix_based_vm(dist_vec_matrix_based, 0);
00102     DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
00103 
00104     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00105     double Cm  = HeartConfig::Instance()->GetCapacitance();
00106 
00107     if(!(this->mBathSimulation))
00108     {
00109         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00110              index!= dist_vec_matrix_based.End();
00111              ++index)
00112         {
00113             double V = distributed_current_solution_vm[index];
00114             double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00115                        - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00116 
00117             dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00118             dist_vec_matrix_based_phie[index] = 0.0;
00119         }
00120     }
00121     else
00122     {
00123         DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
00124 
00125         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00126              index!= dist_vec_matrix_based.End();
00127              ++index)
00128         {
00129 
00130             if ( !HeartRegionCode::IsRegionBath( this->mpMesh->GetNode(index.Global)->GetRegion() ))
00131             {
00132                 double V = distributed_current_solution_vm[index];
00133                 double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00134                            - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00135 
00136                 dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00137             }
00138             else
00139             {
00140                 dist_vec_matrix_based_vm[index] = 0.0;
00141             }
00142 
00143             dist_vec_matrix_based_phie[index] = 0.0;
00144 
00145         }
00146     }
00147 
00148     dist_vec_matrix_based.Restore();
00149 
00151     // b = Mz
00153     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00154 
00155     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00156     // the event will be begun again inside mpBidomainAssembler->AssembleVector();
00157     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00158 
00159 
00161     // apply Neumann boundary conditions
00163     mpBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions); // as the BCC can change
00164     mpBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00165     mpBidomainNeumannSurfaceTermAssembler->AssembleVector();
00166 
00167 
00169     // apply correction term
00171     if(mpBidomainCorrectionTermAssembler)
00172     {
00173         mpBidomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00174         // don't need to set current solution
00175         mpBidomainCorrectionTermAssembler->AssembleVector();
00176     }
00177 
00178     this->mpLinearSystem->FinaliseRhsVector();
00179 
00180     this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
00181 
00182     if(this->mBathSimulation)
00183     {
00184         this->mpLinearSystem->FinaliseLhsMatrix();
00185         this->FinaliseForBath(computeMatrix,true);
00186     }
00187 
00188     if(computeMatrix)
00189     {
00190         this->mpLinearSystem->FinaliseLhsMatrix();
00191     }
00192     this->mpLinearSystem->FinaliseRhsVector();
00193 }
00194 
00195 
00196 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00197 BidomainSolver<ELEMENT_DIM,SPACE_DIM>::BidomainSolver(
00198         bool bathSimulation,
00199         AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00200         BidomainTissue<SPACE_DIM>* pTissue,
00201         BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,2>* pBoundaryConditions,
00202         unsigned numQuadPoints)
00203     : AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
00204 {
00205     // Tell tissue there's no need to replicate ionic caches
00206     pTissue->SetCacheReplication(false);
00207     mVecForConstructingRhs = NULL;
00208 
00209     // create assembler
00210     if(bathSimulation)
00211     {
00212         mpBidomainAssembler = new BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mNumQuadPoints);
00213     }
00214     else
00215     {
00216         mpBidomainAssembler = new BidomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mNumQuadPoints);
00217     }
00218 
00219 
00220     mpBidomainNeumannSurfaceTermAssembler = new BidomainNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM>(pMesh,pBoundaryConditions);
00221 
00222     if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00223     {
00224         mpBidomainCorrectionTermAssembler
00225             = new BidomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mNumQuadPoints);
00226         //We are going to need those caches after all
00227         pTissue->SetCacheReplication(true);
00228     }
00229     else
00230     {
00231         mpBidomainCorrectionTermAssembler = NULL;
00232     }
00233 }
00234 
00235 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00236 BidomainSolver<ELEMENT_DIM,SPACE_DIM>::~BidomainSolver()
00237 {
00238     delete mpBidomainAssembler;
00239     delete mpBidomainNeumannSurfaceTermAssembler;
00240 
00241     if(mVecForConstructingRhs)
00242     {
00243         VecDestroy(mVecForConstructingRhs);
00244         MatDestroy(mMassMatrix);
00245     }
00246 
00247     if(mpBidomainCorrectionTermAssembler)
00248     {
00249         delete mpBidomainCorrectionTermAssembler;
00250     }
00251 }
00252 
00254 // explicit instantiation
00256 
00257 template class BidomainSolver<1,1>;
00258 template class BidomainSolver<2,2>;
00259 template class BidomainSolver<3,3>;
00260 
Generated on Thu Dec 22 13:00:06 2011 for Chaste by  doxygen 1.6.3