24 #ifdef FLOW123D_HAVE_BDDCML 27 #endif // FLOW123D_HAVE_BDDCML 39 return it::Record(
"Bddc",
"BDDCML (Balancing Domain Decomposition by Constraints - Multi-Level) solver settings.")
42 "If not, we use the value 1.0e-7."),
43 "Residual tolerance relative to the initial error.")
45 "If not, we use the value 1000."),
46 "Maximum number of outer iterations of the linear solver.")
49 "Maximum number of iterations of the linear solver with non-decreasing residual.")
51 "Number of levels in the multilevel method (=2 for the standard BDDC).")
53 "Use adaptive selection of constraints in BDDCML.")
55 "Level of verbosity of the BDDCML library:\n\n - 0 - no output,\n - 1 - mild output,\n - 2 - detailed output.")
65 const bool swap_sign )
69 #ifdef FLOW123D_HAVE_BDDCML 73 throw ExcInputMessage << EI_Message(
"Unsupported solver BDDC. Compiled without support for the BDDCML solver.")
'' 74 #endif // FLOW123D_HAVE_BDDCML 91 const int nDim,
const int numNodes,
FMT_UNUSED const int numDofs,
102 #ifdef FLOW123D_HAVE_BDDCML 104 uint num_of_local_subdomains = 1;
109 num_of_local_subdomains);
116 isngn_.resize(isngn.size());
117 std::copy( isngn.begin(), isngn.end(),
isngn_.begin() );
118 OLD_ASSERT( numDofs == static_cast<int>(
size_),
"Global problem size mismatch!" );
120 bddcml_ -> loadRawMesh( nDim, numNodes, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
126 PetscInt numDofsSubInt =
static_cast<int>(
isngn_.size( ) );
129 ISLocalToGlobalMapping subdomainMapping;
130 ierr = ISLocalToGlobalMappingCreate(
comm_, 1, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
132 IS subdomainIndexSet;
134 ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
135 ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
140 chkerr(ISDestroy( &subdomainIndexSet ));
141 chkerr(ISDestroy( &from ));
151 #endif // FLOW123D_HAVE_BDDCML 156 #ifdef FLOW123D_HAVE_BDDCML 157 namespace ublas = boost::numeric::ublas;
161 ublas::matrix< double >
mat( nrow, ncol );
163 std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
164 std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
166 for (
int i = 0; i < nrow; i++ ) {
167 for (
int j = 0; j < ncol; j++ ) {
168 mat( i, j ) = vals[i*ncol + j];
175 bddcml_ -> insertToMatrix( mat, myRows, myCols );
176 #endif // FLOW123D_HAVE_BDDCML 181 #ifdef FLOW123D_HAVE_BDDCML 182 namespace ublas = boost::numeric::ublas;
185 ublas::vector< double >
vec( nrow );
187 std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
189 for (
int i = 0; i < nrow; i++ ) {
196 bddcml_ -> insertToRhs( vec, myRows );
197 #endif // FLOW123D_HAVE_BDDCML 202 #ifdef FLOW123D_HAVE_BDDCML 203 bddcml_ -> insertToDiagonalWeights( global_index, value );
204 #endif // FLOW123D_HAVE_BDDCML 209 #ifdef FLOW123D_HAVE_BDDCML 211 #endif // FLOW123D_HAVE_BDDCML 217 #ifdef FLOW123D_HAVE_BDDCML 219 #endif // FLOW123D_HAVE_BDDCML 225 #ifdef FLOW123D_HAVE_BDDCML 226 bddcml_ -> finishMatAssembly( );
227 #endif // FLOW123D_HAVE_BDDCML 232 #ifdef FLOW123D_HAVE_BDDCML 234 #endif // FLOW123D_HAVE_BDDCML 239 #ifdef FLOW123D_HAVE_BDDCML 249 LogOut().fmt(
"BDDCML converged reason: {} ( 0 means OK ) \n",
bddcml_ -> giveConvergedReason() );
250 LogOut().fmt(
"BDDCML converged in {} iterations. \n",
bddcml_ -> giveNumIterations() );
251 LogOut().fmt(
"BDDCML estimated condition number is {} \n",
bddcml_ -> giveCondNumber() );
259 double * locSolVecArray;
262 VecRestoreArray(
locSolVec_, &locSolVecArray );
274 #endif // FLOW123D_HAVE_BDDCML 291 #ifdef FLOW123D_HAVE_BDDCML 300 #endif // FLOW123D_HAVE_BDDCML 375 <<
"zzz(:,1:2) = zzz(:,1:2) + ones(size(zzz),2);\n" 376 <<
"matrix_bddc = spconvert(zzz);\n";
int bddcml_verbosity_level_
SetValuesMode status_
Set value status of the linear system.
LinSys_BDDC(const Distribution *rows_ds, const bool swap_sign=false)
double normRhs()
Get norm of right-hand side.
double get_solution_precision() override
VecScatter VSpetscToSubScatter_
scatter from solution_ to locSolVec_
std::vector< int > isngn_
indices of subdomain nodes in global numbering
PetscErrorCode rhs_zero_entries() override
virtual void set_from_input(const Input::Record in_rec)
Wrappers for linear systems based on MPIAIJ and MATIS format.
ConstraintVec_ constraints_
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
void print_matrix(std::ostream &out)
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
unsigned size_
global number of matrix rows, i.e. problem size
#define LogOut()
Macro defining 'log' record of log.
static const Input::Type::Record & get_input_type()
bool use_adaptive_bddc_
should adaptive BDDC be used?
static constexpr bool value
void set_from_input(const Input::Record in_rec) override
void apply_constrains(double scalar=1.) override
ArmaMat< double, N, M > mat
unsigned int max_it_
maximum number of iterations of linear solver
Vec locSolVec_
local solution PETSc vector - sequential
#define START_TIMER(tag)
Starts a timer with specified tag.
const Distribution * rows_ds_
final distribution of rows of MH matrix
int max_nondecr_it_
parameters expected from input file:
double a_tol_
absolute tolerance of linear solver
static const int registrar
Registrar of class to factory.
void diagonal_weights_set_value(int global_index, double value)
double residual_norm_
local solution array pointing into Vec solution_
const bool swap_sign_
swap sign of matrix and rhs entries, e.g. to make the matrix SPD
int number_of_levels_
number of levels in the multilevel method
la::BddcmlWrapper Bddcml_
la::BddcmlWrapper::MatrixType BDDCMatrixType
double r_tol_
relative tolerance of linear solver
Bddcml_ * bddcml_
BDDCML wrapper.
PetscErrorCode mat_zero_entries() override
MPI_Comm get_comm() const
Returns communicator.
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
void finish_assembly() override
static Input::Type::Abstract & get_input_type()
std::vector< double > locSolution_
subdomain solution
void writeMatrix(std::ostream &out)
Outputs.
void load_mesh(BDDCMatrixType matrix_type, 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)
Vec solution_
PETSc vector constructed with vb array.
LinSys::SolveInfo solve() override
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Sets tolerances. Note that BDDC does not use a_tol.
void rhs_set_values(int nrow, int *rows, double *vals) override
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.