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