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