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