Flow123d
JS_before_hm-1626-gde32303
|
Multilevel BDDC based linear system solver. More...
#include <bddcml_wrapper.hh>
Public Types | |
enum | matrixTypeEnum { GENERAL, SPD, SYMMETRICGENERAL, SPD_VIA_SYMMETRICGENERAL } |
typedef enum matrixTypeEnum | MatrixType |
Public Member Functions | |
BddcmlWrapper (const unsigned numDofs, const unsigned numDofsSub, const MatrixType matrixType=GENERAL, const MPI_Comm comm=MPI_COMM_WORLD, int numSubLoc=1) | |
~BddcmlWrapper () | |
Destructor frees the BDDCML structures. More... | |
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. More... | |
void | prepareMatAssembly (unsigned numElements, unsigned elMatSize) |
Prepare assembly of matrix - reserve space in underlying coordinate matrix object. More... | |
void | insertToMatrix (const SubMat_ &eStiff, const VecUint_ &rowIndices, const VecUint_ &colIndices) |
insert submatrix to system matrix - into coordinate matrix object More... | |
void | finishMatAssembly () |
Finalize assembly of matrix. More... | |
void | insertToRhs (const SubVec_ &eForce, const VecUint_ &dofIndices) |
insert subvector to system rhs vector More... | |
void | writeMatrix (std::ostream &out) |
Outputs. More... | |
void | finishAssembly () |
Compatibility call which does absolutely nothing. More... | |
void | applyConstraints (ConstraintVec_ &constraints, const double factor, const double scalar) |
void | fixDOF (const unsigned index, const double scalar=1.) |
void | insertToDiagonalWeights (const int &index, const double &value) |
void | finishDiagAssembly () |
Finalize assembly of diagonal. More... | |
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. More... | |
double | normRhs () |
Get norm of right-hand side. More... | |
double | normSol () |
Get norm of solution. More... | |
int | giveNumIterations () |
Give number of iteration. More... | |
int | giveConvergedReason () |
Give reason of convergence. More... | |
double | giveCondNumber () |
Give estimate of condition number. More... | |
void | clearMatrix () |
Fill matrix with zeros. More... | |
void | clearRhs () |
Fill RHS with zeros. More... | |
void | clearBC () |
Blank arrays for boundary conditions. More... | |
void | clearSol () |
Fill solution with zeros. More... | |
void | destroyMatrix () |
Destroy matrix. More... | |
void | destroyRhs () |
Destroy RHS. More... | |
void | destroyBC () |
Destroy boundary conditions. More... | |
void | destroySol () |
Destroy solution. More... | |
template<typename VEC1 , typename VEC2 > | |
void | giveSolution (const VEC1 &dofIndices, VEC2 &result) const |
Give access to the solution vector. More... | |
Private Types | |
typedef ublas::matrix< double > | SubMat_ |
typedef ublas::vector< double > | SubVec_ |
typedef std::vector< unsigned > | VecUint_ |
typedef std::vector< std::pair< unsigned, double > > | ConstraintVec_ |
typedef std::map< unsigned, unsigned > | Global2LocalMap_ |
Private Attributes | |
const int | numDofs_ |
number of total dofs More... | |
const int | numDofsSub_ |
number of subdomain dofs More... | |
const MatrixType | matrixType_ |
type of matrix More... | |
const MPI_Comm | comm_ |
communicator More... | |
int | numSubLoc_ |
number of subdomains on process More... | |
int | nProc_ |
number of processors in communicator More... | |
int | rank_ |
ID of processors in communicator. More... | |
int | nDim_ |
space dimension More... | |
int | meshDim_ |
mesh dimension, may differ from space dimension - e.g. 2 for shells in 3D More... | |
int | numElemSub_ |
number of elements in subdomain More... | |
int | numElem_ |
total number of elements More... | |
int | numNodesSub_ |
number of nodes in subdomain More... | |
int | numNodes_ |
total number of nodes More... | |
bool | meshLoaded_ |
std::vector< int > | inet_ |
safety state variable to warn if solver is called without loading subdomain mesh More... | |
std::vector< int > | nnet_ |
std::vector< int > | nndf_ |
std::vector< double > | xyz_ |
std::vector< double > | element_data_ |
array with data on elements, here used for permeability for each element More... | |
std::vector< int > | isngn_ |
Indices of Subdomain Nodes in Global Numbering. More... | |
std::vector< int > | isvgvn_ |
Indices of Subdomain Variables in Global Variable Numbering ( i.e. dof mapping ) More... | |
std::vector< int > | isegn_ |
Indices of Subdomain Elements in Global Numbering. More... | |
Global2LocalMap_ | global2LocalDofMap_ |
type for storage of global to local map More... | |
std::vector< int > | ifix_ |
indices of fixed boundary conditions - ( 0 for free variable, 1 for constrained dof ) More... | |
std::vector< double > | fixv_ |
values of fixed boundary conditions at places of ifix_ - values outside fixed dofs are ignored More... | |
la::MatrixCoo< int, double > | coo_ |
matrix in coordinate format (COO) More... | |
bool | isMatAssembled_ |
true if matrix is assembled More... | |
la::MatrixCoo< int, double > | diagWeightsCoo_ |
diagonal of subdomain matrix in sparse format and global numbers More... | |
bool | isDiagAssembled_ |
true if diagonal is assembled More... | |
std::vector< double > | rhsVec_ |
vector with RHS values restricted to subdomain, i.e. values are repeated at shared nodes More... | |
double | normRhs_ |
norm of the global right-hand side More... | |
std::vector< double > | sol_ |
distributed solution vector More... | |
double | normSol_ |
norm of the global solution More... | |
int | numIter_ |
required number of iterations More... | |
int | convergedReason_ |
double | condNumber_ |
Multilevel BDDC based linear system solver.
This class implements solution of a system of linear equations by multilevel BDDC implementation in BDDCML solver.
The class uses assembly of the matrix into the coordinate matrix format implemented in matrix_coo.hpp file.
prepareMatAssembly - (optional) - reserves memory for sparse coordinate matrix insertToMatrix - called many times, inserts values into coordinate matrix insertToRhs - called many times, inserts values into RHS. If matrix is created and assembled, structure of the vector is deduced from it. Otherwise, it is created as empty. applyConstraints - currently only passes list of boundary conditions into the BDDCML solver, rather then modifying the matrix and rhs
Definition at line 80 of file bddcml_wrapper.hh.
|
private |
Definition at line 96 of file bddcml_wrapper.hh.
|
private |
Definition at line 271 of file bddcml_wrapper.hh.
typedef enum matrixTypeEnum la::BddcmlWrapper::MatrixType |
Definition at line 90 of file bddcml_wrapper.hh.
|
private |
Definition at line 93 of file bddcml_wrapper.hh.
|
private |
Definition at line 94 of file bddcml_wrapper.hh.
|
private |
Definition at line 95 of file bddcml_wrapper.hh.
Definition at line 83 of file bddcml_wrapper.hh.
la::BddcmlWrapper::BddcmlWrapper | ( | const unsigned | numDofs, |
const unsigned | numDofsSub, | ||
const MatrixType | matrixType = GENERAL , |
||
const MPI_Comm | comm = MPI_COMM_WORLD , |
||
int | numSubLoc = 1 |
||
) |
Constructor with total number of dofs, subdomain number of dofs, (optional) type of matrix, (optional) communicator ( MPI_COMM_WORLD is the default ), (optional) number of subdomains ( if not given, use number of processes ),
Definition at line 9 of file bddcml_wrapper.cc.
la::BddcmlWrapper::~BddcmlWrapper | ( | ) |
Destructor frees the BDDCML structures.
Destructor frees the BDDC solver structures
Definition at line 46 of file bddcml_wrapper.cc.
void la::BddcmlWrapper::applyConstraints | ( | ConstraintVec_ & | constraints, |
const double | factor, | ||
const double | scalar | ||
) |
Apply constraints Currently, scalar is not used in the function and it is kept only for compatibility with constraint functors. Proper scalar is enforced within the BDDC solver.
Apply boundary conditions Since BDDCML handles constraints, just fill proper BC arrays. Scalar is not used, it is kept for compatibility with constraint functors.
Definition at line 248 of file bddcml_wrapper.cc.
|
inline |
Blank arrays for boundary conditions.
Definition at line 209 of file bddcml_wrapper.hh.
|
inline |
Fill matrix with zeros.
Definition at line 196 of file bddcml_wrapper.hh.
|
inline |
Fill RHS with zeros.
Definition at line 203 of file bddcml_wrapper.hh.
|
inline |
Fill solution with zeros.
Definition at line 215 of file bddcml_wrapper.hh.
|
inline |
Destroy boundary conditions.
Definition at line 227 of file bddcml_wrapper.hh.
|
inline |
Destroy matrix.
Definition at line 221 of file bddcml_wrapper.hh.
|
inline |
Destroy RHS.
Definition at line 225 of file bddcml_wrapper.hh.
|
inline |
Destroy solution.
Definition at line 232 of file bddcml_wrapper.hh.
|
inline |
Compatibility call which does absolutely nothing.
Definition at line 144 of file bddcml_wrapper.hh.
void la::BddcmlWrapper::finishDiagAssembly | ( | ) |
Finalize assembly of diagonal.
Routine for diagonal finalization.
Definition at line 302 of file bddcml_wrapper.cc.
void la::BddcmlWrapper::finishMatAssembly | ( | ) |
Finalize assembly of matrix.
Routine for triplet finalization.
Definition at line 205 of file bddcml_wrapper.cc.
void la::BddcmlWrapper::fixDOF | ( | const unsigned | index, |
const double | scalar = 1. |
||
) |
Set any DOF to zero
[in] | index | Global index of DOF to fix |
[in] | scalar | Value to put on diagonal (default 1.) |
Simply set the dof with index 'index' to zero
Definition at line 275 of file bddcml_wrapper.cc.
|
inline |
Give estimate of condition number.
Definition at line 193 of file bddcml_wrapper.hh.
|
inline |
Give reason of convergence.
Definition at line 190 of file bddcml_wrapper.hh.
|
inline |
Give number of iteration.
Definition at line 187 of file bddcml_wrapper.hh.
void la::BddcmlWrapper::giveSolution | ( | const VEC1 & | dofIndices, |
VEC2 & | result | ||
) | const |
Give access to the solution vector.
Definition at line 303 of file bddcml_wrapper.hh.
void la::BddcmlWrapper::insertToDiagonalWeights | ( | const int & | index, |
const double & | value | ||
) |
Insert a value to diagonal for creating weights in BDDC
[in] | index | Global index of DOF |
[in] | value | Value to put on diagonal |
Insert elements on diagonal for weights.
Definition at line 287 of file bddcml_wrapper.cc.
void la::BddcmlWrapper::insertToMatrix | ( | const SubMat_ & | subMat, |
const VecUint_ & | rowIndices, | ||
const VecUint_ & | colIndices | ||
) |
insert submatrix to system matrix - into coordinate matrix object
Insert element stiffness matrix into global system matrix
Definition at line 182 of file bddcml_wrapper.cc.
insert subvector to system rhs vector
Insert subvector to the system's RHS vector
Definition at line 221 of file bddcml_wrapper.cc.
void la::BddcmlWrapper::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.
Routine for loading raw mesh data into the solver - for cases of strange meshes, where these are not Mesh or Grid objects, user can create own raw description to exploit the flexibility of mesh format underlaying BDDCML.
Definition at line 80 of file bddcml_wrapper.cc.
|
inline |
Get norm of right-hand side.
Definition at line 181 of file bddcml_wrapper.hh.
|
inline |
Get norm of solution.
Definition at line 184 of file bddcml_wrapper.hh.
void la::BddcmlWrapper::prepareMatAssembly | ( | unsigned | numElements, |
unsigned | elMatSize | ||
) |
Prepare assembly of matrix - reserve space in underlying coordinate matrix object.
Routine for preparing MATRIX assembly. It reserves memory in triplet.
Definition at line 172 of file bddcml_wrapper.cc.
void la::BddcmlWrapper::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.
Solve the system by BDDCML
< number of subdomains on first level
< numbers of subdomains on levels ( should be monotonically decreasing, e.g. for 3 levels 1024, 32, 1 )
tol | tolerance on relative residual ||res||/||rhs|| |
numLevels | number of levels |
numSubAtLevels | number of subdomains at levels |
verboseLevel | level of verbosity of BDDCML library ( 0 - only fatal errors reported, 1 - mild output, 2 - detailed output ) |
maxIt | maximum number of iterations |
ndecrMax | maximum number of iterations with non-decreasing residual ( used to stop diverging process ) |
use_adaptive | if adaptive constraints should be used by BDDCML |
Definition at line 319 of file bddcml_wrapper.cc.
|
inline |
Outputs.
Definition at line 141 of file bddcml_wrapper.hh.
|
private |
communicator
Definition at line 242 of file bddcml_wrapper.hh.
|
private |
estimated condition number of the preconditioned operator only provided for SPD problems solved by PCG ( -1.0 if not meaningful )
Definition at line 294 of file bddcml_wrapper.hh.
|
private |
reason of convergence: 0 - converged relative residual, desirable state -1 - reached maximal number of iterations -2 - reached maximal number of iterations with non-decreasing residual
Definition at line 290 of file bddcml_wrapper.hh.
|
private |
matrix in coordinate format (COO)
Definition at line 278 of file bddcml_wrapper.hh.
|
private |
diagonal of subdomain matrix in sparse format and global numbers
Definition at line 280 of file bddcml_wrapper.hh.
|
private |
array with data on elements, here used for permeability for each element
Definition at line 265 of file bddcml_wrapper.hh.
|
private |
values of fixed boundary conditions at places of ifix_ - values outside fixed dofs are ignored
Definition at line 276 of file bddcml_wrapper.hh.
|
private |
type for storage of global to local map
map from global dof indices to subdomain local indices
Definition at line 272 of file bddcml_wrapper.hh.
|
private |
indices of fixed boundary conditions - ( 0 for free variable, 1 for constrained dof )
Definition at line 275 of file bddcml_wrapper.hh.
|
private |
safety state variable to warn if solver is called without loading subdomain mesh
indices of nodes on elements (linearized connectivity) - numbering w.r.t. local subdomain numbering
Definition at line 257 of file bddcml_wrapper.hh.
|
private |
true if diagonal is assembled
Definition at line 281 of file bddcml_wrapper.hh.
|
private |
Indices of Subdomain Elements in Global Numbering.
Definition at line 269 of file bddcml_wrapper.hh.
|
private |
true if matrix is assembled
Definition at line 279 of file bddcml_wrapper.hh.
|
private |
Indices of Subdomain Nodes in Global Numbering.
Definition at line 267 of file bddcml_wrapper.hh.
|
private |
Indices of Subdomain Variables in Global Variable Numbering ( i.e. dof mapping )
Definition at line 268 of file bddcml_wrapper.hh.
|
private |
type of matrix
Definition at line 241 of file bddcml_wrapper.hh.
|
private |
mesh dimension, may differ from space dimension - e.g. 2 for shells in 3D
Definition at line 249 of file bddcml_wrapper.hh.
|
private |
Definition at line 256 of file bddcml_wrapper.hh.
|
private |
space dimension
Definition at line 248 of file bddcml_wrapper.hh.
|
private |
number of nodal degrees of freedom ( entry for each node ), e.g. [ 3, 3, 3, ... ] for 3D solids
Definition at line 261 of file bddcml_wrapper.hh.
|
private |
numbers of nodes on elements (determines chunks of nodes for each element in inet), e.g. [ 8, 8, 8, ... ] for linear hexahedra
Definition at line 258 of file bddcml_wrapper.hh.
|
private |
norm of the global right-hand side
Definition at line 284 of file bddcml_wrapper.hh.
|
private |
norm of the global solution
Definition at line 287 of file bddcml_wrapper.hh.
|
private |
number of processors in communicator
Definition at line 245 of file bddcml_wrapper.hh.
|
private |
number of total dofs
Definition at line 239 of file bddcml_wrapper.hh.
|
private |
number of subdomain dofs
Definition at line 240 of file bddcml_wrapper.hh.
|
private |
total number of elements
Definition at line 252 of file bddcml_wrapper.hh.
|
private |
number of elements in subdomain
Definition at line 251 of file bddcml_wrapper.hh.
|
private |
required number of iterations
Definition at line 289 of file bddcml_wrapper.hh.
|
private |
total number of nodes
Definition at line 254 of file bddcml_wrapper.hh.
|
private |
number of nodes in subdomain
Definition at line 253 of file bddcml_wrapper.hh.
|
private |
number of subdomains on process
Definition at line 243 of file bddcml_wrapper.hh.
|
private |
ID of processors in communicator.
Definition at line 246 of file bddcml_wrapper.hh.
|
private |
vector with RHS values restricted to subdomain, i.e. values are repeated at shared nodes
Definition at line 283 of file bddcml_wrapper.hh.
|
private |
distributed solution vector
Definition at line 286 of file bddcml_wrapper.hh.
|
private |
linearized array of coordinates - all x values, all y values, all z values ( if present ) 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) ]
Definition at line 263 of file bddcml_wrapper.hh.