Flow123d  jenkins-Flow123d-windows-release-multijob-285
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 #include "system/sys_profiler.hh"
43 
44 
45 
46 namespace it = Input::Type;
47 
48 
49 it::Record LinSys_BDDC::input_type = it::Record("Bddc", "Solver setting.")
51  .declare_key("max_nondecr_it", it::Integer(0), it::Default("30"),
52  "Maximum number of iterations of the linear solver with non-decreasing residual.")
53  .declare_key("number_of_levels", it::Integer(0), it::Default("2"),
54  "Number of levels in the multilevel method (=2 for the standard BDDC).")
55  .declare_key("use_adaptive_bddc", it::Bool(), it::Default("false"),
56  "Use adaptive selection of constraints in BDDCML.")
57  .declare_key("bddcml_verbosity_level", it::Integer(0,2), it::Default("0"),
58  "Level of verbosity of the BDDCML library: 0 - no output, 1 - mild output, 2 - detailed output.");
59 
60 
61 
62 LinSys_BDDC::LinSys_BDDC( const unsigned numDofsSub,
63  const Distribution * rows_ds,
64  const int matrixTypeInt,
65  const int numSubLoc,
66  const bool swap_sign )
67  : LinSys( rows_ds ),
68  swap_sign_(swap_sign)
69 {
70 #ifdef HAVE_BDDCML
71  // from the point of view of assembly, BDDC linsys is in the ADD state
73 
75  switch ( matrixTypeInt ) {
76  case 0:
77  matrixType = la::BddcmlWrapper::GENERAL;
78  break;
79  case 1:
80  matrixType = la::BddcmlWrapper::SPD;
81  break;
82  case 2:
84  break;
85  case 3:
87  break;
88  default:
89  matrixType = la::BddcmlWrapper::GENERAL;
90  ASSERT( true, "Unknown matrix type %d", matrixTypeInt );
91  }
92 
93  bddcml_ = new Bddcml_( size_,
94  numDofsSub,
95  matrixType,
96  rows_ds->get_comm(),
97  numSubLoc );
98 
99  // prepare space for local solution
100  locSolution_.resize( numDofsSub );
101 
102  // identify it with PETSc vector
103  PetscErrorCode ierr;
104  PetscInt numDofsSubInt = static_cast<PetscInt>( numDofsSub );
105  //ierr = VecCreateSeq( PETSC_COMM_SELF, numDofsSubInt, &locSolVec_ );
106  ierr = VecCreateSeqWithArray( PETSC_COMM_SELF, 1, numDofsSub, &(locSolution_[0]), &locSolVec_ );
107  CHKERRV( ierr );
108 #else
109  xprintf(UsrErr, "Compiled without support for BDDCML solver.\n");
110 #endif // HAVE_BDDCML
111 }
112 
113 void LinSys_BDDC::load_mesh( const int nDim, const int numNodes, const int numDofs,
114  const std::vector<int> & inet,
115  const std::vector<int> & nnet,
116  const std::vector<int> & nndf,
117  const std::vector<int> & isegn,
118  const std::vector<int> & isngn,
119  const std::vector<int> & isvgvn,
120  const std::vector<double> & xyz,
121  const std::vector<double> & element_permeability,
122  const int meshDim )
123 {
124 #ifdef HAVE_BDDCML
125  // simply pass the data to BDDCML solver
126  isngn_.resize(isngn.size());
127  std::copy( isngn.begin(), isngn.end(), isngn_.begin() );
128  ASSERT( numDofs == size_, "Global problem size mismatch!" );
129 
130  bddcml_ -> loadRawMesh( nDim, numNodes, numDofs, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
131 
132  // create a map for BDDCML to PETSc vector
133  PetscErrorCode ierr;
134 
135  // local index set
136  PetscInt numDofsSubInt = static_cast<int>( isngn_.size( ) );
138 
139  ISLocalToGlobalMapping subdomainMapping;
140  ierr = ISLocalToGlobalMappingCreate( comm_, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
141 
142  IS subdomainIndexSet;
143  IS from;
144  ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
145  ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
146 
147 
148  // create scatter from parallel PETSc vector to local indices of subdomain
149  ierr = VecScatterCreate( solution_, from, locSolVec_, subdomainIndexSet, &VSpetscToSubScatter_ ); CHKERRV( ierr );
150  ierr = ISDestroy( &subdomainIndexSet ); CHKERRV( ierr );
151  ierr = ISDestroy( &from ); CHKERRV( ierr );
152 
153  //double * locSolVecArray;
154  //ierr = VecGetArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
155  //std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
156  //ierr = VecRestoreArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
157 
158  // scatter local solutions back to global one
159  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
160  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
161 #endif // HAVE_BDDCML
162 }
163 
164 void LinSys_BDDC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
165 {
166 #ifdef HAVE_BDDCML
167  namespace ublas = boost::numeric::ublas;
168 
169  std::vector< unsigned > myRows( nrow );
170  std::vector< unsigned > myCols( ncol );
171  ublas::matrix< double > mat( nrow, ncol );
172 
173  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
174  std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
175 
176  for ( int i = 0; i < nrow; i++ ) {
177  for ( int j = 0; j < ncol; j++ ) {
178  mat( i, j ) = vals[i*ncol + j];
179  }
180  }
181  if (swap_sign_) {
182  mat = -mat;
183  }
184 
185  bddcml_ -> insertToMatrix( mat, myRows, myCols );
186 #endif // HAVE_BDDCML
187 }
188 
189 void LinSys_BDDC::rhs_set_values( int nrow, int *rows, double *vals)
190 {
191 #ifdef HAVE_BDDCML
192  namespace ublas = boost::numeric::ublas;
193 
194  std::vector< unsigned > myRows( nrow );
195  ublas::vector< double > vec( nrow );
196 
197  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
198 
199  for ( int i = 0; i < nrow; i++ ) {
200  vec( i ) = vals[i];
201  }
202  if (swap_sign_) {
203  vec = -vec;
204  }
205 
206  bddcml_ -> insertToRhs( vec, myRows );
207 #endif // HAVE_BDDCML
208 }
209 
210 void LinSys_BDDC::diagonal_weights_set_value( int global_index, double value )
211 {
212 #ifdef HAVE_BDDCML
213  bddcml_ -> insertToDiagonalWeights( global_index, value );
214 #endif // HAVE_BDDCML
215 }
216 
218 {
219 #ifdef HAVE_BDDCML
220  bddcml_ -> clearMatrix( );
221 #endif // HAVE_BDDCML
222  return 0;
223 }
224 
226 {
227 #ifdef HAVE_BDDCML
228  bddcml_ -> clearRhs( );
229 #endif // HAVE_BDDCML
230  return 0;
231 }
232 
234 {
235 #ifdef HAVE_BDDCML
236  bddcml_ -> finishMatAssembly( );
237 #endif // HAVE_BDDCML
238 }
239 
240 void LinSys_BDDC::apply_constrains( double scalar )
241 {
242 #ifdef HAVE_BDDCML
243  bddcml_ -> applyConstraints( constraints_, 1., scalar );
244 #endif // HAVE_BDDCML
245 }
246 
247 int LinSys_BDDC::solve() // ! params are not currently used
248 {
249 #ifdef HAVE_BDDCML
250  std::vector<int> * numSubAtLevels = NULL; //!< number of subdomains at levels
251 
252  {
253  START_TIMER("BDDC linear solver");
254  START_TIMER("BDDC linear iteration");
256 
257 
258  xprintf(MsgLog,"BDDCML converged reason: %d ( 0 means OK ) \n", bddcml_ -> giveConvergedReason() );
259  xprintf(MsgLog,"BDDCML converged in %d iterations. \n", bddcml_ -> giveNumIterations() );
260  xprintf(MsgLog,"BDDCML estimated condition number is %f \n", bddcml_ -> giveCondNumber() );
261  ADD_CALLS(bddcml_ -> giveNumIterations());
262  }
263 
264  // download local solution
265  bddcml_ -> giveSolution( isngn_, locSolution_ );
266 
267 
268  double * locSolVecArray;
269  VecGetArray( locSolVec_, &locSolVecArray );
270  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
271  VecRestoreArray( locSolVec_, &locSolVecArray );
272 
273  // scatter local solutions back to global one
274  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
275  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
276 
277  // upper bound on the residual error
279 
280  return bddcml_ -> giveConvergedReason();
281 #else
282  return 0;
283 #endif // HAVE_BDDCML
284 
285 }
286 
288 {
289  // common values
290  r_tol_ = in_rec.val<double>("r_tol");
291  max_it_ = in_rec.val<int>("max_it");
292 
293  // BDDCML specific parameters
294  max_nondecr_it_ = in_rec.val<int>("max_nondecr_it");
295  number_of_levels_ = in_rec.val<int>("number_of_levels");
296  use_adaptive_bddc_ = in_rec.val<bool>("use_adaptive_bddc");
297  bddcml_verbosity_level_ = in_rec.val<int>("bddcml_verbosity_level");
298 }
299 
301 {
302 #ifdef HAVE_BDDCML
303  isngn_.clear();
304  locSolution_.clear();
305 
306  PetscErrorCode ierr;
307  ierr = VecDestroy( &locSolVec_ ); CHKERRV( ierr );
308 
309  ierr = VecScatterDestroy( &VSpetscToSubScatter_ ); CHKERRV( ierr );
310 
311  delete bddcml_;
312 #endif // HAVE_BDDCML
313 };
314 
315 // construct global solution
316 //void LinSys_BDDC::gatherSolution_( )
317 //{
318 //#ifdef HAVE_BDDCML
319 // int ierr;
320 //
321 // // merge solution on root
322 // int rank;
323 // MPI_Comm_rank( comm_, &rank );
324 // int nProc;
325 // MPI_Comm_size( comm_, &nProc );
326 //
327 // globalSolution_.resize( size_ );
328 // std::vector<double> locSolutionNeib;
329 // if ( rank == 0 ) {
330 // // merge my own data
331 // for ( unsigned int i = 0; i < isngn_.size(); i++ ) {
332 // int ind = isngn_[i];
333 // globalSolution_[ind] = locSolution_[i];
334 // }
335 // for ( int iProc = 1; iProc < nProc; iProc++ ) {
336 // // receive length
337 // int length;
338 // MPI_Status status;
339 // ierr = MPI_Recv( &length, 1, MPI_INT, iProc, iProc, comm_, &status );
340 //
341 // // receive indices
342 // std::vector<int> inds( length );
343 // ierr = MPI_Recv( &(inds[0]), length, MPI_INT, iProc, iProc, comm_, &status );
344 //
345 // // receive local solution
346 // locSolutionNeib.resize( length );
347 // ierr = MPI_Recv( &(locSolutionNeib[0]), length, MPI_DOUBLE, iProc, iProc, comm_, &status );
348 //
349 // // merge data
350 // for ( int i = 0; i < length; i++ ) {
351 // int ind = inds[i];
352 // globalSolution_[ind] = locSolutionNeib[i];
353 // }
354 // }
355 // }
356 // else {
357 // // send my solution to root
358 // int length = isngn_.size();
359 // ierr = MPI_Send( &length, 1, MPI_INT, 0, rank, comm_ );
360 // ierr = MPI_Send( &(isngn_[0]), length, MPI_INT, 0, rank, comm_ );
361 // ierr = MPI_Send( &(locSolution_[0]), length, MPI_DOUBLE, 0, rank, comm_ );
362 // }
363 // // broadcast global solution from root
364 // ierr = MPI_Bcast( &(globalSolution_[0]), globalSolution_.size(), MPI_DOUBLE, 0, comm_ );
365 //#endif // HAVE_BDDCML
366 //}
367 
369 {
370  double bnorm=0.0;
371  VecNorm(locSolVec_, NORM_2, &bnorm);
372 
373  return max(a_tol_, r_tol_*bnorm);
374 }
375 
376 
int bddcml_verbosity_level_
Definition: linsys_BDDC.hh:102
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:576
double normRhs()
Get norm of right-hand side.
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)
Definition: linsys_BDDC.cc:164
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:225
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:332
Wrappers for linear systems based on MPIAIJ and MATIS format.
ConstraintVec_ constraints_
Definition: linsys.hh:597
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:197
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:579
Class for declaration of the integral input data.
Definition: type_base.hh:355
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:113
Definition: system.hh:72
void rhs_set_values(int nrow, int *rows, double *vals)
Definition: linsys_BDDC.cc:189
bool use_adaptive_bddc_
should adaptive BDDC be used?
Definition: linsys_BDDC.hh:101
#define ASSERT(...)
Definition: global_defs.h:121
MPI_Comm comm_
Definition: linsys.hh:575
static Input::Type::AbstractRecord input_type
Definition: linsys.hh:90
int max_it_
Definition: linsys.hh:573
Vec locSolVec_
local solution PETSc vector - sequential
Definition: linsys_BDDC.hh:111
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:100
#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_
Definition: linsys.hh:572
void diagonal_weights_set_value(int global_index, double value)
Definition: linsys_BDDC.cc:210
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:595
const bool swap_sign_
swap sign of matrix and rhs entries, e.g. to make the matrix SPD
Definition: linsys_BDDC.hh:107
double get_solution_precision()
Definition: linsys_BDDC.cc:368
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:72
symmetric positive definite,
general symmetric ( e.g. saddle point ),
static Input::Type::Record input_type
Definition: linsys_BDDC.hh:51
double r_tol_
Definition: linsys.hh:571
Bddcml_ * bddcml_
BDDCML wrapper.
Definition: linsys_BDDC.hh:115
PetscErrorCode mat_zero_entries() override
Definition: linsys_BDDC.cc:217
MPI_Comm get_comm() const
Returns communicator.
Abstract linear system class.
void apply_constrains(double scalar=1.)
Definition: linsys_BDDC.cc:240
Record type proxy class.
Definition: type_record.hh:169
void finish_assembly()
Definition: linsys_BDDC.cc:233
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:62
void set_from_input(const Input::Record in_rec)
Definition: linsys_BDDC.cc:287
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:591
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:430