[petsc-users] Solving a simpler PDE results in error while solving the original equation works in snes/tutorials/ex13.c

Tony Wu wutony0408 at gmail.com
Mon Jul 10 15:56:51 CDT 2023


In snes/tutorials/ex13.c,
there is a function SetupPrimalProblem(), which sets up the f_0 and f_1 in
PetscDSSetResidual() as described in
https://petsc.org/release/manualpages/DT/PetscDSSetResidual/#petscdssetresidual
,
and also the boundary conditions using the exact solution
trig_inhomogeneous_u().

Now, I am trying to modify f_0, f_1, and the boundary conditions in order
to solve the 2D equation
u+\Delta u=f,
where f(x,y)=x and the boundary conditions are determined by the exact
solution u(x,y)=x.

In order to do this, I defined the following point-wise functions:

static void f00(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt
uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar
u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt
aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar
a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const
PetscScalar constants[], PetscScalar f0[])
{
  f0[0] -= u[0];
  f0[0] += x[0];
}

static void f11(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt
uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar
u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt
aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar
a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const
PetscScalar constants[], PetscScalar f1[])
{
  PetscInt d;
  for (d = 0; d < dim; ++d)
  {
    f1[d] = u_x[d];
  }
}

static void u_exact(PetscInt dim, PetscReal time, const PetscReal x[],
PetscInt Nc, PetscScalar *u, void *ctx)
{
  *u = 0.0;
  *u += x[0];
  return PETSC_SUCCESS;
}

I then replaced the f0, f1_u and ex in SetupPrimalProblem()
by the corresponding f00, f11, and u_exact.
.

However, I got the following error message:

[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[0]PETSC ERROR:
[0]PETSC ERROR: See https://petsc.org/release/faq/ for trouble shooting.
[0]PETSC ERROR: Petsc Release Version 3.19.3, unknown
[0]PETSC ERROR: ./ex13 on a arch-linux-c-debug named LAPTOP by user Tue Jul
11 04:06:13 2023
[0]PETSC ERROR: Configure options --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3
--FOPTFLAGS=-O3 --with-shared-libraries --with-cxx-dialect=C++11
--with-make-np=8 --download-hdf5=yes --download-mumps=yes
--download-parmetis=yes --download-metis=yes --download-superlu=yes
--download-hypre=yes --download-superlu_dist=yes --download-scalapack=yes
--download-triangle=yes
[0]PETSC ERROR: #1 DMProjectPoint_Func_Private() at
/home/user/petsc/src/dm/impls/plex/plexproject.c:128
[0]PETSC ERROR: #2 DMProjectPoint_Private() at
/home/user/petsc/src/dm/impls/plex/plexproject.c:395
[0]PETSC ERROR: #3 DMProjectLocal_Generic_Plex() at
/home/user/petsc/src/dm/impls/plex/plexproject.c:852
[0]PETSC ERROR: #4 DMProjectFunctionLabelLocal_Plex() at
/home/user/petsc/src/dm/impls/plex/plexproject.c:921
[0]PETSC ERROR: #5 DMProjectFunctionLabelLocal() at
/home/user/petsc/src/dm/interface/dm.c:8071
[0]PETSC ERROR: #6 DMPlexInsertBoundaryValuesEssential() at
/home/user/petsc/src/dm/impls/plex/plexfem.c:795
[0]PETSC ERROR: #7 DMPlexInsertBoundaryValues_Plex() at
/home/user/petsc/src/dm/impls/plex/plexfem.c:1046
[0]PETSC ERROR: #8 DMPlexInsertBoundaryValues() at
/home/user/petsc/src/dm/impls/plex/plexfem.c:1154
[0]PETSC ERROR: #9 DMPlexSNESComputeBoundaryFEM() at
/home/user/petsc/src/snes/utils/dmplexsnes.c:1372
[0]PETSC ERROR: #10 SNESComputeFunction_DMLocal() at
/home/user/petsc/src/snes/utils/dmlocalsnes.c:63
[0]PETSC ERROR: #11 SNESComputeFunction() at
/home/user/petsc/src/snes/interface/snes.c:2389
[0]PETSC ERROR: #12 SNESSolve_NEWTONLS() at
/home/user/petsc/src/snes/impls/ls/ls.c:172
[0]PETSC ERROR: #13 SNESSolve() at
/home/user/petsc/src/snes/interface/snes.c:4644
[0]PETSC ERROR: #14 PetscConvEstGetConvRateSNES_Private() at
/home/user/petsc/src/snes/utils/convest.c:369
[0]PETSC ERROR: #15 PetscConvEstGetConvRate() at
/home/user/petsc/src/snes/utils/convest.c:468
[0]PETSC ERROR: #16 SNESSolve() at
/home/user/petsc/src/snes/interface/snes.c:4561
[0]PETSC ERROR: #17 main() at ex13.c:431
[0]PETSC ERROR: Reached the main program with an out-of-range error code
584475600. This should never happen
[0]PETSC ERROR: PETSc Option Table entries:
[0]PETSC ERROR: -convest_num_refine 2 (source: command line)
[0]PETSC ERROR: -potential_petscspace_degree 1 (source: command line)
[0]PETSC ERROR: -snes_convergence_estimate (source: command line)
[0]PETSC ERROR: -snes_monitor (source: command line)
[0]PETSC ERROR: ----------------End of Error Message -------send entire
error message to petsc-maint at mcs.anl.gov----------
--------------------------------------------------------------------------
MPI_ABORT was invoked on rank 0 in communicator MPI_COMM_SELF
with errorcode 584475600.

NOTE: invoking MPI_ABORT causes Open MPI to kill all MPI processes.
You may or may not see output from other processes, depending on
exactly when Open MPI kills them.
--------------------------------------------------------------------------


I have no idea why the program solves the original, more complicated,
equation successfully while not being able to solve a simpler one.
Any help or reference would be much appreciated!
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230711/56ea8e87/attachment.html>


More information about the petsc-users mailing list