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