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