1ce78bad3SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2ce78bad3SBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
3ce78bad3SBarry Smith #include <petscds.h>
4ce78bad3SBarry Smith #include <petsc/private/petscimpl.h>
5ce78bad3SBarry Smith #include <petsc/private/petscfeimpl.h>
6ce78bad3SBarry Smith #include <petscdmceed.h>
7ce78bad3SBarry Smith #include <petscdmplexceed.h>
8ce78bad3SBarry Smith
DMPlexSNESComputeResidualCEED(DM dm,Vec locX,Vec locF,PetscCtx ctx)9*2a8381b2SBarry Smith PetscErrorCode DMPlexSNESComputeResidualCEED(DM dm, Vec locX, Vec locF, PetscCtx ctx)
10ce78bad3SBarry Smith {
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 PetscCall(DMCeedComputeGeometry(dm, sd));
19ce78bad3SBarry Smith
20ce78bad3SBarry Smith PetscCall(VecGetCeedVectorRead(locX, ceed, &clocX));
21ce78bad3SBarry Smith PetscCall(VecGetCeedVector(locF, ceed, &clocF));
22ce78bad3SBarry Smith PetscCallCEED(CeedOperatorApplyAdd(sd->op, clocX, clocF, CEED_REQUEST_IMMEDIATE));
23ce78bad3SBarry Smith PetscCall(VecRestoreCeedVectorRead(locX, &clocX));
24ce78bad3SBarry Smith PetscCall(VecRestoreCeedVector(locF, &clocF));
25ce78bad3SBarry Smith
26ce78bad3SBarry Smith {
27ce78bad3SBarry Smith DM_Plex *mesh = (DM_Plex *)dm->data;
28ce78bad3SBarry Smith
29ce78bad3SBarry Smith if (mesh->printFEM) {
30ce78bad3SBarry Smith PetscSection section;
31ce78bad3SBarry Smith Vec locFbc;
32ce78bad3SBarry Smith PetscInt pStart, pEnd, p, maxDof;
33ce78bad3SBarry Smith PetscScalar *zeroes;
34ce78bad3SBarry Smith
35ce78bad3SBarry Smith PetscCall(DMGetLocalSection(dm, §ion));
36ce78bad3SBarry Smith PetscCall(VecDuplicate(locF, &locFbc));
37ce78bad3SBarry Smith PetscCall(VecCopy(locF, locFbc));
38ce78bad3SBarry Smith PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
39ce78bad3SBarry Smith PetscCall(PetscSectionGetMaxDof(section, &maxDof));
40ce78bad3SBarry Smith PetscCall(PetscCalloc1(maxDof, &zeroes));
41ce78bad3SBarry Smith for (p = pStart; p < pEnd; ++p) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
42ce78bad3SBarry Smith PetscCall(PetscFree(zeroes));
43ce78bad3SBarry Smith PetscCall(DMPrintLocalVec(dm, "Residual", mesh->printTol, locFbc));
44ce78bad3SBarry Smith PetscCall(VecDestroy(&locFbc));
45ce78bad3SBarry Smith }
46ce78bad3SBarry Smith }
47ce78bad3SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
48ce78bad3SBarry Smith }
49