Flow123d
release_3.0.0-968-gc87a28e79
|
Go to the documentation of this file.
23 #ifdef FLOW123D_HAVE_BDDCML
26 #endif // FLOW123D_HAVE_BDDCML
38 return it::Record(
"Bddc",
"BDDCML (Balancing Domain Decomposition by Constraints - Multi-Level) solver settings.")
41 "If not, we use the value 1.0e-7."),
42 "Residual tolerance relative to the initial error.")
44 "If not, we use the value 1000."),
45 "Maximum number of outer iterations of the linear solver.")
48 "Maximum number of iterations of the linear solver with non-decreasing residual.")
50 "Number of levels in the multilevel method (=2 for the standard BDDC).")
52 "Use adaptive selection of constraints in BDDCML.")
54 "Level of verbosity of the BDDCML library:\n\n - 0 - no output,\n - 1 - mild output,\n - 2 - detailed output.")
65 const int matrixTypeInt,
67 const bool swap_sign )
71 #ifdef FLOW123D_HAVE_BDDCML
76 switch ( matrixTypeInt ) {
91 OLD_ASSERT(
true,
"Unknown matrix type %d", matrixTypeInt );
110 xprintf(
UsrErr,
"Compiled without support for BDDCML solver.\n");
111 #endif // FLOW123D_HAVE_BDDCML
138 #ifdef FLOW123D_HAVE_BDDCML
140 isngn_.resize(isngn.size());
141 std::copy( isngn.begin(), isngn.end(),
isngn_.begin() );
144 bddcml_ -> loadRawMesh( nDim, numNodes, numDofs, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
150 PetscInt numDofsSubInt =
static_cast<int>(
isngn_.size( ) );
153 ISLocalToGlobalMapping subdomainMapping;
154 ierr = ISLocalToGlobalMappingCreate(
comm_, 1, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
156 IS subdomainIndexSet;
158 ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
159 ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
164 chkerr(ISDestroy( &subdomainIndexSet ));
165 chkerr(ISDestroy( &from ));
175 #endif // FLOW123D_HAVE_BDDCML
180 #ifdef FLOW123D_HAVE_BDDCML
181 namespace ublas = boost::numeric::ublas;
185 ublas::matrix< double > mat( nrow, ncol );
187 std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
188 std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
190 for (
int i = 0; i < nrow; i++ ) {
191 for (
int j = 0; j < ncol; j++ ) {
192 mat( i, j ) = vals[i*ncol + j];
199 bddcml_ -> insertToMatrix( mat, myRows, myCols );
200 #endif // FLOW123D_HAVE_BDDCML
205 #ifdef FLOW123D_HAVE_BDDCML
206 namespace ublas = boost::numeric::ublas;
209 ublas::vector< double > vec( nrow );
211 std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
213 for (
int i = 0; i < nrow; i++ ) {
220 bddcml_ -> insertToRhs( vec, myRows );
221 #endif // FLOW123D_HAVE_BDDCML
226 #ifdef FLOW123D_HAVE_BDDCML
227 bddcml_ -> insertToDiagonalWeights( global_index,
value );
228 #endif // FLOW123D_HAVE_BDDCML
233 #ifdef FLOW123D_HAVE_BDDCML
235 #endif // FLOW123D_HAVE_BDDCML
241 #ifdef FLOW123D_HAVE_BDDCML
243 #endif // FLOW123D_HAVE_BDDCML
249 #ifdef FLOW123D_HAVE_BDDCML
250 bddcml_ -> finishMatAssembly( );
251 #endif // FLOW123D_HAVE_BDDCML
256 #ifdef FLOW123D_HAVE_BDDCML
258 #endif // FLOW123D_HAVE_BDDCML
263 #ifdef FLOW123D_HAVE_BDDCML
273 LogOut().fmt(
"BDDCML converged reason: {} ( 0 means OK ) \n",
bddcml_ -> giveConvergedReason() );
274 LogOut().fmt(
"BDDCML converged in {} iterations. \n",
bddcml_ -> giveNumIterations() );
275 LogOut().fmt(
"BDDCML estimated condition number is {} \n",
bddcml_ -> giveCondNumber() );
283 double * locSolVecArray;
286 VecRestoreArray(
locSolVec_, &locSolVecArray );
298 #endif // FLOW123D_HAVE_BDDCML
315 #ifdef FLOW123D_HAVE_BDDCML
325 #endif // FLOW123D_HAVE_BDDCML
400 <<
"zzz(:,1:2) = zzz(:,1:2) + ones(size(zzz),2);\n"
401 <<
"matrix_bddc = spconvert(zzz);\n";
Vec locSolVec_
local solution PETSc vector - sequential
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
std::vector< int > isngn_
indices of subdomain nodes in global numbering
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
double residual_norm_
local solution array pointing into Vec solution_
int number_of_levels_
number of levels in the multilevel method
bool use_adaptive_bddc_
should adaptive BDDC be used?
static const Input::Type::Record & get_input_type()
VecScatter VSpetscToSubScatter_
scatter from solution_ to locSolVec_
static constexpr bool value
@ SYMMETRICGENERAL
general symmetric ( e.g. saddle point ),
Wrappers for linear systems based on MPIAIJ and MATIS format.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
int bddcml_verbosity_level_
double normRhs()
Get norm of right-hand side.
void finish_assembly() override
Vec solution_
PETSc vector constructed with vb array.
#define LogOut()
Macro defining 'log' record of log.
unsigned size_
global number of matrix rows, i.e. problem size
void rhs_set_values(int nrow, int *rows, double *vals) override
double a_tol_
absolute tolerance of linear solver
la::BddcmlWrapper Bddcml_
enum matrixTypeEnum MatrixType
std::vector< double > locSolution_
subdomain solution
SetValuesMode status_
Set value status of the linear system.
@ SPD_VIA_SYMMETRICGENERAL
@ GENERAL
general (default),
double get_solution_precision() override
PetscErrorCode rhs_zero_entries() override
void set_from_input(const Input::Record in_rec) override
virtual void set_from_input(const Input::Record in_rec)
double r_tol_
relative tolerance of linear solver
static const int registrar
Registrar of class to factory.
ConstraintVec_ constraints_
int max_nondecr_it_
parameters expected from input file:
LinSys_BDDC(const unsigned numDofsSub, const Distribution *rows_ds, const int matrixTypeInt=0, const int numSubLoc=1, const bool swap_sign=false)
LinSys::SolveInfo solve() override
const bool swap_sign_
swap sign of matrix and rhs entries, e.g. to make the matrix SPD
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
Bddcml_ * bddcml_
BDDCML wrapper.
static Input::Type::Abstract & get_input_type()
void writeMatrix(std::ostream &out)
Outputs.
PetscErrorCode mat_zero_entries() override
void load_mesh(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_permeability, const int meshDim)
void apply_constrains(double scalar=1.) override
#define START_TIMER(tag)
Starts a timer with specified tag.
MPI_Comm get_comm() const
Returns communicator.
void print_matrix(std::ostream &out)
unsigned int max_it_
maximum number of iterations of linear solver
void diagonal_weights_set_value(int global_index, double value)
@ SPD
symmetric positive definite,