xref: /petsc/src/dm/impls/swarm/tests/ex12.c (revision 98d129c30f3ee9fdddc40fdbc5a989b7be64f888)
1 static char help[] = "Test periodic DMDA for DMSwarm point location.\n";
2 
3 #include <petscdmda.h>
4 #include <petscdmswarm.h>
5 
6 typedef struct {
7   PetscInt dim; // Mesh dimension
8   PetscInt Np;  // Number of particles along each dimension
9 } UserContext;
10 
11 static PetscErrorCode ProcessOptions(UserContext *options)
12 {
13   PetscFunctionBeginUser;
14   options->dim = 3;
15   options->Np  = -1;
16 
17   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &options->dim, NULL));
18   PetscCall(PetscOptionsGetInt(NULL, NULL, "-np", &options->Np, NULL));
19   PetscFunctionReturn(PETSC_SUCCESS);
20 }
21 
22 static PetscErrorCode CreateMesh(DM *da, UserContext *user)
23 {
24   PetscReal gmin[3] = {0, 0., 0.}, gmax[3] = {0, 0., 0.};
25 
26   PetscFunctionBeginUser;
27   PetscCall(DMDACreate(PETSC_COMM_WORLD, da));
28   PetscCall(DMSetDimension(*da, user->dim));
29   PetscCall(DMDASetSizes(*da, 7, 7, 7));
30   PetscCall(DMDASetBoundaryType(*da, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED));
31   PetscCall(DMDASetDof(*da, 2));
32   PetscCall(DMDASetStencilType(*da, DMDA_STENCIL_BOX));
33   PetscCall(DMDASetStencilWidth(*da, 1));
34   PetscCall(DMDASetElementType(*da, DMDA_ELEMENT_Q1));
35   PetscCall(DMSetFromOptions(*da));
36   PetscCall(DMSetUp(*da));
37   PetscCall(DMDASetUniformCoordinates(*da, 0., 1., 0., 1., 0., 1.));
38   PetscCall(DMSetApplicationContext(*da, user));
39   PetscCall(DMView(*da, PETSC_VIEWER_STDOUT_WORLD));
40   PetscCall(DMGetBoundingBox(*da, gmin, gmax));
41   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "min: (%g, %g, %g) max: (%g, %g, %g)\n", gmin[0], gmin[1], gmin[2], gmax[0], gmax[1], gmax[2]));
42   PetscFunctionReturn(PETSC_SUCCESS);
43 }
44 
45 static PetscErrorCode CreateSwarm(DM mesh, DM *swarm, UserContext *user)
46 {
47   MPI_Comm    comm;
48   PetscMPIInt size;
49   PetscInt    dim;
50 
51   PetscFunctionBeginUser;
52   PetscCall(PetscObjectGetComm((PetscObject)mesh, &comm));
53   PetscCallMPI(MPI_Comm_size(comm, &size));
54   PetscCall(DMCreate(comm, swarm));
55   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*swarm), "pic_"));
56   PetscCall(DMSetType(*swarm, DMSWARM));
57   PetscCall(PetscObjectSetName((PetscObject)*swarm, "ions"));
58   PetscCall(DMGetDimension(mesh, &dim));
59   PetscCall(DMSetDimension(*swarm, dim));
60   PetscCall(DMSwarmSetType(*swarm, DMSWARM_PIC));
61   PetscCall(DMSwarmSetCellDM(*swarm, mesh));
62   PetscCall(DMSwarmInitializeFieldRegister(*swarm));
63   PetscCall(DMSwarmFinalizeFieldRegister(*swarm));
64   PetscCall(DMSwarmSetLocalSizes(*swarm, user->Np / size, 0));
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 static PetscErrorCode InitializeParticles(DM sw, UserContext *user)
69 {
70   DM        da;
71   PetscReal gmin[3], gmax[3];
72   PetscInt  ndir[3];
73 
74   PetscFunctionBeginUser;
75   PetscCall(DMSwarmGetCellDM(sw, &da));
76   PetscCall(DMGetBoundingBox(da, gmin, gmax));
77   ndir[0] = user->Np;
78   ndir[1] = user->Np;
79   ndir[2] = user->Np;
80   PetscCall(DMSwarmSetPointsUniformCoordinates(sw, gmin, gmax, ndir, INSERT_VALUES));
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83 
84 int main(int argc, char **args)
85 {
86   DM          dm, sw;
87   UserContext user;
88 
89   PetscFunctionBeginUser;
90   PetscCall(PetscInitialize(&argc, &args, NULL, help));
91   PetscCall(ProcessOptions(&user));
92   PetscCall(CreateMesh(&dm, &user));
93   PetscCall(CreateSwarm(dm, &sw, &user));
94 
95   PetscCall(InitializeParticles(sw, &user));
96   PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
97   PetscCall(DMViewFromOptions(sw, NULL, "-sw_view"));
98 
99   PetscCall(DMDestroy(&dm));
100   PetscCall(DMDestroy(&sw));
101   PetscCall(PetscFinalize());
102   return 0;
103 }
104 
105 // ./dmswarm-coords -nx 8 -ny 8 -np 4 -pic_sw_view -dim 2 -periodic 0
106 // ./dmswarm-coords -nx 8 -ny 8 -np 4 -pic_sw_view -dim 2 -periodic 1
107 // ./dmswarm-coords -nx 8 -ny 8 -nz 8 -np 4 -pic_sw_view -dim 3 -periodic 0
108 // ./dmswarm-coords -nx 8 -ny 8 -nz 8 -np 4 -pic_sw_view -dim 3 -periodic 1
109 
110 /*TEST
111 
112   build:
113     requires: double !complex
114 
115   testset:
116     suffix: 2d
117     args: -dim 2 -da_grid_x 8 -da_grid_y 8 -np 4 -da_bd_all {{none periodic}} -pic_sw_view
118 
119     test:
120       suffix: p1
121 
122     test:
123       suffix: p2
124       nsize: 2
125 
126     test:
127       suffix: p4
128       nsize: 4
129 
130   testset:
131     suffix: 3d
132     args: -dim 3 -da_grid_x 8 -da_grid_y 8 -da_grid_z 8 -np 4 -da_bd_all {{none periodic}} -pic_sw_view
133 
134     test:
135       suffix: p1
136 
137     test:
138       suffix: p2
139       nsize: 2
140 
141     test:
142       suffix: p4
143       nsize: 4
144 
145 TEST*/
146