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