[petsc-users] DGMRES Right preconditioner
Desire NUENTSA WAKAM
desire.nuentsa_wakam at inria.fr
Wed Oct 5 03:18:41 CDT 2011
Yes Barry, that is what the code is doing.
Thanks
Desire
On 10/05/2011 05:09 AM, Barry Smith wrote:
> 1) It seems this is what the code is already doing?
>
> if (ksp->pc_side == PC_LEFT) {
> /* Apply the first preconditioner */
> ierr = KSP_PCApplyBAorAB (ksp,VEC_VV (it), VEC_TEMP,VEC_TEMP_MATOP);
> CHKERRQ (ierr);
> /* Then apply Deflation as a preconditioner */
> ierr=KSPDGMRESApplyDeflation (ksp, VEC_TEMP, VEC_VV (1+it));
> CHKERRQ (ierr);
> } else if (ksp->pc_side == PC_RIGHT) {
> ierr=KSPDGMRESApplyDeflation (ksp, VEC_VV (it), VEC_TEMP);
> CHKERRQ (ierr);
> ierr=KSP_PCApplyBAorAB (ksp, VEC_TEMP, VEC_VV (1+it), VEC_TEMP_MATOP);
> CHKERRQ (ierr);
>
> For example, with right preconditioning it is applying M^-1(the deflation thingy) to VEC_VV(it) putting the result in VEC_TEMP then it is applying the static preconditioner (for example ILU) M_D^-1 then applying the operator A.
>
> 2) you are certainly free to edit dgmres.c and do a different algorithm and see how it goes; if you show that it works well or better than the implemented code then people are likely to accept your suggestion. But you need to take the initiative with your ideas to try them out and demonstrate they are good ideas.
>
>
> Barry
>
> On Oct 4, 2011, at 9:47 PM, Gong Ding wrote:
>
>
>> The order of the preconditioner is the problem.
>> I would like to use my own preconditioner i.e. M_D first, then the deflation preconditioner M.
>> That is solving
>> (A*M_D^-1*M^-1) (M M_D) x = b
>>
>> M_D is a static preconditioner build from A (i.e ILU or some preconditioner based on domain decomposition). And M is the deflation preconditioner build from AM_D^-1.
>> Anyway, since M is a dynamic preconditioner built from GMRES Arnoldi process, it is easy to be constructed from AM_D^-1.
>> However build preconditioner such as ILU from AM^-1 is more difficult.
>>
>> I think it is reasonable to change the order of preconditioner in DGMRES code. Or at least gives an option to do this.
>>
>>
>>
>>> Hi Gong,
>>> You do not need to modify the source code for the right preconditioning.
>>> The additional preconditioner (say M_D^{-1}) is built from the
>>> approximate eigenvalues of M^{-1}A or AM^{-1}.
>>> Hence you are solving either M_D^{-1}M_{-1}Ax = M_D^{-1}M_{-1}b or
>>> AM^{-1}M^{-1}_Dy = b for x = M^{-1}M^{-1}_Dy
>>>
>>> Thanks
>>> Desire
>>>
>>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20111005/80f6c925/attachment.htm>
More information about the petsc-users
mailing list