xref: /petsc/src/ts/tests/ex27.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1d0c080abSJoseph Pusztay static char help[] = "Particle Basis Landau Example using nonlinear solve + Implicit Midpoint-like time stepping.";
2d0c080abSJoseph Pusztay 
3d0c080abSJoseph Pusztay /* TODO
4d0c080abSJoseph Pusztay 
5d0c080abSJoseph Pusztay 1) SNES is sensitive to epsilon (but not to h). Should we do continuation in it?
6d0c080abSJoseph Pusztay 
7d0c080abSJoseph Pusztay 2) Put this timestepper in library, maybe by changing DG
8d0c080abSJoseph Pusztay 
9d0c080abSJoseph Pusztay 3) Add monitor to visualize distributions
10d0c080abSJoseph Pusztay 
11d0c080abSJoseph Pusztay */
12d0c080abSJoseph Pusztay 
13d0c080abSJoseph Pusztay /* References
14d0c080abSJoseph Pusztay   [1] https://arxiv.org/abs/1910.03080v2
15d0c080abSJoseph Pusztay */
16d0c080abSJoseph Pusztay 
17d0c080abSJoseph Pusztay #include <petscdmplex.h>
18d0c080abSJoseph Pusztay #include <petscdmswarm.h>
19d0c080abSJoseph Pusztay #include <petscts.h>
20d0c080abSJoseph Pusztay #include <petscviewer.h>
21d0c080abSJoseph Pusztay #include <petscmath.h>
22d0c080abSJoseph Pusztay 
23d0c080abSJoseph Pusztay typedef struct {
24d0c080abSJoseph Pusztay   /* Velocity space grid and functions */
25d0c080abSJoseph Pusztay   PetscInt  N;         /* The number of partices per spatial cell */
26d0c080abSJoseph Pusztay   PetscReal L;         /* Velocity space is [-L, L]^d */
27d0c080abSJoseph Pusztay   PetscReal h;         /* Spacing for grid 2L / N^{1/d} */
28d0c080abSJoseph Pusztay   PetscReal epsilon;   /* gaussian regularization parameter */
29d0c080abSJoseph Pusztay   PetscReal momentTol; /* Tolerance for checking moment conservation */
30d0c080abSJoseph Pusztay } AppCtx;
31d0c080abSJoseph Pusztay 
32*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33*d71ae5a4SJacob Faibussowitsch {
34d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
35d0c080abSJoseph Pusztay   options->N         = 1;
36d0c080abSJoseph Pusztay   options->momentTol = 100.0 * PETSC_MACHINE_EPSILON;
37d0c080abSJoseph Pusztay   options->L         = 1.0;
38d0c080abSJoseph Pusztay   options->h         = -1.0;
39d0c080abSJoseph Pusztay   options->epsilon   = -1.0;
40d0c080abSJoseph Pusztay 
41d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-N", "Number of particles per spatial cell", "ex27.c", options->N, &options->N, NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-L", "Velocity-space extent", "ex27.c", options->L, &options->L, NULL));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-h", "Velocity-space resolution", "ex27.c", options->h, &options->h, NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-epsilon", "Mollifier regularization parameter", "ex27.c", options->epsilon, &options->epsilon, NULL));
46d0609cedSBarry Smith   PetscOptionsEnd();
47d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
48d0c080abSJoseph Pusztay }
49d0c080abSJoseph Pusztay 
50*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
51*d71ae5a4SJacob Faibussowitsch {
52d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
549566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
559566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
569566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
57d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
58d0c080abSJoseph Pusztay }
59d0c080abSJoseph Pusztay 
60*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialCoordinates(DM sw)
61*d71ae5a4SJacob Faibussowitsch {
62d0c080abSJoseph Pusztay   AppCtx        *user;
63d0c080abSJoseph Pusztay   PetscRandom    rnd, rndv;
64d0c080abSJoseph Pusztay   DM             dm;
65d0c080abSJoseph Pusztay   DMPolytopeType ct;
66d0c080abSJoseph Pusztay   PetscBool      simplex;
67d0c080abSJoseph Pusztay   PetscReal     *centroid, *coords, *velocity, *xi0, *v0, *J, *invJ, detJ, *vals;
68d0c080abSJoseph Pusztay   PetscInt       dim, d, cStart, cEnd, c, Np, p;
69d0c080abSJoseph Pusztay 
70d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
71d0c080abSJoseph Pusztay   /* Randomization for coordinates */
729566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
739566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
749566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
759566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rndv));
769566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rndv, -1., 1.));
779566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rndv));
789566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(sw, &user));
79d0c080abSJoseph Pusztay   Np = user->N;
809566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
819566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dm));
829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
84d0c080abSJoseph Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
86d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
879566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
889566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity));
899566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals));
90d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
91d0c080abSJoseph Pusztay     if (Np == 1) {
929566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
93ad540459SPierre Jolivet       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
94d0c080abSJoseph Pusztay       vals[c] = 1.0;
95d0c080abSJoseph Pusztay     } else {
969566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
97d0c080abSJoseph Pusztay       for (p = 0; p < Np; ++p) {
98d0c080abSJoseph Pusztay         const PetscInt n   = c * Np + p;
99d0c080abSJoseph Pusztay         PetscReal      sum = 0.0, refcoords[3];
100d0c080abSJoseph Pusztay 
101d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) {
1029566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
103d0c080abSJoseph Pusztay           sum += refcoords[d];
104d0c080abSJoseph Pusztay         }
1059371c9d4SSatish Balay         if (simplex && sum > 0.0)
1069371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
107d0c080abSJoseph Pusztay         vals[n] = 1.0;
1089566063dSJacob Faibussowitsch         PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim]));
109d0c080abSJoseph Pusztay       }
110d0c080abSJoseph Pusztay     }
111d0c080abSJoseph Pusztay   }
112d0c080abSJoseph Pusztay   /* Random velocity IC */
113d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
114d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
115d0c080abSJoseph Pusztay       for (d = 0; d < dim; ++d) {
116d0c080abSJoseph Pusztay         PetscReal v_val;
117d0c080abSJoseph Pusztay 
1189566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rndv, &v_val));
119d0c080abSJoseph Pusztay         velocity[p * dim + d] = v_val;
120d0c080abSJoseph Pusztay       }
121d0c080abSJoseph Pusztay     }
122d0c080abSJoseph Pusztay   }
1239566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1249566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity));
1259566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals));
1269566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
1279566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1289566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rndv));
129d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
130d0c080abSJoseph Pusztay }
131d0c080abSJoseph Pusztay 
132d0c080abSJoseph Pusztay /* Get velocities from swarm and place in solution vector */
133*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
134*d71ae5a4SJacob Faibussowitsch {
135d0c080abSJoseph Pusztay   DM           dm;
136d0c080abSJoseph Pusztay   AppCtx      *user;
137d0c080abSJoseph Pusztay   PetscReal   *velocity;
138d0c080abSJoseph Pusztay   PetscScalar *initialConditions;
139d0c080abSJoseph Pusztay   PetscInt     dim, d, cStart, cEnd, c, Np, p, n;
140d0c080abSJoseph Pusztay 
141d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
1429566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &n));
1439566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmSw, &user));
144d0c080abSJoseph Pusztay   Np = user->N;
1459566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dmSw, &dm));
1469566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1489566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **)&velocity));
1499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &initialConditions));
150d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
151d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
152d0c080abSJoseph Pusztay       const PetscInt n = c * Np + p;
153ad540459SPierre Jolivet       for (d = 0; d < dim; d++) initialConditions[n * dim + d] = velocity[n * dim + d];
154d0c080abSJoseph Pusztay     }
155d0c080abSJoseph Pusztay   }
1569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &initialConditions));
1579566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **)&velocity));
158d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
159d0c080abSJoseph Pusztay }
160d0c080abSJoseph Pusztay 
161*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
162*d71ae5a4SJacob Faibussowitsch {
163d0c080abSJoseph Pusztay   PetscInt *cellid;
164d0c080abSJoseph Pusztay   PetscInt  dim, cStart, cEnd, c, Np = user->N, p;
165d0c080abSJoseph Pusztay   PetscBool view = PETSC_FALSE;
166d0c080abSJoseph Pusztay 
167d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
1689566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1699566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1709566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1719566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
172d0c080abSJoseph Pusztay   /* h = 2L/n and N = n^d */
173d0c080abSJoseph Pusztay   if (user->h < 0.) user->h = 2. * user->L / PetscPowReal(user->N, 1. / dim);
174d0c080abSJoseph Pusztay   /* From Section 4 in [1], \epsilon = 0.64 h^.98 */
175d0c080abSJoseph Pusztay   if (user->epsilon < 0.) user->epsilon = 0.64 * pow(user->h, 1.98);
1769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL));
17748a46eb9SPierre Jolivet   if (view) PetscCall(PetscPrintf(PETSC_COMM_SELF, "N: %" PetscInt_FMT " L: %g h: %g eps: %g\n", user->N, (double)user->L, (double)user->h, (double)user->epsilon));
1789566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1799566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
1809566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
1819566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1829566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
1839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1849566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
1859566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
1869566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
187d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
188d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
189d0c080abSJoseph Pusztay       const PetscInt n = c * Np + p;
190d0c080abSJoseph Pusztay       cellid[n]        = c;
191d0c080abSJoseph Pusztay     }
192d0c080abSJoseph Pusztay   }
1939566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
1949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
1959566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
196d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
197d0c080abSJoseph Pusztay }
198d0c080abSJoseph Pusztay 
199d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
200*d71ae5a4SJacob Faibussowitsch static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
201*d71ae5a4SJacob Faibussowitsch {
202d0c080abSJoseph Pusztay   PetscInt d;
203d0c080abSJoseph Pusztay 
204d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d];
205d0c080abSJoseph Pusztay }
206d0c080abSJoseph Pusztay 
207d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
208*d71ae5a4SJacob Faibussowitsch static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
209*d71ae5a4SJacob Faibussowitsch {
210d0c080abSJoseph Pusztay   PetscReal sum = 0.0;
211d0c080abSJoseph Pusztay   PetscInt  d;
212d0c080abSJoseph Pusztay 
213d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d];
214d0c080abSJoseph Pusztay   return sum;
215d0c080abSJoseph Pusztay }
216d0c080abSJoseph Pusztay 
217d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
218*d71ae5a4SJacob Faibussowitsch static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
219*d71ae5a4SJacob Faibussowitsch {
220d0c080abSJoseph Pusztay   PetscScalar z[2];
2219371c9d4SSatish Balay   z[0] = x[0];
2229371c9d4SSatish Balay   z[1] = x[ldx];
223d0c080abSJoseph Pusztay   y[0] += A[0] * z[0] + A[1] * z[1];
224d0c080abSJoseph Pusztay   y[ldx] += A[2] * z[0] + A[3] * z[1];
225d0c080abSJoseph Pusztay   (void)PetscLogFlops(6.0);
226d0c080abSJoseph Pusztay }
227d0c080abSJoseph Pusztay 
228d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */
229*d71ae5a4SJacob Faibussowitsch static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
230*d71ae5a4SJacob Faibussowitsch {
231d0c080abSJoseph Pusztay   PetscScalar z[3];
2329371c9d4SSatish Balay   z[0] = x[0];
2339371c9d4SSatish Balay   z[1] = x[ldx];
2349371c9d4SSatish Balay   z[2] = x[ldx * 2];
235d0c080abSJoseph Pusztay   y[0] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
236d0c080abSJoseph Pusztay   y[ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
237d0c080abSJoseph Pusztay   y[ldx * 2] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
238d0c080abSJoseph Pusztay   (void)PetscLogFlops(15.0);
239d0c080abSJoseph Pusztay }
240d0c080abSJoseph Pusztay 
241d0c080abSJoseph Pusztay /*
242d0c080abSJoseph Pusztay   Gaussian - The Gaussian function G(x)
243d0c080abSJoseph Pusztay 
244d0c080abSJoseph Pusztay   Input Parameters:
245d0c080abSJoseph Pusztay +  dim   - The number of dimensions, or size of x
246d0c080abSJoseph Pusztay .  mu    - The mean, or center
247d0c080abSJoseph Pusztay .  sigma - The standard deviation, or width
248d0c080abSJoseph Pusztay -  x     - The evaluation point of the function
249d0c080abSJoseph Pusztay 
250d0c080abSJoseph Pusztay   Output Parameter:
251d0c080abSJoseph Pusztay . ret - The value G(x)
252d0c080abSJoseph Pusztay */
253*d71ae5a4SJacob Faibussowitsch static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[])
254*d71ae5a4SJacob Faibussowitsch {
255d0c080abSJoseph Pusztay   PetscReal arg = 0.0;
256d0c080abSJoseph Pusztay   PetscInt  d;
257d0c080abSJoseph Pusztay 
258d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]);
259d0c080abSJoseph Pusztay   return PetscPowReal(2.0 * PETSC_PI * sigma, -dim / 2.0) * PetscExpReal(-arg / (2.0 * sigma));
260d0c080abSJoseph Pusztay }
261d0c080abSJoseph Pusztay 
262d0c080abSJoseph Pusztay /*
263d0c080abSJoseph Pusztay   ComputeGradS - Compute grad_v dS_eps/df
264d0c080abSJoseph Pusztay 
265d0c080abSJoseph Pusztay   Input Parameters:
266d0c080abSJoseph Pusztay + dim      - The dimension
267d0c080abSJoseph Pusztay . Np       - The number of particles
268d0c080abSJoseph Pusztay . vp       - The velocity v_p of the particle at which we evaluate
269d0c080abSJoseph Pusztay . velocity - The velocity field for all particles
270d0c080abSJoseph Pusztay . epsilon  - The regularization strength
271d0c080abSJoseph Pusztay 
272d0c080abSJoseph Pusztay   Output Parameter:
273d0c080abSJoseph Pusztay . integral - The output grad_v dS_eps/df (v_p)
274d0c080abSJoseph Pusztay 
275d0c080abSJoseph Pusztay   Note:
276d0c080abSJoseph Pusztay   This comes from (3.6) in [1], and we are computing
277d0c080abSJoseph Pusztay $   \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q)
278d0c080abSJoseph Pusztay   which is discretized by using a one-point quadrature in each box l at its center v^c_l
279d0c080abSJoseph Pusztay $   \sum_l h^d \nabla\psi_\epsilon(v_p - v^c_l) \log\left( \sum_q w_q \psi_\epsilon(v^c_l - v_q) \right)
280d0c080abSJoseph Pusztay   where h^d is the volume of each box.
281d0c080abSJoseph Pusztay */
282*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx)
283*d71ae5a4SJacob Faibussowitsch {
284d0c080abSJoseph Pusztay   PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5 * h - L;
285d0c080abSJoseph Pusztay   PetscInt  nx = roundf(2. * L / h);
286d0c080abSJoseph Pusztay   PetscInt  ny = dim > 1 ? nx : 1;
287d0c080abSJoseph Pusztay   PetscInt  nz = dim > 2 ? nx : 1;
288d0c080abSJoseph Pusztay   PetscInt  i, j, k, d, q, dbg = 0;
289d0c080abSJoseph Pusztay 
290d0c080abSJoseph Pusztay   PetscFunctionBeginHot;
291d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) integral[d] = 0.0;
292d0c080abSJoseph Pusztay   for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) {
293d0c080abSJoseph Pusztay     for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) {
294d0c080abSJoseph Pusztay       for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) {
295d0c080abSJoseph Pusztay         PetscReal sum = 0.0;
296d0c080abSJoseph Pusztay 
29763a3b9bcSJacob Faibussowitsch         if (dbg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(%" PetscInt_FMT " %" PetscInt_FMT ") vc_l: %g %g\n", i, j, (double)vc_l[0], (double)vc_l[1]));
298d0c080abSJoseph Pusztay         /* \log \sum_k \psi(v - v_k)  */
299d0c080abSJoseph Pusztay         for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q * dim], epsilon, vc_l);
300d0c080abSJoseph Pusztay         sum = PetscLogReal(sum);
301d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) integral[d] += (-1. / (epsilon)) * PetscAbsReal(vp[d] - vc_l[d]) * (Gaussian(dim, vp, epsilon, vc_l)) * sum;
302d0c080abSJoseph Pusztay       }
303d0c080abSJoseph Pusztay     }
304d0c080abSJoseph Pusztay   }
305d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
306d0c080abSJoseph Pusztay }
307d0c080abSJoseph Pusztay 
308d0c080abSJoseph Pusztay /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */
309*d71ae5a4SJacob Faibussowitsch static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[])
310*d71ae5a4SJacob Faibussowitsch {
311d0c080abSJoseph Pusztay   PetscReal xi[3], xi2, xi3, mag;
312d0c080abSJoseph Pusztay   PetscInt  d, e;
313d0c080abSJoseph Pusztay 
314d0c080abSJoseph Pusztay   PetscFunctionBeginHot;
315d0c080abSJoseph Pusztay   DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi);
316d0c080abSJoseph Pusztay   xi2 = DMPlex_DotD_Internal(dim, xi, xi);
317d0c080abSJoseph Pusztay   mag = PetscSqrtReal(xi2);
318d0c080abSJoseph Pusztay   xi3 = xi2 * mag;
319d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) {
320ad540459SPierre Jolivet     for (e = 0; e < dim; ++e) Q[d * dim + e] = -xi[d] * xi[e] / xi3;
321d0c080abSJoseph Pusztay     Q[d * dim + d] += 1. / mag;
322d0c080abSJoseph Pusztay   }
323d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
324d0c080abSJoseph Pusztay }
325d0c080abSJoseph Pusztay 
326*d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
327*d71ae5a4SJacob Faibussowitsch {
328d0c080abSJoseph Pusztay   AppCtx            *user = (AppCtx *)ctx;
329d0c080abSJoseph Pusztay   PetscInt           dbg  = 0;
330d0c080abSJoseph Pusztay   DM                 sw;  /* Particles */
331d0c080abSJoseph Pusztay   Vec                sol; /* Solution vector at current time */
332d0c080abSJoseph Pusztay   const PetscScalar *u;   /* input solution vector */
333d0c080abSJoseph Pusztay   PetscScalar       *r;
334d0c080abSJoseph Pusztay   PetscReal         *velocity;
335d0c080abSJoseph Pusztay   PetscInt           dim, Np, p, q;
336d0c080abSJoseph Pusztay 
337d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
3389566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(R));
3399566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
3409566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
3419566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
3429566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &sol));
3439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(sol, &velocity));
3449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
3459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
346d0c080abSJoseph Pusztay   Np /= dim;
3479566063dSJacob Faibussowitsch   if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Part  ppr     x        y\n"));
348d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
349d0c080abSJoseph Pusztay     PetscReal gradS_p[3] = {0., 0., 0.};
350d0c080abSJoseph Pusztay 
3519566063dSJacob Faibussowitsch     PetscCall(ComputeGradS(dim, Np, &velocity[p * dim], velocity, gradS_p, user));
352d0c080abSJoseph Pusztay     for (q = 0; q < Np; ++q) {
353d0c080abSJoseph Pusztay       PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9];
354d0c080abSJoseph Pusztay 
355d0c080abSJoseph Pusztay       if (q == p) continue;
3569566063dSJacob Faibussowitsch       PetscCall(ComputeGradS(dim, Np, &velocity[q * dim], velocity, gradS_q, user));
357d0c080abSJoseph Pusztay       DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS);
3589566063dSJacob Faibussowitsch       PetscCall(QCompute(dim, &u[p * dim], &u[q * dim], Q));
359d0c080abSJoseph Pusztay       switch (dim) {
360*d71ae5a4SJacob Faibussowitsch       case 2:
361*d71ae5a4SJacob Faibussowitsch         DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p * dim]);
362*d71ae5a4SJacob Faibussowitsch         break;
363*d71ae5a4SJacob Faibussowitsch       case 3:
364*d71ae5a4SJacob Faibussowitsch         DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p * dim]);
365*d71ae5a4SJacob Faibussowitsch         break;
366*d71ae5a4SJacob Faibussowitsch       default:
367*d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %" PetscInt_FMT, dim);
368d0c080abSJoseph Pusztay       }
369d0c080abSJoseph Pusztay     }
37063a3b9bcSJacob Faibussowitsch     if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final %4" PetscInt_FMT " %10.8lf %10.8lf\n", p, r[p * dim + 0], r[p * dim + 1]));
371d0c080abSJoseph Pusztay   }
3729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
3739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
3749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(sol, &velocity));
3759566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(R, NULL, "-residual_view"));
376d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
377d0c080abSJoseph Pusztay }
378d0c080abSJoseph Pusztay 
379d0c080abSJoseph Pusztay /*
380d0c080abSJoseph Pusztay  TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform
381d0c080abSJoseph Pusztay  the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid
382d0c080abSJoseph Pusztay  to migrate between.
383d0c080abSJoseph Pusztay */
384*d71ae5a4SJacob Faibussowitsch static PetscErrorCode UpdateSwarm(TS ts)
385*d71ae5a4SJacob Faibussowitsch {
386d0c080abSJoseph Pusztay   PetscInt           idx, n;
387d0c080abSJoseph Pusztay   const PetscScalar *u;
388d0c080abSJoseph Pusztay   PetscScalar       *velocity;
389d0c080abSJoseph Pusztay   DM                 sw;
390d0c080abSJoseph Pusztay   Vec                sol;
391d0c080abSJoseph Pusztay 
392d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
3939566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
3949566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity));
3959566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &sol));
3969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(sol, &u));
3979566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(sol, &n));
398d0c080abSJoseph Pusztay   for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx];
3999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(sol, &u));
4009566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity));
401d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
402d0c080abSJoseph Pusztay }
403d0c080abSJoseph Pusztay 
404*d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitializeSolve(TS ts, Vec u)
405*d71ae5a4SJacob Faibussowitsch {
406d0c080abSJoseph Pusztay   DM      dm;
407d0c080abSJoseph Pusztay   AppCtx *user;
408d0c080abSJoseph Pusztay 
409d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
4109566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
4119566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &user));
4129566063dSJacob Faibussowitsch   PetscCall(SetInitialCoordinates(dm));
4139566063dSJacob Faibussowitsch   PetscCall(SetInitialConditions(dm, u));
414d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
415d0c080abSJoseph Pusztay }
416d0c080abSJoseph Pusztay 
417*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
418*d71ae5a4SJacob Faibussowitsch {
419d0c080abSJoseph Pusztay   TS       ts;     /* nonlinear solver */
420d0c080abSJoseph Pusztay   DM       dm, sw; /* Velocity space mesh and Particle Swarm */
421d0c080abSJoseph Pusztay   Vec      u, v;   /* problem vector */
422d0c080abSJoseph Pusztay   MPI_Comm comm;
423d0c080abSJoseph Pusztay   AppCtx   user;
424d0c080abSJoseph Pusztay 
425327415f7SBarry Smith   PetscFunctionBeginUser;
4269566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
427d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
4289566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
429d0c080abSJoseph Pusztay   /* Initialize objects and set initial conditions */
4309566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
4319566063dSJacob Faibussowitsch   PetscCall(CreateParticles(dm, &sw, &user));
4329566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(sw, &user));
4339566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(sw, "velocity"));
4349566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
4359566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, sw));
4369566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10.0));
4379566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.1));
4389566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1));
4399566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
4409566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user));
4419566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
4429566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
4439566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v));
4449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v, &u));
4459566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v));
4469566063dSJacob Faibussowitsch   PetscCall(TSComputeInitialCondition(ts, u));
4479566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, UpdateSwarm));
4489566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
4499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
4509566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
4519566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
4529566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
454b122ec5aSJacob Faibussowitsch   return 0;
455d0c080abSJoseph Pusztay }
456d0c080abSJoseph Pusztay 
457d0c080abSJoseph Pusztay /*TEST
458d0c080abSJoseph Pusztay    build:
459d0c080abSJoseph Pusztay      requires: triangle !single !complex
460d0c080abSJoseph Pusztay    test:
461d0c080abSJoseph Pusztay      suffix: midpoint
46230602db0SMatthew G. Knepley      args: -N 3 -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1 -dm_view \
463d0c080abSJoseph Pusztay            -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_monitor_frequency 1 -snes_fd
464d0c080abSJoseph Pusztay TEST*/
465