Chaste Release::3.1
BidomainSolver.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 
00037 #include "BidomainSolver.hpp"
00038 #include "BidomainAssembler.hpp"
00039 #include "BidomainWithBathAssembler.hpp"
00040 #include "PetscMatTools.hpp"
00041 
00042 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00043 void BidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00044 {
00045     if (this->mpLinearSystem != NULL)
00046     {
00047         return;
00048     }
00049     AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(initialSolution);
00050 
00051     // initialise matrix-based RHS vector and matrix, and use the linear
00052     // system rhs as a template
00053     Vec& r_template = this->mpLinearSystem->rGetRhsVector();
00054     VecDuplicate(r_template, &mVecForConstructingRhs);
00055     PetscInt ownership_range_lo;
00056     PetscInt ownership_range_hi;
00057     VecGetOwnershipRange(r_template, &ownership_range_lo, &ownership_range_hi);
00058     PetscInt local_size = ownership_range_hi - ownership_range_lo;
00059     PetscTools::SetupMat(mMassMatrix, 2*this->mpMesh->GetNumNodes(), 2*this->mpMesh->GetNumNodes(),
00060                          2*this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00061                          local_size, local_size);
00062 }
00063 
00064 
00065 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00066 void BidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetupLinearSystem(
00067         Vec currentSolution,
00068         bool computeMatrix)
00069 {
00070     assert(this->mpLinearSystem->rGetLhsMatrix() != NULL);
00071     assert(this->mpLinearSystem->rGetRhsVector() != NULL);
00072     assert(currentSolution != NULL);
00073 
00074 
00076     // set up LHS matrix (and mass matrix)
00078     if (computeMatrix)
00079     {
00080         mpBidomainAssembler->SetMatrixToAssemble(this->mpLinearSystem->rGetLhsMatrix());
00081         mpBidomainAssembler->AssembleMatrix();
00082 
00083         // the BidomainMassMatrixAssembler deals with the mass matrix
00084         // for both bath and nonbath problems
00085         assert(SPACE_DIM==ELEMENT_DIM);
00086         BidomainMassMatrixAssembler<SPACE_DIM> mass_matrix_assembler(this->mpMesh);
00087         mass_matrix_assembler.SetMatrixToAssemble(mMassMatrix);
00088         mass_matrix_assembler.Assemble();
00089 
00090         this->mpLinearSystem->SwitchWriteModeLhsMatrix();
00091         PetscMatTools::Finalise(mMassMatrix);
00092     }
00093 
00094 
00095     HeartEventHandler::BeginEvent(HeartEventHandler::ASSEMBLE_RHS);
00096 
00098     // Set up z in b=Mz
00100     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00101 
00102     // dist stripe for the current Voltage
00103     DistributedVector distributed_current_solution = p_factory->CreateDistributedVector(currentSolution);
00104     DistributedVector::Stripe distributed_current_solution_vm(distributed_current_solution, 0);
00105 
00106     // dist stripe for z
00107     DistributedVector dist_vec_matrix_based = p_factory->CreateDistributedVector(mVecForConstructingRhs);
00108     DistributedVector::Stripe dist_vec_matrix_based_vm(dist_vec_matrix_based, 0);
00109     DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
00110 
00111     double Am = HeartConfig::Instance()->GetSurfaceAreaToVolumeRatio();
00112     double Cm  = HeartConfig::Instance()->GetCapacitance();
00113 
00114     if(!(this->mBathSimulation))
00115     {
00116         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00117              index!= dist_vec_matrix_based.End();
00118              ++index)
00119         {
00120             double V = distributed_current_solution_vm[index];
00121             double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00122                        - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00123 
00124             dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00125             dist_vec_matrix_based_phie[index] = 0.0;
00126         }
00127     }
00128     else
00129     {
00130         DistributedVector::Stripe dist_vec_matrix_based_phie(dist_vec_matrix_based, 1);
00131 
00132         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00133              index!= dist_vec_matrix_based.End();
00134              ++index)
00135         {
00136 
00137             if ( !HeartRegionCode::IsRegionBath( this->mpMesh->GetNode(index.Global)->GetRegion() ))
00138             {
00139                 double V = distributed_current_solution_vm[index];
00140                 double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00141                            - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00142 
00143                 dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00144             }
00145             else
00146             {
00147                 dist_vec_matrix_based_vm[index] = 0.0;
00148             }
00149 
00150             dist_vec_matrix_based_phie[index] = 0.0;
00151 
00152         }
00153     }
00154 
00155     dist_vec_matrix_based.Restore();
00156 
00158     // b = Mz
00160     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00161 
00162     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00163     // the event will be begun again inside mpBidomainAssembler->AssembleVector();
00164     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00165 
00166 
00168     // apply Neumann boundary conditions
00170     mpBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions); // as the BCC can change
00171     mpBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00172     mpBidomainNeumannSurfaceTermAssembler->AssembleVector();
00173 
00174 
00176     // apply correction term
00178     if(mpBidomainCorrectionTermAssembler)
00179     {
00180         mpBidomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00181         // don't need to set current solution
00182         mpBidomainCorrectionTermAssembler->AssembleVector();
00183     }
00184 
00185     this->mpLinearSystem->FinaliseRhsVector();
00186 
00187     this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
00188 
00189     if(this->mBathSimulation)
00190     {
00191         this->mpLinearSystem->FinaliseLhsMatrix();
00192         this->FinaliseForBath(computeMatrix,true);
00193     }
00194 
00195     if(computeMatrix)
00196     {
00197         this->mpLinearSystem->FinaliseLhsMatrix();
00198     }
00199     this->mpLinearSystem->FinaliseRhsVector();
00200 }
00201 
00202 
00203 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00204 BidomainSolver<ELEMENT_DIM,SPACE_DIM>::BidomainSolver(
00205         bool bathSimulation,
00206         AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00207         BidomainTissue<SPACE_DIM>* pTissue,
00208         BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,2>* pBoundaryConditions,
00209         unsigned numQuadPoints)
00210     : AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
00211 {
00212     // Tell tissue there's no need to replicate ionic caches
00213     pTissue->SetCacheReplication(false);
00214     mVecForConstructingRhs = NULL;
00215 
00216     // create assembler
00217     if(bathSimulation)
00218     {
00219         mpBidomainAssembler = new BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mNumQuadPoints);
00220     }
00221     else
00222     {
00223         mpBidomainAssembler = new BidomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mNumQuadPoints);
00224     }
00225 
00226 
00227     mpBidomainNeumannSurfaceTermAssembler = new BidomainNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM>(pMesh,pBoundaryConditions);
00228 
00229     if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00230     {
00231         mpBidomainCorrectionTermAssembler
00232             = new BidomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue,this->mNumQuadPoints);
00233         //We are going to need those caches after all
00234         pTissue->SetCacheReplication(true);
00235     }
00236     else
00237     {
00238         mpBidomainCorrectionTermAssembler = NULL;
00239     }
00240 }
00241 
00242 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00243 BidomainSolver<ELEMENT_DIM,SPACE_DIM>::~BidomainSolver()
00244 {
00245     delete mpBidomainAssembler;
00246     delete mpBidomainNeumannSurfaceTermAssembler;
00247 
00248     if(mVecForConstructingRhs)
00249     {
00250         PetscTools::Destroy(mVecForConstructingRhs);
00251         PetscTools::Destroy(mMassMatrix);
00252     }
00253 
00254     if(mpBidomainCorrectionTermAssembler)
00255     {
00256         delete mpBidomainCorrectionTermAssembler;
00257     }
00258 }
00259 
00261 // explicit instantiation
00263 
00264 template class BidomainSolver<1,1>;
00265 template class BidomainSolver<2,2>;
00266 template class BidomainSolver<3,3>;
00267