[petsc-users] PCMG with PCREDISTRIBUTE

Matthew Knepley knepley at gmail.com
Fri Jun 30 12:28:21 CDT 2023


On Fri, Jun 30, 2023 at 12:08 PM Carl-Johan Thore <carl-johan.thore at liu.se>
wrote:

> “Possibly, but if you are doing FD, then there is built-in topology in
> DMDA that is not present in Plex, so
>
> finding the neighbors in the right order is harder (possible, but harder,
> we address this in some new work that is not yet merged). There is also
> structured adaptive support with DMForest, but this also does not
> preserve the stencil.
>
>>
>
>
> I’m using an FEM which doesn’t utilize such neighbor info, so perhaps Plex
> or DMForest could be easier to use then
>

Oh yes. Then I would recommend at least looking at Plex and Forest. FEM is
much better supported than in DMDA. By the end of the year, I should have
everything in place to seamlessly support libCEED for assembly as well.

    Matt


> “The efficiency of active set VI solvers in PETSc demonstrates to me that
> solving reduced systems can be done efficiently with geometric multigrid
> using a structured grid so I would not suggest giving up on what you
> started.
>
>
>
>     You can do it in two steps
>
>
>
> 1) Use PCREDISTRIBUTE but hack the code in redistribute.c to not move dof
> between MPI ranks, just have it remove the locked rows/columns (to start
> just run on one MPI rank since then nothing is moved) Then in your  code
> you just need to pull out the appropriate rows and columns of the
> interpolation that correspond to the dof you have kept and pass this
> smaller interpolation to the inner KSP PCMG. This is straightforward and
> like what is in DMSetVI.  The MG convergence should be just as good as on
> the full system.
>
>
>
> 2) the only problem with 1 is it is likely to be poorly load balanced (but
> you can make some runs to see how imbalanced it is, that will depend
> exactly on what parts are locked and what MPI processes they are on).  So
> if it is poorly balanced then you would need to get out of redistribute.c a
> mapping for each kept dof to what MPI rank it is moved to and use that to
> move the entries in the reduced interpolation you have created.
>
>
>
>   If you do succeed it would actually be useful code that we could add to
> PCREDISTRIBUTE for more general use by others.
>
>
>
>   Barry
>
>>
>
>
> Thanks, that seems doable, if not super easy. I might try that
>
>
>
> Kind regards
>
> /Carl-Johan
>
>
>
> *From:* Barry Smith <bsmith at petsc.dev>
> *Sent:* Friday, June 30, 2023 5:21 PM
> *To:* Matthew Knepley <knepley at gmail.com>
> *Cc:* Carl-Johan Thore <carl-johan.thore at liu.se>; petsc-users at mcs.anl.gov
> *Subject:* Re: [petsc-users] PCMG with PCREDISTRIBUTE
>
>
>
>
>
>
>
> On Jun 30, 2023, at 10:22 AM, Matthew Knepley <knepley at gmail.com> wrote:
>
>
>
> On Fri, Jun 30, 2023 at 10:16 AM Carl-Johan Thore via petsc-users <
> petsc-users at mcs.anl.gov> wrote:
>
> Thanks for the quick reply and the suggestions!
>
>
>
> “ … you should first check that the PCMG works quite well “
>
>
>
> Yes, the PCMG works very well for the full system.
>
>
>
> “I am guessing that your code is slightly different than ex42.c because
> you take the interpolation matrix provided by the DM
>
> and give it to the inner KSP PCMG?. So you solve problem 2 but not problem
> 1.”
>
>
>
> Yes, it’s slightly different so problem 2 should be solved.
>
>
>
> It looked somewhat complicated to get PCMG to work with redistribute, so
> I’ll try with PCGAMG first
>
> (it ran immediately with redistribute, but was slower than PCMG on my,
> very small, test problem. I’ll try
>
> to tune the settings).
>
>
>
> A related question: I’m here using a DMDA for a structured grid but I’m
> locking so many DOFs that for many of the elements
>
> all DOFs are locked. In such a case could it make sense to switch/convert
> the DMDA to a DMPlex containing only those
>
> elements that actually have DOFs?
>
>
>
> Possibly, but if you are doing FD, then there is built-in topology in DMDA
> that is not present in Plex, so
>
> finding the neighbors in the right order is harder (possible, but harder,
> we address this in some new work that is not yet merged). There is also
> structured adaptive support with DMForest, but this also does not
> preserve the stencil.
>
>
>
>    The efficiency of active set VI solvers in PETSc demonstrates to me
> that solving reduced systems can be done efficiently with geometric
> multigrid using a structured grid so I would not suggest giving up on what
> you started.
>
>
>
>     You can do it in two steps
>
>
>
> 1) Use PCREDISTRIBUTE but hack the code in redistribute.c to not move dof
> between MPI ranks, just have it remove the locked rows/columns (to start
> just run on one MPI rank since then nothing is moved) Then in your  code
> you just need to pull out the appropriate rows and columns of the
> interpolation that correspond to the dof you have kept and pass this
> smaller interpolation to the inner KSP PCMG. This is straightforward and
> like what is in DMSetVI.  The MG convergence should be just as good as on
> the full system.
>
>
>
> 2) the only problem with 1 is it is likely to be poorly load balanced (but
> you can make some runs to see how imbalanced it is, that will depend
> exactly on what parts are locked and what MPI processes they are on).  So
> if it is poorly balanced then you would need to get out of redistribute.c a
> mapping for each kept dof to what MPI rank it is moved to and use that to
> move the entries in the reduced interpolation you have created.
>
>
>
>   If you do succeed it would actually be useful code that we could add to
> PCREDISTRIBUTE for more general use by others.
>
>
>
>   Barry
>
>
>
>
>
>
>
>
>
>   Thanks,
>
>
>
>     Matt
>
>
>
> *From:* Barry Smith <bsmith at petsc.dev>
> *Sent:* Friday, June 30, 2023 3:57 PM
> *To:* Carl-Johan Thore <carl-johan.thore at liu.se>
> *Cc:* petsc-users at mcs.anl.gov
> *Subject:* Re: [petsc-users] PCMG with PCREDISTRIBUTE
>
>
>
>
>
>    Oh, I forgot to mention you should first check that the PCMG works
> quite well for the full system (without the PCREDISTRIBUTE); the convergence
>
> on the redistributed system (assuming you did all the work to get PCMG to
> work for you) should be very similar to (but not measurably better) than
> the convergence on the full system.
>
>
>
>
>
>
>
> On Jun 30, 2023, at 9:17 AM, Barry Smith <bsmith at petsc.dev> wrote:
>
>
>
>
>
>    ex42.c provides directly the interpolation/restriction needed to move
> between levels in the loop
>
>
>
>  for (k = 1; k < nlevels; k++) {
>
>     PetscCall(DMCreateInterpolation(da_list[k - 1], da_list[k], &R, NULL));
>
>     PetscCall(PCMGSetInterpolation(pc, k, R));
>
>     PetscCall(MatDestroy(&R));
>
>   }
>
>
>
> The more standard alternative to this is to call KSPSetDM() and have the
> PCMG setup use the DM
>
> to construct the interpolations (I don't know why ex42.c does this
> construction itself instead of having the KSPSetDM() process handle it but
> that doesn't matter). The end result is the same in both cases.
>
>
>
> Since PCREDISTRIBUTE  builds its own  new matrix (by using only certain
> rows and columns of the original matrix) the original interpolation
>
> cannot be used for two reasons
>
>
>
> 1) (since it is for the full system) It is for the wrong problem.
>
>
>
> 2) In addition, if you ran with ex42.c the inner KSP does not have access
> to the interpolation that was constructed so you could not get PCMG to to
> work as indicated below.
>
>
>
> I am guessing that your code is slightly different than ex42.c because you
> take the interpolation matrix provided by the DM
>
> and give it to the inner KSP PCMG?. So you solve problem 2 but not problem
> 1.
>
>
>
> So the short answer is that there is no "canned" way to use the PCMG
> process trivially with PCDISTRIBUTE.
>
>
>
> To do what you want requires two additional steps
>
>
>
> 1) after you construct the full interpolation matrix  (by using the DM)
> you need to remove the rows associated with the dof that have been removed
> by the "locked" variables (and the columns that are associated with coarse
> grid points that live on the removed points) so that the interpolation is
> the correct "size" for the smaller problem
>
>
>
> 2) since PCREDISTRIBUTE actually moves dof of freedom between MPI
> processes for load balancing after it has removed the locked variables you
> would need to do the exact same movement for the rows of the interpolation
> matrix that you have constructed (after you have removed the "locked" rows
> of the interpolation.
>
>
>
> Lots of bookkeeping to acheive 1 and 2 but conceptually simple.
>
>
>
> As an experiment you can try using PCGAMG on the redistributed matrix
> -redistribute_pc_type gamg to use algebraic multigrid just to see the time
> and convergence rates. Since GAMG creates its own interpolation based on
> the matrix and it will be built on the smaller redistributed matrix there
> will no issue with the wrong "sized" interpolation. Of course you have the
> overhead of algebraic multigrid and cannot take advantage of geometric
> multigrid.  The GAMG approach may be satisfactory to your needs.
>
>
>
> If you are game for looking more closely at using redistribute with
> geometric multigrid and PETSc (which will require digging into PETSc source
> code and using internal information in the PETSc source code) you can start
> by looking at how we solve variational problems with SNES using reduced
> space active set methods. SNESVINEWTONRSLS /src/snes/impls/vi/rs/virs.c
> This code solves problem 1 see() it builds the entire interpolation and
> then pulls out the required non-locked part. Reduced space active set
> methods essentially lock the constrained dof and solve a smaller system
> without those dof at each iteration.
>
>
>
> But it does not solve problem 2. Moving the rows of the "smaller"
> interpolation to the correct MPI process based on where PCREDISTRIBUTE
> moved rows. To do this would requring looking at the PCREDISTRUBUTE code
> and extracting the information of where each row was moving and performing
> the process for the interpolation matrix.
>
> src/ksp/pc/impls/redistribute/redistribute.c
>
>
>
>   Barry
>
>
>
>
>
>
>
>
>
>
>
>
>
>
>
>
>
>
>
> On Jun 30, 2023, at 8:21 AM, Carl-Johan Thore via petsc-users <
> petsc-users at mcs.anl.gov> wrote:
>
>
>
> Hi,
>
>
>
> I'm trying to run an iterative solver (FGMRES for example) with PCMG as
> preconditioner. The setup of PCMG
>
> is done roughly as in ex42 of the PETSc-tutorials (
> https://petsc.org/main/src/ksp/ksp/tutorials/ex42.c.html).
>
> Since I have many locked  degrees-of-freedom I would like to use
> PCREDISTRIBUTE. However, this
>
> results in (30039 is the number of DOFs after redistribute and 55539 the
> number before):
>
>
>
> [0]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
>
> [0]PETSC ERROR: Nonconforming object sizes
>
> [0]PETSC ERROR: Matrix dimensions of A and P are incompatible for
> MatProductType PtAP: A 30039x30039, P 55539x7803
>
> [0]PETSC ERROR: See https://petsc.org/release/faq/ for trouble shooting.
>
> [0]PETSC ERROR: Petsc Development GIT revision: v3.19.0-238-g512d1ae6db4
> GIT Date: 2023-04-24 16:37:00 +0200
>
> [0]PETSC ERROR: topopt on a arch-linux-c-opt Fri Jun 30 13:28:41 2023
>
> [0]PETSC ERROR: Configure options COPTFLAGS="-O3 -march=native"
> CXXOPTFLAGS="-O3 -march=native" FOPTFLAGS="-O3 -march=native"
> CUDAOPTFLAGS=-O3 --with-cuda --with-cusp --with-debugging=0
> --download-scalapack --download-hdf5 --download-zlib --download-mumps
> --download-parmetis --download-metis --download-ptscotch --download-hypre
> --download-spai
>
> [0]PETSC ERROR: #1 MatProductSetFromOptions_Private() at
> /mnt/c/mathware/petsc/src/mat/interface/matproduct.c:420
>
> [0]PETSC ERROR: #2 MatProductSetFromOptions() at
> /mnt/c/mathware/petsc/src/mat/interface/matproduct.c:541
>
> [0]PETSC ERROR: #3 MatPtAP() at
> /mnt/c/mathware/petsc/src/mat/interface/matrix.c:9868
>
> [0]PETSC ERROR: #4 MatGalerkin() at
> /mnt/c/mathware/petsc/src/mat/interface/matrix.c:10899
>
> [0]PETSC ERROR: #5 PCSetUp_MG() at
> /mnt/c/mathware/petsc/src/ksp/pc/impls/mg/mg.c:1029
>
> [0]PETSC ERROR: #6 PCSetUp() at
> /mnt/c/mathware/petsc/src/ksp/pc/interface/precon.c:994
>
> [0]PETSC ERROR: #7 KSPSetUp() at
> /mnt/c/mathware/petsc/src/ksp/ksp/interface/itfunc.c:406
>
> [0]PETSC ERROR: #8 PCSetUp_Redistribute() at
> /mnt/c/mathware/petsc/src/ksp/pc/impls/redistribute/redistribute.c:327
>
> [0]PETSC ERROR: #9 PCSetUp() at
> /mnt/c/mathware/petsc/src/ksp/pc/interface/precon.c:994
>
> [0]PETSC ERROR: #10 KSPSetUp() at
> /mnt/c/mathware/petsc/src/ksp/ksp/interface/itfunc.c:406
>
> [0]PETSC ERROR: #11 KSPSolve_Private() at
> /mnt/c/mathware/petsc/src/ksp/ksp/interface/itfunc.c:824
>
> [0]PETSC ERROR: #12 KSPSolve() at
> /mnt/c/mathware/petsc/src/ksp/ksp/interface/itfunc.c:1070
>
>
>
> It’s clear what happens I think, and it kind of make since not all levels
> are redistributed as they should (?).
>
> Is it possible to use PCMG with PCREDISTRIBUTE in an easy way?
>
>
>
> Kind regards,
>
> Carl-Johan
>
>
>
>
>
>
>
>
> --
>
> 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
>
>
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>
>
>


-- 
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

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230630/ef28ac01/attachment-0001.html>


More information about the petsc-users mailing list