MatrixBasedBidomainSolver.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
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 "MatrixBasedBidomainSolver.hpp"
00031 #include "BidomainAssembler.hpp"
00032 #include "BidomainWithBathAssembler.hpp"
00033 
00034 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00035 void MatrixBasedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00036 {
00037     if (this->mpLinearSystem != NULL)
00038     {
00039         return;
00040     }
00041     AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00042 
00043     // initialise matrix-based RHS vector and matrix, and use the linear
00044     // system rhs as a template
00045     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00046     VecDuplicate(r_template, &mVecForConstructingRhs);
00047     PetscInt ownership_range_lo;
00048     PetscInt ownership_range_hi;
00049     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00050     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00051     PetscTools::SetupMat(mMassMatrix, 2*this->mpMesh->GetNumNodes(), 2*this->mpMesh->GetNumNodes(), 
00052                          2*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00053                          local_size, local_size); 
00054 }
00055 
00056 
00057 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00058 void MatrixBasedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(
00059         Vec currentSolution,
00060         bool computeMatrix)
00061 {
00062     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00063     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00064     assert(currentSolution != NULL);
00065 
00066     // create assembler
00067     if(!this->mpBidomainAssembler)
00068     {
00069         if(this->mBathSimulation)
00070         {
00071             this->mpBidomainAssembler = new BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mDt,this->mNumQuadPoints);
00072         }
00073         else
00074         {
00075             this->mpBidomainAssembler = new BidomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mDt,this->mNumQuadPoints);
00076         }
00077     }    
00078 
00079 
00081     // set up LHS matrix (and mass matrix)
00083     if(computeMatrix)
00084     {
00085         this->mpBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00086         this->mpBidomainAssembler->AssembleMatrix();
00087 
00088         // the BidomainMassMatrixAssembler deals with the mass matrix
00089         // for both bath and nonbath problems 
00090         assert(SPACE_DIM==ELEMENT_DIM);
00091         BidomainMassMatrixAssembler<SPACE_DIM> mass_matrix_assembler(this->mpMesh);
00092         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00093         mass_matrix_assembler.Assemble();
00094 
00095         this->mpLinearSystem->AssembleIntermediateLhsMatrix();
00096         MatAssemblyBegin(mMassMatrix, MAT_FINAL_ASSEMBLY);
00097         MatAssemblyEnd(mMassMatrix, MAT_FINAL_ASSEMBLY);
00098     }
00099 
00100 
00101     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00102 
00104     // Set up z in b=Mz
00106     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00107 
00108     // dist stripe for the current Voltage
00109     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00110     DistributedVector::Stripe distributed_current_solution_vm(distributed_current_solution, 0);
00111 
00112     // dist stripe for z
00113     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00114     DistributedVector::Stripe dist_vec_matrix_based_vm(dist_vec_matrix_based, 0);
00115     DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
00116 
00117     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00118     double Cm  = HeartConfig::Instance()->GetCapacitance();
00119 
00120     if(!(this->mBathSimulation))
00121     {
00122         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00123              index!= dist_vec_matrix_based.End();
00124              ++index)
00125         {
00126             double V = distributed_current_solution_vm[index];
00127             double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00128                        - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00129 
00130             dist_vec_matrix_based_vm[index] = Am*Cm*V*this->mDtInverse + F;
00131             dist_vec_matrix_based_phie[index] = 0.0;
00132         }
00133     }
00134     else
00135     {
00136         DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
00137 
00138         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00139              index!= dist_vec_matrix_based.End();
00140              ++index)
00141         {
00142     
00143             if (this->mpMesh->GetNode(index.Global)->GetRegion() != HeartRegionCode::BATH)
00144             {
00145                 double V = distributed_current_solution_vm[index];
00146                 double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00147                            - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00148     
00149                 dist_vec_matrix_based_vm[index] = Am*Cm*V*this->mDtInverse + F;
00150             }
00151             else
00152             {
00153                 dist_vec_matrix_based_vm[index] = 0.0;
00154             }
00155     
00156             dist_vec_matrix_based_phie[index] = 0.0;
00157     
00158         }
00159     }
00160 
00161     dist_vec_matrix_based.Restore();
00162 
00164     // b = Mz
00166     this->mpLinearSystem->ZeroRhsVector();
00167 
00168     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00169 
00170     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00171     // the event will be begun again inside this->mpBidomainAssembler->AssembleVector(); 
00172     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00173 
00174 
00176     // apply Neumann boundary conditions
00178 
00179     this->mpBidomainAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00180 
00181     // MUST be called every timestep in case the bcc has been reset, see comment in
00182     // this->ResetBoundaryConditionsContainer()
00183     this->mpBidomainAssembler->SetApplyNeummanBoundaryConditionsToVector(this->mpBoundaryConditions);
00184 
00185     this->mpBidomainAssembler->OnlyAssembleOnSurfaceElements();
00186     // note: don't need this for neumann bcs, would introduce parallel replication overhead
00187     //mpBidomainAssembler->SetCurrentSolution(currentSolution);
00188     this->mpBidomainAssembler->AssembleVector();
00189 
00190     this->mpLinearSystem->AssembleRhsVector();
00191 
00192     this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
00193 
00194     if(this->mBathSimulation)
00195     {
00196         this->FinaliseForBath(computeMatrix,true);
00197     }
00198 
00199     if(computeMatrix)
00200     {
00201         this->mpLinearSystem->AssembleFinalLhsMatrix();
00202     }
00203     this->mpLinearSystem->AssembleRhsVector();
00204 }
00205 
00206 
00207 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00208 MatrixBasedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::MatrixBasedBidomainSolver(
00209         bool bathSimulation,
00210         AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00211         BidomainTissue<SPACE_DIM>* pTissue,
00212         BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,2>* pBoundaryConditions,
00213         unsigned numQuadPoints)
00214     : AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
00215 {
00216     // Tell tissue there's no need to replicate ionic caches
00217     pTissue->SetCacheReplication(false);
00218     mVecForConstructingRhs = NULL;
00219 }
00220 
00221 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00222 MatrixBasedBidomainSolver<ELEMENT_DIM,SPACE_DIM>::~MatrixBasedBidomainSolver()
00223 {
00224     if(mVecForConstructingRhs)
00225     {
00226         VecDestroy(mVecForConstructingRhs);
00227         MatDestroy(mMassMatrix);
00228     }    
00229 }
00230 
00232 // explicit instantiation
00234 
00235 template class MatrixBasedBidomainSolver<1,1>;
00236 template class MatrixBasedBidomainSolver<2,2>;
00237 template class MatrixBasedBidomainSolver<3,3>;
00238 

Generated on Mon Nov 1 12:35:17 2010 for Chaste by  doxygen 1.5.5