1dc2eb70dSMatthew G. Knepley static const char help[] = "Tests for injecting basis functions";
2dc2eb70dSMatthew G. Knepley
3dc2eb70dSMatthew G. Knepley #include <petscdmplex.h>
4dc2eb70dSMatthew G. Knepley #include <petscfe.h>
5dc2eb70dSMatthew G. Knepley #include <petscds.h>
6dc2eb70dSMatthew G. Knepley
7dc2eb70dSMatthew G. Knepley typedef struct {
8dc2eb70dSMatthew G. Knepley PetscInt its; /* Number of replications for timing */
9dc2eb70dSMatthew G. Knepley } AppCtx;
10dc2eb70dSMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12d71ae5a4SJacob Faibussowitsch {
13dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
14dc2eb70dSMatthew G. Knepley options->its = 1;
15dc2eb70dSMatthew G. Knepley
16d0609cedSBarry Smith PetscOptionsBegin(comm, "", "FE Injection Options", "PETSCFE");
179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL));
18d0609cedSBarry Smith PetscOptionsEnd();
193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
20dc2eb70dSMatthew G. Knepley }
21dc2eb70dSMatthew G. Knepley
trig_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)22*2a8381b2SBarry Smith static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
23d71ae5a4SJacob Faibussowitsch {
24dc2eb70dSMatthew G. Knepley PetscInt d;
25dc2eb70dSMatthew G. Knepley *u = 0.0;
26dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
273ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
28dc2eb70dSMatthew G. Knepley }
29dc2eb70dSMatthew G. Knepley
f0_trig_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])30d71ae5a4SJacob Faibussowitsch static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
31d71ae5a4SJacob Faibussowitsch {
32dc2eb70dSMatthew G. Knepley PetscInt d;
33dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
34dc2eb70dSMatthew G. Knepley }
35dc2eb70dSMatthew G. Knepley
f1_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])36d71ae5a4SJacob Faibussowitsch static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
37d71ae5a4SJacob Faibussowitsch {
38dc2eb70dSMatthew G. Knepley PetscInt d;
39dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d];
40dc2eb70dSMatthew G. Knepley }
41dc2eb70dSMatthew G. Knepley
g3_uu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g3[])42d71ae5a4SJacob Faibussowitsch static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
43d71ae5a4SJacob Faibussowitsch {
44dc2eb70dSMatthew G. Knepley PetscInt d;
45dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
46dc2eb70dSMatthew G. Knepley }
47dc2eb70dSMatthew G. Knepley
SetupPrimalProblem(DM dm,AppCtx * user)48d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
49d71ae5a4SJacob Faibussowitsch {
50dc2eb70dSMatthew G. Knepley PetscDS ds;
51dc2eb70dSMatthew G. Knepley DMLabel label;
52dc2eb70dSMatthew G. Knepley const PetscInt id = 1;
53dc2eb70dSMatthew G. Knepley
54dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
559566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
569566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
579566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
589566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
599566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
6057d50842SBarry Smith if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)trig_u, NULL, user, NULL));
613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
62dc2eb70dSMatthew G. Knepley }
63dc2eb70dSMatthew G. Knepley
SetupDiscretization(DM dm,const char name[],PetscErrorCode (* setup)(DM,AppCtx *),AppCtx * user)64d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
65d71ae5a4SJacob Faibussowitsch {
66dc2eb70dSMatthew G. Knepley DM cdm = dm;
67dc2eb70dSMatthew G. Knepley PetscFE fe;
68dc2eb70dSMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN];
69dc2eb70dSMatthew G. Knepley PetscInt dim;
70dc2eb70dSMatthew G. Knepley
71dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
739566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
749566063dSJacob Faibussowitsch PetscCall(DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe));
759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name));
76dc2eb70dSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */
779566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
789566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
799566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user));
80dc2eb70dSMatthew G. Knepley while (cdm) {
819566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
829566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
83dc2eb70dSMatthew G. Knepley }
849566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
86dc2eb70dSMatthew G. Knepley }
87dc2eb70dSMatthew G. Knepley
8849abdd8aSBarry Smith /* PetscObjectContainerCompose() compose requires void ** signature on destructor */
PetscFEGeomDestroy_Void(PetscCtxRt ctx)89*2a8381b2SBarry Smith static PetscErrorCode PetscFEGeomDestroy_Void(PetscCtxRt ctx)
90d71ae5a4SJacob Faibussowitsch {
9149abdd8aSBarry Smith return PetscFEGeomDestroy((PetscFEGeom **)ctx);
92dc2eb70dSMatthew G. Knepley }
93dc2eb70dSMatthew G. Knepley
CellRangeGetFEGeom(IS cellIS,DMField coordField,PetscQuadrature quad,PetscFEGeomMode mode,PetscFEGeom ** geom)94ac9d17c7SMatthew G. Knepley PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscFEGeomMode mode, PetscFEGeom **geom)
95d71ae5a4SJacob Faibussowitsch {
96dc2eb70dSMatthew G. Knepley char composeStr[33] = {0};
97dc2eb70dSMatthew G. Knepley PetscObjectId id;
98dc2eb70dSMatthew G. Knepley PetscContainer container;
99dc2eb70dSMatthew G. Knepley
100dc2eb70dSMatthew G. Knepley PetscFunctionBegin;
1019566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)quad, &id));
10263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%" PetscInt64_FMT "\n", id));
1039566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)cellIS, composeStr, (PetscObject *)&container));
104dc2eb70dSMatthew G. Knepley if (container) {
105*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(container, geom));
106dc2eb70dSMatthew G. Knepley } else {
107ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, mode, geom));
10849abdd8aSBarry Smith PetscCall(PetscObjectContainerCompose((PetscObject)cellIS, composeStr, *geom, PetscFEGeomDestroy_Void));
109dc2eb70dSMatthew G. Knepley }
1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
111dc2eb70dSMatthew G. Knepley }
112dc2eb70dSMatthew G. Knepley
CellRangeRestoreFEGeom(IS cellIS,DMField coordField,PetscQuadrature quad,PetscBool faceData,PetscFEGeom ** geom)113d71ae5a4SJacob Faibussowitsch PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
114d71ae5a4SJacob Faibussowitsch {
115dc2eb70dSMatthew G. Knepley PetscFunctionBegin;
116dc2eb70dSMatthew G. Knepley *geom = NULL;
1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
118dc2eb70dSMatthew G. Knepley }
119dc2eb70dSMatthew G. Knepley
CreateFEGeometry(DM dm,PetscDS ds,IS cellIS,PetscQuadrature * affineQuad,PetscFEGeom ** affineGeom,PetscQuadrature ** quads,PetscFEGeom *** geoms)120d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
121d71ae5a4SJacob Faibussowitsch {
122dc2eb70dSMatthew G. Knepley DMField coordField;
123dc2eb70dSMatthew G. Knepley PetscInt Nf, f, maxDegree;
124dc2eb70dSMatthew G. Knepley
125dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
126dc2eb70dSMatthew G. Knepley *affineQuad = NULL;
127dc2eb70dSMatthew G. Knepley *affineGeom = NULL;
128dc2eb70dSMatthew G. Knepley *quads = NULL;
129dc2eb70dSMatthew G. Knepley *geoms = NULL;
1309566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
1319566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField));
1329566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
133dc2eb70dSMatthew G. Knepley if (maxDegree <= 1) {
1349566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad));
135ac9d17c7SMatthew G. Knepley if (*affineQuad) PetscCall(CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FEGEOM_BASIC, affineGeom));
136dc2eb70dSMatthew G. Knepley } else {
1379566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nf, quads, Nf, geoms));
138dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) {
139dc2eb70dSMatthew G. Knepley PetscFE fe;
140dc2eb70dSMatthew G. Knepley
1419566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1429566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &(*quads)[f]));
1439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*quads)[f]));
144ac9d17c7SMatthew G. Knepley PetscCall(CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FEGEOM_BASIC, &(*geoms)[f]));
145dc2eb70dSMatthew G. Knepley }
146dc2eb70dSMatthew G. Knepley }
1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
148dc2eb70dSMatthew G. Knepley }
149dc2eb70dSMatthew G. Knepley
DestroyFEGeometry(DM dm,PetscDS ds,IS cellIS,PetscQuadrature * affineQuad,PetscFEGeom ** affineGeom,PetscQuadrature ** quads,PetscFEGeom *** geoms)150d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms)
151d71ae5a4SJacob Faibussowitsch {
152dc2eb70dSMatthew G. Knepley DMField coordField;
153dc2eb70dSMatthew G. Knepley PetscInt Nf, f;
154dc2eb70dSMatthew G. Knepley
155dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
1569566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
1579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(dm, &coordField));
158dc2eb70dSMatthew G. Knepley if (*affineQuad) {
1599566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom));
1609566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(affineQuad));
161dc2eb70dSMatthew G. Knepley } else {
162dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) {
1639566063dSJacob Faibussowitsch PetscCall(CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]));
1649566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(*quads)[f]));
165dc2eb70dSMatthew G. Knepley }
1669566063dSJacob Faibussowitsch PetscCall(PetscFree2(*quads, *geoms));
167dc2eb70dSMatthew G. Knepley }
1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
169dc2eb70dSMatthew G. Knepley }
170dc2eb70dSMatthew G. Knepley
TestEvaluation(DM dm)171d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestEvaluation(DM dm)
172d71ae5a4SJacob Faibussowitsch {
173dc2eb70dSMatthew G. Knepley PetscFE fe;
174dc2eb70dSMatthew G. Knepley PetscSpace sp;
175dc2eb70dSMatthew G. Knepley PetscReal *points;
176dc2eb70dSMatthew G. Knepley PetscReal *B, *D, *H;
177dc2eb70dSMatthew G. Knepley PetscInt dim, Nb, b, Nc, c, Np, p;
178dc2eb70dSMatthew G. Knepley
179dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
1809566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
1819566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
182dc2eb70dSMatthew G. Knepley Np = 6;
1839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Np * dim, &points));
184dc2eb70dSMatthew G. Knepley if (dim == 3) {
1859371c9d4SSatish Balay points[0] = -1.0;
1869371c9d4SSatish Balay points[1] = -1.0;
1879371c9d4SSatish Balay points[2] = -1.0;
1889371c9d4SSatish Balay points[3] = 1.0;
1899371c9d4SSatish Balay points[4] = -1.0;
1909371c9d4SSatish Balay points[5] = -1.0;
1919371c9d4SSatish Balay points[6] = -1.0;
1929371c9d4SSatish Balay points[7] = 1.0;
1939371c9d4SSatish Balay points[8] = -1.0;
1949371c9d4SSatish Balay points[9] = -1.0;
1959371c9d4SSatish Balay points[10] = -1.0;
1969371c9d4SSatish Balay points[11] = 1.0;
1979371c9d4SSatish Balay points[12] = 1.0;
1989371c9d4SSatish Balay points[13] = -1.0;
1999371c9d4SSatish Balay points[14] = 1.0;
2009371c9d4SSatish Balay points[15] = -1.0;
2019371c9d4SSatish Balay points[16] = 1.0;
2029371c9d4SSatish Balay points[17] = 1.0;
203dc2eb70dSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now");
2049566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp));
2059566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(sp, &Nb));
2069566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
2079566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np * Nb * Nc, MPIU_REAL, &B));
2089566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim, MPIU_REAL, &D));
2099566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Np * Nb * Nc * dim * dim, MPIU_REAL, &H));
2109566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/));
211dc2eb70dSMatthew G. Knepley for (p = 0; p < Np; ++p) {
2129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p));
213dc2eb70dSMatthew G. Knepley for (b = 0; b < Nb; ++b) {
2149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b));
2159566063dSJacob Faibussowitsch for (c = 0; c < Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)B[(p * Nb + b) * Nc + c]));
2169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
217dc2eb70dSMatthew G. Knepley #if 0
218dc2eb70dSMatthew G. Knepley for (c = 0; c < Nc; ++c) {
2199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c));
2209566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)]));
2219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
222dc2eb70dSMatthew G. Knepley }
223dc2eb70dSMatthew G. Knepley #endif
224dc2eb70dSMatthew G. Knepley }
225dc2eb70dSMatthew G. Knepley }
2269566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np * Nb, MPIU_REAL, &B));
2279566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim, MPIU_REAL, &D));
2289566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Np * Nb * dim * dim, MPIU_REAL, &H));
2299566063dSJacob Faibussowitsch PetscCall(PetscFree(points));
2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
231dc2eb70dSMatthew G. Knepley }
232dc2eb70dSMatthew G. Knepley
TestIntegration(DM dm,PetscInt cbs,PetscInt its)233d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its)
234d71ae5a4SJacob Faibussowitsch {
235dc2eb70dSMatthew G. Knepley PetscDS ds;
236dc2eb70dSMatthew G. Knepley PetscFEGeom *chunkGeom = NULL;
237dc2eb70dSMatthew G. Knepley PetscQuadrature affineQuad, *quads = NULL;
238dc2eb70dSMatthew G. Knepley PetscFEGeom *affineGeom, **geoms = NULL;
239dc2eb70dSMatthew G. Knepley PetscScalar *u, *elemVec;
240dc2eb70dSMatthew G. Knepley IS cellIS;
241dc2eb70dSMatthew G. Knepley PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k;
242dc2eb70dSMatthew G. Knepley
243dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
2449566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth));
2459566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "depth", depth, &cellIS));
2469566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
24707218a29SMatthew G. Knepley PetscCall(DMGetCellDS(dm, cStart, &ds, NULL));
2489566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
2499566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2509566063dSJacob Faibussowitsch PetscCall(CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2519566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(chunkSize * totDim, &u, chunkSize * totDim, &elemVec));
252dc2eb70dSMatthew G. Knepley /* Assumptions:
253dc2eb70dSMatthew G. Knepley - Single field
254dc2eb70dSMatthew G. Knepley - No input data
255dc2eb70dSMatthew G. Knepley - No auxiliary data
256dc2eb70dSMatthew G. Knepley - No time-dependence
257dc2eb70dSMatthew G. Knepley */
258dc2eb70dSMatthew G. Knepley for (i = 0; i < its; ++i) {
259dc2eb70dSMatthew G. Knepley for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) {
260dc2eb70dSMatthew G. Knepley const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS;
261dc2eb70dSMatthew G. Knepley
2629566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elemVec, chunkSize * totDim));
263dc2eb70dSMatthew G. Knepley /* TODO Replace with DMPlexGetCellFields() */
264dc2eb70dSMatthew G. Knepley for (k = 0; k < chunkSize * totDim; ++k) u[k] = 1.0;
265dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) {
266dc2eb70dSMatthew G. Knepley PetscFormKey key;
267dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
268dc2eb70dSMatthew G. Knepley /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */
269dc2eb70dSMatthew G. Knepley
2709371c9d4SSatish Balay key.label = NULL;
2719371c9d4SSatish Balay key.value = 0;
2729371c9d4SSatish Balay key.field = f;
2739371c9d4SSatish Balay key.part = 0;
2749566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom));
2759566063dSJacob Faibussowitsch PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec));
276dc2eb70dSMatthew G. Knepley }
277dc2eb70dSMatthew G. Knepley }
278dc2eb70dSMatthew G. Knepley }
2799566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom));
2809566063dSJacob Faibussowitsch PetscCall(DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms));
2819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS));
2829566063dSJacob Faibussowitsch PetscCall(PetscFree2(u, elemVec));
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
284dc2eb70dSMatthew G. Knepley }
285dc2eb70dSMatthew G. Knepley
TestUnisolvence(DM dm)286d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestUnisolvence(DM dm)
287d71ae5a4SJacob Faibussowitsch {
288dc2eb70dSMatthew G. Knepley Mat M;
289dc2eb70dSMatthew G. Knepley Vec v;
290dc2eb70dSMatthew G. Knepley
291dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser;
2929566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v));
2939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v));
2949566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dm, dm, &M));
2959566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-mass_view"));
2969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M));
2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
298dc2eb70dSMatthew G. Knepley }
299dc2eb70dSMatthew G. Knepley
main(int argc,char ** argv)300d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
301d71ae5a4SJacob Faibussowitsch {
302dc2eb70dSMatthew G. Knepley DM dm;
303dc2eb70dSMatthew G. Knepley AppCtx ctx;
304dc2eb70dSMatthew G. Knepley PetscMPIInt size;
305dc2eb70dSMatthew G. Knepley
306327415f7SBarry Smith PetscFunctionBeginUser;
3079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
309bd158744SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only.");
3109566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3119566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
3129566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX));
3139566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
3159566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)dm, NULL, "-dm_view"));
3169566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx));
3179566063dSJacob Faibussowitsch PetscCall(TestEvaluation(dm));
3189566063dSJacob Faibussowitsch PetscCall(TestIntegration(dm, 1, ctx.its));
3199566063dSJacob Faibussowitsch PetscCall(TestUnisolvence(dm));
3209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
3219566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
322b122ec5aSJacob Faibussowitsch return 0;
323dc2eb70dSMatthew G. Knepley }
324dc2eb70dSMatthew G. Knepley
325dc2eb70dSMatthew G. Knepley /*TEST
326dc2eb70dSMatthew G. Knepley test:
327dc2eb70dSMatthew G. Knepley suffix: 0
328dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0
329dc2eb70dSMatthew G. Knepley
330dc2eb70dSMatthew G. Knepley test:
331dc2eb70dSMatthew G. Knepley suffix: 2
332dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
333dc2eb70dSMatthew G. Knepley -field_petscspace_type sum \
334dc2eb70dSMatthew G. Knepley -field_petscspace_variables 3 \
335dc2eb70dSMatthew G. Knepley -field_petscspace_components 3 \
336dc2eb70dSMatthew G. Knepley -field_petscspace_sum_spaces 2 \
337dc2eb70dSMatthew G. Knepley -field_petscspace_sum_concatenate false \
338dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_variables 3 \
339dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_components 3 \
340dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_degree 1 \
341dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_variables 3 \
342dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_components 3 \
343dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_type wxy \
344dc2eb70dSMatthew G. Knepley -field_petscdualspace_form_degree 0 \
345dc2eb70dSMatthew G. Knepley -field_petscdualspace_order 1 \
346dc2eb70dSMatthew G. Knepley -field_petscdualspace_components 3
347dc2eb70dSMatthew G. Knepley
348dc2eb70dSMatthew G. Knepley TEST*/
349