xref: /petsc/src/ts/tests/ex28.c (revision 2cf5aabc426dcedfaf8346422880b98b77869311)
1d0c080abSJoseph Pusztay static char help[] = "Example application of the Bhatnagar-Gross-Krook (BGK) collision operator.\n\
2d0c080abSJoseph Pusztay This example is a 0D-1V setting for the kinetic equation\n\
3d0c080abSJoseph Pusztay https://en.wikipedia.org/wiki/Bhatnagar%E2%80%93Gross%E2%80%93Krook_operator\n";
4d0c080abSJoseph Pusztay 
5d0c080abSJoseph Pusztay #include <petscdmplex.h>
6d0c080abSJoseph Pusztay #include <petscdmswarm.h>
7d0c080abSJoseph Pusztay #include <petscts.h>
8d0c080abSJoseph Pusztay #include <petscdraw.h>
9d0c080abSJoseph Pusztay #include <petscviewer.h>
10d0c080abSJoseph Pusztay 
11d0c080abSJoseph Pusztay typedef struct {
12d0c080abSJoseph Pusztay   PetscInt    particlesPerCell; /* The number of partices per cell */
13d0c080abSJoseph Pusztay   PetscReal   momentTol;        /* Tolerance for checking moment conservation */
14d0c080abSJoseph Pusztay   PetscBool   monitorhg;        /* Flag for using the TS histogram monitor */
15d0c080abSJoseph Pusztay   PetscBool   monitorsp;        /* Flag for using the TS scatter monitor */
16d0c080abSJoseph Pusztay   PetscBool   monitorks;        /* Monitor to perform KS test to the maxwellian */
17d0c080abSJoseph Pusztay   PetscBool   error;            /* Flag for printing the error */
18d0c080abSJoseph Pusztay   PetscInt    ostep;            /* print the energy at each ostep time steps */
19d0c080abSJoseph Pusztay   PetscDraw   draw;             /* The draw object for histogram monitoring */
20d0c080abSJoseph Pusztay   PetscDrawHG drawhg;           /* The histogram draw context for monitoring */
21d0c080abSJoseph Pusztay   PetscDrawSP drawsp;           /* The scatter plot draw context for the monitor */
22d0c080abSJoseph Pusztay   PetscDrawSP drawks;           /* Scatterplot draw context for KS test */
23d0c080abSJoseph Pusztay } AppCtx;
24d0c080abSJoseph Pusztay 
25d0c080abSJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
26d0c080abSJoseph Pusztay {
27d0c080abSJoseph Pusztay   PetscErrorCode ierr;
28d0c080abSJoseph Pusztay 
29d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
30d0c080abSJoseph Pusztay   options->monitorhg        = PETSC_FALSE;
31d0c080abSJoseph Pusztay   options->monitorsp        = PETSC_FALSE;
32d0c080abSJoseph Pusztay   options->monitorks        = PETSC_FALSE;
33d0c080abSJoseph Pusztay   options->particlesPerCell = 1;
34d0c080abSJoseph Pusztay   options->momentTol        = 100.0*PETSC_MACHINE_EPSILON;
35d0c080abSJoseph Pusztay   options->ostep            = 100;
36d0c080abSJoseph Pusztay 
37d0c080abSJoseph Pusztay   ierr = PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");CHKERRQ(ierr);
38d0c080abSJoseph Pusztay   ierr = PetscOptionsBool("-monitorhg", "Flag to use the TS histogram monitor", "ex28.c", options->monitorhg, &options->monitorhg, NULL);CHKERRQ(ierr);
39d0c080abSJoseph Pusztay   ierr = PetscOptionsBool("-monitorsp", "Flag to use the TS scatter plot monitor", "ex28.c", options->monitorsp, &options->monitorsp, NULL);CHKERRQ(ierr);
40d0c080abSJoseph Pusztay   ierr = PetscOptionsBool("-monitorks", "Flag to plot KS test results", "ex28.c", options->monitorks, &options->monitorks, NULL);CHKERRQ(ierr);
41d0c080abSJoseph Pusztay   ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex28.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr);
42d0c080abSJoseph Pusztay   ierr = PetscOptionsInt("-output_step", "Number of time steps between output", "ex28.c", options->ostep, &options->ostep, PETSC_NULL);CHKERRQ(ierr);
43d0c080abSJoseph Pusztay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
44d0c080abSJoseph Pusztay 
45d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
46d0c080abSJoseph Pusztay }
47d0c080abSJoseph Pusztay 
48d0c080abSJoseph Pusztay /* Create the mesh for velocity space */
49d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
50d0c080abSJoseph Pusztay {
51d0c080abSJoseph Pusztay   PetscErrorCode ierr;
52d0c080abSJoseph Pusztay 
53d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
5430602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
5530602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
56d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
57d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
58d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
59d0c080abSJoseph Pusztay }
60d0c080abSJoseph Pusztay 
61d0c080abSJoseph Pusztay /* Since we are putting the same number of particles in each cell, this amounts to a uniform distribution of v */
62d0c080abSJoseph Pusztay static PetscErrorCode SetInitialCoordinates(DM sw)
63d0c080abSJoseph Pusztay {
64d0c080abSJoseph Pusztay   AppCtx        *user;
65d0c080abSJoseph Pusztay   PetscRandom    rnd;
66d0c080abSJoseph Pusztay   DM             dm;
67d0c080abSJoseph Pusztay   DMPolytopeType ct;
68d0c080abSJoseph Pusztay   PetscBool      simplex;
69d0c080abSJoseph Pusztay   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *vals;
70d0c080abSJoseph Pusztay   PetscInt       dim, d, cStart, cEnd, c, Np, p;
71d0c080abSJoseph Pusztay   PetscErrorCode ierr;
72d0c080abSJoseph Pusztay 
73d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
74d0c080abSJoseph Pusztay   ierr = PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rnd);CHKERRQ(ierr);
75d0c080abSJoseph Pusztay   ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr);
76d0c080abSJoseph Pusztay   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
77d0c080abSJoseph Pusztay 
783ec1f749SStefano Zampini   ierr = DMGetApplicationContext(sw, &user);CHKERRQ(ierr);
79d0c080abSJoseph Pusztay   Np   = user->particlesPerCell;
80d0c080abSJoseph Pusztay   ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr);
81d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellDM(sw, &dm);CHKERRQ(ierr);
82d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
83d0c080abSJoseph Pusztay   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
84d0c080abSJoseph Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
85d0c080abSJoseph Pusztay   ierr = PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr);
86d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
87d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
88d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr);
89d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
90d0c080abSJoseph Pusztay     if (Np == 1) {
91d0c080abSJoseph Pusztay       ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
92d0c080abSJoseph Pusztay       for (d = 0; d < dim; ++d) {
93d0c080abSJoseph Pusztay         coords[c*dim+d] = centroid[d];
94d0c080abSJoseph Pusztay         if ((coords[c*dim+d] >= -1) && (coords[c*dim+d] <= 1)) {
95d0c080abSJoseph Pusztay           vals[c] = 1.0;
96d0c080abSJoseph Pusztay         } else {
97d0c080abSJoseph Pusztay           vals[c] = 0.;
98d0c080abSJoseph Pusztay         }
99d0c080abSJoseph Pusztay       }
100d0c080abSJoseph Pusztay     } else {
101d0c080abSJoseph Pusztay       ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */
102d0c080abSJoseph Pusztay       for (p = 0; p < Np; ++p) {
103d0c080abSJoseph Pusztay         const PetscInt n   = c*Np + p;
104d0c080abSJoseph Pusztay         PetscReal      sum = 0.0, refcoords[3];
105d0c080abSJoseph Pusztay 
106d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) {
107d0c080abSJoseph Pusztay           ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr);
108d0c080abSJoseph Pusztay           sum += refcoords[d];
109d0c080abSJoseph Pusztay         }
110d0c080abSJoseph Pusztay         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
111d0c080abSJoseph Pusztay         vals[n] = 1.0;
112d0c080abSJoseph Pusztay         ierr = DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n*dim]);CHKERRQ(ierr);
113d0c080abSJoseph Pusztay       }
114d0c080abSJoseph Pusztay     }
115d0c080abSJoseph Pusztay   }
116d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
117d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr);
118d0c080abSJoseph Pusztay   ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr);
119d0c080abSJoseph Pusztay   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
120d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
121d0c080abSJoseph Pusztay }
122d0c080abSJoseph Pusztay 
123d0c080abSJoseph Pusztay /* The intiial conditions are just the initial particle weights */
124d0c080abSJoseph Pusztay static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
125d0c080abSJoseph Pusztay {
126d0c080abSJoseph Pusztay   DM             dm;
127d0c080abSJoseph Pusztay   AppCtx        *user;
128d0c080abSJoseph Pusztay   PetscReal     *vals;
129d0c080abSJoseph Pusztay   PetscScalar   *initialConditions;
130d0c080abSJoseph Pusztay   PetscInt       dim, d, cStart, cEnd, c, Np, p, n;
131d0c080abSJoseph Pusztay   PetscErrorCode ierr;
132d0c080abSJoseph Pusztay 
133d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
134d0c080abSJoseph Pusztay   ierr = VecGetLocalSize(u, &n);CHKERRQ(ierr);
1353ec1f749SStefano Zampini   ierr = DMGetApplicationContext(dmSw, &user);CHKERRQ(ierr);
136d0c080abSJoseph Pusztay   Np   = user->particlesPerCell;
137d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
138d0c080abSJoseph Pusztay   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
139d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
140d0c080abSJoseph Pusztay   if (n != (cEnd-cStart)*Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "TS solution local size %D != %D nm particles", n, (cEnd-cStart)*Np);
141d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(dmSw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr);
142d0c080abSJoseph Pusztay   ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr);
143d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
144d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
145d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
146d0c080abSJoseph Pusztay       for (d = 0; d < dim; d++) {
147d0c080abSJoseph Pusztay         initialConditions[n] = vals[n];
148d0c080abSJoseph Pusztay       }
149d0c080abSJoseph Pusztay     }
150d0c080abSJoseph Pusztay   }
151d0c080abSJoseph Pusztay   ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr);
152d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(dmSw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr);
153d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
154d0c080abSJoseph Pusztay }
155d0c080abSJoseph Pusztay 
156d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
157d0c080abSJoseph Pusztay {
158d0c080abSJoseph Pusztay   PetscInt      *cellid;
159d0c080abSJoseph Pusztay   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
160d0c080abSJoseph Pusztay   PetscErrorCode ierr;
161d0c080abSJoseph Pusztay 
162d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
163d0c080abSJoseph Pusztay   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
164d0c080abSJoseph Pusztay   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
165d0c080abSJoseph Pusztay   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
166d0c080abSJoseph Pusztay   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
167d0c080abSJoseph Pusztay   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
168d0c080abSJoseph Pusztay   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
169d0c080abSJoseph Pusztay   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL);CHKERRQ(ierr);
170d0c080abSJoseph Pusztay   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR);CHKERRQ(ierr);
171d0c080abSJoseph Pusztay   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
172d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
173d0c080abSJoseph Pusztay   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
174d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
175d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
176d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
177d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
178d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
179d0c080abSJoseph Pusztay       cellid[n] = c;
180d0c080abSJoseph Pusztay     }
181d0c080abSJoseph Pusztay   }
182d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
183d0c080abSJoseph Pusztay   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
184d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
185d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
186d0c080abSJoseph Pusztay }
187d0c080abSJoseph Pusztay 
188d0c080abSJoseph Pusztay /*
189d0c080abSJoseph Pusztay   f_t = 1/\tau \left( f_eq - f \right)
190d0c080abSJoseph Pusztay   n_t = 1/\tau \left( \int f_eq - \int f \right) = 1/\tau (n - n) = 0
191d0c080abSJoseph Pusztay   v_t = 1/\tau \left( \int v f_eq - \int v f \right) = 1/\tau (v - v) = 0
192d0c080abSJoseph Pusztay   E_t = 1/\tau \left( \int v^2 f_eq - \int v^2 f \right) = 1/\tau (T - T) = 0
193d0c080abSJoseph Pusztay 
194d0c080abSJoseph Pusztay   Let's look at a single cell:
195d0c080abSJoseph Pusztay 
196d0c080abSJoseph Pusztay     \int_C f_t             = 1/\tau \left( \int_C f_eq - \int_C f \right)
197d0c080abSJoseph Pusztay     \sum_{x_i \in C} w^i_t = 1/\tau \left( F_eq(C) - \sum_{x_i \in C} w_i \right)
198d0c080abSJoseph Pusztay */
199d0c080abSJoseph Pusztay 
200d0c080abSJoseph Pusztay /* This computes the 1D Maxwellian distribution for given mass n, velocity v, and temperature T */
201d0c080abSJoseph Pusztay static PetscReal ComputePDF(PetscReal m, PetscReal n, PetscReal T, PetscReal v[])
202d0c080abSJoseph Pusztay {
203d0c080abSJoseph Pusztay   return (n/PetscSqrtReal(2.0*PETSC_PI*T/m)) * PetscExpReal(-0.5*m*PetscSqr(v[0])/T);
204d0c080abSJoseph Pusztay }
205d0c080abSJoseph Pusztay 
206d0c080abSJoseph Pusztay /*
207d0c080abSJoseph Pusztay   erf z = \frac{2}{\sqrt\pi} \int^z_0 dt e^{-t^2} and erf \infty = 1
208d0c080abSJoseph Pusztay 
209d0c080abSJoseph Pusztay   We begin with our distribution
210d0c080abSJoseph Pusztay 
211d0c080abSJoseph Pusztay     \sqrt{\frac{m}{2 \pi T}} e^{-m v^2/2T}
212d0c080abSJoseph Pusztay 
213d0c080abSJoseph Pusztay   Let t = \sqrt{m/2T} v, z = \sqrt{m/2T} w, so that we now have
214d0c080abSJoseph Pusztay 
215d0c080abSJoseph Pusztay       \sqrt{\frac{m}{2 \pi T}} \int^w_0 dv e^{-m v^2/2T}
216d0c080abSJoseph Pusztay     = \sqrt{\frac{m}{2 \pi T}} \int^{\sqrt{m/2T} w}_0 \sqrt{2T/m} dt e^{-t^2}
217d0c080abSJoseph Pusztay     = 1/\sqrt{\pi} \int^{\sqrt{m/2T} w}_0 dt e^{-t^2}
218d0c080abSJoseph Pusztay     = 1/2 erf(\sqrt{m/2T} w)
219d0c080abSJoseph Pusztay */
220d0c080abSJoseph Pusztay static PetscReal ComputeCDF(PetscReal m, PetscReal n, PetscReal T, PetscReal va, PetscReal vb)
221d0c080abSJoseph Pusztay {
222d0c080abSJoseph Pusztay   PetscReal alpha = PetscSqrtReal(0.5*m/T);
223d0c080abSJoseph Pusztay   PetscReal za    = alpha*va;
224d0c080abSJoseph Pusztay   PetscReal zb    = alpha*vb;
225d0c080abSJoseph Pusztay   PetscReal sum   = 0.0;
226d0c080abSJoseph Pusztay 
227d0c080abSJoseph Pusztay   sum += zb >= 0 ? erf(zb) : -erf(-zb);
228d0c080abSJoseph Pusztay   sum -= za >= 0 ? erf(za) : -erf(-za);
229d0c080abSJoseph Pusztay   return 0.5 * n * sum;
230d0c080abSJoseph Pusztay }
231d0c080abSJoseph Pusztay 
232d0c080abSJoseph Pusztay static PetscErrorCode CheckDistribution(DM dm, PetscReal m, PetscReal n, PetscReal T, PetscReal v[])
233d0c080abSJoseph Pusztay {
234d0c080abSJoseph Pusztay   PetscSection   coordSection;
235d0c080abSJoseph Pusztay   Vec            coordsLocal;
236d0c080abSJoseph Pusztay   PetscReal     *xq, *wq;
237d0c080abSJoseph Pusztay   PetscReal      vmin, vmax, neq, veq, Teq;
238d0c080abSJoseph Pusztay   PetscInt       Nq = 100, q, cStart, cEnd, c;
239d0c080abSJoseph Pusztay   PetscErrorCode ierr;
240d0c080abSJoseph Pusztay 
241d0c080abSJoseph Pusztay   PetscFunctionBegin;
242d0c080abSJoseph Pusztay   ierr = DMGetBoundingBox(dm, &vmin, &vmax);CHKERRQ(ierr);
243d0c080abSJoseph Pusztay   /* Check analytic over entire line */
244d0c080abSJoseph Pusztay   neq  = ComputeCDF(m, n, T, vmin, vmax);
245d0c080abSJoseph Pusztay   if (PetscAbsReal(neq - n) > PETSC_SMALL) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", neq, n, neq-n);
246d0c080abSJoseph Pusztay   /* Check analytic over cells */
247d0c080abSJoseph Pusztay   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
248d0c080abSJoseph Pusztay   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
249d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
250d0c080abSJoseph Pusztay   neq  = 0.0;
251d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
252d0c080abSJoseph Pusztay     PetscScalar *vcoords = NULL;
253d0c080abSJoseph Pusztay 
254d0c080abSJoseph Pusztay     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords);CHKERRQ(ierr);
255d0c080abSJoseph Pusztay     neq += ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
256d0c080abSJoseph Pusztay     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords);CHKERRQ(ierr);
257d0c080abSJoseph Pusztay   }
258d0c080abSJoseph Pusztay   if (PetscAbsReal(neq - n) > PETSC_SMALL) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell Int f %g != %g mass (%g)", neq, n, neq-n);
259d0c080abSJoseph Pusztay   /* Check quadrature over entire line */
260d0c080abSJoseph Pusztay   ierr = PetscMalloc2(Nq, &xq, Nq, &wq);CHKERRQ(ierr);
261d0c080abSJoseph Pusztay   ierr = PetscDTGaussQuadrature(100, vmin, vmax, xq, wq);CHKERRQ(ierr);
262d0c080abSJoseph Pusztay   neq  = 0.0;
263d0c080abSJoseph Pusztay   for (q = 0; q < Nq; ++q) {
264d0c080abSJoseph Pusztay     neq += ComputePDF(m, n, T, &xq[q])*wq[q];
265d0c080abSJoseph Pusztay   }
266d0c080abSJoseph Pusztay   if (PetscAbsReal(neq - n) > PETSC_SMALL) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", neq, n, neq-n);
267d0c080abSJoseph Pusztay   /* Check omemnts with quadrature */
268d0c080abSJoseph Pusztay   veq  = 0.0;
269d0c080abSJoseph Pusztay   for (q = 0; q < Nq; ++q) {
270d0c080abSJoseph Pusztay     veq += xq[q]*ComputePDF(m, n, T, &xq[q])*wq[q];
271d0c080abSJoseph Pusztay   }
272d0c080abSJoseph Pusztay   veq /= neq;
273d0c080abSJoseph Pusztay   if (PetscAbsReal(veq - v[0]) > PETSC_SMALL) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v f %g != %g velocity (%g)", veq, v[0], veq-v[0]);
274d0c080abSJoseph Pusztay   Teq  = 0.0;
275d0c080abSJoseph Pusztay   for (q = 0; q < Nq; ++q) {
276d0c080abSJoseph Pusztay     Teq += PetscSqr(xq[q])*ComputePDF(m, n, T, &xq[q])*wq[q];
277d0c080abSJoseph Pusztay   }
278d0c080abSJoseph Pusztay   Teq = Teq * m/neq - PetscSqr(veq);
279d0c080abSJoseph Pusztay   if (PetscAbsReal(Teq - T) > PETSC_SMALL) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v^2 f %g != %g temperature (%g)", Teq, T, Teq-T);
280d0c080abSJoseph Pusztay   ierr = PetscFree2(xq, wq);CHKERRQ(ierr);
281d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
282d0c080abSJoseph Pusztay }
283d0c080abSJoseph Pusztay 
284d0c080abSJoseph Pusztay static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
285d0c080abSJoseph Pusztay {
286d0c080abSJoseph Pusztay   const PetscScalar *u;
287d0c080abSJoseph Pusztay   PetscSection       coordSection;
288d0c080abSJoseph Pusztay   Vec                coordsLocal;
289d0c080abSJoseph Pusztay   PetscScalar       *r, *coords;
290d0c080abSJoseph Pusztay   PetscReal          n = 0.0, v = 0.0, E = 0.0, T = 0.0, m = 1.0, cn = 0.0, cv = 0.0, cE = 0.0, pE = 0.0, eqE = 0.0;
291d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, Ncp, p, cStart, cEnd, c;
292d0c080abSJoseph Pusztay   DM                 dmSw, plex;
293d0c080abSJoseph Pusztay   PetscErrorCode     ierr;
294d0c080abSJoseph Pusztay 
295d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
296d0c080abSJoseph Pusztay   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
297d0c080abSJoseph Pusztay   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
298d0c080abSJoseph Pusztay   ierr = VecGetArray(R, &r);CHKERRQ(ierr);
299d0c080abSJoseph Pusztay   ierr = TSGetDM(ts, &dmSw);CHKERRQ(ierr);
300d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellDM(dmSw, &plex);CHKERRQ(ierr);
301d0c080abSJoseph Pusztay   ierr = DMGetDimension(dmSw, &dim);CHKERRQ(ierr);
302d0c080abSJoseph Pusztay   ierr = DMGetCoordinatesLocal(plex, &coordsLocal);CHKERRQ(ierr);
303d0c080abSJoseph Pusztay   ierr = DMGetCoordinateSection(plex, &coordSection);CHKERRQ(ierr);
304d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
305d0c080abSJoseph Pusztay   Np  /= dim;
306d0c080abSJoseph Pusztay   Ncp  = Np / (cEnd - cStart);
307d0c080abSJoseph Pusztay   /* Calculate moments of particle distribution, note that velocity is in the coordinate */
308d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
309d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
310d0c080abSJoseph Pusztay     PetscReal m1 = 0.0, m2 = 0.0;
311d0c080abSJoseph Pusztay 
312d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {m1 += PetscRealPart(coords[p*dim+d]); m2 += PetscSqr(PetscRealPart(coords[p*dim+d]));}
313d0c080abSJoseph Pusztay     n += u[p];
314d0c080abSJoseph Pusztay     v += u[p]*m1;
315d0c080abSJoseph Pusztay     E += u[p]*m2;
316d0c080abSJoseph Pusztay   }
317d0c080abSJoseph Pusztay   v /= n;
318d0c080abSJoseph Pusztay   T  = E*m/n - v*v;
319d0c080abSJoseph Pusztay   ierr = PetscInfo4(ts, "Time %.2f: mass %.4f velocity: %+.4f temperature: %.4f\n", t, n, v, T);CHKERRQ(ierr);
320d0c080abSJoseph Pusztay   ierr = CheckDistribution(plex, m, n, T, &v);CHKERRQ(ierr);
321d0c080abSJoseph Pusztay   /*
322d0c080abSJoseph Pusztay      Begin cellwise evaluation of the collision operator. Essentially, penalize the weights of the particles
323d0c080abSJoseph Pusztay      in that cell against the maxwellian for the number of particles expected to be in that cell
324d0c080abSJoseph Pusztay   */
325d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
326d0c080abSJoseph Pusztay     PetscScalar *vcoords = NULL;
327d0c080abSJoseph Pusztay     PetscReal    relaxation = 1.0, neq;
328d0c080abSJoseph Pusztay     PetscInt     sp      = c*Ncp, q;
329d0c080abSJoseph Pusztay 
330d0c080abSJoseph Pusztay     /* Calculate equilibrium occupation for this velocity cell */
331d0c080abSJoseph Pusztay     ierr = DMPlexVecGetClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords);CHKERRQ(ierr);
332d0c080abSJoseph Pusztay     neq  = ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
333d0c080abSJoseph Pusztay     ierr = DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords);CHKERRQ(ierr);
334d0c080abSJoseph Pusztay     for (q = 0; q < Ncp; ++q) r[sp+q] = (1.0/relaxation)*(neq - u[sp+q]);
335d0c080abSJoseph Pusztay   }
336d0c080abSJoseph Pusztay   /* Check update */
337d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
338d0c080abSJoseph Pusztay     PetscReal m1 = 0.0, m2 = 0.0;
339d0c080abSJoseph Pusztay     PetscScalar *vcoords = NULL;
340d0c080abSJoseph Pusztay 
341d0c080abSJoseph Pusztay     for (d = 0; d < dim; ++d) {m1 += PetscRealPart(coords[p*dim+d]); m2 += PetscSqr(PetscRealPart(coords[p*dim+d]));}
342d0c080abSJoseph Pusztay     cn += r[p];
343d0c080abSJoseph Pusztay     cv += r[p]*m1;
344d0c080abSJoseph Pusztay     cE += r[p]*m2;
345d0c080abSJoseph Pusztay     pE  += u[p]*m2;
346d0c080abSJoseph Pusztay     ierr = DMPlexVecGetClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords);CHKERRQ(ierr);
347d0c080abSJoseph Pusztay     eqE += ComputeCDF(m, n, T, vcoords[0], vcoords[1])*m2;
348d0c080abSJoseph Pusztay     ierr = DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords);CHKERRQ(ierr);
349d0c080abSJoseph Pusztay   }
350d0c080abSJoseph Pusztay   ierr = PetscInfo6(ts, "Time %.2f: mass update %.8f velocity update: %+.8f energy update: %.8f (%.8f, %.8f)\n", t, cn, cv, cE, pE, eqE);CHKERRQ(ierr);
351d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
352d0c080abSJoseph Pusztay   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
353d0c080abSJoseph Pusztay   ierr = VecRestoreArray(R, &r);CHKERRQ(ierr);
354d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
355d0c080abSJoseph Pusztay }
356d0c080abSJoseph Pusztay 
357d0c080abSJoseph Pusztay static PetscErrorCode HGMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
358d0c080abSJoseph Pusztay {
359d0c080abSJoseph Pusztay   AppCtx            *user  = (AppCtx *) ctx;
360d0c080abSJoseph Pusztay   const PetscScalar *u;
361d0c080abSJoseph Pusztay   DM                 sw, dm;
362d0c080abSJoseph Pusztay   PetscInt           dim, Np, p;
363d0c080abSJoseph Pusztay   PetscErrorCode     ierr;
364d0c080abSJoseph Pusztay 
365d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
366d0c080abSJoseph Pusztay   if (step < 0) PetscFunctionReturn(0);
367d0c080abSJoseph Pusztay   if (((user->ostep > 0) && (!(step % user->ostep)))) {
368d0c080abSJoseph Pusztay     PetscDrawAxis axis;
369d0c080abSJoseph Pusztay 
370d0c080abSJoseph Pusztay     ierr = TSGetDM(ts, &sw);CHKERRQ(ierr);
371d0c080abSJoseph Pusztay     ierr = DMSwarmGetCellDM(sw, &dm);CHKERRQ(ierr);
372d0c080abSJoseph Pusztay     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
373d0c080abSJoseph Pusztay     ierr = PetscDrawHGReset(user->drawhg);CHKERRQ(ierr);
374d0c080abSJoseph Pusztay     ierr = PetscDrawHGGetAxis(user->drawhg,&axis);CHKERRQ(ierr);
375d0c080abSJoseph Pusztay     ierr = PetscDrawAxisSetLabels(axis,"Particles","V","f(V)");CHKERRQ(ierr);
376d0c080abSJoseph Pusztay     ierr = PetscDrawAxisSetLimits(axis, -3, 3, 0, 100);CHKERRQ(ierr);
377d0c080abSJoseph Pusztay     ierr = PetscDrawHGSetLimits(user->drawhg,-3.0, 3.0, 0, 10);CHKERRQ(ierr);
378d0c080abSJoseph Pusztay     ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
379d0c080abSJoseph Pusztay     Np  /= dim;
380d0c080abSJoseph Pusztay     ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
381d0c080abSJoseph Pusztay     /* get points from solution vector */
382d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {ierr = PetscDrawHGAddValue(user->drawhg,u[p]);CHKERRQ(ierr);}
383d0c080abSJoseph Pusztay     ierr = PetscDrawHGDraw(user->drawhg);CHKERRQ(ierr);
384d0c080abSJoseph Pusztay     ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
385d0c080abSJoseph Pusztay   }
386d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
387d0c080abSJoseph Pusztay }
388d0c080abSJoseph Pusztay 
389d0c080abSJoseph Pusztay static PetscErrorCode SPMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
390d0c080abSJoseph Pusztay {
391d0c080abSJoseph Pusztay   AppCtx            *user  = (AppCtx *) ctx;
392d0c080abSJoseph Pusztay   const PetscScalar *u;
393d0c080abSJoseph Pusztay   PetscReal         *v, *coords;
394d0c080abSJoseph Pusztay   PetscInt           Np, p;
395d0c080abSJoseph Pusztay   DM                 dmSw;
396d0c080abSJoseph Pusztay   PetscErrorCode     ierr;
397d0c080abSJoseph Pusztay 
398d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
399d0c080abSJoseph Pusztay 
400d0c080abSJoseph Pusztay   if (step < 0) PetscFunctionReturn(0);
401d0c080abSJoseph Pusztay   if (((user->ostep > 0) && (!(step % user->ostep)))) {
402d0c080abSJoseph Pusztay     PetscDrawAxis axis;
403d0c080abSJoseph Pusztay 
404d0c080abSJoseph Pusztay     ierr = TSGetDM(ts, &dmSw);CHKERRQ(ierr);
405d0c080abSJoseph Pusztay     ierr = PetscDrawSPReset(user->drawsp);CHKERRQ(ierr);
406d0c080abSJoseph Pusztay     ierr = PetscDrawSPGetAxis(user->drawsp,&axis);CHKERRQ(ierr);
407d0c080abSJoseph Pusztay     ierr = PetscDrawAxisSetLabels(axis,"Particles","V","w");CHKERRQ(ierr);
408d0c080abSJoseph Pusztay     ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
409d0c080abSJoseph Pusztay     ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
410d0c080abSJoseph Pusztay     /* get points from solution vector */
411d0c080abSJoseph Pusztay     ierr = PetscMalloc1(Np, &v);CHKERRQ(ierr);
412d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
413d0c080abSJoseph Pusztay     ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
414d0c080abSJoseph Pusztay     for (p = 0; p < Np-1; ++p) {ierr = PetscDrawSPAddPoint(user->drawsp, &coords[p], &v[p]);CHKERRQ(ierr);}
415d0c080abSJoseph Pusztay     ierr = PetscDrawSPDraw(user->drawsp, PETSC_TRUE);CHKERRQ(ierr);
416d0c080abSJoseph Pusztay     ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
417d0c080abSJoseph Pusztay     ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
418d0c080abSJoseph Pusztay     ierr = PetscFree(v);CHKERRQ(ierr);
419d0c080abSJoseph Pusztay   }
420d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
421d0c080abSJoseph Pusztay }
422d0c080abSJoseph Pusztay 
423d0c080abSJoseph Pusztay static PetscErrorCode KSConv(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
424d0c080abSJoseph Pusztay {
425d0c080abSJoseph Pusztay   AppCtx            *user  = (AppCtx *) ctx;
426d0c080abSJoseph Pusztay   const PetscScalar *u;
427d0c080abSJoseph Pusztay   PetscScalar       sup;
428d0c080abSJoseph Pusztay   PetscReal         *v, *coords, T=0., vel=0., step_cast, w_sum;
429d0c080abSJoseph Pusztay   PetscInt           dim, Np, p, cStart, cEnd;
430d0c080abSJoseph Pusztay   DM                 sw, plex;
431d0c080abSJoseph Pusztay   PetscErrorCode     ierr;
432d0c080abSJoseph Pusztay 
433d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
434d0c080abSJoseph Pusztay   if (step < 0) PetscFunctionReturn(0);
435d0c080abSJoseph Pusztay   if (((user->ostep > 0) && (!(step % user->ostep)))) {
436d0c080abSJoseph Pusztay     PetscDrawAxis axis;
437d0c080abSJoseph Pusztay     ierr = PetscDrawSPGetAxis(user->drawks,&axis);CHKERRQ(ierr);
438d0c080abSJoseph Pusztay     ierr = PetscDrawAxisSetLabels(axis,"Particles","t","D_n");CHKERRQ(ierr);
439d0c080abSJoseph Pusztay     ierr = PetscDrawSPSetLimits(user->drawks,0.,100,0.,3.5);CHKERRQ(ierr);
440d0c080abSJoseph Pusztay     ierr = TSGetDM(ts, &sw);CHKERRQ(ierr);
441d0c080abSJoseph Pusztay     ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
442d0c080abSJoseph Pusztay     ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
443d0c080abSJoseph Pusztay     /* get points from solution vector */
444d0c080abSJoseph Pusztay     ierr = PetscMalloc1(Np, &v);CHKERRQ(ierr);
445d0c080abSJoseph Pusztay     ierr = DMSwarmGetCellDM(sw, &plex);CHKERRQ(ierr);
446d0c080abSJoseph Pusztay     ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr);
447d0c080abSJoseph Pusztay     ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
448d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
449d0c080abSJoseph Pusztay     ierr = DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
450d0c080abSJoseph Pusztay     w_sum = 0.;
451d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
452d0c080abSJoseph Pusztay       w_sum += u[p];
453d0c080abSJoseph Pusztay       T += u[p]*coords[p]*coords[p];
454d0c080abSJoseph Pusztay       vel += u[p]*coords[p];
455d0c080abSJoseph Pusztay     }
456d0c080abSJoseph Pusztay     vel /= w_sum;
457d0c080abSJoseph Pusztay     T = T/w_sum - vel*vel;
458d0c080abSJoseph Pusztay     sup = 0.0;
459d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
460d0c080abSJoseph Pusztay         PetscReal tmp = 0.;
461d0c080abSJoseph Pusztay         tmp = PetscAbs(u[p]-ComputePDF(1.0, w_sum, T, &coords[p*dim]));
462d0c080abSJoseph Pusztay         if (tmp > sup) sup = tmp;
463d0c080abSJoseph Pusztay     }
464d0c080abSJoseph Pusztay     step_cast = (PetscReal)step;
465d0c080abSJoseph Pusztay     ierr = PetscDrawSPAddPoint(user->drawks, &step_cast, &sup);CHKERRQ(ierr);
466d0c080abSJoseph Pusztay     ierr = PetscDrawSPDraw(user->drawks, PETSC_TRUE);CHKERRQ(ierr);
467d0c080abSJoseph Pusztay     ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
468d0c080abSJoseph Pusztay     ierr = DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
469d0c080abSJoseph Pusztay     ierr = PetscFree(v);CHKERRQ(ierr);
470d0c080abSJoseph Pusztay   }
471d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
472d0c080abSJoseph Pusztay }
473d0c080abSJoseph Pusztay 
474d0c080abSJoseph Pusztay static PetscErrorCode InitializeSolve(TS ts, Vec u)
475d0c080abSJoseph Pusztay {
476d0c080abSJoseph Pusztay   DM             dm;
477d0c080abSJoseph Pusztay   AppCtx        *user;
478d0c080abSJoseph Pusztay   PetscErrorCode ierr;
479d0c080abSJoseph Pusztay 
480d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
481d0c080abSJoseph Pusztay   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
4823ec1f749SStefano Zampini   ierr = DMGetApplicationContext(dm, &user);CHKERRQ(ierr);
483d0c080abSJoseph Pusztay   ierr = SetInitialCoordinates(dm);CHKERRQ(ierr);
484d0c080abSJoseph Pusztay   ierr = SetInitialConditions(dm, u);CHKERRQ(ierr);
485d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
486d0c080abSJoseph Pusztay }
487d0c080abSJoseph Pusztay   /*
488d0c080abSJoseph Pusztay      A single particle is generated for each velocity space cell using the dmswarmpicfield_coor and is used to evaluate collisions in that cell.
489d0c080abSJoseph Pusztay      0 weight ghost particles are initialized outside of a small velocity domain to ensure the tails of the amxwellian are resolved.
490d0c080abSJoseph Pusztay    */
491d0c080abSJoseph Pusztay int main(int argc,char **argv)
492d0c080abSJoseph Pusztay {
493d0c080abSJoseph Pusztay   TS             ts;     /* nonlinear solver */
494d0c080abSJoseph Pusztay   DM             dm, sw; /* Velocity space mesh and Particle Swarm */
495d0c080abSJoseph Pusztay   Vec            u, w;   /* swarm vector */
496d0c080abSJoseph Pusztay   MPI_Comm       comm;
497d0c080abSJoseph Pusztay   AppCtx         user;
498d0c080abSJoseph Pusztay   PetscErrorCode ierr;
499d0c080abSJoseph Pusztay 
500d0c080abSJoseph Pusztay   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
501d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
502d0c080abSJoseph Pusztay   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
503d0c080abSJoseph Pusztay   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
504d0c080abSJoseph Pusztay   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
505d0c080abSJoseph Pusztay   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
506d0c080abSJoseph Pusztay   ierr = TSCreate(comm, &ts);CHKERRQ(ierr);
507d0c080abSJoseph Pusztay   ierr = TSSetDM(ts, sw);CHKERRQ(ierr);
508d0c080abSJoseph Pusztay   ierr = TSSetMaxTime(ts, 10.0);CHKERRQ(ierr);
509d0c080abSJoseph Pusztay   ierr = TSSetTimeStep(ts, 0.01);CHKERRQ(ierr);
510d0c080abSJoseph Pusztay   ierr = TSSetMaxSteps(ts, 100000);CHKERRQ(ierr);
511d0c080abSJoseph Pusztay   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
512d0c080abSJoseph Pusztay   if (user.monitorhg) {
513d0c080abSJoseph Pusztay     ierr = PetscDrawCreate(comm, NULL, "monitor", 0,0,400,300, &user.draw);CHKERRQ(ierr);
514d0c080abSJoseph Pusztay     ierr = PetscDrawSetFromOptions(user.draw);CHKERRQ(ierr);
515d0c080abSJoseph Pusztay     ierr = PetscDrawHGCreate(user.draw, 20, &user.drawhg);CHKERRQ(ierr);
516d0c080abSJoseph Pusztay     ierr = PetscDrawHGSetColor(user.drawhg,3);CHKERRQ(ierr);
517d0c080abSJoseph Pusztay     ierr = TSMonitorSet(ts, HGMonitor, &user, NULL);CHKERRQ(ierr);
518d0c080abSJoseph Pusztay   }
519d0c080abSJoseph Pusztay   else if (user.monitorsp) {
520d0c080abSJoseph Pusztay     ierr = PetscDrawCreate(comm, NULL, "monitor", 0,0,400,300, &user.draw);CHKERRQ(ierr);
521d0c080abSJoseph Pusztay     ierr = PetscDrawSetFromOptions(user.draw);CHKERRQ(ierr);
522*2cf5aabcSBarry Smith     ierr = PetscDrawSPCreate(user.draw, 1, &user.drawsp);CHKERRQ(ierr);
523d0c080abSJoseph Pusztay     ierr = TSMonitorSet(ts, SPMonitor, &user, NULL);CHKERRQ(ierr);
524d0c080abSJoseph Pusztay   }
525d0c080abSJoseph Pusztay   else if (user.monitorks) {
526d0c080abSJoseph Pusztay     ierr = PetscDrawCreate(comm, NULL, "monitor", 0,0,400,300, &user.draw);CHKERRQ(ierr);
527d0c080abSJoseph Pusztay     ierr = PetscDrawSetFromOptions(user.draw);CHKERRQ(ierr);
528*2cf5aabcSBarry Smith     ierr = PetscDrawSPCreate(user.draw, 1, &user.drawks);CHKERRQ(ierr);
529d0c080abSJoseph Pusztay     ierr = TSMonitorSet(ts, KSConv, &user, NULL);CHKERRQ(ierr);
530d0c080abSJoseph Pusztay   }
531d0c080abSJoseph Pusztay   ierr = TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user);CHKERRQ(ierr);
532d0c080abSJoseph Pusztay   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
533d0c080abSJoseph Pusztay   ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr);
534d0c080abSJoseph Pusztay   ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w);CHKERRQ(ierr);
535d0c080abSJoseph Pusztay   ierr = VecDuplicate(w, &u);CHKERRQ(ierr);
536d0c080abSJoseph Pusztay   ierr = VecCopy(w, u);CHKERRQ(ierr);
537d0c080abSJoseph Pusztay   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w);CHKERRQ(ierr);
538d0c080abSJoseph Pusztay   ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr);
539d0c080abSJoseph Pusztay   ierr = TSSolve(ts, u);CHKERRQ(ierr);
540d0c080abSJoseph Pusztay   if (user.monitorhg) {
5411e1ea65dSPierre Jolivet     ierr = PetscDrawSave(user.draw);CHKERRQ(ierr);
542d0c080abSJoseph Pusztay     ierr = PetscDrawHGDestroy(&user.drawhg);CHKERRQ(ierr);
5431e1ea65dSPierre Jolivet     ierr = PetscDrawDestroy(&user.draw);CHKERRQ(ierr);
544d0c080abSJoseph Pusztay   }
545d0c080abSJoseph Pusztay   if (user.monitorsp) {
5461e1ea65dSPierre Jolivet     ierr = PetscDrawSave(user.draw);CHKERRQ(ierr);
547d0c080abSJoseph Pusztay     ierr = PetscDrawSPDestroy(&user.drawsp);CHKERRQ(ierr);
5481e1ea65dSPierre Jolivet     ierr = PetscDrawDestroy(&user.draw);CHKERRQ(ierr);
549d0c080abSJoseph Pusztay   }
550d0c080abSJoseph Pusztay   if (user.monitorks) {
5511e1ea65dSPierre Jolivet     ierr = PetscDrawSave(user.draw);CHKERRQ(ierr);
552d0c080abSJoseph Pusztay     ierr = PetscDrawSPDestroy(&user.drawks);CHKERRQ(ierr);
5531e1ea65dSPierre Jolivet     ierr = PetscDrawDestroy(&user.draw);CHKERRQ(ierr);
554d0c080abSJoseph Pusztay   }
555d0c080abSJoseph Pusztay   ierr = VecDestroy(&u);CHKERRQ(ierr);
556d0c080abSJoseph Pusztay   ierr = TSDestroy(&ts);CHKERRQ(ierr);
557d0c080abSJoseph Pusztay   ierr = DMDestroy(&sw);CHKERRQ(ierr);
558d0c080abSJoseph Pusztay   ierr = DMDestroy(&dm);CHKERRQ(ierr);
559d0c080abSJoseph Pusztay   ierr = PetscFinalize();
560d0c080abSJoseph Pusztay   return ierr;
561d0c080abSJoseph Pusztay }
562d0c080abSJoseph Pusztay 
563d0c080abSJoseph Pusztay /*TEST
564d0c080abSJoseph Pusztay    build:
56530602db0SMatthew G. Knepley      requires: double !complex
566d0c080abSJoseph Pusztay    test:
567d0c080abSJoseph Pusztay      suffix: 1
56830602db0SMatthew G. Knepley      args: -particles_per_cell 1 -output_step 10 -ts_type euler -dm_plex_dim 1 -dm_plex_box_faces 200 -dm_plex_box_lower -10 -dm_plex_box_upper 10 -dm_view -monitorsp
569d0c080abSJoseph Pusztay    test:
570d0c080abSJoseph Pusztay      suffix: 2
57130602db0SMatthew G. Knepley      args: -particles_per_cell 1 -output_step 50 -ts_type euler -dm_plex_dim 1 -dm_plex_box_faces 200 -dm_plex_box_lower -10 -dm_plex_box_upper 10 -dm_view -monitorks
572d0c080abSJoseph Pusztay TEST*/
573