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

Buesing, Henrik hbuesing at eonerc.rwth-aachen.de
Tue Jan 9 09:32:59 CST 2018


[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

It looks like you give it the number of vertices, not cells, so +1 to your grid size.

[Buesing, Henrik] Yes! But, I do not want that. I want to give the cells…

That is the interface. You will get m-1 cells so give it m = num_cells + 1.

[Buesing, Henrik] Ok, I can do that. But, does the space discretization matter for the coarsening and refining. With DMDA_Q0, I tell him it is cell-centered two-point flux approximation, right?


[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_run/000023/000000_copy_sub_launch/work/shem_fw64gnu_const.x on a gnu_openmpi named linuxbmc0008.rz.RWTH-Aachen.DE<http://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/precon.c
[0]PETSC ERROR: #8 KSPSetUp() line 381 in /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: #9 KSPSolve() line 612 in /rwthfs/rz/cluster/work/hb111949/Code/petsc/src/ksp/ksp/interface/itfunc.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<tel:%2B49%20%280%29241%2080%2049907>
52074 Aachen, Germany     |    Fax +49 (0)241 80 49889<tel:%2B49%20%280%29241%2080%2049889>
------------------------------------------------------
http://www.eonerc.rwth-aachen.de/GGE
hbuesing at eonerc.rwth-aachen.de<mailto:hbuesing at eonerc.rwth-aachen.de>
------------------------------------------------------

> -----Ursprüngliche Nachricht-----
> Von: Smith, Barry F. [mailto:bsmith at mcs.anl.gov<mailto:bsmith at mcs.anl.gov>]
> Gesendet: 09 January 2018 15:10
> An: Buesing, Henrik <hbuesing at eonerc.rwth-aachen.de<mailto:hbuesing at eonerc.rwth-aachen.de>>
> Cc: petsc-users <petsc-users at mcs.anl.gov<mailto: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-
<mailto:hbuesing at eonerc.rwth-%0b>> aachen.de<http://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<tel:+49%20241%208049907>
> > 52074 Aachen, Germany     |    Fax +49 (0)241 80 49889<tel:+49%20241%208049889>
> > ------------------------------------------------------
> > http://www.eonerc.rwth-aachen.de/GGE
> > hbuesing at eonerc.rwth-aachen.de<mailto:hbuesing at eonerc.rwth-aachen.de>
> > ------------------------------------------------------


-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20180109/92d250bb/attachment-0001.html>


More information about the petsc-users mailing list