Flow123d  release_2.2.0-914-gf1a3a4f
Public Types | Public Member Functions | Private Types | Private Attributes | List of all members
la::BddcmlWrapper Class Reference

Multilevel BDDC based linear system solver. More...

#include <bddcml_wrapper.hpp>

Collaboration diagram for la::BddcmlWrapper:
Collaboration graph
[legend]

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 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. 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_
 

Detailed Description

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.hpp.

Member Typedef Documentation

typedef std::vector< std::pair<unsigned,double> > la::BddcmlWrapper::ConstraintVec_
private

Definition at line 96 of file bddcml_wrapper.hpp.

typedef std::map<unsigned,unsigned> la::BddcmlWrapper::Global2LocalMap_
private

Definition at line 271 of file bddcml_wrapper.hpp.

Definition at line 90 of file bddcml_wrapper.hpp.

typedef ublas::matrix< double > la::BddcmlWrapper::SubMat_
private

Definition at line 93 of file bddcml_wrapper.hpp.

typedef ublas::vector< double > la::BddcmlWrapper::SubVec_
private

Definition at line 94 of file bddcml_wrapper.hpp.

typedef std::vector< unsigned > la::BddcmlWrapper::VecUint_
private

Definition at line 95 of file bddcml_wrapper.hpp.

Member Enumeration Documentation

Enumerator
GENERAL 

general (default),

SPD 

symmetric positive definite,

SYMMETRICGENERAL 

general symmetric ( e.g. saddle point ),

SPD_VIA_SYMMETRICGENERAL 

(advanced) interface problem is symm. pos. def., but subdomain problems are general symmetric ( e.g. saddle point )

Definition at line 83 of file bddcml_wrapper.hpp.

Constructor & Destructor Documentation

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 ),

la::BddcmlWrapper::~BddcmlWrapper ( )

Destructor frees the BDDCML structures.

Member Function Documentation

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.

void la::BddcmlWrapper::clearBC ( )
inline

Blank arrays for boundary conditions.

Definition at line 209 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::clearMatrix ( )
inline

Fill matrix with zeros.

Definition at line 196 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::clearRhs ( )
inline

Fill RHS with zeros.

Definition at line 203 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::clearSol ( )
inline

Fill solution with zeros.

Definition at line 215 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::destroyBC ( )
inline

Destroy boundary conditions.

Definition at line 227 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::destroyMatrix ( )
inline

Destroy matrix.

Definition at line 221 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::destroyRhs ( )
inline

Destroy RHS.

Definition at line 225 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::destroySol ( )
inline

Destroy solution.

Definition at line 232 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::finishAssembly ( )
inline

Compatibility call which does absolutely nothing.

Definition at line 144 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::finishDiagAssembly ( )

Finalize assembly of diagonal.

void la::BddcmlWrapper::finishMatAssembly ( )

Finalize assembly of matrix.

void la::BddcmlWrapper::fixDOF ( const unsigned  index,
const double  scalar = 1. 
)

Set any DOF to zero

Parameters
[in]indexGlobal index of DOF to fix
[in]scalarValue to put on diagonal (default 1.)
double la::BddcmlWrapper::giveCondNumber ( )
inline

Give estimate of condition number.

Definition at line 193 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::giveConvergedReason ( )
inline

Give reason of convergence.

Definition at line 190 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::giveNumIterations ( )
inline

Give number of iteration.

Definition at line 187 of file bddcml_wrapper.hpp.

template<typename VEC1 , typename VEC2 >
void la::BddcmlWrapper::giveSolution ( const VEC1 &  dofIndices,
VEC2 &  result 
) const

Give access to the solution vector.

void la::BddcmlWrapper::insertToDiagonalWeights ( const int &  index,
const double &  value 
)

Insert a value to diagonal for creating weights in BDDC

Parameters
[in]indexGlobal index of DOF
[in]valueValue to put on diagonal
void la::BddcmlWrapper::insertToMatrix ( const SubMat_ eStiff,
const VecUint_ rowIndices,
const VecUint_ colIndices 
)

insert submatrix to system matrix - into coordinate matrix object

void la::BddcmlWrapper::insertToRhs ( const SubVec_ eForce,
const VecUint_ dofIndices 
)

insert subvector to system rhs vector

void la::BddcmlWrapper::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.

double la::BddcmlWrapper::normRhs ( )
inline

Get norm of right-hand side.

Definition at line 181 of file bddcml_wrapper.hpp.

Here is the caller graph for this function:

double la::BddcmlWrapper::normSol ( )
inline

Get norm of solution.

Definition at line 184 of file bddcml_wrapper.hpp.

void la::BddcmlWrapper::prepareMatAssembly ( unsigned  numElements,
unsigned  elMatSize 
)

Prepare assembly of matrix - reserve space in underlying coordinate matrix object.

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.

