xref: /petsc/src/ts/tests/ex35.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
22*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Test of colorized scatter plot", "");PetscCall(ierr);
23*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-Np", "Number of particles", "ex35.c", options->Np, &options->Np, PETSC_NULL));
24*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dim", "Number of dimensions", "ex35.c", options->dim_inp, &options->dim_inp, PETSC_NULL));
25*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(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.;
355f80ce2aSJacob Faibussowitsch   PetscInt   maxIter = 100;
368c87cf4dSdanfinn 
37*9566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(maxIter,&ck));
388c87cf4dSdanfinn   ck[0] = 1;
398c87cf4dSdanfinn   r = ck[0]*((PetscSqrtReal(PETSC_PI)/2.)*x);
405f80ce2aSJacob Faibussowitsch   for (PetscInt k = 1; k < maxIter; ++k){
415f80ce2aSJacob 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*9566063dSJacob Faibussowitsch   PetscCallAbort(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 
66*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
678c87cf4dSdanfinn   comm = PETSC_COMM_WORLD;
68*9566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
698c87cf4dSdanfinn 
708c87cf4dSdanfinn   Np = user.Np;
718c87cf4dSdanfinn   dim = user.dim;
728c87cf4dSdanfinn 
73*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Np*dim, &x, Np*dim, &v));
748c87cf4dSdanfinn 
75*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &rngx));
76*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rngx, 0., 1.));
77*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rngx));
78*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetSeed(rngx, 1034));
79*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSeed(rngx));
808c87cf4dSdanfinn 
81*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &rng1));
82*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rng1, 0., 1.));
83*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rng1));
84*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetSeed(rng1, 3084));
85*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSeed(rng1));
868c87cf4dSdanfinn 
87*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &rng2));
88*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rng2, 0., 1.));
89*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rng2));
90*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetSeed(rng2, 2397));
91*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSeed(rng2));
928c87cf4dSdanfinn 
938c87cf4dSdanfinn   /* Set particle positions and velocities */
948c87cf4dSdanfinn   if (user.dim_inp == 1) {
958c87cf4dSdanfinn     for (p = 0; p < Np; ++p) {
968c87cf4dSdanfinn       PetscReal temp;
97*9566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rngx, &value));
988c87cf4dSdanfinn       x[p*dim] = value;
998c87cf4dSdanfinn       x[p*dim+1] = 0.;
1008c87cf4dSdanfinn       temp = erfinv(2*value-1);
1018c87cf4dSdanfinn       v[p*dim] = temp;
1028c87cf4dSdanfinn       v[p*dim+1] = 0.;
1038c87cf4dSdanfinn     }
1048c87cf4dSdanfinn   } else if (user.dim_inp == 2) {
1058c87cf4dSdanfinn     /*
1068c87cf4dSdanfinn      Use Box-Muller to sample a distribution of velocities for the maxwellian.
1078c87cf4dSdanfinn      https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
1088c87cf4dSdanfinn     */
109*9566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm,&randVec));
110*9566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(randVec,PETSC_DECIDE, Np*dim));
111*9566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(randVec));
1128c87cf4dSdanfinn 
113*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, Np*dim/2, 0, 2, &isvx));
114*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, Np*dim/2, 1, 2, &isvy));
115*9566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(randVec, isvx, &subvecvx));
116*9566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(randVec, isvy, &subvecvy));
117*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(subvecvx, rng1));
118*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(subvecvy, rng2));
119*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(randVec, isvx, &subvecvx));
120*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(randVec, isvy, &subvecvy));
121*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(randVec, &randVecNums));
1228c87cf4dSdanfinn 
1238c87cf4dSdanfinn     for (p = 0; p < Np; ++p) {
1248c87cf4dSdanfinn       PetscReal u1, u2, mag, zx, zy;
1258c87cf4dSdanfinn 
1268c87cf4dSdanfinn       u1 = PetscRealPart(randVecNums[p*dim]);
1278c87cf4dSdanfinn       u2 = PetscRealPart(randVecNums[p*dim+1]);
1288c87cf4dSdanfinn 
1298c87cf4dSdanfinn       x[p*dim] = u1;
1308c87cf4dSdanfinn       x[p*dim+1] = u2;
1318c87cf4dSdanfinn 
1328c87cf4dSdanfinn       mag = PetscSqrtReal(-2.0 * PetscLogReal(u1));
1338c87cf4dSdanfinn 
1348c87cf4dSdanfinn       zx = mag * PetscCosReal(2*PETSC_PI*u2) + 0;
1358c87cf4dSdanfinn       zy = mag * PetscSinReal(2*PETSC_PI*u2) + 0;
1368c87cf4dSdanfinn 
1378c87cf4dSdanfinn       v[p*dim] = zx;
1388c87cf4dSdanfinn       v[p*dim+1] = zy;
1398c87cf4dSdanfinn     }
140*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isvx));
141*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isvy));
142*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&subvecvx));
143*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&subvecvy));
144*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&randVec));
1458c87cf4dSdanfinn   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %D", dim);
1468c87cf4dSdanfinn 
147*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawCreate(comm, NULL, "monitor_particle_positions", 0,0,400,300, &positionDraw));
148*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetFromOptions(positionDraw));
149*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPCreate(positionDraw, 10, &positionDrawSP));
150*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPSetDimension(positionDrawSP,1));
151*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPGetAxis(positionDrawSP, &axis));
152*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPReset(positionDrawSP));
153*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "x", "y"));
154*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetSave(positionDraw, "ex35_pos.ppm"));
155*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPReset(positionDrawSP));
156*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPSetLimits(positionDrawSP, 0, 1, 0, 1));
1578c87cf4dSdanfinn   for (p = 0; p < Np; ++p) {
1588c87cf4dSdanfinn     speed = PetscSqrtReal(PetscSqr(v[p*dim]) + PetscSqr(v[p*dim+1]));
159*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPAddPointColorized(positionDrawSP, &x[p*dim], &x[p*dim+1], &speed));
1608c87cf4dSdanfinn   }
161*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDraw(positionDrawSP, PETSC_TRUE));
162*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(positionDraw));
1638c87cf4dSdanfinn 
164*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, v));
165*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rngx));
166*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rng1));
167*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rng2));
1688c87cf4dSdanfinn 
169*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSPDestroy(&positionDrawSP));
170*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawDestroy(&positionDraw));
171*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
172b122ec5aSJacob Faibussowitsch   return 0;
1738c87cf4dSdanfinn }
1748c87cf4dSdanfinn 
1758c87cf4dSdanfinn /*TEST
1768c87cf4dSdanfinn    test:
1778c87cf4dSdanfinn      suffix: 1D
1788c87cf4dSdanfinn      args: -Np 50\
1798c87cf4dSdanfinn      -dim 1
1808c87cf4dSdanfinn    test:
1818c87cf4dSdanfinn      suffix: 2D
1828c87cf4dSdanfinn      args: -Np 50\
1838c87cf4dSdanfinn      -dim 2
1848c87cf4dSdanfinn TEST*/
185