Chaste  Release::3.4
LinearSystem.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, 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 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF 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 
54 LinearSystem::LinearSystem(PetscInt lhsVectorSize, unsigned rowPreallocation)
55  :mPrecondMatrix(NULL),
56  mSize(lhsVectorSize),
57  mMatNullSpace(NULL),
58  mDestroyMatAndVec(true),
59  mKspIsSetup(false),
60  mNonZerosUsed(0.0),
61  mMatrixIsConstant(false),
62  mTolerance(1e-6),
63  mUseAbsoluteTolerance(false),
64  mDirichletBoundaryConditionsVector(NULL),
65  mpBlockDiagonalPC(NULL),
66  mpLDUFactorisationPC(NULL),
67  mpTwoLevelsBlockDiagonalPC(NULL),
68  mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
69  mPrecondMatrixIsNotLhs(false),
70  mRowPreallocation(rowPreallocation),
71  mUseFixedNumberIterations(false),
72  mEvaluateNumItsEveryNSolves(UINT_MAX),
73  mpConvergenceTestContext(NULL),
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 
109 LinearSystem::LinearSystem(PetscInt lhsVectorSize, Mat lhsMatrix, Vec rhsVector)
110  :mPrecondMatrix(NULL),
111  mSize(lhsVectorSize),
112  mMatNullSpace(NULL),
113  mDestroyMatAndVec(true),
114  mKspIsSetup(false),
115  mNonZerosUsed(0.0),
116  mMatrixIsConstant(false),
117  mTolerance(1e-6),
118  mUseAbsoluteTolerance(false),
119  mDirichletBoundaryConditionsVector(NULL),
120  mpBlockDiagonalPC(NULL),
121  mpLDUFactorisationPC(NULL),
122  mpTwoLevelsBlockDiagonalPC(NULL),
123  mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
124  mPrecondMatrixIsNotLhs(false),
125  mUseFixedNumberIterations(false),
126  mEvaluateNumItsEveryNSolves(UINT_MAX),
127  mpConvergenceTestContext(NULL),
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 
146 LinearSystem::LinearSystem(Vec templateVector, unsigned rowPreallocation, bool newAllocationError)
147  :mPrecondMatrix(NULL),
148  mMatNullSpace(NULL),
149  mDestroyMatAndVec(true),
150  mKspIsSetup(false),
151  mMatrixIsConstant(false),
152  mTolerance(1e-6),
153  mUseAbsoluteTolerance(false),
154  mDirichletBoundaryConditionsVector(NULL),
155  mpBlockDiagonalPC(NULL),
156  mpLDUFactorisationPC(NULL),
157  mpTwoLevelsBlockDiagonalPC(NULL),
158  mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
159  mPrecondMatrixIsNotLhs(false),
160  mRowPreallocation(rowPreallocation),
161  mUseFixedNumberIterations(false),
162  mEvaluateNumItsEveryNSolves(UINT_MAX),
163  mpConvergenceTestContext(NULL),
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 
187 LinearSystem::LinearSystem(Vec residualVector, Mat jacobianMatrix)
188  :mPrecondMatrix(NULL),
189  mMatNullSpace(NULL),
190  mDestroyMatAndVec(false),
191  mKspIsSetup(false),
192  mMatrixIsConstant(false),
193  mTolerance(1e-6),
194  mUseAbsoluteTolerance(false),
195  mDirichletBoundaryConditionsVector(NULL),
196  mpBlockDiagonalPC(NULL),
197  mpLDUFactorisationPC(NULL),
198  mpTwoLevelsBlockDiagonalPC(NULL),
199  mpBathNodes( boost::shared_ptr<std::vector<PetscInt> >() ),
200  mPrecondMatrixIsNotLhs(false),
201  mRowPreallocation(UINT_MAX),
202  mUseFixedNumberIterations(false),
203  mEvaluateNumItsEveryNSolves(UINT_MAX),
204  mpConvergenceTestContext(NULL),
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;
254  delete mpLDUFactorisationPC;
256 
257  if (mDestroyMatAndVec)
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
296  if (PetscTools::AmMaster())
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 
321 {
324 }
325 
327 {
330 }
331 
333 {
335 }
336 
338 {
340 }
341 
343 {
345 }
346 
348 {
350 }
351 
353 {
355 }
356 
358 {
360 }
361 
363 {
365 }
366 
368 {
370 }
371 
372 void LinearSystem::SetMatrixRow(PetscInt row, double value)
373 {
374  PetscMatTools::SetRow(mLhsMatrix, row, value);
375 }
376 
378 {
380 }
381 
382 void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal(std::vector<unsigned>& rRows, double diagonalValue)
383 {
385 }
386 
387 void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector<unsigned>& rRowColIndices, double diagonalValue)
388 {
390 }
391 
393 {
395 }
396 
398 {
400 }
401 
403 {
405 }
406 
408 {
409  ZeroRhsVector();
410  ZeroLhsMatrix();
411 }
412 
413 unsigned LinearSystem::GetSize() const
414 {
415  return (unsigned) mSize;
416 }
417 
418 void 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, NULL, &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  */
472  ResetKspSolver();
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 
484 {
485  lo = mOwnershipRangeLo;
486  hi = mOwnershipRangeHi;
487 }
488 
490 {
491  return PetscMatTools::GetElement(mLhsMatrix, row, col);
492 }
493 
495 {
497 }
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 
510 {
511  return mRhsVector;
512 }
513 
515 {
516  return mRhsVector;
517 }
518 
520 {
521  return mLhsMatrix;
522 }
523 
525 {
526  return mLhsMatrix;
527 }
528 
530 {
532  {
533  EXCEPTION("LHS matrix used for preconditioner construction");
534  }
535 
536  return mPrecondMatrix;
537 }
538 
540 {
542 }
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 
579 void LinearSystem::SetMatrixIsConstant(bool matrixIsConstant)
580 {
581  mMatrixIsConstant = matrixIsConstant;
582 }
583 
584 void 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 
594 void LinearSystem::SetAbsoluteTolerance(double absoluteTolerance)
595 {
596  mTolerance = absoluteTolerance;
597  mUseAbsoluteTolerance = true;
598  if (mKspIsSetup)
599  {
600  KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
601  }
602 }
603 
604 void LinearSystem::SetKspType(const char *kspType)
605 {
606  mKspType = kspType;
607  if (mKspIsSetup)
608  {
609  KSPSetType(mKspSolver, kspType);
610  KSPSetFromOptions(mKspSolver);
611  }
612 }
613 
614 void 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 = NULL;
627  delete mpLDUFactorisationPC;
628  mpLDUFactorisationPC = NULL;
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 = NULL;
640  delete mpLDUFactorisationPC;
641  mpLDUFactorisationPC = NULL;
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 = NULL;
653  delete mpLDUFactorisationPC;
654  mpLDUFactorisationPC = NULL;
657 
658  if (!mpBathNodes)
659  {
660  TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
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 = NULL;
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  const bool is_small = (mSize <= 6);
699 
700 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR >= 5 ) //PETSc 3.5 or later
701  if (mMatrixIsConstant && (!is_small))
702  {
703  // Attempt to emulate SAME_PRECONDITIONER below
704  KSPSetReusePreconditioner(mKspSolver, PETSC_TRUE);
705  }
706 #else
707  /*
708  * The preconditioner flag (last argument) in the following calls says
709  * how to reuse the preconditioner on subsequent iterations.
710  */
711  MatStructure preconditioner_over_successive_calls;
712 
713  if (mMatrixIsConstant)
714  {
715  preconditioner_over_successive_calls = SAME_PRECONDITIONER;
716  }
717  else
718  {
719  preconditioner_over_successive_calls = SAME_NONZERO_PATTERN;
720  }
721 #endif
722 
723  if (mMatNullSpace) // Adding null-space to the matrix (new style) has to happen *before* KSPSetOperators
724  {
725 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
726  // Setting null space in the KSP was deprecated in PETSc 3.6, but setting the null space
727  // for the matrix appeared in PETSc 3.3 so 3.3, 3.4, 3.5 can do either
728 
729  PETSCEXCEPT( MatSetNullSpace(mLhsMatrix, mMatNullSpace) );
730 #else
731  PETSCEXCEPT( KSPSetNullSpace(mKspSolver, mMatNullSpace) );
732 #endif
733  }
734 
736  {
737 #if ( PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5 ) //PETSc 3.5 or later
738  KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix);
739 #else
740  KSPSetOperators(mKspSolver, mLhsMatrix, mPrecondMatrix, preconditioner_over_successive_calls);
741 #endif
742  }
743  else
744  {
745 #if ( PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5 ) //PETSc 3.5 or later
746  KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix);
747 #else
748  KSPSetOperators(mKspSolver, mLhsMatrix, mLhsMatrix, preconditioner_over_successive_calls);
749 #endif
750  }
751 
752  // Set either absolute or relative tolerance of the KSP solver.
753  // The default is to use relative tolerance (1e-6)
755  {
756  KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
757  }
758  else
759  {
760  KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
761  }
762 
763  // Set ksp and pc types
764 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
765  if (mKspType == "chebychev")
766  {
767  KSPSetType(mKspSolver, "chebyshev");
768  }
769  else
770  {
771  KSPSetType(mKspSolver, mKspType.c_str());
772  }
773 #else
774  KSPSetType(mKspSolver, mKspType.c_str());
775 #endif
776  KSPGetPC(mKspSolver, &prec);
777 
778  // Turn off pre-conditioning if the system size is very small
779  if (is_small)
780  {
781  PCSetType(prec, PCNONE);
782  }
783  else
784  {
785 #ifdef TRACE_KSP
786  Timer::Reset();
787 #endif
788  if (mPcType == "blockdiagonal")
789  {
791 #ifdef TRACE_KSP
792  if (PetscTools::AmMaster())
793  {
794  Timer::Print("Purpose-build preconditioner creation");
795  }
796 #endif
797  }
798  else if (mPcType == "ldufactorisation")
799  {
801 #ifdef TRACE_KSP
802  if (PetscTools::AmMaster())
803  {
804  Timer::Print("Purpose-build preconditioner creation");
805  }
806 #endif
807  }
808  else if (mPcType == "twolevelsblockdiagonal")
809  {
810  if (!mpBathNodes)
811  {
812  TERMINATE("You must provide a list of bath nodes when using TwoLevelsBlockDiagonalPC");
813  }
815 #ifdef TRACE_KSP
816  if (PetscTools::AmMaster())
817  {
818  Timer::Print("Purpose-build preconditioner creation");
819  }
820 #endif
821 
822  }
823  else
824  {
825  PCSetType(prec, mPcType.c_str());
826  }
827  }
828 
829  KSPSetFromOptions(mKspSolver);
830 
831  if (lhsGuess)
832  {
833  // Assume that the user of this method will always be kind enough to give us a reasonable guess.
834  KSPSetInitialGuessNonzero(mKspSolver,PETSC_TRUE);
835  }
836  /*
837  * Non-adaptive Chebyshev: the required spectrum approximation is computed just once
838  * at the beginning of the simulation. This is done with two extra CG solves.
839  */
840  if (mKspType == "chebychev" && !mUseFixedNumberIterations)
841  {
842 #ifdef TRACE_KSP
843  Timer::Reset();
844 #endif
845  // Preconditioned matrix spectrum is approximated with CG
846  KSPSetType(mKspSolver,"cg");
847  KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
848  KSPSetUp(mKspSolver);
849 
850  VecDuplicate(mRhsVector, &chebyshev_lhs_vector);
851  if (lhsGuess)
852  {
853  VecCopy(lhsGuess, chebyshev_lhs_vector);
854  }
855 
856  // Smallest eigenvalue is approximated to default tolerance
857  KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
858 
859  PetscReal *r_eig = new PetscReal[mSize];
860  PetscReal *c_eig = new PetscReal[mSize];
861  PetscInt eigs_computed;
862  KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
863 
864  mEigMin = r_eig[0];
865 
866  // Largest eigenvalue is approximated to machine precision
867  KSPSetTolerances(mKspSolver, DBL_EPSILON, DBL_EPSILON, PETSC_DEFAULT, PETSC_DEFAULT);
868  KSPSetUp(mKspSolver);
869  if (lhsGuess)
870  {
871  VecCopy(lhsGuess, chebyshev_lhs_vector);
872  }
873 
874  KSPSolve(mKspSolver, mRhsVector, chebyshev_lhs_vector);
875  KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
876 
877  mEigMax = r_eig[eigs_computed-1];
878 
879 #ifdef TRACE_KSP
880  /*
881  * Under certain circumstances (see Golub&Overton 1988), underestimating
882  * the spectrum of the preconditioned operator improves convergence rate.
883  * See publication for a discussion and for definition of alpha and sigma_one.
884  */
885  if (PetscTools::AmMaster())
886  {
887  std::cout << "EIGS: ";
888  for (int index=0; index<eigs_computed; index++)
889  {
890  std::cout << r_eig[index] << ", ";
891  }
892  std::cout << std::endl;
893  }
894 
895  if (PetscTools::AmMaster()) std::cout << "EIGS "<< mEigMax << " " << mEigMin <<std::endl;
896  double alpha = 2/(mEigMax+mEigMin);
897  double sigma_one = 1 - alpha*mEigMin;
898  if (PetscTools::AmMaster()) std::cout << "sigma_1 = 1 - alpha*mEigMin = "<< sigma_one <<std::endl;
899 #endif
900 
901  // Set Chebyshev solver and max/min eigenvalues
902  assert(mKspType == "chebychev");
903 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
904  KSPSetType(mKspSolver, "chebyshev");
905  KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
906 #else
907  KSPSetType(mKspSolver, mKspType.c_str());
908  KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
909 #endif
910  KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
912  {
913  KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
914  }
915  else
916  {
917  KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
918  }
919 
920  delete[] r_eig;
921  delete[] c_eig;
922 
923 #ifdef TRACE_KSP
924  if (PetscTools::AmMaster())
925  {
926  Timer::Print("Computing extremal eigenvalues");
927  }
928 #endif
929  }
930 
931 #ifdef TRACE_KSP
932  Timer::Reset();
933 #endif
934 
935  KSPSetUp(mKspSolver);
936 
937  if (chebyshev_lhs_vector)
938  {
939  PetscTools::Destroy(chebyshev_lhs_vector);
940  }
941 
942 #ifdef TRACE_KSP
943  if (PetscTools::AmMaster())
944  {
945  Timer::Print("KSPSetUP (contains preconditioner creation for PETSc preconditioners)");
946  }
947 #endif
948 
949  mKspIsSetup = true;
950 
951  HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
952  }
953  else
954  {
955  if (mNonZerosUsed!=mat_info.nz_used)
956  {
957  WARNING("LinearSystem doesn't like the non-zero pattern of a matrix to change. (I think you changed it).");
958  mNonZerosUsed = mat_info.nz_used;
959  }
960 // PetscScalar norm;
961 // MatNorm(mLhsMatrix, NORM_FROBENIUS, &norm);
962 // if (fabs(norm - mMatrixNorm) > 0)
963 // {
964 // EXCEPTION("LinearSystem doesn't allow the matrix norm to change");
965 // }
966  }
967 
968  HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
969  // Create solution vector
971  Vec lhs_vector;
972  VecDuplicate(mRhsVector, &lhs_vector); // Sets the same size (doesn't copy)
973  if (lhsGuess)
974  {
975  VecCopy(lhsGuess, lhs_vector);
976  // If this wasn't done at construction time then it may be too late for this:
977  // KSPSetInitialGuessNonzero(mKspSolver, PETSC_TRUE);
978  // Is it possible to warn the user?
979  }
980 
981  // Check if the right hand side is small (but non-zero), PETSc can diverge immediately
982  // with a non-zero initial guess. Here we check for this and alter the initial guess to zero.
983  PetscReal rhs_norm;
984  VecNorm(mRhsVector, NORM_1, &rhs_norm);
985  double rhs_zero_tol = 1e-15;
986 
987  if (rhs_norm < rhs_zero_tol && rhs_norm > 0.0)
988  {
989  WARNING("Using zero initial guess due to small right hand side vector");
990  PetscVecTools::Zero(lhs_vector);
991  }
992 
993  HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
994 // // Double check that the mRhsVector contains sensible values
995 // double* p_rhs,* p_guess;
996 // VecGetArray(mRhsVector, &p_rhs);
997 // VecGetArray(lhsGuess, &p_guess);
998 // for (int global_index=mOwnershipRangeLo; global_index<mOwnershipRangeHi; global_index++)
999 // {
1000 // int local_index = global_index - mOwnershipRangeLo;
1001 // assert(!std::isnan(p_rhs[local_index]));
1002 // assert(!std::isnan(p_guess[local_index]));
1003 // if (p_rhs[local_index] != p_rhs[local_index])
1004 // {
1005 // std::cout << "********* PETSc nan\n";
1006 // assert(0);
1007 // }
1008 // }
1009 // std::cout << "b[0] = " << p_rhs[0] << ", guess[0] = " << p_guess[0] << "\n";
1010 // VecRestoreArray(mRhsVector, &p_rhs);
1011 // VecRestoreArray(lhsGuess, &p_guess);
1012 // // Test A*guess
1013 // Vec temp;
1014 // VecDuplicate(mRhsVector, &temp);
1015 // MatMult(mLhsMatrix, lhs_vector, temp);
1016 // double* p_temp;
1017 // VecGetArray(temp, &p_temp);
1018 // std::cout << "temp[0] = " << p_temp[0] << "\n";
1019 // VecRestoreArray(temp, &p_temp);
1020 // PetscTools::Destroy(temp);
1021 // // Dump the matrix to file
1022 // PetscViewer viewer;
1023 // PetscViewerASCIIOpen(PETSC_COMM_WORLD,"mat.output",&viewer);
1024 // MatView(mLhsMatrix, viewer);
1025 // PetscViewerFlush(viewer);
1026 // PetscViewerDestroy(PETSC_DESTROY_PARAM(viewer);
1027 // // Dump the rhs vector to file
1028 // PetscViewerASCIIOpen(PETSC_COMM_WORLD,"vec.output",&viewer);
1029 // PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1030 // VecView(mRhsVector, viewer);
1031 // PetscViewerFlush(viewer);
1032 // PetscViewerDestroy(PETSC_DESTROY_PARAM(viewer);
1033 
1034  try
1035  {
1036  HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
1037 
1038 #ifdef TRACE_KSP
1039  Timer::Reset();
1040 #endif
1041 
1042  // Current solve has to be done with tolerance-based stop criteria in order to record iterations taken
1044  {
1045 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
1046  KSPSetNormType(mKspSolver, KSP_PRECONDITIONED_NORM);
1047 #else
1048  KSPSetNormType(mKspSolver, KSP_NORM_PRECONDITIONED);
1049 #endif
1050 
1051 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
1053  {
1054  #if ( PETSC_VERSION_MINOR >= 5 ) //PETSc 3.5 or later
1055  KSPConvergedDefaultCreate(&mpConvergenceTestContext);
1056  #else
1057  KSPDefaultConvergedCreate(&mpConvergenceTestContext);
1058  #endif
1059  }
1060  #if ( PETSC_VERSION_MINOR >= 5 ) //PETSc 3.5 or later
1061  // From PETSc 3.5, KSPDefaultConverged became KSPConvergedDefault.
1062  KSPSetConvergenceTest(mKspSolver, KSPConvergedDefault, &mpConvergenceTestContext, PETSC_NULL);
1063  #else
1064  KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, &mpConvergenceTestContext, PETSC_NULL);
1065  #endif
1066 
1067 #else
1068  KSPSetConvergenceTest(mKspSolver, KSPDefaultConverged, PETSC_NULL);
1069 #endif
1070 
1072  {
1073  KSPSetTolerances(mKspSolver, DBL_EPSILON, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT);
1074  }
1075  else
1076  {
1077  KSPSetTolerances(mKspSolver, mTolerance, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
1078  }
1079 
1081  std::stringstream num_it_str;
1082  num_it_str << 1000;
1083  PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
1084 
1085  // Adaptive Chebyshev: reevaluate spectrum with cg
1086  if (mKspType == "chebychev")
1087  {
1088  // You can estimate preconditioned matrix spectrum with CG
1089  KSPSetType(mKspSolver,"cg");
1090  KSPSetComputeEigenvalues(mKspSolver, PETSC_TRUE);
1091  }
1092 
1093  KSPSetFromOptions(mKspSolver);
1094  KSPSetUp(mKspSolver);
1095  }
1096 
1097  PETSCEXCEPT(KSPSolve(mKspSolver, mRhsVector, lhs_vector));
1098  HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_LINEAR_SYSTEM);
1099 
1100 #ifdef TRACE_KSP
1101  PetscInt num_it;
1102  KSPGetIterationNumber(mKspSolver, &num_it);
1103  if (PetscTools::AmMaster())
1104  {
1105  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)
1106  Timer::Print("Solve");
1107  }
1108 
1109  mTotalNumIterations += num_it;
1110  if ((unsigned) num_it > mMaxNumIterations)
1111  {
1112  mMaxNumIterations = num_it;
1113  }
1114 #endif
1115 
1116  // Check that solver converged and throw if not
1117  KSPConvergedReason reason;
1118  KSPGetConvergedReason(mKspSolver, &reason);
1119 
1120  if (mUseFixedNumberIterations && PETSC_VERSION_MAJOR < 3)
1121  {
1122  WARNING("Not explicitly checking convergence reason when using fixed number of iterations and PETSc 2");
1123  }
1124  else
1125  {
1126  KSPEXCEPT(reason);
1127  }
1128 
1130  {
1131  // Adaptive Chebyshev: reevaluate spectrum with cg
1132  if (mKspType == "chebychev")
1133  {
1134  PetscReal *r_eig = new PetscReal[mSize];
1135  PetscReal *c_eig = new PetscReal[mSize];
1136  PetscInt eigs_computed;
1137  KSPComputeEigenvalues(mKspSolver, mSize, r_eig, c_eig, &eigs_computed);
1138 
1139  mEigMin = r_eig[0];
1140  /*
1141  * Using max() covers a borderline case found in TestChasteBenchmarksForPreDiCT where there's a big
1142  * gap in the spectrum between ~1.2 and ~2.5. Some reevaluations pick up 2.5 and others don't. If it
1143  * is not picked up, Chebyshev will diverge after 10 solves or so.
1144  */
1145  mEigMax = std::max(mEigMax,r_eig[eigs_computed-1]);
1146 
1147  delete[] r_eig;
1148  delete[] c_eig;
1149 
1150  assert(mKspType == "chebychev");
1151 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 3) //PETSc 3.3 or later
1152  KSPSetType(mKspSolver, "chebyshev");
1153  KSPChebyshevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
1154 #else
1155  KSPSetType(mKspSolver, mKspType.c_str());
1156  KSPChebychevSetEigenvalues(mKspSolver, mEigMax, mEigMin);
1157 #endif
1158  KSPSetComputeEigenvalues(mKspSolver, PETSC_FALSE);
1159  }
1160 
1161 #if ((PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR <= 2))
1162  if (mKspType == "chebychev")
1163  {
1164  // See #1695 for more details.
1165  EXCEPTION("Chebyshev with fixed number of iterations is known to be broken in PETSc <= 2.3.2");
1166  }
1167 
1168  KSPSetNormType(mKspSolver, KSP_NO_NORM);
1169 #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
1170  KSPSetNormType(mKspSolver, KSP_NORM_NONE);
1171 #else
1172  KSPSetNormType(mKspSolver, KSP_NORM_NO);
1173 #endif
1174 
1175 #if (PETSC_VERSION_MAJOR == 2)
1176  KSPSetConvergenceTest(mKspSolver, KSPSkipConverged, PETSC_NULL);
1177 #endif
1178 
1179  PetscInt num_it;
1180  KSPGetIterationNumber(mKspSolver, &num_it);
1181  std::stringstream num_it_str;
1182  num_it_str << num_it;
1183  PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
1184 
1185  KSPSetFromOptions(mKspSolver);
1186  KSPSetUp(mKspSolver);
1187 
1189  }
1190 
1191  mNumSolves++;
1192 
1193  }
1194  catch (const Exception& e)
1195  {
1196  // Destroy solution vector on error to avoid memory leaks
1197  PetscTools::Destroy(lhs_vector);
1198  throw e;
1199  }
1200 
1201  return lhs_vector;
1202 }
1203 
1205 {
1206  mPrecondMatrixIsNotLhs = precondIsDifferent;
1207 
1209  {
1210  if (mRowPreallocation == UINT_MAX)
1211  {
1212  /*
1213  * At the time of writing, this line will be reached if the constructor
1214  * with signature LinearSystem(Vec residualVector, Mat jacobianMatrix) is
1215  * called with jacobianMatrix=NULL and preconditioning matrix different
1216  * from lhs is used.
1217  *
1218  * If this combination is ever required you will need to work out
1219  * matrix allocation (mRowPreallocation) here.
1220  */
1221  NEVER_REACHED;
1222  }
1223 
1225  PetscTools::SetupMat(mPrecondMatrix, mSize, mSize, mRowPreallocation, local_size, local_size);
1226  }
1227 }
1228 
1229 void LinearSystem::SetUseFixedNumberIterations(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
1230 {
1231 
1232  mUseFixedNumberIterations = useFixedNumberIterations;
1233  mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
1234 }
1235 
1237 {
1238  if (mKspIsSetup)
1239  {
1240  KSPDestroy(PETSC_DESTROY_PARAM(mKspSolver));
1241  }
1242 
1243  mKspIsSetup = false;
1245 
1246  /*
1247  * Reset max number of iterations. This option is stored in the configuration database and
1248  * explicitely read in with KSPSetFromOptions() everytime a KSP object is created. Therefore,
1249  * destroying the KSP object will not ensure that it is set back to default.
1250  */
1252  std::stringstream num_it_str;
1253  num_it_str << 1000;
1254  PetscOptionsSetValue("-ksp_max_it", num_it_str.str().c_str());
1255 }
1256 
1257 // Serialization for Boost >= 1.36
void AssembleFinalLinearSystem()
Mat & rGetPrecondMatrix()
void FinalisePrecondMatrix()
unsigned mRowPreallocation
void FinaliseLhsMatrix()
void AddToRhsVectorElement(PetscInt row, double value)
void GetOwnershipRange(PetscInt &lo, PetscInt &hi)
bool mUseAbsoluteTolerance
void * mpConvergenceTestContext
void ZeroLinearSystem()
static double GetElement(Vec vector, PetscInt row)
void SetKspType(const char *kspType)
#define PETSC_DESTROY_PARAM(x)
Definition: PetscTools.hpp:69
void ZeroMatrixRowsWithValueOnDiagonal(std::vector< unsigned > &rRows, double diagonalValue)
unsigned GetSize() const
void ZeroRhsVector()
unsigned mEvaluateNumItsEveryNSolves
#define EXCEPTION(message)
Definition: Exception.hpp:143
void SetRelativeTolerance(double relativeTolerance)
Vec Solve(Vec lhsGuess=NULL)
static void SetElement(Vec vector, PetscInt row, double value)
void SetMatrixElement(PetscInt row, PetscInt col, double value)
void AddToMatrixElement(PetscInt row, PetscInt col, double value)
double mTolerance
Vec mDirichletBoundaryConditionsVector
static bool AmMaster()
Definition: PetscTools.cpp:120
Vec & rGetDirichletBoundaryConditionsVector()
PCTwoLevelsBlockDiagonal * mpTwoLevelsBlockDiagonalPC
bool mMatrixIsConstant
static Vec CreateVec(int size, int localSize=PETSC_DECIDE, bool ignoreOffProcEntries=true)
Definition: PetscTools.cpp:214
#define TERMINATE(message)
Definition: Exception.hpp:168
void ZeroMatrixRowsAndColumnsWithValueOnDiagonal(std::vector< unsigned > &rRowColIndices, double diagonalValue)
double GetRhsVectorElement(PetscInt row)
boost::shared_ptr< std::vector< PetscInt > > mpBathNodes
Mat & rGetLhsMatrix()
double GetMatrixElement(PetscInt row, PetscInt col)
unsigned mNumSolves
void RemoveNullSpace()
void ZeroMatrixColumn(PetscInt col)
#define NEVER_REACHED
Definition: Exception.hpp:206
void SetPrecondMatrixIsDifferentFromLhs(bool precondIsDifferent=true)
PetscInt mOwnershipRangeHi
std::string mKspType
Vec & rGetRhsVector()
static void AddToElement(Mat matrix, PetscInt row, PetscInt col, double value)
static void Print(std::string message)
Definition: Timer.cpp:59
void SetMatrixIsSymmetric(bool isSymmetric=true)
LinearSystem(PetscInt lhsVectorSize, unsigned rowPreallocation=UINT_MAX)
static Vec GetMatrixRowDistributed(Mat matrix, unsigned rowIndex)
static void Zero(Mat matrix)
PCBlockDiagonal * mpBlockDiagonalPC
void SetNullBasis(Vec nullbasis[], unsigned numberOfBases)
static void Display(Mat matrix)
static double GetElement(Mat matrix, PetscInt row, PetscInt col)
MatNullSpace mMatNullSpace
static void SetElement(Mat matrix, PetscInt row, PetscInt col, double value)
PCLDUFactorisation * mpLDUFactorisationPC
static void ZeroRowsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > &rRows, double diagonalValue)
static void ZeroRowsAndColumnsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > rowColIndices, double diagonalValue)
static void SwitchWriteMode(Mat matrix)
void SetPcType(const char *pcType, boost::shared_ptr< std::vector< PetscInt > > pBathNodes=boost::shared_ptr< std::vector< PetscInt > >())
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)
Definition: PetscTools.cpp:268
bool mDestroyMatAndVec
void SetAbsoluteTolerance(double absoluteTolerance)
void ZeroLhsMatrix()
Vec GetRhsVector() const
static void Zero(Vec vector)
void FinaliseRhsVector()
bool mForceSpectrumReevaluation
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
void SetRhsVectorElement(PetscInt row, double value)
Mat GetLhsMatrix() const
std::string mPcType
PetscReal mEigMax
static void Display(Vec vector)
void AssembleIntermediateLinearSystem()
void ResetKspSolver()
double mNonZerosUsed
static void SetOption(Mat matrix, MatOption option)
bool IsMatrixSymmetric()
static void Finalise(Mat matrix)
#define CHASTE_CLASS_EXPORT(T)
void SwitchWriteModeLhsMatrix()
PetscInt mSize
unsigned GetNumIterations() const
static void Reset()
Definition: Timer.cpp:44
static void AddToElement(Vec vector, PetscInt row, double value)
PetscInt mOwnershipRangeLo
void SetUseFixedNumberIterations(bool useFixedNumberIterations=true, unsigned evaluateNumItsEveryNSolves=UINT_MAX)
bool mPrecondMatrixIsNotLhs
void SetMatrixRow(PetscInt row, double value)
PetscReal mEigMin
Vec GetMatrixRowDistributed(unsigned rowIndex)
static void SetRow(Mat matrix, PetscInt row, double value)
void SetMatrixIsConstant(bool matrixIsConstant)
bool mUseFixedNumberIterations
static void ZeroColumn(Mat matrix, PetscInt col)
static void Finalise(Vec vector)
void DisplayMatrix()