18c87cf4dSdanfinn static char help[] = "Test of Colorized Scatter Plot.\n"; 28c87cf4dSdanfinn 38c87cf4dSdanfinn #include "petscdraw.h" 48c87cf4dSdanfinn #include "petscvec.h" 58c87cf4dSdanfinn #include "petscis.h" 68c87cf4dSdanfinn 78c87cf4dSdanfinn typedef struct { 88c87cf4dSdanfinn PetscInt Np; /* total number of particles */ 98c87cf4dSdanfinn PetscInt dim; 108c87cf4dSdanfinn PetscInt dim_inp; 118c87cf4dSdanfinn } AppCtx; 128c87cf4dSdanfinn 138c87cf4dSdanfinn static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 148c87cf4dSdanfinn { 158c87cf4dSdanfinn PetscErrorCode ierr; 168c87cf4dSdanfinn 178c87cf4dSdanfinn PetscFunctionBeginUser; 188c87cf4dSdanfinn options->dim = 2; 198c87cf4dSdanfinn options->dim_inp = 2; 208c87cf4dSdanfinn options->Np = 100; 218c87cf4dSdanfinn 228c87cf4dSdanfinn ierr = PetscOptionsBegin(comm, "", "Test of colorized scatter plot", "");CHKERRQ(ierr); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-Np", "Number of particles", "ex35.c", options->Np, &options->Np, PETSC_NULL)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-dim", "Number of dimensions", "ex35.c", options->dim_inp, &options->dim_inp, PETSC_NULL)); 258c87cf4dSdanfinn ierr = PetscOptionsEnd();CHKERRQ(ierr); 268c87cf4dSdanfinn PetscFunctionReturn(0); 278c87cf4dSdanfinn } 288c87cf4dSdanfinn 298c87cf4dSdanfinn /* 308c87cf4dSdanfinn ref: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/ 318c87cf4dSdanfinn */ 328c87cf4dSdanfinn PetscReal erfinv(PetscReal x) 338c87cf4dSdanfinn { 348c87cf4dSdanfinn PetscReal *ck, r = 0.; 35*5f80ce2aSJacob Faibussowitsch PetscInt maxIter = 100; 368c87cf4dSdanfinn 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(maxIter,&ck)); 388c87cf4dSdanfinn ck[0] = 1; 398c87cf4dSdanfinn r = ck[0]*((PetscSqrtReal(PETSC_PI)/2.)*x); 40*5f80ce2aSJacob Faibussowitsch for (PetscInt k = 1; k < maxIter; ++k){ 41*5f80ce2aSJacob Faibussowitsch for (PetscInt m = 0; m <= k-1; ++m){ 428c87cf4dSdanfinn PetscReal denom = (m+1.)*(2.*m+1.); 438c87cf4dSdanfinn ck[k] += (ck[m]*ck[k-1-m])/denom; 448c87cf4dSdanfinn } 458c87cf4dSdanfinn PetscReal temp = 2.*k+1.; 468c87cf4dSdanfinn r += (ck[k]/temp)*PetscPowReal((PetscSqrtReal(PETSC_PI)/2.)*x,2.*k+1.); 478c87cf4dSdanfinn } 48*5f80ce2aSJacob Faibussowitsch CHKERRABORT(PETSC_COMM_SELF,PetscFree(ck)); 498c87cf4dSdanfinn return r; 508c87cf4dSdanfinn } 518c87cf4dSdanfinn 528c87cf4dSdanfinn int main(int argc, char **argv) 538c87cf4dSdanfinn { 548c87cf4dSdanfinn PetscInt p, dim, Np; 558c87cf4dSdanfinn PetscScalar *randVecNums; 568c87cf4dSdanfinn PetscReal speed, value, *x, *v; 578c87cf4dSdanfinn PetscRandom rngx, rng1, rng2; 588c87cf4dSdanfinn Vec randVec, subvecvx, subvecvy; 598c87cf4dSdanfinn IS isvx, isvy; 608c87cf4dSdanfinn AppCtx user; 618c87cf4dSdanfinn PetscDrawAxis axis; 628c87cf4dSdanfinn PetscDraw positionDraw; 638c87cf4dSdanfinn PetscDrawSP positionDrawSP; 648c87cf4dSdanfinn MPI_Comm comm; 658c87cf4dSdanfinn PetscErrorCode ierr; 668c87cf4dSdanfinn 678c87cf4dSdanfinn ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 688c87cf4dSdanfinn comm = PETSC_COMM_WORLD; 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 708c87cf4dSdanfinn 718c87cf4dSdanfinn Np = user.Np; 728c87cf4dSdanfinn dim = user.dim; 738c87cf4dSdanfinn 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Np*dim, &x, Np*dim, &v)); 758c87cf4dSdanfinn 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(comm, &rngx)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rngx, 0., 1.)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rngx)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetSeed(rngx, 1034)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSeed(rngx)); 818c87cf4dSdanfinn 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(comm, &rng1)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rng1, 0., 1.)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rng1)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetSeed(rng1, 3084)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSeed(rng1)); 878c87cf4dSdanfinn 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(comm, &rng2)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rng2, 0., 1.)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rng2)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetSeed(rng2, 2397)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSeed(rng2)); 938c87cf4dSdanfinn 948c87cf4dSdanfinn /* Set particle positions and velocities */ 958c87cf4dSdanfinn if (user.dim_inp == 1) { 968c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 978c87cf4dSdanfinn PetscReal temp; 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rngx, &value)); 998c87cf4dSdanfinn x[p*dim] = value; 1008c87cf4dSdanfinn x[p*dim+1] = 0.; 1018c87cf4dSdanfinn temp = erfinv(2*value-1); 1028c87cf4dSdanfinn v[p*dim] = temp; 1038c87cf4dSdanfinn v[p*dim+1] = 0.; 1048c87cf4dSdanfinn } 1058c87cf4dSdanfinn } else if (user.dim_inp == 2) { 1068c87cf4dSdanfinn /* 1078c87cf4dSdanfinn Use Box-Muller to sample a distribution of velocities for the maxwellian. 1088c87cf4dSdanfinn https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform 1098c87cf4dSdanfinn */ 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(comm,&randVec)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(randVec,PETSC_DECIDE, Np*dim)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(randVec)); 1138c87cf4dSdanfinn 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(comm, Np*dim/2, 0, 2, &isvx)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(comm, Np*dim/2, 1, 2, &isvy)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSubVector(randVec, isvx, &subvecvx)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSubVector(randVec, isvy, &subvecvy)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(subvecvx, rng1)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(subvecvy, rng2)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreSubVector(randVec, isvx, &subvecvx)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreSubVector(randVec, isvy, &subvecvy)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(randVec, &randVecNums)); 1238c87cf4dSdanfinn 1248c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 1258c87cf4dSdanfinn PetscReal u1, u2, mag, zx, zy; 1268c87cf4dSdanfinn 1278c87cf4dSdanfinn u1 = PetscRealPart(randVecNums[p*dim]); 1288c87cf4dSdanfinn u2 = PetscRealPart(randVecNums[p*dim+1]); 1298c87cf4dSdanfinn 1308c87cf4dSdanfinn x[p*dim] = u1; 1318c87cf4dSdanfinn x[p*dim+1] = u2; 1328c87cf4dSdanfinn 1338c87cf4dSdanfinn mag = PetscSqrtReal(-2.0 * PetscLogReal(u1)); 1348c87cf4dSdanfinn 1358c87cf4dSdanfinn zx = mag * PetscCosReal(2*PETSC_PI*u2) + 0; 1368c87cf4dSdanfinn zy = mag * PetscSinReal(2*PETSC_PI*u2) + 0; 1378c87cf4dSdanfinn 1388c87cf4dSdanfinn v[p*dim] = zx; 1398c87cf4dSdanfinn v[p*dim+1] = zy; 1408c87cf4dSdanfinn } 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isvx)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isvy)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&subvecvx)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&subvecvy)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&randVec)); 1468c87cf4dSdanfinn } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %D", dim); 1478c87cf4dSdanfinn 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawCreate(comm, NULL, "monitor_particle_positions", 0,0,400,300, &positionDraw)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetFromOptions(positionDraw)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPCreate(positionDraw, 10, &positionDrawSP)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPSetDimension(positionDrawSP,1)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPGetAxis(positionDrawSP, &axis)); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPReset(positionDrawSP)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis, "Particles", "x", "y")); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetSave(positionDraw, "ex35_pos.ppm")); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPReset(positionDrawSP)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPSetLimits(positionDrawSP, 0, 1, 0, 1)); 1588c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 1598c87cf4dSdanfinn speed = PetscSqrtReal(PetscSqr(v[p*dim]) + PetscSqr(v[p*dim+1])); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPAddPointColorized(positionDrawSP, &x[p*dim], &x[p*dim+1], &speed)); 1618c87cf4dSdanfinn } 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPDraw(positionDrawSP, PETSC_TRUE)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSave(positionDraw)); 1648c87cf4dSdanfinn 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(x, v)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rngx)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rng1)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rng2)); 1698c87cf4dSdanfinn 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPDestroy(&positionDrawSP)); 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawDestroy(&positionDraw)); 1728c87cf4dSdanfinn ierr = PetscFinalize(); 1738c87cf4dSdanfinn return ierr; 1748c87cf4dSdanfinn } 1758c87cf4dSdanfinn 1768c87cf4dSdanfinn /*TEST 1778c87cf4dSdanfinn test: 1788c87cf4dSdanfinn suffix: 1D 1798c87cf4dSdanfinn args: -Np 50\ 1808c87cf4dSdanfinn -dim 1 1818c87cf4dSdanfinn test: 1828c87cf4dSdanfinn suffix: 2D 1838c87cf4dSdanfinn args: -Np 50\ 1848c87cf4dSdanfinn -dim 2 1858c87cf4dSdanfinn TEST*/ 186