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 #if ( PETSC_VERSION_MINOR>=5 )
00288 KSPConvergedDefaultDestroy(mpConvergenceTestContext);
00289 #else
00290 KSPDefaultConvergedDestroy(mpConvergenceTestContext);
00291 #endif
00292 }
00293 #endif
00294
00295 #ifdef TRACE_KSP
00296 if (PetscTools::AmMaster())
00297 {
00298 if (mNumSolves > 0)
00299 {
00300 double ave_num_iterations = mTotalNumIterations/(double)mNumSolves;
00301
00302 std::cout << std::endl << "KSP iterations report:" << std::endl;
00303 std::cout << "mNumSolves" << "\t" << "mTotalNumIterations" << "\t" << "mMaxNumIterations" << "\t" << "mAveNumIterations" << std::endl;
00304 std::cout << mNumSolves << "\t" << mTotalNumIterations << "\t" << mMaxNumIterations << "\t" << ave_num_iterations << std::endl;
00305 }
00306 }
00307 #endif
00308 }
00309
00310 void LinearSystem::SetMatrixElement(PetscInt row, PetscInt col, double value)
00311 {
00312 PetscMatTools::SetElement(mLhsMatrix, row, col, value);
00313 }
00314
00315 void LinearSystem::AddToMatrixElement(PetscInt row, PetscInt col, double value)
00316 {
00317 PetscMatTools::AddToElement(mLhsMatrix, row, col, value);
00318 }
00319
00320 void LinearSystem::AssembleFinalLinearSystem()
00321 {
00322 FinaliseLhsMatrix();
00323 FinaliseRhsVector();
00324 }
00325
00326 void LinearSystem::AssembleIntermediateLinearSystem()
00327 {
00328 SwitchWriteModeLhsMatrix();
00329 FinaliseRhsVector();
00330 }
00331
00332 void LinearSystem::FinaliseLhsMatrix()
00333 {
00334 PetscMatTools::Finalise(mLhsMatrix);
00335 }
00336
00337 void LinearSystem::SwitchWriteModeLhsMatrix()
00338 {
00339 PetscMatTools::SwitchWriteMode(mLhsMatrix);
00340 }
00341
00342 void LinearSystem::FinalisePrecondMatrix()
00343 {
00344 PetscMatTools::Finalise(mPrecondMatrix);
00345 }
00346
00347 void LinearSystem::FinaliseRhsVector()
00348 {
00349 PetscVecTools::Finalise(mRhsVector);
00350 }
00351
00352 void LinearSystem::SetRhsVectorElement(PetscInt row, double value)
00353 {
00354 PetscVecTools::SetElement(mRhsVector, row, value);
00355 }
00356
00357 void LinearSystem::AddToRhsVectorElement(PetscInt row, double value)
00358 {
00359 PetscVecTools::AddToElement(mRhsVector, row, value);
00360 }
00361
00362 void LinearSystem::DisplayMatrix()
00363 {
00364 PetscMatTools::Display(mLhsMatrix);
00365 }
00366
00367 void LinearSystem::DisplayRhs()
00368 {
00369 PetscVecTools::Display(mRhsVector);
00370 }
00371
00372 void LinearSystem::SetMatrixRow(PetscInt row, double value)
00373 {
00374 PetscMatTools::SetRow(mLhsMatrix, row, value);
00375 }
00376
00377 Vec LinearSystem::GetMatrixRowDistributed(unsigned rowIndex)
00378 {
00379 return PetscMatTools::GetMatrixRowDistributed(mLhsMatrix, rowIndex);
00380 }
00381
00382 void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal(std::vector<unsigned>& rRows, double diagonalValue)
00383 {
00384 PetscMatTools::ZeroRowsWithValueOnDiagonal(mLhsMatrix, rRows, diagonalValue);
00385 }
00386
00387 void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector<unsigned>& rRowColIndices, double diagonalValue)
00388 {
00389 PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mLhsMatrix, rRowColIndices, diagonalValue);
00390 }
00391
00392 void LinearSystem::ZeroMatrixColumn(PetscInt col)
00393 {
00394 PetscMatTools::ZeroColumn(mLhsMatrix, col);
00395 }
00396
00397 void LinearSystem::ZeroRhsVector()
00398 {
00399 PetscVecTools::Zero(mRhsVector);
00400 }
00401
00402 void LinearSystem::ZeroLhsMatrix()
00403 {
00404 PetscMatTools::Zero(mLhsMatrix);
00405 }
00406
00407 void LinearSystem::ZeroLinearSystem()
00408 {
00409 ZeroRhsVector();
00410 ZeroLhsMatrix();
00411 }
00412
00413 unsigned LinearSystem::GetSize() const
00414 {
00415 return (unsigned) mSize;
00416 }
00417
00418 void LinearSystem::SetNullBasis(Vec nullBasis[], unsigned numberOfBases)
00419 {
00420 #ifndef NDEBUG
00421
00422 for (unsigned vec_index=0; vec_index<numberOfBases; vec_index++)
00423 {
00424 PetscReal l2_norm;
00425 VecNorm(nullBasis[vec_index], NORM_2, &l2_norm);
00426 if (fabs(l2_norm-1.0) > 1e-08)
00427 {
00428 EXCEPTION("One of the vectors in the null space is not normalised");
00429 }
00430 }
00431
00432
00433 for (unsigned vec_index=1; vec_index<numberOfBases; vec_index++)
00434 {
00435
00436 unsigned num_vectors_ahead = numberOfBases-vec_index;
00437 boost::scoped_array<PetscScalar> dot_products(new PetscScalar[num_vectors_ahead]);
00438 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00439 VecMDot(num_vectors_ahead, nullBasis[vec_index-1], &nullBasis[vec_index], dot_products.get());
00440 #else
00441 VecMDot(nullBasis[vec_index-1], num_vectors_ahead, &nullBasis[vec_index], dot_products.get());
00442 #endif
00443 for (unsigned index=0; index<num_vectors_ahead; index++)
00444 {
00445 if (fabs(dot_products[index]) > 1e-08 )
00446 {
00447 EXCEPTION("The null space is not orthogonal.");
00448 }
00449 }
00450 }
00451
00452 #endif
00453
00454 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, numberOfBases, nullBasis, &mMatNullSpace) );
00455 }
00456
00457 void LinearSystem::RemoveNullSpace()
00458 {
00459
00460 if (mMatNullSpace)
00461 {
00462 PETSCEXCEPT( MatNullSpaceDestroy(PETSC_DESTROY_PARAM(mMatNullSpace)) );
00463 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, 0, NULL, &mMatNullSpace) );
00464 if (mKspIsSetup)
00465 {
00466 PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
00467 }
00468
00469 }
00470 }
00471
00472 void LinearSystem::GetOwnershipRange(PetscInt& lo, PetscInt& hi)
00473 {
00474 lo = mOwnershipRangeLo;
00475 hi = mOwnershipRangeHi;
00476 }
00477
00478 double LinearSystem::GetMatrixElement(PetscInt row, PetscInt col)
00479 {
00480 return PetscMatTools::GetElement(mLhsMatrix, row, col);
00481 }
00482
00483 double LinearSystem::GetRhsVectorElement(PetscInt row)
00484 {
00485 return PetscVecTools::GetElement(mRhsVector, row);
00486 }
00487
00488 unsigned LinearSystem::GetNumIterations() const
00489 {
00490 assert(this->mKspIsSetup);
00491
00492 PetscInt num_its;
00493 KSPGetIterationNumber(this->mKspSolver, &num_its);
00494
00495 return (unsigned) num_its;
00496 }
00497
00498 Vec& LinearSystem::rGetRhsVector()
00499 {
00500 return mRhsVector;
00501 }
00502
00503 Vec LinearSystem::GetRhsVector() const
00504 {
00505 return mRhsVector;
00506 }
00507
00508 Mat& LinearSystem::rGetLhsMatrix()
00509 {
00510 return mLhsMatrix;
00511 }
00512
00513 Mat LinearSystem::GetLhsMatrix() const
00514 {
00515 return mLhsMatrix;
00516 }
00517
00518 Mat& LinearSystem::rGetPrecondMatrix()
00519 {
00520 if (!mPrecondMatrixIsNotLhs)
00521 {
00522 EXCEPTION("LHS matrix used for preconditioner construction");
00523 }
00524
00525 return mPrecondMatrix;
00526 }
00527
00528 Vec& LinearSystem::rGetDirichletBoundaryConditionsVector()
00529 {
00530 return mDirichletBoundaryConditionsVector;
00531 }
00532
00533 void LinearSystem::SetMatrixIsSymmetric(bool isSymmetric)
00534 {
00536
00537 if (isSymmetric)
00538 {
00539 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRIC);
00540 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL);
00541 }
00542 else
00543 {
00544
00545 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
00546 MatSetOption(mLhsMatrix, MAT_SYMMETRIC, PETSC_FALSE);
00547 MatSetOption(mLhsMatrix, MAT_STRUCTURALLY_SYMMETRIC, PETSC_FALSE);
00548 MatSetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL, PETSC_FALSE);
00549 #else
00550 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRIC);
00551 MatSetOption(mLhsMatrix, MAT_NOT_STRUCTURALLY_SYMMETRIC);
00552 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRY_ETERNAL);
00553 #endif
00554 }
00555 }
00556
00557 bool LinearSystem::IsMatrixSymmetric()
00558 {
00559 PetscBool symmetry_flag_is_set;
00560 PetscBool symmetry_flag;
00561
00562 MatIsSymmetricKnown(mLhsMatrix, &symmetry_flag_is_set, &symmetry_flag);
00563
00564
00565 return symmetry_flag_is_set && symmetry_flag;
00566 }
00567
00568 void LinearSystem::SetMatrixIsConstant(bool matrixIsConstant)
00569 {
00570 mMatrixIsConstant = matrixIsConstant;
00571 }
00572
00573 void LinearSystem::SetRelativeTolerance(double relativeTolerance)
00574 {
00575 mTolerance = relativeTolerance;
00576 mUseAbsoluteTolerance = false;
00577 if (mKspIsSetup)
00578 {
00579 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00580 }
00581 }
00582
00583 void LinearSystem::SetAbsoluteTolerance(double absoluteTolerance)
00584 {
00585 mTolerance = absoluteTolerance;
00586 mUseAbsoluteTolerance = true;
00587 if (mKspIsSetup)
00588 {
00589 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00590 }
00591 }
00592
00593 void LinearSystem::SetKspType(const char *kspType)
00594 {
00595 mKspType = kspType;
00596 if (mKspIsSetup)
00597 {
00598 KSPSetType(mKspSolver, kspType);
00599 KSPSetFromOptions(mKspSolver);
00600 }
00601 }
00602
00603 void LinearSystem::SetPcType(const char* pcType, boost::shared_ptr<std::vector<PetscInt> > pBathNodes)
00604 {
00605 mPcType = pcType;
00606 mpBathNodes = pBathNodes;
00607
00608 if (mKspIsSetup)
00609 {
00610 if (mPcType == "blockdiagonal")
00611 {
00612
00614 delete mpBlockDiagonalPC;
00615 mpBlockDiagonalPC = NULL;
00616 delete mpLDUFactorisationPC;
00617 mpLDUFactorisationPC = NULL;
00618 delete mpTwoLevelsBlockDiagonalPC;
00619 mpTwoLevelsBlockDiagonalPC = NULL;
00620
00621 mpBlockDiagonalPC = new PCBlockDiagonal(mKspSolver);
00622 }
00623 else if (mPcType == "ldufactorisation")
00624 {
00625
00627 delete mpBlockDiagonalPC;
00628 mpBlockDiagonalPC = NULL;
00629 delete mpLDUFactorisationPC;
00630 mpLDUFactorisationPC = NULL;
00631 delete mpTwoLevelsBlockDiagonalPC;
00632 mpTwoLevelsBlockDiagonalPC = NULL;
00633
00634 mpLDUFactorisationPC = new PCLDUFactorisation(mKspSolver);
00635 }
00636 else if (mPcType == "twolevelsblockdiagonal")
00637 {
00638
00640 delete mpBlockDiagonalPC;
00641 mpBlockDiagonalPC = NULL;
00642 delete mpLDUFactorisationPC;
00643 mpLDUFactorisationPC = NULL;
00644 delete mpTwoLevelsBlockDiagonalPC;
00645 mpTwoLevelsBlockDiagonalPC = NULL;
00646
00647 if (!mpBathNodes)
00648 {
00649 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
00650 }
00651 mpTwoLevelsBlockDiagonalPC = new PCTwoLevelsBlockDiagonal(mKspSolver, *mpBathNodes);
00652 }
00653 else
00654 {
00655 PC prec;
00656 KSPGetPC(mKspSolver, &prec);
00657 PCSetType(prec, pcType);
00658 }
00659 KSPSetFromOptions(mKspSolver);
00660 }
00661 }
00662
00663 Vec LinearSystem::Solve(Vec lhsGuess)
00664 {
00665
00666
00667
00668
00669
00670
00671
00672 MatInfo mat_info;
00673 MatGetInfo(mLhsMatrix, MAT_GLOBAL_SUM, &mat_info);
00674
00675 if (!mKspIsSetup)
00676 {
00677
00678 Vec chebyshev_lhs_vector = NULL;
00679
00680 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00681 mNonZerosUsed=mat_info.nz_used;
00682
00683 PC prec;
00684
00685 KSPCreate(PETSC_COMM_WORLD, &mKspSolver);
00686
00687 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00688 if (mMatrixIsConstant)
00689 {
00690
00691 KSPSetReusePreconditioner(mKspSolver, PETSC_TRUE);
00692 }
00693 #else
00694
00695
00696
00697
00698 MatStructure preconditioner_over_successive_calls;
00699
00700 if (mMatrixIsConstant)
00701 {
00702 preconditioner_over_successive_calls = SAME_PRECONDITIONER;
00703 }
00704 else
00705 {
00706 preconditioner_over_successive_calls = SAME_NONZERO_PATTERN;
00707 }
00708 #endif
00709
00710 if (mPrecondMatrixIsNotLhs)
00711 {
00712 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00713 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix);
00714 #else
00715 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix, preconditioner_over_successive_calls);
00716 #endif
00717 }
00718 else
00719 {
00720 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00721 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix);
00722 #else
00723 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix, preconditioner_over_successive_calls);
00724 #endif
00725 }
00726
00727
00728
00729 if (mUseAbsoluteTolerance)
00730 {
00731 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00732 }
00733 else
00734 {
00735 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00736 }
00737
00738
00739 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
00740 if (mKspType == "chebychev")
00741 {
00742 KSPSetType(mKspSolver, "chebyshev");
00743 }
00744 else
00745 {
00746 KSPSetType(mKspSolver, mKspType.c_str());
00747 }
00748 #else
00749 KSPSetType(mKspSolver, mKspType.c_str());
00750 #endif
00751 KSPGetPC(mKspSolver, &prec);
00752
00753
00754 if (mSize <= 6)
00755 {
00756 PCSetType(prec, PCNONE);
00757 }
00758 else
00759 {
00760 #ifdef TRACE_KSP
00761 Timer::Reset();
00762 #endif
00763 if (mPcType == "blockdiagonal")
00764 {
00765 mpBlockDiagonalPC = new PCBlockDiagonal(mKspSolver);
00766 #ifdef TRACE_KSP
00767 if (PetscTools::AmMaster())
00768 {
00769 Timer::Print("Purpose-build preconditioner creation");
00770 }
00771 #endif
00772 }
00773 else if (mPcType == "ldufactorisation")
00774 {
00775 mpLDUFactorisationPC = new PCLDUFactorisation(mKspSolver);
00776 #ifdef TRACE_KSP
00777 if (PetscTools::AmMaster())
00778 {
00779 Timer::Print("Purpose-build preconditioner creation");
00780 }
00781 #endif
00782 }
00783 else if (mPcType == "twolevelsblockdiagonal")
00784 {
00785 if (!mpBathNodes)
00786 {
00787 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
00788 }
00789 mpTwoLevelsBlockDiagonalPC = new PCTwoLevelsBlockDiagonal(mKspSolver, *mpBathNodes);
00790 #ifdef TRACE_KSP
00791 if (PetscTools::AmMaster())
00792 {
00793 Timer::Print("Purpose-build preconditioner creation");
00794 }
00795 #endif
00796
00797 }
00798 else
00799 {
00800 PCSetType(prec, mPcType.c_str());
00801 }
00802 }
00803
00804 if (mMatNullSpace)
00805 {
00807 PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
00808 }
00809
00810 KSPSetFromOptions(mKspSolver);
00811 if (lhsGuess)
00812 {
00813
00814 KSPSetInitialGuessNonzero(mKspSolver,PETSC_TRUE);
00815 }
00816
00817
00818
00819
00820 if (mKspType == "chebychev" && !mUseFixedNumberIterations)
00821 {
00822 #ifdef TRACE_KSP
00823 Timer::Reset();
00824 #endif
00825
00826 KSPSetType(mKspSolver,"cg");
00827 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
00828 KSPSetUp(mKspSolver);
00829
00830 VecDuplicate(mRhsVector, &chebyshev_lhs_vector);
00831 if (lhsGuess)
00832 {
00833 VecCopy(lhsGuess, chebyshev_lhs_vector);
00834 }
00835
00836
00837 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
00838
00839 PetscReal *r_eig = new PetscReal[mSize];
00840 PetscReal *c_eig = new PetscReal[mSize];
00841 PetscInt eigs_computed;
00842 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
00843
00844 mEigMin = r_eig[0];
00845
00846
00847 KSPSetTolerances(mKspSolver, DBL_EPSILON, DBL_EPSILON, PETSC_DEFAULT, PETSC_DEFAULT);
00848 KSPSetUp(mKspSolver);
00849 if (lhsGuess)
00850 {
00851 VecCopy(lhsGuess, chebyshev_lhs_vector);
00852 }
00853
00854 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
00855 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
00856
00857 mEigMax = r_eig[eigs_computed-1];
00858
00859 #ifdef TRACE_KSP
00860
00861
00862
00863
00864
00865 if (PetscTools::AmMaster())
00866 {
00867 std::cout << "EIGS: ";
00868 for (int index=0; index<eigs_computed; index++)
00869 {
00870 std::cout << r_eig[index] << ", ";
00871 }
00872 std::cout << std::endl;
00873 }
00874
00875 if (PetscTools::AmMaster()) std::cout << "EIGS "<< mEigMax << " " << mEigMin <<std::endl;
00876 double alpha = 2/(mEigMax+mEigMin);
00877 double sigma_one = 1 - alpha*mEigMin;
00878 if (PetscTools::AmMaster()) std::cout << "sigma_1 = 1 - alpha*mEigMin = "<< sigma_one <<std::endl;
00879 #endif
00880
00881
00882 assert(mKspType == "chebychev");
00883 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
00884 KSPSetType(mKspSolver, "chebyshev");
00885 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
00886 #else
00887 KSPSetType(mKspSolver, mKspType.c_str());
00888 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
00889 #endif
00890 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
00891 if (mUseAbsoluteTolerance)
00892 {
00893 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
00894 }
00895 else
00896 {
00897 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
00898 }
00899
00900 delete[] r_eig;
00901 delete[] c_eig;
00902
00903 #ifdef TRACE_KSP
00904 if (PetscTools::AmMaster())
00905 {
00906 Timer::Print("Computing extremal eigenvalues");
00907 }
00908 #endif
00909 }
00910
00911 #ifdef TRACE_KSP
00912 Timer::Reset();
00913 #endif
00914
00915 KSPSetUp(mKspSolver);
00916
00917 if (chebyshev_lhs_vector)
00918 {
00919 PetscTools::Destroy(chebyshev_lhs_vector);
00920 }
00921
00922 #ifdef TRACE_KSP
00923 if (PetscTools::AmMaster())
00924 {
00925 Timer::Print("KSPSetUP (contains preconditioner creation for PETSc preconditioners)");
00926 }
00927 #endif
00928
00929 mKspIsSetup = true;
00930
00931 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00932 }
00933 else
00934 {
00935 if (mNonZerosUsed!=mat_info.nz_used)
00936 {
00937 WARNING("LinearSystem doesn't like the non-zero pattern of a matrix to change. (I think you changed it).");
00938 mNonZerosUsed = mat_info.nz_used;
00939 }
00940
00941
00942
00943
00944
00945
00946 }
00947
00948 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00949
00951 Vec lhs_vector;
00952 VecDuplicate(mRhsVector, &lhs_vector);
00953 if (lhsGuess)
00954 {
00955 VecCopy(lhsGuess, lhs_vector);
00956
00957
00958
00959 }
00960
00961
00962
00963 PetscReal rhs_norm;
00964 VecNorm(mRhsVector, NORM_1, &rhs_norm);
00965 double rhs_zero_tol = 1e-15;
00966
00967 if (rhs_norm < rhs_zero_tol && rhs_norm > 0.0)
00968 {
00969 WARNING("Using zero initial guess due to small right hand side vector");
00970 PetscVecTools::Zero(lhs_vector);
00971 }
00972
00973 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
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
01000
01001
01002
01003
01004
01005
01006
01007
01008
01009
01010
01011
01012
01013
01014 try
01015 {
01016 HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
01017
01018 #ifdef TRACE_KSP
01019 Timer::Reset();
01020 #endif
01021
01022
01023 if (mUseFixedNumberIterations && (mNumSolves%mEvaluateNumItsEveryNSolves==0 || mForceSpectrumReevaluation))
01024 {
01025 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
01026 KSPSetNormType(mKspSolver, KSP_PRECONDITIONED_NORM);
01027 #else
01028 KSPSetNormType(mKspSolver, KSP_NORM_PRECONDITIONED);
01029 #endif
01030
01031 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
01032 if (!mpConvergenceTestContext)
01033 {
01034 #if ( PETSC_VERSION_MINOR>=5 )
01035 KSPConvergedDefaultCreate(&mpConvergenceTestContext);
01036 #else
01037 KSPDefaultConvergedCreate(&mpConvergenceTestContext);
01038 #endif
01039 }
01040 #if ( PETSC_VERSION_MINOR>=5 )
01041
01042 KSPSetConvergenceTest(mKspSolver, KSPConvergedDefault, &mpConvergenceTestContext, PETSC_NULL);
01043 #else
01044 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, &mpConvergenceTestContext, PETSC_NULL);
01045 #endif
01046
01047 #else
01048 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, PETSC_NULL);
01049 #endif
01050
01051 if (mUseAbsoluteTolerance)
01052 {
01053 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
01054 }
01055 else
01056 {
01057 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
01058 }
01059
01061 std::stringstream num_it_str;
01062 num_it_str << 1000;
01063 PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01064
01065
01066 if (mKspType == "chebychev")
01067 {
01068
01069 KSPSetType(mKspSolver,"cg");
01070 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
01071 }
01072
01073 KSPSetFromOptions(mKspSolver);
01074 KSPSetUp(mKspSolver);
01075 }
01076
01077 PETSCEXCEPT(KSPSolve(mKspSolver, mRhsVector, lhs_vector));
01078 HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
01079
01080 #ifdef TRACE_KSP
01081 PetscInt num_it;
01082 KSPGetIterationNumber(mKspSolver, &num_it);
01083 if (PetscTools::AmMaster())
01084 {
01085 std::cout << "++ Solve: " << mNumSolves << " NumIterations: " << num_it << " ";
01086 Timer::Print("Solve");
01087 }
01088
01089 mTotalNumIterations += num_it;
01090 if ((unsigned) num_it > mMaxNumIterations)
01091 {
01092 mMaxNumIterations = num_it;
01093 }
01094 #endif
01095
01096
01097 KSPConvergedReason reason;
01098 KSPGetConvergedReason(mKspSolver, &reason);
01099
01100 if (mUseFixedNumberIterations && PETSC_VERSION_MAJOR < 3)
01101 {
01102 WARNING("Not explicitly checking convergence reason when using fixed number of iterations and PETSc 2");
01103 }
01104 else
01105 {
01106 KSPEXCEPT(reason);
01107 }
01108
01109 if (mUseFixedNumberIterations && (mNumSolves%mEvaluateNumItsEveryNSolves==0 || mForceSpectrumReevaluation))
01110 {
01111
01112 if (mKspType == "chebychev")
01113 {
01114 PetscReal *r_eig = new PetscReal[mSize];
01115 PetscReal *c_eig = new PetscReal[mSize];
01116 PetscInt eigs_computed;
01117 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
01118
01119 mEigMin = r_eig[0];
01120
01121
01122
01123
01124
01125 mEigMax = std::max(mEigMax,r_eig[eigs_computed-1]);
01126
01127 delete[] r_eig;
01128 delete[] c_eig;
01129
01130 assert(mKspType == "chebychev");
01131 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
01132 KSPSetType(mKspSolver, "chebyshev");
01133 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
01134 #else
01135 KSPSetType(mKspSolver, mKspType.c_str());
01136 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
01137 #endif
01138 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
01139 }
01140
01141 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
01142 if (mKspType == "chebychev")
01143 {
01144
01145 EXCEPTION("Chebyshev with fixed number of iterations is known to be broken in PETSc <= 2.3.2");
01146 }
01147
01148 KSPSetNormType(mKspSolver, KSP_NO_NORM);
01149 #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
01150 KSPSetNormType(mKspSolver, KSP_NORM_NONE);
01151 #else
01152 KSPSetNormType(mKspSolver, KSP_NORM_NO);
01153 #endif
01154
01155 #if (PETSC_VERSION_MAJOR == 2)
01156 KSPSetConvergenceTest(mKspSolver, KSPSkipConverged, PETSC_NULL);
01157 #endif
01158
01159 PetscInt num_it;
01160 KSPGetIterationNumber(mKspSolver, &num_it);
01161 std::stringstream num_it_str;
01162 num_it_str << num_it;
01163 PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01164
01165 KSPSetFromOptions(mKspSolver);
01166 KSPSetUp(mKspSolver);
01167
01168 mForceSpectrumReevaluation=false;
01169 }
01170
01171 mNumSolves++;
01172
01173 }
01174 catch (const Exception& e)
01175 {
01176
01177 PetscTools::Destroy(lhs_vector);
01178 throw e;
01179 }
01180
01181 return lhs_vector;
01182 }
01183
01184 void LinearSystem::SetPrecondMatrixIsDifferentFromLhs(bool precondIsDifferent)
01185 {
01186 mPrecondMatrixIsNotLhs = precondIsDifferent;
01187
01188 if (mPrecondMatrixIsNotLhs)
01189 {
01190 if (mRowPreallocation == UINT_MAX)
01191 {
01192
01193
01194
01195
01196
01197
01198
01199
01200
01201 NEVER_REACHED;
01202 }
01203
01204 PetscInt local_size = mOwnershipRangeHi - mOwnershipRangeLo;
01205 PetscTools::SetupMat(mPrecondMatrix, mSize, mSize, mRowPreallocation, local_size, local_size);
01206 }
01207 }
01208
01209 void LinearSystem::SetUseFixedNumberIterations(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
01210 {
01211
01212 mUseFixedNumberIterations = useFixedNumberIterations;
01213 mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
01214 }
01215
01216 void LinearSystem::ResetKspSolver()
01217 {
01218 if (mKspIsSetup)
01219 {
01220 KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
01221 }
01222
01223 mKspIsSetup = false;
01224 mForceSpectrumReevaluation = true;
01225
01226
01227
01228
01229
01230
01232 std::stringstream num_it_str;
01233 num_it_str << 1000;
01234 PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
01235 }
01236
01237
01238 #include "SerializationExportWrapperForCpp.hpp"
01239 CHASTE_CLASS_EXPORT(LinearSystem)