Flow123d
JS_before_hm-2199-g85616a0ec
|
Go to the documentation of this file.
24 #ifdef FLOW123D_HAVE_BDDCML
27 #endif // FLOW123D_HAVE_BDDCML
40 return it::Record(
"Bddc",
"BDDCML (Balancing Domain Decomposition by Constraints - Multi-Level) solver settings.")
43 "If not, we use the value 1.0e-7."),
44 "Residual tolerance relative to the initial error.")
46 "If not, we use the value 1000."),
47 "Maximum number of outer iterations of the linear solver.")
50 "Maximum number of iterations of the linear solver with non-decreasing residual.")
52 "Number of levels in the multilevel method (=2 for the standard BDDC).")
54 "Use adaptive selection of constraints in BDDCML.")
56 "Level of verbosity of the BDDCML library:\n\n - 0 - no output,\n - 1 - mild output,\n - 2 - detailed output.")
66 const bool swap_sign )
70 #ifdef FLOW123D_HAVE_BDDCML
74 throw Input::ExcInputMessage() << EI_Message(
"Unsupported solver BDDC. Compiled without support for the BDDCML solver.");
75 #endif // FLOW123D_HAVE_BDDCML
92 const int nDim,
const int numNodes,
FMT_UNUSED const int numDofs,
103 #ifdef FLOW123D_HAVE_BDDCML
105 uint num_of_local_subdomains = 1;
110 num_of_local_subdomains);
117 isngn_.resize(isngn.size());
118 std::copy( isngn.begin(), isngn.end(),
isngn_.begin() );
119 OLD_ASSERT( numDofs ==
static_cast<int>(
size_),
"Global problem size mismatch!" );
121 bddcml_ -> loadRawMesh( nDim, numNodes, inet, nnet, nndf, isegn, isngn, isvgvn, xyz, element_permeability, meshDim );
127 PetscInt numDofsSubInt =
static_cast<int>(
isngn_.size( ) );
130 ISLocalToGlobalMapping subdomainMapping;
131 ierr = ISLocalToGlobalMappingCreate(
comm_, 1, numDofsSubInt, &(idx[0]), PETSC_COPY_VALUES, &subdomainMapping ); CHKERRV( ierr );
133 IS subdomainIndexSet;
135 ierr = ISCreateStride( PETSC_COMM_SELF, numDofsSubInt, 0, 1, &subdomainIndexSet );
136 ierr = ISLocalToGlobalMappingApplyIS( subdomainMapping, subdomainIndexSet, &from );
141 chkerr(ISDestroy( &subdomainIndexSet ));
142 chkerr(ISDestroy( &from ));
152 #endif // FLOW123D_HAVE_BDDCML
157 #ifdef FLOW123D_HAVE_BDDCML
158 namespace ublas = boost::numeric::ublas;
162 ublas::matrix< double >
mat( nrow, ncol );
164 std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
165 std::copy( &(cols[0]), &(cols[ncol]), myCols.begin() );
167 for (
int i = 0; i < nrow; i++ ) {
168 for (
int j = 0; j < ncol; j++ ) {
169 mat( i, j ) = vals[i*ncol + j];
176 bddcml_ -> insertToMatrix(
mat, myRows, myCols );
177 #endif // FLOW123D_HAVE_BDDCML
182 #ifdef FLOW123D_HAVE_BDDCML
183 namespace ublas = boost::numeric::ublas;
186 ublas::vector< double >
vec( nrow );
188 std::copy( &(rows[0]), &(rows[nrow]), myRows.begin() );
190 for (
int i = 0; i < nrow; i++ ) {
198 #endif // FLOW123D_HAVE_BDDCML
203 #ifdef FLOW123D_HAVE_BDDCML
204 bddcml_ -> insertToDiagonalWeights( global_index,
value );
205 #endif // FLOW123D_HAVE_BDDCML
210 #ifdef FLOW123D_HAVE_BDDCML
212 #endif // FLOW123D_HAVE_BDDCML
218 #ifdef FLOW123D_HAVE_BDDCML
220 #endif // FLOW123D_HAVE_BDDCML
226 #ifdef FLOW123D_HAVE_BDDCML
227 bddcml_ -> finishMatAssembly( );
228 #endif // FLOW123D_HAVE_BDDCML
233 #ifdef FLOW123D_HAVE_BDDCML
235 #endif // FLOW123D_HAVE_BDDCML
240 #ifdef FLOW123D_HAVE_BDDCML
250 LogOut().fmt(
"BDDCML converged reason: {} ( 0 means OK ) \n",
bddcml_ -> giveConvergedReason() );
251 LogOut().fmt(
"BDDCML converged in {} iterations. \n",
bddcml_ -> giveNumIterations() );
252 LogOut().fmt(
"BDDCML estimated condition number is {} \n",
bddcml_ -> giveCondNumber() );
260 double * locSolVecArray;
263 VecRestoreArray(
locSolVec_, &locSolVecArray );
275 #endif // FLOW123D_HAVE_BDDCML
292 #ifdef FLOW123D_HAVE_BDDCML
301 #endif // FLOW123D_HAVE_BDDCML
376 <<
"zzz(:,1:2) = zzz(:,1:2) + ones(size(zzz),2);\n"
377 <<
"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
Sets tolerances. Note that BDDC does not use a_tol.
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
Wrappers for linear systems based on MPIAIJ and MATIS format.
la::BddcmlWrapper::MatrixType BDDCMatrixType
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
ArmaMat< double, N, M > mat
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.
#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_
std::vector< double > locSolution_
subdomain solution
SetValuesMode status_
Set value status of the linear system.
double get_solution_precision() override
PetscErrorCode rhs_zero_entries() override
void set_from_input(const Input::Record in_rec) override
const Distribution * rows_ds_
final distribution of rows of MH matrix
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::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
LinSys_BDDC(const Distribution *rows_ds, const bool swap_sign=false)
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)