xref: /petsc/src/dm/impls/swarm/tests/ex1.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1ee25b4d0Sjosephpu static const char help[] = "Test initialization and migration with swarm.\n";
2461fcd32Sjosephpu 
3461fcd32Sjosephpu #include <petscdmplex.h>
4461fcd32Sjosephpu #include <petscdmswarm.h>
5461fcd32Sjosephpu 
6461fcd32Sjosephpu typedef struct {
7461fcd32Sjosephpu   PetscReal L[3]; /* Dimensions of the mesh bounding box */
88ea0bac9SZach Atkins   PetscBool setClosurePermutation;
9461fcd32Sjosephpu } AppCtx;
10461fcd32Sjosephpu 
ProcessOptions(MPI_Comm comm,AppCtx * options)11461fcd32Sjosephpu static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12461fcd32Sjosephpu {
13461fcd32Sjosephpu   PetscFunctionBeginUser;
14461fcd32Sjosephpu   PetscOptionsBegin(comm, "", "Swarm configuration options", "DMSWARM");
158ea0bac9SZach Atkins   options->setClosurePermutation = PETSC_FALSE;
168ea0bac9SZach Atkins   PetscCall(PetscOptionsBool("-set_closure_permutation", "Set the closure permutation to tensor ordering", NULL, options->setClosurePermutation, &options->setClosurePermutation, NULL));
17461fcd32Sjosephpu   PetscOptionsEnd();
18461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
19461fcd32Sjosephpu }
20461fcd32Sjosephpu 
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * user)21461fcd32Sjosephpu static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
22461fcd32Sjosephpu {
23461fcd32Sjosephpu   PetscReal low[3], high[3];
24461fcd32Sjosephpu   PetscInt  cdim, d;
25461fcd32Sjosephpu 
26461fcd32Sjosephpu   PetscFunctionBeginUser;
27461fcd32Sjosephpu   PetscCall(DMCreate(comm, dm));
28461fcd32Sjosephpu   PetscCall(DMSetType(*dm, DMPLEX));
29461fcd32Sjosephpu   PetscCall(DMSetFromOptions(*dm));
30461fcd32Sjosephpu   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
318ea0bac9SZach Atkins   if (user->setClosurePermutation) {
328ea0bac9SZach Atkins     DM cdm;
338ea0bac9SZach Atkins 
348ea0bac9SZach Atkins     // -- Set tensor permutation
358ea0bac9SZach Atkins     PetscCall(DMGetCoordinateDM(*dm, &cdm));
368ea0bac9SZach Atkins     PetscCall(DMPlexSetClosurePermutationTensor(*dm, PETSC_DETERMINE, NULL));
378ea0bac9SZach Atkins     PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
388ea0bac9SZach Atkins   }
39461fcd32Sjosephpu   PetscCall(DMGetCoordinateDim(*dm, &cdim));
40461fcd32Sjosephpu   PetscCall(DMGetBoundingBox(*dm, low, high));
41ee25b4d0Sjosephpu   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "dim: %" PetscInt_FMT "\n", cdim));
42461fcd32Sjosephpu   for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d];
43461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
44461fcd32Sjosephpu }
45461fcd32Sjosephpu 
46461fcd32Sjosephpu /*
47461fcd32Sjosephpu   This function initializes all particles on rank 0.
48461fcd32Sjosephpu   They are sent to other ranks to test migration across non nearest neighbors
49461fcd32Sjosephpu */
CreateSwarm(DM dm,DM * sw,AppCtx * user)50461fcd32Sjosephpu static PetscErrorCode CreateSwarm(DM dm, DM *sw, AppCtx *user)
51461fcd32Sjosephpu {
5219307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
53461fcd32Sjosephpu   PetscInt      particleInitSize = 10;
54461fcd32Sjosephpu   PetscReal    *coords, upper[3], lower[3];
5519307e5cSMatthew G. Knepley   PetscInt     *swarm_cellid, Np, dim, Nfc;
56ee25b4d0Sjosephpu   PetscMPIInt   rank, size;
57461fcd32Sjosephpu   MPI_Comm      comm;
5819307e5cSMatthew G. Knepley   const char  **coordFields, *cellid;
59461fcd32Sjosephpu 
60461fcd32Sjosephpu   PetscFunctionBegin;
61461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
62ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_rank(comm, &rank));
63ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_size(comm, &size));
64461fcd32Sjosephpu   PetscCall(DMGetBoundingBox(dm, lower, upper));
65461fcd32Sjosephpu   PetscCall(DMCreate(PETSC_COMM_WORLD, sw));
66461fcd32Sjosephpu   PetscCall(DMGetDimension(dm, &dim));
67461fcd32Sjosephpu   PetscCall(DMSetType(*sw, DMSWARM));
68461fcd32Sjosephpu   PetscCall(DMSetDimension(*sw, dim));
69461fcd32Sjosephpu   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
70461fcd32Sjosephpu   PetscCall(DMSwarmSetCellDM(*sw, dm));
71461fcd32Sjosephpu   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
72461fcd32Sjosephpu   PetscCall(DMSwarmSetLocalSizes(*sw, rank == 0 ? particleInitSize : 0, 0));
73461fcd32Sjosephpu   PetscCall(DMSetFromOptions(*sw));
7419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(*sw, &celldm));
7519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
7619307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
7719307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
7819307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, coordFields[0], NULL, NULL, (void **)&coords));
7919307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
80461fcd32Sjosephpu   PetscCall(DMSwarmGetLocalSize(*sw, &Np));
81461fcd32Sjosephpu   for (PetscInt p = 0; p < Np; ++p) {
82*ac530a7eSPierre Jolivet     for (PetscInt d = 0; d < dim; ++d) coords[p * dim + d] = 0.5 * (upper[d] - lower[d]);
83461fcd32Sjosephpu     coords[p * dim + 1] = (upper[1] - lower[1]) / particleInitSize * p + lower[1];
8419307e5cSMatthew G. Knepley     swarm_cellid[p]     = 0;
85461fcd32Sjosephpu   }
8619307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, coordFields[0], NULL, NULL, (void **)&coords));
8719307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
88461fcd32Sjosephpu   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
89461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
90461fcd32Sjosephpu }
91461fcd32Sjosephpu 
92461fcd32Sjosephpu /*
93461fcd32Sjosephpu   Configure the swarm on rank 0 with all particles
94461fcd32Sjosephpu   located there, then migrate where they need to be
95461fcd32Sjosephpu */
CheckMigrate(DM sw)96461fcd32Sjosephpu static PetscErrorCode CheckMigrate(DM sw)
97461fcd32Sjosephpu {
98461fcd32Sjosephpu   Vec       preMigrate, postMigrate, tmp;
99461fcd32Sjosephpu   PetscInt  preSize, postSize;
100461fcd32Sjosephpu   PetscReal prenorm, postnorm;
101461fcd32Sjosephpu 
102461fcd32Sjosephpu   PetscFunctionBeginUser;
103461fcd32Sjosephpu   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
104461fcd32Sjosephpu   PetscCall(VecDuplicate(tmp, &preMigrate));
105461fcd32Sjosephpu   PetscCall(VecCopy(tmp, preMigrate));
106461fcd32Sjosephpu   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
107461fcd32Sjosephpu   PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
108461fcd32Sjosephpu   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
109461fcd32Sjosephpu   PetscCall(VecDuplicate(tmp, &postMigrate));
110461fcd32Sjosephpu   PetscCall(VecCopy(tmp, postMigrate));
111461fcd32Sjosephpu   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
112461fcd32Sjosephpu   PetscCall(VecGetSize(preMigrate, &preSize));
113461fcd32Sjosephpu   PetscCall(VecGetSize(postMigrate, &postSize));
114e978a55eSPierre Jolivet   PetscCheck(preSize == postSize, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Particles either lost or duplicated. Pre migrate global size %" PetscInt_FMT " != Post migrate size %" PetscInt_FMT, preSize, postSize);
115461fcd32Sjosephpu   PetscCall(VecNorm(preMigrate, NORM_2, &prenorm));
116461fcd32Sjosephpu   PetscCall(VecNorm(postMigrate, NORM_2, &postnorm));
117ee25b4d0Sjosephpu   PetscCheck(PetscAbsReal(prenorm - postnorm) < 100. * PETSC_MACHINE_EPSILON, PETSC_COMM_SELF, PETSC_ERR_COR, "Particle coordinates corrupted in migrate with abs(norm(pre) - norm(post)) = %.16g", PetscAbsReal(prenorm - postnorm));
118461fcd32Sjosephpu   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Migrate check passes.\n"));
119461fcd32Sjosephpu   PetscCall(VecDestroy(&preMigrate));
120461fcd32Sjosephpu   PetscCall(VecDestroy(&postMigrate));
121461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
122461fcd32Sjosephpu }
123461fcd32Sjosephpu 
124461fcd32Sjosephpu /*
125461fcd32Sjosephpu   Checks added points persist on migrate
126461fcd32Sjosephpu */
CheckPointInsertion(DM sw)127461fcd32Sjosephpu static PetscErrorCode CheckPointInsertion(DM sw)
128461fcd32Sjosephpu {
129ee25b4d0Sjosephpu   PetscInt    Np_pre, Np_post;
130ee25b4d0Sjosephpu   PetscMPIInt rank, size;
131461fcd32Sjosephpu   MPI_Comm    comm;
132461fcd32Sjosephpu 
133461fcd32Sjosephpu   PetscFunctionBeginUser;
134461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
135ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_rank(comm, &rank));
136ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_size(comm, &size));
137461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Basic point insertion check...\n"));
138461fcd32Sjosephpu   PetscCall(DMSwarmGetSize(sw, &Np_pre));
139461fcd32Sjosephpu   if (rank == 0) PetscCall(DMSwarmAddPoint(sw));
140461fcd32Sjosephpu   PetscCall(DMSwarmGetSize(sw, &Np_post));
141461fcd32Sjosephpu   PetscCheck(Np_post == (Np_pre + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Particle insertion failed. Global size pre insertion: %" PetscInt_FMT " global size post insertion: %" PetscInt_FMT, Np_pre, Np_post);
142461fcd32Sjosephpu   PetscCall(CheckMigrate(sw));
143461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Basic point insertion check passes.\n"));
144461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
145461fcd32Sjosephpu }
146461fcd32Sjosephpu 
147461fcd32Sjosephpu /*
148461fcd32Sjosephpu   Checks tie breaking works properly when a particle
149461fcd32Sjosephpu   is located at a shared boundary. The higher rank should
150baca6076SPierre Jolivet   receive the particle while the lower rank deletes it.
151461fcd32Sjosephpu 
152461fcd32Sjosephpu   TODO: Currently only works for 2 procs.
153461fcd32Sjosephpu */
CheckPointInsertion_Boundary(DM sw)154461fcd32Sjosephpu static PetscErrorCode CheckPointInsertion_Boundary(DM sw)
155461fcd32Sjosephpu {
156ee25b4d0Sjosephpu   PetscInt    Np_loc_pre, Np_loc_post, dim;
157ee25b4d0Sjosephpu   PetscMPIInt rank, size;
158461fcd32Sjosephpu   PetscReal   lbox_low[3], lbox_high[3], gbox_low[3], gbox_high[3];
159461fcd32Sjosephpu   MPI_Comm    comm;
160461fcd32Sjosephpu   DM          cdm;
161461fcd32Sjosephpu 
162461fcd32Sjosephpu   PetscFunctionBeginUser;
163461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
164ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_rank(comm, &rank));
165ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_size(comm, &size));
166461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Rank boundary point insertion check...\n"));
167461fcd32Sjosephpu   PetscCall(DMSwarmGetCellDM(sw, &cdm));
168461fcd32Sjosephpu   PetscCall(DMGetDimension(cdm, &dim));
169461fcd32Sjosephpu   PetscCall(DMGetBoundingBox(cdm, gbox_low, gbox_high));
170461fcd32Sjosephpu   if (rank == 0) {
171461fcd32Sjosephpu     PetscReal *coords;
172ee25b4d0Sjosephpu     PetscInt   adjacentdim = 0, Np;
173461fcd32Sjosephpu 
174461fcd32Sjosephpu     PetscCall(DMGetLocalBoundingBox(cdm, lbox_low, lbox_high));
175461fcd32Sjosephpu     // find a face that belongs to the neighbor.
176461fcd32Sjosephpu     for (PetscInt d = 0; d < dim; ++d) {
177461fcd32Sjosephpu       if ((gbox_high[d] - lbox_high[d]) != 0.) adjacentdim = d;
178461fcd32Sjosephpu     }
179461fcd32Sjosephpu     PetscCall(DMSwarmAddPoint(sw));
180461fcd32Sjosephpu     PetscCall(DMSwarmGetLocalSize(sw, &Np));
181461fcd32Sjosephpu     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
182461fcd32Sjosephpu     for (PetscInt d = 0; d < dim; ++d) coords[(Np - 1) * dim + d] = 0.5 * (lbox_high[d] - lbox_low[d]);
183461fcd32Sjosephpu     coords[(Np - 1) * dim + adjacentdim] = lbox_high[adjacentdim];
184461fcd32Sjosephpu     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
185461fcd32Sjosephpu   }
186461fcd32Sjosephpu   PetscCall(DMSwarmGetLocalSize(sw, &Np_loc_pre));
187461fcd32Sjosephpu   PetscCall(CheckMigrate(sw));
188461fcd32Sjosephpu   PetscCall(DMSwarmGetLocalSize(sw, &Np_loc_post));
189ee25b4d0Sjosephpu   if (rank == 0) PetscCheck(Np_loc_pre == (Np_loc_post + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Migration tie breaking failed on rank %d. Particle on boundary not sent.", rank);
190baca6076SPierre Jolivet   if (rank == 1) PetscCheck(Np_loc_pre == (Np_loc_post - 1), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Migration tie breaking failed on rank %d. Particle on boundary not received.", rank);
191461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Rank boundary point insertion check passes.\n"));
192461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
193461fcd32Sjosephpu }
194461fcd32Sjosephpu 
main(int argc,char ** argv)195461fcd32Sjosephpu int main(int argc, char **argv)
196461fcd32Sjosephpu {
197461fcd32Sjosephpu   DM       dm, sw;
198461fcd32Sjosephpu   MPI_Comm comm;
199461fcd32Sjosephpu   AppCtx   user;
200461fcd32Sjosephpu 
201461fcd32Sjosephpu   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
202461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
203461fcd32Sjosephpu   PetscCall(ProcessOptions(comm, &user));
204461fcd32Sjosephpu   PetscCall(CreateMesh(comm, &dm, &user));
205461fcd32Sjosephpu   PetscCall(CreateSwarm(dm, &sw, &user));
206461fcd32Sjosephpu   PetscCall(CheckMigrate(sw));
207461fcd32Sjosephpu   PetscCall(CheckPointInsertion(sw));
208461fcd32Sjosephpu   PetscCall(CheckPointInsertion_Boundary(sw));
209461fcd32Sjosephpu   PetscCall(DMDestroy(&sw));
210461fcd32Sjosephpu   PetscCall(DMDestroy(&dm));
211ee25b4d0Sjosephpu   PetscCall(PetscFinalize());
212461fcd32Sjosephpu   return 0;
213461fcd32Sjosephpu }
214461fcd32Sjosephpu 
215461fcd32Sjosephpu /*TEST
216ee25b4d0Sjosephpu 
217461fcd32Sjosephpu   # Swarm does not handle complex or quad
218461fcd32Sjosephpu   build:
219461fcd32Sjosephpu     requires: !complex double
220ee25b4d0Sjosephpu   # swarm_migrate_hash and swarm_migrate_scan test swarm migration against point location types
221ee25b4d0Sjosephpu   # with a distributed mesh where ranks overlap by 1. Points in the shared boundary should
222ee25b4d0Sjosephpu   # be sent to the process which has the highest rank that has that portion of the domain.
223461fcd32Sjosephpu   test:
224953fc092SRylanor     suffix: swarm_migrate_vec_hip_scan
225953fc092SRylanor     nsize: 2
226953fc092SRylanor     requires: hip
227953fc092SRylanor     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
228953fc092SRylanor           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
229953fc092SRylanor           -dm_plex_hash_location false -dm_vec_type hip -dm_plex_hash_location false
230953fc092SRylanor   test:
231953fc092SRylanor     suffix: swarm_migrate_vec_hip_hash
232953fc092SRylanor     nsize: 2
233953fc092SRylanor     requires: hip
234953fc092SRylanor     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
235953fc092SRylanor           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
236953fc092SRylanor           -dm_plex_hash_location false -dm_vec_type hip -dm_plex_hash_location true
237953fc092SRylanor   test:
238953fc092SRylanor     suffix: swarm_migrate_vec_hip_hash_tensor_permutation
239953fc092SRylanor     nsize: 2
240953fc092SRylanor     requires: hip
241953fc092SRylanor     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
242953fc092SRylanor           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
243953fc092SRylanor           -dm_plex_hash_location false -dm_vec_type hip -dm_plex_hash_location true\
244953fc092SRylanor           -set_closure_permutation
245953fc092SRylanor   test:
246461fcd32Sjosephpu     suffix: swarm_migrate_hash
247461fcd32Sjosephpu     nsize: 2
248ee25b4d0Sjosephpu     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
249ee25b4d0Sjosephpu           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
250ee25b4d0Sjosephpu           -dm_plex_hash_location true
251461fcd32Sjosephpu     filter: grep -v marker | grep -v atomic | grep -v usage
252461fcd32Sjosephpu   test:
2538ea0bac9SZach Atkins     suffix: swarm_migrate_hash_tensor_permutation
2548ea0bac9SZach Atkins     nsize: 2
2558ea0bac9SZach Atkins     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
2568ea0bac9SZach Atkins           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
2578ea0bac9SZach Atkins           -dm_plex_hash_location true -set_closure_permutation
2588ea0bac9SZach Atkins     filter: grep -v marker | grep -v atomic | grep -v usage
2598ea0bac9SZach Atkins   test:
260461fcd32Sjosephpu     suffix: swarm_migrate_scan
261461fcd32Sjosephpu     nsize: 2
262ee25b4d0Sjosephpu     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
263ee25b4d0Sjosephpu           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
264ee25b4d0Sjosephpu           -dm_plex_hash_location false
265461fcd32Sjosephpu     filter: grep -v marker | grep -v atomic | grep -v usage
2668ea0bac9SZach Atkins   test:
2678ea0bac9SZach Atkins     suffix: swarm_migrate_scan_tensor_permutation
2688ea0bac9SZach Atkins     nsize: 2
2698ea0bac9SZach Atkins     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
2708ea0bac9SZach Atkins           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
2718ea0bac9SZach Atkins           -dm_plex_hash_location false -set_closure_permutation
2728ea0bac9SZach Atkins     filter: grep -v marker | grep -v atomic | grep -v usage
273461fcd32Sjosephpu TEST*/
274