Flow123d
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 namespace it = Input::Type;
45 
46 
47 it::Record LinSys_BDDC::input_type = it::Record("Bddc", "Solver setting.")
49  .declare_key("max_nondecr_it", it::Integer(0), it::Default("30"),
50  "Maximum number of iterations of the linear solver with non-decreasing residual.")
51  .declare_key("number_of_levels", it::Integer(0), it::Default("2"),
52  "Number of levels in the multilevel method (=2 for the standard BDDC).")
53  .declare_key("use_adaptive_bddc", it::Bool(), it::Default("false"),
54  "Use adaptive selection of constraints in BDDCML.")
55  .declare_key("bddcml_verbosity_level", it::Integer(0,2), it::Default("0"),
56  "Level of verbosity of the BDDCML library: 0 - no output, 1 - mild output, 2 - detailed output.");
57 
58 
59 
60 LinSys_BDDC::LinSys_BDDC( const unsigned numDofsSub,
61  const Distribution * rows_ds,
62  const int matrixTypeInt,
63  const int numSubLoc,
64  const bool swap_sign )
65  : LinSys( rows_ds ),
66  swap_sign_(swap_sign)
67 {
68 #ifdef HAVE_BDDCML
69  // set type
70  //type = LinSys::BDDC;
71 
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  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  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  //std::cout << "IDX: \n";
139  //std::copy( idx.begin(), idx.end(), std::ostream_iterator<PetscInt>( std::cout, " " ) );
140 
141  ISLocalToGlobalMapping subdomainMapping;
142  ierr = ISLocalToGlobalMappingCreate( comm_, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
143 
144  IS subdomainIndexSet;
145  IS from;
146  ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
147  ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
148 
149  //ierr = ISCreateGeneral( comm_, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
150  //ISView( subdomainIndexSet, PETSC_VIEWER_STDOUT_WORLD);
151 
152 
153  // create scatter from parallel PETSc vector to local indices of subdomain
154  ierr = VecScatterCreate( solution_, from, locSolVec_, subdomainIndexSet, &VSpetscToSubScatter_ ); CHKERRV( ierr );
155  ierr = ISDestroy( &subdomainIndexSet ); CHKERRV( ierr );
156  ierr = ISDestroy( &from ); CHKERRV( ierr );
157 
158  //VecScatterView(VSpetscToSubScatter_,PETSC_VIEWER_STDOUT_SELF);
159 
160  double * locSolVecArray;
161  ierr = VecGetArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
162  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
163  ierr = VecRestoreArray( locSolVec_, &locSolVecArray ); CHKERRV( ierr );
164 
165  // scatter local solutions back to global one
166  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
167  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
168 #endif // HAVE_BDDCML
169 }
170 
172 {
173 #ifdef HAVE_BDDCML
174  // simply pass the data to BDDCML solver
175  bddcml_ -> loadDiagonal( diag );
176 #endif // HAVE_BDDCML
177 }
178 
179 void LinSys_BDDC::mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals )
180 {
181 #ifdef HAVE_BDDCML
182  namespace ublas = boost::numeric::ublas;
183 
184  std::vector< unsigned > myRows( nrow );
185  std::vector< unsigned > myCols( ncol );
186  ublas::matrix< double > mat( nrow, ncol );
187 
188  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
189  std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
190 
191  for ( unsigned i = 0; i < nrow; i++ ) {
192  for ( unsigned j = 0; j < ncol; j++ ) {
193  mat( i, j ) = vals[i*ncol + j];
194  }
195  }
196  if (swap_sign_) {
197  mat = -mat;
198  }
199 
200  bddcml_ -> insertToMatrix( mat, myRows, myCols );
201 #endif // HAVE_BDDCML
202 }
203 
204 void LinSys_BDDC::rhs_set_values( int nrow, int *rows, double *vals)
205 {
206 #ifdef HAVE_BDDCML
207  namespace ublas = boost::numeric::ublas;
208 
209  std::vector< unsigned > myRows( nrow );
210  ublas::vector< double > vec( nrow );
211 
212  std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
213 
214  for ( unsigned i = 0; i < nrow; i++ ) {
215  vec( i ) = vals[i];
216  }
217  if (swap_sign_) {
218  vec = -vec;
219  }
220 
221  bddcml_ -> insertToRhs( vec, myRows );
222 #endif // HAVE_BDDCML
223 }
224 
226 {
227 #ifdef HAVE_BDDCML
228  bddcml_ -> finishMatAssembly( );
229 #endif // HAVE_BDDCML
230 }
231 
232 void LinSys_BDDC::apply_constrains( double scalar )
233 {
234 #ifdef HAVE_BDDCML
235  bddcml_ -> applyConstraints( constraints_, 1., scalar );
236 #endif // HAVE_BDDCML
237 }
238 
239 int LinSys_BDDC::solve() // ! params are not currently used
240 {
241 #ifdef HAVE_BDDCML
242  std::vector<int> * numSubAtLevels = NULL; //!< number of subdomains at levels
243 
245 
246  DBGMSG("BDDCML converged reason: %d ( 0 means OK ) \n", bddcml_ -> giveConvergedReason() );
247  DBGMSG("BDDCML converged in %d iterations. \n", bddcml_ -> giveNumIterations() );
248  DBGMSG("BDDCML estimated condition number is %f \n", bddcml_ -> giveCondNumber() );
249 
250  // download local solution
251  bddcml_ -> giveSolution( isngn_, locSolution_ );
252 
253 
254  double * locSolVecArray;
255  PetscErrorCode ierr;
256  ierr = VecGetArray( locSolVec_, &locSolVecArray );
257  std::copy( locSolution_.begin(), locSolution_.end(), locSolVecArray );
258  ierr = VecRestoreArray( locSolVec_, &locSolVecArray );
259 
260  // scatter local solutions back to global one
261  VecScatterBegin( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
262  VecScatterEnd( VSpetscToSubScatter_, locSolVec_, solution_, INSERT_VALUES, SCATTER_REVERSE );
263 
264  // upper bound on the residual error
266 
267  return bddcml_ -> giveConvergedReason();
268 #else
269  return 0;
270 #endif // HAVE_BDDCML
271 
272 }
273 
275 {
276 #ifdef HAVE_BDDCML
277  this -> gatherSolution_( );
278  globalSolution.resize( globalSolution_.size( ) );
279  std::copy( globalSolution_.begin(), globalSolution_.end(), globalSolution.begin() );
280 #endif // HAVE_BDDCML
281 }
282 
284 {
285 #ifdef HAVE_BDDCML
286  globalSolution_.resize( globalSolution.size( ) );
287  std::copy( globalSolution.begin(), globalSolution.end(), globalSolution_.begin() );
288 #endif // HAVE_BDDCML
289 }
290 
292 {
293  // common values
294  r_tol_ = in_rec.val<double>("r_tol");
295  max_it_ = in_rec.val<int>("max_it");
296 
297  // BDDCML specific parameters
298  max_nondecr_it_ = in_rec.val<int>("max_nondecr_it");
299  number_of_levels_ = in_rec.val<int>("number_of_levels");
300  use_adaptive_bddc_ = in_rec.val<bool>("use_adaptive_bddc");
301  bddcml_verbosity_level_ = in_rec.val<int>("bddcml_verbosity_level");
302 }
303 
305 {
306 #ifdef HAVE_BDDCML
307  isngn_.clear();
308  locSolution_.clear();
309 
310  PetscErrorCode ierr;
311  ierr = VecDestroy( &locSolVec_ ); CHKERRV( ierr );
312 
313  ierr = VecScatterDestroy( &VSpetscToSubScatter_ ); CHKERRV( ierr );
314 
315  delete bddcml_;
316 #endif // HAVE_BDDCML
317 };
318 
319 // construct global solution
321 {
322 #ifdef HAVE_BDDCML
323  int ierr;
324 
325  // merge solution on root
326  int rank;
327  MPI_Comm_rank( comm_, &rank );
328  int nProc;
329  MPI_Comm_size( comm_, &nProc );
330 
331  globalSolution_.resize( size_ );
332  std::vector<double> locSolutionNeib;
333  if ( rank == 0 ) {
334  // merge my own data
335  for ( int i = 0; i < isngn_.size(); i++ ) {
336  int ind = isngn_[i];
337  globalSolution_[ind] = locSolution_[i];
338  }
339  for ( int iProc = 1; iProc < nProc; iProc++ ) {
340  // receive length
341  int length;
343  ierr = MPI_Recv( &length, 1, MPI_INT, iProc, iProc, comm_, &status );
344 
345  // receive indices
346  std::vector<int> inds( length );
347  ierr = MPI_Recv( &(inds[0]), length, MPI_INT, iProc, iProc, comm_, &status );
348 
349  // receive local solution
350  locSolutionNeib.resize( length );
351  ierr = MPI_Recv( &(locSolutionNeib[0]), length, MPI_DOUBLE, iProc, iProc, comm_, &status );
352 
353  // merge data
354  for ( int i = 0; i < length; i++ ) {
355  int ind = inds[i];
356  globalSolution_[ind] = locSolutionNeib[i];
357  }
358  }
359  }
360  else {
361  // send my solution to root
362  int length = isngn_.size();
363  ierr = MPI_Send( &length, 1, MPI_INT, 0, rank, comm_ );
364  ierr = MPI_Send( &(isngn_[0]), length, MPI_INT, 0, rank, comm_ );
365  ierr = MPI_Send( &(locSolution_[0]), length, MPI_DOUBLE, 0, rank, comm_ );
366  }
367  // broadcast global solution from root
368  ierr = MPI_Bcast( &(globalSolution_[0]), globalSolution_.size(), MPI_DOUBLE, 0, comm_ );
369 #endif // HAVE_BDDCML
370 }
371 
373 {
374  double bnorm=0.0;
375  VecNorm(locSolVec_, NORM_2, &bnorm);
376 
377  return max(a_tol_, r_tol_*bnorm);
378 }
379 
380