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 PetscFunctionBeginUser; 168c87cf4dSdanfinn options->dim = 2; 178c87cf4dSdanfinn options->dim_inp = 2; 188c87cf4dSdanfinn options->Np = 100; 19d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Test of colorized scatter plot", ""); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-Np", "Number of particles", "ex35.c", options->Np, &options->Np, PETSC_NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dim", "Number of dimensions", "ex35.c", options->dim_inp, &options->dim_inp, PETSC_NULL)); 22d0609cedSBarry Smith PetscOptionsEnd(); 238c87cf4dSdanfinn PetscFunctionReturn(0); 248c87cf4dSdanfinn } 258c87cf4dSdanfinn 268c87cf4dSdanfinn /* 278c87cf4dSdanfinn ref: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/ 288c87cf4dSdanfinn */ 298c87cf4dSdanfinn PetscReal erfinv(PetscReal x) 308c87cf4dSdanfinn { 318c87cf4dSdanfinn PetscReal *ck, r = 0.; 325f80ce2aSJacob Faibussowitsch PetscInt maxIter = 100; 338c87cf4dSdanfinn 349566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(maxIter,&ck)); 358c87cf4dSdanfinn ck[0] = 1; 368c87cf4dSdanfinn r = ck[0]*((PetscSqrtReal(PETSC_PI)/2.)*x); 375f80ce2aSJacob Faibussowitsch for (PetscInt k = 1; k < maxIter; ++k){ 385f80ce2aSJacob Faibussowitsch for (PetscInt m = 0; m <= k-1; ++m){ 398c87cf4dSdanfinn PetscReal denom = (m+1.)*(2.*m+1.); 408c87cf4dSdanfinn ck[k] += (ck[m]*ck[k-1-m])/denom; 418c87cf4dSdanfinn } 428c87cf4dSdanfinn PetscReal temp = 2.*k+1.; 438c87cf4dSdanfinn r += (ck[k]/temp)*PetscPowReal((PetscSqrtReal(PETSC_PI)/2.)*x,2.*k+1.); 448c87cf4dSdanfinn } 459566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,PetscFree(ck)); 468c87cf4dSdanfinn return r; 478c87cf4dSdanfinn } 488c87cf4dSdanfinn 498c87cf4dSdanfinn int main(int argc, char **argv) 508c87cf4dSdanfinn { 518c87cf4dSdanfinn PetscInt p, dim, Np; 528c87cf4dSdanfinn PetscScalar *randVecNums; 538c87cf4dSdanfinn PetscReal speed, value, *x, *v; 548c87cf4dSdanfinn PetscRandom rngx, rng1, rng2; 558c87cf4dSdanfinn Vec randVec, subvecvx, subvecvy; 568c87cf4dSdanfinn IS isvx, isvy; 578c87cf4dSdanfinn AppCtx user; 588c87cf4dSdanfinn PetscDrawAxis axis; 598c87cf4dSdanfinn PetscDraw positionDraw; 608c87cf4dSdanfinn PetscDrawSP positionDrawSP; 618c87cf4dSdanfinn MPI_Comm comm; 628c87cf4dSdanfinn 639566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 648c87cf4dSdanfinn comm = PETSC_COMM_WORLD; 659566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 668c87cf4dSdanfinn 678c87cf4dSdanfinn Np = user.Np; 688c87cf4dSdanfinn dim = user.dim; 698c87cf4dSdanfinn 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Np*dim, &x, Np*dim, &v)); 718c87cf4dSdanfinn 729566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rngx)); 739566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rngx, 0., 1.)); 749566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rngx)); 759566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rngx, 1034)); 769566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(rngx)); 778c87cf4dSdanfinn 789566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rng1)); 799566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rng1, 0., 1.)); 809566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rng1)); 819566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rng1, 3084)); 829566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(rng1)); 838c87cf4dSdanfinn 849566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rng2)); 859566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rng2, 0., 1.)); 869566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rng2)); 879566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rng2, 2397)); 889566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(rng2)); 898c87cf4dSdanfinn 908c87cf4dSdanfinn /* Set particle positions and velocities */ 918c87cf4dSdanfinn if (user.dim_inp == 1) { 928c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 938c87cf4dSdanfinn PetscReal temp; 949566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rngx, &value)); 958c87cf4dSdanfinn x[p*dim] = value; 968c87cf4dSdanfinn x[p*dim+1] = 0.; 978c87cf4dSdanfinn temp = erfinv(2*value-1); 988c87cf4dSdanfinn v[p*dim] = temp; 998c87cf4dSdanfinn v[p*dim+1] = 0.; 1008c87cf4dSdanfinn } 1018c87cf4dSdanfinn } else if (user.dim_inp == 2) { 1028c87cf4dSdanfinn /* 1038c87cf4dSdanfinn Use Box-Muller to sample a distribution of velocities for the maxwellian. 1048c87cf4dSdanfinn https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform 1058c87cf4dSdanfinn */ 1069566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&randVec)); 1079566063dSJacob Faibussowitsch PetscCall(VecSetSizes(randVec,PETSC_DECIDE, Np*dim)); 1089566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(randVec)); 1098c87cf4dSdanfinn 1109566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Np*dim/2, 0, 2, &isvx)); 1119566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Np*dim/2, 1, 2, &isvy)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(randVec, isvx, &subvecvx)); 1139566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(randVec, isvy, &subvecvy)); 1149566063dSJacob Faibussowitsch PetscCall(VecSetRandom(subvecvx, rng1)); 1159566063dSJacob Faibussowitsch PetscCall(VecSetRandom(subvecvy, rng2)); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(randVec, isvx, &subvecvx)); 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(randVec, isvy, &subvecvy)); 1189566063dSJacob Faibussowitsch PetscCall(VecGetArray(randVec, &randVecNums)); 1198c87cf4dSdanfinn 1208c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 1218c87cf4dSdanfinn PetscReal u1, u2, mag, zx, zy; 1228c87cf4dSdanfinn 1238c87cf4dSdanfinn u1 = PetscRealPart(randVecNums[p*dim]); 1248c87cf4dSdanfinn u2 = PetscRealPart(randVecNums[p*dim+1]); 1258c87cf4dSdanfinn 1268c87cf4dSdanfinn x[p*dim] = u1; 1278c87cf4dSdanfinn x[p*dim+1] = u2; 1288c87cf4dSdanfinn 1298c87cf4dSdanfinn mag = PetscSqrtReal(-2.0 * PetscLogReal(u1)); 1308c87cf4dSdanfinn 1318c87cf4dSdanfinn zx = mag * PetscCosReal(2*PETSC_PI*u2) + 0; 1328c87cf4dSdanfinn zy = mag * PetscSinReal(2*PETSC_PI*u2) + 0; 1338c87cf4dSdanfinn 1348c87cf4dSdanfinn v[p*dim] = zx; 1358c87cf4dSdanfinn v[p*dim+1] = zy; 1368c87cf4dSdanfinn } 1379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isvx)); 1389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isvy)); 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&subvecvx)); 1409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&subvecvy)); 1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&randVec)); 142*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %" PetscInt_FMT, dim); 1438c87cf4dSdanfinn 1449566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, NULL, "monitor_particle_positions", 0,0,400,300, &positionDraw)); 1459566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(positionDraw)); 1469566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(positionDraw, 10, &positionDrawSP)); 1479566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetDimension(positionDrawSP,1)); 1489566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(positionDrawSP, &axis)); 1499566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(positionDrawSP)); 1509566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "x", "y")); 1519566063dSJacob Faibussowitsch PetscCall(PetscDrawSetSave(positionDraw, "ex35_pos.ppm")); 1529566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(positionDrawSP)); 1539566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetLimits(positionDrawSP, 0, 1, 0, 1)); 1548c87cf4dSdanfinn for (p = 0; p < Np; ++p) { 1558c87cf4dSdanfinn speed = PetscSqrtReal(PetscSqr(v[p*dim]) + PetscSqr(v[p*dim+1])); 1569566063dSJacob Faibussowitsch PetscCall(PetscDrawSPAddPointColorized(positionDrawSP, &x[p*dim], &x[p*dim+1], &speed)); 1578c87cf4dSdanfinn } 1589566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(positionDrawSP, PETSC_TRUE)); 1599566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(positionDraw)); 1608c87cf4dSdanfinn 1619566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, v)); 1629566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rngx)); 1639566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rng1)); 1649566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rng2)); 1658c87cf4dSdanfinn 1669566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&positionDrawSP)); 1679566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&positionDraw)); 1689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 169b122ec5aSJacob Faibussowitsch return 0; 1708c87cf4dSdanfinn } 1718c87cf4dSdanfinn 1728c87cf4dSdanfinn /*TEST 1738c87cf4dSdanfinn test: 1748c87cf4dSdanfinn suffix: 1D 1758c87cf4dSdanfinn args: -Np 50\ 1768c87cf4dSdanfinn -dim 1 1778c87cf4dSdanfinn test: 1788c87cf4dSdanfinn suffix: 2D 1798c87cf4dSdanfinn args: -Np 50\ 1808c87cf4dSdanfinn -dim 2 1818c87cf4dSdanfinn TEST*/ 182