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