<div dir="ltr"><div dir="ltr">On Mon, Jul 10, 2023 at 5:02 PM Tony Wu <<a href="mailto:wutony0408@gmail.com">wutony0408@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">In snes/tutorials/ex13.c,<br>there is a function SetupPrimalProblem(), which sets up the f_0 and f_1 in PetscDSSetResidual() as described in <a href="https://petsc.org/release/manualpages/DT/PetscDSSetResidual/#petscdssetresidual" target="_blank">https://petsc.org/release/manualpages/DT/PetscDSSetResidual/#petscdssetresidual</a>,<br>and also the boundary conditions using the exact solution trig_inhomogeneous_u().<br><br>Now, I am trying to modify f_0, f_1, and the boundary conditions in order to solve the 2D equation<br>u+\Delta u=f,<br>where f(x,y)=x and the boundary conditions are determined by the exact solution u(x,y)=x.<br><br>In order to do this, I defined the following point-wise functions:<br><br>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[])<br>{<br> f0[0] -= u[0];<br> f0[0] += x[0];<br>}<br><br>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[])<br>{<br> PetscInt d;<br> for (d = 0; d < dim; ++d)<br> {<br> f1[d] = u_x[d];<br> }<br>}<br><br>static void u_exact(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)<br></div></blockquote><div><br></div><div>This function is declared "void" but should be "PetscErrorCode". There must have been a warning on compile.</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">{<br> *u = 0.0;<br> *u += x[0];<br> return PETSC_SUCCESS;<br>}<br><br>I then replaced the f0, f1_u and ex in SetupPrimalProblem()<br>by the corresponding f00, f11, and u_exact.<br>.<br><br>However, I got the following error message:<div><br>[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[0]PETSC ERROR: <br>[0]PETSC ERROR: See <a href="https://petsc.org/release/faq/" target="_blank">https://petsc.org/release/faq/</a> for trouble shooting.<br>[0]PETSC ERROR: Petsc Release Version 3.19.3, unknown <br>[0]PETSC ERROR: ./ex13 on a arch-linux-c-debug named LAPTOP by user Tue Jul 11 04:06:13 2023<br>[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<br>[0]PETSC ERROR: #1 DMProjectPoint_Func_Private() at /home/user/petsc/src/dm/impls/plex/plexproject.c:128<br>[0]PETSC ERROR: #2 DMProjectPoint_Private() at /home/user/petsc/src/dm/impls/plex/plexproject.c:395<br>[0]PETSC ERROR: #3 DMProjectLocal_Generic_Plex() at /home/user/petsc/src/dm/impls/plex/plexproject.c:852<br>[0]PETSC ERROR: #4 DMProjectFunctionLabelLocal_Plex() at /home/user/petsc/src/dm/impls/plex/plexproject.c:921<br>[0]PETSC ERROR: #5 DMProjectFunctionLabelLocal() at /home/user/petsc/src/dm/interface/dm.c:8071<br>[0]PETSC ERROR: #6 DMPlexInsertBoundaryValuesEssential() at /home/user/petsc/src/dm/impls/plex/plexfem.c:795<br>[0]PETSC ERROR: #7 DMPlexInsertBoundaryValues_Plex() at /home/user/petsc/src/dm/impls/plex/plexfem.c:1046<br>[0]PETSC ERROR: #8 DMPlexInsertBoundaryValues() at /home/user/petsc/src/dm/impls/plex/plexfem.c:1154<br>[0]PETSC ERROR: #9 DMPlexSNESComputeBoundaryFEM() at /home/user/petsc/src/snes/utils/dmplexsnes.c:1372<br>[0]PETSC ERROR: #10 SNESComputeFunction_DMLocal() at /home/user/petsc/src/snes/utils/dmlocalsnes.c:63<br>[0]PETSC ERROR: #11 SNESComputeFunction() at /home/user/petsc/src/snes/interface/snes.c:2389<br>[0]PETSC ERROR: #12 SNESSolve_NEWTONLS() at /home/user/petsc/src/snes/impls/ls/ls.c:172<br>[0]PETSC ERROR: #13 SNESSolve() at /home/user/petsc/src/snes/interface/snes.c:4644<br>[0]PETSC ERROR: #14 PetscConvEstGetConvRateSNES_Private() at /home/user/petsc/src/snes/utils/convest.c:369<br>[0]PETSC ERROR: #15 PetscConvEstGetConvRate() at /home/user/petsc/src/snes/utils/convest.c:468<br>[0]PETSC ERROR: #16 SNESSolve() at /home/user/petsc/src/snes/interface/snes.c:4561<br>[0]PETSC ERROR: #17 main() at ex13.c:431<br>[0]PETSC ERROR: Reached the main program with an out-of-range error code 584475600. This should never happen<br>[0]PETSC ERROR: PETSc Option Table entries:<br>[0]PETSC ERROR: -convest_num_refine 2 (source: command line)<br>[0]PETSC ERROR: -potential_petscspace_degree 1 (source: command line)<br>[0]PETSC ERROR: -snes_convergence_estimate (source: command line)<br>[0]PETSC ERROR: -snes_monitor (source: command line)<br>[0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint@mcs.anl.gov----------<br>--------------------------------------------------------------------------<br>MPI_ABORT was invoked on rank 0 in communicator MPI_COMM_SELF<br>with errorcode 584475600.<br><br>NOTE: invoking MPI_ABORT causes Open MPI to kill all MPI processes.<br>You may or may not see output from other processes, depending on<br>exactly when Open MPI kills them.<br>--------------------------------------------------------------------------<br><br><br>I have no idea why the program solves the original, more complicated, equation successfully while not being able to solve a simpler one.<br>Any help or reference would be much appreciated!<br></div></div>
</blockquote></div><br clear="all"><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>