Flow123d  jenkins-Flow123d-windows-release-multijob-285
linsys_PETSC.cc
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: la_linsys.hh 1299 2011-08-23 21:42:50Z jakub.sistek $
21  * $Revision: 1299 $
22  * $LastChangedBy: jakub.sistek $
23  * $LastChangedDate: 2011-08-23 23:42:50 +0200 (Tue, 23 Aug 2011) $
24  *
25  * @file
26  * @brief Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construction
27  * @author Jakub Sistek
28  *
29  *
30  */
31 
32 // derived from base linsys
33 #include "la/linsys_PETSC.hh"
34 #include "petscvec.h"
35 #include "petscksp.h"
36 #include "petscmat.h"
37 #include "system/sys_profiler.hh"
38 
39 
40 //#include <boost/bind.hpp>
41 
42 namespace it = Input::Type;
43 
44 it::Record LinSys_PETSC::input_type = it::Record("Petsc", "Solver setting.")
46  .declare_key("a_tol", it::Double(0.0), it::Default("1.0e-9"), "Absolute residual tolerance.")
47  .declare_key("options", it::String(), it::Default(""), "Options passed to PETSC before creating KSP instead of default setting.");
48 
49 
51  : LinSys( rows_ds ),
52  init_guess_nonzero(false),
53  matrix_(0)
54 {
55  // create PETSC vectors:
56  PetscErrorCode ierr;
57  // rhs
58  v_rhs_= new double[ rows_ds_->lsize() + 1 ];
59  ierr = VecCreateMPIWithArray( comm_, 1, rows_ds_->lsize(), PETSC_DECIDE, v_rhs_, &rhs_ ); CHKERRV( ierr );
60  ierr = VecZeroEntries( rhs_ ); CHKERRV( ierr );
61 
62  params_ = "";
63  matrix_ = NULL;
64  solution_precision_ = std::numeric_limits<double>::infinity();
65  matrix_changed_ = true;
66  rhs_changed_ = true;
67 }
68 
70  : LinSys(other), params_(other.params_), v_rhs_(NULL), solution_precision_(other.solution_precision_)
71 {
72  MatCopy(other.matrix_, matrix_, DIFFERENT_NONZERO_PATTERN);
73  VecCopy(other.rhs_, rhs_);
74  VecCopy(other.on_vec_, on_vec_);
75  VecCopy(other.off_vec_, off_vec_);
76 }
77 
79 {
80  PetscErrorCode ierr;
81 
82  ierr = VecCreateMPI( comm_, rows_ds_->lsize(), PETSC_DECIDE, &(on_vec_) ); CHKERRV( ierr );
83  ierr = VecDuplicate( on_vec_, &(off_vec_) ); CHKERRV( ierr );
84  status_ = ALLOCATE;
85 }
86 
88 {
89  switch ( status_ ) {
90  case ALLOCATE:
91  this->preallocate_matrix( );
92  break;
93  case INSERT:
94  this->finish_assembly( MAT_FLUSH_ASSEMBLY );
95  break;
96  case ADD:
97  case DONE:
98  break;
99  default:
100  ASSERT( false, "Can not set values. Matrix is not preallocated.\n");
101  }
102  status_ = ADD;
103 }
104 
106 {
107  switch ( status_ ) {
108  case ALLOCATE:
109  this->preallocate_matrix();
110  break;
111  case ADD:
112  this->finish_assembly( MAT_FLUSH_ASSEMBLY );
113  break;
114  case INSERT:
115  case DONE:
116  break;
117  default:
118  ASSERT( false, "Can not set values. Matrix is not preallocated.\n");
119  }
120  status_ = INSERT;
121 }
122 
123 void LinSys_PETSC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
124 {
125  PetscErrorCode ierr;
126 
127  // here vals would need to be converted from double to PetscScalar if it was ever something else than double :-)
128  switch (status_) {
129  case INSERT:
130  case ADD:
131  ierr = MatSetValues(matrix_,nrow,rows,ncol,cols,vals,(InsertMode)status_); CHKERRV( ierr );
132  break;
133  case ALLOCATE:
134  this->preallocate_values(nrow,rows,ncol,cols);
135  break;
136  default: DBGMSG("LS SetValues with non allowed insert mode.\n");
137  }
138 
139  matrix_changed_ = true;
140 }
141 
142 void LinSys_PETSC::rhs_set_values( int nrow, int *rows, double *vals )
143 {
144  PetscErrorCode ierr;
145 
146  switch (status_) {
147  case INSERT:
148  case ADD:
149  ierr = VecSetValues(rhs_,nrow,rows,vals,(InsertMode)status_); CHKERRV( ierr );
150  break;
151  case ALLOCATE:
152  break;
153  default: ASSERT(false, "LinSys's status disallow set values.\n");
154  }
155 
156  rhs_changed_ = true;
157 }
158 
159 void LinSys_PETSC::preallocate_values(int nrow,int *rows,int ncol,int *cols)
160 {
161  int i,j;
162  int col;
163  PetscInt row;
164 
165  for (i=0; i<nrow; i++) {
166  row=rows[i];
167  for(j=0; j<ncol; j++) {
168  col = cols[j];
169  if (rows_ds_->get_proc(row) == rows_ds_->get_proc(col))
170  VecSetValue(on_vec_,row,1.0,ADD_VALUES);
171  else
172  VecSetValue(off_vec_,row,1.0,ADD_VALUES);
173  }
174  }
175 }
176 
178 {
179  ASSERT(status_ == ALLOCATE, "Linear system has to be in ALLOCATE status.");
180 
181  PetscErrorCode ierr;
182  PetscInt *on_nz, *off_nz;
183  PetscScalar *on_array, *off_array;
184 
185  // assembly and get values from counting vectors, destroy them
186  VecAssemblyBegin(on_vec_);
187  VecAssemblyBegin(off_vec_);
188 
189  on_nz = new PetscInt[ rows_ds_->lsize() ];
190  off_nz = new PetscInt[ rows_ds_->lsize() ];
191 
192  VecAssemblyEnd(on_vec_);
193  VecAssemblyEnd(off_vec_);
194 
195  VecGetArray( on_vec_, &on_array );
196  VecGetArray( off_vec_, &off_array );
197 
198  for ( unsigned int i=0; i<rows_ds_->lsize(); i++ ) {
199  on_nz[i] = static_cast<PetscInt>( on_array[i]+0.1 ); // small fraction to ensure correct rounding
200  off_nz[i] = static_cast<PetscInt>( off_array[i]+0.1 );
201  }
202 
203  VecRestoreArray(on_vec_,&on_array);
204  VecRestoreArray(off_vec_,&off_array);
205  VecDestroy(&on_vec_);
206  VecDestroy(&off_vec_);
207 
208  // create PETSC matrix with preallocation
209  if (matrix_ != NULL)
210  {
211  ierr = MatDestroy(&matrix_); CHKERRV( ierr );
212  }
213  ierr = MatCreateAIJ(PETSC_COMM_WORLD, rows_ds_->lsize(), rows_ds_->lsize(), PETSC_DETERMINE, PETSC_DETERMINE,
214  0, on_nz, 0, off_nz, &matrix_); CHKERRV( ierr );
215 
216  if (symmetric_) MatSetOption(matrix_, MAT_SYMMETRIC, PETSC_TRUE);
217  MatSetOption(matrix_, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
218 
219  delete[] on_nz;
220  delete[] off_nz;
221 }
222 
224 {
225  MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
226  this->finish_assembly( assemblyType );
227 }
228 
229 void LinSys_PETSC::finish_assembly( MatAssemblyType assembly_type )
230 {
231  PetscErrorCode ierr;
232 
233  if (status_ == ALLOCATE) {
234  xprintf(Warn, "Finalizing linear system without setting values.\n");
235  this->preallocate_matrix();
236  }
237  ierr = MatAssemblyBegin(matrix_, assembly_type); CHKERRV( ierr );
238  ierr = VecAssemblyBegin(rhs_); CHKERRV( ierr );
239  ierr = MatAssemblyEnd(matrix_, assembly_type); CHKERRV( ierr );
240  ierr = VecAssemblyEnd(rhs_); CHKERRV( ierr );
241 
242  if (assembly_type == MAT_FINAL_ASSEMBLY) status_ = DONE;
243 
244  matrix_changed_ = true;
245  rhs_changed_ = true;
246 }
247 
248 void LinSys_PETSC::apply_constrains( double scalar )
249 {
250  PetscErrorCode ierr;
251 
252  // check that system matrix is assembled
253  ASSERT ( status_ == DONE, "System matrix and right-hand side are not assembled when applying constraints." );
254 
255  // number of constraints
256  PetscInt numConstraints = static_cast<PetscInt>( constraints_.size() );
257 
258  // Additional multiplier for numerical reasons (criterion to be established)
259  const PetscScalar diagScalar = static_cast<PetscScalar>( scalar );
260 
261  std::vector<PetscInt> globalDofs;
263 
264  // Constraint iterators
265  ConstraintVec_::const_iterator cIter = constraints_.begin( );
266  ConstraintVec_::const_iterator cEnd = constraints_.end( );
267  // collect global dof indices and the correpsonding values
268  for ( ; cIter != cEnd; ++cIter ) {
269  globalDofs.push_back( static_cast<PetscInt>( cIter -> first ) );
270  values.push_back( static_cast<PetscScalar>( cIter -> second ) * diagScalar );
271  }
272 
273  // prepare pointers to be passed to PETSc
274  PetscInt * globalDofPtr = this->makePetscPointer_( globalDofs );
275  PetscScalar * valuePtr = this->makePetscPointer_( values );
276 
277  // set matrix rows to zero
278  ierr = MatZeroRows( matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
279  matrix_changed_ = true;
280 
281  // set RHS entries to values (crashes if called with NULL pointers)
282  if ( numConstraints ) {
283  ierr = VecSetValues( rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
284  rhs_changed_ = true;
285  }
286 
287  // perform communication in the rhs vector
288  ierr = VecAssemblyBegin( rhs_ ); CHKERRV( ierr );
289  ierr = VecAssemblyEnd( rhs_ ); CHKERRV( ierr );
290 }
291 
292 
294 {
295  init_guess_nonzero = set_nonzero;
296 }
297 
298 
300 {
301  KSP system;
302  KSPConvergedReason reason;
303 
304  const char *petsc_dflt_opt;
305  int nits;
306 
307  // -mat_no_inode ... inodes are usefull only for
308  // vector problems e.g. MH without Schur complement reduction
309  if (rows_ds_->np() > 1) {
310  // parallel setting
311  if (this->is_positive_definite())
312  //petsc_dflt_opt="-ksp_type cg -ksp_diagonal_scale_fix -pc_type asm -pc_asm_overlap 4 -sub_pc_type icc -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
313  petsc_dflt_opt="-ksp_type bcgs -ksp_diagonal_scale_fix -pc_type asm -pc_asm_overlap 4 -sub_pc_type ilu -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
314  else
315  petsc_dflt_opt="-ksp_type bcgs -ksp_diagonal_scale_fix -pc_type asm -pc_asm_overlap 4 -sub_pc_type ilu -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
316 
317  }
318  else {
319  // serial setting
320  if (this->is_positive_definite())
321  //petsc_dflt_opt="-ksp_type cg -pc_type icc -pc_factor_levels 3 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
322  petsc_dflt_opt="-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
323  else
324  petsc_dflt_opt="-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
325  }
326 
327  if (params_ == "") params_ = petsc_dflt_opt;
328  xprintf(MsgLog,"inserting petsc options: %s\n",params_.c_str());
329  PetscOptionsInsertString(params_.c_str()); // overwrites previous options values
330 
331  MatSetOption( matrix_, MAT_USE_INODES, PETSC_FALSE );
332 
333  KSPCreate( comm_, &system );
334  KSPSetOperators(system, matrix_, matrix_, DIFFERENT_NONZERO_PATTERN);
335 
336  // TODO take care of tolerances - shall we support both input file and command line petsc setting
337  KSPSetTolerances(system, r_tol_, a_tol_, PETSC_DEFAULT,PETSC_DEFAULT);
338  KSPSetFromOptions(system);
339  // We set the KSP flag set_initial_guess_nonzero
340  // unless KSP type is preonly.
341  // In such case PETSc fails (version 3.4.1)
342  if (init_guess_nonzero)
343  {
344  KSPType type;
345  KSPGetType(system, &type);
346  if (strcmp(type, KSPPREONLY) != 0)
347  KSPSetInitialGuessNonzero(system, PETSC_TRUE);
348  }
349 
350  {
351  START_TIMER("PETSC linear solver");
352  START_TIMER("PETSC linear iteration");
353  KSPSolve(system, rhs_, solution_ );
354  KSPGetConvergedReason(system,&reason);
355  KSPGetIterationNumber(system,&nits);
356  ADD_CALLS(nits);
357  }
358  // substitute by PETSc call for residual
359  VecNorm(rhs_, NORM_2, &residual_norm_);
360 
361  xprintf(MsgLog,"convergence reason %d, number of iterations is %d\n", reason, nits);
362 
363  // get residual norm
364  KSPGetResidualNorm(system, &solution_precision_);
365 
366  // TODO: I do not understand this
367  //Profiler::instance()->set_timer_subframes("SOLVING MH SYSTEM", nits);
368 
369  KSPDestroy(&system);
370 
371  return static_cast<int>(reason);
372 
373 }
374 
376 {
377  std::string matFileName = "flow123d_matrix.m";
378  std::string rhsFileName = "flow123d_rhs.m";
379  std::string solFileName = "flow123d_sol.m";
380 
381  PetscViewer myViewer;
382 
383  PetscViewerASCIIOpen(comm_,matFileName.c_str(),&myViewer);
384  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
385  MatView( matrix_, myViewer );
386  PetscViewerDestroy(&myViewer);
387 
388  PetscViewerASCIIOpen(comm_,rhsFileName.c_str(),&myViewer);
389  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
390  VecView( rhs_, myViewer );
391  PetscViewerDestroy(&myViewer);
392 
393  if ( solution_ != NULL ) {
394  PetscViewerASCIIOpen(comm_,solFileName.c_str(),&myViewer);
395  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
396  VecView( solution_, myViewer );
397  PetscViewerDestroy(&myViewer);
398  }
399 }
400 
402 {
403  PetscErrorCode ierr;
404 
405  if (matrix_ != NULL) { ierr = MatDestroy(&matrix_); CHKERRV( ierr ); }
406  ierr = VecDestroy(&rhs_); CHKERRV( ierr );
407 
408  if (v_rhs_ != NULL) delete[] v_rhs_;
409 }
410 
412 {
413  if (! in_rec.is_empty()) {
414  // common values
415  LinSys::set_from_input( in_rec );
416 
417  // PETSc specific setting
418  a_tol_ = in_rec.val<double>("a_tol");
419  params_ = in_rec.val<string>("options");
420  }
421 }
422 
423 
425 {
426  return solution_precision_;
427 }
428 
unsigned int get_proc(unsigned int idx) const
get processor of the given index
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:576
#define DBGMSG(...)
Definition: global_defs.h:196
void start_insert_assembly()
void apply_constrains(double scalar=1.)
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:548
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:588
double * v_rhs_
local RHS array pointing to Vec rhs_
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:39
LinSys_PETSC(const Distribution *rows_ds)
Definition: linsys_PETSC.cc:50
bool symmetric_
Definition: linsys.hh:583
ConstraintVec_ constraints_
Definition: linsys.hh:597
Vec rhs_
PETSc vector constructed with vx array.
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:197
static Input::Type::Record input_type
Definition: linsys_PETSC.hh:46
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
void set_from_input(const Input::Record in_rec)
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
void rhs_set_values(int nrow, int *rows, double *vals)
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:589
Definition: system.hh:72
#define ASSERT(...)
Definition: global_defs.h:121
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:392
MPI_Comm comm_
Definition: linsys.hh:575
Definition: system.hh:72
static Input::Type::AbstractRecord input_type
Definition: linsys.hh:90
T * makePetscPointer_(std::vector< T > &array)
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:100
#define START_TIMER(tag)
Starts a timer with specified tag.
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:581
double a_tol_
Definition: linsys.hh:572
void finish_assembly()
bool is_empty() const
Definition: accessors.hh:401
unsigned int np() const
get num of processors
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:595
Vec on_vec_
Vectors for counting non-zero entries in diagonal block.
void preallocate_matrix()
void start_allocation()
Definition: linsys_PETSC.cc:78
std::string params_
command-line-like options for the PETSc solver
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)
void start_add_assembly()
Definition: linsys_PETSC.cc:87
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
double r_tol_
Definition: linsys.hh:571
double get_solution_precision()
Abstract linear system class.
bool is_positive_definite()
Definition: linsys.hh:504
void set_initial_guess_nonzero(bool set_nonzero=true)
Mat matrix_
Petsc matrix of the problem.
Record type proxy class.
Definition: type_record.hh:169
bool init_guess_nonzero
flag for starting from nonzero guess
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:591
double solution_precision_
unsigned int lsize(int proc) const
get local size
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:430