[petsc-dev] Right-preconditioned GMRES
Pierre Jolivet
pierre.jolivet at enseeiht.fr
Thu Feb 9 15:52:08 CST 2017
There, attached.
Thanks.
PS: here is the patch also, just in case
diff --git a/src/ksp/ksp/examples/tutorials/ex1.c b/src/ksp/ksp/examples/tutorials/ex1.c
index 6382c09..d9dc69a 100644
--- a/src/ksp/ksp/examples/tutorials/ex1.c
+++ b/src/ksp/ksp/examples/tutorials/ex1.c
@@ -123,6 +123,11 @@ int main(int argc,char **args)
routines.
*/
ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
+ ierr = KSPSetUp(ksp);CHKERRQ(ierr);
+ char right_options[] = "-ksp_norm_type unpreconditioned -ksp_pc_side right";
+ PetscOptionsInsertString(NULL,right_options);
+ ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
+ ierr = KSPSetUp(ksp);CHKERRQ(ierr);
if (nonzeroguess) {
PetscScalar p = .5;
> On Feb 9, 2017, at 2:56 PM, Matthew Knepley <knepley at gmail.com> wrote:
>
> On Thu, Feb 9, 2017 at 7:48 AM, Pierre Jolivet <pierre.jolivet at enseeiht.fr <mailto:pierre.jolivet at enseeiht.fr>> wrote:
>> On Feb 9, 2017, at 2:31 PM, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>
>> On Thu, Feb 9, 2017 at 7:29 AM, Pierre Jolivet <pierre.jolivet at enseeiht.fr <mailto:pierre.jolivet at enseeiht.fr>> wrote:
>>
>>> On Feb 9, 2017, at 2:17 PM, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>>
>>> On Thu, Feb 9, 2017 at 7:06 AM, Pierre Jolivet <pierre.jolivet at enseeiht.fr <mailto:pierre.jolivet at enseeiht.fr>> wrote:
>>>
>>>> On Feb 9, 2017, at 1:37 PM, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>>>
>>>> On Thu, Feb 9, 2017 at 4:56 AM, Pierre Jolivet <pierre.jolivet at enseeiht.fr <mailto:pierre.jolivet at enseeiht.fr>> wrote:
>>>> Hello,
>>>> I _might_ be doing something wrong when calling KSPSetFromOptions/KSPSetUp, but still, am I supposed to get this kind of error?
>>>> [0]PETSC ERROR: No support for this operation for this object type
>>>> [0]PETSC ERROR: KSP gmres does not support PRECONDITIONED with RIGHT
>>>>
>>>> If not, I’ll try to send a MWE (it basically depends on where I set the option -ksp_pc_side right).
>>>>
>>>> I believe its telling you that it cannot display the preconditioned residual with right preconditioning, which is true.
>>>> Are you requesting both?
>>>
>>> No, I’m just using -ksp_view and -ksp_pc_side right. In the “standard” case, it obviously work.
>>> However, the error is triggered if I set -ksp_pc_side right after a first call to KSPSetOperators()/KSPSetUp() without destroying the KSP.
>>> Should I destroy my KSP before changing from left- (default) to right-preconditioning if KSPSetOperators()/KSPSetUp has already been called?
>>>
>>> Its likely you have to use
>>>
>>> http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetNormType.html <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetNormType.html>
>>>
>>> or
>>>
>>> -ksp_norm_type unpreconditioned
>>
>> If I use this after my first call to KSPSetUp, I now end up with:
>> [0]PETSC ERROR: KSP gmres does not support UNPRECONDITIONED with LEFT
>>
>> Yes, you would need to use it at the same time that you change the preconditioning side.
>
> Yes, I’m using PetscOptionsInsert, and whether I let PETSc parse "-ksp_pc_side right -ksp_norm_type unpreconditioned” or "-ksp_norm_type unpreconditioned -ksp_pc_side right” I end up with the same error afterwards.
>
> Can you make an example do it (SNES ex5)? Or can you send a small example that does it?
>
> Thanks,
>
> Matt
>
> Thanks.
>
>> Thanks,
>>
>> Matt
>>> since we would not have a chance to change it since its been set already.
>>>
>>> Thanks,
>>>
>>> Matt
>>>
>>> Thanks.
>>>
>>>> Thanks,
>>>>
>>>> Matt
>>>>
>>>> Thanks,
>>>> Pierre
>>>>
>>>>
>>>>
>>>> --
>>>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>>>> -- Norbert Wiener
>>>
>>>
>>>
>>>
>>> --
>>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>>> -- Norbert Wiener
>>
>>
>>
>>
>> --
>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>> -- Norbert Wiener
>
>
>
>
> --
> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
> -- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20170209/b165cb01/attachment.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: ex1.c
Type: application/octet-stream
Size: 6777 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20170209/b165cb01/attachment.obj>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20170209/b165cb01/attachment-0001.html>
More information about the petsc-dev
mailing list