> On May 17, 2016, at 1:14 PM, Pratapa, Phanisri P <[email protected]> wrote:
>
> Barry,
>
> What I mean is that: I want to implement multigrid preconditioning on a new
> Krylov method that is not currently in PETSc. For this, I was hoping that I
> could replace the KSPFGMRES smoother (default) with the "new solver" I have.
> But I do not have the new solver as a PETSc KSP yet.
Understood. We don't provide a KSPSHELL because we consider it so easy to
implement a new KSP in PETSc that having a KSPSHELL is unnecessary. If your new
Krylov method is "stand-alone" and not, for example, a modification of GMRES
here is how to proceed. Say your KSP is called mykrylov
Make a directory src/ksp/ksp/impls/mykrylov Copy src/ksp/ksp/impls/cg/cg.c
and src/ksp/ksp/impls/cg/cgimpl.c and src/ksp/ksp/impls/cg/makefile to that
directory changing the .c and .h file names. Then edit the three files to
reflect the new names. Edit
src/ksp/ksp/impls/makefile and add the mykrylov directory to the list of
directories (variable DIRS). Then code your new Krylov method inside the two
new files you copied over following the directions given in the file.
If your new Krylov method is an extension/modification of GMRES then it is
possible to reuse most of the GMRES implementation in PETSc but a bit more
involved. Under the src/ksp/ksp/impls/gmres directory are several variants
fgmres, pgmres, lgmres, pipefgmres, agmres, and dgmres. I would recommend
picking the one that most closely resembles your new method and copying it as
above and modifying it to match your algorithm.
Barry
>
> Thank you,
>
> Regards,
>
> Pradeep
> ________________________________________
> From: Barry Smith <[email protected]>
> Sent: Tuesday, May 17, 2016 2:08:18 PM
> To: Pratapa, Phanisri P
> Cc: [email protected]
> Subject: Re: [petsc-users] User defined solver for PCMG
>
>> On May 17, 2016, at 12:21 PM, Pratapa, Phanisri P <[email protected]>
>> wrote:
>>
>> Hi,
>>
>> I am trying to find out if one can use a user defined linear solver function
>> in PCMG (instead of the default GMRES). According to the petsc manual, I can
>> change the solver/smoother through the KSP context and the available
>> solvers, but I am interested in using my own function (solver).
>
> What do you mean here by "solver"? Do you want to provide a new Krylov
> method that is not currently in PETSc or a new preconditioner that is
> specific to your problem and cannot be written as a composition of
> preconditioners and Krylov methods already in PETSc?
>
> An example of your own custom preconditioner could be an SOR iteration that
> you hand code based on the stencil and doesn't use a stored matrix. In that
> case you would access the PC on the level and use PCSetType(subpc,PCSHELL)
> and PCShellSetApply(subpc, your custom function).
>
> Barry
>
>>
>> Thank you,
>>
>> Regards,
>>
>> Pradeep
>