00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #include "LinearSystem.hpp"
00037
00038 #include <iostream>
00039 #include <cassert>
00040 #include <algorithm>
00041 #include <boost/scoped_array.hpp>
00042
00043 #include "PetscException.hpp"
00044 #include "OutputFileHandler.hpp"
00045 #include "PetscTools.hpp"
00046 #include "HeartEventHandler.hpp"
00047 #include "Timer.hpp"
00048 #include "Warnings.hpp"
00049
00051
00053
00054 LinearSystem::LinearSystem(PetscInt lhsVectorSize, unsigned rowPreallocation)
00055 :mPrecondMatrix(NULL),
00056 mSize(lhsVectorSize),
00057 mMatNullSpace(NULL),
00058 mDestroyMatAndVec(true),
00059 mKspIsSetup(false),
00060 mNonZerosUsed(0.0),
00061 mMatrixIsConstant(false),
00062 mTolerance(1e-6),
00063 mUseAbsoluteTolerance(false),
00064 mDirichletBoundaryConditionsVector(NULL),
00065 mpBlockDiagonalPC(NULL),
00066 mpLDUFactorisationPC(NULL),
00067 mpTwoLevelsBlockDiagonalPC(NULL),
00068 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00069 mPrecondMatrixIsNotLhs(false),
00070 mRowPreallocation(rowPreallocation),
00071 mUseFixedNumberIterations(false),
00072 mEvaluateNumItsEveryNSolves(UINT_MAX),
00073 mpConvergenceTestContext(NULL),
00074 mEigMin(DBL_MAX),
00075 mEigMax(DBL_MIN),
00076 mForceSpectrumReevaluation(false)
00077 {
00078 assert(lhsVectorSize > 0);
00079 if (mRowPreallocation == UINT_MAX)
00080 {
00081
00082 if (lhsVectorSize<15)
00083 {
00084 mRowPreallocation=lhsVectorSize;
00085 }
00086 else
00087 {
00088 EXCEPTION("You must provide a rowPreallocation argument for a large sparse system");
00089 }
00090 }
00091
00092 mRhsVector = PetscTools::CreateVec(mSize);
00093 PetscTools::SetupMat(mLhsMatrix, mSize, mSize, mRowPreallocation, PETSC_DECIDE, PETSC_DECIDE);
00094
00095 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00096
00099 mKspType = "gmres";
00100 mPcType = "jacobi";
00101
00102 mNumSolves = 0;
00103 #ifdef TRACE_KSP
00104 mTotalNumIterations = 0;
00105 mMaxNumIterations = 0;
00106 #endif
00107 }
00108
00109 LinearSystem::LinearSystem(PetscInt lhsVectorSize, Mat lhsMatrix, Vec rhsVector)
00110 :mPrecondMatrix(NULL),
00111 mSize(lhsVectorSize),
00112 mMatNullSpace(NULL),
00113 mDestroyMatAndVec(true),
00114 mKspIsSetup(false),
00115 mNonZerosUsed(0.0),
00116 mMatrixIsConstant(false),
00117 mTolerance(1e-6),
00118 mUseAbsoluteTolerance(false),
00119 mDirichletBoundaryConditionsVector(NULL),
00120 mpBlockDiagonalPC(NULL),
00121 mpLDUFactorisationPC(NULL),
00122 mpTwoLevelsBlockDiagonalPC(NULL),
00123 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00124 mPrecondMatrixIsNotLhs(false),
00125 mUseFixedNumberIterations(false),
00126 mEvaluateNumItsEveryNSolves(UINT_MAX),
00127 mpConvergenceTestContext(NULL),
00128 mEigMin(DBL_MAX),
00129 mEigMax(DBL_MIN),
00130 mForceSpectrumReevaluation(false)
00131 {
00132 assert(lhsVectorSize > 0);
00133
00134 mLhsMatrix = lhsMatrix;
00135 mRhsVector = rhsVector;
00136
00137 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00138
00139 mNumSolves = 0;
00140 #ifdef TRACE_KSP
00141 mTotalNumIterations = 0;
00142 mMaxNumIterations = 0;
00143 #endif
00144 }
00145
00146 LinearSystem::LinearSystem(Vec templateVector, unsigned rowPreallocation, bool newAllocationError)
00147 :mPrecondMatrix(NULL),
00148 mMatNullSpace(NULL),
00149 mDestroyMatAndVec(true),
00150 mKspIsSetup(false),
00151 mMatrixIsConstant(false),
00152 mTolerance(1e-6),
00153 mUseAbsoluteTolerance(false),
00154 mDirichletBoundaryConditionsVector(NULL),
00155 mpBlockDiagonalPC(NULL),
00156 mpLDUFactorisationPC(NULL),
00157 mpTwoLevelsBlockDiagonalPC(NULL),
00158 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00159 mPrecondMatrixIsNotLhs(false),
00160 mRowPreallocation(rowPreallocation),
00161 mUseFixedNumberIterations(false),
00162 mEvaluateNumItsEveryNSolves(UINT_MAX),
00163 mpConvergenceTestContext(NULL),
00164 mEigMin(DBL_MAX),
00165 mEigMax(DBL_MIN),
00166 mForceSpectrumReevaluation(false)
00167 {
00168 VecDuplicate(templateVector, &mRhsVector);
00169 VecGetSize(mRhsVector, &mSize);
00170 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00171 PetscInt local_size = mOwnershipRangeHi - mOwnershipRangeLo;
00172
00173 PetscTools::SetupMat(mLhsMatrix, mSize, mSize, mRowPreallocation, local_size, local_size, true, newAllocationError);
00174
00177 mKspType = "gmres";
00178 mPcType = "jacobi";
00179
00180 mNumSolves = 0;
00181 #ifdef TRACE_KSP
00182 mTotalNumIterations = 0;
00183 mMaxNumIterations = 0;
00184 #endif
00185 }
00186
00187 LinearSystem::LinearSystem(Vec residualVector, Mat jacobianMatrix)
00188 :mPrecondMatrix(NULL),
00189 mMatNullSpace(NULL),
00190 mDestroyMatAndVec(false),
00191 mKspIsSetup(false),
00192 mMatrixIsConstant(false),
00193 mTolerance(1e-6),
00194 mUseAbsoluteTolerance(false),
00195 mDirichletBoundaryConditionsVector(NULL),
00196 mpBlockDiagonalPC(NULL),
00197 mpLDUFactorisationPC(NULL),
00198 mpTwoLevelsBlockDiagonalPC(NULL),
00199 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
00200 mPrecondMatrixIsNotLhs(false),
00201 mRowPreallocation(UINT_MAX),
00202 mUseFixedNumberIterations(false),
00203 mEvaluateNumItsEveryNSolves(UINT_MAX),
00204 mpConvergenceTestContext(NULL),
00205 mEigMin(DBL_MAX),
00206 mEigMax(DBL_MIN),
00207 mForceSpectrumReevaluation(false)
00208 {
00209 assert(residualVector || jacobianMatrix);
00210 mRhsVector = residualVector;
00211 mLhsMatrix = jacobianMatrix;
00212
00213 PetscInt mat_size=0, vec_size=0;
00214 if (mRhsVector)
00215 {
00216 VecGetSize(mRhsVector, &vec_size);
00217 mSize = (unsigned)vec_size;
00218 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
00219 }
00220 if (mLhsMatrix)
00221 {
00222 PetscInt mat_cols;
00223 MatGetSize(mLhsMatrix, &mat_size, &mat_cols);
00224 assert(mat_size == mat_cols);
00225 mSize = (unsigned)mat_size;
00226 MatGetOwnershipRange(mLhsMatrix, &mOwnershipRangeLo, &mOwnershipRangeHi);
00227
00228 MatInfo matrix_info;
00229 MatGetInfo(mLhsMatrix, MAT_GLOBAL_MAX, &matrix_info);
00230
00231
00232
00233
00234
00235 mRowPreallocation = (unsigned) matrix_info.nz_allocated / mSize;
00236 }
00237 assert(!mRhsVector || !mLhsMatrix || vec_size == mat_size);
00238
00241 mKspType = "gmres";
00242 mPcType = "jacobi";
00243
00244 mNumSolves = 0;
00245 #ifdef TRACE_KSP
00246 mTotalNumIterations = 0;
00247 mMaxNumIterations = 0;
00248 #endif
00249 }
00250
00251 LinearSystem::~LinearSystem()
00252 {
00253 delete mpBlockDiagonalPC;
00254 delete mpLDUFactorisationPC;
00255 delete mpTwoLevelsBlockDiagonalPC;
00256
00257 if (mDestroyMatAndVec)
00258 {
00259 PetscTools::Destroy(mRhsVector);
00260 PetscTools::Destroy(mLhsMatrix);
00261 }
00262
00263 if (mPrecondMatrixIsNotLhs)
00264 {
00265 PetscTools::Destroy(mPrecondMatrix);
00266 }
00267
00268 if (mMatNullSpace)
00269 {
00270 MatNullSpaceDestroy(PETSC_DESTROY_PARAM(mMatNullSpace));
00271 }
00272
00273 if (mKspIsSetup)
00274 {
00275 KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
00276 }
00277
00278 if (mDirichletBoundaryConditionsVector)
00279 {
00281 PetscTools::Destroy(mDirichletBoundaryConditionsVector);
00282 }
00283
00284 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
00285 if (mpConvergenceTestContext)
00286 {
00287 KSPDefaultConvergedDestroy(mpConvergenceTestContext);
00288 }
00289 #endif
00290
00291 #ifdef TRACE_KSP
00292 if (PetscTools::AmMaster())
00293 {
00294 if (mNumSolves > 0)
00295 {
00296 double ave_num_iterations = mTotalNumIterations/(double)mNumSolves;
00297
00298 std::cout << std::endl << "KSP iterations report:" << std::endl;
00299 std::cout << "mNumSolves" << "\t" << "mTotalNumIterations" << "\t" << "mMaxNumIterations" << "\t" << "mAveNumIterations" << std::endl;
00300 std::cout << mNumSolves << "\t" << mTotalNumIterations << "\t" << mMaxNumIterations << "\t" << ave_num_iterations << std::endl;
00301 }
00302 }
00303 #endif
00304 }
00305
00306 void LinearSystem::SetMatrixElement(PetscInt row, PetscInt col, double value)
00307 {
00308 PetscMatTools::SetElement(mLhsMatrix, row, col, value);
00309 }
00310
00311 void LinearSystem::AddToMatrixElement(PetscInt row, PetscInt col, double value)
00312 {
00313 PetscMatTools::AddToElement(mLhsMatrix, row, col, value);
00314 }
00315
00316 void LinearSystem::AssembleFinalLinearSystem()
00317 {
00318 FinaliseLhsMatrix();
00319 FinaliseRhsVector();
00320 }
00321
00322 void LinearSystem::AssembleIntermediateLinearSystem()
00323 {
00324 SwitchWriteModeLhsMatrix();
00325 FinaliseRhsVector();
00326 }
00327
00328 void LinearSystem::FinaliseLhsMatrix()
00329 {
00330 PetscMatTools::Finalise(mLhsMatrix);
00331 }
00332
00333 void LinearSystem::SwitchWriteModeLhsMatrix()
00334 {
00335 PetscMatTools::SwitchWriteMode(mLhsMatrix);
00336 }
00337
00338 void LinearSystem::FinalisePrecondMatrix()
00339 {
00340 PetscMatTools::Finalise(mPrecondMatrix);
00341 }
00342
00343 void LinearSystem::FinaliseRhsVector()
00344 {
00345 PetscVecTools::Finalise(mRhsVector);
00346 }
00347
00348 void LinearSystem::SetRhsVectorElement(PetscInt row, double value)
00349 {
00350 PetscVecTools::SetElement(mRhsVector, row, value);
00351 }
00352
00353 void LinearSystem::AddToRhsVectorElement(PetscInt row, double value)
00354 {
00355 PetscVecTools::AddToElement(mRhsVector, row, value);
00356 }
00357
00358 void LinearSystem::DisplayMatrix()
00359 {
00360 PetscMatTools::Display(mLhsMatrix);
00361 }
00362
00363 void LinearSystem::DisplayRhs()
00364 {
00365 PetscVecTools::Display(mRhsVector);
00366 }
00367
00368 void LinearSystem::SetMatrixRow(PetscInt row, double value)
00369 {
00370 PetscMatTools::SetRow(mLhsMatrix, row, value);
00371 }
00372
00373 Vec LinearSystem::GetMatrixRowDistributed(unsigned rowIndex)
00374 {
00375 return PetscMatTools::GetMatrixRowDistributed(mLhsMatrix, rowIndex);
00376 }
00377
00378 void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal(std::vector<unsigned>& rRows, double diagonalValue)
00379 {
00380 PetscMatTools::ZeroRowsWithValueOnDiagonal(mLhsMatrix, rRows, diagonalValue);
00381 }
00382
00383 void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector<unsigned>& rRowColIndices, double diagonalValue)
00384 {
00385 PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mLhsMatrix, rRowColIndices, diagonalValue);
00386 }
00387
00388 void LinearSystem::ZeroMatrixColumn(PetscInt col)
00389 {
00390 PetscMatTools::ZeroColumn(mLhsMatrix, col);
00391 }
00392
00393 void LinearSystem::ZeroRhsVector()
00394 {
00395 PetscVecTools::Zero(mRhsVector);
00396 }
00397
00398 void LinearSystem::ZeroLhsMatrix()
00399 {
00400 PetscMatTools::Zero(mLhsMatrix);
00401 }
00402
00403 void LinearSystem::ZeroLinearSystem()
00404 {
00405 ZeroRhsVector();
00406 ZeroLhsMatrix();
00407 }
00408
00409 unsigned LinearSystem::GetSize() const
00410 {
00411 return (unsigned) mSize;
00412 }
00413
00414 void LinearSystem::SetNullBasis(Vec nullBasis[], unsigned numberOfBases)
00415 {
00416 #ifndef NDEBUG
00417
00418 for (unsigned vec_index=0; vec_index<numberOfBases; vec_index++)
00419 {
00420 PetscReal l2_norm;
00421 VecNorm(nullBasis[vec_index], NORM_2, &l2_norm);
00422 if (fabs(l2_norm-1.0) > 1e-08)
00423 {
00424 EXCEPTION("One of the vectors in the null space is not normalised");
00425 }
00426 }
00427
00428
00429 for (unsigned vec_index=1; vec_index<numberOfBases; vec_index++)
00430 {
00431
00432 unsigned num_vectors_ahead = numberOfBases-vec_index;
00433 boost::scoped_array<PetscScalar> dot_products(new PetscScalar[num_vectors_ahead]);
00434 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00435 VecMDot(num_vectors_ahead, nullBasis[vec_index-1], &nullBasis[vec_index], dot_products.get());
00436 #else
00437 VecMDot(nullBasis[vec_index-1], num_vectors_ahead, &nullBasis[vec_index], dot_products.get());
00438 #endif
00439 for (unsigned index=0; index<num_vectors_ahead; index++)
00440 {
00441 if (fabs(dot_products[index]) > 1e-08 )
00442 {
00443 EXCEPTION("The null space is not orthogonal.");
00444 }
00445 }
00446 }
00447
00448 #endif
00449
00450 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, numberOfBases, nullBasis, &mMatNullSpace) );
00451 }
00452
00453 void LinearSystem::RemoveNullSpace()
00454 {
00455
00456 if (mMatNullSpace)
00457 {
00458 PETSCEXCEPT( MatNullSpaceDestroy(PETSC_DESTROY_PARAM(mMatNullSpace)) );
00459 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, 0, NULL, &mMatNullSpace) );
00460 if (mKspIsSetup)
00461 {
00462 PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
00463 }
00464
00465 }
00466 }
00467
00468 void LinearSystem::GetOwnershipRange(PetscInt& lo, PetscInt& hi)
00469 {
00470 lo = mOwnershipRangeLo;
00471 hi = mOwnershipRangeHi;
00472 }
00473
00474 double LinearSystem::GetMatrixElement(PetscInt row, PetscInt col)
00475 {
00476 return PetscMatTools::GetElement(mLhsMatrix, row, col);
00477 }
00478
00479 double LinearSystem::GetRhsVectorElement(PetscInt row)
00480 {
00481 return PetscVecTools::GetElement(mRhsVector, row);
00482 }
00483
00484 unsigned LinearSystem::GetNumIterations() const
00485 {
00486 assert(this->mKspIsSetup);
00487
00488 PetscInt num_its;
00489 KSPGetIterationNumber(this->mKspSolver, &num_its);
00490
00491 return (unsigned) num_its;
00492 }
00493
00494 Vec& LinearSystem::rGetRhsVector()
00495 {
00496 return mRhsVector;
00497 }
00498
00499 Vec LinearSystem::GetRhsVector() const
00500 {
00501 return mRhsVector;
00502 }
00503
00504 Mat& LinearSystem::rGetLhsMatrix()
00505 {
00506 return mLhsMatrix;
00507 }
00508
00509 Mat LinearSystem::GetLhsMatrix() const
00510 {
00511 return mLhsMatrix;
00512 }
00513
00514 Mat& LinearSystem::rGetPrecondMatrix()
00515 {
00516 if (!mPrecondMatrixIsNotLhs)
00517 {
00518 EXCEPTION("LHS matrix used for preconditioner construction");
00519 }
00520
00521 return mPrecondMatrix;
00522 }
00523
00524 Vec& LinearSystem::rGetDirichletBoundaryConditionsVector()
00525 {
00526 return mDirichletBoundaryConditionsVector;
00527 }
00528
00529 void LinearSystem::SetMatrixIsSymmetric(bool isSymmetric)
00530 {
00532
00533 if (isSymmetric)
00534 {
00535 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRIC);
00536 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL);
00537 }
00538 else
00539 {
00540
00541 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
00542 MatSetOption(mLhsMatrix, MAT_SYMMETRIC, PETSC_FALSE);
00543 MatSetOption(mLhsMatrix, MAT_STRUCTURALLY_SYMMETRIC, PETSC_FALSE);
00544 MatSetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL, PETSC_FALSE);
00545 #else
00546 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRIC);
00547 MatSetOption(mLhsMatrix, MAT_NOT_STRUCTURALLY_SYMMETRIC);
00548 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRY_ETERNAL);
00549 #endif
00550 }
00551 }
00552
00553 bool LinearSystem::IsMatrixSymmetric()
00554 {
00555 PetscTruth symmetry_flag_is_set;
00556 PetscTruth symmetry_flag;
00557
00558 MatIsSymmetricKnown(mLhsMatrix, &symmetry_flag_is_set, &symmetry_flag);
00559
00560
00561 return symmetry_flag_is_set && symmetry_flag;
00562 }
00563
00564 void LinearSystem::SetMatrixIsConstant(bool matrixIsConstant)
00565 {
00566 mMatrixIsConstant = matrixIsConstant;
00567 }
00568
00569 void LinearSystem::SetRelativeTolerance(double relativeTolerance)
00570 {
00571 mTolerance = relativeTolerance;
00572 mUseAbsoluteTolerance = false;
00573 if (mKspIsSetup)
00574 {
00575 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00576 }
00577 }
00578
00579 void LinearSystem::SetAbsoluteTolerance(double absoluteTolerance)
00580 {
00581 mTolerance = absoluteTolerance;
00582 mUseAbsoluteTolerance = true;
00583 if (mKspIsSetup)
00584 {
00585 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00586 }
00587 }
00588
00589 void LinearSystem::SetKspType(const char *kspType)
00590 {
00591 mKspType = kspType;
00592 if (mKspIsSetup)
00593 {
00594 KSPSetType(mKspSolver, kspType);
00595 KSPSetFromOptions(mKspSolver);
00596 }
00597 }
00598
00599 void LinearSystem::SetPcType(const char* pcType, boost::shared_ptr<std::vector<PetscInt> > pBathNodes)
00600 {
00601 mPcType = pcType;
00602 mpBathNodes = pBathNodes;
00603
00604 if (mKspIsSetup)
00605 {
00606 if (mPcType == "blockdiagonal")
00607 {
00608
00610 delete mpBlockDiagonalPC;
00611 mpBlockDiagonalPC = NULL;
00612 delete mpLDUFactorisationPC;
00613 mpLDUFactorisationPC = NULL;
00614 delete mpTwoLevelsBlockDiagonalPC;
00615 mpTwoLevelsBlockDiagonalPC = NULL;
00616
00617 mpBlockDiagonalPC = new PCBlockDiagonal(mKspSolver);
00618 }
00619 else if (mPcType == "ldufactorisation")
00620 {
00621
00623 delete mpBlockDiagonalPC;
00624 mpBlockDiagonalPC = NULL;
00625 delete mpLDUFactorisationPC;
00626 mpLDUFactorisationPC = NULL;
00627 delete mpTwoLevelsBlockDiagonalPC;
00628 mpTwoLevelsBlockDiagonalPC = NULL;
00629
00630 mpLDUFactorisationPC = new PCLDUFactorisation(mKspSolver);
00631 }
00632 else if (mPcType == "twolevelsblockdiagonal")
00633 {
00634
00636 delete mpBlockDiagonalPC;
00637 mpBlockDiagonalPC = NULL;
00638 delete mpLDUFactorisationPC;
00639 mpLDUFactorisationPC = NULL;
00640 delete mpTwoLevelsBlockDiagonalPC;
00641 mpTwoLevelsBlockDiagonalPC = NULL;
00642
00643 if (!mpBathNodes)
00644 {
00645 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
00646 }
00647 mpTwoLevelsBlockDiagonalPC = new PCTwoLevelsBlockDiagonal(mKspSolver, *mpBathNodes);
00648 }
00649 else
00650 {
00651 PC prec;
00652 KSPGetPC(mKspSolver, &prec);
00653 PCSetType(prec, pcType);
00654 }
00655 KSPSetFromOptions(mKspSolver);
00656 }
00657 }
00658
00659 Vec LinearSystem::Solve(Vec lhsGuess)
00660 {
00661
00662
00663
00664
00665
00666
00667
00668 MatInfo mat_info;
00669 MatGetInfo(mLhsMatrix, MAT_GLOBAL_SUM, &mat_info);
00670
00671 if (!mKspIsSetup)
00672 {
00673
00674 Vec chebyshev_lhs_vector = NULL;
00675
00676 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00677 mNonZerosUsed=mat_info.nz_used;
00678
00679 PC prec;
00680
00681 KSPCreate(PETSC_COMM_WORLD, &mKspSolver);
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692 MatStructure preconditioner_over_successive_calls;
00693
00694 if (mMatrixIsConstant)
00695 {
00696 preconditioner_over_successive_calls = SAME_PRECONDITIONER;
00697 }
00698 else
00699 {
00700 preconditioner_over_successive_calls = SAME_NONZERO_PATTERN;
00701 }
00702
00703 if (mPrecondMatrixIsNotLhs)
00704 {
00705 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix, preconditioner_over_successive_calls);
00706 }
00707 else
00708 {
00709 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix, preconditioner_over_successive_calls);
00710 }
00711
00712
00713
00714 if (mUseAbsoluteTolerance)
00715 {
00716 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00717 }
00718 else
00719 {
00720 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00721 }
00722
00723
00724 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
00725 if (mKspType == "chebychev")
00726 {
00727 KSPSetType(mKspSolver, "chebyshev");
00728 }
00729 else
00730 {
00731 KSPSetType(mKspSolver, mKspType.c_str());
00732 }
00733 #else
00734 KSPSetType(mKspSolver, mKspType.c_str());
00735 #endif
00736 KSPGetPC(mKspSolver, &prec);
00737
00738
00739 if (mSize <= 6)
00740 {
00741 PCSetType(prec, PCNONE);
00742 }
00743 else
00744 {
00745 #ifdef TRACE_KSP
00746 Timer::Reset();
00747 #endif
00748 if (mPcType == "blockdiagonal")
00749 {
00750 mpBlockDiagonalPC = new PCBlockDiagonal(mKspSolver);
00751 #ifdef TRACE_KSP
00752 if (PetscTools::AmMaster())
00753 {
00754 Timer::Print("Purpose-build preconditioner creation");
00755 }
00756 #endif
00757 }
00758 else if (mPcType == "ldufactorisation")
00759 {
00760 mpLDUFactorisationPC = new PCLDUFactorisation(mKspSolver);
00761 #ifdef TRACE_KSP
00762 if (PetscTools::AmMaster())
00763 {
00764 Timer::Print("Purpose-build preconditioner creation");
00765 }
00766 #endif
00767 }
00768 else if (mPcType == "twolevelsblockdiagonal")
00769 {
00770 if (!mpBathNodes)
00771 {
00772 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
00773 }
00774 mpTwoLevelsBlockDiagonalPC = new PCTwoLevelsBlockDiagonal(mKspSolver, *mpBathNodes);
00775 #ifdef TRACE_KSP
00776 if (PetscTools::AmMaster())
00777 {
00778 Timer::Print("Purpose-build preconditioner creation");
00779 }
00780 #endif
00781
00782 }
00783 else
00784 {
00785 PCSetType(prec, mPcType.c_str());
00786 }
00787 }
00788
00789 if (mMatNullSpace)
00790 {
00792 PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
00793 }
00794
00795 KSPSetFromOptions(mKspSolver);
00796 if (lhsGuess)
00797 {
00798
00799 KSPSetInitialGuessNonzero(mKspSolver,PETSC_TRUE);
00800 }
00801
00802
00803
00804
00805 if (mKspType == "chebychev" && !mUseFixedNumberIterations)
00806 {
00807 #ifdef TRACE_KSP
00808 Timer::Reset();
00809 #endif
00810
00811 KSPSetType(mKspSolver,"cg");
00812 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
00813 KSPSetUp(mKspSolver);
00814
00815 VecDuplicate(mRhsVector, &chebyshev_lhs_vector);
00816 if (lhsGuess)
00817 {
00818 VecCopy(lhsGuess, chebyshev_lhs_vector);
00819 }
00820
00821
00822 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
00823
00824 PetscReal *r_eig = new PetscReal[mSize];
00825 PetscReal *c_eig = new PetscReal[mSize];
00826 PetscInt eigs_computed;
00827 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
00828
00829 mEigMin = r_eig[0];
00830
00831
00832 KSPSetTolerances(mKspSolver, DBL_EPSILON, DBL_EPSILON, PETSC_DEFAULT, PETSC_DEFAULT);
00833 KSPSetUp(mKspSolver);
00834 if (lhsGuess)
00835 {
00836 VecCopy(lhsGuess, chebyshev_lhs_vector);
00837 }
00838
00839 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
00840 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
00841
00842 mEigMax = r_eig[eigs_computed-1];
00843
00844 #ifdef TRACE_KSP
00845
00846
00847
00848
00849
00850 if (PetscTools::AmMaster())
00851 {
00852 std::cout << "EIGS: ";
00853 for (int index=0; index<eigs_computed; index++)
00854 {
00855 std::cout << r_eig[index] << ", ";
00856 }
00857 std::cout << std::endl;
00858 }
00859
00860 if (PetscTools::AmMaster()) std::cout << "EIGS "<< mEigMax << " " << mEigMin <<std::endl;
00861 double alpha = 2/(mEigMax+mEigMin);
00862 double sigma_one = 1 - alpha*mEigMin;
00863 if (PetscTools::AmMaster()) std::cout << "sigma_1 = 1 - alpha*mEigMin = "<< sigma_one <<std::endl;
00864 #endif
00865
00866
00867 assert(mKspType == "chebychev");
00868 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
00869 KSPSetType(mKspSolver, "chebyshev");
00870 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
00871 #else
00872 KSPSetType(mKspSolver, mKspType.c_str());
00873 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
00874 #endif
00875 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
00876 if (mUseAbsoluteTolerance)
00877 {
00878 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00879 }
00880 else
00881 {
00882 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00883 }
00884
00885 delete[] r_eig;
00886 delete[] c_eig;
00887
00888 #ifdef TRACE_KSP
00889 if (PetscTools::AmMaster())
00890 {
00891 Timer::Print("Computing extremal eigenvalues");
00892 }
00893 #endif
00894 }
00895
00896 #ifdef TRACE_KSP
00897 Timer::Reset();
00898 #endif
00899
00900 KSPSetUp(mKspSolver);
00901
00902 if (chebyshev_lhs_vector)
00903 {
00904 PetscTools::Destroy(chebyshev_lhs_vector);
00905 }
00906
00907 #ifdef TRACE_KSP
00908 if (PetscTools::AmMaster())
00909 {
00910 Timer::Print("KSPSetUP (contains preconditioner creation for PETSc preconditioners)");
00911 }
00912 #endif
00913
00914 mKspIsSetup = true;
00915
00916 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00917 }
00918 else
00919 {
00920 if (mNonZerosUsed!=mat_info.nz_used)
00921 {
00922 WARNING("LinearSystem doesn't like the non-zero pattern of a matrix to change. (I think you changed it).");
00923 mNonZerosUsed = mat_info.nz_used;
00924 }
00925
00926
00927
00928
00929
00930
00931 }
00932
00933 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00934
00936 Vec lhs_vector;
00937 VecDuplicate(mRhsVector, &lhs_vector);
00938 if (lhsGuess)
00939 {
00940 VecCopy(lhsGuess, lhs_vector);
00941
00942
00943
00944 }
00945
00946
00947
00948 PetscReal rhs_norm;
00949 VecNorm(mRhsVector, NORM_1, &rhs_norm);
00950 double rhs_zero_tol = 1e-15;
00951
00952 if (rhs_norm < rhs_zero_tol && rhs_norm > 0.0)
00953 {
00954 WARNING("Using zero initial guess due to small right hand side vector");
00955 PetscVecTools::Zero(lhs_vector);
00956 }
00957
00958 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00959
00960
00961
00962
00963
00964
00965
00966
00967
00968
00969
00970
00971
00972
00973
00974
00975
00976
00977
00978
00979
00980
00981
00982
00983
00984
00985
00986
00987
00988
00989
00990
00991
00992
00993
00994
00995
00996
00997
00998
00999 try
01000 {
01001 HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
01002
01003 #ifdef TRACE_KSP
01004 Timer::Reset();
01005 #endif
01006
01007
01008 if (mUseFixedNumberIterations && (mNumSolves%mEvaluateNumItsEveryNSolves==0 || mForceSpectrumReevaluation))
01009 {
01010 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
01011 KSPSetNormType(mKspSolver, KSP_PRECONDITIONED_NORM);
01012 #else
01013 KSPSetNormType(mKspSolver, KSP_NORM_PRECONDITIONED);
01014 #endif
01015
01016 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
01017 if (!mpConvergenceTestContext)
01018 {
01019 KSPDefaultConvergedCreate(&mpConvergenceTestContext);
01020 }
01021 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, &mpConvergenceTestContext, PETSC_NULL);
01022 #else
01023 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, PETSC_NULL);
01024 #endif
01025
01026 if (mUseAbsoluteTolerance)
01027 {
01028 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
01029 }
01030 else
01031 {
01032 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
01033 }
01034
01036 std::stringstream num_it_str;
01037 num_it_str << 1000;
01038 PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01039
01040
01041 if (mKspType == "chebychev")
01042 {
01043
01044 KSPSetType(mKspSolver,"cg");
01045 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
01046 }
01047
01048 KSPSetFromOptions(mKspSolver);
01049 KSPSetUp(mKspSolver);
01050 }
01051
01052 PETSCEXCEPT(KSPSolve(mKspSolver, mRhsVector, lhs_vector));
01053 HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
01054
01055 #ifdef TRACE_KSP
01056 PetscInt num_it;
01057 KSPGetIterationNumber(mKspSolver, &num_it);
01058 if (PetscTools::AmMaster())
01059 {
01060 std::cout << "++ Solve: " << mNumSolves << " NumIterations: " << num_it << " ";
01061 Timer::Print("Solve");
01062 }
01063
01064 mTotalNumIterations += num_it;
01065 if ((unsigned) num_it > mMaxNumIterations)
01066 {
01067 mMaxNumIterations = num_it;
01068 }
01069 #endif
01070
01071
01072 KSPConvergedReason reason;
01073 KSPGetConvergedReason(mKspSolver, &reason);
01074
01075 if (mUseFixedNumberIterations && PETSC_VERSION_MAJOR < 3)
01076 {
01077 WARNING("Not explicitly checking convergence reason when using fixed number of iterations and PETSc 2");
01078 }
01079 else
01080 {
01081 KSPEXCEPT(reason);
01082 }
01083
01084 if (mUseFixedNumberIterations && (mNumSolves%mEvaluateNumItsEveryNSolves==0 || mForceSpectrumReevaluation))
01085 {
01086
01087 if (mKspType == "chebychev")
01088 {
01089 PetscReal *r_eig = new PetscReal[mSize];
01090 PetscReal *c_eig = new PetscReal[mSize];
01091 PetscInt eigs_computed;
01092 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
01093
01094 mEigMin = r_eig[0];
01095
01096
01097
01098
01099
01100 mEigMax = std::max(mEigMax,r_eig[eigs_computed-1]);
01101
01102 delete[] r_eig;
01103 delete[] c_eig;
01104
01105 assert(mKspType == "chebychev");
01106 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
01107 KSPSetType(mKspSolver, "chebyshev");
01108 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
01109 #else
01110 KSPSetType(mKspSolver, mKspType.c_str());
01111 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
01112 #endif
01113 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
01114 }
01115
01116 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
01117 if (mKspType == "chebychev")
01118 {
01119
01120 EXCEPTION("Chebyshev with fixed number of iterations is known to be broken in PETSc <= 2.3.2");
01121 }
01122
01123 KSPSetNormType(mKspSolver, KSP_NO_NORM);
01124 #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
01125 KSPSetNormType(mKspSolver, KSP_NORM_NONE);
01126 #else
01127 KSPSetNormType(mKspSolver, KSP_NORM_NO);
01128 #endif
01129
01130 #if (PETSC_VERSION_MAJOR == 2)
01131 KSPSetConvergenceTest(mKspSolver, KSPSkipConverged, PETSC_NULL);
01132 #endif
01133
01134 PetscInt num_it;
01135 KSPGetIterationNumber(mKspSolver, &num_it);
01136 std::stringstream num_it_str;
01137 num_it_str << num_it;
01138 PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01139
01140 KSPSetFromOptions(mKspSolver);
01141 KSPSetUp(mKspSolver);
01142
01143 mForceSpectrumReevaluation=false;
01144 }
01145
01146 mNumSolves++;
01147
01148 }
01149 catch (const Exception& e)
01150 {
01151
01152 PetscTools::Destroy(lhs_vector);
01153 throw e;
01154 }
01155
01156 return lhs_vector;
01157 }
01158
01159 void LinearSystem::SetPrecondMatrixIsDifferentFromLhs(bool precondIsDifferent)
01160 {
01161 mPrecondMatrixIsNotLhs = precondIsDifferent;
01162
01163 if (mPrecondMatrixIsNotLhs)
01164 {
01165 if (mRowPreallocation == UINT_MAX)
01166 {
01167
01168
01169
01170
01171
01172
01173
01174
01175
01176 NEVER_REACHED;
01177 }
01178
01179 PetscInt local_size = mOwnershipRangeHi - mOwnershipRangeLo;
01180 PetscTools::SetupMat(mPrecondMatrix, mSize, mSize, mRowPreallocation, local_size, local_size);
01181 }
01182 }
01183
01184 void LinearSystem::SetUseFixedNumberIterations(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
01185 {
01186
01187 mUseFixedNumberIterations = useFixedNumberIterations;
01188 mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
01189 }
01190
01191 void LinearSystem::ResetKspSolver()
01192 {
01193 if (mKspIsSetup)
01194 {
01195 KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
01196 }
01197
01198 mKspIsSetup = false;
01199 mForceSpectrumReevaluation = true;
01200
01201
01202
01203
01204
01205
01207 std::stringstream num_it_str;
01208 num_it_str << 1000;
01209 PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01210 }
01211
01212
01213 #include "SerializationExportWrapperForCpp.hpp"
01214 CHASTE_CLASS_EXPORT(LinearSystem)