Flow123d  release_2.2.0-20-gb8056ca
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  own_solution_(false)
108  {
109  int lsizeInt = static_cast<int>( rows_ds->lsize() );
110  int sizeInt;
111  MPI_Allreduce ( &lsizeInt, &sizeInt, 1, MPI_INT, MPI_SUM, comm_ );
112  size_ = static_cast<unsigned>( sizeInt );
113 
117  };
118 
119  /**
120  * Copy constructor.
121  */
122  LinSys(LinSys &other)
123  : r_tol_(other.r_tol_), a_tol_(other.a_tol_), max_it_(other.max_it_), comm_(other.comm_), status_(other.status_),
124  lsize_( other.rows_ds_->lsize() ), size_(other.size_), rows_ds_(other.rows_ds_), symmetric_(other.symmetric_),
129 
130  {
131  OLD_ASSERT( false, "Using copy constructor of LinSys is not allowed!");
132  set_solution(other.v_solution_);
133  };
134 
135  /**
136  * Returns global system size.
137  */
138  inline unsigned int size()
139  {
140  return size_;
141  }
142 
143  /**
144  * Returns local system size. (for partitioning of solution vectors)
145  * for PETSC_MPIAIJ it is also partitioning of the matrix
146  */
147  inline unsigned int vec_lsize()
148  {
149  return lsize_;
150  }
151 
152  /**
153  * Returns PETSC matrix (only for PETSC solvers)
154  *
155  * If matrix is changed, method set_matrix_changed() must be called.
156  * Example:
157  * @CODE
158  * MatDiagonalSet(schur->get_matrix(), new_diagonal, ADD_VALUES);
159  * schur->set_matrix_changed();
160  * @ENDCODE
161  */
162  virtual const Mat *get_matrix()
163  {
164  OLD_ASSERT( false, "Function get_matrix is not implemented for linsys type %s \n.", typeid(*this).name() );
165  return NULL;
166  }
167 
168  /**
169  * Returns RHS vector (only for PETSC solvers)
170  *
171  * If vector is changed, method set_rhs_changed() must be called.
172  * Example:
173  * @CODE
174  * VecScale(schur->get_rhs(), -1.0);
175  * schur->set_rhs_changed();
176  * @ENDCODE
177  */
178  virtual const Vec *get_rhs()
179  {
180  OLD_ASSERT( false, "Function get_rhs is not implemented for linsys type %s \n.", typeid(*this).name() );
181  return NULL;
182  }
183 
184  /**
185  * Sets matrix changed flag.
186  */
188  { matrix_changed_ = true;}
189 
190  /**
191  * Sets rhs changed flag (only for PETSC solvers)
192  */
194  { rhs_changed_ = true; }
195 
196  /**
197  * Set relative tolerance, absolute tolerance, and maximum number of iterations of the linear solver.
198  *
199  * For each of these three parameters we first look for the value at user input
200  * if not set we use the value provided to this method and finally the default values are
201  * set by the call of this method in the constructor.
202  */
203  virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it) = 0;
204 
205  /**
206  * Returns true if the system matrix has changed since the last solve.
207  */
209  { return matrix_changed_;}
210 
211  /**
212  * Returns true if the system RHS has changed since the last solve.
213  */
215  { return rhs_changed_;}
216 
217 
218  /**
219  * Sets PETSC matrix (only for PETSC solvers)
220  */
221  virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
222  {
223  OLD_ASSERT( false, "Function set_matrix is not implemented for linsys type %s \n.", typeid(*this).name() );
224  return 0;
225  }
226 
227  /**
228  * Sets RHS vector (only for PETSC solvers)
229  */
230  virtual PetscErrorCode set_rhs(Vec &rhs)
231  {
232  OLD_ASSERT( false, "Function set_rhs is not implemented for linsys type %s \n.", typeid(*this).name() );
233  return 0;
234  }
235 
236  /**
237  * Clears entries of the matrix
238  */
239  virtual PetscErrorCode mat_zero_entries()
240  {
241  OLD_ASSERT( false, "Function mat_zero_entries is not implemented for linsys type %s \n.", typeid(*this).name() );
242  return 0;
243  }
244 
245  /**
246  * Clears entries of the right-hand side
247  */
248  virtual PetscErrorCode rhs_zero_entries()
249  {
250  OLD_ASSERT( false, "Function vec_zero_entries is not implemented for linsys type %s \n.", typeid(*this).name() );
251  return 0;
252  }
253 
254  /**
255  * Returns PETSC vector with solution. Underlying array can be provided on construction.
256  */
257  const Vec &get_solution()
258  {
259  return solution_;
260  }
261 
262  /**
263  * Create PETSc solution
264  */
265  void set_solution(double *sol_array) {
266  if (sol_array == NULL) {
267  v_solution_ = new double[ rows_ds_->lsize() + 1 ];
268  own_solution_ = true;
269  }
270  else {
271  v_solution_ = sol_array;
272  own_solution_ = false;
273  }
274  PetscErrorCode ierr;
275  ierr = VecCreateMPIWithArray( comm_,1, rows_ds_->lsize(), PETSC_DECIDE, v_solution_, &solution_ ); CHKERRV( ierr );
276  }
277 
278  /**
279  * Returns PETSC subarray with solution. Underlying array can be provided on construction.
280  */
282  {
283  return v_solution_;
284  }
285 
286  /**
287  * Switch linear system into allocating assembly. (only for PETSC_MPIAIJ_preallocate_by_assembly)
288  */
289  virtual void start_allocation()
290  {
291  OLD_ASSERT( false, "Function start_allocation is not implemented for linsys type %s \n.", typeid(*this).name() );
292  }
293 
294  /**
295  * Switch linear system into adding assembly. (the only one supported by triplets ??)
296  */
297  virtual void start_add_assembly()
298  {
299  OLD_ASSERT( false, "Function start_add_assembly is not implemented for linsys type %s \n.", typeid(*this).name() );
300  }
301 
302  /**
303  * Switch linear system into insert assembly. (not currently used)
304  */
305  virtual void start_insert_assembly()
306  {
307  OLD_ASSERT( false, "Function start_insert_assembly is not implemented for linsys type %s \n.", typeid(*this).name() );
308  }
309 
310  /**
311  * Finish assembly of the whole system. For PETSC this should call MatEndAssembly with MAT_FINAL_ASSEMBLY
312  */
313  virtual void finish_assembly( )=0;
314 
315  /**
316  * Assembly full rectangular submatrix into the system matrix.
317  * Should be virtual, implemented differently in particular solvers.
318  */
319  virtual void mat_set_values(int nrow,int *rows,int ncol,int *cols,double *vals)=0;
320 
321  /**
322  * Shortcut for assembling just one element into the matrix.
323  * Similarly we can provide method accepting armadillo matrices.
324  */
325  void mat_set_value(int row,int col,double val)
326  { mat_set_values(1,&row,1,&col,&val); }
327 
328  /**
329  * Set values of the system right-hand side.
330  * Should be virtual, implemented differently in particular solvers.
331  */
332  virtual void rhs_set_values(int nrow,int *rows,double *vals)=0;
333 
334  /**
335  * Shorcut for assembling just one element into RHS vector.
336  */
337  void rhs_set_value(int row,double val)
338  { rhs_set_values(1,&row,&val); }
339 
340 
341  /// Set values in the system matrix and values in the right-hand side vector on corresponding rows.
342  inline void set_values(int nrow,int *rows,int ncol,int *cols,PetscScalar *mat_vals, PetscScalar *rhs_vals)
343  {
344  mat_set_values(nrow, rows, ncol, cols, mat_vals);
345  rhs_set_values(nrow, rows, rhs_vals);
346  }
347 
348  /**
349  * Shortcut to assembly into matrix and RHS in one call, possibly apply Dirichlet boundary conditions.
350  * @p row_dofs - are global indices of rows of dense @p matrix and rows of dense vector @rhs in global system
351  * @p col_dofs - are global indices of columns of the matrix, and possibly
352  *
353  * Application of Dirichlet conditions:
354  * 1) Rows with negative dofs are set to zero.
355  * 2) Cols with negative dofs are eliminated.
356  * 3) If there are entries on global diagonal. We determine value K either from diagonal of local matrix, or (if it is zero) from
357  * diagonal average.
358  *
359  * Caveats:
360  * - can not set dirichlet condition on zero dof
361  * - Armadillo stores matrix in column first form (Fortran like) which makes it not well suited
362  * for passing local matrices.
363  *
364  */
365  void set_values(std::vector<int> &row_dofs, std::vector<int> &col_dofs,
366  const arma::mat &matrix, const arma::vec &rhs,
367  const arma::vec &row_solution, const arma::vec &col_solution)
368 
369  {
370  arma::mat tmp = matrix.t();
371  arma::vec tmp_rhs = rhs;
372  bool negative_row = false;
373  bool negative_col = false;
374 
375  for(unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
376  if (row_dofs[l_row] < 0) {
377  tmp_rhs(l_row)=0.0;
378  tmp.col(l_row).zeros();
379  negative_row=true;
380  }
381 
382  for(unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
383  if (col_dofs[l_col] < 0) {
384  tmp_rhs -= matrix.col(l_col) * col_solution[l_col];
385  tmp.row(l_col).zeros();
386  negative_col=true;
387  }
388 
389 
390  if (negative_row && negative_col) {
391  // look for diagonal entry
392  for(unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
393  if (row_dofs[l_row] < 0)
394  for(unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
395  if (col_dofs[l_col] < 0 && row_dofs[l_row] == col_dofs[l_col]) {
396  double new_diagonal = fabs(matrix.at(l_row,l_col));
397  if (new_diagonal == 0.0) {
398  if (matrix.is_square()) {
399  new_diagonal = arma::sum( abs(matrix.diag())) / matrix.n_rows;
400  } else {
401  new_diagonal = arma::accu( abs(matrix) ) / matrix.n_elem;
402  }
403  }
404  tmp.at(l_col, l_row) = new_diagonal;
405  tmp_rhs(l_row) = new_diagonal * row_solution[l_row];
406  }
407 
408  }
409 
410  if (negative_row)
411  for(int &row : row_dofs) row=abs(row);
412 
413  if (negative_col)
414  for(int &col : col_dofs) col=abs(col);
415 
416 
417  mat_set_values(row_dofs.size(), const_cast<int *>(&(row_dofs[0])),
418  col_dofs.size(), const_cast<int *>(&(col_dofs[0])), tmp.memptr() );
419  rhs_set_values(row_dofs.size(), const_cast<int *>(&(row_dofs[0])), tmp_rhs.memptr() );
420  }
421 
422  /**
423  * Adds Dirichlet constrain.
424  * @param row - global number of row that should be eliminated.
425  * @param value - solution value at the given row
426  */
427  void add_constraint(int row, double value) {
428 
429  constraints_.push_back( Constraint_( static_cast<unsigned>( row ), value ) );
430  }
431 
432  /**
433  * Apply constrains to assembled matrix. Constrains are given by pairs: global row index, value.
434  * i.e. typedef pair<unsigned int, double> Constrain;
435  *
436  * What is th meaning of ( const double factor ) form Cambridge code?
437  */
438  virtual void apply_constrains( double scalar )=0;
439 
440  /**
441  * Solve the system and return convergence reason.
442  */
443  virtual int solve()=0;
444 
445  /**
446  * Returns norm of the initial right hand side
447  */
449  return residual_norm_;
450  };
451 
452  /**
453  * Explicitly compute residual and its norm for current solution.
454  */
455  virtual double compute_residual() =0;
456 
457  /**
458  * Returns information on relative solver accuracy
459  */
461  return r_tol_;
462  };
463 
464  /**
465  * Returns information on absolute solver accuracy
466  */
467  virtual double get_absolute_accuracy(){
468  return 0.0;
469  };
470 
471  /**
472  * Provides user knowledge about symmetry.
473  */
474  inline void set_symmetric(bool flag = true)
475  {
476  symmetric_ = flag;
477  if (!flag) {
478  set_positive_definite(false);
479  set_negative_definite(false);
480  }
481  }
482 
483  inline bool is_symmetric()
484  { return symmetric_; }
485 
486  /**
487  * Provides user knowledge about positive definiteness.
488  */
489  inline void set_positive_definite(bool flag = true)
490  {
491  positive_definite_ = flag;
492  if (flag) {
493  set_symmetric();
494  set_negative_definite(false);
495  }
496  }
497 
498  /**
499  * Provides user knowledge about negative definiteness.
500  */
501  inline void set_negative_definite(bool flag = true)
502  {
503  negative_definite_ = flag;
504  if (flag) {
505  set_symmetric();
506  set_positive_definite(false);
507  }
508  }
509 
510  inline bool is_positive_definite()
511  { return positive_definite_; }
512 
513  inline bool is_negative_definite()
514  { return negative_definite_; }
515 
516  /// TODO: In fact we want to know if the matrix is already preallocated
517  /// However to do this we need explicit finalisation of preallocating cycle.
518  inline bool is_new() {
519  return ( status_ == NONE );
520  }
521 
522  inline bool is_preallocated() {
523  return ( status_ == INSERT || status_ == ADD);
524  }
525 
526  /**
527  * Provides user knowledge about positive definiteness via symmetric general approach.
528  * This is useful for solving Darcy flow by mixed hybrid method, where blocks on subdomains are saddle point but
529  * interface among subdomains is only at the block of Lagrange multipliers and is symmetric positive definite.
530  * Problem condensed to interface can thus be solved by PCG method, although original problem is saddle point.
531  */
532  inline void set_spd_via_symmetric_general(bool flag = true)
533  {
535  if (flag) set_symmetric();
536  }
537 
539  { return spd_via_symmetric_general_; }
540 
541 
542  /**
543  * Output the system in the Matlab format possibly with given ordering.
544  * Rather we shoud provide output operator <<, since it is more flexible.
545  */
546  virtual void view()
547  {
548  OLD_ASSERT( false, "Function view is not implemented for linsys type %s \n.", typeid(*this).name() );
549  }
550 
551  /**
552  * Sets basic parameters of LinSys defined by user in input file and used to calculate
553  */
554  virtual void set_from_input(const Input::Record in_rec)
555  {
556  in_rec_ = in_rec;
558  }
559 
560  /**
561  * Get precision of solving
562  */
563  virtual double get_solution_precision() = 0;
564 
565  virtual ~LinSys()
566  {
567  PetscErrorCode ierr;
568  if ( solution_ ) { ierr = VecDestroy(&solution_); CHKERRV( ierr ); }
569  if ( own_solution_ ) delete[] v_solution_;
570  }
571 
572 protected:
573  // Default values initialized in constructor
574  static constexpr double default_r_tol_ = 1e-7;
575  static constexpr double default_a_tol_ = 1e-11;
576  static constexpr unsigned int default_max_it_ = 1000;
577 
578  double r_tol_; ///< relative tolerance of linear solver
579  double a_tol_; ///< absolute tolerance of linear solver
580  unsigned int max_it_; ///< maximum number of iterations of linear solver
581 
583  SetValuesMode status_; //!< Set value status of the linear system.
584 
585  const unsigned lsize_; //!< local number of matrix rows (non-overlapping division of rows)
586  unsigned size_; //!< global number of matrix rows, i.e. problem size
587 
588  const Distribution * rows_ds_; //!< final distribution of rows of MH matrix
589 
594 
595  bool matrix_changed_; //!< true if the matrix was changed since the last solve
596  bool rhs_changed_; //!< true if the right hand side was changed since the last solve
597 
598  Vec solution_; //!< PETSc vector constructed with vb array.
599  double *v_solution_; //!< local solution array pointing into Vec solution_
600  bool own_solution_; //!< Indicates if the solution array has been allocated by this class
601 
602  double residual_norm_; //!< local solution array pointing into Vec solution_
603 
604  ConstraintVec_ constraints_;
605 
606  std::vector<double> globalSolution_; //!< global solution in numbering for linear system
607 
608  Input::Record in_rec_; // structure contained parameters of LinSys defined in input file
609 
610 };
611 
612 #endif /* LA_LINSYS_HH_ */
unsigned int size()
Definition: linsys.hh:138
std::vector< double > globalSolution_
global solution in numbering for linear system
Definition: linsys.hh:606
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:583
bool is_new()
Definition: linsys.hh:518
virtual void start_insert_assembly()
Definition: linsys.hh:305
bool is_negative_definite()
Definition: linsys.hh:513
bool spd_via_symmetric_general_
Definition: linsys.hh:593
std::pair< unsigned, double > Constraint_
Definition: linsys.hh:89
double get_residual_norm()
Definition: linsys.hh:448
void set_symmetric(bool flag=true)
Definition: linsys.hh:474
double * v_solution_
local solution array pointing into Vec solution_
Definition: linsys.hh:599
int MPI_Comm
Definition: mpi.h:141
virtual void view()
Definition: linsys.hh:546
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:554
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:595
virtual void start_add_assembly()
Definition: linsys.hh:297
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:239
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:365
const unsigned lsize_
local number of matrix rows (non-overlapping division of rows)
Definition: linsys.hh:585
void set_rhs_changed()
Definition: linsys.hh:193
static constexpr double default_r_tol_
Definition: linsys.hh:574
bool symmetric_
Definition: linsys.hh:590
unsigned int vec_lsize()
Definition: linsys.hh:147
ConstraintVec_ constraints_
Definition: linsys.hh:604
#define MPI_SUM
Definition: mpi.h:196
bool is_matrix_changed()
Definition: linsys.hh:208
void set_negative_definite(bool flag=true)
Definition: linsys.hh:501
virtual void start_allocation()
Definition: linsys.hh:289
virtual void finish_assembly()=0
virtual double get_solution_precision()=0
void add_constraint(int row, double value)
Definition: linsys.hh:427
void set_spd_via_symmetric_general(bool flag=true)
Definition: linsys.hh:532
unsigned size_
global number of matrix rows, i.e. problem size
Definition: linsys.hh:586
Input::Record in_rec_
Definition: linsys.hh:608
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:596
LinSys(const Distribution *rows_ds)
Definition: linsys.hh:103
virtual ~LinSys()
Definition: linsys.hh:565
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
virtual void apply_constrains(double scalar)=0
static constexpr bool value
Definition: json.hpp:87
bool is_preallocated()
Definition: linsys.hh:522
void mat_set_value(int row, int col, double val)
Definition: linsys.hh:325
#define OLD_ASSERT(...)
Definition: global_defs.h:131
Global macros to enhance readability and debugging, general constants.
virtual PetscErrorCode set_rhs(Vec &rhs)
Definition: linsys.hh:230
MPI_Comm comm_
Definition: linsys.hh:582
double get_relative_accuracy()
Definition: linsys.hh:460
unsigned int max_it_
maximum number of iterations of linear solver
Definition: linsys.hh:580
double * get_solution_array()
Definition: linsys.hh:281
const Vec & get_solution()
Definition: linsys.hh:257
static constexpr unsigned int default_max_it_
Definition: linsys.hh:576
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:588
double a_tol_
absolute tolerance of linear solver
Definition: linsys.hh:579
SetValuesMode
Definition: linsys.hh:80
static constexpr double default_a_tol_
Definition: linsys.hh:575
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:602
bool is_rhs_changed()
Definition: linsys.hh:214
bool negative_definite_
Definition: linsys.hh:592
void set_solution(double *sol_array)
Definition: linsys.hh:265
Class for declaration of polymorphic Record.
virtual const Vec * get_rhs()
Definition: linsys.hh:178
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:248
bool is_symmetric()
Definition: linsys.hh:483
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:342
double r_tol_
relative tolerance of linear solver
Definition: linsys.hh:578
#define MPI_INT
Definition: mpi.h:160
Abstract linear system class.
bool is_positive_definite()
Definition: linsys.hh:510
void set_matrix_changed()
Definition: linsys.hh:187
void set_positive_definite(bool flag=true)
Definition: linsys.hh:489
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
Definition: linsys.hh:221
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:591
bool own_solution_
Indicates if the solution array has been allocated by this class.
Definition: linsys.hh:600
virtual const Mat * get_matrix()
Definition: linsys.hh:162
virtual double get_absolute_accuracy()
Definition: linsys.hh:467
virtual double compute_residual()=0
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:598
LinSys(LinSys &other)
Definition: linsys.hh:122
bool is_spd_via_symmetric_general()
Definition: linsys.hh:538
virtual int solve()=0
void rhs_set_value(int row, double val)
Definition: linsys.hh:337
unsigned int lsize(int proc) const
get local size