xref: /petsc/src/dm/impls/plex/tests/ex25.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test DMCreateLocalVector_Plex, DMPlexGetCellFields and DMPlexRestoreCellFields work properly for 0 fields/cells/DS dimension\n\n";
2*c4762a1bSJed Brown static char FILENAME[] = "ex25.c";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscdmplex.h>
5*c4762a1bSJed Brown #include <petscds.h>
6*c4762a1bSJed Brown #include <petscsnes.h>
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown typedef struct {
10*c4762a1bSJed Brown   PetscInt  test;
11*c4762a1bSJed Brown   PetscInt  dim;                          /* The topological mesh dimension */
12*c4762a1bSJed Brown   PetscInt  faces[3];                     /* Number of faces per dimension */
13*c4762a1bSJed Brown   PetscBool simplex;                      /* Use simplices or hexes */
14*c4762a1bSJed Brown   PetscBool interpolate;                  /* Interpolate mesh */
15*c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
16*c4762a1bSJed Brown } AppCtx;
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
19*c4762a1bSJed Brown {
20*c4762a1bSJed Brown   PetscInt dim;
21*c4762a1bSJed Brown   PetscErrorCode ierr;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   PetscFunctionBegin;
24*c4762a1bSJed Brown   options->test         = 0;
25*c4762a1bSJed Brown   options->dim          = 3;
26*c4762a1bSJed Brown   options->simplex      = PETSC_TRUE;
27*c4762a1bSJed Brown   options->interpolate  = PETSC_FALSE;
28*c4762a1bSJed Brown   options->filename[0]  = '\0';
29*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Zero-sized DMPlexGetCellFields Test Options", "DMPLEX");CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-test", "Test to run", FILENAME, options->test, &options->test, NULL,0);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", FILENAME, options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
32*c4762a1bSJed Brown   if (options->dim > 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "dimension set to %d, must be <= 3", options->dim);
33*c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", FILENAME, options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = PetscOptionsBool("-interpolate", "Interpolate the mesh", FILENAME, options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "The mesh file", FILENAME, options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
36*c4762a1bSJed Brown   options->faces[0] = 1; options->faces[1] = 1; options->faces[2] = 1;
37*c4762a1bSJed Brown   dim = options->dim;
38*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", FILENAME, options->faces, &dim, NULL);CHKERRQ(ierr);
39*c4762a1bSJed Brown   if (dim) options->dim = dim;
40*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
41*c4762a1bSJed Brown   PetscFunctionReturn(0);
42*c4762a1bSJed Brown }
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *options, DM *dm)
45*c4762a1bSJed Brown {
46*c4762a1bSJed Brown   PetscInt       dim          = options->dim;
47*c4762a1bSJed Brown   PetscInt      *faces        = options->faces;
48*c4762a1bSJed Brown   PetscBool      simplex      = options->simplex;
49*c4762a1bSJed Brown   PetscBool      interpolate  = options->interpolate;
50*c4762a1bSJed Brown   const char    *filename     = options->filename;
51*c4762a1bSJed Brown   size_t         len;
52*c4762a1bSJed Brown   PetscMPIInt    rank;
53*c4762a1bSJed Brown   PetscErrorCode ierr;
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   PetscFunctionBegin;
56*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
58*c4762a1bSJed Brown   if (len) {
59*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
60*c4762a1bSJed Brown     ierr = DMGetDimension(*dm, &options->dim);CHKERRQ(ierr);
61*c4762a1bSJed Brown   } else {
62*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, interpolate, dm);CHKERRQ(ierr);
63*c4762a1bSJed Brown   }
64*c4762a1bSJed Brown   PetscFunctionReturn(0);
65*c4762a1bSJed Brown }
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown /* no discretization is given so DMGetNumFields yields 0 */
68*c4762a1bSJed Brown static PetscErrorCode test0(DM dm, AppCtx *options)
69*c4762a1bSJed Brown {
70*c4762a1bSJed Brown   Vec            locX;
71*c4762a1bSJed Brown   PetscErrorCode ierr;
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   PetscFunctionBegin;
74*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
76*c4762a1bSJed Brown   PetscFunctionReturn(0);
77*c4762a1bSJed Brown }
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown /* no discretization is given so DMGetNumFields and PetscDSGetTotalDimension yield 0 */
80*c4762a1bSJed Brown static PetscErrorCode test1(DM dm, AppCtx *options)
81*c4762a1bSJed Brown {
82*c4762a1bSJed Brown   IS             cells;
83*c4762a1bSJed Brown   Vec            locX, locX_t, locA;
84*c4762a1bSJed Brown   PetscScalar    *u, *u_t, *a;
85*c4762a1bSJed Brown   PetscErrorCode ierr;
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   PetscFunctionBegin;
88*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &cells);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locA);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = DMPlexGetCellFields(    dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locA);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = ISDestroy(&cells);CHKERRQ(ierr);
98*c4762a1bSJed Brown   PetscFunctionReturn(0);
99*c4762a1bSJed Brown }
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown /* no discretization is given so DMGetNumFields and PetscDSGetTotalDimension yield 0 */
102*c4762a1bSJed Brown static PetscErrorCode test2(DM dm, AppCtx *options)
103*c4762a1bSJed Brown {
104*c4762a1bSJed Brown   IS             cells;
105*c4762a1bSJed Brown   Vec            locX, locX_t, locA;
106*c4762a1bSJed Brown   PetscScalar    *u, *u_t, *a;
107*c4762a1bSJed Brown   PetscMPIInt    rank;
108*c4762a1bSJed Brown   PetscErrorCode ierr;
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   PetscFunctionBegin;
111*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF, rank ? 0 : 1, 0, 1, &cells);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locA);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = DMPlexGetCellFields(    dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locA);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = ISDestroy(&cells);CHKERRQ(ierr);
122*c4762a1bSJed Brown   PetscFunctionReturn(0);
123*c4762a1bSJed Brown }
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown static PetscErrorCode test3(DM dm, AppCtx *options)
126*c4762a1bSJed Brown {
127*c4762a1bSJed Brown   PetscDS        ds;
128*c4762a1bSJed Brown   PetscFE        fe;
129*c4762a1bSJed Brown   PetscErrorCode ierr;
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown   PetscFunctionBegin;
132*c4762a1bSJed Brown   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = PetscDSSetDiscretization(ds, 0, (PetscObject)fe);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
136*c4762a1bSJed Brown   ierr = test1(dm, options);CHKERRQ(ierr);
137*c4762a1bSJed Brown   PetscFunctionReturn(0);
138*c4762a1bSJed Brown }
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown static PetscErrorCode test4(DM dm, AppCtx *options)
141*c4762a1bSJed Brown {
142*c4762a1bSJed Brown   PetscFE        fe;
143*c4762a1bSJed Brown   PetscErrorCode ierr;
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   PetscFunctionBegin;
146*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = test2(dm, options);CHKERRQ(ierr);
151*c4762a1bSJed Brown   PetscFunctionReturn(0);
152*c4762a1bSJed Brown }
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown static PetscErrorCode test5(DM dm, AppCtx *options)
155*c4762a1bSJed Brown {
156*c4762a1bSJed Brown   IS             cells;
157*c4762a1bSJed Brown   Vec            locX, locX_t, locA;
158*c4762a1bSJed Brown   PetscScalar    *u, *u_t, *a;
159*c4762a1bSJed Brown   PetscErrorCode ierr;
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown   PetscFunctionBegin;
162*c4762a1bSJed Brown   locX_t = NULL;
163*c4762a1bSJed Brown   locA = NULL;
164*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &cells);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = DMPlexGetCellFields(    dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = ISDestroy(&cells);CHKERRQ(ierr);
170*c4762a1bSJed Brown   PetscFunctionReturn(0);
171*c4762a1bSJed Brown }
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown static PetscErrorCode test6(DM dm, AppCtx *options)
174*c4762a1bSJed Brown {
175*c4762a1bSJed Brown   IS             cells;
176*c4762a1bSJed Brown   Vec            locX, locX_t, locA;
177*c4762a1bSJed Brown   PetscScalar    *u, *u_t, *a;
178*c4762a1bSJed Brown   PetscMPIInt    rank;
179*c4762a1bSJed Brown   PetscErrorCode ierr;
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   PetscFunctionBegin;
182*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
183*c4762a1bSJed Brown   locX_t = NULL;
184*c4762a1bSJed Brown   locA = NULL;
185*c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF, rank ? 0 : 1, 0, 1, &cells);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = DMPlexGetCellFields(    dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = ISDestroy(&cells);CHKERRQ(ierr);
191*c4762a1bSJed Brown   PetscFunctionReturn(0);
192*c4762a1bSJed Brown }
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown static PetscErrorCode test7(DM dm, AppCtx *options)
195*c4762a1bSJed Brown {
196*c4762a1bSJed Brown   PetscFE        fe;
197*c4762a1bSJed Brown   PetscErrorCode ierr;
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown   PetscFunctionBegin;
200*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = test5(dm, options);CHKERRQ(ierr);
205*c4762a1bSJed Brown   PetscFunctionReturn(0);
206*c4762a1bSJed Brown }
207*c4762a1bSJed Brown 
208*c4762a1bSJed Brown static PetscErrorCode test8(DM dm, AppCtx *options)
209*c4762a1bSJed Brown {
210*c4762a1bSJed Brown   PetscFE        fe;
211*c4762a1bSJed Brown   PetscErrorCode ierr;
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown   PetscFunctionBegin;
214*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr);
215*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = test6(dm, options);CHKERRQ(ierr);
219*c4762a1bSJed Brown   PetscFunctionReturn(0);
220*c4762a1bSJed Brown }
221*c4762a1bSJed Brown 
222*c4762a1bSJed Brown int main(int argc, char **argv)
223*c4762a1bSJed Brown {
224*c4762a1bSJed Brown   MPI_Comm       comm;
225*c4762a1bSJed Brown   DM             dm;
226*c4762a1bSJed Brown   AppCtx         options;
227*c4762a1bSJed Brown   PetscErrorCode ierr;
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
230*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
231*c4762a1bSJed Brown   ierr = ProcessOptions(comm, &options);CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = CreateMesh(comm, &options, &dm);CHKERRQ(ierr);
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown   switch (options.test) {
235*c4762a1bSJed Brown     case 0: ierr = test0(dm, &options);CHKERRQ(ierr); break;
236*c4762a1bSJed Brown     case 1: ierr = test1(dm, &options);CHKERRQ(ierr); break;
237*c4762a1bSJed Brown     case 2: ierr = test2(dm, &options);CHKERRQ(ierr); break;
238*c4762a1bSJed Brown     case 3: ierr = test3(dm, &options);CHKERRQ(ierr); break;
239*c4762a1bSJed Brown     case 4: ierr = test4(dm, &options);CHKERRQ(ierr); break;
240*c4762a1bSJed Brown     case 5: ierr = test5(dm, &options);CHKERRQ(ierr); break;
241*c4762a1bSJed Brown     case 6: ierr = test6(dm, &options);CHKERRQ(ierr); break;
242*c4762a1bSJed Brown     case 7: ierr = test7(dm, &options);CHKERRQ(ierr); break;
243*c4762a1bSJed Brown     case 8: ierr = test8(dm, &options);CHKERRQ(ierr); break;
244*c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No such test: %D", options.test);
245*c4762a1bSJed Brown   }
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
248*c4762a1bSJed Brown   ierr = PetscFinalize();
249*c4762a1bSJed Brown   return ierr;
250*c4762a1bSJed Brown }
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown /*TEST
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   test:
255*c4762a1bSJed Brown     suffix: 0
256*c4762a1bSJed Brown     requires: ctetgen
257*c4762a1bSJed Brown     args: -test 0
258*c4762a1bSJed Brown   test:
259*c4762a1bSJed Brown     suffix: 1
260*c4762a1bSJed Brown     requires: ctetgen
261*c4762a1bSJed Brown     args: -test 1
262*c4762a1bSJed Brown   test:
263*c4762a1bSJed Brown     suffix: 2
264*c4762a1bSJed Brown     requires: ctetgen
265*c4762a1bSJed Brown     args: -test 2
266*c4762a1bSJed Brown   test:
267*c4762a1bSJed Brown     suffix: 3
268*c4762a1bSJed Brown     requires: ctetgen
269*c4762a1bSJed Brown     args: -test 3
270*c4762a1bSJed Brown   test:
271*c4762a1bSJed Brown     suffix: 4
272*c4762a1bSJed Brown     requires: ctetgen
273*c4762a1bSJed Brown     args: -test 4
274*c4762a1bSJed Brown   test:
275*c4762a1bSJed Brown     suffix: 5
276*c4762a1bSJed Brown     requires: ctetgen
277*c4762a1bSJed Brown     args: -test 5
278*c4762a1bSJed Brown   test:
279*c4762a1bSJed Brown     suffix: 6
280*c4762a1bSJed Brown     requires: ctetgen
281*c4762a1bSJed Brown     args: -test 6
282*c4762a1bSJed Brown   test:
283*c4762a1bSJed Brown     suffix: 7
284*c4762a1bSJed Brown     requires: ctetgen
285*c4762a1bSJed Brown     args: -test 7
286*c4762a1bSJed Brown   test:
287*c4762a1bSJed Brown     suffix: 8
288*c4762a1bSJed Brown     requires: ctetgen
289*c4762a1bSJed Brown     args: -test 8
290*c4762a1bSJed Brown   test:
291*c4762a1bSJed Brown     suffix: 9
292*c4762a1bSJed Brown     requires: ctetgen
293*c4762a1bSJed Brown     nsize: 2
294*c4762a1bSJed Brown     args: -test 1
295*c4762a1bSJed Brown   test:
296*c4762a1bSJed Brown     suffix: 10
297*c4762a1bSJed Brown     requires: ctetgen
298*c4762a1bSJed Brown     nsize: 2
299*c4762a1bSJed Brown     args: -test 2
300*c4762a1bSJed Brown   test:
301*c4762a1bSJed Brown     suffix: 11
302*c4762a1bSJed Brown     requires: ctetgen
303*c4762a1bSJed Brown     nsize: 2
304*c4762a1bSJed Brown     args: -test 3
305*c4762a1bSJed Brown   test:
306*c4762a1bSJed Brown     suffix: 12
307*c4762a1bSJed Brown     requires: ctetgen
308*c4762a1bSJed Brown     nsize: 2
309*c4762a1bSJed Brown     args: -test 4
310*c4762a1bSJed Brown 
311*c4762a1bSJed Brown TEST*/
312*c4762a1bSJed Brown 
313