Flow123d  JB_transport-d4c8564
linsys_PETSC.cc
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_PETSC.cc
15  * @brief Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construction
16  * @author Jakub Sistek
17  */
18 
19 // derived from base linsys
20 #include "la/linsys_PETSC.hh"
21 #include "petscvec.h"
22 #include "petscksp.h"
23 #include "petscmat.h"
24 #include "system/sys_profiler.hh"
25 #include "system/system.hh"
26 
27 
28 //#include <boost/bind.hpp>
29 
30 namespace it = Input::Type;
31 
33  return it::Record("Petsc", "PETSc solver settings.\n It provides interface to various PETSc solvers. The convergence criteria is:\n"
34  "```\n"
35  "norm( res_i ) < max( norm( res_0 ) * r_tol, a_tol )\n"
36  "```\n"
37  "where ```res_i``` is the residuum vector after i-th iteration of the solver and ```res_0``` is the estimate of the norm of the initial residual. "
38  "If the initial guess of the solution is provided (usually only for transient equations) the residual of this estimate is used, "
39  "otherwise the norm of preconditioned RHS is used. "
40  "The default norm is (($L_2$)) norm of preconditioned residual: (($ P^{-1}(Ax-b)$)), usage of other norm may be prescribed using the 'option' key. "
41  "See also PETSc documentation for KSPSetNormType.")
43  .declare_key("r_tol", it::Double(0.0, 1.0), it::Default::read_time("Default value is set by the nonlinear solver or the equation. "
44  "If not, we use the value 1.0e-7."),
45  "Residual tolerance relative to the initial error.")
46  .declare_key("a_tol", it::Double(0.0), it::Default::read_time("Default value is set by the nonlinear solver or the equation. "
47  "If not, we use the value 1.0e-11."),
48  "Absolute residual tolerance.")
49  .declare_key("max_it", it::Integer(0), it::Default::read_time("Default value is set by the nonlinear solver or the equation. "
50  "If not, we use the value 1000."),
51  "Maximum number of outer iterations of the linear solver.")
52  .declare_key("options", it::String(), it::Default("\"\""), "This options is passed to PETSC to create a particular KSP (Krylov space method).\n"
53  "If the string is left empty (by default), the internal default options is used.")
54  .close();
55 }
56 
57 
59 
60 
61 LinSys_PETSC::LinSys_PETSC( const Distribution * rows_ds, const std::string &params)
62  : LinSys( rows_ds ),
63  params_(params),
64  init_guess_nonzero(false),
65  matrix_(0)
66 {
67  // create PETSC vectors:
68  PetscErrorCode ierr;
69  // rhs
70  v_rhs_= new double[ rows_ds_->lsize() + 1 ];
71  ierr = VecCreateMPIWithArray( comm_, 1, rows_ds_->lsize(), PETSC_DECIDE, v_rhs_, &rhs_ ); CHKERRV( ierr );
72  ierr = VecZeroEntries( rhs_ ); CHKERRV( ierr );
73  VecDuplicate(rhs_, &residual_);
74 
75  matrix_ = NULL;
76  solution_precision_ = std::numeric_limits<double>::infinity();
77  matrix_changed_ = true;
78  rhs_changed_ = true;
79 }
80 
82  : LinSys(other), params_(other.params_), v_rhs_(NULL), solution_precision_(other.solution_precision_)
83 {
84  MatCopy(other.matrix_, matrix_, DIFFERENT_NONZERO_PATTERN);
85  VecCopy(other.rhs_, rhs_);
86  VecCopy(other.on_vec_, on_vec_);
87  VecCopy(other.off_vec_, off_vec_);
88 }
89 
90 void LinSys_PETSC::set_tolerances(double r_tol, double a_tol, unsigned int max_it)
91 {
92  if (! in_rec_.is_empty()) {
93  // input record is set
94  r_tol_ = in_rec_.val<double>("r_tol", r_tol);
95  a_tol_ = in_rec_.val<double>("a_tol", a_tol);
96  max_it_ = in_rec_.val<unsigned int>("max_it", max_it);
97  } else {
98  r_tol_ = r_tol;
99  a_tol_ = a_tol;
100  max_it_ = max_it;
101 
102  }
103 }
104 
105 
107 {
108  PetscErrorCode ierr;
109 
110  ierr = VecCreateMPI( comm_, rows_ds_->lsize(), PETSC_DECIDE, &(on_vec_) ); CHKERRV( ierr );
111  ierr = VecDuplicate( on_vec_, &(off_vec_) ); CHKERRV( ierr );
112  status_ = ALLOCATE;
113 }
114 
116 {
117  switch ( status_ ) {
118  case ALLOCATE:
119  this->preallocate_matrix( );
120  break;
121  case INSERT:
122  this->finish_assembly( MAT_FLUSH_ASSEMBLY );
123  break;
124  case ADD:
125  case DONE:
126  break;
127  default:
128  ASSERT_PERMANENT(false).error("Can not set values. Matrix is not preallocated.\n");
129  }
130  status_ = ADD;
131 }
132 
134 {
135  switch ( status_ ) {
136  case ALLOCATE:
137  this->preallocate_matrix();
138  break;
139  case ADD:
140  this->finish_assembly( MAT_FLUSH_ASSEMBLY );
141  break;
142  case INSERT:
143  case DONE:
144  break;
145  default:
146  ASSERT_PERMANENT(false).error("Can not set values. Matrix is not preallocated.\n");
147  }
148  status_ = INSERT;
149 }
150 
151 void LinSys_PETSC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
152 {
153  // here vals would need to be converted from double to PetscScalar if it was ever something else than double :-)
154  switch (status_) {
155  case INSERT:
156  case ADD:
157  chkerr(MatSetValues(matrix_,nrow,rows,ncol,cols,vals,(InsertMode)status_));
158  break;
159  case ALLOCATE:
160  this->preallocate_values(nrow,rows,ncol,cols);
161  break;
162  default: DebugOut() << "LS SetValues with non allowed insert mode.\n";
163  }
164 
165  matrix_changed_ = true;
166 }
167 
168 void LinSys_PETSC::rhs_set_values( int nrow, int *rows, double *vals )
169 {
170  PetscErrorCode ierr;
171 
172  switch (status_) {
173  case INSERT:
174  case ADD:
175  ierr = VecSetValues(rhs_,nrow,rows,vals,(InsertMode)status_); CHKERRV( ierr );
176  break;
177  case ALLOCATE:
178  break;
179  default: ASSERT_PERMANENT(false).error("LinSys's status disallow set values.\n");
180  }
181 
182  rhs_changed_ = true;
183 }
184 
185 void LinSys_PETSC::preallocate_values(int nrow,int *rows,int ncol,int *cols)
186 {
187  int i,j;
188  int col;
189  PetscInt row;
190 
191  for (i=0; i<nrow; i++) {
192  row=rows[i];
193  for(j=0; j<ncol; j++) {
194  col = cols[j];
195  if (rows_ds_->get_proc(row) == rows_ds_->get_proc(col))
196  VecSetValue(on_vec_,row,1.0,ADD_VALUES);
197  else
198  VecSetValue(off_vec_,row,1.0,ADD_VALUES);
199  }
200  }
201 }
202 
204 {
205  ASSERT_EQ(status_, ALLOCATE).error("Linear system has to be in ALLOCATE status.");
206 
207  PetscErrorCode ierr;
208  PetscInt *on_nz, *off_nz;
209  PetscScalar *on_array, *off_array;
210 
211  // assembly and get values from counting vectors, destroy them
212  VecAssemblyBegin(on_vec_);
213  VecAssemblyBegin(off_vec_);
214 
215  on_nz = new PetscInt[ rows_ds_->lsize() ];
216  off_nz = new PetscInt[ rows_ds_->lsize() ];
217 
218  VecAssemblyEnd(on_vec_);
219  VecAssemblyEnd(off_vec_);
220 
221  VecGetArray( on_vec_, &on_array );
222  VecGetArray( off_vec_, &off_array );
223 
224  for ( unsigned int i=0; i<rows_ds_->lsize(); i++ ) {
225  on_nz[i] = std::min( rows_ds_->lsize(), static_cast<uint>( on_array[i]+0.1 ) ); // small fraction to ensure correct rounding
226  off_nz[i] = std::min( rows_ds_->size() - rows_ds_->lsize(), static_cast<uint>( off_array[i]+0.1 ) );
227  }
228 
229  VecRestoreArray(on_vec_,&on_array);
230  VecRestoreArray(off_vec_,&off_array);
231  VecDestroy(&on_vec_);
232  VecDestroy(&off_vec_);
233 
234  // create PETSC matrix with preallocation
235  if (matrix_ != NULL)
236  {
237  chkerr(MatDestroy(&matrix_));
238  }
239  ierr = MatCreateAIJ(PETSC_COMM_WORLD, rows_ds_->lsize(), rows_ds_->lsize(), PETSC_DETERMINE, PETSC_DETERMINE,
240  0, on_nz, 0, off_nz, &matrix_); CHKERRV( ierr );
241 
242  if (symmetric_) MatSetOption(matrix_, MAT_SYMMETRIC, PETSC_TRUE);
243  MatSetOption(matrix_, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
244 
245  // This option is used in order to assembly larger local matrices with own non-zero structure.
246  // Zero entries are ignored so we must prevent adding exact zeroes.
247  // Add LocalSystem::almost_zero for entries that should not be eliminated.
248  MatSetOption(matrix_, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
249 
250 
251 
252  delete[] on_nz;
253  delete[] off_nz;
254 }
255 
257 {
258  MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
259  this->finish_assembly( assemblyType );
260 }
261 
262 void LinSys_PETSC::finish_assembly( MatAssemblyType assembly_type )
263 {
264  PetscErrorCode ierr;
265 
266  if (status_ == ALLOCATE) {
267  WarningOut() << "Finalizing linear system without setting values.\n";
268  this->preallocate_matrix();
269  }
270  ierr = MatAssemblyBegin(matrix_, assembly_type); CHKERRV( ierr );
271  ierr = VecAssemblyBegin(rhs_); CHKERRV( ierr );
272  ierr = MatAssemblyEnd(matrix_, assembly_type); CHKERRV( ierr );
273  ierr = VecAssemblyEnd(rhs_); CHKERRV( ierr );
274 
275  if (assembly_type == MAT_FINAL_ASSEMBLY) status_ = DONE;
276 
277  //PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_INDEX);
278  //MatView(matrix_, PETSC_VIEWER_STDOUT_SELF);
279  //VecView(rhs_, PETSC_VIEWER_STDOUT_SELF);
280  //this->view();
281 
282  matrix_changed_ = true;
283  rhs_changed_ = true;
284 }
285 
286 void LinSys_PETSC::apply_constrains( double scalar )
287 {
288  PetscErrorCode ierr;
289 
290  // check that system matrix is assembled
291  ASSERT_EQ(status_, DONE).error("System matrix and right-hand side are not assembled when applying constraints." );
292 
293  // number of constraints
294  PetscInt numConstraints = static_cast<PetscInt>( constraints_.size() );
295 
296  // Additional multiplier for numerical reasons (criterion to be established)
297  const PetscScalar diagScalar = static_cast<PetscScalar>( scalar );
298 
299  std::vector<PetscInt> globalDofs;
301 
302  // Constraint iterators
303  ConstraintVec_::const_iterator cIter = constraints_.begin( );
304  ConstraintVec_::const_iterator cEnd = constraints_.end( );
305  // collect global dof indices and the correpsonding values
306  for ( ; cIter != cEnd; ++cIter ) {
307  globalDofs.push_back( static_cast<PetscInt>( cIter -> first ) );
308  values.push_back( static_cast<PetscScalar>( cIter -> second ) * diagScalar );
309  }
310 
311  // prepare pointers to be passed to PETSc
312  PetscInt * globalDofPtr = this->makePetscPointer_( globalDofs );
313  PetscScalar * valuePtr = this->makePetscPointer_( values );
314 
315  // set matrix rows to zero
316  ierr = MatZeroRows( matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
317  matrix_changed_ = true;
318 
319  // set RHS entries to values (crashes if called with NULL pointers)
320  if ( numConstraints ) {
321  ierr = VecSetValues( rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
322  rhs_changed_ = true;
323  }
324 
325  // perform communication in the rhs vector
326  ierr = VecAssemblyBegin( rhs_ ); CHKERRV( ierr );
327  ierr = VecAssemblyEnd( rhs_ ); CHKERRV( ierr );
328 }
329 
330 
332 {
333  init_guess_nonzero = set_nonzero;
334 }
335 
336 
338 {
339 
340  const char *petsc_dflt_opt;
341  int nits;
342 
343  // -mat_no_inode ... inodes are usefull only for
344  // vector problems e.g. MH without Schur complement reduction
345 
346  /* Comment to PETSc options:
347  *
348  * -ksp_diagonal_scale scales the matrix before solution, while -ksp_diagonal_scale_fix just fixes the scaling after solution
349  * -pc_asm_type basic enforces classical Schwartz method, which seems more stable for positive definite systems.
350  * The default 'restricted' probably violates s.p.d. structure, many tests fail.
351  */
352  if (rows_ds_->np() > 1) {
353  // parallel setting
354  if (this->is_positive_definite())
355  petsc_dflt_opt="-ksp_type cg -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_type asm -pc_asm_type basic -pc_asm_overlap 4 -sub_pc_type icc -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
356  //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";
357  else
358  petsc_dflt_opt="-ksp_type bcgs -ksp_diagonal_scale -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";
359 
360  }
361  else {
362  // serial setting
363  if (this->is_positive_definite())
364  petsc_dflt_opt="-ksp_type cg -pc_type icc -pc_factor_levels 3 -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
365  //petsc_dflt_opt="-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
366  else
367  petsc_dflt_opt="-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
368  }
369 
370  if (params_ == "") params_ = petsc_dflt_opt;
371  LogOut().fmt("inserting petsc options: {}\n",params_.c_str());
372 
373  // now takes an optional PetscOptions object as the first argument
374  // value NULL will preserve previous behaviour previous behavior.
375  PetscOptionsInsertString(NULL, params_.c_str()); // overwrites previous options values
376 
377  MatSetOption( matrix_, MAT_USE_INODES, PETSC_FALSE );
378 
379  chkerr(KSPCreate( comm_, &system ));
380  chkerr(KSPSetOperators(system, matrix_, matrix_));
381 
382 
383  // TODO take care of tolerances - shall we support both input file and command line petsc setting
384  chkerr(KSPSetTolerances(system, r_tol_, a_tol_, PETSC_DEFAULT,PETSC_DEFAULT));
385  chkerr(KSPSetTolerances(system, r_tol_, a_tol_, PETSC_DEFAULT, max_it_));
386  KSPSetFromOptions(system);
387  // We set the KSP flag set_initial_guess_nonzero
388  // unless KSP type is preonly.
389  // In such case PETSc fails (version 3.4.1)
390  if (init_guess_nonzero)
391  {
392  KSPType type;
393  KSPGetType(system, &type);
394  if (strcmp(type, KSPPREONLY) != 0)
395  KSPSetInitialGuessNonzero(system, PETSC_TRUE);
396  }
397 
398  {
399  START_TIMER("PETSC linear solver");
400  START_TIMER("PETSC linear iteration");
401  chkerr(KSPSolve(system, rhs_, solution_ ));
402  KSPGetConvergedReason(system,&reason);
403  KSPGetIterationNumber(system,&nits);
404  ADD_CALLS(nits);
405  }
406  // substitute by PETSc call for residual
407  VecNorm(rhs_, NORM_2, &residual_norm_);
408 
409  LogOut().fmt("convergence reason {}, number of iterations is {}\n", reason, nits);
410 
411  // get residual norm
412  KSPGetResidualNorm(system, &solution_precision_);
413 
414  // TODO: I do not understand this
415  //Profiler::instance()->set_timer_subframes("SOLVING MH SYSTEM", nits);
416 
417  chkerr(KSPDestroy(&system));
418 
419  return LinSys::SolveInfo(static_cast<int>(reason), static_cast<int>(nits));
420 
421 }
422 
423 void LinSys_PETSC::view(string text )
424 {
425  FilePath matFileName(text + "_flow123d_matrix.m",FilePath::FileType::output_file);
426  FilePath rhsFileName(text + "_flow123d_rhs.m",FilePath::FileType::output_file);
427  FilePath solFileName(text + "_flow123d_sol.m",FilePath::FileType::output_file);
428 
429  PetscViewer myViewer;
430 
431  if ( matrix_ != NULL ) {
432  PetscViewerASCIIOpen(comm_,((string)matFileName).c_str(),&myViewer);
433  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
434  MatView( matrix_, myViewer );
435  PetscViewerDestroy(&myViewer);
436  }
437  else
438  WarningOut() << "PetscViewer: the matrix of LinSys is not set.\n";
439 
440  if ( rhs_ != NULL ) {
441  PetscViewerASCIIOpen(comm_,((string)rhsFileName).c_str(),&myViewer);
442  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
443  VecView( rhs_, myViewer );
444  PetscViewerDestroy(&myViewer);
445  }
446  else
447  WarningOut() << "PetscViewer: the rhs of LinSys is not set.\n";
448 
449  if ( solution_ != NULL ) {
450  PetscViewerASCIIOpen(comm_,((string)solFileName).c_str(),&myViewer);
451  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
452  VecView( solution_, myViewer );
453  PetscViewerDestroy(&myViewer);
454  }
455  else
456  WarningOut() << "PetscViewer: the solution of LinSys is not set.\n";
457 }
458 
460 {
461  if (matrix_ != NULL) { chkerr(MatDestroy(&matrix_)); }
462  chkerr(VecDestroy(&rhs_));
463 
464  if (residual_ != NULL) chkerr(VecDestroy(&residual_));
465  if (v_rhs_ != NULL) delete[] v_rhs_;
466 }
467 
468 
469 
471 {
472  LinSys::set_from_input( in_rec );
473 
474  // PETSC specific parameters
475  // If parameters are specified in input file, they are used,
476  // otherwise keep settings provided in constructor of LinSys_PETSC.
477  std::string user_params = in_rec.val<string>("options");
478  if (user_params != "") params_ = user_params;
479 }
480 
481 
483 {
484  return solution_precision_;
485 }
486 
487 
489 {
490  MatMult(matrix_, solution_, residual_);
491  VecAXPY(residual_,-1.0, rhs_);
492  double residual_norm;
493  VecNorm(residual_, NORM_2, &residual_norm);
494  return residual_norm;
495 }
Distribution::np
unsigned int np() const
get num of processors
Definition: distribution.hh:105
LinSys_PETSC::system
KSP system
Definition: linsys_PETSC.hh:187
LinSys_PETSC::residual_
Vec residual_
Definition: linsys_PETSC.hh:177
LinSys
Definition: la_linsys_new.hh:169
Distribution::lsize
unsigned int lsize(int proc) const
get local size
Definition: distribution.hh:115
LinSys_PETSC::start_insert_assembly
void start_insert_assembly() override
Definition: linsys_PETSC.cc:133
LinSys::residual_norm_
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:688
LinSys_PETSC::set_initial_guess_nonzero
void set_initial_guess_nonzero(bool set_nonzero=true)
Definition: linsys_PETSC.cc:331
LinSys::symmetric_
bool symmetric_
Definition: linsys.hh:675
LinSys::SolveInfo
Definition: linsys.hh:104
LinSys_PETSC::params_
std::string params_
command-line-like options for the PETSc solver
Definition: linsys_PETSC.hh:171
LinSys_PETSC::init_guess_nonzero
bool init_guess_nonzero
flag for starting from nonzero guess
Definition: linsys_PETSC.hh:173
LinSys_PETSC::LinSys_PETSC
LinSys_PETSC(const Distribution *rows_ds, const std::string &params="")
Definition: linsys_PETSC.cc:61
Input::Type::Integer
Class for declaration of the integral input data.
Definition: type_base.hh:483
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
LinSys_PETSC::apply_constrains
void apply_constrains(double scalar=1.) override
Definition: linsys_PETSC.cc:286
LinSys_PETSC::reason
KSPConvergedReason reason
Definition: linsys_PETSC.hh:188
FilePath
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
Input::Type::Double
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
chkerr
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:142
LinSys_PETSC
Definition: linsys_PETSC.hh:43
std::vector
Definition: doxy_dummy_defs.hh:7
Input::Type::Default::read_time
static Default read_time(const std::string &description)
The factory function to make an default value that will be specified at the time when a key will be r...
Definition: type_record.hh:97
Distribution::size
unsigned int size() const
get global size
Definition: distribution.hh:118
system.hh
LinSys::ALLOCATE
@ ALLOCATE
Definition: linsys.hh:99
ConstraintType::cols
@ cols
LinSys_PETSC::solve
LinSys::SolveInfo solve() override
Definition: linsys_PETSC.cc:337
uint
unsigned int uint
Definition: mh_dofhandler.hh:101
linsys_PETSC.hh
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
LinSys_PETSC::set_tolerances
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: linsys_PETSC.cc:90
LinSys_PETSC::matrix_
Mat matrix_
Petsc matrix of the problem.
Definition: linsys_PETSC.hh:175
LinSys_PETSC::registrar
static const int registrar
Registrar of class to factory.
Definition: linsys_PETSC.hh:150
LinSys::DONE
@ DONE
Definition: linsys.hh:100
LinSys_PETSC::mat_set_values
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
Definition: linsys_PETSC.cc:151
Distribution::get_proc
unsigned int get_proc(unsigned int idx) const
get processor of the given index
Definition: distribution.cc:130
ConstraintType::rows
@ rows
Input::Type::Record::size
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Input::Type::Record::derive_from
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
Distribution
Definition: distribution.hh:50
LinSys::solution_
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:683
LogOut
#define LogOut()
Macro defining 'log' record of log.
Definition: logger.hh:281
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
LinSys::a_tol_
double a_tol_
absolute tolerance of linear solver
Definition: linsys.hh:664
sys_profiler.hh
ASSERT_PERMANENT
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
LinSys::in_rec_
Input::Record in_rec_
Definition: linsys.hh:694
LinSys_PETSC::v_rhs_
double * v_rhs_
local RHS array pointing to Vec rhs_
Definition: linsys_PETSC.hh:179
LinSys_PETSC::finish_assembly
void finish_assembly() override
Definition: linsys_PETSC.cc:256
LinSys_PETSC::solution_precision_
double solution_precision_
Definition: linsys_PETSC.hh:185
LinSys::rhs_changed_
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:681
LinSys_PETSC::get_solution_precision
double get_solution_precision() override
Definition: linsys_PETSC.cc:482
LinSys_PETSC::rhs_set_values
void rhs_set_values(int nrow, int *rows, double *vals) override
Definition: linsys_PETSC.cc:168
LinSys::status_
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:668
LinSys_PETSC::on_vec_
Vec on_vec_
Vectors for counting non-zero entries in diagonal block.
Definition: linsys_PETSC.hh:181
LinSys_PETSC::preallocate_matrix
void preallocate_matrix()
Definition: linsys_PETSC.cc:203
ASSERT_EQ
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
Input::Type::Record::declare_key
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
LinSys::rows_ds_
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:673
Input::Type
Definition: balance.hh:41
LinSys_PETSC::preallocate_values
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
Definition: linsys_PETSC.cc:185
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
LinSys::set_from_input
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:640
LinSys_PETSC::start_add_assembly
void start_add_assembly() override
Definition: linsys_PETSC.cc:115
LinSys::INSERT
@ INSERT
Definition: linsys.hh:97
LinSys::r_tol_
double r_tol_
relative tolerance of linear solver
Definition: linsys.hh:663
LinSys_PETSC::set_from_input
void set_from_input(const Input::Record in_rec) override
Definition: linsys_PETSC.cc:470
Input::Type::String
Class for declaration of the input data that are in string format.
Definition: type_base.hh:582
LinSys::constraints_
ConstraintVec_ constraints_
Definition: linsys.hh:690
LinSys::matrix_changed_
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:680
WarningOut
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
ADD_CALLS
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
Definition: sys_profiler.hh:186
LinSys_PETSC::makePetscPointer_
T * makePetscPointer_(std::vector< T > &array)
Definition: linsys_PETSC.hh:154
LinSys_PETSC::start_allocation
void start_allocation() override
Definition: linsys_PETSC.cc:106
LinSys::get_input_type
static Input::Type::Abstract & get_input_type()
Definition: linsys.cc:29
LinSys_PETSC::off_vec_
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
Definition: linsys_PETSC.hh:182
LinSys_PETSC::view
void view(string text="") override
Definition: linsys_PETSC.cc:423
LinSys::ADD
@ ADD
Definition: linsys.hh:98
LinSys::is_positive_definite
bool is_positive_definite()
Definition: linsys.hh:596
DebugOut
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
LinSys::comm_
MPI_Comm comm_
Definition: linsys.hh:667
START_TIMER
#define START_TIMER(tag)
Starts a timer with specified tag.
Definition: sys_profiler.hh:115
LinSys_PETSC::rhs_
Vec rhs_
PETSc vector constructed with vx array.
Definition: linsys_PETSC.hh:176
Input::Record::is_empty
bool is_empty() const
Definition: accessors.hh:365
LinSys_PETSC::compute_residual
double compute_residual() override
Definition: linsys_PETSC.cc:488
LinSys_PETSC::~LinSys_PETSC
~LinSys_PETSC()
Definition: linsys_PETSC.cc:459
LinSys_PETSC::get_input_type
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
LinSys::max_it_
unsigned int max_it_
maximum number of iterations of linear solver
Definition: linsys.hh:665