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