[petsc-users] How to set KSP Solver Options at Runtime?

Abhyankar, Shrirang G. abhyshr at mcs.anl.gov
Wed Jul 9 23:04:30 CDT 2014


Shuangshuang,
   You don't need 

PCRegister("ourjacobi", PCCreate_Jacobi);
  KSPGetPC(ksp, &pc);
  PCSetType(pc, "ourjacobi");


in your code if you want to use PETSc's native Jacobi preconditioner. All
you have to do is set it at run-time  -pc_type jacobi.
In PETSc, XXXSetFromOptions routine checks if a particular solver is
specified via a run-time option. If it is then that solver is set otherwise
it uses the default. KSPSetFromOptions sets the linear solver and the
preconditioner (runtime options -ksp_type and -pc_type),
SNESSetFromOptions 
the nonlinear solver (runtime option -snes_type), and TSSetFromOptions the
time-stepping method (runtime option -ts_type).

There are other bunch of options for tuning preconditioners. Run your code
with -help option to see the list of options available.

Shri


From:  "Jin, Shuangshuang" <Shuangshuang.Jin at pnnl.gov>
Date:  Wed, 9 Jul 2014 23:37:19 +0000
To:  "petsc-users at mcs.anl.gov" <petsc-users at mcs.anl.gov>
Subject:  [petsc-users] How to set KSP Solver Options at Runtime?


>Hi, I¹m using the KSP solver as following to solve a linear system Ax=b,
>where A is a 16767x16767 square matrix, b is a 16767 length vector.
> 
>EXTERN_C_BEGIN
>extern PetscErrorCode PCCreate_Jacobi(PC);
>EXTERN_C_END
> 
>static PetscErrorCode solvingAxb(Mat A, Vec b, PetscInt nbus, Vec &x,
>const int me)
>{
>  PetscErrorCode ierr;
> 
>  KSP ksp; // linear solver context
>  Vec u; // exact solution
>  PetscInt its;
>  PetscReal norm; // norm of solution error
>  PetscLogDouble t1, t2;
>  PetscViewer viewer;
>  PC pc; // preconditioner context
>  PetscInt Istart, Iend;
> 
>  ierr = VecDuplicate(b, &u); CHKERRQ(ierr);
>  ierr = VecDuplicate(b, &x); CHKERRQ(ierr);
> 
>  ierr = KSPCreate(PETSC_COMM_WORLD, &ksp); CHKERRQ(ierr);
> 
>  ierr = KSPSetOperators(ksp, A, A, DIFFERENT_NONZERO_PATTERN);
>CHKERRQ(ierr);
> 
>  PCRegister("ourjacobi", PCCreate_Jacobi);
> 
>  KSPGetPC(ksp, &pc);
>  PCSetType(pc, "ourjacobi");
> 
>  ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
> 
>  ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
> 
>  ierr = MatMult(A,x,u); CHKERRQ(ierr);
>  ierr = VecAXPY(u, -1.0, b); CHKERRQ(ierr);
>  ierr = VecNorm(u, NORM_2, &norm); CHKERRQ(ierr);
>  ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
> 
>  ierr = VecDestroy(&b); CHKERRQ(ierr);
>  ierr = VecDestroy(&u); CHKERRQ(ierr);
> 
>  ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
> 
>  PetscFunctionReturn(0);
>}
> 
>I would like to know how to set the solver options at runtime to make it
>run faster, such as ksp_type, pc_type, and etc? It takes very long time
>to solve the system if I use no options.
> 
>Thanks,
>Shuangshuang
> 



More information about the petsc-users mailing list