Chaste Commit::ed021841b4f5a92f33f7f94a19c918edfb2bb29d
LinearSystem.cpp
1/*
2
3Copyright (c) 2005-2026, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, are permitted provided that the following conditions are met:
14 * Redistributions of source code must retain the above copyright notice,
15 this list of conditions and the following disclaimer.
16 * Redistributions in binary form must reproduce the above copyright notice,
17 this list of conditions and the following disclaimer in the documentation
18 and/or other materials provided with the distribution.
19 * Neither the name of the University of Oxford nor the names of its
20 contributors may be used to endorse or promote products derived from this
21 software without specific prior written permission.
22
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#include "LinearSystem.hpp"
37
38#include <iostream>
39#include <cassert>
40#include <algorithm>
41#include <boost/scoped_array.hpp>
42
43#include "PetscException.hpp"
44#include "OutputFileHandler.hpp"
45#include "PetscTools.hpp"
46#include "HeartEventHandler.hpp"
47#include "Timer.hpp"
48#include "Warnings.hpp"
49
51// Implementation
53
54LinearSystem::LinearSystem(PetscInt lhsVectorSize, unsigned rowPreallocation)
55 :mPrecondMatrix(nullptr),
56 mSize(lhsVectorSize),
57 mMatNullSpace(nullptr),
58 mDestroyMatAndVec(true),
59 mKspIsSetup(false),
60 mNonZerosUsed(0.0),
61 mMatrixIsConstant(false),
62 mTolerance(1e-6),
63 mUseAbsoluteTolerance(false),
64 mDirichletBoundaryConditionsVector(nullptr),
65 mpBlockDiagonalPC(nullptr),
66 mpLDUFactorisationPC(nullptr),
67 mpTwoLevelsBlockDiagonalPC(nullptr),
68 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
69 mPrecondMatrixIsNotLhs(false),
70 mRowPreallocation(rowPreallocation),
71 mUseFixedNumberIterations(false),
72 mEvaluateNumItsEveryNSolves(UINT_MAX),
73 mpConvergenceTestContext(nullptr),
74 mEigMin(DBL_MAX),
75 mEigMax(DBL_MIN),
76 mForceSpectrumReevaluation(false)
77{
78 assert(lhsVectorSize > 0);
79 if (mRowPreallocation == UINT_MAX)
80 {
81 // Automatic preallocation if it's a small matrix
82 if (lhsVectorSize<15)
83 {
84 mRowPreallocation=lhsVectorSize;
85 }
86 else
87 {
88 EXCEPTION("You must provide a rowPreallocation argument for a large sparse system");
89 }
90 }
91
93 PetscTools::SetupMat(mLhsMatrix, mSize, mSize, mRowPreallocation, PETSC_DECIDE, PETSC_DECIDE);
94
95 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
96
99 mKspType = "gmres";
100 mPcType = "jacobi";
101
102 mNumSolves = 0;
103#ifdef TRACE_KSP
104 mTotalNumIterations = 0;
105 mMaxNumIterations = 0;
106#endif
107}
108
109LinearSystem::LinearSystem(PetscInt lhsVectorSize, Mat lhsMatrix, Vec rhsVector)
110 :mPrecondMatrix(nullptr),
111 mSize(lhsVectorSize),
112 mMatNullSpace(nullptr),
113 mDestroyMatAndVec(true),
114 mKspIsSetup(false),
115 mNonZerosUsed(0.0),
116 mMatrixIsConstant(false),
117 mTolerance(1e-6),
118 mUseAbsoluteTolerance(false),
119 mDirichletBoundaryConditionsVector(nullptr),
120 mpBlockDiagonalPC(nullptr),
121 mpLDUFactorisationPC(nullptr),
122 mpTwoLevelsBlockDiagonalPC(nullptr),
123 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
124 mPrecondMatrixIsNotLhs(false),
125 mUseFixedNumberIterations(false),
126 mEvaluateNumItsEveryNSolves(UINT_MAX),
127 mpConvergenceTestContext(nullptr),
128 mEigMin(DBL_MAX),
129 mEigMax(DBL_MIN),
130 mForceSpectrumReevaluation(false)
131{
132 assert(lhsVectorSize > 0);
133 // Conveniently, PETSc Mats and Vecs are actually pointers
134 mLhsMatrix = lhsMatrix;
135 mRhsVector = rhsVector;
136
137 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
138
139 mNumSolves = 0;
140#ifdef TRACE_KSP
141 mTotalNumIterations = 0;
142 mMaxNumIterations = 0;
143#endif
144}
145
146LinearSystem::LinearSystem(Vec templateVector, unsigned rowPreallocation, bool newAllocationError)
147 :mPrecondMatrix(nullptr),
148 mMatNullSpace(nullptr),
149 mDestroyMatAndVec(true),
150 mKspIsSetup(false),
151 mMatrixIsConstant(false),
152 mTolerance(1e-6),
153 mUseAbsoluteTolerance(false),
154 mDirichletBoundaryConditionsVector(nullptr),
155 mpBlockDiagonalPC(nullptr),
156 mpLDUFactorisationPC(nullptr),
157 mpTwoLevelsBlockDiagonalPC(nullptr),
158 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
159 mPrecondMatrixIsNotLhs(false),
160 mRowPreallocation(rowPreallocation),
161 mUseFixedNumberIterations(false),
162 mEvaluateNumItsEveryNSolves(UINT_MAX),
163 mpConvergenceTestContext(nullptr),
164 mEigMin(DBL_MAX),
165 mEigMax(DBL_MIN),
166 mForceSpectrumReevaluation(false)
167{
168 VecDuplicate(templateVector, &mRhsVector);
169 VecGetSize(mRhsVector, &mSize);
170 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
172
173 PetscTools::SetupMat(mLhsMatrix, mSize, mSize, mRowPreallocation, local_size, local_size, true, newAllocationError);
174
177 mKspType = "gmres";
178 mPcType = "jacobi";
179
180 mNumSolves = 0;
181#ifdef TRACE_KSP
182 mTotalNumIterations = 0;
183 mMaxNumIterations = 0;
184#endif
185}
186
187LinearSystem::LinearSystem(Vec residualVector, Mat jacobianMatrix)
188 :mPrecondMatrix(nullptr),
189 mMatNullSpace(nullptr),
190 mDestroyMatAndVec(false),
191 mKspIsSetup(false),
192 mMatrixIsConstant(false),
193 mTolerance(1e-6),
194 mUseAbsoluteTolerance(false),
195 mDirichletBoundaryConditionsVector(nullptr),
196 mpBlockDiagonalPC(nullptr),
197 mpLDUFactorisationPC(nullptr),
198 mpTwoLevelsBlockDiagonalPC(nullptr),
199 mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
200 mPrecondMatrixIsNotLhs(false),
201 mRowPreallocation(UINT_MAX),
202 mUseFixedNumberIterations(false),
203 mEvaluateNumItsEveryNSolves(UINT_MAX),
204 mpConvergenceTestContext(nullptr),
205 mEigMin(DBL_MAX),
206 mEigMax(DBL_MIN),
207 mForceSpectrumReevaluation(false)
208{
209 assert(residualVector || jacobianMatrix);
210 mRhsVector = residualVector;
211 mLhsMatrix = jacobianMatrix;
212
213 PetscInt mat_size=0, vec_size=0;
214 if (mRhsVector)
215 {
216 VecGetSize(mRhsVector, &vec_size);
217 mSize = (unsigned)vec_size;
218 VecGetOwnershipRange(mRhsVector, &mOwnershipRangeLo, &mOwnershipRangeHi);
219 }
220 if (mLhsMatrix)
221 {
222 PetscInt mat_cols;
223 MatGetSize(mLhsMatrix, &mat_size, &mat_cols);
224 assert(mat_size == mat_cols);
225 mSize = (unsigned)mat_size;
226 MatGetOwnershipRange(mLhsMatrix, &mOwnershipRangeLo, &mOwnershipRangeHi);
227
228 MatInfo matrix_info;
229 MatGetInfo(mLhsMatrix, MAT_GLOBAL_MAX, &matrix_info);
230
231 /*
232 * Assuming that mLhsMatrix was created with PetscTools::SetupMat, the value
233 * below should be equivalent to what was used as preallocation in that call.
234 */
235 mRowPreallocation = (unsigned) matrix_info.nz_allocated / mSize;
236 }
237 assert(!mRhsVector || !mLhsMatrix || vec_size == mat_size);
238
241 mKspType = "gmres";
242 mPcType = "jacobi";
243
244 mNumSolves = 0;
245#ifdef TRACE_KSP
246 mTotalNumIterations = 0;
247 mMaxNumIterations = 0;
248#endif
249}
250
252{
253 delete mpBlockDiagonalPC;
256
258 {
261 }
262
264 {
266 }
267
268 if (mMatNullSpace)
269 {
270 MatNullSpaceDestroy(PETSC_DESTROY_PARAM(mMatNullSpace));
271 }
272
273 if (mKspIsSetup)
274 {
275 KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
276 }
277
279 {
282 }
283
284#if (PETSC_VERSION_MAJOR == 3) // PETSc 3.x.x
286 {
287#if (PETSC_VERSION_MINOR >= 24) // PETSc 3.24 or later
288 KSPConvergedDefaultDestroy(&mpConvergenceTestContext);
289#elif (PETSC_VERSION_MINOR >= 5) // PETSc 3.5 to 3.23
290 KSPConvergedDefaultDestroy(mpConvergenceTestContext);
291#else
292 KSPDefaultConvergedDestroy(mpConvergenceTestContext);
293#endif
294 mpConvergenceTestContext = nullptr;
295 }
296#endif
297
298#ifdef TRACE_KSP
300 {
301 if (mNumSolves > 0)
302 {
303 double ave_num_iterations = mTotalNumIterations/(double)mNumSolves;
304
305 std::cout << std::endl << "KSP iterations report:" << std::endl;
306 std::cout << "mNumSolves" << "\t" << "mTotalNumIterations" << "\t" << "mMaxNumIterations" << "\t" << "mAveNumIterations" << std::endl;
307 std::cout << mNumSolves << "\t" << mTotalNumIterations << "\t" << mMaxNumIterations << "\t" << ave_num_iterations << std::endl;
308 }
309 }
310#endif
311}
312
314{
315 PetscMatTools::SetElement(mLhsMatrix, row, col, value);
316}
317
319{
320 PetscMatTools::AddToElement(mLhsMatrix, row, col, value);
321}
322
328
334
339
344
349
354
356{
358}
359
361{
363}
364
369
374
376{
378}
379
384
385void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal(std::vector<unsigned>& rRows, double diagonalValue)
386{
388}
389
390void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector<unsigned>& rRowColIndices, double diagonalValue)
391{
393}
394
399
404
409
415
416unsigned LinearSystem::GetSize() const
417{
418 return (unsigned) mSize;
419}
420
421void LinearSystem::SetNullBasis(Vec nullBasis[], unsigned numberOfBases)
422{
423#ifndef NDEBUG
424 // Check all the vectors of the base are normal
425 for (unsigned vec_index=0; vec_index<numberOfBases; vec_index++)
426 {
427 PetscReal l2_norm;
428 VecNorm(nullBasis[vec_index], NORM_2, &l2_norm);
429 if (fabs(l2_norm-1.0) > 1e-08)
430 {
431 EXCEPTION("One of the vectors in the null space is not normalised");
432 }
433 }
434
435 // Check all the vectors of the base are orthogonal
436 for (unsigned vec_index=1; vec_index<numberOfBases; vec_index++)
437 {
438 // The strategy is to check the (i-1)-th vector against vectors from i to n with VecMDot. This should be the most efficient way of doing it.
439 unsigned num_vectors_ahead = numberOfBases-vec_index;
440 boost::scoped_array<PetscScalar> dot_products(new PetscScalar[num_vectors_ahead]);
441#if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
442 VecMDot(num_vectors_ahead, nullBasis[vec_index-1], &nullBasis[vec_index], dot_products.get());
443#else
444 VecMDot(nullBasis[vec_index-1], num_vectors_ahead, &nullBasis[vec_index], dot_products.get());
445#endif
446 for (unsigned index=0; index<num_vectors_ahead; index++)
447 {
448 if (fabs(dot_products[index]) > 1e-08 )
449 {
450 EXCEPTION("The null space is not orthogonal.");
451 }
452 }
453 }
454
455#endif
456
457 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, numberOfBases, nullBasis, &mMatNullSpace) );
458}
459
461{
462 // Only remove if previously set
463 if (mMatNullSpace)
464 {
465 PETSCEXCEPT( MatNullSpaceDestroy(PETSC_DESTROY_PARAM(mMatNullSpace)) );
466 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, 0, nullptr, &mMatNullSpace) );
467#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
468 // Setting null space in the KSP was deprecated in PETSc 3.6, but setting the null space
469 // for the matrix appeared in PETSc 3.3 so 3.3, 3.4, 3.5 can do either
470 PETSCEXCEPT( MatSetNullSpace(mLhsMatrix, mMatNullSpace) );
471 /* This is heavy-handed:
472 * Changing the Null-space appears to only work if we reset the KSP
473 * Adding null-space to the matrix has to happen *before* KSPSetOperators
474 */
476#else
477 if (mKspIsSetup)
478 {
479 PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
480 }
481 //else: it will be set next time Solve() is called
482#endif
483 }
484}
485
491
496
501
503{
504 assert(this->mKspIsSetup);
505
506 PetscInt num_its;
507 KSPGetIterationNumber(this->mKspSolver, &num_its);
508
509 return (unsigned) num_its;
510}
511
516
518{
519 return mRhsVector;
520}
521
526
528{
529 return mLhsMatrix;
530}
531
533{
535 {
536 EXCEPTION("LHS matrix used for preconditioner construction");
537 }
538
539 return mPrecondMatrix;
540}
541
546
548{
550
551 if (isSymmetric)
552 {
553 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRIC);
554 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL);
555 }
556 else
557 {
558// don't have a PetscMatTools method for setting options to false
559#if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
560 MatSetOption(mLhsMatrix, MAT_SYMMETRIC, PETSC_FALSE);
561 MatSetOption(mLhsMatrix, MAT_STRUCTURALLY_SYMMETRIC, PETSC_FALSE);
562 MatSetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL, PETSC_FALSE);
563#else
564 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRIC);
565 MatSetOption(mLhsMatrix, MAT_NOT_STRUCTURALLY_SYMMETRIC);
566 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRY_ETERNAL);
567#endif
568 }
569}
570
572{
573 PetscBool symmetry_flag_is_set;
574 PetscBool symmetry_flag;
575
576 MatIsSymmetricKnown(mLhsMatrix, &symmetry_flag_is_set, &symmetry_flag);
577
578 // If the flag is not set we assume is a non-symmetric matrix
579 return symmetry_flag_is_set && symmetry_flag;
580}
581
582void LinearSystem::SetMatrixIsConstant(bool matrixIsConstant)
583{
584 mMatrixIsConstant = matrixIsConstant;
585}
586
587void LinearSystem::SetRelativeTolerance(double relativeTolerance)
588{
589 mTolerance = relativeTolerance;
590 mUseAbsoluteTolerance = false;
591 if (mKspIsSetup)
592 {
593 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
594 }
595}
596
597void LinearSystem::SetAbsoluteTolerance(double absoluteTolerance)
598{
599 mTolerance = absoluteTolerance;
601 if (mKspIsSetup)
602 {
603 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
604 }
605}
606
607void LinearSystem::SetKspType(const char *kspType)
608{
609 mKspType = kspType;
610 if (mKspIsSetup)
611 {
612 KSPSetType(mKspSolver, kspType);
613 KSPSetFromOptions(mKspSolver);
614 }
615}
616
617void LinearSystem::SetPcType(const char* pcType, boost::shared_ptr<std::vector<PetscInt> > pBathNodes)
618{
619 mPcType = pcType;
620 mpBathNodes = pBathNodes;
621
622 if (mKspIsSetup)
623 {
624 if (mPcType == "blockdiagonal")
625 {
626 // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
628 delete mpBlockDiagonalPC;
629 mpBlockDiagonalPC = nullptr;
631 mpLDUFactorisationPC = nullptr;
634
636 }
637 else if (mPcType == "ldufactorisation")
638 {
639 // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
641 delete mpBlockDiagonalPC;
642 mpBlockDiagonalPC = nullptr;
644 mpLDUFactorisationPC = nullptr;
647
649 }
650 else if (mPcType == "twolevelsblockdiagonal")
651 {
652 // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
654 delete mpBlockDiagonalPC;
655 mpBlockDiagonalPC = nullptr;
657 mpLDUFactorisationPC = nullptr;
660
661 if (!mpBathNodes)
662 {
663 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC"); // LCOV_EXCL_LINE
664 }
666 }
667 else
668 {
669 PC prec;
670 KSPGetPC(mKspSolver, &prec);
671 PCSetType(prec, pcType);
672 }
673 KSPSetFromOptions(mKspSolver);
674 }
675}
676
678{
679 /*
680 * The following lines are very useful for debugging:
681 * MatView(mLhsMatrix, PETSC_VIEWER_STDOUT_WORLD);
682 * VecView(mRhsVector, PETSC_VIEWER_STDOUT_WORLD);
683 */
684
685 // Double check that the non-zero pattern hasn't changed
686 MatInfo mat_info;
687 MatGetInfo(mLhsMatrix, MAT_GLOBAL_SUM, &mat_info);
688
689 if (!mKspIsSetup)
690 {
691 // Create PETSc Vec that may be required if we use a Chebyshev solver
692 Vec chebyshev_lhs_vector = nullptr;
693
694 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
695 mNonZerosUsed=mat_info.nz_used;
696 //MatNorm(mLhsMatrix, NORM_FROBENIUS, &mMatrixNorm);
697 PC prec; //Type of pre-conditioner
698
699 KSPCreate(PETSC_COMM_WORLD, &mKspSolver);
700
701 if (mMatNullSpace) // Adding null-space to the matrix (new style) has to happen *before* KSPSetOperators
702 {
703#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
704 // Setting null space in the KSP was deprecated in PETSc 3.6, but setting the null space
705 // for the matrix appeared in PETSc 3.3 so 3.3, 3.4, 3.5 can do either
706
707 PETSCEXCEPT(MatSetNullSpace(mLhsMatrix, mMatNullSpace));
708#else
709 PETSCEXCEPT(KSPSetNullSpace(mKspSolver, mMatNullSpace));
710#endif
711 }
712#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
713 // Do nothing. Note that reusing the pre-conditioner in later PETSc versions is done after the pre-conditioner is formed
714 // (This comment is retained here so that the #if logic is consistent.)
715#else
716 /*
717 * The preconditioner flag (last argument) in the following calls says
718 * how to reuse the preconditioner on subsequent iterations.
719 */
720 MatStructure preconditioner_over_successive_calls;
721
723 {
724 preconditioner_over_successive_calls = SAME_PRECONDITIONER;
725 }
726 else
727 {
728 preconditioner_over_successive_calls = SAME_NONZERO_PATTERN;
729 }
730#endif
731
733 {
734#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
735 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix);
736#else
737 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix, preconditioner_over_successive_calls);
738#endif
739 }
740 else
741 {
742#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
743 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix);
744#else
745 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix, preconditioner_over_successive_calls);
746#endif
747 }
748#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
750 {
751 // Emulate SAME_PRECONDITIONER as above
752 KSPSetReusePreconditioner(mKspSolver, PETSC_TRUE);
753 }
754#endif
755
756 // Set either absolute or relative tolerance of the KSP solver.
757 // The default is to use relative tolerance (1e-6)
759 {
760 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
761 }
762 else
763 {
764 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
765 }
766
767 // Set ksp and pc types
768#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
769 if (mKspType == "chebychev")
770 {
771 KSPSetType(mKspSolver, "chebyshev");
772 }
773 else
774 {
775 KSPSetType(mKspSolver, mKspType.c_str());
776 }
777#else
778 KSPSetType(mKspSolver, mKspType.c_str());
779#endif
780 KSPGetPC(mKspSolver, &prec);
781
782 // Turn off pre-conditioning if the system size is very small
783 const bool is_small = (mSize <= 6);
784 if (is_small)
785 {
786 PCSetType(prec, PCNONE);
787 }
788 else
789 {
790#ifdef TRACE_KSP
791 Timer::Reset();
792#endif
793 if (mPcType == "blockdiagonal")
794 {
796#ifdef TRACE_KSP
798 {
799 Timer::Print("Purpose-build preconditioner creation");
800 }
801#endif
802 }
803 else if (mPcType == "ldufactorisation")
804 {
806#ifdef TRACE_KSP
808 {
809 Timer::Print("Purpose-build preconditioner creation");
810 }
811#endif
812 }
813 else if (mPcType == "twolevelsblockdiagonal")
814 {
815 if (!mpBathNodes)
816 {
817 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC"); // LCOV_EXCL_LINE
818 }
820#ifdef TRACE_KSP
822 {
823 Timer::Print("Purpose-build preconditioner creation");
824 }
825#endif
826
827 }
828 else
829 {
830 PCSetType(prec, mPcType.c_str());
831 }
832 }
833
834 KSPSetFromOptions(mKspSolver);
835
836 if (lhsGuess)
837 {
838 // Assume that the user of this method will always be kind enough to give us a reasonable guess.
839 KSPSetInitialGuessNonzero(mKspSolver,PETSC_TRUE);
840 }
841 /*
842 * Non-adaptive Chebyshev: the required spectrum approximation is computed just once
843 * at the beginning of the simulation. This is done with two extra CG solves.
844 */
845 if (mKspType == "chebychev" && !mUseFixedNumberIterations)
846 {
847#ifdef TRACE_KSP
848 Timer::Reset();
849#endif
850 // Preconditioned matrix spectrum is approximated with CG
851 KSPSetType(mKspSolver,"cg");
852 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
853 KSPSetUp(mKspSolver);
854
855 VecDuplicate(mRhsVector, &chebyshev_lhs_vector);
856 if (lhsGuess)
857 {
858 VecCopy(lhsGuess, chebyshev_lhs_vector);
859 }
860
861 // Smallest eigenvalue is approximated to default tolerance
862 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
863
864 PetscReal *r_eig = new PetscReal[mSize];
865 PetscReal *c_eig = new PetscReal[mSize];
866 PetscInt eigs_computed;
867 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
868
869 mEigMin = r_eig[0];
870
871 // Largest eigenvalue is approximated to machine precision
872 KSPSetTolerances(mKspSolver, DBL_EPSILON, DBL_EPSILON, PETSC_DEFAULT, PETSC_DEFAULT);
873 KSPSetUp(mKspSolver);
874 if (lhsGuess)
875 {
876 VecCopy(lhsGuess, chebyshev_lhs_vector);
877 }
878
879 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
880 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
881
882 mEigMax = r_eig[eigs_computed-1];
883
884#ifdef TRACE_KSP
885 /*
886 * Under certain circumstances (see Golub&Overton 1988), underestimating
887 * the spectrum of the preconditioned operator improves convergence rate.
888 * See publication for a discussion and for definition of alpha and sigma_one.
889 */
891 {
892 std::cout << "EIGS: ";
893 for (int index=0; index<eigs_computed; index++)
894 {
895 std::cout << r_eig[index] << ", ";
896 }
897 std::cout << std::endl;
898 }
899
900 if (PetscTools::AmMaster()) std::cout << "EIGS "<< mEigMax << " " << mEigMin <<std::endl;
901 double alpha = 2/(mEigMax+mEigMin);
902 double sigma_one = 1 - alpha*mEigMin;
903 if (PetscTools::AmMaster()) std::cout << "sigma_1 = 1 - alpha*mEigMin = "<< sigma_one <<std::endl;
904#endif
905
906 // Set Chebyshev solver and max/min eigenvalues
907 assert(mKspType == "chebychev");
908#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
909 KSPSetType(mKspSolver, "chebyshev");
910 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
911#else
912 KSPSetType(mKspSolver, mKspType.c_str());
913 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
914#endif
915 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
917 {
918 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
919 }
920 else
921 {
922 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
923 }
924
925 delete[] r_eig;
926 delete[] c_eig;
927
928#ifdef TRACE_KSP
930 {
931 Timer::Print("Computing extremal eigenvalues");
932 }
933#endif
934 }
935
936#ifdef TRACE_KSP
937 Timer::Reset();
938#endif
939
940 KSPSetUp(mKspSolver);
941
942 if (chebyshev_lhs_vector)
943 {
944 PetscTools::Destroy(chebyshev_lhs_vector);
945 }
946
947#ifdef TRACE_KSP
949 {
950 Timer::Print("KSPSetUP (contains preconditioner creation for PETSc preconditioners)");
951 }
952#endif
953
954 mKspIsSetup = true;
955
956 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
957 }
958 else
959 {
960 if (mNonZerosUsed!=mat_info.nz_used)
961 {
962 WARNING("LinearSystem doesn't like the non-zero pattern of a matrix to change. (I think you changed it).");
963 mNonZerosUsed = mat_info.nz_used;
964 }
965// PetscScalar norm;
966// MatNorm(mLhsMatrix, NORM_FROBENIUS, &norm);
967// if (fabs(norm - mMatrixNorm) > 0)
968// {
969// EXCEPTION("LinearSystem doesn't allow the matrix norm to change");
970// }
971 }
972
973 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
974 // Create solution vector
976 Vec lhs_vector;
977 VecDuplicate(mRhsVector, &lhs_vector); // Sets the same size (doesn't copy)
978 if (lhsGuess)
979 {
980 VecCopy(lhsGuess, lhs_vector);
981 // If this wasn't done at construction time then it may be too late for this:
982 // KSPSetInitialGuessNonzero(mKspSolver, PETSC_TRUE);
983 // Is it possible to warn the user?
984 }
985
986 // Check if the right hand side is small (but non-zero), PETSc can diverge immediately
987 // with a non-zero initial guess. Here we check for this and alter the initial guess to zero.
988 PetscReal rhs_norm;
989 VecNorm(mRhsVector, NORM_1, &rhs_norm);
990 double rhs_zero_tol = 1e-15;
991
992 if (rhs_norm < rhs_zero_tol && rhs_norm > 0.0)
993 {
994 WARNING("Using zero initial guess due to small right hand side vector");
995 PetscVecTools::Zero(lhs_vector);
996 }
997
998 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
999// // Double check that the mRhsVector contains sensible values
1000// double* p_rhs,* p_guess;
1001// VecGetArray(mRhsVector, &p_rhs);
1002// VecGetArray(lhsGuess, &p_guess);
1003// for (int global_index=mOwnershipRangeLo; global_index<mOwnershipRangeHi; global_index++)
1004// {
1005// int local_index = global_index - mOwnershipRangeLo;
1006// assert(!std::isnan(p_rhs[local_index]));
1007// assert(!std::isnan(p_guess[local_index]));
1008// if (p_rhs[local_index] != p_rhs[local_index])
1009// {
1010// std::cout << "********* PETSc nan\n";
1011// assert(0);
1012// }
1013// }
1014// std::cout << "b[0] = " << p_rhs[0] << ", guess[0] = " << p_guess[0] << "\n";
1015// VecRestoreArray(mRhsVector, &p_rhs);
1016// VecRestoreArray(lhsGuess, &p_guess);
1017// // Test A*guess
1018// Vec temp;
1019// VecDuplicate(mRhsVector, &temp);
1020// MatMult(mLhsMatrix, lhs_vector, temp);
1021// double* p_temp;
1022// VecGetArray(temp, &p_temp);
1023// std::cout << "temp[0] = " << p_temp[0] << "\n";
1024// VecRestoreArray(temp, &p_temp);
1025// PetscTools::Destroy(temp);
1026// // Dump the matrix to file
1027// PetscViewer viewer;
1028// PetscViewerASCIIOpen(PETSC_COMM_WORLD,"mat.output",&viewer);
1029// MatView(mLhsMatrix, viewer);
1030// PetscViewerFlush(viewer);
1031// PetscViewerDestroy(PETSC_DESTROY_PARAM(viewer);
1032// // Dump the rhs vector to file
1033// PetscViewerASCIIOpen(PETSC_COMM_WORLD,"vec.output",&viewer);
1034// PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1035// VecView(mRhsVector, viewer);
1036// PetscViewerFlush(viewer);
1037// PetscViewerDestroy(PETSC_DESTROY_PARAM(viewer);
1038
1039 try
1040 {
1041 HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
1042
1043#ifdef TRACE_KSP
1044 Timer::Reset();
1045#endif
1046
1047 // Current solve has to be done with tolerance-based stop criteria in order to record iterations taken
1049 {
1050#if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
1051 KSPSetNormType(mKspSolver, KSP_PRECONDITIONED_NORM);
1052#else
1053 KSPSetNormType(mKspSolver, KSP_NORM_PRECONDITIONED);
1054#endif
1055
1056#if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
1058 {
1059 #if (PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
1060 KSPConvergedDefaultCreate(&mpConvergenceTestContext);
1061 #else
1062 KSPDefaultConvergedCreate(&mpConvergenceTestContext);
1063 #endif
1064 }
1065 #if (PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
1066 // From PETSc 3.5, KSPDefaultConverged became KSPConvergedDefault.
1067 KSPSetConvergenceTest(mKspSolver, KSPConvergedDefault, &mpConvergenceTestContext, CHASTE_PETSC_NULLPTR);
1068 #else
1069 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, &mpConvergenceTestContext, CHASTE_PETSC_NULLPTR);
1070 #endif
1071
1072#else
1073 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, CHASTE_PETSC_NULLPTR);
1074#endif
1075
1077 {
1078 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
1079 }
1080 else
1081 {
1082 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
1083 }
1084
1086 std::stringstream num_it_str;
1087 num_it_str << 1000;
1088 PetscTools::SetOption("-ksp_max_it", num_it_str.str().c_str());
1089
1090 // Adaptive Chebyshev: reevaluate spectrum with cg
1091 if (mKspType == "chebychev")
1092 {
1093 // You can estimate preconditioned matrix spectrum with CG
1094 KSPSetType(mKspSolver,"cg");
1095 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
1096 }
1097
1098 KSPSetFromOptions(mKspSolver);
1099 KSPSetUp(mKspSolver);
1100 }
1101
1102 PETSCEXCEPT(KSPSolve(mKspSolver, mRhsVector, lhs_vector));
1103 HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
1104
1105#ifdef TRACE_KSP
1106 PetscInt num_it;
1107 KSPGetIterationNumber(mKspSolver, &num_it);
1109 {
1110 std::cout << "++ Solve: " << mNumSolves << " NumIterations: " << num_it << " "; // don't add std::endl so we get Timer::Print output in the same line (better for grep-ing)
1111 Timer::Print("Solve");
1112 }
1113
1114 mTotalNumIterations += num_it;
1115 if ((unsigned) num_it > mMaxNumIterations)
1116 {
1117 mMaxNumIterations = num_it;
1118 }
1119#endif
1120
1121 // Check that solver converged and throw if not
1122 KSPConvergedReason reason;
1123 KSPGetConvergedReason(mKspSolver, &reason);
1124
1125 if (mUseFixedNumberIterations && PETSC_VERSION_MAJOR < 3)
1126 {
1127 WARNING("Not explicitly checking convergence reason when using fixed number of iterations and PETSc 2");
1128 }
1129 else
1130 {
1131 KSPEXCEPT(reason);
1132 }
1133
1135 {
1136 // Adaptive Chebyshev: reevaluate spectrum with cg
1137 if (mKspType == "chebychev")
1138 {
1139 PetscReal *r_eig = new PetscReal[mSize];
1140 PetscReal *c_eig = new PetscReal[mSize];
1141 PetscInt eigs_computed;
1142 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
1143
1144 mEigMin = r_eig[0];
1145 /*
1146 * Using max() covers a borderline case found in TestChasteBenchmarksForPreDiCT where there's a big
1147 * gap in the spectrum between ~1.2 and ~2.5. Some reevaluations pick up 2.5 and others don't. If it
1148 * is not picked up, Chebyshev will diverge after 10 solves or so.
1149 */
1150 mEigMax = std::max(mEigMax,r_eig[eigs_computed-1]);
1151
1152 delete[] r_eig;
1153 delete[] c_eig;
1154
1155 assert(mKspType == "chebychev");
1156#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
1157 KSPSetType(mKspSolver, "chebyshev");
1158 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
1159#else
1160 KSPSetType(mKspSolver, mKspType.c_str());
1161 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
1162#endif
1163 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
1164 }
1165
1166#if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
1167 if (mKspType == "chebychev")
1168 {
1169 // See #1695 for more details.
1170 EXCEPTION("Chebyshev with fixed number of iterations is known to be broken in PETSc <= 2.3.2");
1171 }
1172
1173 KSPSetNormType(mKspSolver, KSP_NO_NORM);
1174#elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
1175 KSPSetNormType(mKspSolver, KSP_NORM_NONE);
1176 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 7) //PETSc 3.7 or later
1177 /*
1178 * Up to PETSc 3.7.2 the above call also turned off the default convergence test.
1179 * However, in PETSc 3.7.3 (subminor release) this behaviour was removed and so, here,
1180 * we explicitly add it back again.
1181 * See
1182 * https://bitbucket.org/petsc/petsc/commits/eb70c44be3430b039effa3de7e1ca2fab9f75a57
1183 * This following line of code is actually valid from PETSc 3.5.
1184 */
1185 KSPSetConvergenceTest(mKspSolver, KSPConvergedSkip, CHASTE_PETSC_NULLPTR, CHASTE_PETSC_NULLPTR);
1186 #endif
1187#else
1188 KSPSetNormType(mKspSolver, KSP_NORM_NO);
1189#endif
1190
1191#if (PETSC_VERSION_MAJOR == 2)
1192 KSPSetConvergenceTest(mKspSolver, KSPSkipConverged, CHASTE_PETSC_NULLPTR);
1193#endif
1194
1195 PetscInt num_it;
1196 KSPGetIterationNumber(mKspSolver, &num_it);
1197 std::stringstream num_it_str;
1198 num_it_str << num_it;
1199 PetscTools::SetOption("-ksp_max_it", num_it_str.str().c_str());
1200
1201 KSPSetFromOptions(mKspSolver);
1202 KSPSetUp(mKspSolver);
1203
1205 }
1206
1207 mNumSolves++;
1208
1209 }
1210 catch (const Exception& e)
1211 {
1212 // Destroy solution vector on error to avoid memory leaks
1213 PetscTools::Destroy(lhs_vector);
1214 throw e;
1215 }
1216
1217 return lhs_vector;
1218}
1219
1221{
1222 mPrecondMatrixIsNotLhs = precondIsDifferent;
1223
1225 {
1226 if (mRowPreallocation == UINT_MAX)
1227 {
1228 /*
1229 * At the time of writing, this line will be reached if the constructor
1230 * with signature LinearSystem(Vec residualVector, Mat jacobianMatrix) is
1231 * called with jacobianMatrix=NULL and preconditioning matrix different
1232 * from lhs is used.
1233 *
1234 * If this combination is ever required you will need to work out
1235 * matrix allocation (mRowPreallocation) here.
1236 */
1238 }
1239
1242 }
1243}
1244
1245void LinearSystem::SetUseFixedNumberIterations(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
1246{
1247
1248 mUseFixedNumberIterations = useFixedNumberIterations;
1249 mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
1250}
1251
1253{
1254 if (mKspIsSetup)
1255 {
1256 KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
1257 }
1258
1259 mKspIsSetup = false;
1261
1262 /*
1263 * Reset max number of iterations. This option is stored in the configuration database and
1264 * explicitely read in with KSPSetFromOptions() everytime a KSP object is created. Therefore,
1265 * destroying the KSP object will not ensure that it is set back to default.
1266 */
1268 std::stringstream num_it_str;
1269 num_it_str << 1000;
1270 PetscTools::SetOption("-ksp_max_it", num_it_str.str().c_str());
1271}
1272
1273// Serialization for Boost >= 1.36
#define TERMINATE(message)
#define EXCEPTION(message)
#define NEVER_REACHED
#define CHASTE_PETSC_NULLPTR
A macro to define PETSc null pointer based on the PETSc version.
#define PETSC_DESTROY_PARAM(x)
#define CHASTE_CLASS_EXPORT(T)
Vec & rGetDirichletBoundaryConditionsVector()
void SetNullBasis(Vec nullbasis[], unsigned numberOfBases)
void SetPrecondMatrixIsDifferentFromLhs(bool precondIsDifferent=true)
boost::shared_ptr< std::vector< PetscInt > > mpBathNodes
Vec Solve(Vec lhsGuess=nullptr)
void FinaliseLhsMatrix()
Vec mDirichletBoundaryConditionsVector
PetscReal mEigMin
void SetMatrixIsSymmetric(bool isSymmetric=true)
void AssembleFinalLinearSystem()
Mat GetLhsMatrix() const
bool mUseAbsoluteTolerance
void ZeroMatrixColumn(PetscInt col)
bool mForceSpectrumReevaluation
void SetMatrixIsConstant(bool matrixIsConstant)
unsigned mRowPreallocation
void FinalisePrecondMatrix()
void GetOwnershipRange(PetscInt &lo, PetscInt &hi)
Mat & rGetLhsMatrix()
PCLDUFactorisation * mpLDUFactorisationPC
void SwitchWriteModeLhsMatrix()
void AddToMatrixElement(PetscInt row, PetscInt col, double value)
double GetMatrixElement(PetscInt row, PetscInt col)
PetscInt mOwnershipRangeLo
void SetKspType(const char *kspType)
void SetMatrixElement(PetscInt row, PetscInt col, double value)
unsigned GetNumIterations() const
std::string mKspType
bool IsMatrixSymmetric()
PetscInt mOwnershipRangeHi
Vec GetRhsVector() const
void SetUseFixedNumberIterations(bool useFixedNumberIterations=true, unsigned evaluateNumItsEveryNSolves=UINT_MAX)
void FinaliseRhsVector()
std::string mPcType
MatNullSpace mMatNullSpace
void ZeroMatrixRowsWithValueOnDiagonal(std::vector< unsigned > &rRows, double diagonalValue)
unsigned GetSize() const
void ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector< unsigned > &rRowColIndices, double diagonalValue)
void SetAbsoluteTolerance(double absoluteTolerance)
void AssembleIntermediateLinearSystem()
void SetRhsVectorElement(PetscInt row, double value)
void SetMatrixRow(PetscInt row, double value)
void * mpConvergenceTestContext
PCTwoLevelsBlockDiagonal * mpTwoLevelsBlockDiagonalPC
Mat & rGetPrecondMatrix()
void ZeroLinearSystem()
Vec GetMatrixRowDistributed(unsigned rowIndex)
PetscReal mEigMax
PCBlockDiagonal * mpBlockDiagonalPC
void AddToRhsVectorElement(PetscInt row, double value)
unsigned mEvaluateNumItsEveryNSolves
Vec & rGetRhsVector()
bool mPrecondMatrixIsNotLhs
void SetRelativeTolerance(double relativeTolerance)
void RemoveNullSpace()
bool mUseFixedNumberIterations
double GetRhsVectorElement(PetscInt row)
void SetPcType(const char *pcType, boost::shared_ptr< std::vector< PetscInt > > pBathNodes=boost::shared_ptr< std::vector< PetscInt > >())
LinearSystem(PetscInt lhsVectorSize, unsigned rowPreallocation=UINT_MAX)
unsigned mNumSolves
double mNonZerosUsed
static double GetElement(Mat matrix, PetscInt row, PetscInt col)
static void AddToElement(Mat matrix, PetscInt row, PetscInt col, double value)
static void Zero(Mat matrix)
static void SetRow(Mat matrix, PetscInt row, double value)
static void ZeroRowsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > &rRows, double diagonalValue)
static void SetOption(Mat matrix, MatOption option)
static void ZeroRowsAndColumnsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > rowColIndices, double diagonalValue)
static void SwitchWriteMode(Mat matrix)
static Vec GetMatrixRowDistributed(Mat matrix, unsigned rowIndex)
static void Finalise(Mat matrix)
static void ZeroColumn(Mat matrix, PetscInt col)
static void SetElement(Mat matrix, PetscInt row, PetscInt col, double value)
static void Display(Mat matrix)
static void Destroy(Vec &rVec)
static bool AmMaster()
static Vec CreateVec(int size, int localSize=PETSC_DECIDE, bool ignoreOffProcEntries=true)
static void SetOption(const char *pOptionName, const char *pOptionValue)
static void SetupMat(Mat &rMat, int numRows, int numColumns, unsigned rowPreallocation, int numLocalRows=PETSC_DECIDE, int numLocalColumns=PETSC_DECIDE, bool ignoreOffProcEntries=true, bool newAllocationError=true)
static double GetElement(Vec vector, PetscInt row)
static void AddToElement(Vec vector, PetscInt row, double value)
static void Finalise(Vec vector)
static void Display(Vec vector)
static void SetElement(Vec vector, PetscInt row, double value)
static void Zero(Vec vector)
static void Print(std::string message)
Definition Timer.cpp:59
static void Reset()
Definition Timer.cpp:44