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