Actual source code: dmplextsceed.c
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/tsimpl.h>
3: #include <petsc/private/snesimpl.h>
4: #include <petscds.h>
5: #include <petscfv.h>
7: PetscErrorCode DMPlexTSComputeRHSFunctionFVMCEED(DM dm, PetscReal time, Vec locX, Vec F, PetscCtx ctx)
8: {
9: PetscFV fv;
10: Vec locF;
11: Ceed ceed;
12: DMCeed sd = dm->dmceed;
13: CeedVector clocX, clocF;
15: PetscFunctionBegin;
16: PetscCall(DMGetCeed(dm, &ceed));
17: PetscCheck(sd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DM has no CEED data. Call DMCeedCreate() before computing the residual.");
18: if (time == 0.) PetscCall(DMCeedComputeGeometry(dm, sd));
19: PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fv));
20: PetscCall(DMPlexInsertBoundaryValuesFVM(dm, fv, locX, time, NULL));
21: PetscCall(DMGetLocalVector(dm, &locF));
22: PetscCall(VecZeroEntries(locF));
23: PetscCall(VecGetCeedVectorRead(locX, ceed, &clocX));
24: PetscCall(VecGetCeedVector(locF, ceed, &clocF));
25: PetscCallCEED(CeedOperatorApplyAdd(sd->op, clocX, clocF, CEED_REQUEST_IMMEDIATE));
26: PetscCall(VecRestoreCeedVectorRead(locX, &clocX));
27: PetscCall(VecRestoreCeedVector(locF, &clocF));
28: PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
29: PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
30: PetscCall(DMRestoreLocalVector(dm, &locF));
31: PetscCall(VecViewFromOptions(F, NULL, "-fv_rhs_view"));
32: PetscFunctionReturn(PETSC_SUCCESS);
33: }