Flow123d  release_1.8.2-1603-g0109a2b
linsys.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file linsys.hh
15  * @brief Wrappers for linear systems based on MPIAIJ and MATIS format.
16  * @author Jan Brezina
17  *
18  * Linear system only keep together matrix, RHS and the solution vector.
19  *
20  */
21 
22 //=============================================================================
23 //
24 // LINER SYSTEM ROUTINES - linear system use wrapers of PETSc assembling routins
25 // in order to allow counting of allocation and filling of matrix to be done by the same code
26 //
27 // we want to allow allocations of nonlocal rows, to this end we count on- and off-processor
28 // members in an parallel Vector
29 //
30 //=============================================================================
31 
32 #ifndef LA_LINSYS_HH_
33 #define LA_LINSYS_HH_
34 
35 /**
36  * @brief Abstract linear system class.
37  *
38  * Linear system consists of Matrix, RHS and solution.
39  * It provides general methods for:
40  * - matrix preallocation
41  * - assembling matrix and RHS
42  * - application of linear constrains (Dirichlet conditions) either during assembly or
43  * on assembled system
44  * - solution of the system
45  * - output in Matlab format
46  *
47  * Method operates on the system as single object. But some methods for matrix only manipulation
48  * can be provided until we have matrix as separate class.
49  *
50  * TODO:
51  * - simplify constructors
52  * - introduce set_solution_vector() and set_rhs() in order to solve multiple systems with same matrix
53  * - simplify constructor, make common interface, rather introduce particular seters for particular solvers
54  * - which parameters expose through Input::Record
55  * - why we need lsize in constructors if we have Distribution
56  */
57 
58 #include "system/global_defs.h"
59 #include "la/distribution.hh"
61 #include "input/accessors.hh"
62 
63 
64 #include <mpi.h>
65 
66 #include <vector>
67 #include <armadillo>
68 
69 // PETSc includes
70 #include "petscmat.h"
71 
72 
73 class LinSys
74 {
75 friend class SchurComplement;
76 public:
77  // Abstract Input Record for LinSys initialization
79 
80  typedef enum {
81  INSERT=INSERT_VALUES,
82  ADD=ADD_VALUES,
86  } SetValuesMode;
87 
88 protected:
89  typedef std::pair<unsigned,double> Constraint_;
91 
92 public:
93  /**
94  * Constructor.
95  * Constructor of abstract class should not be called directly, but is used for initialization of member common
96  * to all particular solvers.
97  *
98  * @param comm - MPI communicator
99  *
100  * TODO: Vector solution_ is now initialized to NULL, but it should be rather allocated
101  * in the constructor instead of the method set_solution().
102  */
103  LinSys(const Distribution *rows_ds)
104  : comm_( rows_ds->get_comm() ), status_( NONE ), lsize_( rows_ds->lsize() ), rows_ds_(rows_ds),
105  symmetric_( false ), positive_definite_( false ), negative_definite_( false ),
106  spd_via_symmetric_general_( false ), solution_(NULL), v_solution_(NULL)
107  {
108  int lsizeInt = static_cast<int>( rows_ds->lsize() );
109  int sizeInt;
110  MPI_Allreduce ( &lsizeInt, &sizeInt, 1, MPI_INT, MPI_SUM, comm_ );
111  size_ = static_cast<unsigned>( sizeInt );
112 
116  };
117 
118  /**
119  * Copy constructor.
120  */
121  LinSys(LinSys &other)
122  : r_tol_(other.r_tol_), a_tol_(other.a_tol_), max_it_(other.max_it_), comm_(other.comm_), status_(other.status_),
123  lsize_( other.rows_ds_->lsize() ), size_(other.size_), rows_ds_(other.rows_ds_), symmetric_(other.symmetric_),
128 
129  {
130  OLD_ASSERT( false, "Using copy constructor of LinSys is not allowed!");
131  set_solution(other.v_solution_);
132  };
133 
134  /**
135  * Returns global system size.
136  */
137  inline unsigned int size()
138  {
139  return size_;
140  }
141 
142  /**
143  * Returns local system size. (for partitioning of solution vectors)
144  * for PETSC_MPIAIJ it is also partitioning of the matrix
145  */
146  inline unsigned int vec_lsize()
147  {
148  return lsize_;
149  }
150 
151  /**
152  * Returns PETSC matrix (only for PETSC solvers)
153  *
154  * If matrix is changed, method set_matrix_changed() must be called.
155  * Example:
156  * @CODE
157  * MatDiagonalSet(schur->get_matrix(), new_diagonal, ADD_VALUES);
158  * schur->set_matrix_changed();
159  * @ENDCODE
160  */
161  virtual const Mat *get_matrix()
162  {
163  OLD_ASSERT( false, "Function get_matrix is not implemented for linsys type %s \n.", typeid(*this).name() );
164  return NULL;
165  }
166 
167  /**
168  * Returns RHS vector (only for PETSC solvers)
169  *
170  * If vector is changed, method set_rhs_changed() must be called.
171  * Example:
172  * @CODE
173  * VecScale(schur->get_rhs(), -1.0);
174  * schur->set_rhs_changed();
175  * @ENDCODE
176  */
177  virtual const Vec *get_rhs()
178  {
179  OLD_ASSERT( false, "Function get_rhs is not implemented for linsys type %s \n.", typeid(*this).name() );
180  return NULL;
181  }
182 
183  /**
184  * Sets matrix changed flag.
185  */
187  { matrix_changed_ = true;}
188 
189  /**
190  * Sets rhs changed flag (only for PETSC solvers)
191  */
193  { rhs_changed_ = true; }
194 
195  /**
196  * Set relative tolerance, absolute tolerance, and maximum number of iterations of the linear solver.
197  *
198  * For each of these three parameters we first look for the value at user input
199  * if not set we use the value provided to this method and finally the default values are
200  * set by the call of this method in the constructor.
201  */
202  virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it) = 0;
203 
204  /**
205  * Returns true if the system matrix has changed since the last solve.
206  */
208  { return matrix_changed_;}
209 
210  /**
211  * Returns true if the system RHS has changed since the last solve.
212  */
214  { return rhs_changed_;}
215 
216 
217  /**
218  * Sets PETSC matrix (only for PETSC solvers)
219  */
220  virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
221  {
222  OLD_ASSERT( false, "Function set_matrix is not implemented for linsys type %s \n.", typeid(*this).name() );
223  return 0;
224  }
225 
226  /**
227  * Sets RHS vector (only for PETSC solvers)
228  */
229  virtual PetscErrorCode set_rhs(Vec &rhs)
230  {
231  OLD_ASSERT( false, "Function set_rhs is not implemented for linsys type %s \n.", typeid(*this).name() );
232  return 0;
233  }
234 
235  /**
236  * Clears entries of the matrix
237  */
238  virtual PetscErrorCode mat_zero_entries()
239  {
240  OLD_ASSERT( false, "Function mat_zero_entries is not implemented for linsys type %s \n.", typeid(*this).name() );
241  return 0;
242  }
243 
244  /**
245  * Clears entries of the right-hand side
246  */
247  virtual PetscErrorCode rhs_zero_entries()
248  {
249  OLD_ASSERT( false, "Function vec_zero_entries is not implemented for linsys type %s \n.", typeid(*this).name() );
250  return 0;
251  }
252 
253  /**
254  * Returns PETSC vector with solution. Underlying array can be provided on construction.
255  */
256  const Vec &get_solution()
257  {
258  return solution_;
259  }
260 
261  /**
262  * Create PETSc solution
263  */
264  void set_solution(double *sol_array) {
265  if (sol_array == NULL) {
266  v_solution_ = new double[ rows_ds_->lsize() + 1 ];
267  own_solution_ = true;
268  }
269  else {
270  v_solution_ = sol_array;
271  own_solution_ = false;
272  }
273  PetscErrorCode ierr;
274  ierr = VecCreateMPIWithArray( comm_,1, rows_ds_->lsize(), PETSC_DECIDE, v_solution_, &solution_ ); CHKERRV( ierr );
275  }
276 
277  /**
278  * Returns PETSC subarray with solution. Underlying array can be provided on construction.
279  */
281  {
282  return v_solution_;
283  }
284 
285  /**
286  * Switch linear system into allocating assembly. (only for PETSC_MPIAIJ_preallocate_by_assembly)
287  */
288  virtual void start_allocation()
289  {
290  OLD_ASSERT( false, "Function start_allocation is not implemented for linsys type %s \n.", typeid(*this).name() );
291  }
292 
293  /**
294  * Switch linear system into adding assembly. (the only one supported by triplets ??)
295  */
296  virtual void start_add_assembly()
297  {
298  OLD_ASSERT( false, "Function start_add_assembly is not implemented for linsys type %s \n.", typeid(*this).name() );
299  }
300 
301  /**
302  * Switch linear system into insert assembly. (not currently used)
303  */
304  virtual void start_insert_assembly()
305  {
306  OLD_ASSERT( false, "Function start_insert_assembly is not implemented for linsys type %s \n.", typeid(*this).name() );
307  }
308 
309  /**
310  * Finish assembly of the whole system. For PETSC this should call MatEndAssembly with MAT_FINAL_ASSEMBLY
311  */
312  virtual void finish_assembly( )=0;
313 
314  /**
315  * Assembly full rectangular submatrix into the system matrix.
316  * Should be virtual, implemented differently in particular solvers.
317  */
318  virtual void mat_set_values(int nrow,int *rows,int ncol,int *cols,double *vals)=0;
319 
320  /**
321  * Shortcut for assembling just one element into the matrix.
322  * Similarly we can provide method accepting armadillo matrices.
323  */
324  void mat_set_value(int row,int col,double val)
325  { mat_set_values(1,&row,1,&col,&val); }
326 
327  /**
328  * Set values of the system right-hand side.
329  * Should be virtual, implemented differently in particular solvers.
330  */
331  virtual void rhs_set_values(int nrow,int *rows,double *vals)=0;
332 
333  /**
334  * Shorcut for assembling just one element into RHS vector.
335  */
336  void rhs_set_value(int row,double val)
337  { rhs_set_values(1,&row,&val); }
338 
339 
340  /// Set values in the system matrix and values in the right-hand side vector on corresponding rows.
341  inline void set_values(int nrow,int *rows,int ncol,int *cols,PetscScalar *mat_vals, PetscScalar *rhs_vals)
342  {
343  mat_set_values(nrow, rows, ncol, cols, mat_vals);
344  rhs_set_values(nrow, rows, rhs_vals);
345  }
346 
347  /**
348  * Shortcut to assembly into matrix and RHS in one call, possibly apply Dirichlet boundary conditions.
349  * @p row_dofs - are global indices of rows of dense @p matrix and rows of dense vector @rhs in global system
350  * @p col_dofs - are global indices of columns of the matrix, and possibly
351  *
352  * Application of Dirichlet conditions:
353  * 1) Rows with negative dofs are set to zero.
354  * 2) Cols with negative dofs are eliminated.
355  * 3) If there are entries on global diagonal. We determine value K either from diagonal of local matrix, or (if it is zero) from
356  * diagonal average.
357  *
358  * Caveats:
359  * - can not set dirichlet condition on zero dof
360  * - Armadillo stores matrix in column first form (Fortran like) which makes it not well suited
361  * for passing local matrices.
362  *
363  */
364  void set_values(std::vector<int> &row_dofs, std::vector<int> &col_dofs,
365  const arma::mat &matrix, const arma::vec &rhs,
366  const arma::vec &row_solution, const arma::vec &col_solution)
367 
368  {
369  arma::mat tmp = matrix.t();
370  arma::vec tmp_rhs = rhs;
371  bool negative_row = false;
372  bool negative_col = false;
373 
374  for(unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
375  if (row_dofs[l_row] < 0) {
376  tmp_rhs(l_row)=0.0;
377  tmp.col(l_row).zeros();
378  negative_row=true;
379  }
380 
381  for(unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
382  if (col_dofs[l_col] < 0) {
383  tmp_rhs -= matrix.col(l_col) * col_solution[l_col];
384  tmp.row(l_col).zeros();
385  negative_col=true;
386  }
387 
388 
389  if (negative_row && negative_col) {
390  // look for diagonal entry
391  for(unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
392  if (row_dofs[l_row] < 0)
393  for(unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
394  if (col_dofs[l_col] < 0 && row_dofs[l_row] == col_dofs[l_col]) {
395  double new_diagonal = fabs(matrix.at(l_row,l_col));
396  if (new_diagonal == 0.0) {
397  if (matrix.is_square()) {
398  new_diagonal = arma::sum( abs(matrix.diag())) / matrix.n_rows;
399  } else {
400  new_diagonal = arma::accu( abs(matrix) ) / matrix.n_elem;
401  }
402  }
403  tmp.at(l_col, l_row) = new_diagonal;
404  tmp_rhs(l_row) = new_diagonal * row_solution[l_row];
405  }
406 
407  }
408 
409  if (negative_row)
410  for(int &row : row_dofs) row=abs(row);
411 
412  if (negative_col)
413  for(int &col : col_dofs) col=abs(col);
414 
415 
416  mat_set_values(row_dofs.size(), const_cast<int *>(&(row_dofs[0])),
417  col_dofs.size(), const_cast<int *>(&(col_dofs[0])), tmp.memptr() );
418  rhs_set_values(row_dofs.size(), const_cast<int *>(&(row_dofs[0])), tmp_rhs.memptr() );
419  }
420 
421  /**
422  * Adds Dirichlet constrain.
423  * @param row - global number of row that should be eliminated.
424  * @param value - solution value at the given row
425  */
426  void add_constraint(int row, double value) {
427 
428  constraints_.push_back( Constraint_( static_cast<unsigned>( row ), value ) );
429  }
430 
431  /**
432  * Apply constrains to assembled matrix. Constrains are given by pairs: global row index, value.
433  * i.e. typedef pair<unsigned int, double> Constrain;
434  *
435  * What is th meaning of ( const double factor ) form Cambridge code?
436  */
437  virtual void apply_constrains( double scalar )=0;
438 
439  /**
440  * Solve the system and return convergence reason.
441  */
442  virtual int solve()=0;
443 
444  /**
445  * Returns norm of the initial right hand side
446  */
448  return residual_norm_;
449  };
450 
451  /**
452  * Explicitly compute residual and its norm for current solution.
453  */
454  virtual double compute_residual() =0;
455 
456  /**
457  * Returns information on relative solver accuracy
458  */
460  return r_tol_;
461  };
462 
463  /**
464  * Returns information on absolute solver accuracy
465  */
466  virtual double get_absolute_accuracy(){
467  return 0.0;
468  };
469 
470  /**
471  * Provides user knowledge about symmetry.
472  */
473  inline void set_symmetric(bool flag = true)
474  {
475  symmetric_ = flag;
476  if (!flag) {
477  set_positive_definite(false);
478  set_negative_definite(false);
479  }
480  }
481 
482  inline bool is_symmetric()
483  { return symmetric_; }
484 
485  /**
486  * Provides user knowledge about positive definiteness.
487  */
488  inline void set_positive_definite(bool flag = true)
489  {
490  positive_definite_ = flag;
491  if (flag) {
492  set_symmetric();
493  set_negative_definite(false);
494  }
495  }
496 
497  /**
498  * Provides user knowledge about negative definiteness.
499  */
500  inline void set_negative_definite(bool flag = true)
501  {
502  negative_definite_ = flag;
503  if (flag) {
504  set_symmetric();
505  set_positive_definite(false);
506  }
507  }
508 
509  inline bool is_positive_definite()
510  { return positive_definite_; }
511 
512  inline bool is_negative_definite()
513  { return negative_definite_; }
514 
515  /// TODO: In fact we want to know if the matrix is already preallocated
516  /// However to do this we need explicit finalisation of preallocating cycle.
517  inline bool is_new() {
518  return ( status_ == NONE );
519  }
520 
521  inline bool is_preallocated() {
522  return ( status_ == INSERT || status_ == ADD);
523  }
524 
525  /**
526  * Provides user knowledge about positive definiteness via symmetric general approach.
527  * This is useful for solving Darcy flow by mixed hybrid method, where blocks on subdomains are saddle point but
528  * interface among subdomains is only at the block of Lagrange multipliers and is symmetric positive definite.
529  * Problem condensed to interface can thus be solved by PCG method, although original problem is saddle point.
530  */
531  inline void set_spd_via_symmetric_general(bool flag = true)
532  {
534  if (flag) set_symmetric();
535  }
536 
538  { return spd_via_symmetric_general_; }
539 
540 
541  /**
542  * Output the system in the Matlab format possibly with given ordering.
543  * Rather we shoud provide output operator <<, since it is more flexible.
544  */
545  virtual void view()
546  {
547  OLD_ASSERT( false, "Function view is not implemented for linsys type %s \n.", typeid(*this).name() );
548  }
549 
550  /**
551  * Sets basic parameters of LinSys defined by user in input file and used to calculate
552  */
553  virtual void set_from_input(const Input::Record in_rec)
554  {
555  in_rec_ = in_rec;
557  }
558 
559  /**
560  * Get precision of solving
561  */
562  virtual double get_solution_precision() = 0;
563 
564  virtual ~LinSys()
565  {
566  PetscErrorCode ierr;
567  if ( solution_ ) { ierr = VecDestroy(&solution_); CHKERRV( ierr ); }
568  if ( own_solution_ ) delete[] v_solution_;
569  }
570 
571 protected:
572  // Default values initialized in constructor
573  static constexpr double default_r_tol_ = 1e-7;
574  static constexpr double default_a_tol_ = 1e-11;
575  static constexpr unsigned int default_max_it_ = 1000;
576 
577  double r_tol_; ///< relative tolerance of linear solver
578  double a_tol_; ///< absolute tolerance of linear solver
579  unsigned int max_it_; ///< maximum number of iterations of linear solver
580 
582  SetValuesMode status_; //!< Set value status of the linear system.
583 
584  const unsigned lsize_; //!< local number of matrix rows (non-overlapping division of rows)
585  unsigned size_; //!< global number of matrix rows, i.e. problem size
586 
587  const Distribution * rows_ds_; //!< final distribution of rows of MH matrix
588 
593 
594  bool matrix_changed_; //!< true if the matrix was changed since the last solve
595  bool rhs_changed_; //!< true if the right hand side was changed since the last solve
596 
597  Vec solution_; //!< PETSc vector constructed with vb array.
598  double *v_solution_; //!< local solution array pointing into Vec solution_
599  bool own_solution_; //!< Indicates if the solution array has been allocated by this class
600 
601  double residual_norm_; //!< local solution array pointing into Vec solution_
602 
603  ConstraintVec_ constraints_;
604 
605  std::vector<double> globalSolution_; //!< global solution in numbering for linear system
606 
607  Input::Record in_rec_; // structure contained parameters of LinSys defined in input file
608 
609 };
610 
611 #endif /* LA_LINSYS_HH_ */
unsigned int size()
Definition: linsys.hh:137
std::vector< double > globalSolution_
global solution in numbering for linear system
Definition: linsys.hh:605
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:582
bool is_new()
Definition: linsys.hh:517
virtual void start_insert_assembly()
Definition: linsys.hh:304
bool is_negative_definite()
Definition: linsys.hh:512
bool spd_via_symmetric_general_
Definition: linsys.hh:592
std::pair< unsigned, double > Constraint_
Definition: linsys.hh:89
double get_residual_norm()
Definition: linsys.hh:447
void set_symmetric(bool flag=true)
Definition: linsys.hh:473
double * v_solution_
local solution array pointing into Vec solution_
Definition: linsys.hh:598
int MPI_Comm
Definition: mpi.h:141
virtual void view()
Definition: linsys.hh:545
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:553
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:594
virtual void start_add_assembly()
Definition: linsys.hh:296
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:238
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
void set_values(std::vector< int > &row_dofs, std::vector< int > &col_dofs, const arma::mat &matrix, const arma::vec &rhs, const arma::vec &row_solution, const arma::vec &col_solution)
Definition: linsys.hh:364
const unsigned lsize_
local number of matrix rows (non-overlapping division of rows)
Definition: linsys.hh:584
void set_rhs_changed()
Definition: linsys.hh:192
static constexpr double default_r_tol_
Definition: linsys.hh:573
bool symmetric_
Definition: linsys.hh:589
unsigned int vec_lsize()
Definition: linsys.hh:146
ConstraintVec_ constraints_
Definition: linsys.hh:603
#define MPI_SUM
Definition: mpi.h:196
bool is_matrix_changed()
Definition: linsys.hh:207
void set_negative_definite(bool flag=true)
Definition: linsys.hh:500
virtual void start_allocation()
Definition: linsys.hh:288
virtual void finish_assembly()=0
virtual double get_solution_precision()=0
void add_constraint(int row, double value)
Definition: linsys.hh:426
void set_spd_via_symmetric_general(bool flag=true)
Definition: linsys.hh:531
unsigned size_
global number of matrix rows, i.e. problem size
Definition: linsys.hh:585
Input::Record in_rec_
Definition: linsys.hh:607
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:595
LinSys(const Distribution *rows_ds)
Definition: linsys.hh:103
virtual ~LinSys()
Definition: linsys.hh:564
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
virtual void apply_constrains(double scalar)=0
bool is_preallocated()
Definition: linsys.hh:521
void mat_set_value(int row, int col, double val)
Definition: linsys.hh:324
#define OLD_ASSERT(...)
Definition: global_defs.h:128
Global macros to enhance readability and debugging, general constants.
virtual PetscErrorCode set_rhs(Vec &rhs)
Definition: linsys.hh:229
MPI_Comm comm_
Definition: linsys.hh:581
double get_relative_accuracy()
Definition: linsys.hh:459
unsigned int max_it_
maximum number of iterations of linear solver
Definition: linsys.hh:579
double * get_solution_array()
Definition: linsys.hh:280
const Vec & get_solution()
Definition: linsys.hh:256
static constexpr unsigned int default_max_it_
Definition: linsys.hh:575
Accessor to the data with type Type::Record.
Definition: accessors.hh:277
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:587
double a_tol_
absolute tolerance of linear solver
Definition: linsys.hh:578
SetValuesMode
Definition: linsys.hh:80
static constexpr double default_a_tol_
Definition: linsys.hh:574
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:601
bool is_rhs_changed()
Definition: linsys.hh:213
bool negative_definite_
Definition: linsys.hh:591
void set_solution(double *sol_array)
Definition: linsys.hh:264
Class for declaration of polymorphic Record.
virtual const Vec * get_rhs()
Definition: linsys.hh:177
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:247
bool is_symmetric()
Definition: linsys.hh:482
Support classes for parallel programing.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
std::vector< Constraint_ > ConstraintVec_
Definition: linsys.hh:90
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows...
Definition: linsys.hh:341
double r_tol_
relative tolerance of linear solver
Definition: linsys.hh:577
#define MPI_INT
Definition: mpi.h:160
Abstract linear system class.
bool is_positive_definite()
Definition: linsys.hh:509
void set_matrix_changed()
Definition: linsys.hh:186
void set_positive_definite(bool flag=true)
Definition: linsys.hh:488
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
Definition: linsys.hh:220
static Input::Type::Abstract & get_input_type()
Definition: linsys.cc:29
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
bool positive_definite_
Definition: linsys.hh:590
bool own_solution_
Indicates if the solution array has been allocated by this class.
Definition: linsys.hh:599
virtual const Mat * get_matrix()
Definition: linsys.hh:161
virtual double get_absolute_accuracy()
Definition: linsys.hh:466
virtual double compute_residual()=0
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:597
LinSys(LinSys &other)
Definition: linsys.hh:121
bool is_spd_via_symmetric_general()
Definition: linsys.hh:537
virtual int solve()=0
void rhs_set_value(int row, double val)
Definition: linsys.hh:336
unsigned int lsize(int proc) const
get local size