xref: /petsc/src/dm/field/tutorials/ex1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Demonstration of creating and viewing DMFields objects.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmfield.h>
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown 
7*9371c9d4SSatish Balay static PetscErrorCode ViewResults(PetscViewer viewer, PetscInt N, PetscInt dim, PetscScalar *B, PetscScalar *D, PetscScalar *H, PetscReal *rB, PetscReal *rD, PetscReal *rH) {
8c4762a1bSJed Brown   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "B:\n"));
109566063dSJacob Faibussowitsch   PetscCall(PetscScalarView(N, B, viewer));
119566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "D:\n"));
129566063dSJacob Faibussowitsch   PetscCall(PetscScalarView(N * dim, D, viewer));
139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "H:\n"));
149566063dSJacob Faibussowitsch   PetscCall(PetscScalarView(N * dim * dim, H, viewer));
15c4762a1bSJed Brown 
169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "rB:\n"));
179566063dSJacob Faibussowitsch   PetscCall(PetscRealView(N, rB, viewer));
189566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "rD:\n"));
199566063dSJacob Faibussowitsch   PetscCall(PetscRealView(N * dim, rD, viewer));
209566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "rH:\n"));
219566063dSJacob Faibussowitsch   PetscCall(PetscRealView(N * dim * dim, rH, viewer));
22c4762a1bSJed Brown   PetscFunctionReturn(0);
23c4762a1bSJed Brown }
24c4762a1bSJed Brown 
25*9371c9d4SSatish Balay static PetscErrorCode TestEvaluate(DMField field, PetscInt n, PetscRandom rand) {
26c4762a1bSJed Brown   DM           dm;
27c4762a1bSJed Brown   PetscInt     dim, i, nc;
28c4762a1bSJed Brown   PetscScalar *B, *D, *H;
29c4762a1bSJed Brown   PetscReal   *rB, *rD, *rH;
30c4762a1bSJed Brown   Vec          points;
31c4762a1bSJed Brown   PetscScalar *array;
32c4762a1bSJed Brown   PetscViewer  viewer;
33c4762a1bSJed Brown   MPI_Comm     comm;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBegin;
36c4762a1bSJed Brown   comm = PetscObjectComm((PetscObject)field);
379566063dSJacob Faibussowitsch   PetscCall(DMFieldGetNumComponents(field, &nc));
389566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDM(field, &dm));
399566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
409566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)field), n * dim, PETSC_DETERMINE, &points));
419566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(points, dim));
429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(points, &array));
439566063dSJacob Faibussowitsch   for (i = 0; i < n * dim; i++) PetscCall(PetscRandomGetValue(rand, &array[i]));
449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(points, &array));
459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc6(n * nc, &B, n * nc, &rB, n * nc * dim, &D, n * nc * dim, &rD, n * nc * dim * dim, &H, n * nc * dim * dim, &rH));
469566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluate(field, points, PETSC_SCALAR, B, D, H));
479566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluate(field, points, PETSC_REAL, rB, rD, rH));
48c4762a1bSJed Brown   viewer = PETSC_VIEWER_STDOUT_(comm);
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)points, "Test Points"));
519566063dSJacob Faibussowitsch   PetscCall(VecView(points, viewer));
529566063dSJacob Faibussowitsch   PetscCall(ViewResults(viewer, n * nc, dim, B, D, H, rB, rD, rH));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PetscFree6(B, rB, D, rD, H, rH));
559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&points));
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59*9371c9d4SSatish Balay static PetscErrorCode TestEvaluateFE(DMField field, PetscInt n, PetscInt cStart, PetscInt cEnd, PetscQuadrature quad, PetscRandom rand) {
60c4762a1bSJed Brown   DM           dm;
61c4762a1bSJed Brown   PetscInt     dim, i, nc, nq;
62c4762a1bSJed Brown   PetscInt     N;
63c4762a1bSJed Brown   PetscScalar *B, *D, *H;
64c4762a1bSJed Brown   PetscReal   *rB, *rD, *rH;
65c4762a1bSJed Brown   PetscInt    *cells;
66c4762a1bSJed Brown   IS           cellIS;
67c4762a1bSJed Brown   PetscViewer  viewer;
68c4762a1bSJed Brown   MPI_Comm     comm;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBegin;
71c4762a1bSJed Brown   comm = PetscObjectComm((PetscObject)field);
729566063dSJacob Faibussowitsch   PetscCall(DMFieldGetNumComponents(field, &nc));
739566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDM(field, &dm));
749566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
759566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand, (PetscScalar)cStart, (PetscScalar)cEnd));
769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &cells));
77c4762a1bSJed Brown   for (i = 0; i < n; i++) {
78c4762a1bSJed Brown     PetscReal rc;
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(rand, &rc));
81c4762a1bSJed Brown     cells[i] = PetscFloorReal(rc);
82c4762a1bSJed Brown   }
839566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, n, cells, PETSC_OWN_POINTER, &cellIS));
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)cellIS, "FE Test Cells"));
859566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &nq, NULL, NULL));
86c4762a1bSJed Brown   N = n * nq * nc;
879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc6(N, &B, N, &rB, N * dim, &D, N * dim, &rD, N * dim * dim, &H, N * dim * dim, &rH));
889566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFE(field, cellIS, quad, PETSC_SCALAR, B, D, H));
899566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFE(field, cellIS, quad, PETSC_REAL, rB, rD, rH));
90c4762a1bSJed Brown   viewer = PETSC_VIEWER_STDOUT_(comm);
91c4762a1bSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)quad, "Test quadrature"));
939566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureView(quad, viewer));
949566063dSJacob Faibussowitsch   PetscCall(ISView(cellIS, viewer));
959566063dSJacob Faibussowitsch   PetscCall(ViewResults(viewer, N, dim, B, D, H, rB, rD, rH));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(PetscFree6(B, rB, D, rD, H, rH));
989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102*9371c9d4SSatish Balay static PetscErrorCode TestEvaluateFV(DMField field, PetscInt n, PetscInt cStart, PetscInt cEnd, PetscRandom rand) {
103c4762a1bSJed Brown   DM           dm;
104c4762a1bSJed Brown   PetscInt     dim, i, nc;
105c4762a1bSJed Brown   PetscInt     N;
106c4762a1bSJed Brown   PetscScalar *B, *D, *H;
107c4762a1bSJed Brown   PetscReal   *rB, *rD, *rH;
108c4762a1bSJed Brown   PetscInt    *cells;
109c4762a1bSJed Brown   IS           cellIS;
110c4762a1bSJed Brown   PetscViewer  viewer;
111c4762a1bSJed Brown   MPI_Comm     comm;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBegin;
114c4762a1bSJed Brown   comm = PetscObjectComm((PetscObject)field);
1159566063dSJacob Faibussowitsch   PetscCall(DMFieldGetNumComponents(field, &nc));
1169566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDM(field, &dm));
1179566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1189566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand, (PetscScalar)cStart, (PetscScalar)cEnd));
1199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &cells));
120c4762a1bSJed Brown   for (i = 0; i < n; i++) {
121c4762a1bSJed Brown     PetscReal rc;
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(rand, &rc));
124c4762a1bSJed Brown     cells[i] = PetscFloorReal(rc);
125c4762a1bSJed Brown   }
1269566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, n, cells, PETSC_OWN_POINTER, &cellIS));
1279566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)cellIS, "FV Test Cells"));
128c4762a1bSJed Brown   N = n * nc;
1299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc6(N, &B, N, &rB, N * dim, &D, N * dim, &rD, N * dim * dim, &H, N * dim * dim, &rH));
1309566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFV(field, cellIS, PETSC_SCALAR, B, D, H));
1319566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFV(field, cellIS, PETSC_REAL, rB, rD, rH));
132c4762a1bSJed Brown   viewer = PETSC_VIEWER_STDOUT_(comm);
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(ISView(cellIS, viewer));
1359566063dSJacob Faibussowitsch   PetscCall(ViewResults(viewer, N, dim, B, D, H, rB, rD, rH));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCall(PetscFree6(B, rB, D, rD, H, rH));
1389566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
139c4762a1bSJed Brown   PetscFunctionReturn(0);
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142*9371c9d4SSatish Balay static PetscErrorCode radiusSquared(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) {
143c4762a1bSJed Brown   PetscInt  i;
144c4762a1bSJed Brown   PetscReal r2 = 0.;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBegin;
147c4762a1bSJed Brown   for (i = 0; i < dim; i++) { r2 += PetscSqr(x[i]); }
148*9371c9d4SSatish Balay   for (i = 0; i < Nf; i++) { u[i] = (i + 1) * r2; }
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152*9371c9d4SSatish Balay static PetscErrorCode TestShellEvaluate(DMField field, Vec points, PetscDataType type, void *B, void *D, void *H) {
153c4762a1bSJed Brown   Vec                ctxVec = NULL;
154c4762a1bSJed Brown   const PetscScalar *mult;
155c4762a1bSJed Brown   PetscInt           dim;
156c4762a1bSJed Brown   const PetscScalar *x;
157c4762a1bSJed Brown   PetscInt           Nc, n, i, j, k, l;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBegin;
1609566063dSJacob Faibussowitsch   PetscCall(DMFieldGetNumComponents(field, &Nc));
1619566063dSJacob Faibussowitsch   PetscCall(DMFieldShellGetContext(field, &ctxVec));
1629566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(points, &dim));
1639566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(points, &n));
164c4762a1bSJed Brown   n /= Nc;
1659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctxVec, &mult));
1669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(points, &x));
167c4762a1bSJed Brown   for (i = 0; i < n; i++) {
168c4762a1bSJed Brown     PetscReal r2 = 0.;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown     for (j = 0; j < dim; j++) { r2 += PetscSqr(PetscRealPart(x[i * dim + j])); }
171c4762a1bSJed Brown     for (j = 0; j < Nc; j++) {
172c4762a1bSJed Brown       PetscReal m = PetscRealPart(mult[j]);
173c4762a1bSJed Brown       if (B) {
174c4762a1bSJed Brown         if (type == PETSC_SCALAR) {
175c4762a1bSJed Brown           ((PetscScalar *)B)[i * Nc + j] = m * r2;
176c4762a1bSJed Brown         } else {
177c4762a1bSJed Brown           ((PetscReal *)B)[i * Nc + j] = m * r2;
178c4762a1bSJed Brown         }
179c4762a1bSJed Brown       }
180c4762a1bSJed Brown       if (D) {
181c4762a1bSJed Brown         if (type == PETSC_SCALAR) {
182c4762a1bSJed Brown           for (k = 0; k < dim; k++) ((PetscScalar *)D)[(i * Nc + j) * dim + k] = 2. * m * x[i * dim + k];
183c4762a1bSJed Brown         } else {
184c4762a1bSJed Brown           for (k = 0; k < dim; k++) ((PetscReal *)D)[(i * Nc + j) * dim + k] = 2. * m * PetscRealPart(x[i * dim + k]);
185c4762a1bSJed Brown         }
186c4762a1bSJed Brown       }
187c4762a1bSJed Brown       if (H) {
188c4762a1bSJed Brown         if (type == PETSC_SCALAR) {
189*9371c9d4SSatish Balay           for (k = 0; k < dim; k++)
190*9371c9d4SSatish Balay             for (l = 0; l < dim; l++) ((PetscScalar *)H)[((i * Nc + j) * dim + k) * dim + l] = (k == l) ? 2. * m : 0.;
191c4762a1bSJed Brown         } else {
192*9371c9d4SSatish Balay           for (k = 0; k < dim; k++)
193*9371c9d4SSatish Balay             for (l = 0; l < dim; l++) ((PetscReal *)H)[((i * Nc + j) * dim + k) * dim + l] = (k == l) ? 2. * m : 0.;
194c4762a1bSJed Brown         }
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown     }
197c4762a1bSJed Brown   }
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(points, &x));
1999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctxVec, &mult));
200c4762a1bSJed Brown   PetscFunctionReturn(0);
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203*9371c9d4SSatish Balay static PetscErrorCode TestShellDestroy(DMField field) {
204c4762a1bSJed Brown   Vec ctxVec = NULL;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   PetscFunctionBegin;
2079566063dSJacob Faibussowitsch   PetscCall(DMFieldShellGetContext(field, &ctxVec));
2089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctxVec));
209c4762a1bSJed Brown   PetscFunctionReturn(0);
210c4762a1bSJed Brown }
211c4762a1bSJed Brown 
212*9371c9d4SSatish Balay int main(int argc, char **argv) {
213c4762a1bSJed Brown   DM              dm = NULL;
214c4762a1bSJed Brown   MPI_Comm        comm;
215c4762a1bSJed Brown   char            type[256] = DMPLEX;
216c4762a1bSJed Brown   PetscBool       isda, isplex;
217c4762a1bSJed Brown   PetscInt        dim    = 2;
218c4762a1bSJed Brown   DMField         field  = NULL;
219c4762a1bSJed Brown   PetscInt        nc     = 1;
220c4762a1bSJed Brown   PetscInt        cStart = -1, cEnd = -1;
221c4762a1bSJed Brown   PetscRandom     rand;
222c4762a1bSJed Brown   PetscQuadrature quad          = NULL;
223c4762a1bSJed Brown   PetscInt        pointsPerEdge = 2;
224c4762a1bSJed Brown   PetscInt        numPoint = 0, numFE = 0, numFV = 0;
225c4762a1bSJed Brown   PetscBool       testShell = PETSC_FALSE;
226c4762a1bSJed Brown 
227327415f7SBarry Smith   PetscFunctionBeginUser;
2289566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
229c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
230d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "DMField Tutorial Options", "DM");
2319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-dm_type", "DM implementation on which to define field", "ex1.c", DMList, type, type, 256, NULL));
2329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "DM intrinsic dimension", "ex1.c", dim, &dim, NULL, 1, 3));
2339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_components", "Number of components in field", "ex1.c", nc, &nc, NULL, 1));
2349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_quad_points", "Number of quadrature points per dimension", "ex1.c", pointsPerEdge, &pointsPerEdge, NULL, 1));
2359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_point_tests", "Number of test points for DMFieldEvaluate()", "ex1.c", numPoint, &numPoint, NULL, 0));
2369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fe_tests", "Number of test cells for DMFieldEvaluateFE()", "ex1.c", numFE, &numFE, NULL, 0));
2379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fv_tests", "Number of test cells for DMFieldEvaluateFV()", "ex1.c", numFV, &numFV, NULL, 0));
2389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_shell", "Test the DMFIELDSHELL implementation of DMField", "ex1.c", testShell, &testShell, NULL));
239d0609cedSBarry Smith   PetscOptionsEnd();
240c4762a1bSJed Brown 
24163a3b9bcSJacob Faibussowitsch   PetscCheck(dim <= 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "This examples works for dim <= 3, not %" PetscInt_FMT, dim);
2429566063dSJacob Faibussowitsch   PetscCall(PetscStrncmp(type, DMPLEX, 256, &isplex));
2439566063dSJacob Faibussowitsch   PetscCall(PetscStrncmp(type, DMDA, 256, &isda));
244c4762a1bSJed Brown 
2459566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
2469566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
247c4762a1bSJed Brown   if (isplex) {
248c4762a1bSJed Brown     PetscInt  overlap = 0;
249c4762a1bSJed Brown     Vec       fieldvec;
250c4762a1bSJed Brown     PetscInt  cells[3] = {3, 3, 3};
25130602db0SMatthew G. Knepley     PetscBool simplex;
252c4762a1bSJed Brown     PetscFE   fe;
253c4762a1bSJed Brown 
254d0609cedSBarry Smith     PetscOptionsBegin(comm, "", "DMField DMPlex Options", "DM");
2559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBoundedInt("-overlap", "DMPlex parallel overlap", "ex1.c", overlap, &overlap, NULL, 0));
256d0609cedSBarry Smith     PetscOptionsEnd();
25730602db0SMatthew G. Knepley     if (0) {
2589566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm));
25930602db0SMatthew G. Knepley     } else {
2609566063dSJacob Faibussowitsch       PetscCall(DMCreate(comm, &dm));
2619566063dSJacob Faibussowitsch       PetscCall(DMSetType(dm, DMPLEX));
2629566063dSJacob Faibussowitsch       PetscCall(DMSetFromOptions(dm));
26330602db0SMatthew G. Knepley       CHKMEMQ;
264c4762a1bSJed Brown     }
2659566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
2669566063dSJacob Faibussowitsch     PetscCall(DMPlexIsSimplex(dm, &simplex));
2671baa6e33SBarry Smith     if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad));
2681baa6e33SBarry Smith     else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad));
2699566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
270c4762a1bSJed Brown     if (testShell) {
271c4762a1bSJed Brown       Vec          ctxVec;
272c4762a1bSJed Brown       PetscInt     i;
273c4762a1bSJed Brown       PetscScalar *array;
274c4762a1bSJed Brown 
2759566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, nc, &ctxVec));
2769566063dSJacob Faibussowitsch       PetscCall(VecSetUp(ctxVec));
2779566063dSJacob Faibussowitsch       PetscCall(VecGetArray(ctxVec, &array));
278c4762a1bSJed Brown       for (i = 0; i < nc; i++) array[i] = i + 1.;
2799566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(ctxVec, &array));
2809566063dSJacob Faibussowitsch       PetscCall(DMFieldCreateShell(dm, nc, DMFIELD_VERTEX, (void *)ctxVec, &field));
2819566063dSJacob Faibussowitsch       PetscCall(DMFieldShellSetEvaluate(field, TestShellEvaluate));
2829566063dSJacob Faibussowitsch       PetscCall(DMFieldShellSetDestroy(field, TestShellDestroy));
283c4762a1bSJed Brown     } else {
2849566063dSJacob Faibussowitsch       PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, nc, simplex, NULL, PETSC_DEFAULT, &fe));
2859566063dSJacob Faibussowitsch       PetscCall(PetscFESetName(fe, "MyPetscFE"));
2869566063dSJacob Faibussowitsch       PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
2879566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&fe));
2889566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(dm));
2899566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(dm, &fieldvec));
290c4762a1bSJed Brown       {
291c4762a1bSJed Brown         PetscErrorCode (*func[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
292c4762a1bSJed Brown         void *ctxs[1];
293c4762a1bSJed Brown 
294c4762a1bSJed Brown         func[0] = radiusSquared;
295c4762a1bSJed Brown         ctxs[0] = NULL;
296c4762a1bSJed Brown 
2979566063dSJacob Faibussowitsch         PetscCall(DMProjectFunctionLocal(dm, 0.0, func, ctxs, INSERT_ALL_VALUES, fieldvec));
298c4762a1bSJed Brown       }
2999566063dSJacob Faibussowitsch       PetscCall(DMFieldCreateDS(dm, 0, fieldvec, &field));
3009566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&fieldvec));
301c4762a1bSJed Brown     }
302c4762a1bSJed Brown   } else if (isda) {
303c4762a1bSJed Brown     PetscInt     i;
304c4762a1bSJed Brown     PetscScalar *cv;
305c4762a1bSJed Brown 
306c4762a1bSJed Brown     switch (dim) {
307*9371c9d4SSatish Balay     case 1: PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, 3, 1, 1, NULL, &dm)); break;
308*9371c9d4SSatish Balay     case 2: PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 3, 3, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, &dm)); break;
309*9371c9d4SSatish Balay     default: PetscCall(DMDACreate3d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 3, 3, 3, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, NULL, &dm)); break;
310c4762a1bSJed Brown     }
3119566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
3129566063dSJacob Faibussowitsch     PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd));
3139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nc * (1 << dim), &cv));
314c4762a1bSJed Brown     for (i = 0; i < nc * (1 << dim); i++) {
315c4762a1bSJed Brown       PetscReal rv;
316c4762a1bSJed Brown 
3179566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &rv));
318c4762a1bSJed Brown       cv[i] = rv;
319c4762a1bSJed Brown     }
3209566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateDA(dm, nc, cv, &field));
3219566063dSJacob Faibussowitsch     PetscCall(PetscFree(cv));
3229566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad));
32398921bdaSJacob Faibussowitsch   } else SETERRQ(comm, PETSC_ERR_SUP, "This test does not run for DM type %s", type);
324c4762a1bSJed Brown 
3259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "mesh"));
3269566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
3279566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3289566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)field, "field"));
3299566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)field, NULL, "-dmfield_view"));
3309566063dSJacob Faibussowitsch   if (numPoint) PetscCall(TestEvaluate(field, numPoint, rand));
3319566063dSJacob Faibussowitsch   if (numFE) PetscCall(TestEvaluateFE(field, numFE, cStart, cEnd, quad, rand));
3329566063dSJacob Faibussowitsch   if (numFV) PetscCall(TestEvaluateFV(field, numFV, cStart, cEnd, rand));
3339566063dSJacob Faibussowitsch   PetscCall(DMFieldDestroy(&field));
3349566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quad));
3359566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
3369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
337b122ec5aSJacob Faibussowitsch   return 0;
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown /*TEST
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   test:
343c4762a1bSJed Brown     suffix: da
344c4762a1bSJed Brown     requires: !complex
345c4762a1bSJed Brown     args: -dm_type da -dim 2 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   test:
348c4762a1bSJed Brown     suffix: da_1
349c4762a1bSJed Brown     requires: !complex
350c4762a1bSJed Brown     args: -dm_type da -dim 1  -num_fe_tests 2
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   test:
353c4762a1bSJed Brown     suffix: da_2
354c4762a1bSJed Brown     requires: !complex
355c4762a1bSJed Brown     args: -dm_type da -dim 2  -num_fe_tests 2
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   test:
358c4762a1bSJed Brown     suffix: da_3
359c4762a1bSJed Brown     requires: !complex
360c4762a1bSJed Brown     args: -dm_type da -dim 3  -num_fe_tests 2
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   test:
363c4762a1bSJed Brown     suffix: ds
364c4762a1bSJed Brown     requires: !complex triangle
36530602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view -petscspace_degree 2 -num_quad_points 1
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   test:
368c4762a1bSJed Brown     suffix: ds_simplex_0
369c4762a1bSJed Brown     requires: !complex triangle
37030602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3  -num_fe_tests 2  -petscspace_degree 0
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   test:
373c4762a1bSJed Brown     suffix: ds_simplex_1
374c4762a1bSJed Brown     requires: !complex triangle
37530602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3  -num_fe_tests 2  -petscspace_degree 1
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   test:
378c4762a1bSJed Brown     suffix: ds_simplex_2
379c4762a1bSJed Brown     requires: !complex triangle
38030602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3  -num_fe_tests 2  -petscspace_degree 2
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   test:
383c4762a1bSJed Brown     suffix: ds_tensor_2_0
384c4762a1bSJed Brown     requires: !complex
38530602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3  -num_fe_tests 2  -petscspace_poly_tensor 1 -petscspace_degree 0 -dm_plex_simplex 0
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   test:
388c4762a1bSJed Brown     suffix: ds_tensor_2_1
389c4762a1bSJed Brown     requires: !complex
39030602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3  -num_fe_tests 2  -petscspace_poly_tensor 1 -petscspace_degree 1 -dm_plex_simplex 0
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   test:
393c4762a1bSJed Brown     suffix: ds_tensor_2_2
394c4762a1bSJed Brown     requires: !complex
39530602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3  -num_fe_tests 2  -petscspace_poly_tensor 1 -petscspace_degree 2 -dm_plex_simplex 0
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   test:
398c4762a1bSJed Brown     suffix: ds_tensor_3_0
399c4762a1bSJed Brown     requires: !complex
40030602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 3 -dm_plex_box_faces 3,3,3  -num_fe_tests 2  -petscspace_poly_tensor 1 -petscspace_degree 0 -dm_plex_simplex 0
401c4762a1bSJed Brown 
402c4762a1bSJed Brown   test:
403c4762a1bSJed Brown     suffix: ds_tensor_3_1
404c4762a1bSJed Brown     requires: !complex
40530602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 3 -dm_plex_box_faces 3,3,3  -num_fe_tests 2  -petscspace_poly_tensor 1 -petscspace_degree 1 -dm_plex_simplex 0
406c4762a1bSJed Brown 
407c4762a1bSJed Brown   test:
408c4762a1bSJed Brown     suffix: ds_tensor_3_2
409c4762a1bSJed Brown     requires: !complex
41030602db0SMatthew G. Knepley     args: -dm_type plex -dm_plex_dim 3 -dm_plex_box_faces 3,3,3  -num_fe_tests 2  -petscspace_poly_tensor 1 -petscspace_degree 2 -dm_plex_simplex 0
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   test:
413c4762a1bSJed Brown     suffix: shell
414c4762a1bSJed Brown     requires: !complex triangle
41530602db0SMatthew G. Knepley     args: -dm_coord_space 0 -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view -num_quad_points 1 -test_shell
416c4762a1bSJed Brown 
417c4762a1bSJed Brown TEST*/
418