xref: /petsc/src/dm/dt/dualspace/impls/lagrange/tests/ex1.c (revision 7e1a0bbe36d2be40a00a95404ece00db4857f70d)
1 #include <petscfe.h>
2 #include <petscdmplex.h>
3 #include <petsc/private/hashmap.h>
4 #include <petsc/private/dmpleximpl.h>
5 #include <petsc/private/petscfeimpl.h>
6 
7 const char help[] = "Test PETSCDUALSPACELAGRANGE\n";
8 
9 typedef struct _PetscHashLagKey {
10   PetscInt  dim;
11   PetscInt  order;
12   PetscInt  formDegree;
13   PetscBool trimmed;
14   PetscInt  tensor;
15   PetscBool continuous;
16 } PetscHashLagKey;
17 
18 #define PetscHashLagKeyHash(key) \
19   PetscHashCombine(PetscHashCombine(PetscHashCombine(PetscHashInt((key).dim), PetscHashInt((key).order)), PetscHashInt((key).formDegree)), PetscHashCombine(PetscHashCombine(PetscHashInt((key).trimmed), PetscHashInt((key).tensor)), PetscHashInt((key).continuous)))
20 
21 #define PetscHashLagKeyEqual(k1, k2) \
22   (((k1).dim == (k2).dim) ? ((k1).order == (k2).order) ? ((k1).formDegree == (k2).formDegree) ? ((k1).trimmed == (k2).trimmed) ? ((k1).tensor == (k2).tensor) ? ((k1).continuous == (k2).continuous) : 0 : 0 : 0 : 0 : 0)
23 
24 PETSC_HASH_MAP(HashLag, PetscHashLagKey, PetscInt, PetscHashLagKeyHash, PetscHashLagKeyEqual, 0)
25 
26 static PetscErrorCode ExpectedNumDofs_Total(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs);
27 static PetscErrorCode ExpectedNumDofs_Interior(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs);
28 
29 static PetscErrorCode ExpectedNumDofs_Total(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs)
30 {
31   PetscFunctionBegin;
32   formDegree = PetscAbsInt(formDegree);
33   /* see femtable.org for the source of most of these values */
34   *nDofs = -1;
35   if (tensor == 0) { /* simplex spaces */
36     if (trimmed) {
37       PetscInt rnchooserk;
38       PetscInt rkm1choosek;
39 
40       PetscCall(PetscDTBinomialInt(order + dim, order + formDegree, &rnchooserk));
41       PetscCall(PetscDTBinomialInt(order + formDegree - 1, formDegree, &rkm1choosek));
42       *nDofs = rnchooserk * rkm1choosek * nCopies;
43     } else {
44       PetscInt rnchooserk;
45       PetscInt rkchoosek;
46 
47       PetscCall(PetscDTBinomialInt(order + dim, order + formDegree, &rnchooserk));
48       PetscCall(PetscDTBinomialInt(order + formDegree, formDegree, &rkchoosek));
49       *nDofs = rnchooserk * rkchoosek * nCopies;
50     }
51   } else if (tensor == 1) { /* hypercubes */
52     if (trimmed) {
53       PetscInt nchoosek;
54       PetscInt rpowk, rp1pownmk;
55 
56       PetscCall(PetscDTBinomialInt(dim, formDegree, &nchoosek));
57       rpowk     = PetscPowInt(order, formDegree);
58       rp1pownmk = PetscPowInt(order + 1, dim - formDegree);
59       *nDofs    = nchoosek * rpowk * rp1pownmk * nCopies;
60     } else {
61       PetscInt nchoosek;
62       PetscInt rp1pown;
63 
64       PetscCall(PetscDTBinomialInt(dim, formDegree, &nchoosek));
65       rp1pown = PetscPowInt(order + 1, dim);
66       *nDofs  = nchoosek * rp1pown * nCopies;
67     }
68   } else { /* prism */
69     PetscInt tracek   = 0;
70     PetscInt tracekm1 = 0;
71     PetscInt fiber0   = 0;
72     PetscInt fiber1   = 0;
73 
74     if (formDegree < dim) {
75       PetscCall(ExpectedNumDofs_Total(dim - 1, order, formDegree, trimmed, 0, 1, &tracek));
76       PetscCall(ExpectedNumDofs_Total(1, order, 0, trimmed, 0, 1, &fiber0));
77     }
78     if (formDegree > 0) {
79       PetscCall(ExpectedNumDofs_Total(dim - 1, order, formDegree - 1, trimmed, 0, 1, &tracekm1));
80       PetscCall(ExpectedNumDofs_Total(1, order, 1, trimmed, 0, 1, &fiber1));
81     }
82     *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies;
83   }
84   PetscFunctionReturn(PETSC_SUCCESS);
85 }
86 
87 static PetscErrorCode ExpectedNumDofs_Interior(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs)
88 {
89   PetscFunctionBegin;
90   formDegree = PetscAbsInt(formDegree);
91   /* see femtable.org for the source of most of these values */
92   *nDofs = -1;
93   if (tensor == 0) { /* simplex spaces */
94     if (trimmed) {
95       if (order + formDegree > dim) {
96         PetscInt eorder      = order + formDegree - dim - 1;
97         PetscInt eformDegree = dim - formDegree;
98         PetscInt rnchooserk;
99         PetscInt rkchoosek;
100 
101         PetscCall(PetscDTBinomialInt(eorder + dim, eorder + eformDegree, &rnchooserk));
102         PetscCall(PetscDTBinomialInt(eorder + eformDegree, eformDegree, &rkchoosek));
103         *nDofs = rnchooserk * rkchoosek * nCopies;
104       } else {
105         *nDofs = 0;
106       }
107 
108     } else {
109       if (order + formDegree > dim) {
110         PetscInt eorder      = order + formDegree - dim;
111         PetscInt eformDegree = dim - formDegree;
112         PetscInt rnchooserk;
113         PetscInt rkm1choosek;
114 
115         PetscCall(PetscDTBinomialInt(eorder + dim, eorder + eformDegree, &rnchooserk));
116         PetscCall(PetscDTBinomialInt(eorder + eformDegree - 1, eformDegree, &rkm1choosek));
117         *nDofs = rnchooserk * rkm1choosek * nCopies;
118       } else {
119         *nDofs = 0;
120       }
121     }
122   } else if (tensor == 1) { /* hypercubes */
123     if (dim < 2) {
124       PetscCall(ExpectedNumDofs_Interior(dim, order, formDegree, trimmed, 0, nCopies, nDofs));
125     } else {
126       PetscInt tracek   = 0;
127       PetscInt tracekm1 = 0;
128       PetscInt fiber0   = 0;
129       PetscInt fiber1   = 0;
130 
131       if (formDegree < dim) {
132         PetscCall(ExpectedNumDofs_Interior(dim - 1, order, formDegree, trimmed, dim > 2 ? 1 : 0, 1, &tracek));
133         PetscCall(ExpectedNumDofs_Interior(1, order, 0, trimmed, 0, 1, &fiber0));
134       }
135       if (formDegree > 0) {
136         PetscCall(ExpectedNumDofs_Interior(dim - 1, order, formDegree - 1, trimmed, dim > 2 ? 1 : 0, 1, &tracekm1));
137         PetscCall(ExpectedNumDofs_Interior(1, order, 1, trimmed, 0, 1, &fiber1));
138       }
139       *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies;
140     }
141   } else { /* prism */
142     PetscInt tracek   = 0;
143     PetscInt tracekm1 = 0;
144     PetscInt fiber0   = 0;
145     PetscInt fiber1   = 0;
146 
147     if (formDegree < dim) {
148       PetscCall(ExpectedNumDofs_Interior(dim - 1, order, formDegree, trimmed, 0, 1, &tracek));
149       PetscCall(ExpectedNumDofs_Interior(1, order, 0, trimmed, 0, 1, &fiber0));
150     }
151     if (formDegree > 0) {
152       PetscCall(ExpectedNumDofs_Interior(dim - 1, order, formDegree - 1, trimmed, 0, 1, &tracekm1));
153       PetscCall(ExpectedNumDofs_Interior(1, order, 1, trimmed, 0, 1, &fiber1));
154     }
155     *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies;
156   }
157   PetscFunctionReturn(PETSC_SUCCESS);
158 }
159 
160 PetscErrorCode testLagrange(PetscHashLag lagTable, DM K, PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensorCell, PetscBool continuous, PetscInt nCopies)
161 {
162   PetscDualSpace  sp;
163   MPI_Comm        comm = PETSC_COMM_SELF;
164   PetscInt        Nk;
165   PetscHashLagKey key;
166   PetscHashIter   iter;
167   PetscBool       missing;
168   PetscInt        spdim, spintdim, exspdim, exspintdim;
169 
170   PetscFunctionBegin;
171   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk));
172   PetscCall(PetscDualSpaceCreate(comm, &sp));
173   PetscCall(PetscDualSpaceSetType(sp, PETSCDUALSPACELAGRANGE));
174   PetscCall(PetscDualSpaceSetDM(sp, K));
175   PetscCall(PetscDualSpaceSetOrder(sp, order));
176   PetscCall(PetscDualSpaceSetFormDegree(sp, formDegree));
177   PetscCall(PetscDualSpaceSetNumComponents(sp, nCopies * Nk));
178   PetscCall(PetscDualSpaceLagrangeSetContinuity(sp, continuous));
179   PetscCall(PetscDualSpaceLagrangeSetTensor(sp, (PetscBool)(tensorCell > 0)));
180   PetscCall(PetscDualSpaceLagrangeSetTrimmed(sp, trimmed));
181   PetscCall(PetscInfo(NULL, "Input: dim %" PetscInt_FMT ", order %" PetscInt_FMT ", trimmed %" PetscInt_FMT ", tensorCell %" PetscInt_FMT ", continuous %" PetscInt_FMT ", formDegree %" PetscInt_FMT ", nCopies %" PetscInt_FMT "\n", dim, order, (PetscInt)trimmed, tensorCell, (PetscInt)continuous, formDegree, nCopies));
182   PetscCall(ExpectedNumDofs_Total(dim, order, formDegree, trimmed, tensorCell, nCopies, &exspdim));
183   if (continuous && dim > 0 && order > 0) {
184     PetscCall(ExpectedNumDofs_Interior(dim, order, formDegree, trimmed, tensorCell, nCopies, &exspintdim));
185   } else {
186     exspintdim = exspdim;
187   }
188   PetscCall(PetscDualSpaceSetUp(sp));
189   PetscCall(PetscDualSpaceGetDimension(sp, &spdim));
190   PetscCall(PetscDualSpaceGetInteriorDimension(sp, &spintdim));
191   PetscCheck(spdim == exspdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dual space dimension %" PetscInt_FMT ", got %" PetscInt_FMT, exspdim, spdim);
192   PetscCheck(spintdim == exspintdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dual space interior dimension %" PetscInt_FMT ", got %" PetscInt_FMT, exspintdim, spintdim);
193   key.dim        = dim;
194   key.formDegree = formDegree;
195   PetscCall(PetscDualSpaceGetOrder(sp, &key.order));
196   PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &key.continuous));
197   if (tensorCell == 2) {
198     key.tensor = 2;
199   } else {
200     PetscBool bTensor;
201 
202     PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &bTensor));
203     key.tensor = bTensor;
204   }
205   PetscCall(PetscDualSpaceLagrangeGetTrimmed(sp, &key.trimmed));
206   PetscCall(PetscInfo(NULL, "After setup:  order %" PetscInt_FMT ", trimmed %" PetscInt_FMT ", tensor %" PetscInt_FMT ", continuous %" PetscInt_FMT "\n", key.order, (PetscInt)key.trimmed, key.tensor, (PetscInt)key.continuous));
207   PetscCall(PetscHashLagPut(lagTable, key, &iter, &missing));
208   if (missing) {
209     DMPolytopeType type;
210 
211     PetscCall(DMPlexGetCellType(K, 0, &type));
212     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "New space: %s, order %" PetscInt_FMT ", trimmed %" PetscInt_FMT ", tensor %" PetscInt_FMT ", continuous %" PetscInt_FMT ", form degree %" PetscInt_FMT "\n", DMPolytopeTypes[type], order, (PetscInt)trimmed, tensorCell, (PetscInt)continuous, formDegree));
213     PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF));
214     {
215       PetscQuadrature  intNodes, allNodes;
216       Mat              intMat, allMat;
217       MatInfo          info;
218       PetscInt         i, j, nodeIdxDim, nodeVecDim, nNodes;
219       const PetscInt  *nodeIdx;
220       const PetscReal *nodeVec;
221 
222       PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
223 
224       PetscCall(PetscLagNodeIndicesGetData_Internal(lag->allNodeIndices, &nodeIdxDim, &nodeVecDim, &nNodes, &nodeIdx, &nodeVec));
225       PetscCheck(nodeVecDim == Nk, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect nodeVecDim");
226       PetscCheck(nNodes == spdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect nNodes");
227 
228       PetscCall(PetscDualSpaceGetAllData(sp, &allNodes, &allMat));
229 
230       PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All nodes:\n"));
231       PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF));
232       PetscCall(PetscQuadratureView(allNodes, PETSC_VIEWER_STDOUT_SELF));
233       PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF));
234       PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All node indices:\n"));
235       for (i = 0; i < spdim; i++) {
236         PetscCall(PetscPrintf(PETSC_COMM_SELF, "("));
237         for (j = 0; j < nodeIdxDim; j++) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT ",", nodeIdx[i * nodeIdxDim + j]));
238         PetscCall(PetscPrintf(PETSC_COMM_SELF, "): ["));
239         for (j = 0; j < nodeVecDim; j++) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g,", (double)nodeVec[i * nodeVecDim + j]));
240         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
241       }
242 
243       PetscCall(MatGetInfo(allMat, MAT_LOCAL, &info));
244       PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All matrix: %" PetscInt_FMT " nonzeros\n", (PetscInt)info.nz_used));
245 
246       PetscCall(PetscDualSpaceGetInteriorData(sp, &intNodes, &intMat));
247       if (intMat && intMat != allMat) {
248         PetscInt         intNodeIdxDim, intNodeVecDim, intNnodes;
249         const PetscInt  *intNodeIdx;
250         const PetscReal *intNodeVec;
251         PetscBool        same;
252 
253         PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior nodes:\n"));
254         PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF));
255         PetscCall(PetscQuadratureView(intNodes, PETSC_VIEWER_STDOUT_SELF));
256         PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF));
257 
258         PetscCall(MatGetInfo(intMat, MAT_LOCAL, &info));
259         PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior matrix: %" PetscInt_FMT " nonzeros\n", (PetscInt)info.nz_used));
260         PetscCall(PetscLagNodeIndicesGetData_Internal(lag->intNodeIndices, &intNodeIdxDim, &intNodeVecDim, &intNnodes, &intNodeIdx, &intNodeVec));
261         PetscCheck(intNodeIdxDim == nodeIdxDim && intNodeVecDim == nodeVecDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices not the same shale as all node indices");
262         PetscCheck(intNnodes == spintdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect interior nNodes");
263         PetscCall(PetscArraycmp(intNodeIdx, nodeIdx, nodeIdxDim * intNnodes, &same));
264         PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices not the same as start of all node indices");
265         PetscCall(PetscArraycmp(intNodeVec, nodeVec, nodeVecDim * intNnodes, &same));
266         PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node vectors not the same as start of all node vectors");
267       } else if (intMat) {
268         PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior data is the same as all data\n"));
269         PetscCheck(intNodes == allNodes, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior nodes should be the same as all nodes");
270         PetscCheck(lag->intNodeIndices == lag->allNodeIndices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices should be the same as all node indices");
271       }
272     }
273     if (dim <= 2 && spintdim) {
274       PetscInt numFaces, o;
275 
276       {
277         DMPolytopeType ct;
278         /* The number of arrangements is no longer based on the number of faces */
279         PetscCall(DMPlexGetCellType(K, 0, &ct));
280         numFaces = DMPolytopeTypeGetNumArrangements(ct) / 2;
281       }
282       for (o = -numFaces; o < numFaces; ++o) {
283         Mat symMat;
284 
285         PetscCall(PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, o, &symMat));
286         PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior node symmetry matrix for orientation %" PetscInt_FMT ":\n", o));
287         PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF));
288         PetscCall(MatView(symMat, PETSC_VIEWER_STDOUT_SELF));
289         PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF));
290         PetscCall(MatDestroy(&symMat));
291       }
292     }
293     PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF));
294   }
295   PetscCall(PetscDualSpaceDestroy(&sp));
296   PetscFunctionReturn(PETSC_SUCCESS);
297 }
298 
299 int main(int argc, char **argv)
300 {
301   PetscInt     dim;
302   PetscHashLag lagTable;
303   PetscInt     tensorCell;
304   PetscInt     order, ordermin, ordermax;
305   PetscBool    continuous;
306   PetscBool    trimmed;
307   DM           dm;
308 
309   PetscFunctionBeginUser;
310   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
311   dim        = 3;
312   tensorCell = 0;
313   continuous = PETSC_FALSE;
314   trimmed    = PETSC_FALSE;
315   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PETSCDUALSPACELAGRANGE test", "none");
316   PetscCall(PetscOptionsRangeInt("-dim", "The spatial dimension", "ex1.c", dim, &dim, NULL, 0, 3));
317   PetscCall(PetscOptionsRangeInt("-tensor", "(0) simplex (1) hypercube (2) wedge", "ex1.c", tensorCell, &tensorCell, NULL, 0, 2));
318   PetscCall(PetscOptionsBool("-continuous", "Whether the dual space has continuity", "ex1.c", continuous, &continuous, NULL));
319   PetscCall(PetscOptionsBool("-trimmed", "Whether the dual space matches a trimmed polynomial space", "ex1.c", trimmed, &trimmed, NULL));
320   PetscOptionsEnd();
321   PetscCall(PetscHashLagCreate(&lagTable));
322 
323   if (tensorCell < 2) {
324     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, (PetscBool)(tensorCell == 0)), &dm));
325   } else {
326     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DM_POLYTOPE_TRI_PRISM, &dm));
327   }
328   ordermin = trimmed ? 1 : 0;
329   ordermax = tensorCell == 2 ? 4 : tensorCell == 1 ? 3 : dim + 2;
330   for (order = ordermin; order <= ordermax; order++) {
331     PetscInt formDegree;
332 
333     for (formDegree = PetscMin(0, -dim + 1); formDegree <= dim; formDegree++) {
334       PetscInt nCopies;
335 
336       for (nCopies = 1; nCopies <= 3; nCopies++) PetscCall(testLagrange(lagTable, dm, dim, order, formDegree, trimmed, tensorCell, continuous, nCopies));
337     }
338   }
339   PetscCall(DMDestroy(&dm));
340   PetscCall(PetscHashLagDestroy(&lagTable));
341   PetscCall(PetscFinalize());
342   return 0;
343 }
344 
345 /*TEST
346 
347  test:
348    suffix: 0
349    args: -dim 0
350 
351  test:
352    suffix: 1_discontinuous_full
353    args: -dim 1 -continuous 0 -trimmed 0
354 
355  test:
356    suffix: 1_continuous_full
357    args: -dim 1 -continuous 1 -trimmed 0
358 
359  test:
360    suffix: 2_simplex_discontinuous_full
361    args: -dim 2 -tensor 0 -continuous 0 -trimmed 0
362 
363  test:
364    suffix: 2_simplex_continuous_full
365    args: -dim 2 -tensor 0 -continuous 1 -trimmed 0
366 
367  test:
368    suffix: 2_tensor_discontinuous_full
369    args: -dim 2 -tensor 1 -continuous 0 -trimmed 0
370 
371  test:
372    suffix: 2_tensor_continuous_full
373    args: -dim 2 -tensor 1 -continuous 1 -trimmed 0
374 
375  test:
376    suffix: 3_simplex_discontinuous_full
377    args: -dim 3 -tensor 0 -continuous 0 -trimmed 0
378 
379  test:
380    suffix: 3_simplex_continuous_full
381    args: -dim 3 -tensor 0 -continuous 1 -trimmed 0
382 
383  test:
384    suffix: 3_tensor_discontinuous_full
385    args: -dim 3 -tensor 1 -continuous 0 -trimmed 0
386 
387  test:
388    suffix: 3_tensor_continuous_full
389    args: -dim 3 -tensor 1 -continuous 1 -trimmed 0
390 
391  test:
392    suffix: 3_wedge_discontinuous_full
393    args: -dim 3 -tensor 2 -continuous 0 -trimmed 0
394 
395  test:
396    suffix: 3_wedge_continuous_full
397    args: -dim 3 -tensor 2 -continuous 1 -trimmed 0
398 
399  test:
400    suffix: 1_discontinuous_trimmed
401    args: -dim 1 -continuous 0 -trimmed 1
402 
403  test:
404    suffix: 1_continuous_trimmed
405    args: -dim 1 -continuous 1 -trimmed 1
406 
407  test:
408    suffix: 2_simplex_discontinuous_trimmed
409    args: -dim 2 -tensor 0 -continuous 0 -trimmed 1
410 
411  test:
412    suffix: 2_simplex_continuous_trimmed
413    args: -dim 2 -tensor 0 -continuous 1 -trimmed 1
414 
415  test:
416    suffix: 2_tensor_discontinuous_trimmed
417    args: -dim 2 -tensor 1 -continuous 0 -trimmed 1
418 
419  test:
420    suffix: 2_tensor_continuous_trimmed
421    args: -dim 2 -tensor 1 -continuous 1 -trimmed 1
422 
423  test:
424    suffix: 3_simplex_discontinuous_trimmed
425    args: -dim 3 -tensor 0 -continuous 0 -trimmed 1
426 
427  test:
428    suffix: 3_simplex_continuous_trimmed
429    args: -dim 3 -tensor 0 -continuous 1 -trimmed 1
430 
431  test:
432    suffix: 3_tensor_discontinuous_trimmed
433    args: -dim 3 -tensor 1 -continuous 0 -trimmed 1
434 
435  test:
436    suffix: 3_tensor_continuous_trimmed
437    args: -dim 3 -tensor 1 -continuous 1 -trimmed 1
438 
439  test:
440    suffix: 3_wedge_discontinuous_trimmed
441    args: -dim 3 -tensor 2 -continuous 0 -trimmed 1
442 
443  test:
444    suffix: 3_wedge_continuous_trimmed
445    args: -dim 3 -tensor 2 -continuous 1 -trimmed 1
446 
447 TEST*/
448