32 return it::Record(
"Petsc",
"Interface to PETSc solvers. Convergence criteria is:\n" 34 "norm( res_n ) < max( norm( res_0 ) * r_tol, a_tol )\n" 36 "where res_i is the residuum vector after i-th iteration of the solver and res_0 is an estimate of the norm of initial residual. " 37 "If the initial guess of the solution is provided (usually only for transient equations) the residual of this estimate is used, " 38 "otherwise the norm of preconditioned RHS is used. " 39 "The default norm is L2 norm of preconditioned residual: (($ P^{-1}(Ax-b)$)), usage of other norm may be prescribed using the 'option' key. " 40 "See also PETSc documentation for KSPSetNormType.")
43 "Relative residual tolerance, (to initial error).")
45 "Absolute residual tolerance.")
47 "Maximum number of outer iterations of the linear solver.")
66 ierr = VecZeroEntries(
rhs_ ); CHKERRV( ierr );
123 OLD_ASSERT(
false,
"Can not set values. Matrix is not preallocated.\n");
141 OLD_ASSERT(
false,
"Can not set values. Matrix is not preallocated.\n");
154 ierr = MatSetValues(
matrix_,nrow,rows,ncol,cols,vals,(InsertMode)
status_); CHKERRV( ierr );
159 default:
DebugOut() <<
"LS SetValues with non allowed insert mode.\n";
172 ierr = VecSetValues(
rhs_,nrow,rows,vals,(InsertMode)
status_); CHKERRV( ierr );
176 default:
OLD_ASSERT(
false,
"LinSys's status disallow set values.\n");
188 for (i=0; i<nrow; i++) {
190 for(j=0; j<ncol; j++) {
193 VecSetValue(
on_vec_,row,1.0,ADD_VALUES);
195 VecSetValue(
off_vec_,row,1.0,ADD_VALUES);
205 PetscInt *on_nz, *off_nz;
206 PetscScalar *on_array, *off_array;
218 VecGetArray(
on_vec_, &on_array );
219 VecGetArray(
off_vec_, &off_array );
222 on_nz[i] =
static_cast<PetscInt
>( on_array[i]+0.1 );
223 off_nz[i] =
static_cast<PetscInt
>( off_array[i]+0.1 );
226 VecRestoreArray(
on_vec_,&on_array);
227 VecRestoreArray(
off_vec_,&off_array);
234 ierr = MatDestroy(&
matrix_); CHKERRV( ierr );
237 0, on_nz, 0, off_nz, &
matrix_); CHKERRV( ierr );
240 MatSetOption(
matrix_, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
248 MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
257 WarningOut() <<
"Finalizing linear system without setting values.\n";
260 ierr = MatAssemblyBegin(
matrix_, assembly_type); CHKERRV( ierr );
261 ierr = VecAssemblyBegin(
rhs_); CHKERRV( ierr );
262 ierr = MatAssemblyEnd(
matrix_, assembly_type); CHKERRV( ierr );
263 ierr = VecAssemblyEnd(
rhs_); CHKERRV( ierr );
265 if (assembly_type == MAT_FINAL_ASSEMBLY)
status_ =
DONE;
276 OLD_ASSERT (
status_ ==
DONE,
"System matrix and right-hand side are not assembled when applying constraints." );
279 PetscInt numConstraints =
static_cast<PetscInt
>(
constraints_.size() );
282 const PetscScalar diagScalar =
static_cast<PetscScalar
>( scalar );
288 ConstraintVec_::const_iterator cIter =
constraints_.begin( );
289 ConstraintVec_::const_iterator cEnd =
constraints_.end( );
291 for ( ; cIter != cEnd; ++cIter ) {
292 globalDofs.push_back( static_cast<PetscInt>( cIter -> first ) );
293 values.push_back( static_cast<PetscScalar>( cIter -> second ) * diagScalar );
301 ierr = MatZeroRows(
matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
305 if ( numConstraints ) {
306 ierr = VecSetValues(
rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
311 ierr = VecAssemblyBegin(
rhs_ ); CHKERRV( ierr );
312 ierr = VecAssemblyEnd(
rhs_ ); CHKERRV( ierr );
325 const char *petsc_dflt_opt;
334 petsc_dflt_opt=
"-ksp_type bcgs -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";
336 petsc_dflt_opt=
"-ksp_type bcgs -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";
343 petsc_dflt_opt=
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
345 petsc_dflt_opt=
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
349 LogOut().fmt(
"inserting petsc options: {}\n",
params_.c_str());
350 PetscOptionsInsertString(
params_.c_str());
352 MatSetOption(
matrix_, MAT_USE_INODES, PETSC_FALSE );
361 KSPSetFromOptions(
system);
368 KSPGetType(
system, &type);
369 if (strcmp(type, KSPPREONLY) != 0)
370 KSPSetInitialGuessNonzero(
system, PETSC_TRUE);
378 KSPGetIterationNumber(
system,&nits);
384 LogOut().fmt(
"convergence reason {}, number of iterations is {}\n",
reason, nits);
394 return static_cast<int>(
reason);
400 std::string matFileName =
"flow123d_matrix.m";
401 std::string rhsFileName =
"flow123d_rhs.m";
402 std::string solFileName =
"flow123d_sol.m";
404 PetscViewer myViewer;
406 PetscViewerASCIIOpen(
comm_,matFileName.c_str(),&myViewer);
407 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
409 PetscViewerDestroy(&myViewer);
411 PetscViewerASCIIOpen(
comm_,rhsFileName.c_str(),&myViewer);
412 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
413 VecView(
rhs_, myViewer );
414 PetscViewerDestroy(&myViewer);
417 PetscViewerASCIIOpen(
comm_,solFileName.c_str(),&myViewer);
418 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
420 PetscViewerDestroy(&myViewer);
428 if (
matrix_ != NULL) { ierr = MatDestroy(&
matrix_); CHKERRV( ierr ); }
429 ierr = VecDestroy(&
rhs_); CHKERRV( ierr );
455 double residual_norm;
456 VecNorm(
residual_, NORM_2, &residual_norm);
457 return residual_norm;
unsigned int get_proc(unsigned int idx) const
get processor of the given index
SetValuesMode status_
Set value status of the linear system.
void apply_constrains(double scalar=1.) override
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
virtual void set_from_input(const Input::Record in_rec)
bool matrix_changed_
true if the matrix was changed since the last solve
double * v_rhs_
local RHS array pointing to Vec rhs_
void set_from_input(const Input::Record in_rec) override
LinSys_PETSC(const Distribution *rows_ds)
ConstraintVec_ constraints_
Vec rhs_
PETSc vector constructed with vx array.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
void start_allocation() override
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
void start_add_assembly() override
#define LogOut()
Macro defining 'log' record of log.
bool rhs_changed_
true if the right hand side was changed since the last solve
double get_solution_precision() override
void start_insert_assembly() override
unsigned int max_it_
maximum number of iterations of linear solver
T * makePetscPointer_(std::vector< T > &array)
#define START_TIMER(tag)
Starts a timer with specified tag.
const Distribution * rows_ds_
final distribution of rows of MH matrix
double a_tol_
absolute tolerance of linear solver
unsigned int np() const
get num of processors
double residual_norm_
local solution array pointing into Vec solution_
Vec on_vec_
Vectors for counting non-zero entries in diagonal block.
void preallocate_matrix()
double compute_residual() override
std::string params_
command-line-like options for the PETSc solver
KSPConvergedReason reason
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
double r_tol_
relative tolerance of linear solver
void finish_assembly() override
void rhs_set_values(int nrow, int *rows, double *vals) override
Abstract linear system class.
bool is_positive_definite()
#define WarningOut()
Macro defining 'warning' record of log.
void set_initial_guess_nonzero(bool set_nonzero=true)
static const Input::Type::Record & get_input_type()
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
static Input::Type::Abstract & get_input_type()
bool init_guess_nonzero
flag for starting from nonzero guess
#define DebugOut()
Macro defining 'debug' record of log.
Vec solution_
PETSc vector constructed with vb array.
double solution_precision_
unsigned int lsize(int proc) const
get local size