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