Parameters
toltolerance on relative residual ||res||/||rhs||
numLevelsnumber of levels
numSubAtLevelsnumber of subdomains at levels
verboseLevellevel of verbosity of BDDCML library ( 0 - only fatal errors reported, 1 - mild output, 2 - detailed output )
maxItmaximum number of iterations
ndecrMaxmaximum number of iterations with non-decreasing residual ( used to stop diverging process )
use_adaptiveif adaptive constraints should be used by BDDCML
void la::BddcmlWrapper::writeMatrix ( std::ostream &  out)
inline

Outputs.

Definition at line 141 of file bddcml_wrapper.hpp.

Here is the caller graph for this function:

Member Data Documentation

const MPI_Comm la::BddcmlWrapper::comm_
private

communicator

Definition at line 242 of file bddcml_wrapper.hpp.

double la::BddcmlWrapper::condNumber_
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.hpp.

int la::BddcmlWrapper::convergedReason_
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.hpp.

la::MatrixCoo<int,double> la::BddcmlWrapper::coo_
private

matrix in coordinate format (COO)

Definition at line 278 of file bddcml_wrapper.hpp.

la::MatrixCoo<int,double> la::BddcmlWrapper::diagWeightsCoo_
private

diagonal of subdomain matrix in sparse format and global numbers

Definition at line 280 of file bddcml_wrapper.hpp.

std::vector<double> la::BddcmlWrapper::element_data_
private

array with data on elements, here used for permeability for each element

Definition at line 265 of file bddcml_wrapper.hpp.

std::vector<double> la::BddcmlWrapper::fixv_
private

values of fixed boundary conditions at places of ifix_ - values outside fixed dofs are ignored

Definition at line 276 of file bddcml_wrapper.hpp.

Global2LocalMap_ la::BddcmlWrapper::global2LocalDofMap_
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.hpp.

std::vector<int> la::BddcmlWrapper::ifix_
private

indices of fixed boundary conditions - ( 0 for free variable, 1 for constrained dof )

Definition at line 275 of file bddcml_wrapper.hpp.

std::vector<int> la::BddcmlWrapper::inet_
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.hpp.

bool la::BddcmlWrapper::isDiagAssembled_
private

true if diagonal is assembled

Definition at line 281 of file bddcml_wrapper.hpp.

std::vector<int> la::BddcmlWrapper::isegn_
private

Indices of Subdomain Elements in Global Numbering.

Definition at line 269 of file bddcml_wrapper.hpp.

bool la::BddcmlWrapper::isMatAssembled_
private

true if matrix is assembled

Definition at line 279 of file bddcml_wrapper.hpp.

std::vector<int> la::BddcmlWrapper::isngn_
private

Indices of Subdomain Nodes in Global Numbering.

Definition at line 267 of file bddcml_wrapper.hpp.

std::vector<int> la::BddcmlWrapper::isvgvn_
private

Indices of Subdomain Variables in Global Variable Numbering ( i.e. dof mapping )

Definition at line 268 of file bddcml_wrapper.hpp.

const MatrixType la::BddcmlWrapper::matrixType_
private

type of matrix

Definition at line 241 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::meshDim_
private

mesh dimension, may differ from space dimension - e.g. 2 for shells in 3D

Definition at line 249 of file bddcml_wrapper.hpp.

bool la::BddcmlWrapper::meshLoaded_
private

Definition at line 256 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::nDim_
private

space dimension

Definition at line 248 of file bddcml_wrapper.hpp.

std::vector<int> la::BddcmlWrapper::nndf_
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.hpp.

std::vector<int> la::BddcmlWrapper::nnet_
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.hpp.

double la::BddcmlWrapper::normRhs_
private

norm of the global right-hand side

Definition at line 284 of file bddcml_wrapper.hpp.

double la::BddcmlWrapper::normSol_
private

norm of the global solution

Definition at line 287 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::nProc_
private

number of processors in communicator

Definition at line 245 of file bddcml_wrapper.hpp.

const int la::BddcmlWrapper::numDofs_
private

number of total dofs

Definition at line 239 of file bddcml_wrapper.hpp.

const int la::BddcmlWrapper::numDofsSub_
private

number of subdomain dofs

Definition at line 240 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::numElem_
private

total number of elements

Definition at line 252 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::numElemSub_
private

number of elements in subdomain

Definition at line 251 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::numIter_
private

required number of iterations

Definition at line 289 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::numNodes_
private

total number of nodes

Definition at line 254 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::numNodesSub_
private

number of nodes in subdomain

Definition at line 253 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::numSubLoc_
private

number of subdomains on process

Definition at line 243 of file bddcml_wrapper.hpp.

int la::BddcmlWrapper::rank_
private

ID of processors in communicator.

Definition at line 246 of file bddcml_wrapper.hpp.

std::vector<double> la::BddcmlWrapper::rhsVec_
private

vector with RHS values restricted to subdomain, i.e. values are repeated at shared nodes

Definition at line 283 of file bddcml_wrapper.hpp.

std::vector<double> la::BddcmlWrapper::sol_
private

distributed solution vector

Definition at line 286 of file bddcml_wrapper.hpp.

std::vector<double> la::BddcmlWrapper::xyz_
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.hpp.


The documentation for this class was generated from the following file: