Chaste Release::3.1
AbstractBidomainSolver.cpp
00001 
00002 /*
00003 
00004 Copyright (c) 2005-2012, University of Oxford.
00005 All rights reserved.
00006 
00007 University of Oxford means the Chancellor, Masters and Scholars of the
00008 University of Oxford, having an administrative office at Wellington
00009 Square, Oxford OX1 2JD, UK.
00010 
00011 This file is part of Chaste.
00012 
00013 Redistribution and use in source and binary forms, with or without
00014 modification, are permitted provided that the following conditions are met:
00015  * Redistributions of source code must retain the above copyright notice,
00016    this list of conditions and the following disclaimer.
00017  * Redistributions in binary form must reproduce the above copyright notice,
00018    this list of conditions and the following disclaimer in the documentation
00019    and/or other materials provided with the distribution.
00020  * Neither the name of the University of Oxford nor the names of its
00021    contributors may be used to endorse or promote products derived from this
00022    software without specific prior written permission.
00023 
00024 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00025 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00026 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00027 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00028 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00029 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00030 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00031 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00032 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00033 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00034 
00035 */
00036 
00037 
00038 #include "AbstractBidomainSolver.hpp"
00039 #include "TetrahedralMesh.hpp"
00040 #include "PetscMatTools.hpp"
00041 #include "PetscVecTools.hpp"
00042 
00043 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00044 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::InitialiseForSolve(Vec initialSolution)
00045 {
00046     if (this->mpLinearSystem != NULL)
00047     {
00048         return;
00049     }
00050 
00051     // linear system created here
00052     AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,2>::InitialiseForSolve(initialSolution);
00053 
00054     if (HeartConfig::Instance()->GetUseAbsoluteTolerance())
00055     {
00056 #ifdef TRACE_KSP
00057         if (PetscTools::AmMaster())
00058         {
00059             std::cout << "Using absolute tolerance: " << mpConfig->GetAbsoluteTolerance() <<"\n";
00060         }
00061 #endif
00062         this->mpLinearSystem->SetAbsoluteTolerance(mpConfig->GetAbsoluteTolerance());
00063     }
00064     else
00065     {
00066 #ifdef TRACE_KSP
00067         if (PetscTools::AmMaster())
00068         {
00069             std::cout << "Using relative tolerance: " << mpConfig->GetRelativeTolerance() <<"\n";
00070         }
00071 #endif
00072         this->mpLinearSystem->SetRelativeTolerance(mpConfig->GetRelativeTolerance());
00073     }
00074 
00075     this->mpLinearSystem->SetKspType(HeartConfig::Instance()->GetKSPSolver());
00076 
00078     if(std::string("twolevelsblockdiagonal") == std::string(HeartConfig::Instance()->GetKSPPreconditioner()))
00079     {
00081         TetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* p_mesh = dynamic_cast<TetrahedralMesh<ELEMENT_DIM,SPACE_DIM>*>(this->mpMesh);
00082         if (p_mesh && PetscTools::IsSequential())
00083         {
00084             boost::shared_ptr<std::vector<PetscInt> > p_bath_nodes(new std::vector<PetscInt>());
00085 
00086             for(unsigned node_index=0; node_index<this->mpMesh->GetNumNodes(); node_index++)
00087             {
00088                 if (HeartRegionCode::IsRegionBath( this->mpMesh->GetNode(node_index)->GetRegion() ))
00089                 {
00090                     p_bath_nodes->push_back(node_index);
00091                 }
00092             }
00093 
00094             this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner(), p_bath_nodes);
00095         }
00096         else
00097         {
00098             TERMINATE("Two levels block diagonal only works with TetrahedralMesh and p=1");
00099         }
00100     }
00101     else
00102     {
00103         this->mpLinearSystem->SetPcType(HeartConfig::Instance()->GetKSPPreconditioner());
00104     }
00105 
00106     if (mRowForAverageOfPhiZeroed==INT_MAX)
00107     {
00108         // not applying average(phi)=0 constraint, so matrix is symmetric
00109         this->mpLinearSystem->SetMatrixIsSymmetric(true);
00110     }
00111     else
00112     {
00113         // applying average(phi)=0 constraint, so matrix is not symmetric
00114         this->mpLinearSystem->SetMatrixIsSymmetric(false);
00115     }
00116 
00117     this->mpLinearSystem->SetUseFixedNumberIterations(HeartConfig::Instance()->GetUseFixedNumberIterationsLinearSolver(), HeartConfig::Instance()->GetEvaluateNumItsEveryNSolves());
00118 }
00119 
00120 
00121 
00122 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00123 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::PrepareForSetupLinearSystem(Vec existingSolution)
00124 {
00125     double time = PdeSimulationTime::GetTime();
00126     double dt = PdeSimulationTime::GetPdeTimeStep();
00127     mpBidomainTissue->SolveCellSystems(existingSolution, time, time+dt);
00128 }
00129 
00130 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00131 Vec AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::GenerateNullBasis() const
00132 {
00133     double sqrt_num_nodes = sqrt((double) this->mpMesh->GetNumNodes());
00134 
00135     Vec null_basis;
00136     DistributedVectorFactory* p_factory = this->mpMesh->GetDistributedVectorFactory();
00137     null_basis = p_factory->CreateVec(2);
00138 
00139     DistributedVector dist_null_basis = p_factory->CreateDistributedVector(null_basis);
00140     DistributedVector::Stripe null_basis_stripe_0(dist_null_basis,0);
00141     DistributedVector::Stripe null_basis_stripe_1(dist_null_basis,1);
00142     for (DistributedVector::Iterator index = dist_null_basis.Begin();
00143          index != dist_null_basis.End();
00144          ++index)
00145     {
00146         null_basis_stripe_0[index] = 0.0;
00147         null_basis_stripe_1[index] = 1.0/sqrt_num_nodes; // normalised vector
00148     }
00149     dist_null_basis.Restore();
00150 
00151     return null_basis;
00152 }
00153 
00154 
00155 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00156 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::FinaliseLinearSystem(Vec existingSolution)
00157 {
00158     // If no dirichlet boundary conditions
00159     //  (i) Check compatibility condition to check we are solving
00160     //      a linear system that can be solved
00161     // Then either
00162     //  (a) If not setting average(phi)=0, we are solving a singular system,
00163     //      so set up a null space
00164     //  (b) Apply average(phi)=0 constraint by altering the last row, to
00165     //      get a non-singular system
00166     //
00167     if (!(GetBoundaryConditions()->HasDirichletBoundaryConditions()))
00168     {
00169         // first check compatibility condition
00170         CheckCompatibilityCondition();
00171 
00172         // Check whether applying average(phi_e)=0 constraint
00173         if (mRowForAverageOfPhiZeroed==INT_MAX)
00174         {
00175             // We're not using the constraint, hence use a null space
00176             if (!mNullSpaceCreated)
00177             {
00178                 // No null space set up, so create one and pass it to the linear system
00179                 Vec null_basis[] = {GenerateNullBasis()};
00180 
00181                 this->mpLinearSystem->SetNullBasis(null_basis, 1);
00182 
00183                 PetscTools::Destroy(null_basis[0]);
00184                 mNullSpaceCreated = true;
00185             }
00186         }
00187         else
00188         {
00189             // Using the average(phi_e)=0 constraint
00190 
00191             // CG (default solver) won't work since the system isn't symmetric anymore. Switch to GMRES
00192             this->mpLinearSystem->SetKspType("gmres"); // Switches the solver
00193             mpConfig->SetKSPSolver("gmres"); // Makes sure this change will be reflected in the XML file written to disk at the end of the simulation.
00194 
00195             // Set average phi_e to zero
00196             unsigned matrix_size = this->mpLinearSystem->GetSize();
00197             if (!this->mMatrixIsAssembled)
00198             {
00199 
00200                 // Set the mRowForAverageOfPhiZeroed-th matrix row to 0 1 0 1 ...
00201                 std::vector<unsigned> row_for_average;
00202                 row_for_average.push_back(mRowForAverageOfPhiZeroed);
00203                 this->mpLinearSystem->ZeroMatrixRowsWithValueOnDiagonal(row_for_average, 0.0);
00204                 for (unsigned col_index=0; col_index<matrix_size; col_index++)
00205                 {
00206                     if (col_index%2 == 1)
00207                     {
00208                         this->mpLinearSystem->SetMatrixElement(mRowForAverageOfPhiZeroed, col_index, 1);
00209                     }
00210 
00211                 }
00212                 this->mpLinearSystem->FinaliseLhsMatrix();
00213 
00214             }
00215             // Set the mRowForAverageOfPhiZeroed-th rhs vector row to 0
00216             this->mpLinearSystem->SetRhsVectorElement(mRowForAverageOfPhiZeroed, 0);
00217 
00218             this->mpLinearSystem->FinaliseRhsVector();
00219         }
00220     }
00221 }
00222 
00223 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00224 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::CheckCompatibilityCondition()
00225 {
00226     if (GetBoundaryConditions()->HasDirichletBoundaryConditions() || this->mRowForAverageOfPhiZeroed!=INT_MAX)
00227     {
00228         // not a singular system, no compatibility condition
00229         return;
00230     }
00231 
00232 #ifndef NDEBUG
00233     DistributedVector distributed_rhs=this->mpMesh->GetDistributedVectorFactory()->CreateDistributedVector(this->mpLinearSystem->rGetRhsVector());
00234     DistributedVector::Stripe distributed_rhs_phi_entries(distributed_rhs,1);
00235 
00236     double local_sum=0;
00237     for (DistributedVector::Iterator index = distributed_rhs.Begin();
00238          index!= distributed_rhs.End();
00239          ++index)
00240     {
00241         local_sum += distributed_rhs_phi_entries[index];
00242     }
00243 
00244     double global_sum;
00245     int mpi_ret = MPI_Allreduce(&local_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00246     assert(mpi_ret == MPI_SUCCESS);
00247 
00248     if(fabs(global_sum)>1e-6) // magic number! sum should really be a sum of zeros and exactly zero though anyway (or a-a+b-b+c-c.. etc in the case of electrodes)
00249     {
00250         #define COVERAGE_IGNORE
00251         // shouldn't ever reach this line but useful to have the error printed out if you do
00252         //std::cout << "Sum of b_{2i+1} = " << global_sum << " (should be zero for compatibility)\n";
00253         EXCEPTION("Linear system does not satisfy compatibility constraint!");
00254         #undef COVERAGE_IGNORE
00255     }
00256 #endif
00257 }
00258 
00259 
00260 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00261 AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::AbstractBidomainSolver(
00262             bool bathSimulation,
00263             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh,
00264             BidomainTissue<SPACE_DIM>* pTissue,
00265             BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,2>* pBoundaryConditions,
00266             unsigned numQuadPoints)
00267     : AbstractDynamicLinearPdeSolver<ELEMENT_DIM,SPACE_DIM,2>(pMesh),
00268       mBathSimulation(bathSimulation),
00269       mpBidomainTissue(pTissue),
00270       mpBoundaryConditions(pBoundaryConditions),
00271       mNumQuadPoints(numQuadPoints)
00272 {
00273     assert(pTissue != NULL);
00274     assert(pBoundaryConditions != NULL);
00275 
00276     mNullSpaceCreated = false;
00277 
00278     // important!
00279     this->mMatrixIsConstant = true;
00280 
00281     mRowForAverageOfPhiZeroed = INT_MAX; //this->mpLinearSystem->GetSize() - 1;
00282     mpConfig = HeartConfig::Instance();
00283 }
00284 
00285 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00286 AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::~AbstractBidomainSolver()
00287 {
00288 }
00289 
00290 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00291 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetFixedExtracellularPotentialNodes(
00292             std::vector<unsigned> fixedExtracellularPotentialNodes)
00293 {
00294     for (unsigned i=0; i<fixedExtracellularPotentialNodes.size(); i++)
00295     {
00296         if (fixedExtracellularPotentialNodes[i] >= this->mpMesh->GetNumNodes() )
00297         {
00298             EXCEPTION("Fixed node number must be less than total number nodes");
00299         }
00300     }
00301 
00302     mFixedExtracellularPotentialNodes = fixedExtracellularPotentialNodes;
00303 
00304     // We will need to recalculate this when HasDirichletBoundaryConditions() is called.
00305     GetBoundaryConditions()->ResetDirichletCommunication();
00306 
00307     for (unsigned i=0; i<mFixedExtracellularPotentialNodes.size(); i++)
00308     {
00309         if (this->mpMesh->GetDistributedVectorFactory()->IsGlobalIndexLocal(mFixedExtracellularPotentialNodes[i]))
00310         {
00311             ConstBoundaryCondition<SPACE_DIM>* p_boundary_condition
00312                  = new ConstBoundaryCondition<SPACE_DIM>(0.0);
00313 
00314             //Throws if node is not owned locally
00315             Node<SPACE_DIM>* p_node = this->mpMesh->GetNode(mFixedExtracellularPotentialNodes[i]);
00316 
00317             GetBoundaryConditions()->AddDirichletBoundaryCondition(p_node, p_boundary_condition, 1);
00318 
00319         }
00320     }
00321 }
00322 
00323 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00324 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::SetRowForAverageOfPhiZeroed(unsigned row)
00325 {
00326     // Row should be odd in C++-like indexing
00327     if (row%2 == 0)
00328     {
00329         EXCEPTION("Row for applying the constraint 'Average of phi_e = zero' should be odd in C++ like indexing");
00330     }
00331 
00332     mRowForAverageOfPhiZeroed = row;
00333 }
00334 
00335 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00336 void AbstractBidomainSolver<ELEMENT_DIM,SPACE_DIM>::FinaliseForBath(bool computeMatrix, bool computeVector)
00337 {
00338     assert(mBathSimulation);
00339     PetscTruth is_matrix_assembled;
00340     MatAssembled(this->mpLinearSystem->GetLhsMatrix(), &is_matrix_assembled);
00341     assert(is_matrix_assembled == PETSC_TRUE);
00342 
00343     /*
00344      *  Before revision 6516, we used to zero out i-th row and column here. It seems to be redundant because they are already zero after assembly.
00345      *  When assembling a bath element you get a matrix subblock that looks like (2D example):
00346      *
00347      *     Vm   0 0 0 0 0 0
00348      *     Vm   0 0 0 0 0 0
00349      *     Vm   0 0 0 0 0 0
00350      *     Phie 0 0 0 x x x
00351      *     Phie 0 0 0 x x x  -> the x subblock is assembled from div(grad_phi) = 0
00352      *     Phie 0 0 0 x x x
00353      *
00354      *  Therefore, all the Vm entries of this node are already 0.
00355      *
00356      *  Explicitly checking it in non-production builds.
00357      */
00358 #ifndef NDEBUG
00359     if(computeMatrix)
00360     {
00361         for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator iter=this->mpMesh->GetNodeIteratorBegin();
00362              iter != this->mpMesh->GetNodeIteratorEnd();
00363              ++iter)
00364         {
00365             if (HeartRegionCode::IsRegionBath( (*iter).GetRegion() ))
00366             {
00367                 int num_equation = 2*iter->GetIndex(); // assumes Vm and Phie are interleaved
00368 
00369                 PetscInt local_lo, local_hi;
00370                 this->mpLinearSystem->GetOwnershipRange(local_lo, local_hi);
00371 
00372                 // If this processor owns i-th row, check it.
00373                 if ((local_lo <= (int)num_equation) && ((int)num_equation < local_hi))
00374                 {
00375                     for (unsigned column=0; column < this->mpLinearSystem->GetSize(); column++)
00376                     {
00377                         assert(this->mpLinearSystem->GetMatrixElement(num_equation, column)==0.0);
00378                     }
00379                 }
00380 
00381                 // Check the local entries of the i-th column
00382                 for (int row=local_lo; row<local_hi; row++)
00383                 {
00384                     assert(this->mpLinearSystem->GetMatrixElement(row, num_equation)==0);
00385                 }
00386             }
00387         }
00388     }
00389 #endif
00390 
00391     /*
00392      *  These two loops are decoupled because interleaving calls to GetMatrixElement and MatSetValue
00393      *  require reassembling the matrix before GetMatrixElement which generates massive communication
00394      *  overhead for large models and/or large core counts.
00395      */
00396 
00397     for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::NodeIterator iter=this->mpMesh->GetNodeIteratorBegin();
00398          iter != this->mpMesh->GetNodeIteratorEnd();
00399          ++iter)
00400     {
00401         if (HeartRegionCode::IsRegionBath((*iter).GetRegion() ))
00402         {
00403             PetscInt index = 2*iter->GetIndex(); // assumes Vm and Phie are interleaved
00404 
00405             if(computeMatrix)
00406             {
00407                 // put 1.0 on the diagonal
00408                 PetscMatTools::SetElement(this->mpLinearSystem->rGetLhsMatrix(), index, index, 1.0);
00409             }
00410 
00411             if(computeVector)
00412             {
00413                 // zero rhs vector entry
00414                 PetscVecTools::SetElement(this->mpLinearSystem->rGetRhsVector(), index, 0.0);
00415             }
00416         }
00417     }
00418 }
00419 
00421 // explicit instantiation
00423 
00424 template class AbstractBidomainSolver<1,1>;
00425 template class AbstractBidomainSolver<2,2>;
00426 template class AbstractBidomainSolver<3,3>;