Flow123d
JS_before_hm-1822-gbb48b12e9
|
Go to the documentation of this file.
33 return it::Record(
"Petsc",
"PETSc solver settings.\n It provides interface to various PETSc solvers. The convergence criteria is:\n"
35 "norm( res_i ) < max( norm( res_0 ) * r_tol, a_tol )\n"
37 "where ```res_i``` is the residuum vector after i-th iteration of the solver and ```res_0``` is the estimate of the norm of the initial residual. "
38 "If the initial guess of the solution is provided (usually only for transient equations) the residual of this estimate is used, "
39 "otherwise the norm of preconditioned RHS is used. "
40 "The default norm is (($L_2$)) norm of preconditioned residual: (($ P^{-1}(Ax-b)$)), usage of other norm may be prescribed using the 'option' key. "
41 "See also PETSc documentation for KSPSetNormType.")
44 "If not, we use the value 1.0e-7."),
45 "Residual tolerance relative to the initial error.")
47 "If not, we use the value 1.0e-11."),
48 "Absolute residual tolerance.")
50 "If not, we use the value 1000."),
51 "Maximum number of outer iterations of the linear solver.")
53 "If the string is left empty (by default), the internal default options is used.")
64 init_guess_nonzero(false),
72 ierr = VecZeroEntries(
rhs_ ); CHKERRV( ierr );
82 :
LinSys(other), params_(other.params_), v_rhs_(NULL), solution_precision_(other.solution_precision_)
128 OLD_ASSERT(
false,
"Can not set values. Matrix is not preallocated.\n");
146 OLD_ASSERT(
false,
"Can not set values. Matrix is not preallocated.\n");
162 default:
DebugOut() <<
"LS SetValues with non allowed insert mode.\n";
175 ierr = VecSetValues(
rhs_,nrow,rows,vals,(InsertMode)
status_); CHKERRV( ierr );
179 default:
OLD_ASSERT(
false,
"LinSys's status disallow set values.\n");
191 for (i=0; i<nrow; i++) {
193 for(j=0; j<ncol; j++) {
196 VecSetValue(
on_vec_,row,1.0,ADD_VALUES);
198 VecSetValue(
off_vec_,row,1.0,ADD_VALUES);
208 PetscInt *on_nz, *off_nz;
209 PetscScalar *on_array, *off_array;
221 VecGetArray(
on_vec_, &on_array );
222 VecGetArray(
off_vec_, &off_array );
229 VecRestoreArray(
on_vec_,&on_array);
230 VecRestoreArray(
off_vec_,&off_array);
240 0, on_nz, 0, off_nz, &
matrix_); CHKERRV( ierr );
243 MatSetOption(
matrix_, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
248 MatSetOption(
matrix_, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
258 MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
267 WarningOut() <<
"Finalizing linear system without setting values.\n";
270 ierr = MatAssemblyBegin(
matrix_, assembly_type); CHKERRV( ierr );
271 ierr = VecAssemblyBegin(
rhs_); CHKERRV( ierr );
272 ierr = MatAssemblyEnd(
matrix_, assembly_type); CHKERRV( ierr );
273 ierr = VecAssemblyEnd(
rhs_); CHKERRV( ierr );
275 if (assembly_type == MAT_FINAL_ASSEMBLY)
status_ =
DONE;
291 OLD_ASSERT (
status_ ==
DONE,
"System matrix and right-hand side are not assembled when applying constraints." );
294 PetscInt numConstraints =
static_cast<PetscInt
>(
constraints_.size() );
297 const PetscScalar diagScalar =
static_cast<PetscScalar
>( scalar );
303 ConstraintVec_::const_iterator cIter =
constraints_.begin( );
304 ConstraintVec_::const_iterator cEnd =
constraints_.end( );
306 for ( ; cIter != cEnd; ++cIter ) {
307 globalDofs.push_back(
static_cast<PetscInt
>( cIter -> first ) );
308 values.push_back(
static_cast<PetscScalar
>( cIter -> second ) * diagScalar );
316 ierr = MatZeroRows(
matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
320 if ( numConstraints ) {
321 ierr = VecSetValues(
rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
326 ierr = VecAssemblyBegin(
rhs_ ); CHKERRV( ierr );
327 ierr = VecAssemblyEnd(
rhs_ ); CHKERRV( ierr );
340 const char *petsc_dflt_opt;
355 petsc_dflt_opt=
"-ksp_type cg -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_type asm -pc_asm_type basic -pc_asm_overlap 4 -sub_pc_type icc -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
358 petsc_dflt_opt=
"-ksp_type bcgs -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_type asm -pc_asm_overlap 4 -sub_pc_type ilu -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
364 petsc_dflt_opt=
"-ksp_type cg -pc_type icc -pc_factor_levels 3 -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
367 petsc_dflt_opt=
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
371 LogOut().fmt(
"inserting petsc options: {}\n",
params_.c_str());
375 PetscOptionsInsertString(NULL,
params_.c_str());
377 MatSetOption(
matrix_, MAT_USE_INODES, PETSC_FALSE );
386 KSPSetFromOptions(
system);
393 KSPGetType(
system, &type);
394 if (strcmp(type, KSPPREONLY) != 0)
395 KSPSetInitialGuessNonzero(
system, PETSC_TRUE);
403 KSPGetIterationNumber(
system,&nits);
409 LogOut().fmt(
"convergence reason {}, number of iterations is {}\n",
reason, nits);
425 FilePath matFileName(text +
"_flow123d_matrix.m",FilePath::FileType::output_file);
426 FilePath rhsFileName(text +
"_flow123d_rhs.m",FilePath::FileType::output_file);
427 FilePath solFileName(text +
"_flow123d_sol.m",FilePath::FileType::output_file);
429 PetscViewer myViewer;
432 PetscViewerASCIIOpen(
comm_,((
string)matFileName).c_str(),&myViewer);
433 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
435 PetscViewerDestroy(&myViewer);
438 WarningOut() <<
"PetscViewer: the matrix of LinSys is not set.\n";
440 if (
rhs_ != NULL ) {
441 PetscViewerASCIIOpen(
comm_,((
string)rhsFileName).c_str(),&myViewer);
442 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
443 VecView(
rhs_, myViewer );
444 PetscViewerDestroy(&myViewer);
447 WarningOut() <<
"PetscViewer: the rhs of LinSys is not set.\n";
450 PetscViewerASCIIOpen(
comm_,((
string)solFileName).c_str(),&myViewer);
451 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
453 PetscViewerDestroy(&myViewer);
456 WarningOut() <<
"PetscViewer: the solution of LinSys is not set.\n";
477 std::string user_params = in_rec.
val<
string>(
"options");
478 if (user_params !=
"")
params_ = user_params;
492 double residual_norm;
493 VecNorm(
residual_, NORM_2, &residual_norm);
494 return residual_norm;
unsigned int np() const
get num of processors
unsigned int lsize(int proc) const
get local size
void start_insert_assembly() override
double residual_norm_
local solution array pointing into Vec solution_
void set_initial_guess_nonzero(bool set_nonzero=true)
std::string params_
command-line-like options for the PETSc solver
bool init_guess_nonzero
flag for starting from nonzero guess
LinSys_PETSC(const Distribution *rows_ds, const std::string ¶ms="")
void apply_constrains(double scalar=1.) override
KSPConvergedReason reason
Dedicated class for storing path to input and output files.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
unsigned int size() const
get global size
LinSys::SolveInfo solve() override
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Mat matrix_
Petsc matrix of the problem.
static const int registrar
Registrar of class to factory.
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
unsigned int get_proc(unsigned int idx) const
get processor of the given index
Vec solution_
PETSc vector constructed with vb array.
#define LogOut()
Macro defining 'log' record of log.
double a_tol_
absolute tolerance of linear solver
double * v_rhs_
local RHS array pointing to Vec rhs_
void finish_assembly() override
double solution_precision_
bool rhs_changed_
true if the right hand side was changed since the last solve
double get_solution_precision() override
void rhs_set_values(int nrow, int *rows, double *vals) override
SetValuesMode status_
Set value status of the linear system.
Vec on_vec_
Vectors for counting non-zero entries in diagonal block.
void preallocate_matrix()
const Distribution * rows_ds_
final distribution of rows of MH matrix
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
virtual void set_from_input(const Input::Record in_rec)
void start_add_assembly() override
double r_tol_
relative tolerance of linear solver
void set_from_input(const Input::Record in_rec) override
ConstraintVec_ constraints_
bool matrix_changed_
true if the matrix was changed since the last solve
#define WarningOut()
Macro defining 'warning' record of log.
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
T * makePetscPointer_(std::vector< T > &array)
void start_allocation() override
static Input::Type::Abstract & get_input_type()
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
void view(string text="") override
bool is_positive_definite()
#define DebugOut()
Macro defining 'debug' record of log.
#define START_TIMER(tag)
Starts a timer with specified tag.
Vec rhs_
PETSc vector constructed with vx array.
double compute_residual() override
static const Input::Type::Record & get_input_type()
unsigned int max_it_
maximum number of iterations of linear solver