xref: /petsc/src/dm/impls/swarm/tests/ex1.c (revision 8ea0bac9d10baee98247a9f5030ae4007889014c)
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 */
8*8ea0bac9SZach Atkins   PetscBool setClosurePermutation;
9461fcd32Sjosephpu } AppCtx;
10461fcd32Sjosephpu 
11461fcd32Sjosephpu static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12461fcd32Sjosephpu {
13461fcd32Sjosephpu   PetscFunctionBeginUser;
14461fcd32Sjosephpu   PetscOptionsBegin(comm, "", "Swarm configuration options", "DMSWARM");
15*8ea0bac9SZach Atkins   options->setClosurePermutation = PETSC_FALSE;
16*8ea0bac9SZach 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 
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"));
31*8ea0bac9SZach Atkins   if (user->setClosurePermutation) {
32*8ea0bac9SZach Atkins     DM cdm;
33*8ea0bac9SZach Atkins 
34*8ea0bac9SZach Atkins     // -- Set tensor permutation
35*8ea0bac9SZach Atkins     PetscCall(DMGetCoordinateDM(*dm, &cdm));
36*8ea0bac9SZach Atkins     PetscCall(DMPlexSetClosurePermutationTensor(*dm, PETSC_DETERMINE, NULL));
37*8ea0bac9SZach Atkins     PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
38*8ea0bac9SZach 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 */
50461fcd32Sjosephpu static PetscErrorCode CreateSwarm(DM dm, DM *sw, AppCtx *user)
51461fcd32Sjosephpu {
52461fcd32Sjosephpu   PetscInt    particleInitSize = 10;
53461fcd32Sjosephpu   PetscReal  *coords, upper[3], lower[3];
54ee25b4d0Sjosephpu   PetscInt   *cellid, Np, dim;
55ee25b4d0Sjosephpu   PetscMPIInt rank, size;
56461fcd32Sjosephpu   MPI_Comm    comm;
57461fcd32Sjosephpu 
58461fcd32Sjosephpu   PetscFunctionBegin;
59461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
60ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_rank(comm, &rank));
61ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_size(comm, &size));
62461fcd32Sjosephpu   PetscCall(DMGetBoundingBox(dm, lower, upper));
63461fcd32Sjosephpu   PetscCall(DMCreate(PETSC_COMM_WORLD, sw));
64461fcd32Sjosephpu   PetscCall(DMGetDimension(dm, &dim));
65461fcd32Sjosephpu   PetscCall(DMSetType(*sw, DMSWARM));
66461fcd32Sjosephpu   PetscCall(DMSetDimension(*sw, dim));
67461fcd32Sjosephpu   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
68461fcd32Sjosephpu   PetscCall(DMSwarmSetCellDM(*sw, dm));
69461fcd32Sjosephpu   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
70461fcd32Sjosephpu   PetscCall(DMSwarmSetLocalSizes(*sw, rank == 0 ? particleInitSize : 0, 0));
71461fcd32Sjosephpu   PetscCall(DMSetFromOptions(*sw));
72461fcd32Sjosephpu   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
73461fcd32Sjosephpu   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
74461fcd32Sjosephpu   PetscCall(DMSwarmGetLocalSize(*sw, &Np));
75461fcd32Sjosephpu   for (PetscInt p = 0; p < Np; ++p) {
76461fcd32Sjosephpu     for (PetscInt d = 0; d < dim; ++d) { coords[p * dim + d] = 0.5 * (upper[d] - lower[d]); }
77461fcd32Sjosephpu     coords[p * dim + 1] = (upper[1] - lower[1]) / particleInitSize * p + lower[1];
78461fcd32Sjosephpu     cellid[p]           = 0;
79461fcd32Sjosephpu   }
80461fcd32Sjosephpu   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
81461fcd32Sjosephpu   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
82461fcd32Sjosephpu   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
83461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
84461fcd32Sjosephpu }
85461fcd32Sjosephpu 
86461fcd32Sjosephpu /*
87461fcd32Sjosephpu   Configure the swarm on rank 0 with all particles
88461fcd32Sjosephpu   located there, then migrate where they need to be
89461fcd32Sjosephpu */
90461fcd32Sjosephpu static PetscErrorCode CheckMigrate(DM sw)
91461fcd32Sjosephpu {
92461fcd32Sjosephpu   Vec       preMigrate, postMigrate, tmp;
93461fcd32Sjosephpu   PetscInt  preSize, postSize;
94461fcd32Sjosephpu   PetscReal prenorm, postnorm;
95461fcd32Sjosephpu 
96461fcd32Sjosephpu   PetscFunctionBeginUser;
97461fcd32Sjosephpu   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
98461fcd32Sjosephpu   PetscCall(VecDuplicate(tmp, &preMigrate));
99461fcd32Sjosephpu   PetscCall(VecCopy(tmp, preMigrate));
100461fcd32Sjosephpu   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
101461fcd32Sjosephpu   PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
102461fcd32Sjosephpu   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
103461fcd32Sjosephpu   PetscCall(VecDuplicate(tmp, &postMigrate));
104461fcd32Sjosephpu   PetscCall(VecCopy(tmp, postMigrate));
105461fcd32Sjosephpu   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &tmp));
106461fcd32Sjosephpu   PetscCall(VecGetSize(preMigrate, &preSize));
107461fcd32Sjosephpu   PetscCall(VecGetSize(postMigrate, &postSize));
108461fcd32Sjosephpu   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);
109461fcd32Sjosephpu   PetscCall(VecNorm(preMigrate, NORM_2, &prenorm));
110461fcd32Sjosephpu   PetscCall(VecNorm(postMigrate, NORM_2, &postnorm));
111ee25b4d0Sjosephpu   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));
112461fcd32Sjosephpu   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Migrate check passes.\n"));
113461fcd32Sjosephpu   PetscCall(VecDestroy(&preMigrate));
114461fcd32Sjosephpu   PetscCall(VecDestroy(&postMigrate));
115461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
116461fcd32Sjosephpu }
117461fcd32Sjosephpu 
118461fcd32Sjosephpu /*
119461fcd32Sjosephpu   Checks added points persist on migrate
120461fcd32Sjosephpu */
121461fcd32Sjosephpu static PetscErrorCode CheckPointInsertion(DM sw)
122461fcd32Sjosephpu {
123ee25b4d0Sjosephpu   PetscInt    Np_pre, Np_post;
124ee25b4d0Sjosephpu   PetscMPIInt rank, size;
125461fcd32Sjosephpu   MPI_Comm    comm;
126461fcd32Sjosephpu 
127461fcd32Sjosephpu   PetscFunctionBeginUser;
128461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
129ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_rank(comm, &rank));
130ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_size(comm, &size));
131461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Basic point insertion check...\n"));
132461fcd32Sjosephpu   PetscCall(DMSwarmGetSize(sw, &Np_pre));
133461fcd32Sjosephpu   if (rank == 0) PetscCall(DMSwarmAddPoint(sw));
134461fcd32Sjosephpu   PetscCall(DMSwarmGetSize(sw, &Np_post));
135461fcd32Sjosephpu   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);
136461fcd32Sjosephpu   PetscCall(CheckMigrate(sw));
137461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Basic point insertion check passes.\n"));
138461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
139461fcd32Sjosephpu }
140461fcd32Sjosephpu 
141461fcd32Sjosephpu /*
142461fcd32Sjosephpu   Checks tie breaking works properly when a particle
143461fcd32Sjosephpu   is located at a shared boundary. The higher rank should
144461fcd32Sjosephpu   recieve the particle while the lower rank deletes it.
145461fcd32Sjosephpu 
146461fcd32Sjosephpu   TODO: Currently only works for 2 procs.
147461fcd32Sjosephpu */
148461fcd32Sjosephpu static PetscErrorCode CheckPointInsertion_Boundary(DM sw)
149461fcd32Sjosephpu {
150ee25b4d0Sjosephpu   PetscInt    Np_loc_pre, Np_loc_post, dim;
151ee25b4d0Sjosephpu   PetscMPIInt rank, size;
152461fcd32Sjosephpu   PetscReal   lbox_low[3], lbox_high[3], gbox_low[3], gbox_high[3];
153461fcd32Sjosephpu   MPI_Comm    comm;
154461fcd32Sjosephpu   DM          cdm;
155461fcd32Sjosephpu 
156461fcd32Sjosephpu   PetscFunctionBeginUser;
157461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
158ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_rank(comm, &rank));
159ee25b4d0Sjosephpu   PetscCallMPI(MPI_Comm_size(comm, &size));
160461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Rank boundary point insertion check...\n"));
161461fcd32Sjosephpu   PetscCall(DMSwarmGetCellDM(sw, &cdm));
162461fcd32Sjosephpu   PetscCall(DMGetDimension(cdm, &dim));
163461fcd32Sjosephpu   PetscCall(DMGetBoundingBox(cdm, gbox_low, gbox_high));
164461fcd32Sjosephpu   if (rank == 0) {
165461fcd32Sjosephpu     PetscReal *coords;
166ee25b4d0Sjosephpu     PetscInt   adjacentdim = 0, Np;
167461fcd32Sjosephpu 
168461fcd32Sjosephpu     PetscCall(DMGetLocalBoundingBox(cdm, lbox_low, lbox_high));
169461fcd32Sjosephpu     // find a face that belongs to the neighbor.
170461fcd32Sjosephpu     for (PetscInt d = 0; d < dim; ++d) {
171461fcd32Sjosephpu       if ((gbox_high[d] - lbox_high[d]) != 0.) adjacentdim = d;
172461fcd32Sjosephpu     }
173461fcd32Sjosephpu     PetscCall(DMSwarmAddPoint(sw));
174461fcd32Sjosephpu     PetscCall(DMSwarmGetLocalSize(sw, &Np));
175461fcd32Sjosephpu     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
176461fcd32Sjosephpu     for (PetscInt d = 0; d < dim; ++d) coords[(Np - 1) * dim + d] = 0.5 * (lbox_high[d] - lbox_low[d]);
177461fcd32Sjosephpu     coords[(Np - 1) * dim + adjacentdim] = lbox_high[adjacentdim];
178461fcd32Sjosephpu     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
179461fcd32Sjosephpu   }
180461fcd32Sjosephpu   PetscCall(DMSwarmGetLocalSize(sw, &Np_loc_pre));
181461fcd32Sjosephpu   PetscCall(CheckMigrate(sw));
182461fcd32Sjosephpu   PetscCall(DMSwarmGetLocalSize(sw, &Np_loc_post));
183ee25b4d0Sjosephpu   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);
184ee25b4d0Sjosephpu   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 recieved.", rank);
185461fcd32Sjosephpu   PetscCall(PetscPrintf(comm, "Rank boundary point insertion check passes.\n"));
186461fcd32Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
187461fcd32Sjosephpu }
188461fcd32Sjosephpu 
189461fcd32Sjosephpu int main(int argc, char **argv)
190461fcd32Sjosephpu {
191461fcd32Sjosephpu   DM       dm, sw;
192461fcd32Sjosephpu   MPI_Comm comm;
193461fcd32Sjosephpu   AppCtx   user;
194461fcd32Sjosephpu 
195461fcd32Sjosephpu   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
196461fcd32Sjosephpu   comm = PETSC_COMM_WORLD;
197461fcd32Sjosephpu   PetscCall(ProcessOptions(comm, &user));
198461fcd32Sjosephpu   PetscCall(CreateMesh(comm, &dm, &user));
199461fcd32Sjosephpu   PetscCall(CreateSwarm(dm, &sw, &user));
200461fcd32Sjosephpu   PetscCall(CheckMigrate(sw));
201461fcd32Sjosephpu   PetscCall(CheckPointInsertion(sw));
202461fcd32Sjosephpu   PetscCall(CheckPointInsertion_Boundary(sw));
203461fcd32Sjosephpu   PetscCall(DMDestroy(&sw));
204461fcd32Sjosephpu   PetscCall(DMDestroy(&dm));
205ee25b4d0Sjosephpu   PetscCall(PetscFinalize());
206461fcd32Sjosephpu   return 0;
207461fcd32Sjosephpu }
208461fcd32Sjosephpu 
209461fcd32Sjosephpu /*TEST
210ee25b4d0Sjosephpu 
211461fcd32Sjosephpu   # Swarm does not handle complex or quad
212461fcd32Sjosephpu   build:
213461fcd32Sjosephpu     requires: !complex double
214ee25b4d0Sjosephpu   # swarm_migrate_hash and swarm_migrate_scan test swarm migration against point location types
215ee25b4d0Sjosephpu   # with a distributed mesh where ranks overlap by 1. Points in the shared boundary should
216ee25b4d0Sjosephpu   # be sent to the process which has the highest rank that has that portion of the domain.
217461fcd32Sjosephpu   test:
218461fcd32Sjosephpu     suffix: swarm_migrate_hash
219461fcd32Sjosephpu     nsize: 2
220ee25b4d0Sjosephpu     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
221ee25b4d0Sjosephpu           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
222ee25b4d0Sjosephpu           -dm_plex_hash_location true
223461fcd32Sjosephpu     filter: grep -v marker | grep -v atomic | grep -v usage
224461fcd32Sjosephpu   test:
225*8ea0bac9SZach Atkins     suffix: swarm_migrate_hash_tensor_permutation
226*8ea0bac9SZach Atkins     nsize: 2
227*8ea0bac9SZach Atkins     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
228*8ea0bac9SZach Atkins           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
229*8ea0bac9SZach Atkins           -dm_plex_hash_location true -set_closure_permutation
230*8ea0bac9SZach Atkins     filter: grep -v marker | grep -v atomic | grep -v usage
231*8ea0bac9SZach Atkins   test:
232461fcd32Sjosephpu     suffix: swarm_migrate_scan
233461fcd32Sjosephpu     nsize: 2
234ee25b4d0Sjosephpu     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
235ee25b4d0Sjosephpu           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
236ee25b4d0Sjosephpu           -dm_plex_hash_location false
237461fcd32Sjosephpu     filter: grep -v marker | grep -v atomic | grep -v usage
238*8ea0bac9SZach Atkins   test:
239*8ea0bac9SZach Atkins     suffix: swarm_migrate_scan_tensor_permutation
240*8ea0bac9SZach Atkins     nsize: 2
241*8ea0bac9SZach Atkins     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_distribute_overlap 1 -dm_plex_box_faces 10,10,10\
242*8ea0bac9SZach Atkins           -dm_plex_box_lower 0.,0.,0. -dm_plex_box_upper 1.,1.,10. -dm_plex_box_bd none,none,none\
243*8ea0bac9SZach Atkins           -dm_plex_hash_location false -set_closure_permutation
244*8ea0bac9SZach Atkins     filter: grep -v marker | grep -v atomic | grep -v usage
245461fcd32Sjosephpu TEST*/
246