[petsc-users] multigrid implementation problem
Nakib Haider Protik
nprot048 at uottawa.ca
Sun Jun 10 20:19:37 CDT 2012
With the suggested changes, the code doesn't run. Is there any example
code of KSP solvers using multigrid preconditioning? All the examples of
multigrid I have come across involve dmmg.
Thanks.
> On Sun, Jun 10, 2012 at 7:54 PM, Nakib Haider Protik
> <nprot048 at uottawa.ca>wrote:
>
>> Okay, so here is the code that works. This gives the correct result, but
>> I
>> want to make sure that it is not giving the correct result accidentally.
>>
>
> 1. I don't know how to emphasize this enough: you need to learn C to write
> C code.
>
> 2. Please run in Valgrind, it will show you several memory leaks here.
>
> 3. Please run with -ksp_view so you can see what gets used. If you did
> this, it would be obvious that multigrid was not being used.
>
>
>> Thank you very much for your time.
>>
>> KSPCreate(PETSC_COMM_WORLD, &ksp);
>> KSPCreate(PETSC_COMM_WORLD, &cksp);
>> KSPCreate(PETSC_COMM_WORLD, &uksp);
>> KSPCreate(PETSC_COMM_WORLD, &dksp);
>>
>
> ^^ You only need to create ksp here, delete the other lines
>
>
>> KSPSetType(ksp, KSPGMRES);
>> KSPSetType(cksp, KSPGMRES);
>> KSPSetType(uksp, KSPGMRES);
>> KSPSetType(dksp, KSPGMRES);
>>
>
> ^^ You don't need any of these
>
>
>>
>> KSPGetPC(ksp, &pc);
>> KSPGetPC(cksp, &pc);
>> KSPGetPC(uksp, &pc);
>> KSPGetPC(dksp, &pc);
>>
>
> ^^ delete all but the first of these
>
>
>>
>> KSPSetFromOptions(ksp);
>>
>> PCSetType(pc, PCMG);
>> PCMGSetLevels(pc, 2, PETSC_NULL);
>> PCMGSetType(pc, PC_MG_MULTIPLICATIVE);
>> PCMGSetCycleType(pc, PC_MG_CYCLE_V);
>> MatDuplicate(A, MAT_COPY_VALUES, &P);
>> PCMGSetCyclesOnLevel(pc, 0, 1);
>> PCMGSetCyclesOnLevel(pc, 1, 1);
>>
>> PCMGGetCoarseSolve(pc, &cksp);
>> PCMGGetSmootherDown(pc, 0, &dksp);
>> PCMGGetSmootherUp(pc, 1, &uksp);
>>
>
> ^^^ delete the last three lines here, they aren't doing anything
>
>
>>
>> PCMGSetInterpolation(pc, 1, P);
>> PCMGSetRestriction(pc, 1, P);
>> PCMGSetResidual(pc, 0, PCMGDefaultResidual, P);
>> PCMGSetResidual(pc, 1, PCMGDefaultResidual, P);
>>
>
> There is no way that you want to be using P for interpolation/restriction
> and as the residual on levels and as the preconditioning matrix on the
> line
> below. The matrices are not even the right size. The reason the above runs
> is because nothing you have done with "pc" is being used: it refers to a
> "dksp" which you are leaking without using.
>
>
>>
>> KSPSetOperators(ksp, A, P, SAME_NONZERO_PATTERN);
>> KSPSolve(ksp, breal, xreal_harm);
>> KSPSolve(ksp, bimag, ximag_harm);
>> //////////////////////////////////////////////////////////
>>
>> > On Sun, Jun 10, 2012 at 7:32 PM, Nakib Haider Protik
>> > <nprot048 at uottawa.ca>wrote:
>> >
>> >> Okay. But if I understood correctly, doing the following would also
>> be
>> >> wrong since I will have made three extra outer KSP objects. I am not
>> >> sure
>> >> I understand what these functions want.
>> >>
>> >> KSP cksp, uksp, dksp;
>> >> PCMGGetCoarseSolve(pc, &cksp);/**/
>> >> PCMGGetSmootherDown(pc, 0, &dksp);/**/
>> >> PCMGGetSmootherUp(pc, 1, &uksp);/**/
>> >>
>> >
>> > This is fine, all of these are getting references to inner objects
>> without
>> > overwriting the outer objects. Of course you don't need to get these
>> > references if you aren't going to modify them.
>> >
>> > I recommend that you experiment with algorithms using the command line
>> and
>> > check that it does what you expect using -ksp_view.
>> >
>>
>>
>> --
>> Nakib :)
>>
>
--
Nakib :)
More information about the petsc-users
mailing list