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);
244 MatSetOption(
matrix_, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
252 MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
261 WarningOut() <<
"Finalizing linear system without setting values.\n";
264 ierr = MatAssemblyBegin(
matrix_, assembly_type); CHKERRV( ierr );
265 ierr = VecAssemblyBegin(
rhs_); CHKERRV( ierr );
266 ierr = MatAssemblyEnd(
matrix_, assembly_type); CHKERRV( ierr );
267 ierr = VecAssemblyEnd(
rhs_); CHKERRV( ierr );
269 if (assembly_type == MAT_FINAL_ASSEMBLY)
status_ =
DONE;
285 OLD_ASSERT (
status_ ==
DONE,
"System matrix and right-hand side are not assembled when applying constraints." );
288 PetscInt numConstraints =
static_cast<PetscInt
>(
constraints_.size() );
291 const PetscScalar diagScalar =
static_cast<PetscScalar
>( scalar );
297 ConstraintVec_::const_iterator cIter =
constraints_.begin( );
298 ConstraintVec_::const_iterator cEnd =
constraints_.end( );
300 for ( ; cIter != cEnd; ++cIter ) {
301 globalDofs.push_back(
static_cast<PetscInt
>( cIter -> first ) );
302 values.push_back(
static_cast<PetscScalar
>( cIter -> second ) * diagScalar );
310 ierr = MatZeroRows(
matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
314 if ( numConstraints ) {
315 ierr = VecSetValues(
rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
320 ierr = VecAssemblyBegin(
rhs_ ); CHKERRV( ierr );
321 ierr = VecAssemblyEnd(
rhs_ ); CHKERRV( ierr );
334 const char *petsc_dflt_opt;
349 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";
352 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";
358 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";
361 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";
365 LogOut().fmt(
"inserting petsc options: {}\n",
params_.c_str());
369 PetscOptionsInsertString(NULL,
params_.c_str());
371 MatSetOption(
matrix_, MAT_USE_INODES, PETSC_FALSE );
380 KSPSetFromOptions(
system);
387 KSPGetType(
system, &type);
388 if (strcmp(type, KSPPREONLY) != 0)
389 KSPSetInitialGuessNonzero(
system, PETSC_TRUE);
397 KSPGetIterationNumber(
system,&nits);
403 LogOut().fmt(
"convergence reason {}, number of iterations is {}\n",
reason, nits);
419 std::string matFileName =
"flow123d_matrix.m";
420 std::string rhsFileName =
"flow123d_rhs.m";
421 std::string solFileName =
"flow123d_sol.m";
423 PetscViewer myViewer;
425 PetscViewerASCIIOpen(
comm_,matFileName.c_str(),&myViewer);
426 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
428 PetscViewerDestroy(&myViewer);
430 PetscViewerASCIIOpen(
comm_,rhsFileName.c_str(),&myViewer);
431 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
432 VecView(
rhs_, myViewer );
433 PetscViewerDestroy(&myViewer);
436 PetscViewerASCIIOpen(
comm_,solFileName.c_str(),&myViewer);
437 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
439 PetscViewerDestroy(&myViewer);
461 std::string user_params = in_rec.
val<
string>(
"options");
462 if (user_params !=
"")
params_ = user_params;
476 double residual_norm;
477 VecNorm(
residual_, NORM_2, &residual_norm);
478 return residual_norm;