Flow123d  JS_before_hm-2196-g9bdbf60c0
bddcml_wrapper.cc
Go to the documentation of this file.
1 //! @author Jakub Sistek
2 //! @date 7/2012
3 
4 #include "la/bddcml_wrapper.hh"
5 
6 #ifdef FLOW123D_HAVE_BDDCML
7 
8 //------------------------------------------------------------------------------
9 /** Constructor for the BDDCML solver object
10  */
11 la::BddcmlWrapper::BddcmlWrapper( const unsigned numDofs,
12  const unsigned numDofsSub,
13  const MatrixType matrixType,
14  const MPI_Comm comm,
15  int numSubLoc )
16  : numDofs_( static_cast<int>( numDofs ) ),
17  numDofsSub_( static_cast<int>( numDofsSub ) ),
18  matrixType_( matrixType ),
19  comm_ (comm),
20  numSubLoc_ ( numSubLoc )
21 {
22 
23  // Orient in the communicator
24  MPI_Comm_size( comm_, &nProc_ );
25  MPI_Comm_rank( comm_, &rank_ );
26 
27  // prepare RHS vector
28  rhsVec_.resize( numDofsSub_, 0. );
29 
30  // prepare arrays to mark boundary conditions
31  ifix_.resize( numDofsSub_, 0 );
32  fixv_.resize( numDofsSub_, 0. );
33 
34  // initialize subdomain solution vector
35  sol_.resize( numDofsSub_, 0. );
36  normSol_ = 0.;
37  numIter_ = 0;
38  convergedReason_ = 0;
39  condNumber_ = 0.;
40 
41  // initialize state of mesh and matrix
42  meshLoaded_ = false;
43  isMatAssembled_ = false;
44 }
45 //------------------------------------------------------------------------------
46 /** Destructor frees the BDDC solver structures
47  */
49 {
50  // clear RHS
51  rhsVec_.clear();
52 
53  // clear BC
54  ifix_.clear();
55  fixv_.clear();
56 
57  // clear solution
58  sol_.clear();
59 
60  // clear mesh data
61  inet_.clear();
62  nnet_.clear();
63  nndf_.clear();
64  xyz_.clear();
65  isngn_.clear();
66  isvgvn_.clear();
67  isegn_.clear();
68 
69  // clear global2local map
70  global2LocalDofMap_.clear();
71 
72  // clear triplet
73  coo_.clear();
74 
75 }
76 
77 //------------------------------------------------------------------------------
78 /** Routine for loading raw mesh data into the solver - for cases of strange meshes,
79  * where these are not Mesh or Grid objects, user can create own raw description
80  * to exploit the flexibility of mesh format underlaying BDDCML.
81  */
82 void la::BddcmlWrapper::loadRawMesh( const int nDim, const int numNodes,
83  const std::vector<int> & inet,
84  const std::vector<int> & nnet,
85  const std::vector<int> & nndf,
86  const std::vector<int> & isegn,
87  const std::vector<int> & isngn,
88  const std::vector<int> & isvgvn,
89  const std::vector<double> & xyz,
90  const std::vector<double> & element_data,
91  const int meshDim )
92 {
93  // space dimension
94  nDim_ = nDim;
95 
96  // set topological dimension of mesh
97  if ( meshDim == 0 ) { // if not specified, use space dimension
98  meshDim_ = nDim_;
99  }
100  else {
101  meshDim_ = static_cast<int>( meshDim ); // if specified, use it
102  }
103 
104  // fill internal variables about mesh
105  numElemSub_ = isegn.size( );
106  numNodesSub_ = isngn.size( );
107 
108  // global number of nodes
109  numNodes_ = numNodes;
110 
111  std::vector<int> procElemStarts( nProc_ + 1 );
112  int numElemSubInt = static_cast<int>( numElemSub_ );
113  MPI_Allgather( &numElemSubInt, 1, MPI_INT,
114  &(procElemStarts[0]), 1, MPI_INT, comm_ );
115  // now procElemStarts_ contains counts of local elements at each process
116  // change it to starts
117  for ( unsigned i = 1; i < static_cast<unsigned>( nProc_ ); i++ )
118  procElemStarts[i] = procElemStarts[i-1] + procElemStarts[i];
119  // shift it one back
120  for ( unsigned i = nProc_; i > 0; i-- )
121  procElemStarts[i] = procElemStarts[i-1];
122  // start from zero
123  procElemStarts[0] = 0;
124 
125  // find global number of elements - assumes nonoverlapping division
126  numElem_ = procElemStarts[ nProc_ ];
127 
128  // check sizes of arrays
129  OLD_ASSERT( std::accumulate( nnet.begin(), nnet.end(), 0 ) == (int)(inet.size()),
130  "array inet size mismatch \n " );
131  OLD_ASSERT( std::accumulate( nndf.begin(), nndf.end(), 0 ) == numDofsSub_,
132  "array nndf content mismatch: %d %d \n ", std::accumulate( nndf.begin(), nndf.end(), 0 ), numDofsSub_ );
133  OLD_ASSERT( static_cast<int>(isvgvn.size()) == numDofsSub_,
134  "array isvgvn size mismatch \n " );
135  OLD_ASSERT( static_cast<int>(nnet.size()) == numElemSub_,
136  "array nnet size mismatch \n " );
137  OLD_ASSERT( static_cast<int>(nndf.size()) == numNodesSub_,
138  "array nndf size mismatch \n " );
139  OLD_ASSERT( static_cast<int>(xyz.size()) == numNodesSub_ * nDim_,
140  "array xyz size mismatch: %d %d \n ", static_cast<int>(xyz.size()), numNodesSub_ * nDim_ );
141 
142  // Simply copy input data into the private object equivalents
143  inet_.resize( inet.size() );
144  std::copy( inet.begin(), inet.end(), inet_.begin() );
145  nnet_.resize( nnet.size() );
146  std::copy( nnet.begin(), nnet.end(), nnet_.begin() );
147  nndf_.resize( nndf.size() );
148  std::copy( nndf.begin(), nndf.end(), nndf_.begin() );
149  isegn_.resize( isegn.size() );
150  std::copy( isegn.begin(), isegn.end(), isegn_.begin() );
151  isngn_.resize( isngn.size() );
152  std::copy( isngn.begin(), isngn.end(), isngn_.begin() );
153  isvgvn_.resize( isvgvn.size() );
154  std::copy( isvgvn.begin(), isvgvn.end(), isvgvn_.begin() );
155  xyz_.resize( xyz.size() );
156  std::copy( xyz.begin(), xyz.end(), xyz_.begin() );
157  element_data_.resize( element_data.size() );
158  std::copy( element_data.begin(), element_data.end(), element_data_.begin() );
159 
160  // prepare auxiliary map for renumbering dofs in the global array
161  for ( unsigned ind = 0; ind < isvgvn_.size(); ++ind ) {
162  global2LocalDofMap_.insert( std::make_pair( static_cast<unsigned>( isvgvn_[ind] ), ind ) );
163  }
164 
165  // change the state
166  meshLoaded_ = true;
167 
168  return;
169 }
170 
171 //------------------------------------------------------------------------------
172 /** Routine for preparing MATRIX assembly. It reserves memory in triplet.
173  */
174 void la::BddcmlWrapper::prepareMatAssembly( unsigned numElements, unsigned elMatSize )
175 {
176  unsigned length = numElements * elMatSize;
177  coo_.prepareAssembly( length );
178 
179  return;
180 }
181 //------------------------------------------------------------------------------
182 /** Insert element stiffness matrix into global system matrix
183  */
184 void la::BddcmlWrapper::insertToMatrix( const SubMat_ & subMat,
185  const VecUint_ & rowIndices,
186  const VecUint_ & colIndices )
187 {
188  // copy dense matrix into triplet
189  for ( unsigned i = 0; i < rowIndices.size(); i++) {
190  for ( unsigned j = 0; j < colIndices.size(); j++) {
191 
192  // insert value
193  if (subMat(i,j) != 0.0)
194  coo_.insert( rowIndices[i], colIndices[j], subMat(i,j) );
195  }
196  }
197 
198  // set flag
199  isMatAssembled_ = false;
200 
201  return;
202 }
203 
204 //------------------------------------------------------------------------------
205 /** Routine for triplet finalization.
206  */
208 {
209  // do nothing if already called
210  if ( isMatAssembled_ ) return;
211 
212  // finish assembly of triplet
213  coo_.finishAssembly( );
214 
215  // set flag
216  isMatAssembled_ = true;
217 
218  return;
219 }
220 //------------------------------------------------------------------------------
221 /** Insert subvector to the system's RHS vector
222  */
223 void la::BddcmlWrapper::insertToRhs( const SubVec_ & subVec,
224  const VecUint_ & dofIndices )
225 {
226  for ( unsigned i = 0; i < dofIndices.size( ); i ++ ) {
227 
228  // map global dof index to subdomain dof
229  Global2LocalMap_::iterator pos = global2LocalDofMap_.find( dofIndices[i] );
230  /*
231  OLD_ASSERT( pos != global2LocalDofMap_.end(),
232  "Cannot remap index %d to local indices in right-hand side,"
233  " perhaps missing call to loadMesh() member function. \n ", dofIndices[i] );
234  */
235  if (pos == global2LocalDofMap_.end()) {
236  DebugOut().every_proc() << "RHS, Missing map to local for idx: " << dofIndices[i];
237  }
238  const unsigned indLoc = pos -> second;
239 
240  rhsVec_[ indLoc ] += subVec[i];
241  }
242  return;
243 }
244 
245 //------------------------------------------------------------------------------
246 /** Apply boundary conditions
247  * Since BDDCML handles constraints, just fill proper BC arrays.
248  * Scalar is not used, it is kept for compatibility with constraint functors.
249  */
250 void la::BddcmlWrapper::applyConstraints( ConstraintVec_ & constraints,
251  const double factor, FMT_UNUSED const double scalar )
252 {
253  // Constraint iterators
254  ConstraintVec_::const_iterator cIter = constraints.begin( );
255  ConstraintVec_::const_iterator cEnd = constraints.end( );
256  // collect global dof indices and the correpsonding values
257  for ( ; cIter != cEnd; ++cIter ) {
258  unsigned globalDof = cIter -> first;
259  double value = (cIter -> second );
260 
261  // map global dof index to subdomain dof
262  Global2LocalMap_::iterator pos = global2LocalDofMap_.find( globalDof );
263  if ( pos != global2LocalDofMap_.end() ) {
264  // I contain the dof
265  unsigned indLoc = pos -> second;
266 
267  ifix_[ indLoc ] = 1;
268  fixv_[ indLoc ] = factor * value;
269  }
270  }
271 
272  return;
273 }
274 //------------------------------------------------------------------------------
275 /** Simply set the dof with index 'index' to zero
276  */
277 void la::BddcmlWrapper::fixDOF( const unsigned index, FMT_UNUSED const double scalar )
278 {
280  thisConstraint.push_back( std::make_pair( index, 0. ) );
281  this -> applyConstraints( thisConstraint, 1., scalar );
282  thisConstraint.clear();
283  return;
284 }
285 
286 //------------------------------------------------------------------------------
287 /** Insert elements on diagonal for weights.
288  */
289 void la::BddcmlWrapper::insertToDiagonalWeights( const int & index,
290  const double & value )
291 {
292  // insert value
293  diagWeightsCoo_.insert( index, index, value );
294 
295  // set flag
296  isDiagAssembled_ = false;
297 
298  return;
299 }
300 
301 //------------------------------------------------------------------------------
302 /** Routine for diagonal finalization.
303  */
305 {
306  // do nothing if already called
307  if ( isDiagAssembled_ ) return;
308 
309  // finish assembly of triplet
310  diagWeightsCoo_.finishAssembly( );
311 
312  // set flag
313  isDiagAssembled_ = true;
314 
315  return;
316 }
317 
318 //------------------------------------------------------------------------------
319 /** Solve the system by BDDCML
320  */
321 void la::BddcmlWrapper::solveSystem( double tol, int numLevels, std::vector<int> * numSubAtLevels,
322  int verboseLevel, int maxIt, int ndecrMax, bool use_adaptive )
323 {
324 
325  // check that mesh was loaded
326  OLD_ASSERT( meshLoaded_ ,
327  "Subdomain mesh was not loaded, perhaps missing call to loadMesh() member function. \n " );
328 
329  //============= initialize BDDC
330  int numSub; //!< number of subdomains on first level
331  // Number of subdomains at first level - if not given, deduce it from MPI communicator
332  MPI_Allreduce ( &numSubLoc_, &numSub, 1, MPI_INT, MPI_SUM, comm_ );
333 
334  // Number of subdomains at levels
335  std::vector<int> numSubLev( numLevels ); //!< numbers of subdomains on levels ( should be monotonically decreasing, e.g. for 3 levels 1024, 32, 1 )
336  if ( numSubAtLevels != NULL ) {
337  // numbers are given
338 
339  OLD_ASSERT( numLevels == static_cast<int>((*numSubAtLevels).size()),
340  "Inconsistent size of numbers of subdomains at levels: %d %d \n",
341  numLevels, static_cast<int>((*numSubAtLevels).size()) );
342  OLD_ASSERT( (*numSubAtLevels)[0] == numSub,
343  "Inconsistent number of subdomains at first level: %d %d \n",
344  (*numSubAtLevels)[0], numSub );
345 
346  std::copy( (*numSubAtLevels).begin(), (*numSubAtLevels).end(),
347  numSubLev.begin() );
348  }
349  else {
350  // assume 2-level BDDC
351  // generate default vector for two levels
352  if ( numLevels == 2 ) {
353  numSubLev[0] = numSub;
354  numSubLev[1] = 1;
355  }
356  else if ( numLevels > 2 ) {
357  double coarsening = pow( numSub, 1. / double(numLevels - 1) );
358  // prescribe number of subdomains on levels so that coarsening is fixed between levels
359  numSubLev[0] = numSub;
360  for (int i = 1; i <= numLevels-2; i++) {
361  int ir = numLevels - i - 1;
362  numSubLev[i] = (int) pow( coarsening, ir );
363  if ( numSubLev[i]%2 != 0 ) {
364  numSubLev[i] = numSubLev[i] + 1;
365  }
366  }
367  numSubLev[numLevels-1] = 1;
368  }
369  else {
370  OLD_ASSERT( false, "Missing numbers of subdomains at levels. \n" );
371  }
372  }
373  //std::cout << "numSubLev ";
374  //std::copy( numSubLev.begin(), numSubLev.end(), std::ostream_iterator<int>( std::cout, " " ) );
375  //std::cout << std::endl;
376 
377  int lnumSubLev = numSubLev.size();
378  int commInt = MPI_Comm_c2f( comm_ ); // translate C MPI communicator hanlder to Fortran integer
379  int numBase = 0; // indexing starts with 0 for C++
380  int justDirectSolve = 0;
381 
382  bddcml_init_c( &numLevels, &(numSubLev[0]), &lnumSubLev, &numSubLoc_, &commInt, &verboseLevel, &numBase, &justDirectSolve );
383 
384 
385  //============= uploading data for subdomain
386  // at the moment, number of subdomains equals number of processors
387  int iSub = rank_;
388 
389  int linet = inet_.size();
390  int lnnet = nnet_.size();
391  int lnndf = nndf_.size();
392  int lisngn = isngn_.size();
393  int lisvgvn = isvgvn_.size();
394  int lisegn = isegn_.size();
395 
396  int lxyz1 = numNodesSub_;
397  int lxyz2 = nDim_;
398 
399  int lifix = ifix_.size();
400  int lfixv = fixv_.size();
401 
402  int lrhsVec = rhsVec_.size();
403  int isRhsCompleteInt = 0; // since only local subassembly of right-hand side was performed, it does not contain
404  // repeated entries
405 
406  std::vector<int> i_sparse;
407  std::vector<int> j_sparse;
408  std::vector<double> a_sparse;
409 
410  // specify if matrix is symmetric, if yes, use only upper triangle
411  bool onlyUpperTriangle;
412  if ( matrixType_ == GENERAL ) {
413  onlyUpperTriangle = false;
414  }
415  else {
416  onlyUpperTriangle = true;
417  }
418 
419  // copy out arrays from triplet
420  this -> finishMatAssembly();
421  coo_.extractArrays( i_sparse, j_sparse, a_sparse, onlyUpperTriangle );
422  coo_.clear();
423 
424  // map global row and column indices to local subdomain indices
425  int indRow = -1;
426  int indCol = -1;
427  int indRowLoc = -1;
428  int indColLoc = -1;
429  for ( unsigned inz = 0; inz < a_sparse.size(); inz++ ) {
430  if ( i_sparse[inz] != indRow ) {
431  indRow = i_sparse[inz];
432  Global2LocalMap_::iterator pos = global2LocalDofMap_.find( static_cast<unsigned> ( indRow ) );
433  OLD_ASSERT( pos != global2LocalDofMap_.end(),
434  "Cannot remap index %d to local indices. \n ", indRow );
435  indRowLoc = static_cast<int> ( pos -> second );
436  }
437  if ( j_sparse[inz] != indCol ) {
438  indCol = j_sparse[inz];
439  Global2LocalMap_::iterator pos = global2LocalDofMap_.find( static_cast<unsigned> ( indCol ) );
440  OLD_ASSERT( pos != global2LocalDofMap_.end(),
441  "Cannot remap index %d to local indices. \n ", indCol );
442  indColLoc = static_cast<int> ( pos -> second );
443  }
444 
445  i_sparse[ inz ] = indRowLoc;
446  j_sparse[ inz ] = indColLoc;
447 
448  //std::cout << " row: " << i_sparse[ inz ] << " col: " << j_sparse[ inz ] << std::endl;
449  }
450  // due to remap, it is not guaranteed that the triplet is still correctly ordered
451  int isMatAssembledInt = 0;
452 
453  int la = a_sparse.size();
454 
455  // diagonal weights for BDDC loaded by user
456  OLD_ASSERT( diagWeightsCoo_.nnz() > 0,
457  "It appears that diagonal weights for BDDC are not loaded. This is currently mandatory. \n " );
458 
459  std::vector<int> i_diag_sparse;
460  std::vector<int> j_diag_sparse;
461  std::vector<double> diag_sparse;
462 
463  // copy out arrays from triplet
464  this -> finishDiagAssembly();
465  diagWeightsCoo_.extractArrays( i_diag_sparse, j_diag_sparse, diag_sparse );
466  diagWeightsCoo_.clear();
467 
468  // map global row indices to local subdomain indices
469  indRow = -1;
470  indRowLoc = -1;
471  for ( unsigned inz = 0; inz < diag_sparse.size(); inz++ ) {
472  indRow = i_diag_sparse[inz];
473  Global2LocalMap_::iterator pos = global2LocalDofMap_.find( static_cast<unsigned> ( indRow ) );
474  OLD_ASSERT( pos != global2LocalDofMap_.end(),
475  "Cannot remap index %d to local indices. \n ", indRow );
476  indRowLoc = static_cast<int> ( pos -> second );
477 
478  i_diag_sparse[ inz ] = indRowLoc;
479  }
480 
481  int lsub_diagonal = numDofsSub_;
482  OLD_ASSERT( lsub_diagonal == static_cast<int>(diag_sparse.size()),
483  "Array length mismatch: %d %d . \n ", lsub_diagonal, static_cast<int>(diag_sparse.size()) );
484 
485  std::vector<double> sub_diagonal( lsub_diagonal, -1. );
486  // permute the vector according to subdomain indexing
487  for ( int i = 0; i < lsub_diagonal; i++ ){
488  indRowLoc = i_diag_sparse[i];
489  double value = diag_sparse[i];
490  sub_diagonal[indRowLoc] = value;
491  }
492  OLD_ASSERT( std::find(sub_diagonal.begin(), sub_diagonal.end(), -1. ) == sub_diagonal.end(),
493  "There are missing entries in the diagonal of weights. \n " );
494 
495  // remove const attribute
496  int numDofsInt = static_cast<int> ( numDofs_ );
497  int numDofsSubInt = static_cast<int> ( numDofsSub_ );
498 
499  int matrixTypeInt;
500  if ( matrixType_ == GENERAL ) matrixTypeInt = 0 ;
501  else if ( matrixType_ == SPD ) matrixTypeInt = 1 ;
502  else if ( matrixType_ == SYMMETRICGENERAL ) matrixTypeInt = 2 ;
503  else if ( matrixType_ == SPD_VIA_SYMMETRICGENERAL ) matrixTypeInt = 2 ;
504  else {OLD_ASSERT( false, "Illegal matrixType \n " );}
505 
506  // download local solution
507  int lsol = sol_.size();
508 
509  // user constraints are not used
510  std::vector<double> userConstraints(1,0.);
511  int lUserConstraints1 = 0;
512  int lUserConstraints2 = 0;
513 
514  int lelement_data1 = 1;
515  int lelement_data2 = element_data_.size();
516 
517  // if > 0, BDDCML will check if the mesh is not composed of several components and if
518  // yes, BDDCML will handle them as independent subdomains for selection of coarse
519  // degrees of freedom
520  int find_components_int = 1;
521  int use_dual_mesh_graph_int = 0;
522  int neighbouring = 0; // How many common nodes between elements forms and edge of dual graph, only available for use_dual_mesh_graph_int == 1.
523 
524  // upload local data to the solver
525  bddcml_upload_subdomain_data_c( &numElem_, &numNodes_, &numDofsInt, &nDim_, &meshDim_,
526  &iSub, &numElemSub_, &numNodesSub_, &numDofsSubInt,
527  &(inet_[0]),&linet, &(nnet_[0]),&lnnet, &(nndf_[0]),&lnndf,
528  &(isngn_[0]),&lisngn, &(isvgvn_[0]),&lisvgvn, &(isegn_[0]),&lisegn,
529  &(xyz_[0]),&lxyz1,&lxyz2,
530  &(ifix_[0]),&lifix, &(fixv_[0]),&lfixv,
531  &(rhsVec_[0]),&lrhsVec, &isRhsCompleteInt,
532  &(sol_[0]), &lsol,
533  &matrixTypeInt, &(i_sparse[0]), &(j_sparse[0]), &(a_sparse[0]), &la, &isMatAssembledInt,
534  &(userConstraints[0]),&lUserConstraints1,&lUserConstraints2,
535  &(element_data_[0]),&lelement_data1,&lelement_data2,
536  &(sub_diagonal[0]), &lsub_diagonal,
537  &find_components_int, &use_dual_mesh_graph_int, &neighbouring );
538  i_sparse.clear();
539  j_sparse.clear();
540  a_sparse.clear();
541 
542 
543  //============= setup of BDDC preconditioner
544 
545  int use_defaults_int = 0;
546  int parallel_division_int = 0;
547  int use_arithmetic_int = 1;
548  int use_adaptive_int;
549  if ( use_adaptive ) {
550  use_adaptive_int = 1;
551  }
552  else {
553  use_adaptive_int = 0;
554  }
555  int use_user_constraints_int = 0;
556  int weights_type = 3;
557 
558  // if > 0, BDDCML will use constraints on connectivity at corners
559  int use_corner_constraints_int = 1;
560 
561  // setup multilevel BDDC preconditioner
562  START_TIMER("BDDC preconditioner setup");
563  bddcml_setup_preconditioner_c( & matrixTypeInt, & use_defaults_int,
564  & parallel_division_int, & use_arithmetic_int,
565  & use_corner_constraints_int,
566  & use_adaptive_int, & use_user_constraints_int,
567  & weights_type );
568  END_TIMER("BDDC preconditioner setup");
569 
570  //============= compute the norm on arrays with overlapping entries - apply weights
571  double normSquaredLoc = 0.;
572  if (nProc_ > 1 ) {
573  bddcml_dotprod_subdomain_c( &iSub, &(rhsVec_[0]), &lrhsVec, &(rhsVec_[0]), &lrhsVec, &normSquaredLoc );
574  } else {
575  for ( double elem: rhsVec_ ) normSquaredLoc += elem*elem;
576  }
577 
578  // sum over processes
579  double normSquared = 0.;
580  MPI_Allreduce ( &normSquaredLoc, &normSquared, 1, MPI_DOUBLE, MPI_SUM, comm_ );
581  normRhs_ = std::sqrt( normSquared );
582 
583  //============= solution of the problem
584  int method;
585  if ( matrixType_ == SPD )
586  method = 0; // PCG - set for SPD matrix
587  else if ( matrixType_ == SPD_VIA_SYMMETRICGENERAL)
588  method = 0; // PCG - set for benign matrix
589  else
590  method = 1; // BICGSTAB - set for any other matrix
591 
592  int recycling_int = 0;
593  int max_number_of_stored_vectors = 100;
594 
595  // call iterative solver
596  START_TIMER("BDDC solve");
597  bddcml_solve_c( & commInt, & method, & tol, & maxIt, & ndecrMax,
598  & recycling_int, & max_number_of_stored_vectors,
599  & numIter_, & convergedReason_, & condNumber_);
600  END_TIMER("BDDC solve");
601 
602 
603  //============= downloading solution from BDDCML solver
604 
605  // download local solution
606  lsol = sol_.size();
607 
608  bddcml_download_local_solution_c( &iSub, &(sol_[0]), &lsol );
609 
610  //============= compute the norm on arrays with overlapping entries - apply weights
611  normSquaredLoc = 0.;
612  if (nProc_ > 1 ) {
613  bddcml_dotprod_subdomain_c( &iSub, &(sol_[0]), &lsol, &(sol_[0]), &lsol, &normSquaredLoc );
614  } else {
615  for ( double elem: sol_ ) normSquaredLoc += elem*elem;
616  }
617 
618  // sum over processes
619  normSquared = 0.;
620  MPI_Allreduce ( &normSquaredLoc, &normSquared, 1, MPI_DOUBLE, MPI_SUM, comm_ );
621  normSol_ = std::sqrt( normSquared );
622 
623  //============= erase memory used by BDDCML solver
624  bddcml_finalize_c();
625 
626  return;
627 }
628 
629 #endif
la::BddcmlWrapper::insertToDiagonalWeights
void insertToDiagonalWeights(const int &index, const double &value)
la::BddcmlWrapper::sol_
std::vector< double > sol_
distributed solution vector
Definition: bddcml_wrapper.hh:288
la::BddcmlWrapper::~BddcmlWrapper
~BddcmlWrapper()
Destructor frees the BDDCML structures.
la::BddcmlWrapper::fixv_
std::vector< double > fixv_
values of fixed boundary conditions at places of ifix_ - values outside fixed dofs are ignored
Definition: bddcml_wrapper.hh:278
la::BddcmlWrapper::insertToRhs
void insertToRhs(const SubVec_ &eForce, const VecUint_ &dofIndices)
insert subvector to system rhs vector
la::BddcmlWrapper::rhsVec_
std::vector< double > rhsVec_
vector with RHS values restricted to subdomain, i.e. values are repeated at shared nodes
Definition: bddcml_wrapper.hh:285
value
static constexpr bool value
Definition: json.hpp:87
MPI_DOUBLE
#define MPI_DOUBLE
Definition: mpi.h:156
std::vector< int >
la::BddcmlWrapper::isvgvn_
std::vector< int > isvgvn_
Indices of Subdomain Variables in Global Variable Numbering ( i.e. dof mapping )
Definition: bddcml_wrapper.hh:270
la::BddcmlWrapper::ifix_
std::vector< int > ifix_
indices of fixed boundary conditions - ( 0 for free variable, 1 for constrained dof )
Definition: bddcml_wrapper.hh:277
la::BddcmlWrapper::applyConstraints
void applyConstraints(ConstraintVec_ &constraints, const double factor, const double scalar)
MPI_Comm_rank
#define MPI_Comm_rank
Definition: mpi.h:236
MPI_Comm_size
#define MPI_Comm_size
Definition: mpi.h:235
la::BddcmlWrapper::isngn_
std::vector< int > isngn_
Indices of Subdomain Nodes in Global Numbering.
Definition: bddcml_wrapper.hh:269
la::BddcmlWrapper::isegn_
std::vector< int > isegn_
Indices of Subdomain Elements in Global Numbering.
Definition: bddcml_wrapper.hh:271
la::BddcmlWrapper::fixDOF
void fixDOF(const unsigned index, const double scalar=1.)
la::BddcmlWrapper::prepareMatAssembly
void prepareMatAssembly(unsigned numElements, unsigned elMatSize)
Prepare assembly of matrix - reserve space in underlying coordinate matrix object.
la::BddcmlWrapper::nnet_
std::vector< int > nnet_
Definition: bddcml_wrapper.hh:260
MPI_INT
#define MPI_INT
Definition: mpi.h:160
la::BddcmlWrapper::finishDiagAssembly
void finishDiagAssembly()
Finalize assembly of diagonal.
la::BddcmlWrapper::xyz_
std::vector< double > xyz_
Definition: bddcml_wrapper.hh:265
la
Definition: bddcml_wrapper.hh:41
la::MatrixCoo::clear
void clear()
destroy all contents
Definition: matrix_coo.hpp:325
la::BddcmlWrapper::solveSystem
void solveSystem(double tol=1.e-7, int numLevels=2, std::vector< int > *numSubAtLevels=NULL, int verboseLevel=0, int maxIt=1000, int ndecrMax=30, bool use_adaptive=false)
Solve the system.
la::BddcmlWrapper::nndf_
std::vector< int > nndf_
Definition: bddcml_wrapper.hh:263
MPI_SUM
#define MPI_SUM
Definition: mpi.h:196
MPI_Comm
int MPI_Comm
Definition: mpi.h:141
la::BddcmlWrapper::finishMatAssembly
void finishMatAssembly()
Finalize assembly of matrix.
bddcml_wrapper.hh
OLD_ASSERT
#define OLD_ASSERT(...)
Definition: global_defs.h:108
MPI_Allreduce
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
la::BddcmlWrapper::BddcmlWrapper
BddcmlWrapper(const unsigned numDofs, const unsigned numDofsSub, const MatrixType matrixType=GENERAL, const MPI_Comm comm=MPI_COMM_WORLD, int numSubLoc=1)
MPI_Allgather
#define MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm)
Definition: mpi.h:580
MPI_Comm_c2f
#define MPI_Comm_c2f(comm)
Definition: mpi.h:261
la::BddcmlWrapper::loadRawMesh
void loadRawMesh(const int nDim, const int numNodes, 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_data, const int meshDim=0)
Load raw data about mesh.
la::BddcmlWrapper::coo_
la::MatrixCoo< int, double > coo_
matrix in coordinate format (COO)
Definition: bddcml_wrapper.hh:280
la::BddcmlWrapper::insertToMatrix
void insertToMatrix(const SubMat_ &eStiff, const VecUint_ &rowIndices, const VecUint_ &colIndices)
insert submatrix to system matrix - into coordinate matrix object
DebugOut
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
START_TIMER
#define START_TIMER(tag)
Starts a timer with specified tag.
Definition: sys_profiler.hh:115
END_TIMER
#define END_TIMER(tag)
Ends a timer with specified tag.
Definition: sys_profiler.hh:149
la::BddcmlWrapper::inet_
std::vector< int > inet_
safety state variable to warn if solver is called without loading subdomain mesh
Definition: bddcml_wrapper.hh:259
la::BddcmlWrapper::global2LocalDofMap_
Global2LocalMap_ global2LocalDofMap_
type for storage of global to local map
Definition: bddcml_wrapper.hh:274
FMT_UNUSED
#define FMT_UNUSED
Definition: posix.h:75