42 namespace it = Input::Type;
52 init_guess_nonzero(false),
60 ierr = VecZeroEntries(
rhs_ ); CHKERRV( ierr );
70 :
LinSys(other), params_(other.params_), v_rhs_(NULL), solution_precision_(other.solution_precision_)
100 ASSERT(
false,
"Can not set values. Matrix is not preallocated.\n");
118 ASSERT(
false,
"Can not set values. Matrix is not preallocated.\n");
131 ierr = MatSetValues(
matrix_,nrow,rows,ncol,cols,vals,(InsertMode)
status_); CHKERRV( ierr );
136 default:
DBGMSG(
"LS SetValues with non allowed insert mode.\n");
149 ierr = VecSetValues(
rhs_,nrow,rows,vals,(InsertMode)
status_); CHKERRV( ierr );
153 default:
ASSERT(
false,
"LinSys's status disallow set values.\n");
165 for (i=0; i<nrow; i++) {
167 for(j=0; j<ncol; j++) {
170 VecSetValue(
on_vec_,row,1.0,ADD_VALUES);
172 VecSetValue(
off_vec_,row,1.0,ADD_VALUES);
182 PetscInt *on_nz, *off_nz;
183 PetscScalar *on_array, *off_array;
195 VecGetArray(
on_vec_, &on_array );
196 VecGetArray(
off_vec_, &off_array );
199 on_nz[i] =
static_cast<PetscInt
>( on_array[i]+0.1 );
200 off_nz[i] =
static_cast<PetscInt
>( off_array[i]+0.1 );
203 VecRestoreArray(
on_vec_,&on_array);
204 VecRestoreArray(
off_vec_,&off_array);
211 ierr = MatDestroy(&
matrix_); CHKERRV( ierr );
214 0, on_nz, 0, off_nz, &
matrix_); CHKERRV( ierr );
217 MatSetOption(
matrix_, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
225 MatAssemblyType assemblyType = MAT_FINAL_ASSEMBLY;
234 xprintf(
Warn,
"Finalizing linear system without setting values.\n");
237 ierr = MatAssemblyBegin(
matrix_, assembly_type); CHKERRV( ierr );
238 ierr = VecAssemblyBegin(
rhs_); CHKERRV( ierr );
239 ierr = MatAssemblyEnd(
matrix_, assembly_type); CHKERRV( ierr );
240 ierr = VecAssemblyEnd(
rhs_); CHKERRV( ierr );
242 if (assembly_type == MAT_FINAL_ASSEMBLY)
status_ =
DONE;
253 ASSERT (
status_ ==
DONE,
"System matrix and right-hand side are not assembled when applying constraints." );
256 PetscInt numConstraints =
static_cast<PetscInt
>(
constraints_.size() );
259 const PetscScalar diagScalar =
static_cast<PetscScalar
>( scalar );
265 ConstraintVec_::const_iterator cIter =
constraints_.begin( );
266 ConstraintVec_::const_iterator cEnd =
constraints_.end( );
268 for ( ; cIter != cEnd; ++cIter ) {
269 globalDofs.push_back( static_cast<PetscInt>( cIter -> first ) );
270 values.push_back( static_cast<PetscScalar>( cIter -> second ) * diagScalar );
278 ierr = MatZeroRows(
matrix_, numConstraints, globalDofPtr, diagScalar, PETSC_NULL, PETSC_NULL ); CHKERRV( ierr );
282 if ( numConstraints ) {
283 ierr = VecSetValues(
rhs_, numConstraints, globalDofPtr, valuePtr, INSERT_VALUES ); CHKERRV( ierr );
288 ierr = VecAssemblyBegin(
rhs_ ); CHKERRV( ierr );
289 ierr = VecAssemblyEnd(
rhs_ ); CHKERRV( ierr );
302 KSPConvergedReason reason;
304 const char *petsc_dflt_opt;
313 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";
315 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";
322 petsc_dflt_opt=
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
324 petsc_dflt_opt=
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 5 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
329 PetscOptionsInsertString(
params_.c_str());
331 MatSetOption(
matrix_, MAT_USE_INODES, PETSC_FALSE );
333 KSPCreate(
comm_, &system );
334 KSPSetOperators(system,
matrix_,
matrix_, DIFFERENT_NONZERO_PATTERN);
337 KSPSetTolerances(system,
r_tol_,
a_tol_, PETSC_DEFAULT,PETSC_DEFAULT);
338 KSPSetFromOptions(system);
345 KSPGetType(system, &type);
346 if (strcmp(type, KSPPREONLY) != 0)
347 KSPSetInitialGuessNonzero(system, PETSC_TRUE);
354 KSPGetConvergedReason(system,&reason);
355 KSPGetIterationNumber(system,&nits);
361 xprintf(
MsgLog,
"convergence reason %d, number of iterations is %d\n", reason, nits);
371 return static_cast<int>(reason);
377 std::string matFileName =
"flow123d_matrix.m";
378 std::string rhsFileName =
"flow123d_rhs.m";
379 std::string solFileName =
"flow123d_sol.m";
381 PetscViewer myViewer;
383 PetscViewerASCIIOpen(
comm_,matFileName.c_str(),&myViewer);
384 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
386 PetscViewerDestroy(&myViewer);
388 PetscViewerASCIIOpen(
comm_,rhsFileName.c_str(),&myViewer);
389 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
390 VecView(
rhs_, myViewer );
391 PetscViewerDestroy(&myViewer);
394 PetscViewerASCIIOpen(
comm_,solFileName.c_str(),&myViewer);
395 PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
397 PetscViewerDestroy(&myViewer);
405 if (
matrix_ != NULL) { ierr = MatDestroy(&
matrix_); CHKERRV( ierr ); }
406 ierr = VecDestroy(&
rhs_); CHKERRV( ierr );
unsigned int get_proc(unsigned int idx) const
get processor of the given index
SetValuesMode status_
Set value status of the linear system.
void start_insert_assembly()
void apply_constrains(double scalar=1.)
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_
LinSys_PETSC(const Distribution *rows_ds)
ConstraintVec_ constraints_
Vec rhs_
PETSc vector constructed with vx array.
static Input::Type::Record input_type
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
void set_from_input(const Input::Record in_rec)
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
void rhs_set_values(int nrow, int *rows, double *vals)
bool rhs_changed_
true if the right hand side was changed since the last solve
static Input::Type::AbstractRecord input_type
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
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()
std::string params_
command-line-like options for the PETSc solver
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)
void start_add_assembly()
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
double get_solution_precision()
Abstract linear system class.
bool is_positive_definite()
void set_initial_guess_nonzero(bool set_nonzero=true)
Mat matrix_
Petsc matrix of the problem.
bool init_guess_nonzero
flag for starting from nonzero guess
Vec solution_
PETSc vector constructed with vb array.
double solution_precision_
unsigned int lsize(int proc) const
get local size