xref: /petsc/src/dm/impls/swarm/tests/ex8.c (revision ef46b1a67e276116c83b5d4ce8efc2932ea4fc0a)
1 static char help[] = "Tests for particle initialization using the KS test\n\n";
2 
3 #include <petscdmswarm.h>
4 #include <petscdmplex.h>
5 #include <petsc/private/petscfeimpl.h>
6 
7 /*
8   View the KS test using
9 
10     -ks_monitor draw -draw_size 500,500 -draw_pause 3
11 
12   Set total number to n0 / Mp = 3e14 / 1e12 =  300 using -dm_swarm_num_particles 300
13 
14 */
15 
16 #define BOLTZMANN_K 1.380649e-23 /* J/K */
17 
18 typedef struct {
19   PetscReal mass[2]; /* Electron, Sr+ Mass [kg] */
20   PetscReal T[2];    /* Electron, Ion Temperature [K] */
21   PetscReal v0[2];   /* Species mean velocity in 1D */
22 } AppCtx;
23 
24 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
25 {
26   PetscFunctionBegin;
27   options->mass[0] = 9.10938356e-31; /* Electron Mass [kg] */
28   options->mass[1] = 87.62 * 1.66054e-27; /* Sr+ Mass [kg] */
29   options->T[0]    = 1.; /* Electron Temperature [K] */
30   options->T[1]    = 25.; /* Sr+ Temperature [K] */
31   options->v0[0]   = PetscSqrtReal(BOLTZMANN_K * options->T[0] / options->mass[0]); /* electron mean velocity in 1D */
32   options->v0[1]   = PetscSqrtReal(BOLTZMANN_K * options->T[1] / options->mass[1]); /* ion mean velocity in 1D */
33 
34   PetscOptionsBegin(comm, "", "KS Test Options", "DMPLEX");
35   PetscOptionsEnd();
36   PetscFunctionReturn(0);
37 }
38 
39 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
40 {
41   PetscFunctionBeginUser;
42   PetscCall(DMCreate(comm, dm));
43   PetscCall(DMSetType(*dm, DMPLEX));
44   PetscCall(DMSetFromOptions(*dm));
45   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
46   PetscFunctionReturn(0);
47 }
48 
49 static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
50 {
51   PetscInt       dim;
52 
53   PetscFunctionBeginUser;
54   PetscCall(DMGetDimension(dm, &dim));
55   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
56   PetscCall(DMSetType(*sw, DMSWARM));
57   PetscCall(DMSetDimension(*sw, dim));
58   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
59   PetscCall(DMSwarmSetCellDM(*sw, dm));
60   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
61   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
62   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
63   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
64   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
65   PetscCall(DMSwarmInitializeCoordinates(*sw));
66   PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, user->v0));
67   PetscCall(DMSetFromOptions(*sw));
68   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
69   PetscCall(DMViewFromOptions(*sw, NULL, "-swarm_view"));
70   PetscFunctionReturn(0);
71 }
72 
73 static PetscErrorCode TestDistribution(DM sw, PetscReal confidenceLevel, AppCtx *user)
74 {
75   Vec            locv;
76   PetscProbFunc  cdf;
77   PetscReal      alpha;
78   PetscInt       dim;
79   MPI_Comm       comm;
80 
81   PetscFunctionBeginUser;
82   PetscCall(PetscObjectGetComm((PetscObject) sw, &comm));
83   PetscCall(DMGetDimension(sw, &dim));
84   switch (dim) {
85     case 1: cdf = PetscCDFMaxwellBoltzmann1D;break;
86     case 2: cdf = PetscCDFMaxwellBoltzmann2D;break;
87     case 3: cdf = PetscCDFMaxwellBoltzmann3D;break;
88     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
89   }
90   PetscCall(DMSwarmCreateLocalVectorFromField(sw, "velocity", &locv));
91   PetscCall(PetscProbComputeKSStatistic(locv, cdf, &alpha));
92   PetscCall(DMSwarmDestroyLocalVectorFromField(sw, "velocity", &locv));
93   if (alpha < confidenceLevel) PetscCall(PetscPrintf(comm, "The KS test accepts the null hypothesis at level %.2g\n", (double) confidenceLevel));
94   else                         PetscCall(PetscPrintf(comm, "The KS test rejects the null hypothesis at level %.2g (%.2g)\n", (double) confidenceLevel, (double) alpha));
95   PetscFunctionReturn(0);
96 }
97 
98 int main(int argc, char **argv)
99 {
100   DM             dm, sw;
101   AppCtx         user;
102 
103   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
104   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
105   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
106   PetscCall(CreateSwarm(dm, &user, &sw));
107   PetscCall(TestDistribution(sw, 0.05, &user));
108   PetscCall(DMDestroy(&sw));
109   PetscCall(DMDestroy(&dm));
110   PetscCall(PetscFinalize());
111   return 0;
112 }
113 
114 /*TEST
115 
116   test:
117     suffix: 0
118     requires: ks !complex
119     args: -dm_plex_dim 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1 -dm_swarm_num_particles 375 -dm_swarm_coordinate_density {{constant gaussian}}
120 
121 TEST*/
122