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 24d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 25d71ae5a4SJacob Faibussowitsch { 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(); 36*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37a4a4ae2eSMatthew G. Knepley } 38a4a4ae2eSMatthew G. Knepley 39d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 40d71ae5a4SJacob Faibussowitsch { 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")); 46*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47a4a4ae2eSMatthew G. Knepley } 48a4a4ae2eSMatthew G. Knepley 49d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw) 50d71ae5a4SJacob Faibussowitsch { 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")); 70*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71a4a4ae2eSMatthew G. Knepley } 72a4a4ae2eSMatthew G. Knepley 73d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestDistribution(DM sw, PetscReal confidenceLevel, AppCtx *user) 74d71ae5a4SJacob Faibussowitsch { 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) { 85d71ae5a4SJacob Faibussowitsch case 1: 86d71ae5a4SJacob Faibussowitsch cdf = PetscCDFMaxwellBoltzmann1D; 87d71ae5a4SJacob Faibussowitsch break; 88d71ae5a4SJacob Faibussowitsch case 2: 89d71ae5a4SJacob Faibussowitsch cdf = PetscCDFMaxwellBoltzmann2D; 90d71ae5a4SJacob Faibussowitsch break; 91d71ae5a4SJacob Faibussowitsch case 3: 92d71ae5a4SJacob Faibussowitsch cdf = PetscCDFMaxwellBoltzmann3D; 93d71ae5a4SJacob Faibussowitsch break; 94d71ae5a4SJacob Faibussowitsch default: 95d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported", dim); 96a4a4ae2eSMatthew G. Knepley } 979566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateLocalVectorFromField(sw, "velocity", &locv)); 989566063dSJacob Faibussowitsch PetscCall(PetscProbComputeKSStatistic(locv, cdf, &alpha)); 999566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyLocalVectorFromField(sw, "velocity", &locv)); 1009566063dSJacob Faibussowitsch if (alpha < confidenceLevel) PetscCall(PetscPrintf(comm, "The KS test accepts the null hypothesis at level %.2g\n", (double)confidenceLevel)); 1019566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "The KS test rejects the null hypothesis at level %.2g (%.2g)\n", (double)confidenceLevel, (double)alpha)); 102*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103a4a4ae2eSMatthew G. Knepley } 104a4a4ae2eSMatthew G. Knepley 105d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 106d71ae5a4SJacob Faibussowitsch { 107a4a4ae2eSMatthew G. Knepley DM dm, sw; 108a4a4ae2eSMatthew G. Knepley AppCtx user; 109a4a4ae2eSMatthew G. Knepley 110327415f7SBarry Smith PetscFunctionBeginUser; 1119566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1129566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1139566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 1149566063dSJacob Faibussowitsch PetscCall(CreateSwarm(dm, &user, &sw)); 1159566063dSJacob Faibussowitsch PetscCall(TestDistribution(sw, 0.05, &user)); 1169566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 1179566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 119b122ec5aSJacob Faibussowitsch return 0; 120a4a4ae2eSMatthew G. Knepley } 121a4a4ae2eSMatthew G. Knepley 122a4a4ae2eSMatthew G. Knepley /*TEST 123a4a4ae2eSMatthew G. Knepley 124a4a4ae2eSMatthew G. Knepley test: 125a4a4ae2eSMatthew G. Knepley suffix: 0 126a4a4ae2eSMatthew G. Knepley requires: ks !complex 12788ff9c58SMatthew 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}} 128a4a4ae2eSMatthew G. Knepley 129a4a4ae2eSMatthew G. Knepley TEST*/ 130