xref: /petsc/src/dm/impls/swarm/tests/ex8.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
188ff9c58SMatthew G. Knepley static char help[] = "Tests for particle initialization using the KS test\n\n";
2a4a4ae2eSMatthew G. Knepley 
3a4a4ae2eSMatthew G. Knepley #include <petscdmswarm.h>
4a4a4ae2eSMatthew G. Knepley #include <petscdmplex.h>
5a4a4ae2eSMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
6a4a4ae2eSMatthew G. Knepley 
7a4a4ae2eSMatthew G. Knepley /*
8a4a4ae2eSMatthew G. Knepley   View the KS test using
9a4a4ae2eSMatthew G. Knepley 
10a4a4ae2eSMatthew G. Knepley     -ks_monitor draw -draw_size 500,500 -draw_pause 3
11a4a4ae2eSMatthew G. Knepley 
12a4a4ae2eSMatthew G. Knepley   Set total number to n0 / Mp = 3e14 / 1e12 =  300 using -dm_swarm_num_particles 300
13a4a4ae2eSMatthew G. Knepley 
14a4a4ae2eSMatthew G. Knepley */
15a4a4ae2eSMatthew G. Knepley 
16a4a4ae2eSMatthew G. Knepley #define BOLTZMANN_K 1.380649e-23 /* J/K */
17a4a4ae2eSMatthew G. Knepley 
18a4a4ae2eSMatthew G. Knepley typedef struct {
19a4a4ae2eSMatthew G. Knepley   PetscReal mass[2]; /* Electron, Sr+ Mass [kg] */
20a4a4ae2eSMatthew G. Knepley   PetscReal T[2];    /* Electron, Ion Temperature [K] */
21a4a4ae2eSMatthew G. Knepley   PetscReal v0[2];   /* Species mean velocity in 1D */
22a4a4ae2eSMatthew G. Knepley } AppCtx;
23a4a4ae2eSMatthew G. Knepley 
24a4a4ae2eSMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
25a4a4ae2eSMatthew G. Knepley {
26a4a4ae2eSMatthew G. Knepley   PetscFunctionBegin;
27a4a4ae2eSMatthew G. Knepley   options->mass[0] = 9.10938356e-31; /* Electron Mass [kg] */
28a4a4ae2eSMatthew G. Knepley   options->mass[1] = 87.62 * 1.66054e-27; /* Sr+ Mass [kg] */
29a4a4ae2eSMatthew G. Knepley   options->T[0]    = 1.; /* Electron Temperature [K] */
30a4a4ae2eSMatthew G. Knepley   options->T[1]    = 25.; /* Sr+ Temperature [K] */
31a4a4ae2eSMatthew G. Knepley   options->v0[0]   = PetscSqrtReal(BOLTZMANN_K * options->T[0] / options->mass[0]); /* electron mean velocity in 1D */
32a4a4ae2eSMatthew G. Knepley   options->v0[1]   = PetscSqrtReal(BOLTZMANN_K * options->T[1] / options->mass[1]); /* ion mean velocity in 1D */
33a4a4ae2eSMatthew G. Knepley 
34d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "KS Test Options", "DMPLEX");
35d0609cedSBarry Smith   PetscOptionsEnd();
36a4a4ae2eSMatthew G. Knepley   PetscFunctionReturn(0);
37a4a4ae2eSMatthew G. Knepley }
38a4a4ae2eSMatthew G. Knepley 
39a4a4ae2eSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
40a4a4ae2eSMatthew G. Knepley {
41a4a4ae2eSMatthew G. Knepley   PetscFunctionBeginUser;
429566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
439566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
449566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
459566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
46a4a4ae2eSMatthew G. Knepley   PetscFunctionReturn(0);
47a4a4ae2eSMatthew G. Knepley }
48a4a4ae2eSMatthew G. Knepley 
49a4a4ae2eSMatthew G. Knepley static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
50a4a4ae2eSMatthew G. Knepley {
51a4a4ae2eSMatthew G. Knepley   PetscInt       dim;
52a4a4ae2eSMatthew G. Knepley 
53a4a4ae2eSMatthew G. Knepley   PetscFunctionBeginUser;
549566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
559566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
569566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
579566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
589566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
599566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
609566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
619566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
629566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
639566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
649566063dSJacob Faibussowitsch   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
659566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeCoordinates(*sw));
669566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, user->v0));
679566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
699566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-swarm_view"));
70a4a4ae2eSMatthew G. Knepley   PetscFunctionReturn(0);
71a4a4ae2eSMatthew G. Knepley }
72a4a4ae2eSMatthew G. Knepley 
73a4a4ae2eSMatthew G. Knepley static PetscErrorCode TestDistribution(DM sw, PetscReal confidenceLevel, AppCtx *user)
74a4a4ae2eSMatthew G. Knepley {
75a4a4ae2eSMatthew G. Knepley   Vec            locv;
76a4a4ae2eSMatthew G. Knepley   PetscProbFunc  cdf;
77a4a4ae2eSMatthew G. Knepley   PetscReal      alpha;
78a4a4ae2eSMatthew G. Knepley   PetscInt       dim;
79a4a4ae2eSMatthew G. Knepley   MPI_Comm       comm;
80a4a4ae2eSMatthew G. Knepley 
81a4a4ae2eSMatthew G. Knepley   PetscFunctionBeginUser;
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) sw, &comm));
839566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
84a4a4ae2eSMatthew G. Knepley   switch (dim) {
85a4a4ae2eSMatthew G. Knepley     case 1: cdf = PetscCDFMaxwellBoltzmann1D;break;
86a4a4ae2eSMatthew G. Knepley     case 2: cdf = PetscCDFMaxwellBoltzmann2D;break;
87a4a4ae2eSMatthew G. Knepley     case 3: cdf = PetscCDFMaxwellBoltzmann3D;break;
8863a3b9bcSJacob Faibussowitsch     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported", dim);
89a4a4ae2eSMatthew G. Knepley   }
909566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateLocalVectorFromField(sw, "velocity", &locv));
919566063dSJacob Faibussowitsch   PetscCall(PetscProbComputeKSStatistic(locv, cdf, &alpha));
929566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyLocalVectorFromField(sw, "velocity", &locv));
939566063dSJacob Faibussowitsch   if (alpha < confidenceLevel) PetscCall(PetscPrintf(comm, "The KS test accepts the null hypothesis at level %.2g\n", (double) confidenceLevel));
949566063dSJacob Faibussowitsch   else                         PetscCall(PetscPrintf(comm, "The KS test rejects the null hypothesis at level %.2g (%.2g)\n", (double) confidenceLevel, (double) alpha));
95a4a4ae2eSMatthew G. Knepley   PetscFunctionReturn(0);
96a4a4ae2eSMatthew G. Knepley }
97a4a4ae2eSMatthew G. Knepley 
98a4a4ae2eSMatthew G. Knepley int main(int argc, char **argv)
99a4a4ae2eSMatthew G. Knepley {
100a4a4ae2eSMatthew G. Knepley   DM             dm, sw;
101a4a4ae2eSMatthew G. Knepley   AppCtx         user;
102a4a4ae2eSMatthew G. Knepley 
103*327415f7SBarry Smith   PetscFunctionBeginUser;
1049566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
1059566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1069566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
1079566063dSJacob Faibussowitsch   PetscCall(CreateSwarm(dm, &user, &sw));
1089566063dSJacob Faibussowitsch   PetscCall(TestDistribution(sw, 0.05, &user));
1099566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
1109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1119566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
112b122ec5aSJacob Faibussowitsch   return 0;
113a4a4ae2eSMatthew G. Knepley }
114a4a4ae2eSMatthew G. Knepley 
115a4a4ae2eSMatthew G. Knepley /*TEST
116a4a4ae2eSMatthew G. Knepley 
117a4a4ae2eSMatthew G. Knepley   test:
118a4a4ae2eSMatthew G. Knepley     suffix: 0
119a4a4ae2eSMatthew G. Knepley     requires: ks !complex
12088ff9c58SMatthew G. Knepley     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}}
121a4a4ae2eSMatthew G. Knepley 
122a4a4ae2eSMatthew G. Knepley TEST*/
123