1*8c87cf4dSdanfinn static char help[] = "Test of Colorized Scatter Plot.\n"; 2*8c87cf4dSdanfinn 3*8c87cf4dSdanfinn #include "petscdraw.h" 4*8c87cf4dSdanfinn #include "petscvec.h" 5*8c87cf4dSdanfinn #include "petscis.h" 6*8c87cf4dSdanfinn 7*8c87cf4dSdanfinn typedef struct { 8*8c87cf4dSdanfinn PetscInt Np; /* total number of particles */ 9*8c87cf4dSdanfinn PetscInt dim; 10*8c87cf4dSdanfinn PetscInt dim_inp; 11*8c87cf4dSdanfinn } AppCtx; 12*8c87cf4dSdanfinn 13*8c87cf4dSdanfinn static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 14*8c87cf4dSdanfinn { 15*8c87cf4dSdanfinn PetscErrorCode ierr; 16*8c87cf4dSdanfinn 17*8c87cf4dSdanfinn PetscFunctionBeginUser; 18*8c87cf4dSdanfinn options->dim = 2; 19*8c87cf4dSdanfinn options->dim_inp = 2; 20*8c87cf4dSdanfinn options->Np = 100; 21*8c87cf4dSdanfinn 22*8c87cf4dSdanfinn ierr = PetscOptionsBegin(comm, "", "Test of colorized scatter plot", "");CHKERRQ(ierr); 23*8c87cf4dSdanfinn ierr = PetscOptionsInt("-Np", "Number of particles", "ex35.c", options->Np, &options->Np, PETSC_NULL);CHKERRQ(ierr); 24*8c87cf4dSdanfinn ierr = PetscOptionsInt("-dim", "Number of dimensions", "ex35.c", options->dim_inp, &options->dim_inp, PETSC_NULL);CHKERRQ(ierr); 25*8c87cf4dSdanfinn ierr = PetscOptionsEnd();CHKERRQ(ierr); 26*8c87cf4dSdanfinn PetscFunctionReturn(0); 27*8c87cf4dSdanfinn } 28*8c87cf4dSdanfinn 29*8c87cf4dSdanfinn /* 30*8c87cf4dSdanfinn ref: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/ 31*8c87cf4dSdanfinn */ 32*8c87cf4dSdanfinn PetscReal erfinv(PetscReal x) 33*8c87cf4dSdanfinn { 34*8c87cf4dSdanfinn PetscReal *ck, r = 0.; 35*8c87cf4dSdanfinn PetscInt k, m, maxIter=100; 36*8c87cf4dSdanfinn PetscErrorCode ierr; 37*8c87cf4dSdanfinn 38*8c87cf4dSdanfinn ierr = PetscCalloc1(maxIter,&ck);CHKERRQ(ierr); 39*8c87cf4dSdanfinn ck[0] = 1; 40*8c87cf4dSdanfinn r = ck[0]*((PetscSqrtReal(PETSC_PI)/2.)*x); 41*8c87cf4dSdanfinn for (k = 1; k < maxIter; ++k){ 42*8c87cf4dSdanfinn for (m = 0; m <= k-1; ++m){ 43*8c87cf4dSdanfinn PetscReal denom = (m+1.)*(2.*m+1.); 44*8c87cf4dSdanfinn ck[k] += (ck[m]*ck[k-1-m])/denom; 45*8c87cf4dSdanfinn } 46*8c87cf4dSdanfinn PetscReal temp = 2.*k+1.; 47*8c87cf4dSdanfinn r += (ck[k]/temp)*PetscPowReal((PetscSqrtReal(PETSC_PI)/2.)*x,2.*k+1.); 48*8c87cf4dSdanfinn } 49*8c87cf4dSdanfinn ierr = PetscFree(ck); 50*8c87cf4dSdanfinn return r; 51*8c87cf4dSdanfinn } 52*8c87cf4dSdanfinn 53*8c87cf4dSdanfinn int main(int argc, char **argv) 54*8c87cf4dSdanfinn { 55*8c87cf4dSdanfinn PetscInt p, dim, Np; 56*8c87cf4dSdanfinn PetscScalar *randVecNums; 57*8c87cf4dSdanfinn PetscReal speed, value, *x, *v; 58*8c87cf4dSdanfinn PetscRandom rngx, rng1, rng2; 59*8c87cf4dSdanfinn Vec randVec, subvecvx, subvecvy; 60*8c87cf4dSdanfinn IS isvx, isvy; 61*8c87cf4dSdanfinn AppCtx user; 62*8c87cf4dSdanfinn PetscDrawAxis axis; 63*8c87cf4dSdanfinn PetscDraw positionDraw; 64*8c87cf4dSdanfinn PetscDrawSP positionDrawSP; 65*8c87cf4dSdanfinn MPI_Comm comm; 66*8c87cf4dSdanfinn PetscErrorCode ierr; 67*8c87cf4dSdanfinn 68*8c87cf4dSdanfinn ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 69*8c87cf4dSdanfinn comm = PETSC_COMM_WORLD; 70*8c87cf4dSdanfinn ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 71*8c87cf4dSdanfinn 72*8c87cf4dSdanfinn Np = user.Np; 73*8c87cf4dSdanfinn dim = user.dim; 74*8c87cf4dSdanfinn 75*8c87cf4dSdanfinn ierr = PetscMalloc2(Np*dim, &x, Np*dim, &v);CHKERRQ(ierr); 76*8c87cf4dSdanfinn 77*8c87cf4dSdanfinn ierr = PetscRandomCreate(comm, &rngx);CHKERRQ(ierr); 78*8c87cf4dSdanfinn ierr = PetscRandomSetInterval(rngx, 0., 1.);CHKERRQ(ierr); 79*8c87cf4dSdanfinn ierr = PetscRandomSetFromOptions(rngx);CHKERRQ(ierr); 80*8c87cf4dSdanfinn ierr = PetscRandomSetSeed(rngx, 1034);CHKERRQ(ierr); 81*8c87cf4dSdanfinn ierr = PetscRandomSeed(rngx);CHKERRQ(ierr); 82*8c87cf4dSdanfinn 83*8c87cf4dSdanfinn ierr = PetscRandomCreate(comm, &rng1);CHKERRQ(ierr); 84*8c87cf4dSdanfinn ierr = PetscRandomSetInterval(rng1, 0., 1.);CHKERRQ(ierr); 85*8c87cf4dSdanfinn ierr = PetscRandomSetFromOptions(rng1);CHKERRQ(ierr); 86*8c87cf4dSdanfinn ierr = PetscRandomSetSeed(rng1, 3084);CHKERRQ(ierr); 87*8c87cf4dSdanfinn ierr = PetscRandomSeed(rng1);CHKERRQ(ierr); 88*8c87cf4dSdanfinn 89*8c87cf4dSdanfinn ierr = PetscRandomCreate(comm, &rng2);CHKERRQ(ierr); 90*8c87cf4dSdanfinn ierr = PetscRandomSetInterval(rng2, 0., 1.);CHKERRQ(ierr); 91*8c87cf4dSdanfinn ierr = PetscRandomSetFromOptions(rng2);CHKERRQ(ierr); 92*8c87cf4dSdanfinn ierr = PetscRandomSetSeed(rng2, 2397);CHKERRQ(ierr); 93*8c87cf4dSdanfinn ierr = PetscRandomSeed(rng2);CHKERRQ(ierr); 94*8c87cf4dSdanfinn 95*8c87cf4dSdanfinn /* Set particle positions and velocities */ 96*8c87cf4dSdanfinn if (user.dim_inp == 1) { 97*8c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 98*8c87cf4dSdanfinn PetscReal temp; 99*8c87cf4dSdanfinn ierr = PetscRandomGetValueReal(rngx, &value);CHKERRQ(ierr); 100*8c87cf4dSdanfinn x[p*dim] = value; 101*8c87cf4dSdanfinn x[p*dim+1] = 0.; 102*8c87cf4dSdanfinn temp = erfinv(2*value-1); 103*8c87cf4dSdanfinn v[p*dim] = temp; 104*8c87cf4dSdanfinn v[p*dim+1] = 0.; 105*8c87cf4dSdanfinn } 106*8c87cf4dSdanfinn } else if (user.dim_inp == 2) { 107*8c87cf4dSdanfinn /* 108*8c87cf4dSdanfinn Use Box-Muller to sample a distribution of velocities for the maxwellian. 109*8c87cf4dSdanfinn https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform 110*8c87cf4dSdanfinn */ 111*8c87cf4dSdanfinn ierr = VecCreate(comm,&randVec); 112*8c87cf4dSdanfinn ierr = VecSetSizes(randVec,PETSC_DECIDE, Np*dim); 113*8c87cf4dSdanfinn ierr = VecSetFromOptions(randVec); 114*8c87cf4dSdanfinn 115*8c87cf4dSdanfinn ierr = ISCreateStride(comm, Np*dim/2, 0, 2, &isvx);CHKERRQ(ierr); 116*8c87cf4dSdanfinn ierr = ISCreateStride(comm, Np*dim/2, 1, 2, &isvy);CHKERRQ(ierr); 117*8c87cf4dSdanfinn ierr = VecGetSubVector(randVec, isvx, &subvecvx);CHKERRQ(ierr); 118*8c87cf4dSdanfinn ierr = VecGetSubVector(randVec, isvy, &subvecvy);CHKERRQ(ierr); 119*8c87cf4dSdanfinn ierr = VecSetRandom(subvecvx, rng1);CHKERRQ(ierr); 120*8c87cf4dSdanfinn ierr = VecSetRandom(subvecvy, rng2);CHKERRQ(ierr); 121*8c87cf4dSdanfinn ierr = VecRestoreSubVector(randVec, isvx, &subvecvx);CHKERRQ(ierr); 122*8c87cf4dSdanfinn ierr = VecRestoreSubVector(randVec, isvy, &subvecvy);CHKERRQ(ierr); 123*8c87cf4dSdanfinn ierr = VecGetArray(randVec, &randVecNums);CHKERRQ(ierr); 124*8c87cf4dSdanfinn 125*8c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 126*8c87cf4dSdanfinn PetscReal u1, u2, mag, zx, zy; 127*8c87cf4dSdanfinn 128*8c87cf4dSdanfinn u1 = PetscRealPart(randVecNums[p*dim]); 129*8c87cf4dSdanfinn u2 = PetscRealPart(randVecNums[p*dim+1]); 130*8c87cf4dSdanfinn 131*8c87cf4dSdanfinn x[p*dim] = u1; 132*8c87cf4dSdanfinn x[p*dim+1] = u2; 133*8c87cf4dSdanfinn 134*8c87cf4dSdanfinn mag = PetscSqrtReal(-2.0 * PetscLogReal(u1)); 135*8c87cf4dSdanfinn 136*8c87cf4dSdanfinn zx = mag * PetscCosReal(2*PETSC_PI*u2) + 0; 137*8c87cf4dSdanfinn zy = mag * PetscSinReal(2*PETSC_PI*u2) + 0; 138*8c87cf4dSdanfinn 139*8c87cf4dSdanfinn v[p*dim] = zx; 140*8c87cf4dSdanfinn v[p*dim+1] = zy; 141*8c87cf4dSdanfinn } 142*8c87cf4dSdanfinn ierr = ISDestroy(&isvx);CHKERRQ(ierr); 143*8c87cf4dSdanfinn ierr = ISDestroy(&isvy);CHKERRQ(ierr); 144*8c87cf4dSdanfinn ierr = VecDestroy(&subvecvx);CHKERRQ(ierr); 145*8c87cf4dSdanfinn ierr = VecDestroy(&subvecvy);CHKERRQ(ierr); 146*8c87cf4dSdanfinn ierr = VecDestroy(&randVec);CHKERRQ(ierr); 147*8c87cf4dSdanfinn } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %D", dim); 148*8c87cf4dSdanfinn 149*8c87cf4dSdanfinn ierr = PetscDrawCreate(comm, NULL, "monitor_particle_positions", 0,0,400,300, &positionDraw);CHKERRQ(ierr); 150*8c87cf4dSdanfinn ierr = PetscDrawSetFromOptions(positionDraw);CHKERRQ(ierr); 151*8c87cf4dSdanfinn ierr = PetscDrawSPCreate(positionDraw, 10, &positionDrawSP);CHKERRQ(ierr); 152*8c87cf4dSdanfinn ierr = PetscDrawSPSetDimension(positionDrawSP,1);CHKERRQ(ierr); 153*8c87cf4dSdanfinn ierr = PetscDrawSPGetAxis(positionDrawSP, &axis);CHKERRQ(ierr); 154*8c87cf4dSdanfinn ierr = PetscDrawSPReset(positionDrawSP);CHKERRQ(ierr); 155*8c87cf4dSdanfinn ierr = PetscDrawAxisSetLabels(axis, "Particles", "x", "y");CHKERRQ(ierr); 156*8c87cf4dSdanfinn ierr = PetscDrawSetSave(positionDraw, "ex35_pos.ppm");CHKERRQ(ierr); 157*8c87cf4dSdanfinn ierr = PetscDrawSPReset(positionDrawSP);CHKERRQ(ierr); 158*8c87cf4dSdanfinn ierr = PetscDrawSPSetLimits(positionDrawSP, 0, 1, 0, 1);CHKERRQ(ierr); 159*8c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 160*8c87cf4dSdanfinn speed = PetscSqrtReal(PetscSqr(v[p*dim]) + PetscSqr(v[p*dim+1])); 161*8c87cf4dSdanfinn ierr = PetscDrawSPAddPointColorized(positionDrawSP, &x[p*dim], &x[p*dim+1], &speed);CHKERRQ(ierr); 162*8c87cf4dSdanfinn } 163*8c87cf4dSdanfinn ierr = PetscDrawSPDraw(positionDrawSP, PETSC_TRUE);CHKERRQ(ierr); 164*8c87cf4dSdanfinn ierr = PetscDrawSave(positionDraw);CHKERRQ(ierr); 165*8c87cf4dSdanfinn 166*8c87cf4dSdanfinn ierr = PetscFree2(x, v);CHKERRQ(ierr); 167*8c87cf4dSdanfinn ierr = PetscRandomDestroy(&rngx);CHKERRQ(ierr); 168*8c87cf4dSdanfinn ierr = PetscRandomDestroy(&rng1);CHKERRQ(ierr); 169*8c87cf4dSdanfinn ierr = PetscRandomDestroy(&rng2);CHKERRQ(ierr); 170*8c87cf4dSdanfinn 171*8c87cf4dSdanfinn ierr = PetscDrawSPDestroy(&positionDrawSP);CHKERRQ(ierr); 172*8c87cf4dSdanfinn ierr = PetscDrawDestroy(&positionDraw);CHKERRQ(ierr); 173*8c87cf4dSdanfinn ierr = PetscFinalize(); 174*8c87cf4dSdanfinn return ierr; 175*8c87cf4dSdanfinn } 176*8c87cf4dSdanfinn 177*8c87cf4dSdanfinn /*TEST 178*8c87cf4dSdanfinn test: 179*8c87cf4dSdanfinn suffix: 1D 180*8c87cf4dSdanfinn args: -Np 50\ 181*8c87cf4dSdanfinn -dim 1 182*8c87cf4dSdanfinn test: 183*8c87cf4dSdanfinn suffix: 2D 184*8c87cf4dSdanfinn args: -Np 50\ 185*8c87cf4dSdanfinn -dim 2 186*8c87cf4dSdanfinn TEST*/ 187