Flow123d  jenkins-Flow123d-windows32-release-multijob-51
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  // 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  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 // HAVE_BDDCML
110 }
111 
112 void LinSys_BDDC::load_mesh( const int nDim, const int numNodes, const int numDofs,
113  const std::vector<int> & inet,
114  const std::vector<int> & nnet,
115  const std::vector<int> & nndf,
116  const std::vector<int> & isegn,
117  const std::vector<int> & isngn,
118  const std::vector<int> & isvgvn,
119  const std::vector<double> & xyz,
120  const std::vector<double> & element_permeability,
121  const int meshDim )
122 {
123 #ifdef HAVE_BDDCML
124  // simply pass the data to BDDCML solver
125  isngn_.resize(isngn.size());
126  std::copy( isngn.begin(), isngn.end(), isngn_.begin() );
127  ASSERT( numDofs == size_, "Global problem size mismatch!" );
128 
129  bddcml_ -> loadRawMesh( nDim, numNodes, numDofs, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
130 
131  // create a map for BDDCML to PETSc vector
132  PetscErrorCode ierr;
133 
134  // local index set
135  PetscInt numDofsSubInt = static_cast<int>( isngn_.size( ) );
137 
138  ISLocalToGlobalMapping subdomainMapping;
139  ierr = ISLocalToGlobalMappingCreate( comm_, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
140 
141  IS subdomainIndexSet;
142  IS from;
143  ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
144  ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
145 
146 
147  // create scatter from parallel PETSc vector to local indices of subdomain
148  ierr = VecScatterCreate( solution_, from, locSolVec_, subdomainIndexSet, &VSpetscToSubScatter_ ); CHKERRV( ierr );
149  ierr = ISDestroy( &subdomainIndexSet ); CHKERRV( ierr );
150  ierr = ISDestroy( &from ); CHKERRV( ierr );
151 
152  //double * locSolVecArray;
153  //ierr = VecGetArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
154  //std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
155  //ierr = VecRestoreArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
156 
157  // scatter local solutions back to global one
158  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
159  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
160 #endif // HAVE_BDDCML
161 }
162 
163 void LinSys_BDDC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
164 {
165 #ifdef HAVE_BDDCML
166  namespace ublas = boost::numeric::ublas;
167 
168  std::vector< unsigned > myRows( nrow );
169  std::vector< unsigned > myCols( ncol );
170  ublas::matrix< double > mat( nrow, ncol );
171 
172  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
173  std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
174 
175  for ( int i = 0; i < nrow; i++ ) {
176  for ( int j = 0; j < ncol; j++ ) {
177  mat( i, j ) = vals[i*ncol + j];
178  }
179  }
180  if (swap_sign_) {
181  mat = -mat;
182  }
183 
184  bddcml_ -> insertToMatrix( mat, myRows, myCols );
185 #endif // HAVE_BDDCML
186 }
187 
188 void LinSys_BDDC::rhs_set_values( int nrow, int *rows, double *vals)
189 {
190 #ifdef HAVE_BDDCML
191  namespace ublas = boost::numeric::ublas;
192 
193  std::vector< unsigned > myRows( nrow );
194  ublas::vector< double > vec( nrow );
195 
196  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
197 
198  for ( int i = 0; i < nrow; i++ ) {
199  vec( i ) = vals[i];
200  }
201  if (swap_sign_) {
202  vec = -vec;
203  }
204 
205  bddcml_ -> insertToRhs( vec, myRows );
206 #endif // HAVE_BDDCML
207 }
208 
209 void LinSys_BDDC::diagonal_weights_set_value( int global_index, double value )
210 {
211 #ifdef HAVE_BDDCML
212  bddcml_ -> insertToDiagonalWeights( global_index, value );
213 #endif // HAVE_BDDCML
214 }
215 
217 {
218 #ifdef HAVE_BDDCML
219  bddcml_ -> clearMatrix( );
220 #endif // HAVE_BDDCML
221  return 0;
222 }
223 
225 {
226 #ifdef HAVE_BDDCML
227  bddcml_ -> clearRhs( );
228 #endif // HAVE_BDDCML
229  return 0;
230 }
231 
233 {
234 #ifdef HAVE_BDDCML
235  bddcml_ -> finishMatAssembly( );
236 #endif // HAVE_BDDCML
237 }
238 
239 void LinSys_BDDC::apply_constrains( double scalar )
240 {
241 #ifdef HAVE_BDDCML
242  bddcml_ -> applyConstraints( constraints_, 1., scalar );
243 #endif // HAVE_BDDCML
244 }
245 
246 int LinSys_BDDC::solve() // ! params are not currently used
247 {
248 #ifdef HAVE_BDDCML
249  std::vector<int> * numSubAtLevels = NULL; //!< number of subdomains at levels
250 
252 
253  DBGMSG("BDDCML converged reason: %d ( 0 means OK ) \n", bddcml_ -> giveConvergedReason() );
254  DBGMSG("BDDCML converged in %d iterations. \n", bddcml_ -> giveNumIterations() );
255  DBGMSG("BDDCML estimated condition number is %f \n", bddcml_ -> giveCondNumber() );
256 
257  // download local solution
258  bddcml_ -> giveSolution( isngn_, locSolution_ );
259 
260 
261  double * locSolVecArray;
262  VecGetArray( locSolVec_, &locSolVecArray );
263  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
264  VecRestoreArray( locSolVec_, &locSolVecArray );
265 
266  // scatter local solutions back to global one
267  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
268  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
269 
270  // upper bound on the residual error
272 
273  return bddcml_ -> giveConvergedReason();
274 #else
275  return 0;
276 #endif // HAVE_BDDCML
277 
278 }
279 
281 {
282  // common values
283  r_tol_ = in_rec.val<double>("r_tol");
284  max_it_ = in_rec.val<int>("max_it");
285 
286  // BDDCML specific parameters
287  max_nondecr_it_ = in_rec.val<int>("max_nondecr_it");
288  number_of_levels_ = in_rec.val<int>("number_of_levels");
289  use_adaptive_bddc_ = in_rec.val<bool>("use_adaptive_bddc");
290  bddcml_verbosity_level_ = in_rec.val<int>("bddcml_verbosity_level");
291 }
292 
294 {
295 #ifdef HAVE_BDDCML
296  isngn_.clear();
297  locSolution_.clear();
298 
299  PetscErrorCode ierr;
300  ierr = VecDestroy( &locSolVec_ ); CHKERRV( ierr );
301 
302  ierr = VecScatterDestroy( &VSpetscToSubScatter_ ); CHKERRV( ierr );
303 
304  delete bddcml_;
305 #endif // HAVE_BDDCML
306 };
307 
308 // construct global solution
309 //void LinSys_BDDC::gatherSolution_( )
310 //{
311 //#ifdef HAVE_BDDCML
312 // int ierr;
313 //
314 // // merge solution on root
315 // int rank;
316 // MPI_Comm_rank( comm_, &rank );
317 // int nProc;
318 // MPI_Comm_size( comm_, &nProc );
319 //
320 // globalSolution_.resize( size_ );
321 // std::vector<double> locSolutionNeib;
322 // if ( rank == 0 ) {
323 // // merge my own data
324 // for ( unsigned int i = 0; i < isngn_.size(); i++ ) {
325 // int ind = isngn_[i];
326 // globalSolution_[ind] = locSolution_[i];
327 // }
328 // for ( int iProc = 1; iProc < nProc; iProc++ ) {
329 // // receive length
330 // int length;
331 // MPI_Status status;
332 // ierr = MPI_Recv( &length, 1, MPI_INT, iProc, iProc, comm_, &status );
333 //
334 // // receive indices
335 // std::vector<int> inds( length );
336 // ierr = MPI_Recv( &(inds[0]), length, MPI_INT, iProc, iProc, comm_, &status );
337 //
338 // // receive local solution
339 // locSolutionNeib.resize( length );
340 // ierr = MPI_Recv( &(locSolutionNeib[0]), length, MPI_DOUBLE, iProc, iProc, comm_, &status );
341 //
342 // // merge data
343 // for ( int i = 0; i < length; i++ ) {
344 // int ind = inds[i];
345 // globalSolution_[ind] = locSolutionNeib[i];
346 // }
347 // }
348 // }
349 // else {
350 // // send my solution to root
351 // int length = isngn_.size();
352 // ierr = MPI_Send( &length, 1, MPI_INT, 0, rank, comm_ );
353 // ierr = MPI_Send( &(isngn_[0]), length, MPI_INT, 0, rank, comm_ );
354 // ierr = MPI_Send( &(locSolution_[0]), length, MPI_DOUBLE, 0, rank, comm_ );
355 // }
356 // // broadcast global solution from root
357 // ierr = MPI_Bcast( &(globalSolution_[0]), globalSolution_.size(), MPI_DOUBLE, 0, comm_ );
358 //#endif // HAVE_BDDCML
359 //}
360 
362 {
363  double bnorm=0.0;
364  VecNorm(locSolVec_, NORM_2, &bnorm);
365 
366  return max(a_tol_, r_tol_*bnorm);
367 }
368 
369 
int bddcml_verbosity_level_
Definition: linsys_BDDC.hh:102
SetValuesMode status_
Set value status of the linear system.
Definition: linsys.hh:570
double normRhs()
Get norm of right-hand side.
#define DBGMSG(...)
Definition: global_defs.h:196
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)
Definition: linsys_BDDC.cc:163
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:224
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
Wrappers for linear systems based on MPIAIJ and MATIS format.
ConstraintVec_ constraints_
Definition: linsys.hh:591
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:169
enum matrixTypeEnum MatrixType
unsigned size_
global number of matrix rows, i.e. problem size
Definition: linsys.hh:573
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:112
void rhs_set_values(int nrow, int *rows, double *vals)
Definition: linsys_BDDC.cc:188
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:569
static Input::Type::AbstractRecord input_type
Definition: linsys.hh:90
int max_it_
Definition: linsys.hh:567
Vec locSolVec_
local solution PETSc vector - sequential
Definition: linsys_BDDC.hh:111
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:100
int max_nondecr_it_
parameters expected from input file:
Definition: linsys_BDDC.hh:99
double a_tol_
Definition: linsys.hh:566
void diagonal_weights_set_value(int global_index, double value)
Definition: linsys_BDDC.cc:209
double residual_norm_
local solution array pointing into Vec solution_
Definition: linsys.hh:589
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:361
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:565
Bddcml_ * bddcml_
BDDCML wrapper.
Definition: linsys_BDDC.hh:115
PetscErrorCode mat_zero_entries() override
Definition: linsys_BDDC.cc:216
MPI_Comm get_comm() const
Returns communicator.
Abstract linear system class.
void apply_constrains(double scalar=1.)
Definition: linsys_BDDC.cc:239
Record type proxy class.
Definition: type_record.hh:161
void finish_assembly()
Definition: linsys_BDDC.cc:232
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
void set_from_input(const Input::Record in_rec)
Definition: linsys_BDDC.cc:280
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:585
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:386