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