Flow123d  jenkins-Flow123d-windows32-release-multijob-51
bddcml_wrapper.hpp
Go to the documentation of this file.
1 //! @author Jakub Sistek
2 //! @date 5/2011
3 
4 #ifndef la_bddcml_wrapper_h
5 #define la_bddcml_wrapper_h
6 //------------------------------------------------------------------------------
7 #include <boost/numeric/ublas/matrix.hpp>
8 #include <boost/numeric/ublas/vector.hpp>
9 #include <ostream>
10 #include <vector>
11 #include <set>
12 #include <cmath>
13 #include <la/matrix_coo.hpp>
14 
15 #include <system/global_defs.h>
16 #include <system/xio.h>
17 #include "system/sys_profiler.hh"
18 
19 extern "C" {
20  #include <bddcml_interface_c.h>
21 }
22 
23 //------------------------------------------------------------------------------
24 namespace la{
25  class BddcmlWrapper;
26  namespace ublas = boost::numeric::ublas;
27 
28  namespace detail_{
29 
30  template< typename DOFHOLDER, typename MAPTYPE >
31  int getLocalNumber( DOFHOLDER * np, MAPTYPE & map ) {
32 
33  unsigned globalIndex = np -> giveId( );
34 
35  typename MAPTYPE::iterator pos = map.find( globalIndex );
36  ASSERT( pos != map.end(),
37  "Cannot remap node index %d to local indices. \n ", globalIndex );
38  int localIndexInt = pos -> second;
39 
40  return localIndexInt;
41  }
42  }
43 
44 }
45 
46 //------------------------------------------------------------------------------
47 /** \brief Multilevel BDDC based linear system solver
48  * \details This class implements solution of a system of linear equations
49  * by multilevel BDDC implementation in BDDCML solver.
50  *
51  * The class uses assembly of the matrix into
52  * the coordinate matrix format implemented in matrix_coo.hpp file.
53  *
54  * prepareMatAssembly - (optional) - reserves memory for sparse coordinate matrix
55  * insertToMatrix - called many times, inserts values into coordinate matrix
56  * insertToRhs - called many times, inserts values into RHS.
57  * If matrix is created and assembled, structure of the vector is deduced from it.
58  * Otherwise, it is created as empty.
59  * applyConstraints - currently only passes list of boundary conditions into the BDDCML solver,
60  * rather then modifying the matrix and rhs
61  * - as a result, no nonsymmetry is forced to SPD and symmetric indefinite problems by constraints,
62  * and these can be accelerated and be more memory efficient by specifying proper matrixType
63  * solveSystem - solves the system
64  */
66 {
67 public:
69  GENERAL, //!< general (default),
70  SPD, //!< symmetric positive definite,
71  SYMMETRICGENERAL, //!< general symmetric ( e.g. saddle point ),
72  SPD_VIA_SYMMETRICGENERAL //!< (advanced) interface problem is symm. pos. def.,
73  //!< but subdomain problems are general symmetric ( e.g. saddle point )
74  };
75  typedef enum matrixTypeEnum MatrixType;
76 
77 private:
78  typedef ublas::matrix< double > SubMat_;
79  typedef ublas::vector< double > SubVec_;
82 
83 public:
84  //! Constructor with
85  //! total number of dofs,
86  //! subdomain number of dofs,
87  //! (optional) type of matrix,
88  //! (optional) communicator ( MPI_COMM_WORLD is the default ),
89  //! (optional) number of subdomains ( if not given, use number of processes ),
90  BddcmlWrapper( const unsigned numDofs,
91  const unsigned numDofsSub,
92  const MatrixType matrixType = GENERAL,
93  const MPI_Comm comm = MPI_COMM_WORLD,
94  int numSubLoc = 1 );
95 
96  //! Destructor frees the BDDCML structures
98 
99  //! Load raw data about mesh
100  void loadRawMesh( const int nDim, const int numNodes, const int numDofs,
101  const std::vector<int> & inet,
102  const std::vector<int> & nnet,
103  const std::vector<int> & nndf,
104  const std::vector<int> & isegn,
105  const std::vector<int> & isngn,
106  const std::vector<int> & isvgvn,
107  const std::vector<double> & xyz,
108  const std::vector<double> & element_data,
109  const int meshDim = 0 );
110 
111  //! Prepare assembly of matrix - reserve space in underlying coordinate matrix object
112  void prepareMatAssembly( unsigned numElements, unsigned elMatSize );
113 
114  //! insert submatrix to system matrix - into coordinate matrix object
115  void insertToMatrix( const SubMat_ & eStiff,
116  const VecUint_ & rowIndices,
117  const VecUint_ & colIndices );
118 
119  //! Finalize assembly of matrix
120  void finishMatAssembly( );
121 
122  //! insert subvector to system rhs vector
123  void insertToRhs( const SubVec_ & eForce, const VecUint_ & dofIndices );
124 
125  //! Outputs
126  void writeMatrix( std::ostream & out ) { coo_.write( out ); }
127 
128  //! Compatibility call which does absolutely nothing
129  void finishAssembly( ) { }
130 
131  //! Apply constraints
132  //! Currently, scalar is not used in the function and it is kept only for compatibility with constraint functors.
133  //! Proper scalar is enforced within the BDDC solver.
134  void applyConstraints( ConstraintVec_ & constraints, const double factor, const double scalar );
135 
136  //! Set any DOF to zero
137  //!
138  //! \param[in] index Global index of DOF to fix
139  //! \param[in] scalar Value to put on diagonal (default 1.)
140  void fixDOF( const unsigned index, const double scalar = 1. );
141 
142  //! Insert a value to diagonal for creating weights in BDDC
143  //!
144  //! \param[in] index Global index of DOF
145  //! \param[in] value Value to put on diagonal
146  void insertToDiagonalWeights( const int & index, const double & value );
147 
148  //! Finalize assembly of diagonal
149  void finishDiagAssembly( );
150 
151  //! Solve the system
152  void solveSystem( double tol = 1.e-7, //!< tolerance on relative residual ||res||/||rhs||
153  int numLevels = 2, //!< number of levels
154  std::vector<int> * numSubAtLevels = NULL, //!< number of subdomains at levels
155  int verboseLevel = 0, //!< level of verbosity of BDDCML library
156  //!< ( 0 - only fatal errors reported,
157  //!< 1 - mild output,
158  //!< 2 - detailed output )
159  int maxIt = 1000, //!< maximum number of iterations
160  int ndecrMax = 30, //!< maximum number of iterations with non-decreasing residual
161  //!< ( used to stop diverging process )
162  bool use_adaptive = false //!< if adaptive constraints should be used by BDDCML
163  );
164 
165  //! Get norm of right-hand side
166  double normRhs( ) { return normRhs_; }
167 
168  //! Get norm of solution
169  double normSol( ) { return normSol_; }
170 
171  //! Give number of iteration
172  int giveNumIterations() { return numIter_ ; }
173 
174  //! Give reason of convergence
176 
177  //! Give estimate of condition number
178  double giveCondNumber() { return condNumber_ ; }
179 
180  //! Fill matrix with zeros
181  void clearMatrix( ) {
182  // remember length
183  unsigned length = coo_.nnz( );
184  coo_.clear( );
185  coo_.prepareAssembly( length );
186  }
187  //! Fill RHS with zeros
188  void clearRhs( ) {
189  std::fill( rhsVec_.begin( ), rhsVec_.end( ), 0. );
190  normRhs_ = 0.;
191  }
192 
193  //! Blank arrays for boundary conditions
194  void clearBC( ) {
195  std::fill( ifix_.begin(), ifix_.end(), 0 );
196  std::fill( fixv_.begin(), fixv_.end(), 0. );
197  }
198 
199  //! Fill solution with zeros
200  void clearSol( )
201  { std::fill( sol_.begin( ), sol_.end( ), 0. );
202  normSol_ = 0.;
203  }
204 
205  //! Destroy matrix
206  void destroyMatrix( ) {
207  coo_.clear( );
208  }
209  //! Destroy RHS
210  void destroyRhs( ) { rhsVec_.clear( ); }
211  //! Destroy boundary conditions
212  void destroyBC( ) {
213  ifix_.clear( );
214  fixv_.clear( );
215  }
216  //! Destroy solution
217  void destroySol( ) { sol_.clear( ); }
218 
219  //! Give access to the solution vector
220  template<typename VEC1,typename VEC2>
221  void giveSolution( const VEC1 & dofIndices, VEC2 & result ) const;
222 
223 private:
224  const int numDofs_; //!< number of total dofs
225  const int numDofsSub_; //!< number of subdomain dofs
226  const MatrixType matrixType_; //!< type of matrix
227  const MPI_Comm comm_; //!< communicator
228  int numSubLoc_; //!< number of subdomains on process
229 
230  int nProc_; //!< number of processors in communicator
231  int rank_; //!< ID of processors in communicator
232 
233  int nDim_; //!< space dimension
234  int meshDim_; //!< mesh dimension, may differ from space dimension - e.g. 2 for shells in 3D
235 
236  int numElemSub_; //!< number of elements in subdomain
237  int numElem_; //!< total number of elements
238  int numNodesSub_; //!< number of nodes in subdomain
239  int numNodes_; //!< total number of nodes
240 
241  bool meshLoaded_; //! safety state variable to warn if solver is called without loading subdomain mesh
242  std::vector<int> inet_; //!< indices of nodes on elements (linearized connectivity) - numbering w.r.t. local subdomain numbering
243  std::vector<int> nnet_; //!< numbers of nodes on elements
244  //!< (determines chunks of nodes for each element in inet),
245  //!< e.g. [ 8, 8, 8, ... ] for linear hexahedra
246  std::vector<int> nndf_; //!< number of nodal degrees of freedom ( entry for each node ),
247  //!< e.g. [ 3, 3, 3, ... ] for 3D solids
248  std::vector<double> xyz_; //!< linearized array of coordinates - all x values, all y values, all z values ( if present )
249  //!< i.e. [ x_0 x_1 x_2 ... x_(numNodes - 1) y_0 y_1 y_2 ... y_(numNodes - 1) z_0 z_1 z_2 ... z_(numNodes - 1) ]
250  std::vector<double> element_data_; //!< array with data on elements, here used for permeability for each element
251 
252  std::vector<int> isngn_; //!< Indices of Subdomain Nodes in Global Numbering
253  std::vector<int> isvgvn_; //!< Indices of Subdomain Variables in Global Variable Numbering ( i.e. dof mapping )
254  std::vector<int> isegn_; //!< Indices of Subdomain Elements in Global Numbering
255 
256  typedef std::map<unsigned,unsigned> Global2LocalMap_; //! type for storage of global to local map
257  Global2LocalMap_ global2LocalDofMap_; //!< map from global dof indices to subdomain local indices
258 
259 
260  std::vector<int> ifix_; //!< indices of fixed boundary conditions - ( 0 for free variable, 1 for constrained dof )
261  std::vector<double> fixv_; //!< values of fixed boundary conditions at places of ifix_ - values outside fixed dofs are ignored
262 
263  la::MatrixCoo<int,double> coo_; //!< matrix in coordinate format (COO)
264  bool isMatAssembled_; //!< true if matrix is assembled
265  la::MatrixCoo<int,double> diagWeightsCoo_; //!< diagonal of subdomain matrix in sparse format and global numbers
266  bool isDiagAssembled_; //!< true if diagonal is assembled
267 
268  std::vector<double> rhsVec_; //!< vector with RHS values restricted to subdomain, i.e. values are repeated at shared nodes
269  double normRhs_; //!< norm of the global right-hand side
270 
271  std::vector<double> sol_; //!< distributed solution vector
272  double normSol_; //!< norm of the global solution
273 
274  int numIter_; //!< required number of iterations
275  int convergedReason_; //!< reason of convergence:
276  //!< 0 - converged relative residual, desirable state
277  //!< -1 - reached maximal number of iterations
278  //!< -2 - reached maximal number of iterations with non-decreasing residual
279  double condNumber_; //!< estimated condition number of the preconditioned operator
280  //!< only provided for SPD problems solved by PCG ( -1.0 if not meaningful )
281 
282 };
283 
284 //-----------------------------------------------------------------------------
285 
286 #include "bddcml_wrapper.ipp"
287 
288 //-----------------------------------------------------------------------------
289 
290 #endif
double giveCondNumber()
Give estimate of condition number.
void applyConstraints(ConstraintVec_ &constraints, const double factor, const double scalar)
void destroyRhs()
Destroy RHS.
double normSol()
Get norm of solution.
void prepareAssembly(const IndexType length)
Definition: matrix_coo.hpp:65
double normRhs()
Get norm of right-hand side.
std::vector< double > fixv_
values of fixed boundary conditions at places of ifix_ - values outside fixed dofs are ignored ...
BddcmlWrapper(const unsigned numDofs, const unsigned numDofsSub, const MatrixType matrixType=GENERAL, const MPI_Comm comm=MPI_COMM_WORLD, int numSubLoc=1)
std::ostream & write(std::ostream &out)
debug-output
Definition: matrix_coo.hpp:332
int MPI_Comm
Definition: mpi.h:141
int giveConvergedReason()
Give reason of convergence.
int nProc_
number of processors in communicator
int numIter_
required number of iterations
void finishMatAssembly()
Finalize assembly of matrix.
int rank_
ID of processors in communicator.
void prepareMatAssembly(unsigned numElements, unsigned elMatSize)
Prepare assembly of matrix - reserve space in underlying coordinate matrix object.
std::vector< int > inet_
safety state variable to warn if solver is called without loading subdomain mesh
std::vector< int > nndf_
~BddcmlWrapper()
Destructor frees the BDDCML structures.
std::vector< double > xyz_
int numNodes_
total number of nodes
void destroySol()
Destroy solution.
Global2LocalMap_ global2LocalDofMap_
type for storage of global to local map
enum matrixTypeEnum MatrixType
IndexType nnz() const
Return number of non-zeros.
Definition: matrix_coo.hpp:60
la::MatrixCoo< int, double > coo_
matrix in coordinate format (COO)
void insertToMatrix(const SubMat_ &eStiff, const VecUint_ &rowIndices, const VecUint_ &colIndices)
insert submatrix to system matrix - into coordinate matrix object
double normSol_
norm of the global solution
int nDim_
space dimension
I/O functions with filename storing, able to track current line in opened file. All standard stdio fu...
const MatrixType matrixType_
type of matrix
int numElem_
total number of elements
int numSubLoc_
number of subdomains on process
std::vector< double > sol_
distributed solution vector
Global macros to enhance readability and debugging, general constants.
void giveSolution(const VEC1 &dofIndices, VEC2 &result) const
Give access to the solution vector.
#define ASSERT(...)
Definition: global_defs.h:121
std::vector< int > isvgvn_
Indices of Subdomain Variables in Global Variable Numbering ( i.e. dof mapping )
std::vector< int > isegn_
Indices of Subdomain Elements in Global Numbering.
int numElemSub_
number of elements in subdomain
void clearSol()
Fill solution with zeros.
int getLocalNumber(DOFHOLDER *np, MAPTYPE &map)
std::vector< double > rhsVec_
vector with RHS values restricted to subdomain, i.e. values are repeated at shared nodes ...
bool isMatAssembled_
true if matrix is assembled
void insertToDiagonalWeights(const int &index, const double &value)
std::vector< int > ifix_
indices of fixed boundary conditions - ( 0 for free variable, 1 for constrained dof ) ...
symmetric positive definite,
general symmetric ( e.g. saddle point ),
void clearMatrix()
Fill matrix with zeros.
bool isDiagAssembled_
true if diagonal is assembled
void destroyBC()
Destroy boundary conditions.
void fixDOF(const unsigned index, const double scalar=1.)
std::vector< int > isngn_
Indices of Subdomain Nodes in Global Numbering.
std::vector< unsigned > VecUint_
void loadRawMesh(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_data, const int meshDim=0)
Load raw data about mesh.
int meshDim_
mesh dimension, may differ from space dimension - e.g. 2 for shells in 3D
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.
int numNodesSub_
number of nodes in subdomain
std::map< unsigned, unsigned > Global2LocalMap_
std::vector< int > nnet_
void finishDiagAssembly()
Finalize assembly of diagonal.
#define MPI_COMM_WORLD
Definition: mpi.h:123
int giveNumIterations()
Give number of iteration.
const MPI_Comm comm_
communicator
ublas::vector< double > SubVec_
void clearRhs()
Fill RHS with zeros.
std::vector< double > element_data_
array with data on elements, here used for permeability for each element
void finishAssembly()
Compatibility call which does absolutely nothing.
la::MatrixCoo< int, double > diagWeightsCoo_
diagonal of subdomain matrix in sparse format and global numbers
void clearBC()
Blank arrays for boundary conditions.
void writeMatrix(std::ostream &out)
Outputs.
const int numDofsSub_
number of subdomain dofs
void destroyMatrix()
Destroy matrix.
Multilevel BDDC based linear system solver.
const int numDofs_
number of total dofs
ublas::matrix< double > SubMat_
void clear()
destroy all contents
Definition: matrix_coo.hpp:308
double normRhs_
norm of the global right-hand side
std::vector< std::pair< unsigned, double > > ConstraintVec_
void insertToRhs(const SubVec_ &eForce, const VecUint_ &dofIndices)
insert subvector to system rhs vector