LinearSystem.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 #include "LinearSystem.hpp"
00030 #include "PetscException.hpp"
00031 #include <iostream>
00032 #include "OutputFileHandler.hpp"
00033 #include "PetscTools.hpp"
00034 #include <cassert>
00035 #include "HeartEventHandler.hpp"
00036 #include "Timer.hpp"
00037 #include "Warnings.hpp"
00038 #include <algorithm>
00039 
00040 
00042 // Implementation
00044 
00045 LinearSystem::LinearSystem(PetscInt lhsVectorSize, unsigned rowPreallocation)
00046    :mPrecondMatrix(NULL),
00047     mSize(lhsVectorSize),
00048     mMatNullSpace(NULL),
00049     mDestroyMatAndVec(true),
00050     mKspIsSetup(false),
00051     mNonZerosUsed(0.0),
00052     mMatrixIsConstant(false),
00053     mTolerance(1e-6),
00054     mUseAbsoluteTolerance(false),
00055     mDirichletBoundaryConditionsVector(NULL),
00056     mpBlockDiagonalPC(NULL),
00057     mpLDUFactorisationPC(NULL),
00058     mpTwoLevelsBlockDiagonalPC(NULL),
00059     mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00060     mPrecondMatrixIsNotLhs(false),
00061     mRowPreallocation(rowPreallocation),
00062     mUseFixedNumberIterations(false),
00063     mEvaluateNumItsEveryNSolves(UINT_MAX),
00064     mpConvergenceTestContext(NULL),
00065     mEigMin(DBL_MAX),
00066     mEigMax(DBL_MIN),
00067     mForceSpectrumReevaluation(false)
00068 {
00069     assert(lhsVectorSize>0);
00070     if (mRowPreallocation == UINT_MAX)
00071     {
00072         //Automatic preallocation if it's a small matrix
00073         if (lhsVectorSize<15)
00074         {
00075             mRowPreallocation=lhsVectorSize;
00076         }
00077         else
00078         {
00079             EXCEPTION("You must provide a rowPreallocation argument for a large sparse system");
00080         }
00081     }
00082     
00083     mRhsVector=PetscTools::CreateVec(mSize);
00084     PetscTools::SetupMat(mLhsMatrix, mSize, mSize, mRowPreallocation, PETSC_DECIDE, PETSC_DECIDE);
00085 
00086     VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00087 
00090     mKspType = "gmres";
00091     mPcType = "jacobi";
00092 
00093     mNumSolves = 0;
00094 #ifdef TRACE_KSP
00095     mTotalNumIterations = 0;
00096     mMaxNumIterations = 0;
00097 #endif
00098 }
00099 
00100 LinearSystem::LinearSystem(PetscInt lhsVectorSize, Mat lhsMatrix, Vec rhsVector)
00101    :mPrecondMatrix(NULL),
00102     mSize(lhsVectorSize),
00103     mMatNullSpace(NULL),
00104     mDestroyMatAndVec(true),
00105     mKspIsSetup(false),
00106     mNonZerosUsed(0.0),
00107     mMatrixIsConstant(false),
00108     mTolerance(1e-6),
00109     mUseAbsoluteTolerance(false),
00110     mDirichletBoundaryConditionsVector(NULL),
00111     mpBlockDiagonalPC(NULL),
00112     mpLDUFactorisationPC(NULL),
00113     mpTwoLevelsBlockDiagonalPC(NULL),
00114     mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00115     mPrecondMatrixIsNotLhs(false),
00116     mUseFixedNumberIterations(false),
00117     mEvaluateNumItsEveryNSolves(UINT_MAX),
00118     mpConvergenceTestContext(NULL),
00119     mEigMin(DBL_MAX),
00120     mEigMax(DBL_MIN),
00121     mForceSpectrumReevaluation(false)
00122 {
00123     assert(lhsVectorSize>0);
00124     // Conveniently, PETSc Mats and Vecs are actually pointers
00125     mLhsMatrix = lhsMatrix;
00126     mRhsVector = rhsVector;
00127 
00128     VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00129 
00130     mNumSolves = 0;
00131 #ifdef TRACE_KSP
00132     mTotalNumIterations = 0;
00133     mMaxNumIterations = 0;
00134 #endif
00135 }
00136 
00137 LinearSystem::LinearSystem(Vec templateVector, unsigned rowPreallocation)
00138    :mPrecondMatrix(NULL),
00139     mMatNullSpace(NULL),
00140     mDestroyMatAndVec(true),
00141     mKspIsSetup(false),
00142     mMatrixIsConstant(false),
00143     mTolerance(1e-6),
00144     mUseAbsoluteTolerance(false),
00145     mDirichletBoundaryConditionsVector(NULL),
00146     mpBlockDiagonalPC(NULL),
00147     mpLDUFactorisationPC(NULL),
00148     mpTwoLevelsBlockDiagonalPC(NULL),
00149     mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00150     mPrecondMatrixIsNotLhs(false),
00151     mRowPreallocation(rowPreallocation),
00152     mUseFixedNumberIterations(false),
00153     mEvaluateNumItsEveryNSolves(UINT_MAX),
00154     mpConvergenceTestContext(NULL),
00155     mEigMin(DBL_MAX),
00156     mEigMax(DBL_MIN),
00157     mForceSpectrumReevaluation(false)
00158 {
00159     VecDuplicate(templateVector, &mRhsVector);
00160     VecGetSize(mRhsVector, &mSize);
00161     VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00162     PetscInt local_size = mOwnershipRangeHi - mOwnershipRangeLo;
00163 
00164     PetscTools::SetupMat(mLhsMatrix, mSize, mSize, mRowPreallocation, local_size, local_size);
00165 
00168     mKspType = "gmres";
00169     mPcType = "jacobi";
00170 
00171     mNumSolves = 0;
00172 #ifdef TRACE_KSP
00173     mTotalNumIterations = 0;
00174     mMaxNumIterations = 0;
00175 #endif
00176 }
00177 
00178 LinearSystem::LinearSystem(Vec residualVector, Mat jacobianMatrix)
00179    :mPrecondMatrix(NULL),
00180     mMatNullSpace(NULL),
00181     mDestroyMatAndVec(false),
00182     mKspIsSetup(false),
00183     mMatrixIsConstant(false),
00184     mTolerance(1e-6),
00185     mUseAbsoluteTolerance(false),
00186     mDirichletBoundaryConditionsVector(NULL),
00187     mpBlockDiagonalPC(NULL),
00188     mpLDUFactorisationPC(NULL),
00189     mpTwoLevelsBlockDiagonalPC(NULL),
00190     mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00191     mPrecondMatrixIsNotLhs(false),
00192     mRowPreallocation(UINT_MAX),
00193     mUseFixedNumberIterations(false),
00194     mEvaluateNumItsEveryNSolves(UINT_MAX),
00195     mpConvergenceTestContext(NULL),
00196     mEigMin(DBL_MAX),
00197     mEigMax(DBL_MIN),
00198     mForceSpectrumReevaluation(false)
00199 {
00200     assert(residualVector || jacobianMatrix);
00201     mRhsVector = residualVector;
00202     mLhsMatrix = jacobianMatrix;
00203 
00204     PetscInt mat_size=0, vec_size=0;
00205     if (mRhsVector)
00206     {
00207         VecGetSize(mRhsVector, &vec_size);
00208         mSize = (unsigned)vec_size;
00209         VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00210     }
00211     if (mLhsMatrix)
00212     {
00213         PetscInt mat_cols;
00214         MatGetSize(mLhsMatrix, &mat_size, &mat_cols);
00215         assert(mat_size == mat_cols);
00216         mSize = (unsigned)mat_size;
00217         MatGetOwnershipRange(mLhsMatrix, &mOwnershipRangeLo, &mOwnershipRangeHi);
00218 
00219         MatInfo matrix_info;
00220         MatGetInfo(mLhsMatrix, MAT_GLOBAL_MAX, &matrix_info);
00221 
00222         /*
00223          *  Assuming that mLhsMatrix was created with PetscTools::SetupMat, the value
00224          *  below should be equivalent to what was used as preallocation in that call.
00225          */
00226         mRowPreallocation = (unsigned) matrix_info.nz_allocated / mSize;
00227     }
00228     assert(!mRhsVector || !mLhsMatrix || vec_size == mat_size);
00229 
00232     mKspType = "gmres";
00233     mPcType = "jacobi";
00234 
00235     mNumSolves = 0;
00236 #ifdef TRACE_KSP
00237     mTotalNumIterations = 0;
00238     mMaxNumIterations = 0;
00239 #endif
00240 }
00241 
00242 LinearSystem::~LinearSystem()
00243 {
00244     delete mpBlockDiagonalPC;
00245     delete mpLDUFactorisationPC;
00246     delete mpTwoLevelsBlockDiagonalPC;
00247 
00248     if (mDestroyMatAndVec)
00249     {
00250         VecDestroy(mRhsVector);
00251         MatDestroy(mLhsMatrix);
00252     }
00253 
00254     if(mPrecondMatrixIsNotLhs)
00255     {
00256         MatDestroy(mPrecondMatrix);
00257     }
00258 
00259     if (mMatNullSpace)
00260     {
00261         MatNullSpaceDestroy(mMatNullSpace);
00262     }
00263 
00264     if (mKspIsSetup)
00265     {
00266         KSPDestroy(mKspSolver);
00267     }
00268 
00269     if (mDirichletBoundaryConditionsVector)
00270     {
00272         VecDestroy(mDirichletBoundaryConditionsVector);
00273     }
00274 
00275 #if (PETSC_VERSION_MAJOR == 3)   
00276     if (mpConvergenceTestContext)
00277     {
00278         KSPDefaultConvergedDestroy(mpConvergenceTestContext);
00279     }
00280 #endif    
00281 
00282 #ifdef TRACE_KSP
00283     if (PetscTools::AmMaster())
00284     {
00285         if (mNumSolves > 0)
00286         {
00287             double ave_num_iterations = mTotalNumIterations/(double)mNumSolves;
00288     
00289             std::cout << std::endl << "KSP iterations report:" << std::endl;
00290             std::cout << "mNumSolves" << "\t" << "mTotalNumIterations" << "\t" << "mMaxNumIterations" << "\t" << "mAveNumIterations" << std::endl;
00291             std::cout << mNumSolves << "\t" << mTotalNumIterations << "\t" << mMaxNumIterations << "\t" << ave_num_iterations << std::endl;
00292         }
00293     }
00294 #endif
00295 
00296 }
00297 
00298 
00299 void LinearSystem::SetMatrixElement(PetscInt row, PetscInt col, double value)
00300 {
00301     PetscMatTools::SetElement(mLhsMatrix, row, col, value);
00302 }
00303 
00304 void LinearSystem::AddToMatrixElement(PetscInt row, PetscInt col, double value)
00305 {
00306     PetscMatTools::AddToElement(mLhsMatrix, row, col, value);
00307 }
00308 
00309 void LinearSystem::AssembleFinalLinearSystem()
00310 {
00311     AssembleFinalLhsMatrix();
00312     AssembleRhsVector();
00313 }
00314 
00315 void LinearSystem::AssembleIntermediateLinearSystem()
00316 {
00317     AssembleIntermediateLhsMatrix();
00318     AssembleRhsVector();
00319 }
00320 
00321 void LinearSystem::AssembleFinalLhsMatrix()
00322 {
00323     PetscMatTools::AssembleFinal(mLhsMatrix);
00324 }
00325 
00326 void LinearSystem::AssembleIntermediateLhsMatrix()
00327 {
00328     PetscMatTools::AssembleIntermediate(mLhsMatrix);
00329 }
00330 
00331 void LinearSystem::AssembleFinalPrecondMatrix()
00332 {
00333     PetscMatTools::AssembleFinal(mPrecondMatrix);
00334 }
00335 
00336 void LinearSystem::AssembleRhsVector()
00337 {
00338     PetscVecTools::Assemble(mRhsVector);
00339 }
00340 
00341 void LinearSystem::SetRhsVectorElement(PetscInt row, double value)
00342 {
00343     PetscVecTools::SetElement(mRhsVector, row, value);
00344 }
00345 
00346 void LinearSystem::AddToRhsVectorElement(PetscInt row, double value)
00347 {
00348     PetscVecTools::AddToElement(mRhsVector, row, value);
00349 }
00350 
00351 void LinearSystem::DisplayMatrix()
00352 {
00353     PetscMatTools::Display(mLhsMatrix);
00354 }
00355 
00356 void LinearSystem::DisplayRhs()
00357 {
00358     PetscVecTools::Display(mRhsVector);
00359 }
00360 
00361 void LinearSystem::SetMatrixRow(PetscInt row, double value)
00362 {
00363     PetscMatTools::SetRow(mLhsMatrix, row, value);
00364 }
00365 
00366 Vec LinearSystem::GetMatrixRowDistributed(unsigned rowIndex)
00367 {
00368     return PetscMatTools::GetMatrixRowDistributed(mLhsMatrix, rowIndex);
00369 }
00370 
00371 void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal(std::vector<unsigned>& rRows, double diagonalValue)
00372 {
00373     PetscMatTools::ZeroRowsWithValueOnDiagonal(mLhsMatrix, rRows, diagonalValue);
00374 }
00375 
00376 
00377 void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector<unsigned>& rRowColIndices, double diagonalValue)
00378 {
00379     PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mLhsMatrix, rRowColIndices, diagonalValue);
00380 }
00381 
00382 void LinearSystem::ZeroMatrixColumn(PetscInt col)
00383 {
00384     PetscMatTools::ZeroColumn(mLhsMatrix, col);
00385 }
00386 
00387 void LinearSystem::ZeroRhsVector()
00388 {
00389     PetscVecTools::Zero(mRhsVector);
00390 }
00391 
00392 void LinearSystem::ZeroLhsMatrix()
00393 {
00394     PetscMatTools::Zero(mLhsMatrix);
00395 }
00396 
00397 void LinearSystem::ZeroLinearSystem()
00398 {
00399     ZeroRhsVector();
00400     ZeroLhsMatrix();
00401 }
00402 
00403 unsigned LinearSystem::GetSize() const
00404 {
00405     return (unsigned) mSize;
00406 }
00407 
00408 void LinearSystem::SetNullBasis(Vec nullBasis[], unsigned numberOfBases)
00409 {
00410 #ifndef NDEBUG
00411     // Check all the vectors of the base are normal
00412     for (unsigned vec_index=0; vec_index<numberOfBases; vec_index++)
00413     {
00414         PetscReal l2_norm;
00415         VecNorm(nullBasis[vec_index], NORM_2, &l2_norm);
00416         if (fabs(l2_norm-1.0) > 1e-08)
00417         {
00418             EXCEPTION("One of the vectors in the null space is not normalised");
00419         }
00420     }
00421 
00422     // Check all the vectors of the base are orthogonal
00423     for (unsigned vec_index=1; vec_index<numberOfBases; vec_index++)
00424     {
00425         // The strategy is to check the (i-1)-th vector against vectors from i to n with VecMDot. This should be the most efficient way of doing it.
00426         unsigned num_vectors_ahead = numberOfBases-vec_index;
00427         PetscScalar dot_products[num_vectors_ahead];
00428 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00429         VecMDot(num_vectors_ahead, nullBasis[vec_index-1], &nullBasis[vec_index], dot_products);
00430 #else
00431         VecMDot(nullBasis[vec_index-1], num_vectors_ahead, &nullBasis[vec_index], dot_products);
00432 #endif
00433         for (unsigned index=0; index<num_vectors_ahead; index++)
00434         {
00435             if (fabs(dot_products[index]) > 1e-08 )
00436             {
00437                 EXCEPTION("The null space is not orthogonal.");
00438             }
00439         }
00440 
00441     }
00442 
00443 #endif
00444 
00445     PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, numberOfBases, nullBasis, &mMatNullSpace) );
00446 }
00447 
00448 void LinearSystem::RemoveNullSpace()
00449 {
00450     // Only remove if previously set
00451     if (mMatNullSpace)
00452     {
00453         PETSCEXCEPT( MatNullSpaceDestroy(mMatNullSpace) );
00454         PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, 0, NULL, &mMatNullSpace) );
00455         if (mKspIsSetup)
00456         {
00457             PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
00458         }
00459         //else: it will be set next time Solve() is called
00460     }
00461 }
00462 
00463 
00464 void LinearSystem::GetOwnershipRange(PetscInt& lo, PetscInt& hi)
00465 {
00466     lo = mOwnershipRangeLo;
00467     hi = mOwnershipRangeHi;
00468 }
00469 
00470 double LinearSystem::GetMatrixElement(PetscInt row, PetscInt col)
00471 {
00472     return PetscMatTools::GetElement(mLhsMatrix, row, col);
00473 }
00474 
00475 double LinearSystem::GetRhsVectorElement(PetscInt row)
00476 {
00477     return PetscVecTools::GetElement(mRhsVector, row);
00478 }
00479 
00480 unsigned LinearSystem::GetNumIterations() const
00481 {
00482     assert(this->mKspIsSetup);
00483 
00484     PetscInt num_its;
00485     KSPGetIterationNumber(this->mKspSolver, &num_its);
00486 
00487     return (unsigned) num_its;
00488 }
00489 
00490 
00491 Vec& LinearSystem::rGetRhsVector()
00492 {
00493     return mRhsVector;
00494 }
00495 
00496 Vec LinearSystem::GetRhsVector() const
00497 {
00498     return mRhsVector;
00499 }
00500 
00501 Mat& LinearSystem::rGetLhsMatrix()
00502 {
00503     return mLhsMatrix;
00504 }
00505 
00506 Mat LinearSystem::GetLhsMatrix() const
00507 {
00508     return mLhsMatrix;
00509 }
00510 
00511 Mat& LinearSystem::rGetPrecondMatrix()
00512 {
00513     if(!mPrecondMatrixIsNotLhs)
00514     {
00515         EXCEPTION("LHS matrix used for preconditioner construction");
00516     }
00517     
00518     return mPrecondMatrix;
00519 }
00520 
00521 Vec& LinearSystem::rGetDirichletBoundaryConditionsVector()
00522 {
00523     return mDirichletBoundaryConditionsVector;
00524 }
00525 
00526 void LinearSystem::SetMatrixIsSymmetric(bool isSymmetric)
00527 {
00529 
00530     if (isSymmetric)
00531     {
00532         PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRIC);
00533         PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL);
00534     }
00535     else
00536     {
00537 // don't have a PetscMatTools method for setting options to false        
00538 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
00539         MatSetOption(mLhsMatrix, MAT_SYMMETRIC, PETSC_FALSE);
00540         MatSetOption(mLhsMatrix, MAT_STRUCTURALLY_SYMMETRIC, PETSC_FALSE);
00541         MatSetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL, PETSC_FALSE);
00542 #else
00543         MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRIC);
00544         MatSetOption(mLhsMatrix, MAT_NOT_STRUCTURALLY_SYMMETRIC);
00545         MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRY_ETERNAL);
00546 #endif
00547     }
00548 }
00549 
00550 bool LinearSystem::IsMatrixSymmetric()
00551 {
00552     PetscTruth symmetry_flag_is_set;
00553     PetscTruth symmetry_flag;
00554 
00555     MatIsSymmetricKnown(mLhsMatrix, &symmetry_flag_is_set, &symmetry_flag);
00556 
00557     // If the flag is not set we assume is a non-symmetric matrix
00558     return symmetry_flag_is_set && symmetry_flag;
00559 }
00560 
00561 void LinearSystem::SetMatrixIsConstant(bool matrixIsConstant)
00562 {
00563     mMatrixIsConstant=matrixIsConstant;
00564 }
00565 
00566 void LinearSystem::SetRelativeTolerance(double relativeTolerance)
00567 {
00568     mTolerance=relativeTolerance;
00569     mUseAbsoluteTolerance=false;
00570     if (mKspIsSetup)
00571     {
00572         KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00573     }
00574 }
00575 
00576 void LinearSystem::SetAbsoluteTolerance(double absoluteTolerance)
00577 {
00578     mTolerance=absoluteTolerance;
00579     mUseAbsoluteTolerance=true;
00580     if (mKspIsSetup)
00581     {
00582         KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00583     }
00584 }
00585 
00586 void LinearSystem::SetKspType(const char *kspType)
00587 {
00588     mKspType = kspType;
00589     if (mKspIsSetup)
00590     {
00591         KSPSetType(mKspSolver, kspType);
00592         KSPSetFromOptions(mKspSolver);
00593     }
00594 }
00595 
00596 void LinearSystem::SetPcType(const char* pcType, boost::shared_ptr<std::vector<PetscInt> > pBathNodes)
00597 {
00598     mPcType=pcType;
00599     mpBathNodes = pBathNodes;
00600     
00601     if (mKspIsSetup)
00602     {
00603         if (mPcType == "blockdiagonal")
00604         {
00605             // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
00607             delete mpBlockDiagonalPC;
00608             mpBlockDiagonalPC = NULL;
00609             delete mpLDUFactorisationPC;
00610             mpLDUFactorisationPC = NULL;
00611             delete mpTwoLevelsBlockDiagonalPC;
00612             mpTwoLevelsBlockDiagonalPC = NULL;
00613 
00614             mpBlockDiagonalPC = new PCBlockDiagonal(mKspSolver);
00615         }
00616         else if (mPcType == "ldufactorisation")
00617         {
00618             // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
00620             delete mpBlockDiagonalPC;
00621             mpBlockDiagonalPC = NULL;
00622             delete mpLDUFactorisationPC;
00623             mpLDUFactorisationPC = NULL;
00624             delete mpTwoLevelsBlockDiagonalPC;
00625             mpTwoLevelsBlockDiagonalPC = NULL;
00626 
00627             mpLDUFactorisationPC = new PCLDUFactorisation(mKspSolver);
00628         }
00629         else if (mPcType == "twolevelsblockdiagonal")
00630         {
00631             // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
00633             delete mpBlockDiagonalPC;
00634             mpBlockDiagonalPC = NULL;
00635             delete mpLDUFactorisationPC;
00636             mpLDUFactorisationPC = NULL;
00637             delete mpTwoLevelsBlockDiagonalPC;
00638             mpTwoLevelsBlockDiagonalPC = NULL;
00639 
00640             if (!mpBathNodes)
00641             {
00642                 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
00643             }
00644             mpTwoLevelsBlockDiagonalPC = new PCTwoLevelsBlockDiagonal(mKspSolver, *mpBathNodes);
00645         }        
00646         else
00647         {
00648             PC prec;
00649             KSPGetPC(mKspSolver, &prec);
00650             PCSetType(prec, pcType);
00651         }
00652         KSPSetFromOptions(mKspSolver);
00653     }
00654 }
00655 
00656 Vec LinearSystem::Solve(Vec lhsGuess)
00657 {
00658     /* The following lines are very useful for debugging
00659      *    MatView(mLhsMatrix,    PETSC_VIEWER_STDOUT_WORLD);
00660      *    VecView(mRhsVector,    PETSC_VIEWER_STDOUT_WORLD);
00661      */
00662     //Double check that the non-zero pattern hasn't changed
00663     MatInfo mat_info;
00664     MatGetInfo(mLhsMatrix, MAT_GLOBAL_SUM, &mat_info);
00665 
00666     if (!mKspIsSetup)
00667     {
00668         // Create PETSc Vec that may be required if we use a Chebyshev solver
00669         Vec chebyshev_lhs_vector = NULL;
00670 
00671         HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00672         mNonZerosUsed=mat_info.nz_used;
00673         //MatNorm(mLhsMatrix, NORM_FROBENIUS, &mMatrixNorm);
00674         PC prec; //Type of pre-conditioner
00675 
00676         KSPCreate(PETSC_COMM_WORLD, &mKspSolver);
00677         //See
00678         //http://www-unix.mcs.anl.gov/petsc/petsc-2/snapshots/petsc-current/docs/manualpages/KSP/KSPSetOperators.html
00679         //The preconditioner flag (last argument) in the following calls says
00680         //how to reuse the preconditioner on subsequent iterations
00681 
00682         MatStructure preconditioner_over_successive_calls;
00683 
00684         if (mMatrixIsConstant)
00685         {
00686             preconditioner_over_successive_calls = SAME_PRECONDITIONER;
00687         }
00688         else
00689         {
00690             preconditioner_over_successive_calls = SAME_NONZERO_PATTERN;
00691         }
00692 
00693         if (mPrecondMatrixIsNotLhs)
00694         {
00695             KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix, preconditioner_over_successive_calls);
00696         }
00697         else
00698         {
00699             KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix, preconditioner_over_successive_calls);
00700         }
00701 
00702         // Set either absolute or relative tolerance of the KSP solver.
00703         // The default is to use relative tolerance (1e-6)
00704         if (mUseAbsoluteTolerance)
00705         {
00706             KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00707         }
00708         else
00709         {
00710             KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00711         }
00712 
00713         // set ksp and pc types
00714         KSPSetType(mKspSolver, mKspType.c_str());
00715         KSPGetPC(mKspSolver, &prec);
00716 
00717 
00718         // Turn off pre-conditioning if the system size is very small
00719         if (mSize <= 4)
00720         {
00721             PCSetType(prec, PCNONE);
00722         }
00723         else
00724         {
00725 #ifdef TRACE_KSP
00726             Timer::Reset();
00727 #endif
00728             if (mPcType == "blockdiagonal")
00729             {
00730                 mpBlockDiagonalPC = new PCBlockDiagonal(mKspSolver);
00731 #ifdef TRACE_KSP
00732                 if (PetscTools::AmMaster())
00733                 {
00734                     Timer::Print("Purpose-build preconditioner creation");
00735                 }
00736 #endif
00737             }
00738             else if (mPcType == "ldufactorisation")
00739             {
00740                 mpLDUFactorisationPC = new PCLDUFactorisation(mKspSolver);
00741 #ifdef TRACE_KSP
00742                 if (PetscTools::AmMaster())
00743                 {
00744                     Timer::Print("Purpose-build preconditioner creation");
00745                 }
00746 #endif
00747             }
00748             else if (mPcType == "twolevelsblockdiagonal")
00749             {
00750                 if (!mpBathNodes)
00751                 {
00752                     TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
00753                 }                
00754                 mpTwoLevelsBlockDiagonalPC = new PCTwoLevelsBlockDiagonal(mKspSolver, *mpBathNodes);
00755 #ifdef TRACE_KSP
00756                 if (PetscTools::AmMaster())
00757                 {
00758                     Timer::Print("Purpose-build preconditioner creation");
00759                 }
00760 #endif
00761 
00762             }
00763             else 
00764             {
00765                 PCSetType(prec, mPcType.c_str());
00766             }
00767         }
00768 
00769         if (mMatNullSpace)
00770         {
00772             PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
00773         }
00774 
00775         if (lhsGuess)
00776         {
00777             //Assume that the user of this method will always be kind enough
00778             //to give us a reasonable guess.
00779             KSPSetInitialGuessNonzero(mKspSolver,PETSC_TRUE);
00780         }
00781 
00782         KSPSetFromOptions(mKspSolver);
00783 
00784         /*
00785          * Non-adaptive Chebyshev: the required spectrum approximation is computed just once 
00786          * at the beginning of the simulation. This is done with two extra CG solves. 
00787          */
00788         if(mKspType == "chebychev" && !mUseFixedNumberIterations)        
00789         {            
00790 #ifdef TRACE_KSP
00791             Timer::Reset();
00792 #endif
00793             // Preconditioned matrix spectrum is approximated with CG
00794             KSPSetType(mKspSolver,"cg");
00795             KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
00796             KSPSetUp(mKspSolver);
00797                             
00798             VecDuplicate(mRhsVector, &chebyshev_lhs_vector);
00799             if (lhsGuess)
00800             {
00801                 VecCopy(lhsGuess, chebyshev_lhs_vector);
00802             }
00803 
00804             // Smallest eigenvalue is approximated to default tolerance
00805             KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
00806 
00807             PetscReal *r_eig = new PetscReal[mSize];
00808             PetscReal *c_eig = new PetscReal[mSize];
00809             PetscInt eigs_computed;
00810             KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
00811 
00812             mEigMin = r_eig[0];
00813     
00814             // Largest eigenvalue is approximated to machine precision
00815             KSPSetTolerances(mKspSolver, DBL_EPSILON, DBL_EPSILON, PETSC_DEFAULT, PETSC_DEFAULT);
00816             KSPSetUp(mKspSolver);
00817             if (lhsGuess)
00818             {
00819                 VecCopy(lhsGuess, chebyshev_lhs_vector);
00820             }
00821 
00822             KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
00823             KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
00824 
00825             mEigMax = r_eig[eigs_computed-1];
00826 
00827 #ifdef TRACE_KSP
00828             /*
00829              *  Under certain circumstances (see Golub&Overton 1988), underestimating
00830              * the spectrum of the preconditioned operator improves convergence rate.
00831              * See publication for a discussion and for definition of alpha and sigma_one.
00832              */
00833             if (PetscTools::AmMaster())
00834             {
00835                 std::cout << "EIGS: ";
00836                 for (int index=0; index<eigs_computed; index++)
00837                 {
00838                     std::cout << r_eig[index] << ", ";
00839                 }
00840                 std::cout << std::endl;
00841             }
00842 
00843             if (PetscTools::AmMaster()) std::cout << "EIGS "<< mEigMax << " " << mEigMin <<std::endl;
00844             double alpha = 2/(mEigMax+mEigMin);
00845             double sigma_one = 1 - alpha*mEigMin;
00846             if (PetscTools::AmMaster()) std::cout << "sigma_1 = 1 - alpha*mEigMin = "<< sigma_one <<std::endl;
00847 #endif
00848 
00849             // Set Chebyshev solver and max/min eigenvalues
00850             assert(mKspType == "chebychev");
00851             KSPSetType(mKspSolver, mKspType.c_str());            
00852             KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
00853             KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
00854             if (mUseAbsoluteTolerance)
00855             {
00856                 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00857             }
00858             else
00859             {
00860                 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00861             }
00862 
00863             delete[] r_eig;
00864             delete[] c_eig;
00865 
00866 #ifdef TRACE_KSP
00867             if (PetscTools::AmMaster()) 
00868             {
00869                 Timer::Print("Computing extremal eigenvalues");
00870             }
00871 #endif
00872         }
00873 
00874 #ifdef TRACE_KSP
00875         Timer::Reset();
00876 #endif
00877 
00878         KSPSetUp(mKspSolver);
00879 
00880         if (chebyshev_lhs_vector)
00881         {            
00882             VecDestroy(chebyshev_lhs_vector);
00883         }
00884 
00885 #ifdef TRACE_KSP
00886         if (PetscTools::AmMaster())
00887         {
00888             Timer::Print("KSPSetUP (contains preconditioner creation for PETSc preconditioners)");
00889         }
00890 #endif
00891 
00892         mKspIsSetup = true;
00893 
00894         HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00895     }
00896     else
00897     {
00898         #define COVERAGE_IGNORE
00899         if (mNonZerosUsed!=mat_info.nz_used)
00900         {
00901             EXCEPTION("LinearSystem doesn't allow the non-zero pattern of a matrix to change. (I think you changed it).");
00902         }
00903 //        PetscScalar norm;
00904 //        MatNorm(mLhsMatrix, NORM_FROBENIUS, &norm);
00905 //        if (fabs(norm - mMatrixNorm) > 0)
00906 //        {
00907 //            EXCEPTION("LinearSystem doesn't allow the matrix norm to change");
00908 //        }
00909         #undef COVERAGE_IGNORE
00910     }
00911 
00912     HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00913     // Create solution vector
00915     Vec lhs_vector;
00916     VecDuplicate(mRhsVector, &lhs_vector);//Sets the same size (doesn't copy)
00917     if (lhsGuess)
00918     {
00919         VecCopy(lhsGuess, lhs_vector);
00920     }
00921     HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00922 //    //Double check that the mRhsVector contains sensible values
00923 //    double* p_rhs,* p_guess;
00924 //    VecGetArray(mRhsVector, &p_rhs);
00925 //    VecGetArray(lhsGuess, &p_guess);
00926 //    for (int global_index=mOwnershipRangeLo; global_index<mOwnershipRangeHi; global_index++)
00927 //    {
00928 //        int local_index = global_index - mOwnershipRangeLo;
00929 //        assert(!std::isnan(p_rhs[local_index]));
00930 //        assert(!std::isnan(p_guess[local_index]));
00931 //        if (p_rhs[local_index] != p_rhs[local_index])
00932 //        {
00933 //            std::cout << "********* PETSc nan\n";
00934 //            assert(0);
00935 //        }
00936 //    }
00937 //    std::cout << "b[0] = " << p_rhs[0] << ", guess[0] = " << p_guess[0] << "\n";
00938 //    VecRestoreArray(mRhsVector, &p_rhs);
00939 //    VecRestoreArray(lhsGuess, &p_guess);
00940 //    // Test A*guess
00941 //    Vec temp;
00942 //    VecDuplicate(mRhsVector, &temp);
00943 //    MatMult(mLhsMatrix, lhs_vector, temp);
00944 //    double* p_temp;
00945 //    VecGetArray(temp, &p_temp);
00946 //    std::cout << "temp[0] = " << p_temp[0] << "\n";
00947 //    VecRestoreArray(temp, &p_temp);
00948 //    VecDestroy(temp);
00949 //    // Dump the matrix to file
00950 //    PetscViewer viewer;
00951 //    PetscViewerASCIIOpen(PETSC_COMM_WORLD,"mat.output",&viewer);
00952 //    MatView(mLhsMatrix, viewer);
00953 //    PetscViewerFlush(viewer);
00954 //    PetscViewerDestroy(viewer);
00955 //    // Dump the rhs vector to file
00956 //    PetscViewerASCIIOpen(PETSC_COMM_WORLD,"vec.output",&viewer);
00957 //    PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
00958 //    VecView(mRhsVector, viewer);
00959 //    PetscViewerFlush(viewer);
00960 //    PetscViewerDestroy(viewer);
00961 
00962     try
00963     {
00964         HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
00965 
00966 #ifdef TRACE_KSP
00967         Timer::Reset();
00968 #endif
00969 
00970         // Current solve has to be done with tolerance-based stop criteria in order to record iterations taken
00971         if(mUseFixedNumberIterations && (mNumSolves%mEvaluateNumItsEveryNSolves==0 || mForceSpectrumReevaluation))
00972         {
00973 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
00974             KSPSetNormType(mKspSolver, KSP_PRECONDITIONED_NORM);
00975 #else
00976             KSPSetNormType(mKspSolver, KSP_NORM_PRECONDITIONED);
00977 #endif
00978 
00979 #if (PETSC_VERSION_MAJOR == 3)
00980             if (!mpConvergenceTestContext)
00981             {
00982                 KSPDefaultConvergedCreate(&mpConvergenceTestContext);
00983             }            
00984             KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, &mpConvergenceTestContext, PETSC_NULL); 
00985 #else
00986             KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, PETSC_NULL);
00987 #endif            
00988 
00989             if (mUseAbsoluteTolerance)
00990             {
00991                 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00992             }
00993             else
00994             {
00995                 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00996             }
00997 
00999             std::stringstream num_it_str;
01000             num_it_str << 1000;
01001             PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01002 
01003             // Adaptive Chebyshev: reevaluate spectrum with cg
01004             if (mKspType == "chebychev")
01005             {
01006                 // You can estimate preconditioned matrix spectrum with CG
01007                 KSPSetType(mKspSolver,"cg");
01008                 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
01009             }
01010             
01011             KSPSetFromOptions(mKspSolver);            
01012             KSPSetUp(mKspSolver);
01013         }
01014 
01015         PETSCEXCEPT(KSPSolve(mKspSolver, mRhsVector, lhs_vector));
01016         HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
01017 
01018 #ifdef TRACE_KSP
01019         PetscInt num_it;
01020         KSPGetIterationNumber(mKspSolver, &num_it);
01021         if (PetscTools::AmMaster())
01022         {
01023             std::cout << "++ Solve: " << mNumSolves << " NumIterations: " << num_it << " "; // don't add std::endl so we get Timer::Print output in the same line (better for grep-ing)
01024             Timer::Print("Solve");
01025         }
01026         
01027         mTotalNumIterations += num_it;
01028         if ((unsigned) num_it > mMaxNumIterations)
01029         {
01030             mMaxNumIterations = num_it;
01031         }
01032 #endif
01033 
01034         // Check that solver converged and throw if not
01035         KSPConvergedReason reason;
01036         KSPGetConvergedReason(mKspSolver, &reason);
01037 
01038         if (mUseFixedNumberIterations && PETSC_VERSION_MAJOR < 3)
01039         {
01040             WARNING("Not explicitly checking convergence reason when using fixed number of iterations and PETSc 2");
01041         }
01042         else
01043         {
01044             KSPEXCEPT(reason);
01045         }
01046 
01047         if(mUseFixedNumberIterations && (mNumSolves%mEvaluateNumItsEveryNSolves==0 || mForceSpectrumReevaluation))
01048         {
01049             // Adaptive Chebyshev: reevaluate spectrum with cg            
01050             if (mKspType == "chebychev")
01051             {
01052                 PetscReal *r_eig = new PetscReal[mSize];
01053                 PetscReal *c_eig = new PetscReal[mSize];
01054                 PetscInt eigs_computed;
01055                 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
01056 
01057                 mEigMin = r_eig[0];
01058                 /*
01059                  * Using max() covers a borderline case found in TestChasteBenchmarksForPreDiCT where there's a big
01060                  * gap in the spectrum between ~1.2 and ~2.5. Some reevaluations pick up 2.5 and others don't. If it
01061                  * is not picked up, Chebyshev will diverge after 10 solves or so. 
01062                  */
01063                 mEigMax = std::max(mEigMax,r_eig[eigs_computed-1]);
01064 
01065                 delete[] r_eig;
01066                 delete[] c_eig;
01067 
01068                 assert(mKspType == "chebychev");
01069                 KSPSetType(mKspSolver, mKspType.c_str());            
01070                 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
01071                 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
01072             }
01073 
01074 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))            
01075             if (mKspType == "chebychev")
01076             {
01077                 // See #1695 for more details.
01078                 EXCEPTION("Chebyshev with fixed number of iterations is known to be broken in PETSc <= 2.3.2");
01079             }
01080 
01081             KSPSetNormType(mKspSolver, KSP_NO_NORM);
01082 #else
01083             KSPSetNormType(mKspSolver, KSP_NORM_NO);
01084 #endif                        
01085 
01086 #if (PETSC_VERSION_MAJOR != 3)
01087             KSPSetConvergenceTest(mKspSolver, KSPSkipConverged, PETSC_NULL);
01088 #endif
01089 
01090             PetscInt num_it;
01091             KSPGetIterationNumber(mKspSolver, &num_it);            
01092             std::stringstream num_it_str;
01093             num_it_str << num_it;
01094             PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01095 
01096             KSPSetFromOptions(mKspSolver);
01097             KSPSetUp(mKspSolver);
01098             
01099             mForceSpectrumReevaluation=false;
01100         }
01101 
01102         mNumSolves++;
01103 
01104     }
01105     catch (const Exception& e)
01106     {
01107         // Destroy solution vector on error to avoid memory leaks
01108         VecDestroy(lhs_vector);
01109         throw e;
01110     }
01111 
01112     return lhs_vector;
01113 }
01114 
01115 
01116 void LinearSystem::SetPrecondMatrixIsDifferentFromLhs(bool precondIsDifferent)
01117 {
01118     mPrecondMatrixIsNotLhs = precondIsDifferent;
01119     
01120     if (mPrecondMatrixIsNotLhs)
01121     {
01122         if (mRowPreallocation == UINT_MAX)
01123         {
01124             /*
01125              *  At the time of writing, this line will be reached if the constructor
01126              *  with signature LinearSystem(Vec residualVector, Mat jacobianMatrix) is
01127              *  called with jacobianMatrix=NULL and preconditioning matrix different
01128              *  from lhs is used.
01129              *
01130              *  If this combination is ever required you will need to work out
01131              *  matrix allocation (mRowPreallocation) here.
01132              */
01133             NEVER_REACHED;
01134         }
01135 
01136         PetscInt local_size = mOwnershipRangeHi - mOwnershipRangeLo;                
01137         PetscTools::SetupMat(mPrecondMatrix, mSize, mSize, mRowPreallocation, local_size, local_size);        
01138     }
01139 }
01140 
01141 void LinearSystem::SetUseFixedNumberIterations(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
01142 {
01143 
01144     mUseFixedNumberIterations = useFixedNumberIterations;
01145     mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
01146 }
01147 
01148 void LinearSystem::ResetKspSolver()
01149 {
01150     if (mKspIsSetup)
01151     {
01152         KSPDestroy(mKspSolver);
01153     }
01154 
01155     mKspIsSetup = false;
01156     mForceSpectrumReevaluation = true;
01157     
01158     /*
01159      *  Reset max number of iterations. This option is stored in the configuration database and 
01160      * explicitely read in with KSPSetFromOptions() everytime a KSP object is created. Therefore,
01161      * destroying the KSP object will not ensure that it is set back to default.
01162      */
01164     std::stringstream num_it_str;
01165     num_it_str << 1000;
01166     PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());    
01167 }
01168 
01169 // Serialization for Boost >= 1.36
01170 #include "SerializationExportWrapperForCpp.hpp"
01171 CHASTE_CLASS_EXPORT(LinearSystem)

Generated on Tue May 31 14:31:42 2011 for Chaste by  doxygen 1.5.5