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