[petsc-users] asymmetric errors at X+/Y+ boundary on periodic box mesh
Yann Jobic
yann.jobic at univ-amu.fr
Fri May 22 14:16:49 CDT 2026
Dear all,
Let me be more precise, all the relevant steps that i'm doing are those :
const char *opts =
"-dm_plex_box_faces 32,32 -dm_plex_simplex 0 "
"-dm_plex_box_bd periodic,periodic -dm_distribute_overlap 1 "
"-dm_plex_adj_cone -dm_plex_adj_closure 0";
PetscCall(PetscOptionsInsertString(NULL, opts));
PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
PetscCall(DMSetType(dm, DMPLEX));
PetscCall(DMSetFromOptions(dm));
{
DM g;
PetscCall(DMPlexConstructGhostCells(dm, NULL, NULL, &g));
PetscCall(DMDestroy(&dm));
dm = g;
}
/* FV: leastsquares + vanleer. */
PetscCall(PetscFVCreate(PETSC_COMM_WORLD, &fvm));
PetscCall(PetscFVSetType(fvm, PETSCFVLEASTSQUARES));
PetscCall(PetscFVSetNumComponents(fvm, 1));
PetscCall(PetscFVSetSpatialDimension(fvm, 2));
PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
{
PetscLimiter lim;
PetscCall(PetscFVGetLimiter(fvm, &lim));
PetscCall(PetscLimiterSetType(lim, PETSCLIMITERVANLEER));
}
PetscCall(DMAddField(dm, NULL, (PetscObject)fvm));
PetscCall(DMCreateDS(dm));
PetscCall(DMGetDS(dm, &ds));
PetscCall(PetscDSSetRiemannSolver(ds, 0, Riemann));
/* IC + TS. */
PetscCall(DMCreateGlobalVector(dm, &U));
PetscCall(PetscObjectSetName((PetscObject)U, "phi"));
{
PetscErrorCode (*f[1])(PetscInt, PetscReal, const PetscReal[],
PetscInt, PetscScalar[], void*) = {Sine};
void *fctx[1] = {NULL};
PetscCall(DMProjectFunction(dm, 0.0, f, fctx, INSERT_ALL_VALUES, U));
}
PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
PetscCall(TSSetType(ts, TSSSP));
PetscCall(TSSSPSetType(ts, TSSSPRKS3));
PetscCall(TSSetDM(ts, dm));
PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, NULL));
PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM,
NULL));
PetscCall(TSSetTimeStep(ts, 1e-3));
PetscCall(TSSetMaxTime(ts, 0.01));
PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
PetscCall(TSSetSolution(ts, U));
PetscCall(TSSolve(ts, U));
It's the important parts of the minimal working example that i gave in
the first email.
Maybe i'm doing something wrong for periodic boundaries ?
I may have forgot an important step ?
Maybe there are CI tests that validate the order of accuracy of
DMPlexTSComputeRHSFunctionFVM on a periodic hex DMPlex with a
smooth manufactured solution? I may have missed it.
Many thanks,
Yann
Le 22/05/2026 à 15:50, Yann Jobic a écrit :
> Ce mail provient de l'extérieur, restons vigilants
>
> Dear all,
>
> I'm writing a small finite-volume scalar advection code on a 2D periodic
> DMPlex hex mesh, using the same pipeline as src/ts/tutorials/ex11.c
> (PetscFV + PetscLimiter + DMPlexTSComputeRHSFunctionFVM + TSSSP).
>
> After only 10 time steps from a smooth periodic initial condition
> phi_0(x,y) = sin(2 pi x) * sin(2 pi y), with constant velocity a =
> (1,1), on a 32x32 hex periodic box, the error |phi - phi_exact| is
> strongly asymmetric: it is roughly 0.18 in the cells touching x=1 and
> y=1, but about a tenth of that everywhere else (including near x=0 and
> y=0). The exact solution is just phi_0 translated by a*t and stays smooth.
>
> I'm using PETSc 3.24.4 (release) and 3.25.1 (release).
>
> I'm attaching a 170-line MWE (mwe.c, no command-line options needed) and
> a screenshot of the problem. The MWE prints
> t=1.000e-02 Linf=1.7897e-01
>
> Is there an additional setup step that I'm missing for the periodic hex
> case (e.g. an explicit call about face localization, or about the SF
> after DMPlexConstructGhostCells)?
>
> Thanks,
>
> Yann
More information about the petsc-users
mailing list