BidomainSolver.cpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, 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         for (DistributedVector::Iterator index = dist_vec_matrix_based.Begin();
00131              index!= dist_vec_matrix_based.End();
00132              ++index)
00133         {
00134 
00135             if ( !HeartRegionCode::IsRegionBath( this->mpMesh->GetNode(index.Global)->GetRegion() ))
00136             {
00137                 double V = distributed_current_solution_vm[index];
00138                 double F = - Am*this->mpBidomainTissue->rGetIionicCacheReplicated()[index.Global]
00139                            - this->mpBidomainTissue->rGetIntracellularStimulusCacheReplicated()[index.Global];
00140 
00141                 dist_vec_matrix_based_vm[index] = Am*Cm*V*PdeSimulationTime::GetPdeTimeStepInverse() + F;
00142             }
00143             else
00144             {
00145                 dist_vec_matrix_based_vm[index] = 0.0;
00146             }
00147 
00148             dist_vec_matrix_based_phie[index] = 0.0;
00149 
00150         }
00151     }
00152 
00153     dist_vec_matrix_based.Restore();
00154 
00156     // b = Mz
00158     MatMult(mMassMatrix, mVecForConstructingRhs, this->mpLinearSystem->rGetRhsVector());
00159 
00160     // assembling RHS is not finished yet, as Neumann bcs are added below, but
00161     // the event will be begun again inside mpBidomainAssembler->AssembleVector();
00162     HeartEventHandler::EndEvent(HeartEventHandler::ASSEMBLE_RHS);
00163 
00164 
00166     // apply Neumann boundary conditions
00168     mpBidomainNeumannSurfaceTermAssembler->ResetBoundaryConditionsContainer(this->mpBoundaryConditions); // as the BCC can change
00169     mpBidomainNeumannSurfaceTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00170     mpBidomainNeumannSurfaceTermAssembler->AssembleVector();
00171 
00172 
00174     // apply correction term
00176     if(mpBidomainCorrectionTermAssembler)
00177     {
00178         mpBidomainCorrectionTermAssembler->SetVectorToAssemble(this->mpLinearSystem->rGetRhsVector(), false/*don't zero vector!*/);
00179         // don't need to set current solution
00180         mpBidomainCorrectionTermAssembler->AssembleVector();
00181     }
00182 
00183     this->mpLinearSystem->FinaliseRhsVector();
00184 
00185     this->mpBoundaryConditions->ApplyDirichletToLinearProblem(*(this->mpLinearSystem), computeMatrix);
00186 
00187     if(this->mBathSimulation)
00188     {
00189         this->mpLinearSystem->FinaliseLhsMatrix();
00190         this->FinaliseForBath(computeMatrix,true);
00191     }
00192 
00193     if(computeMatrix)
00194     {
00195         this->mpLinearSystem->FinaliseLhsMatrix();
00196     }
00197     this->mpLinearSystem->FinaliseRhsVector();
00198 }
00199 
00200 
00201 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00202 BidomainSolver<ELEMENT_DIM,SPACE_DIM>::BidomainSolver(
00203         bool bathSimulation,
00204         AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00205         BidomainTissue<SPACE_DIM>* pTissue,
00206         BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,2>* pBoundaryConditions)
00207     : AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>(bathSimulation,pMesh,pTissue,pBoundaryConditions)
00208 {
00209     // Tell tissue there's no need to replicate ionic caches
00210     pTissue->SetCacheReplication(false);
00211     mVecForConstructingRhs = NULL;
00212 
00213     // create assembler
00214     if(bathSimulation)
00215     {
00216         mpBidomainAssembler = new BidomainWithBathAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue);
00217     }
00218     else
00219     {
00220         mpBidomainAssembler = new BidomainAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue);
00221     }
00222 
00223 
00224     mpBidomainNeumannSurfaceTermAssembler = new BidomainNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM>(pMesh,pBoundaryConditions);
00225 
00226     if(HeartConfig::Instance()->GetUseStateVariableInterpolation())
00227     {
00228         mpBidomainCorrectionTermAssembler
00229             = new BidomainCorrectionTermAssembler<ELEMENT_DIM,SPACE_DIM>(this->mpMesh,this->mpBidomainTissue);
00230         //We are going to need those caches after all
00231         pTissue->SetCacheReplication(true);
00232     }
00233     else
00234     {
00235         mpBidomainCorrectionTermAssembler = NULL;
00236     }
00237 }
00238 
00239 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00240 BidomainSolver<ELEMENT_DIM,SPACE_DIM>::~BidomainSolver()
00241 {
00242     delete mpBidomainAssembler;
00243     delete mpBidomainNeumannSurfaceTermAssembler;
00244 
00245     if(mVecForConstructingRhs)
00246     {
00247         PetscTools::Destroy(mVecForConstructingRhs);
00248         PetscTools::Destroy(mMassMatrix);
00249     }
00250 
00251     if(mpBidomainCorrectionTermAssembler)
00252     {
00253         delete mpBidomainCorrectionTermAssembler;
00254     }
00255 }
00256 
00258 // explicit instantiation
00260 
00261 template class BidomainSolver<1,1>;
00262 template class BidomainSolver<2,2>;
00263 template class BidomainSolver<3,3>;
00264 

Generated by  doxygen 1.6.2