[petsc-dev] The problem about running an example in PETSc

Jed Brown jedbrown at mcs.anl.gov
Tue Apr 30 13:20:11 CDT 2013


Lulu Liu <lulu.liu at kaust.edu.sa> writes:

> src/ts/examples/tutorials/ex29.c I want domain partition is 1x2
>
> mpiexec -n 2 ./ex29 -da_refine 4 -ts_max_steps 10 -ts_rtol 1e-3 -ts_atol
> 1e-3   -ts_type arkimex -ts_monitor -snes_monitor -snes_type ngmres
> -npc_snes_type nasm -npc_snes_nasm_type restrict -da_overlap 4
>  -da_processors_x 1 -da_processors_y 2

Thanks for sending the command to reproduce.

> 0 TS dt 0.05 time 0
>     0 SNES Function norm 1.762942085583e+01
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Argument out of range!
> [0]PETSC ERROR: Too many processors in Y direction: 2 1!
> [0]PETSC ERROR: [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Argument out of range!
> [1]PETSC ERROR: Too many processors in Y direction: 2 1!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Development GIT revision:
> bc51b945dfb4f1d1228fb677a132a806066a0f32  GIT Date: 2013-04-29 10:22:29
> -0500
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Development GIT revision:
> bc51b945dfb4f1d1228fb677a132a806066a0f32  GIT Date: 2013-04-29 10:22:29
> -0500
> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [0]PETSC ERROR: See docs/index.html for manual pages.
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: ./ex29 on a arch-basic named kl-12681.local by liul Tue Apr
> 30 21:07:37 2013
> [0]PETSC ERROR: [1]PETSC ERROR: See docs/changes/index.html for recent
> updates.
> [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [1]PETSC ERROR: See docs/index.html for manual pages.
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: ./ex29 on a arch-basic named kl-12681.local by liul Tue Apr
> 30 21:07:37 2013
> [1]PETSC ERROR: Libraries linked from /Users/liul/soft/petsc/arch-basic/lib
> [1]PETSC ERROR: Configure run at Tue Apr 30 15:46:16 2013
> [1]PETSC ERROR: Configure options --download-cmake --download-f-blas-lapack
> --download-mpich --with-cc=gcc --with-fc=gfortran PETSC_ARCH=arch-basic
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Libraries linked from /Users/liul/soft/petsc/arch-basic/lib
> [0]PETSC ERROR: Configure run at Tue Apr 30 15:46:16 2013
> [0]PETSC ERROR: Configure options --download-cmake --download-f-blas-lapack
> --download-mpich --with-cc=gcc --with-fc=gfortran PETSC_ARCH=arch-basic
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: DMSetUp_DA_2D() line 266 in
> /Users/liul/soft/petsc/src/dm/impls/da/da2.c
> DMSetUp_DA_2D() line 266 in /Users/liul/soft/petsc/src/dm/impls/da/da2.c
> [1]PETSC ERROR: DMSetUp_DA() line 27 in
> /Users/liul/soft/petsc/src/dm/impls/da/dareg.c
> [1]PETSC ERROR: DMSetUp() line 469 in
> /Users/liul/soft/petsc/src/dm/interface/dm.c
> [1]PETSC ERROR: DMDACreate2d() line 888 in
> /Users/liul/soft/petsc/src/dm/impls/da/da2.c
> [1]PETSC ERROR: [0]PETSC ERROR: DMSetUp_DA() line 27 in
> /Users/liul/soft/petsc/src/dm/impls/da/dareg.c
> [0]PETSC ERROR: DMSetUp() line 469 in
> /Users/liul/soft/petsc/src/dm/interface/dm.c
> DMDAGetReducedDMDA() line 389 in
> /Users/liul/soft/petsc/src/dm/impls/da/dacorn.c
> [1]PETSC ERROR: CoefficientSubDomainRestrictHook() line 92 in
> src/ts/examples/tutorials/ex29.c
> [1]PETSC ERROR: DMCreateDomainDecomposition() line 1314 in
> /Users/liul/soft/petsc/src/dm/interface/dm.c
> [1]PETSC ERROR: SNESSetUp_NASM() line 112 in
> /Users/liul/soft/petsc/src/snes/impls/nasm/nasm.c
> [1]PETSC ERROR: SNESSetUp() line 2633 in
> /Users/liul/soft/petsc/src/snes/interface/snes.c
> [0]PETSC ERROR: DMDACreate2d() line 888 in
> /Users/liul/soft/petsc/src/dm/impls/da/da2.c
> [0]PETSC ERROR: DMDAGetReducedDMDA() line 389 in
> /Users/liul/soft/petsc/src/dm/impls/da/dacorn.c

  if (dim == 1) {
    ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr);
  } else if (dim == 2) {
    ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr);
  } else if (dim == 3) {
    ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr);
  }

This code should be fixed to use the precise API instead of the
convenience routines DMDACreate*d(), so that it does NOT call
DMSetFromOptions (or if it does, call it after setting a prefix).

Lulu, would you be up for taking a crack at this and submitting a patch?
Just take the sequence of calls found in DMDACreate3d:

  ierr = DMDACreate(comm, da);CHKERRQ(ierr);
  ierr = DMDASetDim(*da, 3);CHKERRQ(ierr);                  // <-- fix this line
  ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
  ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
  ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
  ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
  ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
  ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
  ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
  /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
  ierr = DMSetFromOptions(*da);CHKERRQ(ierr);              // <--- delete this
  ierr = DMSetUp(*da);CHKERRQ(ierr);
  ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr);  // <--- delete or set a prefix

I can help if you get stuck.

https://bitbucket.org/petsc/petsc/wiki/Home



More information about the petsc-dev mailing list