xref: /petsc/src/dm/impls/swarm/tests/ex3.c (revision d2522c19e8fa9bca20aaca277941d9a63e71db6a)
1 static char help[] = "Example usage of extracting single cells with their associated fields from a swarm and putting it in a new swarm object\n";
2 
3 #include <petscdmplex.h>
4 #include <petscdmswarm.h>
5 #include <petscts.h>
6 
7 typedef struct {
8   PetscInt particlesPerCell; /* The number of partices per cell */
9 } AppCtx;
10 
11 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
12   PetscFunctionBeginUser;
13   options->particlesPerCell = 1;
14 
15   PetscOptionsBegin(comm, "", "CellSwarm Options", "DMSWARM");
16   PetscCall(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex3.c", options->particlesPerCell, &options->particlesPerCell, NULL));
17   PetscOptionsEnd();
18   PetscFunctionReturn(0);
19 }
20 
21 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) {
22   PetscFunctionBeginUser;
23   PetscCall(DMCreate(comm, dm));
24   PetscCall(DMSetType(*dm, DMPLEX));
25   PetscCall(DMSetFromOptions(*dm));
26   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
27   PetscFunctionReturn(0);
28 }
29 
30 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) {
31   PetscInt *cellid;
32   PetscInt  dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
33 
34   PetscFunctionBeginUser;
35   PetscCall(DMGetDimension(dm, &dim));
36   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
37   PetscCall(DMSetType(*sw, DMSWARM));
38   PetscCall(DMSetDimension(*sw, dim));
39   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
40   PetscCall(DMSwarmSetCellDM(*sw, dm));
41   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL));
42   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
43   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
44   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
45   PetscCall(DMSetFromOptions(*sw));
46   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
47   for (c = cStart; c < cEnd; ++c) {
48     for (p = 0; p < Np; ++p) {
49       const PetscInt n = c * Np + p;
50 
51       cellid[n] = c;
52     }
53   }
54   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
55   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
56   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
57   PetscFunctionReturn(0);
58 }
59 
60 int main(int argc, char **argv) {
61   DM       dm, sw, cellsw; /* Mesh and particle managers */
62   MPI_Comm comm;
63   AppCtx   user;
64 
65   PetscFunctionBeginUser;
66   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
67   comm = PETSC_COMM_WORLD;
68   PetscCall(ProcessOptions(comm, &user));
69   PetscCall(CreateMesh(comm, &dm, &user));
70   PetscCall(CreateParticles(dm, &sw, &user));
71   PetscCall(DMSetApplicationContext(sw, &user));
72   PetscCall(DMCreate(comm, &cellsw));
73   PetscCall(PetscObjectSetName((PetscObject)cellsw, "SubParticles"));
74   PetscCall(DMSwarmGetCellSwarm(sw, 1, cellsw));
75   PetscCall(DMViewFromOptions(cellsw, NULL, "-subswarm_view"));
76   PetscCall(DMSwarmRestoreCellSwarm(sw, 1, cellsw));
77   PetscCall(DMDestroy(&sw));
78   PetscCall(DMDestroy(&dm));
79   PetscCall(DMDestroy(&cellsw));
80   PetscCall(PetscFinalize());
81   return 0;
82 }
83 
84 /*TEST
85    build:
86      requires: triangle !single !complex
87    test:
88      suffix: 1
89      args: -particles_per_cell 2 -dm_plex_box_faces 2,1 -dm_view -sw_view -subswarm_view
90 TEST*/
91