18c87cf4dSdanfinn static char help[] = "Test of Colorized Scatter Plot.\n";
28c87cf4dSdanfinn
346233b44SBarry Smith #include <petscdraw.h>
446233b44SBarry Smith #include <petscvec.h>
546233b44SBarry Smith #include <petscis.h>
68c87cf4dSdanfinn
78c87cf4dSdanfinn typedef struct {
88c87cf4dSdanfinn PetscInt Np; /* total number of particles */
98c87cf4dSdanfinn PetscInt dim;
108c87cf4dSdanfinn PetscInt dim_inp;
118c87cf4dSdanfinn } AppCtx;
128c87cf4dSdanfinn
ProcessOptions(MPI_Comm comm,AppCtx * options)13d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
14d71ae5a4SJacob Faibussowitsch {
158c87cf4dSdanfinn PetscFunctionBeginUser;
168c87cf4dSdanfinn options->dim = 2;
178c87cf4dSdanfinn options->dim_inp = 2;
188c87cf4dSdanfinn options->Np = 100;
19d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Test of colorized scatter plot", "");
20f3fa974cSJacob Faibussowitsch PetscCall(PetscOptionsInt("-Np", "Number of particles", "ex35.c", options->Np, &options->Np, NULL));
21f3fa974cSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dim", "Number of dimensions", "ex35.c", options->dim_inp, &options->dim_inp, NULL));
22d0609cedSBarry Smith PetscOptionsEnd();
233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
248c87cf4dSdanfinn }
258c87cf4dSdanfinn
268c87cf4dSdanfinn /*
278c87cf4dSdanfinn ref: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/
288c87cf4dSdanfinn */
erfinv(PetscReal x)29d71ae5a4SJacob Faibussowitsch PetscReal erfinv(PetscReal x)
30d71ae5a4SJacob Faibussowitsch {
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
main(int argc,char ** argv)49d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
50d71ae5a4SJacob Faibussowitsch {
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
63327415f7SBarry Smith PetscFunctionBeginUser;
649566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
658c87cf4dSdanfinn comm = PETSC_COMM_WORLD;
669566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user));
678c87cf4dSdanfinn
688c87cf4dSdanfinn Np = user.Np;
698c87cf4dSdanfinn dim = user.dim;
708c87cf4dSdanfinn
719566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Np * dim, &x, Np * dim, &v));
728c87cf4dSdanfinn
739566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rngx));
749566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rngx, 0., 1.));
759566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rngx));
769566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rngx, 1034));
779566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(rngx));
788c87cf4dSdanfinn
799566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rng1));
809566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rng1, 0., 1.));
819566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rng1));
829566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rng1, 3084));
839566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(rng1));
848c87cf4dSdanfinn
859566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rng2));
869566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rng2, 0., 1.));
879566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rng2));
889566063dSJacob Faibussowitsch PetscCall(PetscRandomSetSeed(rng2, 2397));
899566063dSJacob Faibussowitsch PetscCall(PetscRandomSeed(rng2));
908c87cf4dSdanfinn
918c87cf4dSdanfinn /* Set particle positions and velocities */
928c87cf4dSdanfinn if (user.dim_inp == 1) {
938c87cf4dSdanfinn for (p = 0; p < Np; ++p) {
948c87cf4dSdanfinn PetscReal temp;
959566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rngx, &value));
968c87cf4dSdanfinn x[p * dim] = value;
978c87cf4dSdanfinn x[p * dim + 1] = 0.;
988c87cf4dSdanfinn temp = erfinv(2 * value - 1);
998c87cf4dSdanfinn v[p * dim] = temp;
1008c87cf4dSdanfinn v[p * dim + 1] = 0.;
1018c87cf4dSdanfinn }
1028c87cf4dSdanfinn } else if (user.dim_inp == 2) {
1038c87cf4dSdanfinn /*
1048c87cf4dSdanfinn Use Box-Muller to sample a distribution of velocities for the maxwellian.
1058c87cf4dSdanfinn https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
1068c87cf4dSdanfinn */
1079566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &randVec));
1089566063dSJacob Faibussowitsch PetscCall(VecSetSizes(randVec, PETSC_DECIDE, Np * dim));
1099566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(randVec));
1108c87cf4dSdanfinn
1119566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Np * dim / 2, 0, 2, &isvx));
1129566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, Np * dim / 2, 1, 2, &isvy));
1139566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(randVec, isvx, &subvecvx));
1149566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(randVec, isvy, &subvecvy));
1159566063dSJacob Faibussowitsch PetscCall(VecSetRandom(subvecvx, rng1));
1169566063dSJacob Faibussowitsch PetscCall(VecSetRandom(subvecvy, rng2));
1179566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(randVec, isvx, &subvecvx));
1189566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(randVec, isvy, &subvecvy));
1199566063dSJacob Faibussowitsch PetscCall(VecGetArray(randVec, &randVecNums));
1208c87cf4dSdanfinn
1218c87cf4dSdanfinn for (p = 0; p < Np; ++p) {
1228c87cf4dSdanfinn PetscReal u1, u2, mag, zx, zy;
1238c87cf4dSdanfinn
1248c87cf4dSdanfinn u1 = PetscRealPart(randVecNums[p * dim]);
1258c87cf4dSdanfinn u2 = PetscRealPart(randVecNums[p * dim + 1]);
1268c87cf4dSdanfinn
1278c87cf4dSdanfinn x[p * dim] = u1;
1288c87cf4dSdanfinn x[p * dim + 1] = u2;
1298c87cf4dSdanfinn
1308c87cf4dSdanfinn mag = PetscSqrtReal(-2.0 * PetscLogReal(u1));
1318c87cf4dSdanfinn
1328c87cf4dSdanfinn zx = mag * PetscCosReal(2 * PETSC_PI * u2) + 0;
1338c87cf4dSdanfinn zy = mag * PetscSinReal(2 * PETSC_PI * u2) + 0;
1348c87cf4dSdanfinn
1358c87cf4dSdanfinn v[p * dim] = zx;
1368c87cf4dSdanfinn v[p * dim + 1] = zy;
1378c87cf4dSdanfinn }
1389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isvx));
1399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isvy));
1409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&subvecvx));
1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&subvecvy));
1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&randVec));
14363a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %" PetscInt_FMT, dim);
1448c87cf4dSdanfinn
1459566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, NULL, "monitor_particle_positions", 0, 0, 400, 300, &positionDraw));
1469566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(positionDraw));
1479566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(positionDraw, 10, &positionDrawSP));
1489566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetDimension(positionDrawSP, 1));
1499566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(positionDrawSP, &axis));
1509566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(positionDrawSP));
1519566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "x", "y"));
1529566063dSJacob Faibussowitsch PetscCall(PetscDrawSetSave(positionDraw, "ex35_pos.ppm"));
1539566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(positionDrawSP));
1549566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetLimits(positionDrawSP, 0, 1, 0, 1));
1558c87cf4dSdanfinn for (p = 0; p < Np; ++p) {
1568c87cf4dSdanfinn speed = PetscSqrtReal(PetscSqr(v[p * dim]) + PetscSqr(v[p * dim + 1]));
1579566063dSJacob Faibussowitsch PetscCall(PetscDrawSPAddPointColorized(positionDrawSP, &x[p * dim], &x[p * dim + 1], &speed));
1588c87cf4dSdanfinn }
1599566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(positionDrawSP, PETSC_TRUE));
1609566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(positionDraw));
1618c87cf4dSdanfinn
1629566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, v));
1639566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rngx));
1649566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rng1));
1659566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rng2));
1668c87cf4dSdanfinn
1679566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&positionDrawSP));
1689566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&positionDraw));
1699566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
170b122ec5aSJacob Faibussowitsch return 0;
1718c87cf4dSdanfinn }
1728c87cf4dSdanfinn
1738c87cf4dSdanfinn /*TEST
1748c87cf4dSdanfinn test:
1758c87cf4dSdanfinn suffix: 1D
1768c87cf4dSdanfinn args: -Np 50\
1778c87cf4dSdanfinn -dim 1
178*3886731fSPierre Jolivet output_file: output/empty.out
1798c87cf4dSdanfinn test:
1808c87cf4dSdanfinn suffix: 2D
1818c87cf4dSdanfinn args: -Np 50\
1828c87cf4dSdanfinn -dim 2
183*3886731fSPierre Jolivet output_file: output/empty.out
1848c87cf4dSdanfinn TEST*/
185