1068e8908SMatthew G. Knepley #include "petscdm.h" 288ff9c58SMatthew G. Knepley static char help[] = "Tests for particle initialization using the KS test\n\n"; 3a4a4ae2eSMatthew G. Knepley 4a4a4ae2eSMatthew G. Knepley #include <petscdmswarm.h> 5a4a4ae2eSMatthew G. Knepley #include <petscdmplex.h> 6a4a4ae2eSMatthew G. Knepley #include <petsc/private/petscfeimpl.h> 7a4a4ae2eSMatthew G. Knepley 8a4a4ae2eSMatthew G. Knepley /* 9a4a4ae2eSMatthew G. Knepley View the KS test using 10a4a4ae2eSMatthew G. Knepley 11a4a4ae2eSMatthew G. Knepley -ks_monitor draw -draw_size 500,500 -draw_pause 3 12a4a4ae2eSMatthew G. Knepley 13a4a4ae2eSMatthew G. Knepley Set total number to n0 / Mp = 3e14 / 1e12 = 300 using -dm_swarm_num_particles 300 14a4a4ae2eSMatthew G. Knepley 15a4a4ae2eSMatthew G. Knepley */ 16a4a4ae2eSMatthew G. Knepley 17a4a4ae2eSMatthew G. Knepley #define BOLTZMANN_K 1.380649e-23 /* J/K */ 18a4a4ae2eSMatthew G. Knepley 19a4a4ae2eSMatthew G. Knepley typedef struct { 20a4a4ae2eSMatthew G. Knepley PetscReal mass[2]; /* Electron, Sr+ Mass [kg] */ 21a4a4ae2eSMatthew G. Knepley PetscReal T[2]; /* Electron, Ion Temperature [K] */ 22a4a4ae2eSMatthew G. Knepley PetscReal v0[2]; /* Species mean velocity in 1D */ 23a4a4ae2eSMatthew G. Knepley } AppCtx; 24a4a4ae2eSMatthew G. Knepley 25d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 26d71ae5a4SJacob Faibussowitsch { 27a4a4ae2eSMatthew G. Knepley PetscFunctionBegin; 28a4a4ae2eSMatthew G. Knepley options->mass[0] = 9.10938356e-31; /* Electron Mass [kg] */ 29a4a4ae2eSMatthew G. Knepley options->mass[1] = 87.62 * 1.66054e-27; /* Sr+ Mass [kg] */ 30a4a4ae2eSMatthew G. Knepley options->T[0] = 1.; /* Electron Temperature [K] */ 31a4a4ae2eSMatthew G. Knepley options->T[1] = 25.; /* Sr+ Temperature [K] */ 32a4a4ae2eSMatthew G. Knepley options->v0[0] = PetscSqrtReal(BOLTZMANN_K * options->T[0] / options->mass[0]); /* electron mean velocity in 1D */ 33a4a4ae2eSMatthew G. Knepley options->v0[1] = PetscSqrtReal(BOLTZMANN_K * options->T[1] / options->mass[1]); /* ion mean velocity in 1D */ 34a4a4ae2eSMatthew G. Knepley 35d0609cedSBarry Smith PetscOptionsBegin(comm, "", "KS Test Options", "DMPLEX"); 36d0609cedSBarry Smith PetscOptionsEnd(); 373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38a4a4ae2eSMatthew G. Knepley } 39a4a4ae2eSMatthew G. Knepley 40d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 41d71ae5a4SJacob Faibussowitsch { 42a4a4ae2eSMatthew G. Knepley PetscFunctionBeginUser; 439566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 449566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 459566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 469566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48a4a4ae2eSMatthew G. Knepley } 49a4a4ae2eSMatthew G. Knepley 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw) 51d71ae5a4SJacob Faibussowitsch { 52a4a4ae2eSMatthew G. Knepley PetscInt dim; 53a4a4ae2eSMatthew G. Knepley 54a4a4ae2eSMatthew G. Knepley PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 569566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 579566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 589566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 599566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 609566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 619566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 629566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 639566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT)); 649566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 659566063dSJacob Faibussowitsch PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw)); 669566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeCoordinates(*sw)); 679566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, user->v0)); 689566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 709566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-swarm_view")); 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72a4a4ae2eSMatthew G. Knepley } 73a4a4ae2eSMatthew G. Knepley 74d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestDistribution(DM sw, PetscReal confidenceLevel, AppCtx *user) 75d71ae5a4SJacob Faibussowitsch { 76068e8908SMatthew G. Knepley DM dm; 77068e8908SMatthew G. Knepley Vec locx, locv, locw; 78*f8662bd6SBarry Smith PetscProbFn *cdf; 79068e8908SMatthew G. Knepley PetscReal alpha, gmin, gmax; 80a4a4ae2eSMatthew G. Knepley PetscInt dim; 81a4a4ae2eSMatthew G. Knepley MPI_Comm comm; 82a4a4ae2eSMatthew G. Knepley 83a4a4ae2eSMatthew G. Knepley PetscFunctionBeginUser; 849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 85068e8908SMatthew G. Knepley PetscCall(DMSwarmGetCellDM(sw, &dm)); 86068e8908SMatthew G. Knepley PetscCall(DMGetBoundingBox(dm, &gmin, &gmax)); 879566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 88a4a4ae2eSMatthew G. Knepley switch (dim) { 89d71ae5a4SJacob Faibussowitsch case 1: 90d71ae5a4SJacob Faibussowitsch cdf = PetscCDFMaxwellBoltzmann1D; 91d71ae5a4SJacob Faibussowitsch break; 92d71ae5a4SJacob Faibussowitsch case 2: 93d71ae5a4SJacob Faibussowitsch cdf = PetscCDFMaxwellBoltzmann2D; 94d71ae5a4SJacob Faibussowitsch break; 95d71ae5a4SJacob Faibussowitsch case 3: 96d71ae5a4SJacob Faibussowitsch cdf = PetscCDFMaxwellBoltzmann3D; 97d71ae5a4SJacob Faibussowitsch break; 98d71ae5a4SJacob Faibussowitsch default: 99d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported", dim); 100a4a4ae2eSMatthew G. Knepley } 101068e8908SMatthew G. Knepley PetscCall(DMSwarmCreateLocalVectorFromField(sw, DMSwarmPICField_coor, &locx)); 1029566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateLocalVectorFromField(sw, "velocity", &locv)); 103068e8908SMatthew G. Knepley PetscCall(DMSwarmCreateLocalVectorFromField(sw, "w_q", &locw)); 104068e8908SMatthew G. Knepley PetscCall(VecViewFromOptions(locx, NULL, "-x_view")); 105068e8908SMatthew G. Knepley PetscCall(VecViewFromOptions(locv, NULL, "-v_view")); 106068e8908SMatthew G. Knepley PetscCall(VecViewFromOptions(locw, NULL, "-w_view")); 107068e8908SMatthew G. Knepley // Need to rescale coordinates to compare to distribution 108068e8908SMatthew G. Knepley PetscCall(VecScale(locx, 1. / gmax)); 109068e8908SMatthew G. Knepley PetscCall(PetscProbComputeKSStatisticWeighted(locx, locw, PetscCDFConstant1D, &alpha)); 110068e8908SMatthew G. Knepley if (alpha < confidenceLevel) PetscCall(PetscPrintf(comm, "The KS test for X rejects the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 111068e8908SMatthew G. Knepley else PetscCall(PetscPrintf(comm, "The KS test for X accepts the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 112068e8908SMatthew G. Knepley PetscCall(PetscProbComputeKSStatisticWeighted(locv, locw, PetscCDFGaussian1D, &alpha)); 113068e8908SMatthew G. Knepley if (alpha < confidenceLevel) PetscCall(PetscPrintf(comm, "The KS test for V rejects the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 114068e8908SMatthew G. Knepley else PetscCall(PetscPrintf(comm, "The KS test for V accepts the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 115068e8908SMatthew G. Knepley PetscCall(PetscProbComputeKSStatisticMagnitude(locv, cdf, &alpha)); 116068e8908SMatthew G. Knepley if (alpha < confidenceLevel) PetscCall(PetscPrintf(comm, "The KS test for |V| rejects the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 117068e8908SMatthew G. Knepley else PetscCall(PetscPrintf(comm, "The KS test for |V| accepts the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 118068e8908SMatthew G. Knepley PetscCall(DMSwarmDestroyLocalVectorFromField(sw, DMSwarmPICField_coor, &locx)); 1199566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyLocalVectorFromField(sw, "velocity", &locv)); 120068e8908SMatthew G. Knepley PetscCall(DMSwarmDestroyLocalVectorFromField(sw, "w_q", &locw)); 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122a4a4ae2eSMatthew G. Knepley } 123a4a4ae2eSMatthew G. Knepley 124d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 125d71ae5a4SJacob Faibussowitsch { 126a4a4ae2eSMatthew G. Knepley DM dm, sw; 127a4a4ae2eSMatthew G. Knepley AppCtx user; 128a4a4ae2eSMatthew G. Knepley 129327415f7SBarry Smith PetscFunctionBeginUser; 1309566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1319566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1329566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 1339566063dSJacob Faibussowitsch PetscCall(CreateSwarm(dm, &user, &sw)); 134068e8908SMatthew G. Knepley PetscCall(TestDistribution(sw, 0.01, &user)); 1359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 1369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 138b122ec5aSJacob Faibussowitsch return 0; 139a4a4ae2eSMatthew G. Knepley } 140a4a4ae2eSMatthew G. Knepley 141a4a4ae2eSMatthew G. Knepley /*TEST 142a4a4ae2eSMatthew G. Knepley 143366d4293SMatthew G. Knepley build: 144366d4293SMatthew G. Knepley requires: !complex double 145366d4293SMatthew G. Knepley 146068e8908SMatthew G. Knepley testset: 147a4a4ae2eSMatthew G. Knepley requires: ks !complex 148068e8908SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_lower -4 -dm_plex_box_upper 4 -dm_plex_box_faces 40 -dm_swarm_num_particles 10000 149a4a4ae2eSMatthew G. Knepley 150366d4293SMatthew G. Knepley test: 151068e8908SMatthew G. Knepley suffix: 0_constant 152068e8908SMatthew G. Knepley args: -dm_swarm_coordinate_density constant 153068e8908SMatthew G. Knepley 154068e8908SMatthew G. Knepley test: 155068e8908SMatthew G. Knepley suffix: 0_gaussian 156068e8908SMatthew G. Knepley args: -dm_swarm_coordinate_density gaussian 157366d4293SMatthew G. Knepley 158a4a4ae2eSMatthew G. Knepley TEST*/ 159