Flow123d  release_3.0.0-1111-g6aae175
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  OLD_ASSERT( false, "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  OLD_ASSERT( false, "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: OLD_ASSERT(false, "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  OLD_ASSERT(status_ == ALLOCATE, "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  MatSetOption(matrix_, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
245 
246  delete[] on_nz;
247  delete[] off_nz;
248 }
249 
251 {
252  MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
253  this->finish_assembly( assemblyType );
254 }
255 
256 void LinSys_PETSC::finish_assembly( MatAssemblyType assembly_type )
257 {
258  PetscErrorCode ierr;
259 
260  if (status_ == ALLOCATE) {
261  WarningOut() << "Finalizing linear system without setting values.\n";
262  this->preallocate_matrix();
263  }
264  ierr = MatAssemblyBegin(matrix_, assembly_type); CHKERRV( ierr );
265  ierr = VecAssemblyBegin(rhs_); CHKERRV( ierr );
266  ierr = MatAssemblyEnd(matrix_, assembly_type); CHKERRV( ierr );
267  ierr = VecAssemblyEnd(rhs_); CHKERRV( ierr );
268 
269  if (assembly_type == MAT_FINAL_ASSEMBLY) status_ = DONE;
270 
271  //PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_INDEX);
272  //MatView(matrix_, PETSC_VIEWER_STDOUT_SELF);
273  //VecView(rhs_, PETSC_VIEWER_STDOUT_SELF);
274  //this->view();
275 
276  matrix_changed_ = true;
277  rhs_changed_ = true;
278 }
279 
280 void LinSys_PETSC::apply_constrains( double scalar )
281 {
282  PetscErrorCode ierr;
283 
284  // check that system matrix is assembled
285  OLD_ASSERT ( status_ == DONE, "System matrix and right-hand side are not assembled when applying constraints." );
286 
287  // number of constraints
288  PetscInt numConstraints = static_cast<PetscInt>( constraints_.size() );
289 
290  // Additional multiplier for numerical reasons (criterion to be established)
291  const PetscScalar diagScalar = static_cast<PetscScalar>( scalar );
292 
293  std::vector<PetscInt> globalDofs;
295 
296  // Constraint iterators
297  ConstraintVec_::const_iterator cIter = constraints_.begin( );
298  ConstraintVec_::const_iterator cEnd = constraints_.end( );
299  // collect global dof indices and the correpsonding values
300  for ( ; cIter != cEnd; ++cIter ) {
301  globalDofs.push_back( static_cast<PetscInt>( cIter -> first ) );
302  values.push_back( static_cast<PetscScalar>( cIter -> second ) * diagScalar );
303  }
304 
305  // prepare pointers to be passed to PETSc
306  PetscInt * globalDofPtr = this->makePetscPointer_( globalDofs );
307  PetscScalar * valuePtr = this->makePetscPointer_( values );
308 
309  // set matrix rows to zero
310  ierr = MatZeroRows( matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
311  matrix_changed_ = true;
312 
313  // set RHS entries to values (crashes if called with NULL pointers)
314  if ( numConstraints ) {
315  ierr = VecSetValues( rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
316  rhs_changed_ = true;
317  }
318 
319  // perform communication in the rhs vector
320  ierr = VecAssemblyBegin( rhs_ ); CHKERRV( ierr );
321  ierr = VecAssemblyEnd( rhs_ ); CHKERRV( ierr );
322 }
323 
324 
326 {
327  init_guess_nonzero = set_nonzero;
328 }
329 
330 
332 {
333 
334  const char *petsc_dflt_opt;
335  int nits;
336 
337  // -mat_no_inode ... inodes are usefull only for
338  // vector problems e.g. MH without Schur complement reduction
339 
340  /* Comment to PETSc options:
341  *
342  * -ksp_diagonal_scale scales the matrix before solution, while -ksp_diagonal_scale_fix just fixes the scaling after solution
343  * -pc_asm_type basic enforces classical Schwartz method, which seems more stable for positive definite systems.
344  * The default 'restricted' probably violates s.p.d. structure, many tests fail.
345  */
346  if (rows_ds_->np() > 1) {
347  // parallel setting
348  if (this->is_positive_definite())
349  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";
350  //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";
351  else
352  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";
353 
354  }
355  else {
356  // serial setting
357  if (this->is_positive_definite())
358  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";
359  //petsc_dflt_opt="-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
360  else
361  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";
362  }
363 
364  if (params_ == "") params_ = petsc_dflt_opt;
365  LogOut().fmt("inserting petsc options: {}\n",params_.c_str());
366 
367  // now takes an optional PetscOptions object as the first argument
368  // value NULL will preserve previous behaviour previous behavior.
369  PetscOptionsInsertString(NULL, params_.c_str()); // overwrites previous options values
370 
371  MatSetOption( matrix_, MAT_USE_INODES, PETSC_FALSE );
372 
373  chkerr(KSPCreate( comm_, &system ));
374  chkerr(KSPSetOperators(system, matrix_, matrix_));
375 
376 
377  // TODO take care of tolerances - shall we support both input file and command line petsc setting
378  chkerr(KSPSetTolerances(system, r_tol_, a_tol_, PETSC_DEFAULT,PETSC_DEFAULT));
379  chkerr(KSPSetTolerances(system, r_tol_, a_tol_, PETSC_DEFAULT, max_it_));
380  KSPSetFromOptions(system);
381  // We set the KSP flag set_initial_guess_nonzero
382  // unless KSP type is preonly.
383  // In such case PETSc fails (version 3.4.1)
384  if (init_guess_nonzero)
385  {
386  KSPType type;
387  KSPGetType(system, &type);
388  if (strcmp(type, KSPPREONLY) != 0)
389  KSPSetInitialGuessNonzero(system, PETSC_TRUE);
390  }
391 
392  {
393  START_TIMER("PETSC linear solver");
394  START_TIMER("PETSC linear iteration");
395  chkerr(KSPSolve(system, rhs_, solution_ ));
396  KSPGetConvergedReason(system,&reason);
397  KSPGetIterationNumber(system,&nits);
398  ADD_CALLS(nits);
399  }
400  // substitute by PETSc call for residual
401  VecNorm(rhs_, NORM_2, &residual_norm_);
402 
403  LogOut().fmt("convergence reason {}, number of iterations is {}\n", reason, nits);
404 
405  // get residual norm
406  KSPGetResidualNorm(system, &solution_precision_);
407 
408  // TODO: I do not understand this
409  //Profiler::instance()->set_timer_subframes("SOLVING MH SYSTEM", nits);
410 
411  chkerr(KSPDestroy(&system));
412 
413  return LinSys::SolveInfo(static_cast<int>(reason), static_cast<int>(nits));
414 
415 }
416 
418 {
419  std::string matFileName = "flow123d_matrix.m";
420  std::string rhsFileName = "flow123d_rhs.m";
421  std::string solFileName = "flow123d_sol.m";
422 
423  PetscViewer myViewer;
424 
425  PetscViewerASCIIOpen(comm_,matFileName.c_str(),&myViewer);
426  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
427  MatView( matrix_, myViewer );
428  PetscViewerDestroy(&myViewer);
429 
430  PetscViewerASCIIOpen(comm_,rhsFileName.c_str(),&myViewer);
431  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
432  VecView( rhs_, myViewer );
433  PetscViewerDestroy(&myViewer);
434 
435  if ( solution_ != NULL ) {
436  PetscViewerASCIIOpen(comm_,solFileName.c_str(),&myViewer);
437  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
438  VecView( solution_, myViewer );
439  PetscViewerDestroy(&myViewer);
440  }
441 }
442 
444 {
445  if (matrix_ != NULL) { chkerr(MatDestroy(&matrix_)); }
446  chkerr(VecDestroy(&rhs_));
447 
448  if (residual_ != NULL) chkerr(VecDestroy(&residual_));
449  if (v_rhs_ != NULL) delete[] v_rhs_;
450 }
451 
452 
453 
455 {
456  LinSys::set_from_input( in_rec );
457 
458  // PETSC specific parameters
459  // If parameters are specified in input file, they are used,
460  // otherwise keep settings provided in constructor of LinSys_PETSC.
461  std::string user_params = in_rec.val<string>("options");
462  if (user_params != "") params_ = user_params;
463 }
464 
465 
467 {
468  return solution_precision_;
469 }
470 
471 
473 {
474  MatMult(matrix_, solution_, residual_);
475  VecAXPY(residual_,-1.0, rhs_);
476  double residual_norm;
477  VecNorm(residual_, NORM_2, &residual_norm);
478  return residual_norm;
479 }
unsigned int size() const
get global size
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:649
void apply_constrains(double scalar=1.) override
unsigned int uint
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
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:621
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:661
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:61
void set_from_input(const Input::Record in_rec) override
LinSys::SolveInfo solve() override
bool symmetric_
Definition: linsys.hh:656
ConstraintVec_ constraints_
Definition: linsys.hh:671
Vec rhs_
PETSc vector constructed with vx array.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:147
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
void start_allocation() override
Class for declaration of the integral input data.
Definition: type_base.hh:490
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
Input::Record in_rec_
Definition: linsys.hh:675
void start_add_assembly() override
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:249
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:662
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
double get_solution_precision() override
void start_insert_assembly() override
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:195
#define OLD_ASSERT(...)
Definition: global_defs.h:131
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
MPI_Comm comm_
Definition: linsys.hh:648
unsigned int max_it_
maximum number of iterations of linear solver
Definition: linsys.hh:646
T * makePetscPointer_(std::vector< T > &array)
LinSys_PETSC(const Distribution *rows_ds, const std::string &params="")
Definition: linsys_PETSC.cc:61
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
#define START_TIMER(tag)
Starts a timer with specified tag.
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:654
double a_tol_
absolute tolerance of linear solver
Definition: linsys.hh:645
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:501
bool is_empty() const
Definition: accessors.hh:366
unsigned int np() const
get num of processors
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:669
Vec on_vec_
Vectors for counting non-zero entries in diagonal block.
void preallocate_matrix()
double compute_residual() override
std::string params_
command-line-like options for the PETSc solver
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
KSPConvergedReason reason
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
double r_tol_
relative tolerance of linear solver
Definition: linsys.hh:644
void finish_assembly() override
void view() override
void rhs_set_values(int nrow, int *rows, double *vals) override
bool is_positive_definite()
Definition: linsys.hh:577
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
void set_initial_guess_nonzero(bool set_nonzero=true)
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: linsys_PETSC.cc:90
Mat matrix_
Petsc matrix of the problem.
static const int registrar
Registrar of class to factory.
Record type proxy class.
Definition: type_record.hh:182
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
static Input::Type::Abstract & get_input_type()
Definition: linsys.cc:29
bool init_guess_nonzero
flag for starting from nonzero guess
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:252
Class for declaration of the input data that are in string format.
Definition: type_base.hh:589
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:664
double solution_precision_
unsigned int lsize(int proc) const
get local size