Flow123d  release_2.2.0-21-g2806cf6
linsys_BDDC.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_BDDC.cc
15  * @brief Solver based on Multilevel BDDC - using corresponding class of OpenFTL package
16  * @author Jakub Sistek
17  */
18 
19 #include <mpi.h>
20 #include "config.h"
21 
22 // need BDDCML wrapper
23 #ifdef FLOW123D_HAVE_BDDCML
24  #include <map>
25  #include "la/bddcml_wrapper.hpp"
26 #endif // FLOW123D_HAVE_BDDCML
27 
28 #include "la/linsys.hh"
29 #include "la/linsys_BDDC.hh"
30 #include "system/sys_profiler.hh"
31 
32 
33 
34 namespace it = Input::Type;
35 
36 
38  return it::Record("Bddc", "Solver setting.")
40  .declare_key("r_tol", it::Double(0.0, 1.0), it::Default::read_time("Defalut value set by nonlinear solver or equation. If not we use value 1.0e-7."),
41  "Relative residual tolerance, (to initial error).")
42  .declare_key("max_it", it::Integer(0), it::Default::read_time("Defalut value set by nonlinear solver or equation. If not we use value 1000."),
43  "Maximum number of outer iterations of the linear solver.")
44 
45  .declare_key("max_nondecr_it", it::Integer(0), it::Default("30"),
46  "Maximum number of iterations of the linear solver with non-decreasing residual.")
47  .declare_key("number_of_levels", it::Integer(0), it::Default("2"),
48  "Number of levels in the multilevel method (=2 for the standard BDDC).")
49  .declare_key("use_adaptive_bddc", it::Bool(), it::Default("false"),
50  "Use adaptive selection of constraints in BDDCML.")
51  .declare_key("bddcml_verbosity_level", it::Integer(0,2), it::Default("0"),
52  "Level of verbosity of the BDDCML library:\n\n - 0 - no output\n - 1 - mild output\n - 2 - detailed output.")
53  .close();
54 }
55 
56 
58 
59 
60 
61 LinSys_BDDC::LinSys_BDDC( const unsigned numDofsSub,
62  const Distribution * rows_ds,
63  const int matrixTypeInt,
64  const int numSubLoc,
65  const bool swap_sign )
66  : LinSys( rows_ds ),
67  swap_sign_(swap_sign)
68 {
69 #ifdef FLOW123D_HAVE_BDDCML
70  // from the point of view of assembly, BDDC linsys is in the ADD state
72 
74  switch ( matrixTypeInt ) {
75  case 0:
76  matrixType = la::BddcmlWrapper::GENERAL;
77  break;
78  case 1:
79  matrixType = la::BddcmlWrapper::SPD;
80  break;
81  case 2:
83  break;
84  case 3:
86  break;
87  default:
88  matrixType = la::BddcmlWrapper::GENERAL;
89  OLD_ASSERT( true, "Unknown matrix type %d", matrixTypeInt );
90  }
91 
92  bddcml_ = new Bddcml_( size_,
93  numDofsSub,
94  matrixType,
95  rows_ds->get_comm(),
96  numSubLoc );
97 
98  // prepare space for local solution
99  locSolution_.resize( numDofsSub );
100 
101  // identify it with PETSc vector
102  PetscErrorCode ierr;
103  //PetscInt numDofsSubInt = static_cast<PetscInt>( numDofsSub );
104  //ierr = VecCreateSeq( PETSC_COMM_SELF, numDofsSubInt, &locSolVec_ );
105  ierr = VecCreateSeqWithArray( PETSC_COMM_SELF, 1, numDofsSub, &(locSolution_[0]), &locSolVec_ );
106  CHKERRV( ierr );
107 #else
108  xprintf(UsrErr, "Compiled without support for BDDCML solver.\n");
109 #endif // FLOW123D_HAVE_BDDCML
110 }
111 
112 
113 void LinSys_BDDC::set_tolerances(double r_tol, double a_tol, unsigned int max_it)
114 {
115  if (! in_rec_.is_empty()) {
116  // input record is set
117  r_tol_ = in_rec_.val<double>("r_tol", r_tol);
118  // BDDC does not use a_tol_
119  a_tol_ = 0.01 * r_tol_;
120  max_it_ = in_rec_.val<unsigned int>("max_it", max_it);\
121  }
122 }
123 
124 
125 void LinSys_BDDC::load_mesh( const int nDim, const int numNodes, const int numDofs,
126  const std::vector<int> & inet,
127  const std::vector<int> & nnet,
128  const std::vector<int> & nndf,
129  const std::vector<int> & isegn,
130  const std::vector<int> & isngn,
131  const std::vector<int> & isvgvn,
132  const std::vector<double> & xyz,
133  const std::vector<double> & element_permeability,
134  const int meshDim )
135 {
136 #ifdef FLOW123D_HAVE_BDDCML
137  // simply pass the data to BDDCML solver
138  isngn_.resize(isngn.size());
139  std::copy( isngn.begin(), isngn.end(), isngn_.begin() );
140  OLD_ASSERT( numDofs == size_, "Global problem size mismatch!" );
141 
142  bddcml_ -> loadRawMesh( nDim, numNodes, numDofs, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
143 
144  // create a map for BDDCML to PETSc vector
145  PetscErrorCode ierr;
146 
147  // local index set
148  PetscInt numDofsSubInt = static_cast<int>( isngn_.size( ) );
150 
151  ISLocalToGlobalMapping subdomainMapping;
152  ierr = ISLocalToGlobalMappingCreate( comm_, 1, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
153 
154  IS subdomainIndexSet;
155  IS from;
156  ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
157  ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
158 
159 
160  // create scatter from parallel PETSc vector to local indices of subdomain
161  ierr = VecScatterCreate( solution_, from, locSolVec_, subdomainIndexSet, &VSpetscToSubScatter_ ); CHKERRV( ierr );
162  ierr = ISDestroy( &subdomainIndexSet ); CHKERRV( ierr );
163  ierr = ISDestroy( &from ); CHKERRV( ierr );
164 
165  //double * locSolVecArray;
166  //ierr = VecGetArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
167  //std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
168  //ierr = VecRestoreArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
169 
170  // scatter local solutions back to global one
171  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
172  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
173 #endif // FLOW123D_HAVE_BDDCML
174 }
175 
176 void LinSys_BDDC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
177 {
178 #ifdef FLOW123D_HAVE_BDDCML
179  namespace ublas = boost::numeric::ublas;
180 
181  std::vector< unsigned > myRows( nrow );
182  std::vector< unsigned > myCols( ncol );
183  ublas::matrix< double > mat( nrow, ncol );
184 
185  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
186  std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
187 
188  for ( int i = 0; i < nrow; i++ ) {
189  for ( int j = 0; j < ncol; j++ ) {
190  mat( i, j ) = vals[i*ncol + j];
191  }
192  }
193  if (swap_sign_) {
194  mat = -mat;
195  }
196 
197  bddcml_ -> insertToMatrix( mat, myRows, myCols );
198 #endif // FLOW123D_HAVE_BDDCML
199 }
200 
201 void LinSys_BDDC::rhs_set_values( int nrow, int *rows, double *vals)
202 {
203 #ifdef FLOW123D_HAVE_BDDCML
204  namespace ublas = boost::numeric::ublas;
205 
206  std::vector< unsigned > myRows( nrow );
207  ublas::vector< double > vec( nrow );
208 
209  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
210 
211  for ( int i = 0; i < nrow; i++ ) {
212  vec( i ) = vals[i];
213  }
214  if (swap_sign_) {
215  vec = -vec;
216  }
217 
218  bddcml_ -> insertToRhs( vec, myRows );
219 #endif // FLOW123D_HAVE_BDDCML
220 }
221 
222 void LinSys_BDDC::diagonal_weights_set_value( int global_index, double value )
223 {
224 #ifdef FLOW123D_HAVE_BDDCML
225  bddcml_ -> insertToDiagonalWeights( global_index, value );
226 #endif // FLOW123D_HAVE_BDDCML
227 }
228 
230 {
231 #ifdef FLOW123D_HAVE_BDDCML
232  bddcml_ -> clearMatrix( );
233 #endif // FLOW123D_HAVE_BDDCML
234  return 0;
235 }
236 
238 {
239 #ifdef FLOW123D_HAVE_BDDCML
240  bddcml_ -> clearRhs( );
241 #endif // FLOW123D_HAVE_BDDCML
242  return 0;
243 }
244 
246 {
247 #ifdef FLOW123D_HAVE_BDDCML
248  bddcml_ -> finishMatAssembly( );
249 #endif // FLOW123D_HAVE_BDDCML
250 }
251 
252 void LinSys_BDDC::apply_constrains( double scalar )
253 {
254 #ifdef FLOW123D_HAVE_BDDCML
255  bddcml_ -> applyConstraints( constraints_, 1., scalar );
256 #endif // FLOW123D_HAVE_BDDCML
257 }
258 
259 int LinSys_BDDC::solve() // ! params are not currently used
260 {
261 #ifdef FLOW123D_HAVE_BDDCML
262  std::vector<int> * numSubAtLevels = NULL; //!< number of subdomains at levels
263  START_TIMER("BDDC linear solver");
264 
265  {
266 
267  START_TIMER("BDDC linear iteration");
269 
270 
271  LogOut().fmt("BDDCML converged reason: {} ( 0 means OK ) \n", bddcml_ -> giveConvergedReason() );
272  LogOut().fmt("BDDCML converged in {} iterations. \n", bddcml_ -> giveNumIterations() );
273  LogOut().fmt("BDDCML estimated condition number is {} \n", bddcml_ -> giveCondNumber() );
274  ADD_CALLS(bddcml_ -> giveNumIterations());
275  }
276 
277  // download local solution
278  bddcml_ -> giveSolution( isngn_, locSolution_ );
279 
280 
281  double * locSolVecArray;
282  VecGetArray( locSolVec_, &locSolVecArray );
283  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
284  VecRestoreArray( locSolVec_, &locSolVecArray );
285 
286  // scatter local solutions back to global one
287  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
288  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
289 
290  // upper bound on the residual error
292 
293  return bddcml_ -> giveConvergedReason();
294 #else
295  return 0;
296 #endif // FLOW123D_HAVE_BDDCML
297 
298 }
299 
301 {
302  LinSys::set_from_input(in_rec);
303 
304  // BDDCML specific parameters
305  max_nondecr_it_ = in_rec.val<int>("max_nondecr_it");
306  number_of_levels_ = in_rec.val<int>("number_of_levels");
307  use_adaptive_bddc_ = in_rec.val<bool>("use_adaptive_bddc");
308  bddcml_verbosity_level_ = in_rec.val<int>("bddcml_verbosity_level");
309 }
310 
312 {
313 #ifdef FLOW123D_HAVE_BDDCML
314  isngn_.clear();
315  locSolution_.clear();
316 
317  PetscErrorCode ierr;
318  ierr = VecDestroy( &locSolVec_ ); CHKERRV( ierr );
319 
320  ierr = VecScatterDestroy( &VSpetscToSubScatter_ ); CHKERRV( ierr );
321 
322  delete bddcml_;
323 #endif // FLOW123D_HAVE_BDDCML
324 };
325 
326 // construct global solution
327 //void LinSys_BDDC::gatherSolution_( )
328 //{
329 //#ifdef FLOW123D_HAVE_BDDCML
330 // int ierr;
331 //
332 // // merge solution on root
333 // int rank;
334 // MPI_Comm_rank( comm_, &rank );
335 // int nProc;
336 // MPI_Comm_size( comm_, &nProc );
337 //
338 // globalSolution_.resize( size_ );
339 // std::vector<double> locSolutionNeib;
340 // if ( rank == 0 ) {
341 // // merge my own data
342 // for ( unsigned int i = 0; i < isngn_.size(); i++ ) {
343 // int ind = isngn_[i];
344 // globalSolution_[ind] = locSolution_[i];
345 // }
346 // for ( int iProc = 1; iProc < nProc; iProc++ ) {
347 // // receive length
348 // int length;
349 // MPI_Status status;
350 // ierr = MPI_Recv( &length, 1, MPI_INT, iProc, iProc, comm_, &status );
351 //
352 // // receive indices
353 // std::vector<int> inds( length );
354 // ierr = MPI_Recv( &(inds[0]), length, MPI_INT, iProc, iProc, comm_, &status );
355 //
356 // // receive local solution
357 // locSolutionNeib.resize( length );
358 // ierr = MPI_Recv( &(locSolutionNeib[0]), length, MPI_DOUBLE, iProc, iProc, comm_, &status );
359 //
360 // // merge data
361 // for ( int i = 0; i < length; i++ ) {
362 // int ind = inds[i];
363 // globalSolution_[ind] = locSolutionNeib[i];
364 // }
365 // }
366 // }
367 // else {
368 // // send my solution to root
369 // int length = isngn_.size();
370 // ierr = MPI_Send( &length, 1, MPI_INT, 0, rank, comm_ );
371 // ierr = MPI_Send( &(isngn_[0]), length, MPI_INT, 0, rank, comm_ );
372 // ierr = MPI_Send( &(locSolution_[0]), length, MPI_DOUBLE, 0, rank, comm_ );
373 // }
374 // // broadcast global solution from root
375 // ierr = MPI_Bcast( &(globalSolution_[0]), globalSolution_.size(), MPI_DOUBLE, 0, comm_ );
376 //#endif // FLOW123D_HAVE_BDDCML
377 //}
378 
380 {
381  double bnorm=0.0;
382  VecNorm(locSolVec_, NORM_2, &bnorm);
383 
384  return max(a_tol_, r_tol_*bnorm);
385 }
386 
387 
int bddcml_verbosity_level_
Definition: linsys_BDDC.hh:102
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:583
double normRhs()
Get norm of right-hand side.
int solve() override
Definition: linsys_BDDC.cc:259
double get_solution_precision() override
Definition: linsys_BDDC.cc:379
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
VecScatter VSpetscToSubScatter_
scatter from solution_ to locSolVec_
Definition: linsys_BDDC.hh:112
std::vector< int > isngn_
indices of subdomain nodes in global numbering
Definition: linsys_BDDC.hh:109
PetscErrorCode rhs_zero_entries() override
Definition: linsys_BDDC.cc:237
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:554
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Class for declaration of the input of type Bool.
Definition: type_base.hh:458
Wrappers for linear systems based on MPIAIJ and MATIS format.
ConstraintVec_ constraints_
Definition: linsys.hh:604
enum matrixTypeEnum MatrixType
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
unsigned size_
global number of matrix rows, i.e. problem size
Definition: linsys.hh:586
Class for declaration of the integral input data.
Definition: type_base.hh:489
void load_mesh(const int nDim, const int numNodes, const int numDofs, const std::vector< int > &inet, const std::vector< int > &nnet, const std::vector< int > &nndf, const std::vector< int > &isegn, const std::vector< int > &isngn, const std::vector< int > &isvgvn, const std::vector< double > &xyz, const std::vector< double > &element_permeability, const int meshDim)
Definition: linsys_BDDC.cc:125
Input::Record in_rec_
Definition: linsys.hh:608
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:239
static const Input::Type::Record & get_input_type()
Definition: linsys_BDDC.cc:37
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
bool use_adaptive_bddc_
should adaptive BDDC be used?
Definition: linsys_BDDC.hh:101
static constexpr bool value
Definition: json.hpp:87
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
void set_from_input(const Input::Record in_rec) override
Definition: linsys_BDDC.cc:300
void apply_constrains(double scalar=1.) override
Definition: linsys_BDDC.cc:252
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:540
MPI_Comm comm_
Definition: linsys.hh:582
unsigned int max_it_
maximum number of iterations of linear solver
Definition: linsys.hh:580
Vec locSolVec_
local solution PETSc vector - sequential
Definition: linsys_BDDC.hh:111
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:92
#define START_TIMER(tag)
Starts a timer with specified tag.
int max_nondecr_it_
parameters expected from input file:
Definition: linsys_BDDC.hh:99
double a_tol_
absolute tolerance of linear solver
Definition: linsys.hh:579
static const int registrar
Registrar of class to factory.
Definition: linsys_BDDC.hh:96
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:490
void diagonal_weights_set_value(int global_index, double value)
Definition: linsys_BDDC.cc:222
bool is_empty() const
Definition: accessors.hh:366
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:602
const bool swap_sign_
swap sign of matrix and rhs entries, e.g. to make the matrix SPD
Definition: linsys_BDDC.hh:107
int number_of_levels_
number of levels in the multilevel method
Definition: linsys_BDDC.hh:100
la::BddcmlWrapper Bddcml_
Definition: linsys_BDDC.hh:114
Definition: system.hh:64
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
symmetric positive definite,
general symmetric ( e.g. saddle point ),
double r_tol_
relative tolerance of linear solver
Definition: linsys.hh:578
Bddcml_ * bddcml_
BDDCML wrapper.
Definition: linsys_BDDC.hh:115
PetscErrorCode mat_zero_entries() override
Definition: linsys_BDDC.cc:229
MPI_Comm get_comm() const
Returns communicator.
Abstract linear system class.
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
Definition: linsys_BDDC.cc:176
void finish_assembly() override
Definition: linsys_BDDC.cc:245
Record type proxy class.
Definition: type_record.hh:182
static Input::Type::Abstract & get_input_type()
Definition: linsys.cc:29
std::vector< double > locSolution_
subdomain solution
Definition: linsys_BDDC.hh:110
LinSys_BDDC(const unsigned numDofsSub, const Distribution *rows_ds, const int matrixTypeInt=0, const int numSubLoc=1, const bool swap_sign=false)
Definition: linsys_BDDC.cc:61
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:598
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: linsys_BDDC.cc:113
void rhs_set_values(int nrow, int *rows, double *vals) override
Definition: linsys_BDDC.cc:201
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.