1 static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n";
2
3 #include <petscdmplex.h>
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 #include <petscfe.h>
7 #include <petscds.h>
8 #include <petscksp.h>
9 #include <petscsnes.h>
10 #include <petscsf.h>
11
12 typedef struct {
13 /* Domain and mesh definition */
14 PetscBool useDA; /* Flag DMDA tensor product mesh */
15 PetscBool shearCoords; /* Flag for shear transform */
16 PetscBool nonaffineCoords; /* Flag for non-affine transform */
17 /* Element definition */
18 PetscInt qorder; /* Order of the quadrature */
19 PetscInt numComponents; /* Number of field components */
20 PetscFE fe; /* The finite element */
21 /* Testing space */
22 PetscInt porder; /* Order of polynomials to test */
23 PetscBool RT; /* Test for Raviart-Thomas elements */
24 PetscBool convergence; /* Test for order of convergence */
25 PetscBool convRefine; /* Test for convergence using refinement, otherwise use coarsening */
26 PetscBool constraints; /* Test local constraints */
27 PetscBool tree; /* Test tree routines */
28 PetscBool testFEjacobian; /* Test finite element Jacobian assembly */
29 PetscBool testFVgrad; /* Test finite difference gradient routine */
30 PetscBool testInjector; /* Test finite element injection routines */
31 PetscInt treeCell; /* Cell to refine in tree test */
32 PetscReal constants[3]; /* Constant values for each dimension */
33 } AppCtx;
34
35 /*
36 Derivatives are set as n_i \partial u_j / \partial x_i
37 */
38
39 /* u = 1 */
constant(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)40 PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
41 {
42 AppCtx *user = (AppCtx *)ctx;
43 PetscInt d;
44 for (d = 0; d < dim; ++d) u[d] = user->constants[d];
45 return PETSC_SUCCESS;
46 }
constantDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)47 PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
48 {
49 PetscInt d;
50 for (d = 0; d < dim; ++d) u[d] = 0.0;
51 return PETSC_SUCCESS;
52 }
53
54 /* RT_0: u = (1 + x, 1 + y) or (1 + x, 1 + y, 1 + z) */
rt0(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)55 PetscErrorCode rt0(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
56 {
57 PetscInt d;
58 for (d = 0; d < dim; ++d) u[d] = 1.0 + coords[d];
59 return PETSC_SUCCESS;
60 }
61
rt0Der(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)62 PetscErrorCode rt0Der(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
63 {
64 PetscInt d, e;
65 for (d = 0; d < dim; ++d) {
66 u[d] = 0.0;
67 for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
68 }
69 return PETSC_SUCCESS;
70 }
71
72 /* u = (x + y, y + x) or (x + z, 2y, z + x) */
linear(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)73 PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
74 {
75 PetscInt d;
76 for (d = 0; d < dim; ++d) u[d] = coords[d] + coords[dim - d - 1];
77 return PETSC_SUCCESS;
78 }
linearDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)79 PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
80 {
81 PetscInt d, e;
82 for (d = 0; d < dim; ++d) {
83 u[d] = 0.0;
84 for (e = 0; e < dim; ++e) u[d] += ((d == e ? 1. : 0.) + (d == (dim - e - 1) ? 1. : 0.)) * n[e];
85 }
86 return PETSC_SUCCESS;
87 }
88
89 /* RT_1: u = (1 + x + y + x^2 + xy, 1 + x + y + xy + y^2) or (1 + x + y + z + x^2 + xy + xz, 1 + x + y + z + xy + y^2 + yz, 1 + x + y + z + xz + yz + z^2) */
rt1(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)90 PetscErrorCode rt1(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
91 {
92 if (dim > 2) {
93 u[0] = 1.0 + coords[0] + coords[1] + coords[2] + coords[0] * coords[0] + coords[0] * coords[1] + coords[0] * coords[2];
94 u[1] = 1.0 + coords[0] + coords[1] + coords[2] + coords[0] * coords[1] + coords[1] * coords[1] + coords[1] * coords[2];
95 u[2] = 1.0 + coords[0] + coords[1] + coords[2] + coords[0] * coords[2] + coords[1] * coords[2] + coords[2] * coords[2];
96 } else if (dim > 1) {
97 u[0] = 1.0 + coords[0] + coords[1] + coords[0] * coords[0] + coords[0] * coords[1];
98 u[1] = 1.0 + coords[0] + coords[1] + coords[0] * coords[1] + coords[1] * coords[1];
99 }
100 return PETSC_SUCCESS;
101 }
102
rt1Der(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)103 PetscErrorCode rt1Der(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
104 {
105 if (dim > 2) {
106 u[0] = (1.0 + 2.0 * coords[0] + coords[1] + coords[2]) * n[0] + (1.0 + coords[0]) * n[1] + (1.0 + coords[0]) * n[2];
107 u[1] = (1.0 + coords[1]) * n[0] + (1.0 + coords[0] + 2.0 * coords[1] + coords[2]) * n[1] + (1.0 + coords[1]) * n[2];
108 u[2] = (1.0 + coords[2]) * n[0] + (1.0 + coords[2]) * n[1] + (1.0 + coords[0] + coords[1] + 2.0 * coords[2]) * n[2];
109 } else if (dim > 1) {
110 u[0] = (1.0 + 2.0 * coords[0] + coords[1]) * n[0] + (1.0 + coords[0]) * n[1];
111 u[1] = (1.0 + coords[1]) * n[0] + (1.0 + coords[0] + 2.0 * coords[1]) * n[1];
112 }
113 return PETSC_SUCCESS;
114 }
115
116 /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
quadratic(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)117 PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
118 {
119 if (dim > 2) {
120 u[0] = coords[0] * coords[1];
121 u[1] = coords[1] * coords[2];
122 u[2] = coords[2] * coords[0];
123 } else if (dim > 1) {
124 u[0] = coords[0] * coords[0];
125 u[1] = coords[0] * coords[1];
126 } else if (dim > 0) {
127 u[0] = coords[0] * coords[0];
128 }
129 return PETSC_SUCCESS;
130 }
quadraticDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)131 PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
132 {
133 if (dim > 2) {
134 u[0] = coords[1] * n[0] + coords[0] * n[1];
135 u[1] = coords[2] * n[1] + coords[1] * n[2];
136 u[2] = coords[2] * n[0] + coords[0] * n[2];
137 } else if (dim > 1) {
138 u[0] = 2.0 * coords[0] * n[0];
139 u[1] = coords[1] * n[0] + coords[0] * n[1];
140 } else if (dim > 0) {
141 u[0] = 2.0 * coords[0] * n[0];
142 }
143 return PETSC_SUCCESS;
144 }
145
146 /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
cubic(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)147 PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
148 {
149 if (dim > 2) {
150 u[0] = coords[0] * coords[0] * coords[1];
151 u[1] = coords[1] * coords[1] * coords[2];
152 u[2] = coords[2] * coords[2] * coords[0];
153 } else if (dim > 1) {
154 u[0] = coords[0] * coords[0] * coords[0];
155 u[1] = coords[0] * coords[0] * coords[1];
156 } else if (dim > 0) {
157 u[0] = coords[0] * coords[0] * coords[0];
158 }
159 return PETSC_SUCCESS;
160 }
cubicDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)161 PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
162 {
163 if (dim > 2) {
164 u[0] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
165 u[1] = 2.0 * coords[1] * coords[2] * n[1] + coords[1] * coords[1] * n[2];
166 u[2] = 2.0 * coords[2] * coords[0] * n[2] + coords[2] * coords[2] * n[0];
167 } else if (dim > 1) {
168 u[0] = 3.0 * coords[0] * coords[0] * n[0];
169 u[1] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
170 } else if (dim > 0) {
171 u[0] = 3.0 * coords[0] * coords[0] * n[0];
172 }
173 return PETSC_SUCCESS;
174 }
175
176 /* u = tanh(x) */
trig(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)177 PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
178 {
179 PetscInt d;
180 for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
181 return PETSC_SUCCESS;
182 }
trigDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)183 PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
184 {
185 PetscInt d;
186 for (d = 0; d < dim; ++d) u[d] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
187 return PETSC_SUCCESS;
188 }
189
ProcessOptions(MPI_Comm comm,AppCtx * options)190 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
191 {
192 PetscInt n = 3;
193
194 PetscFunctionBeginUser;
195 options->useDA = PETSC_FALSE;
196 options->shearCoords = PETSC_FALSE;
197 options->nonaffineCoords = PETSC_FALSE;
198 options->qorder = 0;
199 options->numComponents = PETSC_DEFAULT;
200 options->porder = 0;
201 options->RT = PETSC_FALSE;
202 options->convergence = PETSC_FALSE;
203 options->convRefine = PETSC_TRUE;
204 options->constraints = PETSC_FALSE;
205 options->tree = PETSC_FALSE;
206 options->treeCell = 0;
207 options->testFEjacobian = PETSC_FALSE;
208 options->testFVgrad = PETSC_FALSE;
209 options->testInjector = PETSC_FALSE;
210 options->constants[0] = 1.0;
211 options->constants[1] = 2.0;
212 options->constants[2] = 3.0;
213
214 PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");
215 PetscCall(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL));
216 PetscCall(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL));
217 PetscCall(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL));
218 PetscCall(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL, 0));
219 PetscCall(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL, PETSC_DEFAULT));
220 PetscCall(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL, 0));
221 PetscCall(PetscOptionsBool("-RT", "Use the Raviart-Thomas elements", "ex3.c", options->RT, &options->RT, NULL));
222 PetscCall(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL));
223 PetscCall(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL));
224 PetscCall(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL));
225 PetscCall(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL));
226 PetscCall(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL, 0));
227 PetscCall(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL));
228 PetscCall(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL));
229 PetscCall(PetscOptionsBool("-test_injector", "Test finite element injection", "ex3.c", options->testInjector, &options->testInjector, NULL));
230 PetscCall(PetscOptionsRealArray("-constants", "Set the constant values", "ex3.c", options->constants, &n, NULL));
231 PetscOptionsEnd();
232 PetscFunctionReturn(PETSC_SUCCESS);
233 }
234
TransformCoordinates(DM dm,AppCtx * user)235 static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user)
236 {
237 PetscSection coordSection;
238 Vec coordinates;
239 PetscScalar *coords;
240 PetscInt vStart, vEnd, v;
241
242 PetscFunctionBeginUser;
243 if (user->nonaffineCoords) {
244 /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */
245 PetscCall(DMGetCoordinateSection(dm, &coordSection));
246 PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
247 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
248 PetscCall(VecGetArray(coordinates, &coords));
249 for (v = vStart; v < vEnd; ++v) {
250 PetscInt dof, off;
251 PetscReal p = 4.0, r;
252
253 PetscCall(PetscSectionGetDof(coordSection, v, &dof));
254 PetscCall(PetscSectionGetOffset(coordSection, v, &off));
255 switch (dof) {
256 case 2:
257 r = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
258 coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
259 coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
260 break;
261 case 3:
262 r = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
263 coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
264 coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
265 coords[off + 2] = coords[off + 2];
266 break;
267 }
268 }
269 PetscCall(VecRestoreArray(coordinates, &coords));
270 }
271 if (user->shearCoords) {
272 /* x' = x + m y + m z, y' = y + m z, z' = z */
273 PetscCall(DMGetCoordinateSection(dm, &coordSection));
274 PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
275 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
276 PetscCall(VecGetArray(coordinates, &coords));
277 for (v = vStart; v < vEnd; ++v) {
278 PetscInt dof, off;
279 PetscReal m = 1.0;
280
281 PetscCall(PetscSectionGetDof(coordSection, v, &dof));
282 PetscCall(PetscSectionGetOffset(coordSection, v, &off));
283 switch (dof) {
284 case 2:
285 coords[off + 0] = coords[off + 0] + m * coords[off + 1];
286 coords[off + 1] = coords[off + 1];
287 break;
288 case 3:
289 coords[off + 0] = coords[off + 0] + m * coords[off + 1] + m * coords[off + 2];
290 coords[off + 1] = coords[off + 1] + m * coords[off + 2];
291 coords[off + 2] = coords[off + 2];
292 break;
293 }
294 }
295 PetscCall(VecRestoreArray(coordinates, &coords));
296 }
297 PetscFunctionReturn(PETSC_SUCCESS);
298 }
299
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)300 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
301 {
302 PetscInt dim = 2;
303 PetscBool simplex;
304
305 PetscFunctionBeginUser;
306 if (user->useDA) {
307 switch (dim) {
308 case 2:
309 PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm));
310 PetscCall(DMSetFromOptions(*dm));
311 PetscCall(DMSetUp(*dm));
312 PetscCall(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
313 break;
314 default:
315 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %" PetscInt_FMT, dim);
316 }
317 PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
318 } else {
319 PetscCall(DMCreate(comm, dm));
320 PetscCall(DMSetType(*dm, DMPLEX));
321 PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
322 PetscCall(DMSetFromOptions(*dm));
323
324 PetscCall(DMGetDimension(*dm, &dim));
325 PetscCall(DMPlexIsSimplex(*dm, &simplex));
326 PetscCallMPI(MPI_Bcast(&simplex, 1, MPI_C_BOOL, 0, comm));
327 if (user->tree) {
328 DM refTree, ncdm = NULL;
329
330 PetscCall(DMPlexCreateDefaultReferenceTree(comm, dim, simplex, &refTree));
331 PetscCall(DMViewFromOptions(refTree, NULL, "-reftree_dm_view"));
332 PetscCall(DMPlexSetReferenceTree(*dm, refTree));
333 PetscCall(DMDestroy(&refTree));
334 PetscCall(DMPlexTreeRefineCell(*dm, user->treeCell, &ncdm));
335 if (ncdm) {
336 PetscCall(DMDestroy(dm));
337 *dm = ncdm;
338 PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
339 }
340 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "tree_"));
341 PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
342 PetscCall(DMSetFromOptions(*dm));
343 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
344 } else {
345 PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
346 }
347 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "dist_"));
348 PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
349 PetscCall(DMSetFromOptions(*dm));
350 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
351 if (simplex) PetscCall(PetscObjectSetName((PetscObject)*dm, "Simplicial Mesh"));
352 else PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
353 }
354 PetscCall(DMSetFromOptions(*dm));
355 PetscCall(TransformCoordinates(*dm, user));
356 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
357 PetscFunctionReturn(PETSC_SUCCESS);
358 }
359
simple_mass(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 g0[])360 static void simple_mass(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 g0[])
361 {
362 PetscInt d, e;
363 for (d = 0, e = 0; d < dim; d++, e += dim + 1) g0[e] = 1.;
364 }
365
366 /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */
symmetric_gradient_inner_product(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 C[])367 static void symmetric_gradient_inner_product(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 C[])
368 {
369 PetscInt compI, compJ, d, e;
370
371 for (compI = 0; compI < dim; ++compI) {
372 for (compJ = 0; compJ < dim; ++compJ) {
373 for (d = 0; d < dim; ++d) {
374 for (e = 0; e < dim; e++) {
375 if (d == e && d == compI && d == compJ) {
376 C[((compI * dim + compJ) * dim + d) * dim + e] = 1.0;
377 } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) {
378 C[((compI * dim + compJ) * dim + d) * dim + e] = 0.5;
379 } else {
380 C[((compI * dim + compJ) * dim + d) * dim + e] = 0.0;
381 }
382 }
383 }
384 }
385 }
386 }
387
SetupSection(DM dm,AppCtx * user)388 static PetscErrorCode SetupSection(DM dm, AppCtx *user)
389 {
390 PetscFunctionBeginUser;
391 if (user->constraints) {
392 /* test local constraints */
393 DM coordDM;
394 PetscInt fStart, fEnd, f, vStart, vEnd, v;
395 PetscInt edgesx = 2, vertsx;
396 PetscInt edgesy = 2, vertsy;
397 PetscMPIInt size;
398 PetscInt numConst;
399 PetscSection aSec;
400 PetscInt *anchors;
401 PetscInt offset;
402 IS aIS;
403 MPI_Comm comm = PetscObjectComm((PetscObject)dm);
404
405 PetscCallMPI(MPI_Comm_size(comm, &size));
406 PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Local constraint test can only be performed in serial");
407
408 /* we are going to test constraints by using them to enforce periodicity
409 * in one direction, and comparing to the existing method of enforcing
410 * periodicity */
411
412 /* first create the coordinate section so that it does not clone the
413 * constraints */
414 PetscCall(DMGetCoordinateDM(dm, &coordDM));
415
416 /* create the constrained-to-anchor section */
417 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
418 PetscCall(DMPlexGetDepthStratum(dm, 1, &fStart, &fEnd));
419 PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &aSec));
420 PetscCall(PetscSectionSetChart(aSec, PetscMin(fStart, vStart), PetscMax(fEnd, vEnd)));
421
422 /* define the constraints */
423 PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &edgesx, NULL));
424 PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &edgesy, NULL));
425 vertsx = edgesx + 1;
426 vertsy = edgesy + 1;
427 numConst = vertsy + edgesy;
428 PetscCall(PetscMalloc1(numConst, &anchors));
429 offset = 0;
430 for (v = vStart + edgesx; v < vEnd; v += vertsx) {
431 PetscCall(PetscSectionSetDof(aSec, v, 1));
432 anchors[offset++] = v - edgesx;
433 }
434 for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) {
435 PetscCall(PetscSectionSetDof(aSec, f, 1));
436 anchors[offset++] = f - edgesx * edgesy;
437 }
438 PetscCall(PetscSectionSetUp(aSec));
439 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numConst, anchors, PETSC_OWN_POINTER, &aIS));
440
441 PetscCall(DMPlexSetAnchors(dm, aSec, aIS));
442 PetscCall(PetscSectionDestroy(&aSec));
443 PetscCall(ISDestroy(&aIS));
444 }
445 PetscCall(DMSetNumFields(dm, 1));
446 PetscCall(DMSetField(dm, 0, NULL, (PetscObject)user->fe));
447 PetscCall(DMCreateDS(dm));
448 if (user->constraints) {
449 /* test getting local constraint matrix that matches section */
450 PetscSection aSec;
451 IS aIS;
452
453 PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
454 if (aSec) {
455 PetscDS ds;
456 PetscSection cSec, section;
457 PetscInt cStart, cEnd, c, numComp;
458 Mat cMat, mass;
459 Vec local;
460 const PetscInt *anchors;
461
462 PetscCall(DMGetLocalSection(dm, §ion));
463 /* this creates the matrix and preallocates the matrix nonzero structure: we
464 * just have to fill in the values */
465 PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
466 PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
467 PetscCall(ISGetIndices(aIS, &anchors));
468 PetscCall(PetscFEGetNumComponents(user->fe, &numComp));
469 for (c = cStart; c < cEnd; c++) {
470 PetscInt cDof;
471
472 /* is this point constrained? (does it have an anchor?) */
473 PetscCall(PetscSectionGetDof(aSec, c, &cDof));
474 if (cDof) {
475 PetscInt cOff, a, aDof, aOff, j;
476 PetscCheck(cDof == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found %" PetscInt_FMT " anchor points: should be just one", cDof);
477
478 /* find the anchor point */
479 PetscCall(PetscSectionGetOffset(aSec, c, &cOff));
480 a = anchors[cOff];
481
482 /* find the constrained dofs (row in constraint matrix) */
483 PetscCall(PetscSectionGetDof(cSec, c, &cDof));
484 PetscCall(PetscSectionGetOffset(cSec, c, &cOff));
485
486 /* find the anchor dofs (column in constraint matrix) */
487 PetscCall(PetscSectionGetDof(section, a, &aDof));
488 PetscCall(PetscSectionGetOffset(section, a, &aOff));
489
490 PetscCheck(cDof == aDof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point and anchor have different number of dofs: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, aDof);
491 PetscCheck(cDof % numComp == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point dofs not divisible by field components: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, numComp);
492
493 /* put in a simple equality constraint */
494 for (j = 0; j < cDof; j++) PetscCall(MatSetValue(cMat, cOff + j, aOff + j, 1., INSERT_VALUES));
495 }
496 }
497 PetscCall(MatAssemblyBegin(cMat, MAT_FINAL_ASSEMBLY));
498 PetscCall(MatAssemblyEnd(cMat, MAT_FINAL_ASSEMBLY));
499 PetscCall(ISRestoreIndices(aIS, &anchors));
500
501 /* Now that we have constructed the constraint matrix, any FE matrix
502 * that we construct will apply the constraints during construction */
503
504 PetscCall(DMCreateMatrix(dm, &mass));
505 /* get a local variable to serve as the solution */
506 PetscCall(DMGetLocalVector(dm, &local));
507 PetscCall(DMGetDS(dm, &ds));
508 /* set the jacobian to be the mass matrix */
509 PetscCall(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL));
510 /* build the mass matrix */
511 PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, mass, mass, NULL));
512 PetscCall(MatView(mass, PETSC_VIEWER_STDOUT_WORLD));
513 PetscCall(MatDestroy(&mass));
514 PetscCall(DMRestoreLocalVector(dm, &local));
515 }
516 }
517 PetscFunctionReturn(PETSC_SUCCESS);
518 }
519
TestFEJacobian(DM dm,AppCtx * user)520 static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user)
521 {
522 PetscFunctionBeginUser;
523 if (!user->useDA) {
524 Vec local;
525 const Vec *vecs;
526 Mat E;
527 MatNullSpace sp;
528 PetscBool isNullSpace, hasConst;
529 PetscInt dim, n, i;
530 Vec res = NULL, localX, localRes;
531 PetscDS ds;
532
533 PetscCall(DMGetDimension(dm, &dim));
534 PetscCheck(user->numComponents == dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "The number of components %" PetscInt_FMT " must be equal to the dimension %" PetscInt_FMT " for this test", user->numComponents, dim);
535 PetscCall(DMGetDS(dm, &ds));
536 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, symmetric_gradient_inner_product));
537 PetscCall(DMCreateMatrix(dm, &E));
538 PetscCall(DMGetLocalVector(dm, &local));
539 PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, E, E, NULL));
540 PetscCall(DMPlexCreateRigidBody(dm, 0, &sp));
541 PetscCall(MatNullSpaceGetVecs(sp, &hasConst, &n, &vecs));
542 if (n) PetscCall(VecDuplicate(vecs[0], &res));
543 PetscCall(DMCreateLocalVector(dm, &localX));
544 PetscCall(DMCreateLocalVector(dm, &localRes));
545 for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */
546 PetscReal resNorm;
547
548 PetscCall(VecSet(localRes, 0.));
549 PetscCall(VecSet(localX, 0.));
550 PetscCall(VecSet(local, 0.));
551 PetscCall(VecSet(res, 0.));
552 PetscCall(DMGlobalToLocalBegin(dm, vecs[i], INSERT_VALUES, localX));
553 PetscCall(DMGlobalToLocalEnd(dm, vecs[i], INSERT_VALUES, localX));
554 PetscCall(DMSNESComputeJacobianAction(dm, local, localX, localRes, NULL));
555 PetscCall(DMLocalToGlobalBegin(dm, localRes, ADD_VALUES, res));
556 PetscCall(DMLocalToGlobalEnd(dm, localRes, ADD_VALUES, res));
557 PetscCall(VecNorm(res, NORM_2, &resNorm));
558 if (resNorm > PETSC_SMALL) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient action null space vector %" PetscInt_FMT " residual: %E\n", i, (double)resNorm));
559 }
560 PetscCall(VecDestroy(&localRes));
561 PetscCall(VecDestroy(&localX));
562 PetscCall(VecDestroy(&res));
563 PetscCall(MatNullSpaceTest(sp, E, &isNullSpace));
564 if (isNullSpace) {
565 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: PASS\n"));
566 } else {
567 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: FAIL\n"));
568 }
569 PetscCall(MatNullSpaceDestroy(&sp));
570 PetscCall(MatDestroy(&E));
571 PetscCall(DMRestoreLocalVector(dm, &local));
572 }
573 PetscFunctionReturn(PETSC_SUCCESS);
574 }
575
TestInjector(DM dm,AppCtx * user)576 static PetscErrorCode TestInjector(DM dm, AppCtx *user)
577 {
578 DM refTree;
579 PetscMPIInt rank;
580
581 PetscFunctionBegin;
582 PetscCall(DMPlexGetReferenceTree(dm, &refTree));
583 if (refTree) {
584 Mat inj;
585
586 PetscCall(DMPlexComputeInjectorReferenceTree(refTree, &inj));
587 PetscCall(PetscObjectSetName((PetscObject)inj, "Reference Tree Injector"));
588 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
589 if (rank == 0) PetscCall(MatView(inj, PETSC_VIEWER_STDOUT_SELF));
590 PetscCall(MatDestroy(&inj));
591 }
592 PetscFunctionReturn(PETSC_SUCCESS);
593 }
594
TestFVGrad(DM dm,AppCtx * user)595 static PetscErrorCode TestFVGrad(DM dm, AppCtx *user)
596 {
597 MPI_Comm comm;
598 DM dmRedist, dmfv, dmgrad, dmCell, refTree;
599 PetscFV fv;
600 PetscInt dim, nvecs, v, cStart, cEnd, cEndInterior;
601 PetscMPIInt size;
602 Vec cellgeom, grad, locGrad;
603 const PetscScalar *cgeom;
604 PetscReal allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON;
605
606 PetscFunctionBeginUser;
607 comm = PetscObjectComm((PetscObject)dm);
608 PetscCall(DMGetDimension(dm, &dim));
609 /* duplicate DM, give dup. a FV discretization */
610 PetscCall(DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE));
611 PetscCallMPI(MPI_Comm_size(comm, &size));
612 dmRedist = NULL;
613 if (size > 1) PetscCall(DMPlexDistributeOverlap(dm, 1, NULL, &dmRedist));
614 if (!dmRedist) {
615 dmRedist = dm;
616 PetscCall(PetscObjectReference((PetscObject)dmRedist));
617 }
618 PetscCall(PetscFVCreate(comm, &fv));
619 PetscCall(PetscFVSetType(fv, PETSCFVLEASTSQUARES));
620 PetscCall(PetscFVSetNumComponents(fv, user->numComponents));
621 PetscCall(PetscFVSetSpatialDimension(fv, dim));
622 PetscCall(PetscFVSetFromOptions(fv));
623 PetscCall(PetscFVSetUp(fv));
624 {
625 PetscSF pointSF;
626 DMLabel label;
627
628 PetscCall(DMCreateLabel(dmRedist, "Face Sets"));
629 PetscCall(DMGetLabel(dmRedist, "Face Sets", &label));
630 PetscCall(DMGetPointSF(dmRedist, &pointSF));
631 PetscCall(PetscObjectReference((PetscObject)pointSF));
632 PetscCall(DMSetPointSF(dmRedist, NULL));
633 PetscCall(DMPlexMarkBoundaryFaces(dmRedist, 1, label));
634 PetscCall(DMSetPointSF(dmRedist, pointSF));
635 PetscCall(PetscSFDestroy(&pointSF));
636 }
637 PetscCall(DMPlexConstructGhostCells(dmRedist, NULL, NULL, &dmfv));
638 PetscCall(DMDestroy(&dmRedist));
639 PetscCall(DMSetNumFields(dmfv, 1));
640 PetscCall(DMSetField(dmfv, 0, NULL, (PetscObject)fv));
641 PetscCall(DMCreateDS(dmfv));
642 PetscCall(DMPlexGetReferenceTree(dm, &refTree));
643 if (refTree) PetscCall(DMCopyDisc(dmfv, refTree));
644 PetscCall(DMPlexGetGradientDM(dmfv, fv, &dmgrad));
645 PetscCall(DMPlexGetHeightStratum(dmfv, 0, &cStart, &cEnd));
646 nvecs = dim * (dim + 1) / 2;
647 PetscCall(DMPlexGetGeometryFVM(dmfv, NULL, &cellgeom, NULL));
648 PetscCall(VecGetDM(cellgeom, &dmCell));
649 PetscCall(VecGetArrayRead(cellgeom, &cgeom));
650 PetscCall(DMGetGlobalVector(dmgrad, &grad));
651 PetscCall(DMGetLocalVector(dmgrad, &locGrad));
652 PetscCall(DMPlexGetCellTypeStratum(dmgrad, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL));
653 cEndInterior = (cEndInterior < 0) ? cEnd : cEndInterior;
654 for (v = 0; v < nvecs; v++) {
655 Vec locX;
656 PetscInt c;
657 PetscScalar trueGrad[3][3] = {{0.}};
658 const PetscScalar *gradArray;
659 PetscReal maxDiff, maxDiffGlob;
660
661 PetscCall(DMGetLocalVector(dmfv, &locX));
662 /* get the local projection of the rigid body mode */
663 for (c = cStart; c < cEnd; c++) {
664 PetscFVCellGeom *cg;
665 PetscScalar cx[3] = {0., 0., 0.};
666
667 PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
668 if (v < dim) {
669 cx[v] = 1.;
670 } else {
671 PetscInt w = v - dim;
672
673 cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim];
674 cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim];
675 }
676 PetscCall(DMPlexVecSetClosure(dmfv, NULL, locX, c, cx, INSERT_ALL_VALUES));
677 }
678 /* TODO: this isn't in any header */
679 PetscCall(DMPlexReconstructGradientsFVM(dmfv, locX, grad));
680 PetscCall(DMGlobalToLocalBegin(dmgrad, grad, INSERT_VALUES, locGrad));
681 PetscCall(DMGlobalToLocalEnd(dmgrad, grad, INSERT_VALUES, locGrad));
682 PetscCall(VecGetArrayRead(locGrad, &gradArray));
683 /* compare computed gradient to exact gradient */
684 if (v >= dim) {
685 PetscInt w = v - dim;
686
687 trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.;
688 trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.;
689 }
690 maxDiff = 0.;
691 for (c = cStart; c < cEndInterior; c++) {
692 PetscScalar *compGrad;
693 PetscInt i, j, k;
694 PetscReal FrobDiff = 0.;
695
696 PetscCall(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad));
697
698 for (i = 0, k = 0; i < dim; i++) {
699 for (j = 0; j < dim; j++, k++) {
700 PetscScalar diff = compGrad[k] - trueGrad[i][j];
701 FrobDiff += PetscRealPart(diff * PetscConj(diff));
702 }
703 }
704 FrobDiff = PetscSqrtReal(FrobDiff);
705 maxDiff = PetscMax(maxDiff, FrobDiff);
706 }
707 PetscCallMPI(MPIU_Allreduce(&maxDiff, &maxDiffGlob, 1, MPIU_REAL, MPIU_MAX, comm));
708 allVecMaxDiff = PetscMax(allVecMaxDiff, maxDiffGlob);
709 PetscCall(VecRestoreArrayRead(locGrad, &gradArray));
710 PetscCall(DMRestoreLocalVector(dmfv, &locX));
711 }
712 if (allVecMaxDiff < fvTol) {
713 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: PASS\n"));
714 } else {
715 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n", (double)fvTol, (double)allVecMaxDiff));
716 }
717 PetscCall(DMRestoreLocalVector(dmgrad, &locGrad));
718 PetscCall(DMRestoreGlobalVector(dmgrad, &grad));
719 PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
720 PetscCall(DMDestroy(&dmfv));
721 PetscCall(PetscFVDestroy(&fv));
722 PetscFunctionReturn(PETSC_SUCCESS);
723 }
724
ComputeError(DM dm,PetscErrorCode (** exactFuncs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),PetscErrorCode (** exactFuncDers)(PetscInt,PetscReal,const PetscReal[],const PetscReal[],PetscInt,PetscScalar *,void *),void ** exactCtxs,PetscReal * error,PetscReal * errorDer,AppCtx * user)725 static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user)
726 {
727 Vec u;
728 PetscReal n[3] = {1.0, 1.0, 1.0};
729
730 PetscFunctionBeginUser;
731 PetscCall(DMGetGlobalVector(dm, &u));
732 /* Project function into FE function space */
733 PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u));
734 PetscCall(VecViewFromOptions(u, NULL, "-projection_view"));
735 /* Compare approximation to exact in L_2 */
736 PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error));
737 PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer));
738 PetscCall(DMRestoreGlobalVector(dm, &u));
739 PetscFunctionReturn(PETSC_SUCCESS);
740 }
741
CheckFunctions(DM dm,PetscInt order,AppCtx * user)742 static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user)
743 {
744 PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
745 PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
746 void *exactCtxs[3];
747 MPI_Comm comm;
748 PetscReal error, errorDer, tol = PETSC_SMALL;
749
750 PetscFunctionBeginUser;
751 exactCtxs[0] = user;
752 exactCtxs[1] = user;
753 exactCtxs[2] = user;
754 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
755 /* Setup functions to approximate */
756 switch (order) {
757 case 0:
758 if (user->RT) {
759 exactFuncs[0] = rt0;
760 exactFuncDers[0] = rt0Der;
761 } else {
762 exactFuncs[0] = constant;
763 exactFuncDers[0] = constantDer;
764 }
765 break;
766 case 1:
767 if (user->RT) {
768 exactFuncs[0] = rt1;
769 exactFuncDers[0] = rt1Der;
770 } else {
771 exactFuncs[0] = linear;
772 exactFuncDers[0] = linearDer;
773 }
774 break;
775 case 2:
776 exactFuncs[0] = quadratic;
777 exactFuncDers[0] = quadraticDer;
778 break;
779 case 3:
780 exactFuncs[0] = cubic;
781 exactFuncDers[0] = cubicDer;
782 break;
783 default:
784 SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %" PetscInt_FMT, order);
785 }
786 PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
787 /* Report result */
788 if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
789 else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
790 if (errorDer > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer));
791 else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
792 PetscFunctionReturn(PETSC_SUCCESS);
793 }
794
CheckInterpolation(DM dm,PetscBool checkRestrict,PetscInt order,AppCtx * user)795 static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user)
796 {
797 PetscErrorCode (*exactFuncs[1])(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, PetscCtx ctx);
798 PetscErrorCode (*exactFuncDers[1])(PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, PetscCtx ctx);
799 PetscReal n[3] = {1.0, 1.0, 1.0};
800 void *exactCtxs[3];
801 DM rdm, idm, fdm;
802 Mat Interp;
803 Vec iu, fu, scaling;
804 MPI_Comm comm;
805 PetscInt dim;
806 PetscReal error, errorDer, tol = PETSC_SMALL;
807
808 PetscFunctionBeginUser;
809 exactCtxs[0] = user;
810 exactCtxs[1] = user;
811 exactCtxs[2] = user;
812 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
813 PetscCall(DMGetDimension(dm, &dim));
814 PetscCall(DMRefine(dm, comm, &rdm));
815 PetscCall(DMSetCoarseDM(rdm, dm));
816 PetscCall(DMGetCoordinatesLocalSetUp(rdm));
817 PetscCall(DMPlexSetRegularRefinement(rdm, user->convRefine));
818 if (user->tree) {
819 DM refTree;
820 PetscCall(DMPlexGetReferenceTree(dm, &refTree));
821 PetscCall(DMPlexSetReferenceTree(rdm, refTree));
822 }
823 if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
824 PetscCall(SetupSection(rdm, user));
825 /* Setup functions to approximate */
826 switch (order) {
827 case 0:
828 exactFuncs[0] = constant;
829 exactFuncDers[0] = constantDer;
830 break;
831 case 1:
832 exactFuncs[0] = linear;
833 exactFuncDers[0] = linearDer;
834 break;
835 case 2:
836 exactFuncs[0] = quadratic;
837 exactFuncDers[0] = quadraticDer;
838 break;
839 case 3:
840 exactFuncs[0] = cubic;
841 exactFuncDers[0] = cubicDer;
842 break;
843 default:
844 SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order);
845 }
846 idm = checkRestrict ? rdm : dm;
847 fdm = checkRestrict ? dm : rdm;
848 PetscCall(DMGetGlobalVector(idm, &iu));
849 PetscCall(DMGetGlobalVector(fdm, &fu));
850 PetscCall(DMSetApplicationContext(dm, user));
851 PetscCall(DMSetApplicationContext(rdm, user));
852 PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling));
853 /* Project function into initial FE function space */
854 PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu));
855 /* Interpolate function into final FE function space */
856 if (checkRestrict) {
857 PetscCall(MatRestrict(Interp, iu, fu));
858 PetscCall(VecPointwiseMult(fu, scaling, fu));
859 } else PetscCall(MatInterpolate(Interp, iu, fu));
860 /* Compare approximation to exact in L_2 */
861 PetscCall(DMGetCoordinatesLocalSetUp(fdm));
862 PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error));
863 PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer));
864 /* Report result */
865 if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
866 else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
867 if (errorDer > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer));
868 else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
869 PetscCall(DMRestoreGlobalVector(idm, &iu));
870 PetscCall(DMRestoreGlobalVector(fdm, &fu));
871 PetscCall(MatDestroy(&Interp));
872 PetscCall(VecDestroy(&scaling));
873 PetscCall(DMDestroy(&rdm));
874 PetscFunctionReturn(PETSC_SUCCESS);
875 }
876
CheckConvergence(DM dm,PetscInt Nr,AppCtx * user)877 static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user)
878 {
879 DM odm = dm, rdm = NULL, cdm = NULL;
880 PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) = {trig};
881 PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) = {trigDer};
882 void *exactCtxs[3];
883 PetscInt r, c, cStart, cEnd;
884 PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld;
885 double p;
886
887 PetscFunctionBeginUser;
888 if (!user->convergence) PetscFunctionReturn(PETSC_SUCCESS);
889 exactCtxs[0] = user;
890 exactCtxs[1] = user;
891 exactCtxs[2] = user;
892 PetscCall(PetscObjectReference((PetscObject)odm));
893 if (!user->convRefine) {
894 for (r = 0; r < Nr; ++r) {
895 PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
896 PetscCall(DMDestroy(&odm));
897 odm = rdm;
898 }
899 PetscCall(SetupSection(odm, user));
900 }
901 PetscCall(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user));
902 if (user->convRefine) {
903 for (r = 0; r < Nr; ++r) {
904 PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
905 if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
906 PetscCall(SetupSection(rdm, user));
907 PetscCall(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
908 p = PetscLog2Real(errorOld / error);
909 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, p));
910 p = PetscLog2Real(errorDerOld / errorDer);
911 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, p));
912 PetscCall(DMDestroy(&odm));
913 odm = rdm;
914 errorOld = error;
915 errorDerOld = errorDer;
916 }
917 } else {
918 /* PetscCall(ComputeLongestEdge(dm, &lenOld)); */
919 PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
920 lenOld = cEnd - cStart;
921 for (c = 0; c < Nr; ++c) {
922 PetscCall(DMCoarsen(odm, PetscObjectComm((PetscObject)dm), &cdm));
923 if (user->useDA) PetscCall(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
924 PetscCall(SetupSection(cdm, user));
925 PetscCall(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
926 /* PetscCall(ComputeLongestEdge(cdm, &len)); */
927 PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
928 len = cEnd - cStart;
929 rel = error / errorOld;
930 p = PetscLogReal(rel) / PetscLogReal(lenOld / len);
931 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, p));
932 rel = errorDer / errorDerOld;
933 p = PetscLogReal(rel) / PetscLogReal(lenOld / len);
934 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, p));
935 PetscCall(DMDestroy(&odm));
936 odm = cdm;
937 errorOld = error;
938 errorDerOld = errorDer;
939 lenOld = len;
940 }
941 }
942 PetscCall(DMDestroy(&odm));
943 PetscFunctionReturn(PETSC_SUCCESS);
944 }
945
main(int argc,char ** argv)946 int main(int argc, char **argv)
947 {
948 DM dm;
949 AppCtx user; /* user-defined work context */
950 PetscInt dim = 2;
951 PetscBool simplex = PETSC_FALSE;
952
953 PetscFunctionBeginUser;
954 PetscCall(PetscInitialize(&argc, &argv, NULL, help));
955 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
956 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
957 if (!user.useDA) {
958 PetscCall(DMGetDimension(dm, &dim));
959 PetscCall(DMPlexIsSimplex(dm, &simplex));
960 }
961 PetscCall(DMPlexMetricSetFromOptions(dm));
962 user.numComponents = user.numComponents < 0 ? dim : user.numComponents;
963 PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe));
964 PetscCall(SetupSection(dm, &user));
965 if (user.testFEjacobian) PetscCall(TestFEJacobian(dm, &user));
966 if (user.testFVgrad) PetscCall(TestFVGrad(dm, &user));
967 if (user.testInjector) PetscCall(TestInjector(dm, &user));
968 PetscCall(CheckFunctions(dm, user.porder, &user));
969 {
970 PetscDualSpace dsp;
971 PetscInt k;
972
973 PetscCall(PetscFEGetDualSpace(user.fe, &dsp));
974 PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
975 if (dim == 2 && user.constraints == PETSC_FALSE && user.tree == PETSC_FALSE && k == 0) {
976 PetscCall(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user));
977 PetscCall(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user));
978 }
979 }
980 PetscCall(CheckConvergence(dm, 3, &user));
981 PetscCall(PetscFEDestroy(&user.fe));
982 PetscCall(DMDestroy(&dm));
983 PetscCall(PetscFinalize());
984 return 0;
985 }
986
987 /*TEST
988
989 test:
990 suffix: 1
991 requires: triangle
992
993 # 2D P_0 on a triangle
994 test:
995 suffix: p0_2d_0
996 requires: triangle
997 args: -petscspace_degree 1 -qorder 1 -convergence
998 test:
999 suffix: p0_2d_1
1000 requires: triangle
1001 args: -petscspace_degree 1 -qorder 1 -porder 1
1002
1003 # 2D P_1 on a triangle
1004 test:
1005 suffix: p1_2d_0
1006 requires: triangle
1007 args: -petscspace_degree 1 -qorder 1 -convergence
1008 test:
1009 suffix: p1_2d_1
1010 requires: triangle
1011 args: -petscspace_degree 1 -qorder 1 -porder 1
1012 test:
1013 suffix: p1_2d_2
1014 requires: triangle
1015 args: -petscspace_degree 1 -qorder 1 -porder 2
1016 test:
1017 suffix: p1_2d_3
1018 requires: triangle mmg
1019 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1020 test:
1021 suffix: p1_2d_4
1022 requires: triangle mmg
1023 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1024 test:
1025 suffix: p1_2d_5
1026 requires: triangle mmg
1027 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1028
1029 # 3D P_1 on a tetrahedron
1030 test:
1031 suffix: p1_3d_0
1032 requires: ctetgen
1033 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence
1034 test:
1035 suffix: p1_3d_1
1036 requires: ctetgen
1037 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1
1038 test:
1039 suffix: p1_3d_2
1040 requires: ctetgen
1041 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2
1042 test:
1043 suffix: p1_3d_3
1044 requires: ctetgen mmg
1045 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1046 test:
1047 suffix: p1_3d_4
1048 requires: ctetgen mmg
1049 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1050 test:
1051 suffix: p1_3d_5
1052 requires: ctetgen mmg
1053 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1054
1055 # 2D P_2 on a triangle
1056 test:
1057 suffix: p2_2d_0
1058 requires: triangle
1059 args: -petscspace_degree 2 -qorder 2 -convergence
1060 test:
1061 suffix: p2_2d_1
1062 requires: triangle
1063 args: -petscspace_degree 2 -qorder 2 -porder 1
1064 test:
1065 suffix: p2_2d_2
1066 requires: triangle
1067 args: -petscspace_degree 2 -qorder 2 -porder 2
1068 test:
1069 suffix: p2_2d_3
1070 requires: triangle mmg
1071 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
1072 test:
1073 suffix: p2_2d_4
1074 requires: triangle mmg
1075 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1076 test:
1077 suffix: p2_2d_5
1078 requires: triangle mmg
1079 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1080
1081 # 3D P_2 on a tetrahedron
1082 test:
1083 suffix: p2_3d_0
1084 requires: ctetgen
1085 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence
1086 test:
1087 suffix: p2_3d_1
1088 requires: ctetgen
1089 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1
1090 test:
1091 suffix: p2_3d_2
1092 requires: ctetgen
1093 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2
1094 test:
1095 suffix: p2_3d_3
1096 requires: ctetgen mmg
1097 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
1098 test:
1099 suffix: p2_3d_4
1100 requires: ctetgen mmg
1101 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1102 test:
1103 suffix: p2_3d_5
1104 requires: ctetgen mmg
1105 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1106
1107 # 2D Q_1 on a quadrilaterial DA
1108 test:
1109 suffix: q1_2d_da_0
1110 TODO: broken
1111 args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence
1112 test:
1113 suffix: q1_2d_da_1
1114 TODO: broken
1115 args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1
1116 test:
1117 suffix: q1_2d_da_2
1118 TODO: broken
1119 args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2
1120
1121 # 2D P_0 on a quadrilaterial Plex
1122 test:
1123 suffix: p0_2d_plex_0
1124 args: -dm_plex_simplex 0 -petscspace_degree 0 -qorder 1 -convergence
1125 test:
1126 suffix: p0_2d_plex_1
1127 args: -dm_plex_simplex 0 -petscspace_degree 0 -qorder 1 -porder 1
1128
1129 # 2D Q_1 on a quadrilaterial Plex
1130 test:
1131 suffix: q1_2d_plex_0
1132 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence
1133 test:
1134 suffix: q1_2d_plex_1
1135 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1
1136 test:
1137 suffix: q1_2d_plex_2
1138 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2
1139 test:
1140 suffix: q1_2d_plex_3
1141 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords
1142 test:
1143 suffix: q1_2d_plex_4
1144 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords
1145 test:
1146 suffix: q1_2d_plex_5
1147 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence
1148 test:
1149 suffix: q1_2d_plex_6
1150 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence
1151 test:
1152 suffix: q1_2d_plex_7
1153 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence
1154 test:
1155 suffix: q1_2d_plex_8
1156 requires: triangle
1157 args: -dist_dm_refine 1 -dist_dm_plex_transform_type refine_tobox -petscspace_degree 1 -qorder 1 -convergence
1158
1159 # 2D Q_2 on a quadrilaterial
1160 # The derivative interpolation test fails in single because we lose precision
1161 test:
1162 suffix: q2_2d_plex_0
1163 requires: !single
1164 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1165 test:
1166 suffix: q2_2d_plex_1
1167 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1168 test:
1169 suffix: q2_2d_plex_2
1170 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1171 test:
1172 suffix: q2_2d_plex_3
1173 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1174 test:
1175 suffix: q2_2d_plex_4
1176 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1177 # The derivative interpolation test fails in single because we lose precision
1178 test:
1179 suffix: q2_2d_plex_5
1180 requires: !single
1181 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence
1182 # The derivative interpolation test fails in single because we lose precision
1183 test:
1184 suffix: q2_2d_plex_6
1185 requires: !single
1186 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence
1187 test:
1188 suffix: q2_2d_plex_7
1189 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence
1190
1191 # 2D P_3 on a triangle
1192 test:
1193 suffix: p3_2d_0
1194 requires: triangle !single
1195 args: -petscspace_degree 3 -qorder 3 -convergence
1196 test:
1197 suffix: p3_2d_1
1198 requires: triangle !single
1199 args: -petscspace_degree 3 -qorder 3 -porder 1
1200 test:
1201 suffix: p3_2d_2
1202 requires: triangle !single
1203 args: -petscspace_degree 3 -qorder 3 -porder 2
1204 test:
1205 suffix: p3_2d_3
1206 requires: triangle !single
1207 args: -petscspace_degree 3 -qorder 3 -porder 3
1208 test:
1209 suffix: p3_2d_4
1210 requires: triangle mmg
1211 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0
1212 test:
1213 suffix: p3_2d_5
1214 requires: triangle mmg
1215 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0
1216 test:
1217 suffix: p3_2d_6
1218 requires: triangle mmg
1219 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0
1220
1221 # 2D Q_3 on a quadrilaterial
1222 test:
1223 suffix: q3_2d_0
1224 requires: !single
1225 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence
1226 test:
1227 suffix: q3_2d_1
1228 requires: !single
1229 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1
1230 test:
1231 suffix: q3_2d_2
1232 requires: !single
1233 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2
1234 test:
1235 suffix: q3_2d_3
1236 requires: !single
1237 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3
1238
1239 # 2D P_1disc on a triangle/quadrilateral
1240 test:
1241 suffix: p1d_2d_0
1242 requires: triangle
1243 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1244 test:
1245 suffix: p1d_2d_1
1246 requires: triangle
1247 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1248 test:
1249 suffix: p1d_2d_2
1250 requires: triangle
1251 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1252 test:
1253 suffix: p1d_2d_3
1254 requires: triangle
1255 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1256 filter: sed -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g"
1257 test:
1258 suffix: p1d_2d_4
1259 requires: triangle
1260 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1261 test:
1262 suffix: p1d_2d_5
1263 requires: triangle
1264 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1265
1266 # 2D BDM_1 on a triangle
1267 test:
1268 suffix: bdm1_2d_0
1269 requires: triangle
1270 args: -petscspace_degree 1 -petscdualspace_type bdm \
1271 -num_comp 2 -qorder 1 -convergence
1272 test:
1273 suffix: bdm1_2d_1
1274 requires: triangle
1275 args: -petscspace_degree 1 -petscdualspace_type bdm \
1276 -num_comp 2 -qorder 1 -porder 1
1277 test:
1278 suffix: bdm1_2d_2
1279 requires: triangle
1280 args: -petscspace_degree 1 -petscdualspace_type bdm \
1281 -num_comp 2 -qorder 1 -porder 2
1282
1283 # 2D BDM_1 on a quadrilateral
1284 test:
1285 suffix: bdm1q_2d_0
1286 requires: triangle
1287 args: -petscspace_degree 1 -petscdualspace_type bdm \
1288 -petscdualspace_lagrange_tensor 1 \
1289 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence
1290 test:
1291 suffix: bdm1q_2d_1
1292 requires: triangle
1293 args: -petscspace_degree 1 -petscdualspace_type bdm \
1294 -petscdualspace_lagrange_tensor 1 \
1295 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1
1296 test:
1297 suffix: bdm1q_2d_2
1298 requires: triangle
1299 args: -petscspace_degree 1 -petscdualspace_type bdm \
1300 -petscdualspace_lagrange_tensor 1 \
1301 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2
1302
1303 # Test high order quadrature
1304 test:
1305 suffix: p1_quad_2
1306 requires: triangle
1307 args: -petscspace_degree 1 -qorder 2 -porder 1
1308 test:
1309 suffix: p1_quad_5
1310 requires: triangle
1311 args: -petscspace_degree 1 -qorder 5 -porder 1
1312 test:
1313 suffix: p2_quad_3
1314 requires: triangle
1315 args: -petscspace_degree 2 -qorder 3 -porder 2
1316 test:
1317 suffix: p2_quad_5
1318 requires: triangle
1319 args: -petscspace_degree 2 -qorder 5 -porder 2
1320 test:
1321 suffix: q1_quad_2
1322 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1
1323 test:
1324 suffix: q1_quad_5
1325 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1
1326 test:
1327 suffix: q2_quad_3
1328 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1
1329 test:
1330 suffix: q2_quad_5
1331 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1
1332
1333 # Nonconforming tests
1334 test:
1335 suffix: constraints
1336 args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints
1337 test:
1338 suffix: nonconforming_tensor_2
1339 nsize: 4
1340 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1341 test:
1342 suffix: nonconforming_tensor_3
1343 nsize: 4
1344 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 1 -qorder 1 -dm_view ascii::ASCII_INFO_DETAIL
1345 test:
1346 suffix: nonconforming_tensor_2_fv
1347 nsize: 4
1348 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2
1349 test:
1350 suffix: nonconforming_tensor_3_fv
1351 nsize: 4
1352 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -num_comp 3
1353 test:
1354 suffix: nonconforming_tensor_2_hi
1355 requires: !single
1356 nsize: 4
1357 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 4 -qorder 4
1358 test:
1359 suffix: nonconforming_tensor_3_hi
1360 requires: !single skip
1361 nsize: 4
1362 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 4 -qorder 4
1363 test:
1364 suffix: nonconforming_simplex_2
1365 requires: triangle
1366 nsize: 4
1367 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1368 test:
1369 suffix: nonconforming_simplex_2_hi
1370 requires: triangle !single
1371 nsize: 4
1372 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4
1373 test:
1374 suffix: nonconforming_simplex_2_fv
1375 requires: triangle
1376 nsize: 4
1377 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2
1378 test:
1379 suffix: nonconforming_simplex_3
1380 requires: ctetgen
1381 nsize: 4
1382 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1383 test:
1384 suffix: nonconforming_simplex_3_hi
1385 requires: ctetgen skip
1386 nsize: 4
1387 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 4 -qorder 4
1388 test:
1389 suffix: nonconforming_simplex_3_fv
1390 requires: ctetgen
1391 nsize: 4
1392 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3
1393
1394 # 3D WXY on a triangular prism
1395 test:
1396 suffix: wxy_0
1397 args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \
1398 -petscspace_type sum \
1399 -petscspace_variables 3 \
1400 -petscspace_components 3 \
1401 -petscspace_sum_spaces 2 \
1402 -petscspace_sum_concatenate false \
1403 -sumcomp_0_petscspace_variables 3 \
1404 -sumcomp_0_petscspace_components 3 \
1405 -sumcomp_0_petscspace_degree 1 \
1406 -sumcomp_1_petscspace_variables 3 \
1407 -sumcomp_1_petscspace_components 3 \
1408 -sumcomp_1_petscspace_type wxy \
1409 -petscdualspace_refcell triangular_prism \
1410 -petscdualspace_form_degree 0 \
1411 -petscdualspace_order 1 \
1412 -petscdualspace_components 3
1413
1414 # 2D RT_0 on a triangle
1415 test:
1416 suffix: rt0_2d_tri
1417 requires: triangle
1418 args: -qorder 1 -porder 0 -RT \
1419 -petscspace_type ptrimmed \
1420 -petscspace_components 2 \
1421 -petscspace_ptrimmed_form_degree -1 \
1422 -petscdualspace_order 1 \
1423 -petscdualspace_form_degree -1 \
1424 -petscdualspace_lagrange_trimmed true
1425
1426 # 2D RT_0 on a quadrilateral
1427 test:
1428 suffix: rt0_2d_quad
1429 requires: triangle
1430 args: -dm_plex_simplex 0 -qorder 1 -porder 0 -RT \
1431 -petscspace_degree 1 \
1432 -petscspace_type sum \
1433 -petscspace_variables 2 \
1434 -petscspace_components 2 \
1435 -petscspace_sum_spaces 2 \
1436 -petscspace_sum_concatenate true \
1437 -sumcomp_0_petscspace_variables 2 \
1438 -sumcomp_0_petscspace_type tensor \
1439 -sumcomp_0_petscspace_tensor_spaces 2 \
1440 -sumcomp_0_petscspace_tensor_uniform false \
1441 -sumcomp_0_tensorcomp_0_petscspace_degree 1 \
1442 -sumcomp_0_tensorcomp_1_petscspace_degree 0 \
1443 -sumcomp_1_petscspace_variables 2 \
1444 -sumcomp_1_petscspace_type tensor \
1445 -sumcomp_1_petscspace_tensor_spaces 2 \
1446 -sumcomp_1_petscspace_tensor_uniform false \
1447 -sumcomp_1_tensorcomp_0_petscspace_degree 0 \
1448 -sumcomp_1_tensorcomp_1_petscspace_degree 1 \
1449 -petscdualspace_form_degree -1 \
1450 -petscdualspace_order 1 \
1451 -petscdualspace_lagrange_trimmed true
1452
1453 TEST*/
1454
1455 /*
1456 # 2D Q_2 on a quadrilaterial Plex
1457 test:
1458 suffix: q2_2d_plex_0
1459 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1460 test:
1461 suffix: q2_2d_plex_1
1462 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1463 test:
1464 suffix: q2_2d_plex_2
1465 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1466 test:
1467 suffix: q2_2d_plex_3
1468 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1469 test:
1470 suffix: q2_2d_plex_4
1471 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1472 test:
1473 suffix: q2_2d_plex_5
1474 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords
1475 test:
1476 suffix: q2_2d_plex_6
1477 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords
1478 test:
1479 suffix: q2_2d_plex_7
1480 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords
1481
1482 test:
1483 suffix: p1d_2d_6
1484 requires: mmg
1485 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1486 test:
1487 suffix: p1d_2d_7
1488 requires: mmg
1489 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1490 test:
1491 suffix: p1d_2d_8
1492 requires: mmg
1493 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1494 */
1495