[petsc-users] Geometric multigrid with a cell-centered discretization

Buesing, Henrik hbuesing at eonerc.rwth-aachen.de
Wed Jan 10 04:22:27 CST 2018


> >> you need to not attach any
> >> coordinates to your DMDA and it should work. See for example
> >> src/ksp/ksp/examples/tutorials/ex32.c
> >
> > You mean, I should omit the call to DMDASetUniformCoordinates, right?
> 
>  I guess.
> 
> > This means, that my Paraview files will have no domain size info
> anymore.
> 
>   Presumably. But you may be able to work around the issue some way, for
> example setting the coordinates after the solvers are complete.

[Buesing, Henrik] Good idea! I moved the DMDASetUniformCoordinates to after the SNESSolve call. This works! 

With this, I stepped over a weird bug with my own convergence test function. It appears only in optimized mode and not in debug mode. Anyhow, I will send another message for this. 

Thank you for your help!
Henrik



> 
>  Barry
> 
> >
> > I will try it tomorrow and come back to you with the results.
> >
> > Thank you!
> > Henrik
> >
> >>
> >>  Barry
> >>
> >> Sorry about this, someone needs to fix the injection for Q0 case.
> >>
> >>
> >>> On Jan 9, 2018, at 9:12 AM, Buesing, Henrik <hbuesing at eonerc.rwth-
> >> aachen.de> wrote:
> >>>
> >>> Sorry (again), for posting only part of the error message...
> >>>
> >>> So, here it comes. I will first [1] show the code snippets I am
> >>> using, then
> >> second [2] I will show the PETSc options, and finally [3] the error
> >> message I get.
> >>> The grid is 1536 x 1 x 256 in x,y,z-direction. Thus, we see the
> >>> first level of
> >> coarsening in x-direction here.
> >>>
> >>> Thank you!
> >>> Henrik
> >>>
> >>> [1]
> >>>       CALL
> >> DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_GHOSTED,DM_BO
> >> UNDARY_GHOSTED,DM_BOUNDARY_GHOSTED,&
> >>>                       DMDA_STENCIL_STAR, I0g,J0g,K0g, PETSC_DECIDE,
> >> PETSC_DECIDE,PETSC_DECIDE, petsc_two, petsc_one, &
> >>>
> >> PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,
> >> petsc_da, petsc_ierr)
> >>>       CALL DMSetUp(petsc_da,petsc_ierr) ...
> >>>       CALL DMSetFromOptions(petsc_da, petsc_ierr)
> >>>       CALL DMDASetInterpolationType(petsc_da, DMDA_Q0, petsc_ierr)
> >>>       CALL SNESSetDM(petsc_snes, petsc_da, petsc_ierr)
> >>>       CALL SNESSetFromOptions(petsc_snes,petsc_ierr)
> >>>
> >>> [2] -snes_monitor -snes_linesearch_monitor -snes_type newtonls
> >>> -snes_linesearch_type basic -ksp_type fgmres -ksp_max_it 1000
> >>> -snes_stol 1e-06 -snes_rtol 1e-06 -snes_atol 1e-6 -ksp_gmres_restart
> >>> 100 -ksp_monitor -pc_type mg -da_refine_y 1 -ksp_converged_reason
> >>> -snes_converged_reason -ksp_monitor -ksp_view -pc_mg_levels 2
> >>> -options_left
> >>>
> >>> [3]
> >>> [0]PETSC ERROR: --------------------- Error Message
> >>> --------------------------------------------------------------
> >>> [0]PETSC ERROR: Arguments are incompatible [0]PETSC ERROR: Ratio
> >>> between levels: (mx - 1)/(Mx - 1) must be integer: mx 1536 Mx 768
> >>> [0]PETSC ERROR: See
> >> http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble
> >> shooting.
> >>> [0]PETSC ERROR: Petsc Development GIT revision: v3.8.2-48-g851ec02
> >>> GIT Date: 2017-12-05 09:52:17 -0600 [0]PETSC ERROR:
> >>>
> >> /rwthfs/rz/cluster/work/hb111949/EoCoE/aixcelerate17/gge_gnu/bench_ru
> >> n
> >>> /000023/000000_copy_sub_launch/work/shem_fw64gnu_const.x on a
> >>> gnu_openmpi named linuxbmc0008.rz.RWTH-Aachen.DE by hb111949 Tue
> >> Jan
> >>> 9 16:04:39 2018 [0]PETSC ERROR: Configure options
> >>> --download-fblaslapack --with-cc=mpicc -with-fc=mpif90
> >>> --with-cxx=mpicxx --download-hypre --download-superlu_dist
> >>> --download-suitesparse --download-scalapack --download-blacs
> >>> --download-hdf5 --download-parmetis --download-metis
> >>> --with-debugging=0 --download-mumps [0]PETSC ERROR: #1
> >>> DMCreateInjection_DA_3D() line 1195 in
> >>>
> >> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/dm/impls/da/dainterp.
> >> c
> >>> [0]PETSC ERROR: #2 DMCreateInjection_DA() line 1283 in
> >>>
> >> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/dm/impls/da/dainterp.
> >> c
> >>> [0]PETSC ERROR: #3 DMCreateInjection() line 1082 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/dm/interface/dm.c
> >>> [0]PETSC ERROR: #4 DMCoarsen_DA() line 1219 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/dm/impls/da/da.c
> >>> [0]PETSC ERROR: #5 DMCoarsen() line 2427 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/dm/interface/dm.c
> >>> [0]PETSC ERROR: #6 PCSetUp_MG() line 618 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/ksp/pc/impls/mg/mg.c
> >>> [0]PETSC ERROR: #7 PCSetUp() line 924 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/ksp/pc/interface/pre
> >>> co n.c [0]PETSC ERROR: #8 KSPSetUp() line 381 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/ksp/ksp/interface/it
> >>> fu nc.c [0]PETSC ERROR: #9 KSPSolve() line 612 in
> >>> /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/ksp/ksp/interface/it
> >>> fu nc.c [0]PETSC ERROR: [1]PETSC ERROR: --------------------- Error
> >>> Message
> >>> --------------------------------------------------------------
> >>>
> >>>
> >>> --
> >>> Dipl.-Math. Henrik Büsing
> >>> Institute for Applied Geophysics and Geothermal Energy E.ON Energy
> >>> Research Center RWTH Aachen University
> >>> ------------------------------------------------------
> >>> Mathieustr. 10            |    Tel +49 (0)241 80 49907
> >>> 52074 Aachen, Germany     |    Fax +49 (0)241 80 49889
> >>> ------------------------------------------------------
> >>> http://www.eonerc.rwth-aachen.de/GGE
> >>> hbuesing at eonerc.rwth-aachen.de
> >>> ------------------------------------------------------
> >>>
> >>>> -----Ursprüngliche Nachricht-----
> >>>> Von: Smith, Barry F. [mailto:bsmith at mcs.anl.gov]
> >>>> Gesendet: 09 January 2018 15:10
> >>>> An: Buesing, Henrik <hbuesing at eonerc.rwth-aachen.de>
> >>>> Cc: petsc-users <petsc-users at mcs.anl.gov>
> >>>> Betreff: Re: [petsc-users] Geometric multigrid with a cell-centered
> >>>> discretization
> >>>>
> >>>>
> >>>> With cell centered elements each element is sliced in half in each
> >>>> of the three directions for refinement (or you can refine only in 1
> >>>> or 2 of the dimensions if you like). This means that the ratio of
> >>>> mx to MX needs to be 2 (or 4).So your numbers are fine.
> >>>>
> >>>>> On Jan 9, 2018, at 7:59 AM, Buesing, Henrik <hbuesing at eonerc.rwth-
> >>>> aachen.de> wrote:
> >>>>>
> >>>>> Dear all,
> >>>>>
> >>>>> I am using DMDACreate3d to generate my domain and a cell-centered
> >>>>> two-point flux approximation as discretization. I use geometric
> >>>>> Multigrid (-pc_type mg) with DMDASetInterpolationType(petsc_da,
> >>>>> DMDA_Q0, petsc_ierr). With DM_BOUNDARY_GHOSTED as boundary
> >> type in
> >>>>> DMDACreate3d, I get
> >>>>>
> >>>>> [0]PETSC ERROR: Arguments are incompatible [0]PETSC ERROR: Ratio
> >>>>> between levels: (mx - 1)/(Mx - 1) must be integer: mx 1536 Mx 768
> >>>>
> >>>> We need ALL the output from the error (always). It should never be
> >>>> doing this check in your case so something is wrong with the order
> >>>> of your code perhaps.
> >>>>
> >>>>  Barry
> >>>>
> >>>>>
> >>>>> With DM_BOUNDARY_PERIODIC, I get invalid accesses when I assemble
> >>>>> the
> >>>> matrix.
> >>>>>
> >>>>> Can I somehow tell him, that I have a cell-centered grid, or do I
> >>>>> have
> >>>> to resort to using N+1 in my dimensions?
> >>>>>
> >>>>> Thank you!
> >>>>> Henrik
> >>>>>
> >>>>>
> >>>>> --
> >>>>> Dipl.-Math. Henrik Büsing
> >>>>> Institute for Applied Geophysics and Geothermal Energy E.ON Energy
> >>>>> Research Center RWTH Aachen University
> >>>>> ------------------------------------------------------
> >>>>> Mathieustr. 10            |    Tel +49 (0)241 80 49907
> >>>>> 52074 Aachen, Germany     |    Fax +49 (0)241 80 49889
> >>>>> ------------------------------------------------------
> >>>>> http://www.eonerc.rwth-aachen.de/GGE
> >>>>> hbuesing at eonerc.rwth-aachen.de
> >>>>> ------------------------------------------------------
> >>>
> >



More information about the petsc-users mailing list