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