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