Flow123d  jenkins-Flow123d-windows32-release-multijob-28
linsys_BDDC.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id: la_linsys.hh 1299 2011-08-23 21:42:50Z jakub.sistek $
21  * $Revision: 1299 $
22  * $LastChangedBy: jakub.sistek $
23  * $LastChangedDate: 2011-08-23 23:42:50 +0200 (Tue, 23 Aug 2011) $
24  *
25  * @file
26  * @brief Solver based on Multilevel BDDC - using corresponding class of OpenFTL package
27  * @author Jakub Sistek
28  *
29  *
30  */
31 
32 #include <mpi.h>
33 
34 // need BDDCML wrapper
35 #ifdef HAVE_BDDCML
36  #include <map>
37  #include "la/bddcml_wrapper.hpp"
38 #endif // HAVE_BDDCML
39 
40 #include "la/linsys.hh"
41 #include "la/linsys_BDDC.hh"
42 
43 
44 
45 namespace it = Input::Type;
46 
47 
48 it::Record LinSys_BDDC::input_type = it::Record("Bddc", "Solver setting.")
50  .declare_key("max_nondecr_it", it::Integer(0), it::Default("30"),
51  "Maximum number of iterations of the linear solver with non-decreasing residual.")
52  .declare_key("number_of_levels", it::Integer(0), it::Default("2"),
53  "Number of levels in the multilevel method (=2 for the standard BDDC).")
54  .declare_key("use_adaptive_bddc", it::Bool(), it::Default("false"),
55  "Use adaptive selection of constraints in BDDCML.")
56  .declare_key("bddcml_verbosity_level", it::Integer(0,2), it::Default("0"),
57  "Level of verbosity of the BDDCML library: 0 - no output, 1 - mild output, 2 - detailed output.");
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 HAVE_BDDCML
70  // set type
71  //type = LinSys::BDDC;
72 
73  // from the point of view of assembly, BDDC linsys is in the ADD state
75 
77  switch ( matrixTypeInt ) {
78  case 0:
79  matrixType = la::BddcmlWrapper::GENERAL;
80  break;
81  case 1:
82  matrixType = la::BddcmlWrapper::SPD;
83  break;
84  case 2:
86  break;
87  case 3:
89  break;
90  default:
91  matrixType = la::BddcmlWrapper::GENERAL;
92  ASSERT( true, "Unknown matrix type %d", matrixTypeInt );
93  }
94 
95  bddcml_ = new Bddcml_( size_,
96  numDofsSub,
97  matrixType,
98  rows_ds->get_comm(),
99  numSubLoc );
100 
101  // prepare space for local solution
102  locSolution_.resize( numDofsSub );
103 
104  // identify it with PETSc vector
105  PetscErrorCode ierr;
106  PetscInt numDofsSubInt = static_cast<PetscInt>( numDofsSub );
107  ierr = VecCreateSeq( PETSC_COMM_SELF, numDofsSubInt, &locSolVec_ );
108  CHKERRV( ierr );
109 #else
110  xprintf(UsrErr, "Compiled without support for BDDCML solver.\n");
111 #endif // HAVE_BDDCML
112 }
113 
114 void LinSys_BDDC::load_mesh( const int nDim, const int numNodes, const int numDofs,
115  const std::vector<int> & inet,
116  const std::vector<int> & nnet,
117  const std::vector<int> & nndf,
118  const std::vector<int> & isegn,
119  const std::vector<int> & isngn,
120  const std::vector<int> & isvgvn,
121  const std::vector<double> & xyz,
122  const std::vector<double> & element_permeability,
123  const int meshDim )
124 {
125 #ifdef HAVE_BDDCML
126  // simply pass the data to BDDCML solver
127  isngn_.resize(isngn.size());
128  std::copy( isngn.begin(), isngn.end(), isngn_.begin() );
129  ASSERT( numDofs == size_, "Global problem size mismatch!" );
130 
131  bddcml_ -> loadRawMesh( nDim, numNodes, numDofs, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
132 
133  // create a map for BDDCML to PETSc vector
134  PetscErrorCode ierr;
135 
136  // local index set
137  PetscInt numDofsSubInt = static_cast<int>( isngn_.size( ) );
139 
140  //std::cout << "IDX: \n";
141  //std::copy( idx.begin(), idx.end(), std::ostream_iterator<PetscInt>( std::cout, " " ) );
142 
143  ISLocalToGlobalMapping subdomainMapping;
144  ierr = ISLocalToGlobalMappingCreate( comm_, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
145 
146  IS subdomainIndexSet;
147  IS from;
148  ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
149  ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
150 
151  //ierr = ISCreateGeneral( comm_, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
152  //ISView( subdomainIndexSet, PETSC_VIEWER_STDOUT_WORLD);
153 
154 
155  // create scatter from parallel PETSc vector to local indices of subdomain
156  ierr = VecScatterCreate( solution_, from, locSolVec_, subdomainIndexSet, &VSpetscToSubScatter_ ); CHKERRV( ierr );
157  ierr = ISDestroy( &subdomainIndexSet ); CHKERRV( ierr );
158  ierr = ISDestroy( &from ); CHKERRV( ierr );
159 
160  //VecScatterView(VSpetscToSubScatter_,PETSC_VIEWER_STDOUT_SELF);
161 
162  double * locSolVecArray;
163  ierr = VecGetArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
164  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
165  ierr = VecRestoreArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
166 
167  // scatter local solutions back to global one
168  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
169  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
170 #endif // HAVE_BDDCML
171 }
172 
173 void LinSys_BDDC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
174 {
175 #ifdef HAVE_BDDCML
176  namespace ublas = boost::numeric::ublas;
177 
178  std::vector< unsigned > myRows( nrow );
179  std::vector< unsigned > myCols( ncol );
180  ublas::matrix< double > mat( nrow, ncol );
181 
182  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
183  std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
184 
185  for ( int i = 0; i < nrow; i++ ) {
186  for ( int j = 0; j < ncol; j++ ) {
187  mat( i, j ) = vals[i*ncol + j];
188  }
189  }
190  if (swap_sign_) {
191  mat = -mat;
192  }
193 
194  bddcml_ -> insertToMatrix( mat, myRows, myCols );
195 #endif // HAVE_BDDCML
196 }
197 
198 void LinSys_BDDC::rhs_set_values( int nrow, int *rows, double *vals)
199 {
200 #ifdef HAVE_BDDCML
201  namespace ublas = boost::numeric::ublas;
202 
203  std::vector< unsigned > myRows( nrow );
204  ublas::vector< double > vec( nrow );
205 
206  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
207 
208  for ( int i = 0; i < nrow; i++ ) {
209  vec( i ) = vals[i];
210  }
211  if (swap_sign_) {
212  vec = -vec;
213  }
214 
215  bddcml_ -> insertToRhs( vec, myRows );
216 #endif // HAVE_BDDCML
217 }
218 
219 void LinSys_BDDC::diagonal_weights_set_value( int global_index, double value )
220 {
221 #ifdef HAVE_BDDCML
222  bddcml_ -> insertToDiagonalWeights( global_index, value );
223 #endif // HAVE_BDDCML
224 }
225 
227 {
228 #ifdef HAVE_BDDCML
229  bddcml_ -> clearMatrix( );
230 #endif // HAVE_BDDCML
231  return 0;
232 }
233 
235 {
236 #ifdef HAVE_BDDCML
237  bddcml_ -> clearRhs( );
238 #endif // HAVE_BDDCML
239  return 0;
240 }
241 
243 {
244 #ifdef HAVE_BDDCML
245  bddcml_ -> finishMatAssembly( );
246 #endif // HAVE_BDDCML
247 }
248 
249 void LinSys_BDDC::apply_constrains( double scalar )
250 {
251 #ifdef HAVE_BDDCML
252  bddcml_ -> applyConstraints( constraints_, 1., scalar );
253 #endif // HAVE_BDDCML
254 }
255 
256 int LinSys_BDDC::solve() // ! params are not currently used
257 {
258 #ifdef HAVE_BDDCML
259  std::vector<int> * numSubAtLevels = NULL; //!< number of subdomains at levels
260 
262 
263  DBGMSG("BDDCML converged reason: %d ( 0 means OK ) \n", bddcml_ -> giveConvergedReason() );
264  DBGMSG("BDDCML converged in %d iterations. \n", bddcml_ -> giveNumIterations() );
265  DBGMSG("BDDCML estimated condition number is %f \n", bddcml_ -> giveCondNumber() );
266 
267  // download local solution
268  bddcml_ -> giveSolution( isngn_, locSolution_ );
269 
270 
271  double * locSolVecArray;
272  PetscErrorCode ierr;
273  ierr = VecGetArray( locSolVec_, &locSolVecArray );
274  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
275  ierr = VecRestoreArray( locSolVec_, &locSolVecArray );
276 
277  // scatter local solutions back to global one
278  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
279  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
280 
281  // upper bound on the residual error
283 
284  return bddcml_ -> giveConvergedReason();
285 #else
286  return 0;
287 #endif // HAVE_BDDCML
288 
289 }
290 
292 {
293 #ifdef HAVE_BDDCML
294  this -> gatherSolution_( );
295  globalSolution.resize( globalSolution_.size( ) );
296  std::copy( globalSolution_.begin(), globalSolution_.end(), globalSolution.begin() );
297 #endif // HAVE_BDDCML
298 }
299 
301 {
302 #ifdef HAVE_BDDCML
303  globalSolution_.resize( globalSolution.size( ) );
304  std::copy( globalSolution.begin(), globalSolution.end(), globalSolution_.begin() );
305 #endif // HAVE_BDDCML
306 }
307 
309 {
310  // common values
311  r_tol_ = in_rec.val<double>("r_tol");
312  max_it_ = in_rec.val<int>("max_it");
313 
314  // BDDCML specific parameters
315  max_nondecr_it_ = in_rec.val<int>("max_nondecr_it");
316  number_of_levels_ = in_rec.val<int>("number_of_levels");
317  use_adaptive_bddc_ = in_rec.val<bool>("use_adaptive_bddc");
318  bddcml_verbosity_level_ = in_rec.val<int>("bddcml_verbosity_level");
319 }
320 
322 {
323 #ifdef HAVE_BDDCML
324  isngn_.clear();
325  locSolution_.clear();
326 
327  PetscErrorCode ierr;
328  ierr = VecDestroy( &locSolVec_ ); CHKERRV( ierr );
329 
330  ierr = VecScatterDestroy( &VSpetscToSubScatter_ ); CHKERRV( ierr );
331 
332  delete bddcml_;
333 #endif // HAVE_BDDCML
334 };
335 
336 // construct global solution
338 {
339 #ifdef HAVE_BDDCML
340  int ierr;
341 
342  // merge solution on root
343  int rank;
344  MPI_Comm_rank( comm_, &rank );
345  int nProc;
346  MPI_Comm_size( comm_, &nProc );
347 
348  globalSolution_.resize( size_ );
349  std::vector<double> locSolutionNeib;
350  if ( rank == 0 ) {
351  // merge my own data
352  for ( unsigned int i = 0; i < isngn_.size(); i++ ) {
353  int ind = isngn_[i];
354  globalSolution_[ind] = locSolution_[i];
355  }
356  for ( int iProc = 1; iProc < nProc; iProc++ ) {
357  // receive length
358  int length;
360  ierr = MPI_Recv( &length, 1, MPI_INT, iProc, iProc, comm_, &status );
361 
362  // receive indices
363  std::vector<int> inds( length );
364  ierr = MPI_Recv( &(inds[0]), length, MPI_INT, iProc, iProc, comm_, &status );
365 
366  // receive local solution
367  locSolutionNeib.resize( length );
368  ierr = MPI_Recv( &(locSolutionNeib[0]), length, MPI_DOUBLE, iProc, iProc, comm_, &status );
369 
370  // merge data
371  for ( int i = 0; i < length; i++ ) {
372  int ind = inds[i];
373  globalSolution_[ind] = locSolutionNeib[i];
374  }
375  }
376  }
377  else {
378  // send my solution to root
379  int length = isngn_.size();
380  ierr = MPI_Send( &length, 1, MPI_INT, 0, rank, comm_ );
381  ierr = MPI_Send( &(isngn_[0]), length, MPI_INT, 0, rank, comm_ );
382  ierr = MPI_Send( &(locSolution_[0]), length, MPI_DOUBLE, 0, rank, comm_ );
383  }
384  // broadcast global solution from root
385  ierr = MPI_Bcast( &(globalSolution_[0]), globalSolution_.size(), MPI_DOUBLE, 0, comm_ );
386 #endif // HAVE_BDDCML
387 }
388 
390 {
391  double bnorm=0.0;
392  VecNorm(locSolVec_, NORM_2, &bnorm);
393 
394  return max(a_tol_, r_tol_*bnorm);
395 }
396 
397 
#define MPI_Recv(buf, count, datatype, source, tag, comm, status)
Definition: mpi.h:271
int bddcml_verbosity_level_
Definition: linsys_BDDC.hh:109
std::vector< double > globalSolution_
global solution in numbering for linear system
Definition: linsys.hh:621
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:598
double normRhs()
Get norm of right-hand side.
#define DBGMSG(...)
Definition: global_defs.h:195
void gatherSolution_()
Definition: linsys_BDDC.cc:337
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)
Definition: linsys_BDDC.cc:173
VecScatter VSpetscToSubScatter_
scatter from solution_ to locSolVec_
Definition: linsys_BDDC.hh:119
std::vector< int > isngn_
indices of subdomain nodes in global numbering
Definition: linsys_BDDC.hh:116
PetscErrorCode rhs_zero_entries() override
Definition: linsys_BDDC.cc:234
void get_whole_solution(std::vector< double > &globalSolution)
Definition: linsys_BDDC.cc:291
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:39
Class for declaration of the input of type Bool.
Definition: type_base.hh:321
SetValuesMode status
Set value status of the linear system.
Wrappers for linear systems based on MPIAIJ and MATIS format.
ConstraintVec_ constraints_
Definition: linsys.hh:619
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:172
enum matrixTypeEnum MatrixType
#define MPI_Send(buf, count, datatype, dest, tag, comm)
Definition: mpi.h:263
unsigned size_
global number of matrix rows, i.e. problem size
Definition: linsys.hh:601
Class for declaration of the integral input data.
Definition: type_base.hh:341
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:114
void rhs_set_values(int nrow, int *rows, double *vals)
Definition: linsys_BDDC.cc:198
bool use_adaptive_bddc_
should adaptive BDDC be used?
Definition: linsys_BDDC.hh:108
#define ASSERT(...)
Definition: global_defs.h:120
MPI_Comm comm_
Definition: linsys.hh:597
static Input::Type::AbstractRecord input_type
Definition: linsys.hh:92
int max_it_
Definition: linsys.hh:595
Vec locSolVec_
local solution PETSc vector - sequential
Definition: linsys_BDDC.hh:118
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:104
int max_nondecr_it_
parameters expected from input file:
Definition: linsys_BDDC.hh:106
double a_tol_
Definition: linsys.hh:594
void diagonal_weights_set_value(int global_index, double value)
Definition: linsys_BDDC.cc:219
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:617
const bool swap_sign_
swap sign of matrix and rhs entries, e.g. to make the matrix SPD
Definition: linsys_BDDC.hh:114
#define MPI_Comm_size
Definition: mpi.h:235
#define MPI_DOUBLE
Definition: mpi.h:156
double get_solution_precision()
Definition: linsys_BDDC.cc:389
int number_of_levels_
number of levels in the multilevel method
Definition: linsys_BDDC.hh:107
#define MPI_Comm_rank
Definition: mpi.h:236
la::BddcmlWrapper Bddcml_
Definition: linsys_BDDC.hh:121
Definition: system.hh:75
symmetric positive definite,
general symmetric ( e.g. saddle point ),
static Input::Type::Record input_type
Definition: linsys_BDDC.hh:54
double r_tol_
Definition: linsys.hh:593
Bddcml_ * bddcml_
BDDCML wrapper.
Definition: linsys_BDDC.hh:122
PetscErrorCode mat_zero_entries() override
Definition: linsys_BDDC.cc:226
#define MPI_INT
Definition: mpi.h:160
MPI_Comm get_comm() const
Returns communicator.
Abstract linear system class.
#define MPI_Bcast(buffer, count, datatype, root, comm)
Definition: mpi.h:534
void apply_constrains(double scalar=1.)
Definition: linsys_BDDC.cc:249
Record type proxy class.
Definition: type_record.hh:161
void finish_assembly()
Definition: linsys_BDDC.cc:242
std::vector< double > locSolution_
subdomain solution
Definition: linsys_BDDC.hh:117
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
void set_whole_solution(std::vector< double > &globalSolution)
Definition: linsys_BDDC.cc:300
void set_from_input(const Input::Record in_rec)
Definition: linsys_BDDC.cc:308
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:613
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:390