Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
LinearSystem.cpp
1/*
2
3Copyright (c) 2005-2024, 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 >= 5) //PETSc 3.5 or later
288 KSPConvergedDefaultDestroy(mpConvergenceTestContext);
289#else
290 KSPDefaultConvergedDestroy(mpConvergenceTestContext);
291#endif
292 }
293#endif
294
295#ifdef TRACE_KSP
297 {
298 if (mNumSolves > 0)
299 {
300 double ave_num_iterations = mTotalNumIterations/(double)mNumSolves;
301
302 std::cout << std::endl << "KSP iterations report:" << std::endl;
303 std::cout << "mNumSolves" << "\t" << "mTotalNumIterations" << "\t" << "mMaxNumIterations" << "\t" << "mAveNumIterations" << std::endl;
304 std::cout << mNumSolves << "\t" << mTotalNumIterations << "\t" << mMaxNumIterations << "\t" << ave_num_iterations << std::endl;
305 }
306 }
307#endif
308}
309
311{
312 PetscMatTools::SetElement(mLhsMatrix, row, col, value);
313}
314
316{
317 PetscMatTools::AddToElement(mLhsMatrix, row, col, value);
318}
319
325
331
336
341
346
351
353{
355}
356
358{
360}
361
366
371
373{
375}
376
381
382void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal(std::vector<unsigned>& rRows, double diagonalValue)
383{
385}
386
387void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector<unsigned>& rRowColIndices, double diagonalValue)
388{
390}
391
396
401
406
412
413unsigned LinearSystem::GetSize() const
414{
415 return (unsigned) mSize;
416}
417
418void LinearSystem::SetNullBasis(Vec nullBasis[], unsigned numberOfBases)
419{
420#ifndef NDEBUG
421 // Check all the vectors of the base are normal
422 for (unsigned vec_index=0; vec_index<numberOfBases; vec_index++)
423 {
424 PetscReal l2_norm;
425 VecNorm(nullBasis[vec_index], NORM_2, &l2_norm);
426 if (fabs(l2_norm-1.0) > 1e-08)
427 {
428 EXCEPTION("One of the vectors in the null space is not normalised");
429 }
430 }
431
432 // Check all the vectors of the base are orthogonal
433 for (unsigned vec_index=1; vec_index<numberOfBases; vec_index++)
434 {
435 // 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.
436 unsigned num_vectors_ahead = numberOfBases-vec_index;
437 boost::scoped_array<PetscScalar> dot_products(new PetscScalar[num_vectors_ahead]);
438#if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
439 VecMDot(num_vectors_ahead, nullBasis[vec_index-1], &nullBasis[vec_index], dot_products.get());
440#else
441 VecMDot(nullBasis[vec_index-1], num_vectors_ahead, &nullBasis[vec_index], dot_products.get());
442#endif
443 for (unsigned index=0; index<num_vectors_ahead; index++)
444 {
445 if (fabs(dot_products[index]) > 1e-08 )
446 {
447 EXCEPTION("The null space is not orthogonal.");
448 }
449 }
450 }
451
452#endif
453
454 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, numberOfBases, nullBasis, &mMatNullSpace) );
455}
456
458{
459 // Only remove if previously set
460 if (mMatNullSpace)
461 {
462 PETSCEXCEPT( MatNullSpaceDestroy(PETSC_DESTROY_PARAM(mMatNullSpace)) );
463 PETSCEXCEPT( MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_FALSE, 0, nullptr, &mMatNullSpace) );
464#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
465 // Setting null space in the KSP was deprecated in PETSc 3.6, but setting the null space
466 // for the matrix appeared in PETSc 3.3 so 3.3, 3.4, 3.5 can do either
467 PETSCEXCEPT( MatSetNullSpace(mLhsMatrix, mMatNullSpace) );
468 /* This is heavy-handed:
469 * Changing the Null-space appears to only work if we reset the KSP
470 * Adding null-space to the matrix has to happen *before* KSPSetOperators
471 */
473#else
474 if (mKspIsSetup)
475 {
476 PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
477 }
478 //else: it will be set next time Solve() is called
479#endif
480 }
481}
482
488
493
498
500{
501 assert(this->mKspIsSetup);
502
503 PetscInt num_its;
504 KSPGetIterationNumber(this->mKspSolver, &num_its);
505
506 return (unsigned) num_its;
507}
508
513
515{
516 return mRhsVector;
517}
518
523
525{
526 return mLhsMatrix;
527}
528
530{
532 {
533 EXCEPTION("LHS matrix used for preconditioner construction");
534 }
535
536 return mPrecondMatrix;
537}
538
543
545{
547
548 if (isSymmetric)
549 {
550 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRIC);
551 PetscMatTools::SetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL);
552 }
553 else
554 {
555// don't have a PetscMatTools method for setting options to false
556#if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
557 MatSetOption(mLhsMatrix, MAT_SYMMETRIC, PETSC_FALSE);
558 MatSetOption(mLhsMatrix, MAT_STRUCTURALLY_SYMMETRIC, PETSC_FALSE);
559 MatSetOption(mLhsMatrix, MAT_SYMMETRY_ETERNAL, PETSC_FALSE);
560#else
561 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRIC);
562 MatSetOption(mLhsMatrix, MAT_NOT_STRUCTURALLY_SYMMETRIC);
563 MatSetOption(mLhsMatrix, MAT_NOT_SYMMETRY_ETERNAL);
564#endif
565 }
566}
567
569{
570 PetscBool symmetry_flag_is_set;
571 PetscBool symmetry_flag;
572
573 MatIsSymmetricKnown(mLhsMatrix, &symmetry_flag_is_set, &symmetry_flag);
574
575 // If the flag is not set we assume is a non-symmetric matrix
576 return symmetry_flag_is_set && symmetry_flag;
577}
578
579void LinearSystem::SetMatrixIsConstant(bool matrixIsConstant)
580{
581 mMatrixIsConstant = matrixIsConstant;
582}
583
584void LinearSystem::SetRelativeTolerance(double relativeTolerance)
585{
586 mTolerance = relativeTolerance;
587 mUseAbsoluteTolerance = false;
588 if (mKspIsSetup)
589 {
590 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
591 }
592}
593
594void LinearSystem::SetAbsoluteTolerance(double absoluteTolerance)
595{
596 mTolerance = absoluteTolerance;
598 if (mKspIsSetup)
599 {
600 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
601 }
602}
603
604void LinearSystem::SetKspType(const char *kspType)
605{
606 mKspType = kspType;
607 if (mKspIsSetup)
608 {
609 KSPSetType(mKspSolver, kspType);
610 KSPSetFromOptions(mKspSolver);
611 }
612}
613
614void LinearSystem::SetPcType(const char* pcType, boost::shared_ptr<std::vector<PetscInt> > pBathNodes)
615{
616 mPcType = pcType;
617 mpBathNodes = pBathNodes;
618
619 if (mKspIsSetup)
620 {
621 if (mPcType == "blockdiagonal")
622 {
623 // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
625 delete mpBlockDiagonalPC;
626 mpBlockDiagonalPC = nullptr;
628 mpLDUFactorisationPC = nullptr;
631
633 }
634 else if (mPcType == "ldufactorisation")
635 {
636 // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
638 delete mpBlockDiagonalPC;
639 mpBlockDiagonalPC = nullptr;
641 mpLDUFactorisationPC = nullptr;
644
646 }
647 else if (mPcType == "twolevelsblockdiagonal")
648 {
649 // If the previous preconditioner was purpose-built we need to free the appropriate pointer.
651 delete mpBlockDiagonalPC;
652 mpBlockDiagonalPC = nullptr;
654 mpLDUFactorisationPC = nullptr;
657
658 if (!mpBathNodes)
659 {
660 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC"); // LCOV_EXCL_LINE
661 }
663 }
664 else
665 {
666 PC prec;
667 KSPGetPC(mKspSolver, &prec);
668 PCSetType(prec, pcType);
669 }
670 KSPSetFromOptions(mKspSolver);
671 }
672}
673
675{
676 /*
677 * The following lines are very useful for debugging:
678 * MatView(mLhsMatrix, PETSC_VIEWER_STDOUT_WORLD);
679 * VecView(mRhsVector, PETSC_VIEWER_STDOUT_WORLD);
680 */
681
682 // Double check that the non-zero pattern hasn't changed
683 MatInfo mat_info;
684 MatGetInfo(mLhsMatrix, MAT_GLOBAL_SUM, &mat_info);
685
686 if (!mKspIsSetup)
687 {
688 // Create PETSc Vec that may be required if we use a Chebyshev solver
689 Vec chebyshev_lhs_vector = nullptr;
690
691 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
692 mNonZerosUsed=mat_info.nz_used;
693 //MatNorm(mLhsMatrix, NORM_FROBENIUS, &mMatrixNorm);
694 PC prec; //Type of pre-conditioner
695
696 KSPCreate(PETSC_COMM_WORLD, &mKspSolver);
697
698 if (mMatNullSpace) // Adding null-space to the matrix (new style) has to happen *before* KSPSetOperators
699 {
700#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
701 // Setting null space in the KSP was deprecated in PETSc 3.6, but setting the null space
702 // for the matrix appeared in PETSc 3.3 so 3.3, 3.4, 3.5 can do either
703
704 PETSCEXCEPT(MatSetNullSpace(mLhsMatrix, mMatNullSpace));
705#else
706 PETSCEXCEPT(KSPSetNullSpace(mKspSolver, mMatNullSpace));
707#endif
708 }
709#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
710 // Do nothing. Note that reusing the pre-conditioner in later PETSc versions is done after the pre-conditioner is formed
711 // (This comment is retained here so that the #if logic is consistent.)
712#else
713 /*
714 * The preconditioner flag (last argument) in the following calls says
715 * how to reuse the preconditioner on subsequent iterations.
716 */
717 MatStructure preconditioner_over_successive_calls;
718
720 {
721 preconditioner_over_successive_calls = SAME_PRECONDITIONER;
722 }
723 else
724 {
725 preconditioner_over_successive_calls = SAME_NONZERO_PATTERN;
726 }
727#endif
728
730 {
731#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
732 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix);
733#else
734 KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix, preconditioner_over_successive_calls);
735#endif
736 }
737 else
738 {
739#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
740 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix);
741#else
742 KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix, preconditioner_over_successive_calls);
743#endif
744 }
745#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
747 {
748 // Emulate SAME_PRECONDITIONER as above
749 KSPSetReusePreconditioner(mKspSolver, PETSC_TRUE);
750 }
751#endif
752
753 // Set either absolute or relative tolerance of the KSP solver.
754 // The default is to use relative tolerance (1e-6)
756 {
757 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
758 }
759 else
760 {
761 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
762 }
763
764 // Set ksp and pc types
765#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
766 if (mKspType == "chebychev")
767 {
768 KSPSetType(mKspSolver, "chebyshev");
769 }
770 else
771 {
772 KSPSetType(mKspSolver, mKspType.c_str());
773 }
774#else
775 KSPSetType(mKspSolver, mKspType.c_str());
776#endif
777 KSPGetPC(mKspSolver, &prec);
778
779 // Turn off pre-conditioning if the system size is very small
780 const bool is_small = (mSize <= 6);
781 if (is_small)
782 {
783 PCSetType(prec, PCNONE);
784 }
785 else
786 {
787#ifdef TRACE_KSP
788 Timer::Reset();
789#endif
790 if (mPcType == "blockdiagonal")
791 {
793#ifdef TRACE_KSP
795 {
796 Timer::Print("Purpose-build preconditioner creation");
797 }
798#endif
799 }
800 else if (mPcType == "ldufactorisation")
801 {
803#ifdef TRACE_KSP
805 {
806 Timer::Print("Purpose-build preconditioner creation");
807 }
808#endif
809 }
810 else if (mPcType == "twolevelsblockdiagonal")
811 {
812 if (!mpBathNodes)
813 {
814 TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC"); // LCOV_EXCL_LINE
815 }
817#ifdef TRACE_KSP
819 {
820 Timer::Print("Purpose-build preconditioner creation");
821 }
822#endif
823
824 }
825 else
826 {
827 PCSetType(prec, mPcType.c_str());
828 }
829 }
830
831 KSPSetFromOptions(mKspSolver);
832
833 if (lhsGuess)
834 {
835 // Assume that the user of this method will always be kind enough to give us a reasonable guess.
836 KSPSetInitialGuessNonzero(mKspSolver,PETSC_TRUE);
837 }
838 /*
839 * Non-adaptive Chebyshev: the required spectrum approximation is computed just once
840 * at the beginning of the simulation. This is done with two extra CG solves.
841 */
842 if (mKspType == "chebychev" && !mUseFixedNumberIterations)
843 {
844#ifdef TRACE_KSP
845 Timer::Reset();
846#endif
847 // Preconditioned matrix spectrum is approximated with CG
848 KSPSetType(mKspSolver,"cg");
849 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
850 KSPSetUp(mKspSolver);
851
852 VecDuplicate(mRhsVector, &chebyshev_lhs_vector);
853 if (lhsGuess)
854 {
855 VecCopy(lhsGuess, chebyshev_lhs_vector);
856 }
857
858 // Smallest eigenvalue is approximated to default tolerance
859 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
860
861 PetscReal *r_eig = new PetscReal[mSize];
862 PetscReal *c_eig = new PetscReal[mSize];
863 PetscInt eigs_computed;
864 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
865
866 mEigMin = r_eig[0];
867
868 // Largest eigenvalue is approximated to machine precision
869 KSPSetTolerances(mKspSolver, DBL_EPSILON, DBL_EPSILON, PETSC_DEFAULT, PETSC_DEFAULT);
870 KSPSetUp(mKspSolver);
871 if (lhsGuess)
872 {
873 VecCopy(lhsGuess, chebyshev_lhs_vector);
874 }
875
876 KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
877 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
878
879 mEigMax = r_eig[eigs_computed-1];
880
881#ifdef TRACE_KSP
882 /*
883 * Under certain circumstances (see Golub&Overton 1988), underestimating
884 * the spectrum of the preconditioned operator improves convergence rate.
885 * See publication for a discussion and for definition of alpha and sigma_one.
886 */
888 {
889 std::cout << "EIGS: ";
890 for (int index=0; index<eigs_computed; index++)
891 {
892 std::cout << r_eig[index] << ", ";
893 }
894 std::cout << std::endl;
895 }
896
897 if (PetscTools::AmMaster()) std::cout << "EIGS "<< mEigMax << " " << mEigMin <<std::endl;
898 double alpha = 2/(mEigMax+mEigMin);
899 double sigma_one = 1 - alpha*mEigMin;
900 if (PetscTools::AmMaster()) std::cout << "sigma_1 = 1 - alpha*mEigMin = "<< sigma_one <<std::endl;
901#endif
902
903 // Set Chebyshev solver and max/min eigenvalues
904 assert(mKspType == "chebychev");
905#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
906 KSPSetType(mKspSolver, "chebyshev");
907 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
908#else
909 KSPSetType(mKspSolver, mKspType.c_str());
910 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
911#endif
912 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
914 {
915 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
916 }
917 else
918 {
919 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
920 }
921
922 delete[] r_eig;
923 delete[] c_eig;
924
925#ifdef TRACE_KSP
927 {
928 Timer::Print("Computing extremal eigenvalues");
929 }
930#endif
931 }
932
933#ifdef TRACE_KSP
934 Timer::Reset();
935#endif
936
937 KSPSetUp(mKspSolver);
938
939 if (chebyshev_lhs_vector)
940 {
941 PetscTools::Destroy(chebyshev_lhs_vector);
942 }
943
944#ifdef TRACE_KSP
946 {
947 Timer::Print("KSPSetUP (contains preconditioner creation for PETSc preconditioners)");
948 }
949#endif
950
951 mKspIsSetup = true;
952
953 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
954 }
955 else
956 {
957 if (mNonZerosUsed!=mat_info.nz_used)
958 {
959 WARNING("LinearSystem doesn't like the non-zero pattern of a matrix to change. (I think you changed it).");
960 mNonZerosUsed = mat_info.nz_used;
961 }
962// PetscScalar norm;
963// MatNorm(mLhsMatrix, NORM_FROBENIUS, &norm);
964// if (fabs(norm - mMatrixNorm) > 0)
965// {
966// EXCEPTION("LinearSystem doesn't allow the matrix norm to change");
967// }
968 }
969
970 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
971 // Create solution vector
973 Vec lhs_vector;
974 VecDuplicate(mRhsVector, &lhs_vector); // Sets the same size (doesn't copy)
975 if (lhsGuess)
976 {
977 VecCopy(lhsGuess, lhs_vector);
978 // If this wasn't done at construction time then it may be too late for this:
979 // KSPSetInitialGuessNonzero(mKspSolver, PETSC_TRUE);
980 // Is it possible to warn the user?
981 }
982
983 // Check if the right hand side is small (but non-zero), PETSc can diverge immediately
984 // with a non-zero initial guess. Here we check for this and alter the initial guess to zero.
985 PetscReal rhs_norm;
986 VecNorm(mRhsVector, NORM_1, &rhs_norm);
987 double rhs_zero_tol = 1e-15;
988
989 if (rhs_norm < rhs_zero_tol && rhs_norm > 0.0)
990 {
991 WARNING("Using zero initial guess due to small right hand side vector");
992 PetscVecTools::Zero(lhs_vector);
993 }
994
995 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
996// // Double check that the mRhsVector contains sensible values
997// double* p_rhs,* p_guess;
998// VecGetArray(mRhsVector, &p_rhs);
999// VecGetArray(lhsGuess, &p_guess);
1000// for (int global_index=mOwnershipRangeLo; global_index<mOwnershipRangeHi; global_index++)
1001// {
1002// int local_index = global_index - mOwnershipRangeLo;
1003// assert(!std::isnan(p_rhs[local_index]));
1004// assert(!std::isnan(p_guess[local_index]));
1005// if (p_rhs[local_index] != p_rhs[local_index])
1006// {
1007// std::cout << "********* PETSc nan\n";
1008// assert(0);
1009// }
1010// }
1011// std::cout << "b[0] = " << p_rhs[0] << ", guess[0] = " << p_guess[0] << "\n";
1012// VecRestoreArray(mRhsVector, &p_rhs);
1013// VecRestoreArray(lhsGuess, &p_guess);
1014// // Test A*guess
1015// Vec temp;
1016// VecDuplicate(mRhsVector, &temp);
1017// MatMult(mLhsMatrix, lhs_vector, temp);
1018// double* p_temp;
1019// VecGetArray(temp, &p_temp);
1020// std::cout << "temp[0] = " << p_temp[0] << "\n";
1021// VecRestoreArray(temp, &p_temp);
1022// PetscTools::Destroy(temp);
1023// // Dump the matrix to file
1024// PetscViewer viewer;
1025// PetscViewerASCIIOpen(PETSC_COMM_WORLD,"mat.output",&viewer);
1026// MatView(mLhsMatrix, viewer);
1027// PetscViewerFlush(viewer);
1028// PetscViewerDestroy(PETSC_DESTROY_PARAM(viewer);
1029// // Dump the rhs vector to file
1030// PetscViewerASCIIOpen(PETSC_COMM_WORLD,"vec.output",&viewer);
1031// PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1032// VecView(mRhsVector, viewer);
1033// PetscViewerFlush(viewer);
1034// PetscViewerDestroy(PETSC_DESTROY_PARAM(viewer);
1035
1036 try
1037 {
1038 HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
1039
1040#ifdef TRACE_KSP
1041 Timer::Reset();
1042#endif
1043
1044 // Current solve has to be done with tolerance-based stop criteria in order to record iterations taken
1046 {
1047#if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
1048 KSPSetNormType(mKspSolver, KSP_PRECONDITIONED_NORM);
1049#else
1050 KSPSetNormType(mKspSolver, KSP_NORM_PRECONDITIONED);
1051#endif
1052
1053#if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
1055 {
1056 #if (PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
1057 KSPConvergedDefaultCreate(&mpConvergenceTestContext);
1058 #else
1059 KSPDefaultConvergedCreate(&mpConvergenceTestContext);
1060 #endif
1061 }
1062 #if (PETSC_VERSION_MINOR >= 5) //PETSc 3.5 or later
1063 // From PETSc 3.5, KSPDefaultConverged became KSPConvergedDefault.
1064 KSPSetConvergenceTest(mKspSolver, KSPConvergedDefault, &mpConvergenceTestContext, PETSC_NULL);
1065 #else
1066 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, &mpConvergenceTestContext, PETSC_NULL);
1067 #endif
1068
1069#else
1070 KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, PETSC_NULL);
1071#endif
1072
1074 {
1075 KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
1076 }
1077 else
1078 {
1079 KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
1080 }
1081
1083 std::stringstream num_it_str;
1084 num_it_str << 1000;
1085 PetscTools::SetOption("-ksp_max_it", num_it_str.str().c_str());
1086
1087 // Adaptive Chebyshev: reevaluate spectrum with cg
1088 if (mKspType == "chebychev")
1089 {
1090 // You can estimate preconditioned matrix spectrum with CG
1091 KSPSetType(mKspSolver,"cg");
1092 KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
1093 }
1094
1095 KSPSetFromOptions(mKspSolver);
1096 KSPSetUp(mKspSolver);
1097 }
1098
1099 PETSCEXCEPT(KSPSolve(mKspSolver, mRhsVector, lhs_vector));
1100 HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
1101
1102#ifdef TRACE_KSP
1103 PetscInt num_it;
1104 KSPGetIterationNumber(mKspSolver, &num_it);
1106 {
1107 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)
1108 Timer::Print("Solve");
1109 }
1110
1111 mTotalNumIterations += num_it;
1112 if ((unsigned) num_it > mMaxNumIterations)
1113 {
1114 mMaxNumIterations = num_it;
1115 }
1116#endif
1117
1118 // Check that solver converged and throw if not
1119 KSPConvergedReason reason;
1120 KSPGetConvergedReason(mKspSolver, &reason);
1121
1122 if (mUseFixedNumberIterations && PETSC_VERSION_MAJOR < 3)
1123 {
1124 WARNING("Not explicitly checking convergence reason when using fixed number of iterations and PETSc 2");
1125 }
1126 else
1127 {
1128 KSPEXCEPT(reason);
1129 }
1130
1132 {
1133 // Adaptive Chebyshev: reevaluate spectrum with cg
1134 if (mKspType == "chebychev")
1135 {
1136 PetscReal *r_eig = new PetscReal[mSize];
1137 PetscReal *c_eig = new PetscReal[mSize];
1138 PetscInt eigs_computed;
1139 KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
1140
1141 mEigMin = r_eig[0];
1142 /*
1143 * Using max() covers a borderline case found in TestChasteBenchmarksForPreDiCT where there's a big
1144 * gap in the spectrum between ~1.2 and ~2.5. Some reevaluations pick up 2.5 and others don't. If it
1145 * is not picked up, Chebyshev will diverge after 10 solves or so.
1146 */
1147 mEigMax = std::max(mEigMax,r_eig[eigs_computed-1]);
1148
1149 delete[] r_eig;
1150 delete[] c_eig;
1151
1152 assert(mKspType == "chebychev");
1153#if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
1154 KSPSetType(mKspSolver, "chebyshev");
1155 KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
1156#else
1157 KSPSetType(mKspSolver, mKspType.c_str());
1158 KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
1159#endif
1160 KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
1161 }
1162
1163#if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
1164 if (mKspType == "chebychev")
1165 {
1166 // See #1695 for more details.
1167 EXCEPTION("Chebyshev with fixed number of iterations is known to be broken in PETSc <= 2.3.2");
1168 }
1169
1170 KSPSetNormType(mKspSolver, KSP_NO_NORM);
1171#elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
1172 KSPSetNormType(mKspSolver, KSP_NORM_NONE);
1173 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 7) //PETSc 3.7 or later
1174 /*
1175 * Up to PETSc 3.7.2 the above call also turned off the default convergence test.
1176 * However, in PETSc 3.7.3 (subminor release) this behaviour was removed and so, here,
1177 * we explicitly add it back again.
1178 * See
1179 * https://bitbucket.org/petsc/petsc/commits/eb70c44be3430b039effa3de7e1ca2fab9f75a57
1180 * This following line of code is actually valid from PETSc 3.5.
1181 */
1182 KSPSetConvergenceTest(mKspSolver, KSPConvergedSkip, PETSC_NULL, PETSC_NULL);
1183 #endif
1184#else
1185 KSPSetNormType(mKspSolver, KSP_NORM_NO);
1186#endif
1187
1188#if (PETSC_VERSION_MAJOR == 2)
1189 KSPSetConvergenceTest(mKspSolver, KSPSkipConverged, PETSC_NULL);
1190#endif
1191
1192 PetscInt num_it;
1193 KSPGetIterationNumber(mKspSolver, &num_it);
1194 std::stringstream num_it_str;
1195 num_it_str << num_it;
1196 PetscTools::SetOption("-ksp_max_it", num_it_str.str().c_str());
1197
1198 KSPSetFromOptions(mKspSolver);
1199 KSPSetUp(mKspSolver);
1200
1202 }
1203
1204 mNumSolves++;
1205
1206 }
1207 catch (const Exception& e)
1208 {
1209 // Destroy solution vector on error to avoid memory leaks
1210 PetscTools::Destroy(lhs_vector);
1211 throw e;
1212 }
1213
1214 return lhs_vector;
1215}
1216
1218{
1219 mPrecondMatrixIsNotLhs = precondIsDifferent;
1220
1222 {
1223 if (mRowPreallocation == UINT_MAX)
1224 {
1225 /*
1226 * At the time of writing, this line will be reached if the constructor
1227 * with signature LinearSystem(Vec residualVector, Mat jacobianMatrix) is
1228 * called with jacobianMatrix=NULL and preconditioning matrix different
1229 * from lhs is used.
1230 *
1231 * If this combination is ever required you will need to work out
1232 * matrix allocation (mRowPreallocation) here.
1233 */
1235 }
1236
1239 }
1240}
1241
1242void LinearSystem::SetUseFixedNumberIterations(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
1243{
1244
1245 mUseFixedNumberIterations = useFixedNumberIterations;
1246 mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
1247}
1248
1250{
1251 if (mKspIsSetup)
1252 {
1253 KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
1254 }
1255
1256 mKspIsSetup = false;
1258
1259 /*
1260 * Reset max number of iterations. This option is stored in the configuration database and
1261 * explicitely read in with KSPSetFromOptions() everytime a KSP object is created. Therefore,
1262 * destroying the KSP object will not ensure that it is set back to default.
1263 */
1265 std::stringstream num_it_str;
1266 num_it_str << 1000;
1267 PetscTools::SetOption("-ksp_max_it", num_it_str.str().c_str());
1268}
1269
1270// Serialization for Boost >= 1.36
#define TERMINATE(message)
#define EXCEPTION(message)
#define NEVER_REACHED
#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