Flow123d  jenkins-Flow123d-linux-release-multijob-198
linsys.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Wrappers for linear systems based on MPIAIJ and MATIS format.
27  * @author Jan Brezina
28  *
29  * Linear system only keep together matrix, RHS and the solution vector.
30  *
31  */
32 
33 
34 //=============================================================================
35 //
36 // LINER SYSTEM ROUTINES - linear system use wrapers of PETSc assembling routins
37 // in order to allow counting of allocation and filling of matrix to be done by the same code
38 //
39 // we want to allow allocations of nonlocal rows, to this end we count on- and off-processor
40 // members in an parallel Vector
41 //
42 //=============================================================================
43 
44 #ifndef LA_LINSYS_HH_
45 #define LA_LINSYS_HH_
46 
47 /**
48  * @brief Abstract linear system class.
49  *
50  * Linear system consists of Matrix, RHS and solution.
51  * It provides general methods for:
52  * - matrix preallocation
53  * - assembling matrix and RHS
54  * - application of linear constrains (Dirichlet conditions) either during assembly or
55  * on assembled system
56  * - solution of the system
57  * - output in Matlab format
58  *
59  * Method operates on the system as single object. But some methods for matrix only manipulation
60  * can be provided until we have matrix as separate class.
61  *
62  * TODO:
63  * - simplify constructors
64  * - introduce set_solution_vector() and set_rhs() in order to solve multiple systems with same matrix
65  * - simplify constructor, make common interface, rather introduce particular seters for particular solvers
66  * - which parameters expose through Input::Record
67  * - why we need lsize in constructors if we have Distribution
68  */
69 
70 #include "system/global_defs.h"
71 #include "la/distribution.hh"
72 #include "input/input_type.hh"
73 #include "input/accessors.hh"
74 
75 
76 #include <mpi.h>
77 
78 #include <vector>
79 #include <armadillo>
80 
81 // PETSc includes
82 #include "petscmat.h"
83 
84 
85 class LinSys
86 {
87 friend class SchurComplement;
88 public:
89  // Abstract Input Record for LinSys initialization
91 
92  typedef enum {
93  INSERT=INSERT_VALUES,
94  ADD=ADD_VALUES,
98  } SetValuesMode;
99 
100 protected:
101  typedef std::pair<unsigned,double> Constraint_;
103 
104 public:
105  /**
106  * Constructor.
107  * Constructor of abstract class should not be called directly, but is used for initialization of member common
108  * to all particular solvers.
109  *
110  * @param comm - MPI communicator
111  *
112  * TODO: Vector solution_ is now initialized to NULL, but it should be rather allocated
113  * in the constructor instead of the method set_solution().
114  */
115  LinSys(const Distribution *rows_ds)
116  : comm_( rows_ds->get_comm() ), status_( NONE ), lsize_( rows_ds->lsize() ), rows_ds_(rows_ds),
117  symmetric_( false ), positive_definite_( false ), negative_definite_( false ),
118  spd_via_symmetric_general_( false ), solution_(NULL), v_solution_(NULL)
119  {
120  int lsizeInt = static_cast<int>( rows_ds->lsize() );
121  int sizeInt;
122  MPI_Allreduce ( &lsizeInt, &sizeInt, 1, MPI_INT, MPI_SUM, comm_ );
123  size_ = static_cast<unsigned>( sizeInt );
124 
125  };
126 
127  /**
128  * Copy constructor.
129  */
130  LinSys(LinSys &other)
131  : r_tol_(other.r_tol_), a_tol_(other.a_tol_), max_it_(other.max_it_), comm_(other.comm_), status_(other.status_),
132  lsize_( other.rows_ds_->lsize() ), size_(other.size_), rows_ds_(other.rows_ds_), symmetric_(other.symmetric_),
137 
138  {
139  ASSERT( false, "Using copy constructor of LinSys is not allowed!");
140  set_solution(other.v_solution_);
141  };
142 
143  /**
144  * Returns global system size.
145  */
146  inline unsigned int size()
147  {
148  return size_;
149  }
150 
151  /**
152  * Returns local system size. (for partitioning of solution vectors)
153  * for PETSC_MPIAIJ it is also partitioning of the matrix
154  */
155  inline unsigned int vec_lsize()
156  {
157  return lsize_;
158  }
159 
160  /**
161  * Returns PETSC matrix (only for PETSC solvers)
162  *
163  * If matrix is changed, method set_matrix_changed() must be called.
164  * Example:
165  * @CODE
166  * MatDiagonalSet(schur->get_matrix(), new_diagonal, ADD_VALUES);
167  * schur->set_matrix_changed();
168  * @ENDCODE
169  */
170  virtual const Mat *get_matrix()
171  {
172  ASSERT( false, "Function get_matrix is not implemented for linsys type %s \n.", typeid(*this).name() );
173  return NULL;
174  }
175 
176  /**
177  * Returns RHS vector (only for PETSC solvers)
178  *
179  * If vector is changed, method set_rhs_changed() must be called.
180  * Example:
181  * @CODE
182  * VecScale(schur->get_rhs(), -1.0);
183  * schur->set_rhs_changed();
184  * @ENDCODE
185  */
186  virtual const Vec *get_rhs()
187  {
188  ASSERT( false, "Function get_rhs is not implemented for linsys type %s \n.", typeid(*this).name() );
189  return NULL;
190  }
191 
192  /**
193  * Sets matrix changed flag.
194  */
196  { matrix_changed_ = true;}
197 
198  /**
199  * Sets rhs changed flag (only for PETSC solvers)
200  */
202  { rhs_changed_ = true; }
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  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  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  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  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  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  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  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  * Returns information on relative solver accuracy
453  */
455  return r_tol_;
456  };
457 
458  /**
459  * Returns information on absolute solver accuracy
460  */
461  virtual double get_absolute_accuracy(){
462  return 0.0;
463  };
464 
465  /**
466  * Provides user knowledge about symmetry.
467  */
468  inline void set_symmetric(bool flag = true)
469  {
470  symmetric_ = flag;
471  if (!flag) {
472  set_positive_definite(false);
473  set_negative_definite(false);
474  }
475  }
476 
477  inline bool is_symmetric()
478  { return symmetric_; }
479 
480  /**
481  * Provides user knowledge about positive definiteness.
482  */
483  inline void set_positive_definite(bool flag = true)
484  {
485  positive_definite_ = flag;
486  if (flag) {
487  set_symmetric();
488  set_negative_definite(false);
489  }
490  }
491 
492  /**
493  * Provides user knowledge about negative definiteness.
494  */
495  inline void set_negative_definite(bool flag = true)
496  {
497  negative_definite_ = flag;
498  if (flag) {
499  set_symmetric();
500  set_positive_definite(false);
501  }
502  }
503 
504  inline bool is_positive_definite()
505  { return positive_definite_; }
506 
507  inline bool is_negative_definite()
508  { return negative_definite_; }
509 
510  /// TODO: In fact we want to know if the matrix is already preallocated
511  /// However to do this we need explicit finalisation of preallocating cycle.
512  inline bool is_new() {
513  return ( status_ == NONE );
514  }
515 
516  inline bool is_preallocated() {
517  return ( status_ == INSERT || status_ == ADD);
518  }
519 
520  /**
521  * Provides user knowledge about positive definiteness via symmetric general approach.
522  * This is useful for solving Darcy flow by mixed hybrid method, where blocks on subdomains are saddle point but
523  * interface among subdomains is only at the block of Lagrange multipliers and is symmetric positive definite.
524  * Problem condensed to interface can thus be solved by PCG method, although original problem is saddle point.
525  */
526  inline void set_spd_via_symmetric_general(bool flag = true)
527  {
529  if (flag) set_symmetric();
530  }
531 
533  { return spd_via_symmetric_general_; }
534 
535 
536  /**
537  * Output the system in the Matlab format possibly with given ordering.
538  * Rather we shoud provide output operator <<, since it is more flexible.
539  */
540  virtual void view()
541  {
542  ASSERT( false, "Function view is not implemented for linsys type %s \n.", typeid(*this).name() );
543  }
544 
545  /**
546  * Sets basic parameters of LinSys defined by user in input file and used to calculate
547  */
548  virtual void set_from_input(const Input::Record in_rec)
549  {
550  if (! in_rec.is_empty()) {
551  in_rec_ = Input::Record( in_rec );
552  r_tol_ = in_rec.val<double>("r_tol");
553  max_it_ = in_rec.val<int>("max_it");
554  a_tol_ = 0.01 * r_tol_;
555  }
556  }
557 
558  /**
559  * Get precision of solving
560  */
561  virtual double get_solution_precision() = 0;
562 
563  virtual ~LinSys()
564  {
565  PetscErrorCode ierr;
566  if ( solution_ ) { ierr = VecDestroy(&solution_); CHKERRV( ierr ); }
567  if ( own_solution_ ) delete[] v_solution_;
568  }
569 
570 protected:
571  double r_tol_; // relative tolerance of linear solver
572  double a_tol_; // absolute tolerance of linear solver
573  int max_it_; // maximum number of iterations of linear solver
574 
576  SetValuesMode status_; //!< Set value status of the linear system.
577 
578  const unsigned lsize_; //!< local number of matrix rows (non-overlapping division of rows)
579  unsigned size_; //!< global number of matrix rows, i.e. problem size
580 
581  const Distribution * rows_ds_; //!< final distribution of rows of MH matrix
582 
587 
588  bool matrix_changed_; //!< true if the matrix was changed since the last solve
589  bool rhs_changed_; //!< true if the right hand side was changed since the last solve
590 
591  Vec solution_; //!< PETSc vector constructed with vb array.
592  double *v_solution_; //!< local solution array pointing into Vec solution_
593  bool own_solution_; //!< Indicates if the solution array has been allocated by this class
594 
595  double residual_norm_; //!< local solution array pointing into Vec solution_
596 
598 
599  std::vector<double> globalSolution_; //!< global solution in numbering for linear system
600 
601  Input::Record in_rec_; // structure contained parameters of LinSys defined in input file
602 
603 };
604 
605 #endif /* LA_LINSYS_HH_ */
unsigned int size()
Definition: linsys.hh:146
std::vector< double > globalSolution_
global solution in numbering for linear system
Definition: linsys.hh:599
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:576
bool is_new()
Definition: linsys.hh:512
virtual void start_insert_assembly()
Definition: linsys.hh:304
bool is_negative_definite()
Definition: linsys.hh:507
bool spd_via_symmetric_general_
Definition: linsys.hh:586
std::pair< unsigned, double > Constraint_
Definition: linsys.hh:101
double get_residual_norm()
Definition: linsys.hh:447
void set_symmetric(bool flag=true)
Definition: linsys.hh:468
double * v_solution_
local solution array pointing into Vec solution_
Definition: linsys.hh:592
int MPI_Comm
Definition: mpi.h:141
virtual void view()
Definition: linsys.hh:540
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:548
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:588
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:578
void set_rhs_changed()
Definition: linsys.hh:201
bool symmetric_
Definition: linsys.hh:583
unsigned int vec_lsize()
Definition: linsys.hh:155
ConstraintVec_ constraints_
Definition: linsys.hh:597
#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:495
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:526
unsigned size_
global number of matrix rows, i.e. problem size
Definition: linsys.hh:579
Input::Record in_rec_
Definition: linsys.hh:601
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:589
LinSys(const Distribution *rows_ds)
Definition: linsys.hh:115
virtual ~LinSys()
Definition: linsys.hh:563
virtual void apply_constrains(double scalar)=0
bool is_preallocated()
Definition: linsys.hh:516
void mat_set_value(int row, int col, double val)
Definition: linsys.hh:324
Global macros to enhance readability and debugging, general constants.
virtual PetscErrorCode set_rhs(Vec &rhs)
Definition: linsys.hh:229
#define ASSERT(...)
Definition: global_defs.h:121
MPI_Comm comm_
Definition: linsys.hh:575
static Input::Type::AbstractRecord input_type
Definition: linsys.hh:90
double get_relative_accuracy()
Definition: linsys.hh:454
double * get_solution_array()
Definition: linsys.hh:280
int max_it_
Definition: linsys.hh:573
const Vec & get_solution()
Definition: linsys.hh:256
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
const Ret val(const string &key) const
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:581
Class for declaration of polymorphic Record.
Definition: type_record.hh:487
double a_tol_
Definition: linsys.hh:572
SetValuesMode
Definition: linsys.hh:92
bool is_empty() const
Definition: accessors.hh:401
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:595
bool is_rhs_changed()
Definition: linsys.hh:213
bool negative_definite_
Definition: linsys.hh:585
void set_solution(double *sol_array)
Definition: linsys.hh:264
virtual const Vec * get_rhs()
Definition: linsys.hh:186
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:247
bool is_symmetric()
Definition: linsys.hh:477
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:102
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_
Definition: linsys.hh:571
#define MPI_INT
Definition: mpi.h:160
Abstract linear system class.
bool is_positive_definite()
Definition: linsys.hh:504
void set_matrix_changed()
Definition: linsys.hh:195
void set_positive_definite(bool flag=true)
Definition: linsys.hh:483
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
Definition: linsys.hh:220
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
bool positive_definite_
Definition: linsys.hh:584
bool own_solution_
Indicates if the solution array has been allocated by this class.
Definition: linsys.hh:593
virtual const Mat * get_matrix()
Definition: linsys.hh:170
virtual double get_absolute_accuracy()
Definition: linsys.hh:461
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:591
LinSys(LinSys &other)
Definition: linsys.hh:130
bool is_spd_via_symmetric_general()
Definition: linsys.hh:532
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