xref: /petsc/src/ts/tests/ex27.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
ProcessOptions(MPI_Comm comm,AppCtx * options)32d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33d71ae5a4SJacob 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();
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48d0c080abSJoseph Pusztay }
49d0c080abSJoseph Pusztay 
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * user)50d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
51d71ae5a4SJacob 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"));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58d0c080abSJoseph Pusztay }
59d0c080abSJoseph Pusztay 
SetInitialCoordinates(DM sw)60d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialCoordinates(DM sw)
61d71ae5a4SJacob 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));
829f4ada15SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
849566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
85d0c080abSJoseph Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
87d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
889566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
899566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity));
909566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals));
91d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
92d0c080abSJoseph Pusztay     if (Np == 1) {
939566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
94ad540459SPierre Jolivet       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
95d0c080abSJoseph Pusztay       vals[c] = 1.0;
96d0c080abSJoseph Pusztay     } else {
979566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
98d0c080abSJoseph Pusztay       for (p = 0; p < Np; ++p) {
99d0c080abSJoseph Pusztay         const PetscInt n   = c * Np + p;
100d0c080abSJoseph Pusztay         PetscReal      sum = 0.0, refcoords[3];
101d0c080abSJoseph Pusztay 
102d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) {
1039566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
104d0c080abSJoseph Pusztay           sum += refcoords[d];
105d0c080abSJoseph Pusztay         }
1069371c9d4SSatish Balay         if (simplex && sum > 0.0)
1079371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
108d0c080abSJoseph Pusztay         vals[n] = 1.0;
1099566063dSJacob Faibussowitsch         PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim]));
110d0c080abSJoseph Pusztay       }
111d0c080abSJoseph Pusztay     }
112d0c080abSJoseph Pusztay   }
113d0c080abSJoseph Pusztay   /* Random velocity IC */
114d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
115d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
116d0c080abSJoseph Pusztay       for (d = 0; d < dim; ++d) {
117d0c080abSJoseph Pusztay         PetscReal v_val;
118d0c080abSJoseph Pusztay 
1199566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rndv, &v_val));
120d0c080abSJoseph Pusztay         velocity[p * dim + d] = v_val;
121d0c080abSJoseph Pusztay       }
122d0c080abSJoseph Pusztay     }
123d0c080abSJoseph Pusztay   }
1249566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1259566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity));
1269566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals));
1279566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
1289566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1299566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rndv));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131d0c080abSJoseph Pusztay }
132d0c080abSJoseph Pusztay 
133d0c080abSJoseph Pusztay /* Get velocities from swarm and place in solution vector */
SetInitialConditions(DM dmSw,Vec u)134d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
135d71ae5a4SJacob Faibussowitsch {
136d0c080abSJoseph Pusztay   DM           dm;
137d0c080abSJoseph Pusztay   AppCtx      *user;
138d0c080abSJoseph Pusztay   PetscReal   *velocity;
139d0c080abSJoseph Pusztay   PetscScalar *initialConditions;
140d0c080abSJoseph Pusztay   PetscInt     dim, d, cStart, cEnd, c, Np, p, n;
141d0c080abSJoseph Pusztay 
142d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
1439566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &n));
1449566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmSw, &user));
145d0c080abSJoseph Pusztay   Np = user->N;
1469566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dmSw, &dm));
1479566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1499566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **)&velocity));
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &initialConditions));
151d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
152d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
153d0c080abSJoseph Pusztay       const PetscInt n = c * Np + p;
154ad540459SPierre Jolivet       for (d = 0; d < dim; d++) initialConditions[n * dim + d] = velocity[n * dim + d];
155d0c080abSJoseph Pusztay     }
156d0c080abSJoseph Pusztay   }
1579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &initialConditions));
1589566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **)&velocity));
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160d0c080abSJoseph Pusztay }
161d0c080abSJoseph Pusztay 
CreateParticles(DM dm,DM * sw,AppCtx * user)162d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
163d71ae5a4SJacob Faibussowitsch {
16419307e5cSMatthew G. Knepley   DMSwarmCellDM celldm;
16519307e5cSMatthew G. Knepley   PetscInt     *swarm_cellid;
166d0c080abSJoseph Pusztay   PetscInt      dim, cStart, cEnd, c, Np = user->N, p;
167d0c080abSJoseph Pusztay   PetscBool     view = PETSC_FALSE;
16819307e5cSMatthew G. Knepley   const char   *cellid;
169d0c080abSJoseph Pusztay 
170d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
1719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1729566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1739566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1749566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
175d0c080abSJoseph Pusztay   /* h = 2L/n and N = n^d */
176d0c080abSJoseph Pusztay   if (user->h < 0.) user->h = 2. * user->L / PetscPowReal(user->N, 1. / dim);
177d0c080abSJoseph Pusztay   /* From Section 4 in [1], \epsilon = 0.64 h^.98 */
178d0c080abSJoseph Pusztay   if (user->epsilon < 0.) user->epsilon = 0.64 * pow(user->h, 1.98);
1799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL));
18048a46eb9SPierre 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));
1819566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1829566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
1839566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
1849566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1859566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
18619307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(*sw, &celldm));
18719307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
1889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1899566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
1909566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
19119307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
192d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
193d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
194d0c080abSJoseph Pusztay       const PetscInt n = c * Np + p;
19519307e5cSMatthew G. Knepley       swarm_cellid[n]  = c;
196d0c080abSJoseph Pusztay     }
197d0c080abSJoseph Pusztay   }
19819307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
1999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
2009566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202d0c080abSJoseph Pusztay }
203d0c080abSJoseph Pusztay 
204d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
DMPlex_WaxpyD_Internal(PetscInt dim,PetscReal a,const PetscReal * x,const PetscReal * y,PetscReal * w)205d71ae5a4SJacob Faibussowitsch static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
206d71ae5a4SJacob Faibussowitsch {
207d0c080abSJoseph Pusztay   PetscInt d;
208d0c080abSJoseph Pusztay 
209d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d];
210d0c080abSJoseph Pusztay }
211d0c080abSJoseph Pusztay 
212d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
DMPlex_DotD_Internal(PetscInt dim,const PetscScalar * x,const PetscReal * y)213d71ae5a4SJacob Faibussowitsch static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
214d71ae5a4SJacob Faibussowitsch {
215d0c080abSJoseph Pusztay   PetscReal sum = 0.0;
216d0c080abSJoseph Pusztay   PetscInt  d;
217d0c080abSJoseph Pusztay 
218d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d];
219d0c080abSJoseph Pusztay   return sum;
220d0c080abSJoseph Pusztay }
221d0c080abSJoseph Pusztay 
222d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
DMPlex_MultAdd2DReal_Internal(const PetscReal A[],PetscInt ldx,const PetscScalar x[],PetscScalar y[])223d71ae5a4SJacob Faibussowitsch static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
224d71ae5a4SJacob Faibussowitsch {
225d0c080abSJoseph Pusztay   PetscScalar z[2];
2269371c9d4SSatish Balay   z[0] = x[0];
2279371c9d4SSatish Balay   z[1] = x[ldx];
228d0c080abSJoseph Pusztay   y[0] += A[0] * z[0] + A[1] * z[1];
229d0c080abSJoseph Pusztay   y[ldx] += A[2] * z[0] + A[3] * z[1];
230d0c080abSJoseph Pusztay   (void)PetscLogFlops(6.0);
231d0c080abSJoseph Pusztay }
232d0c080abSJoseph Pusztay 
233d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */
DMPlex_MultAdd3DReal_Internal(const PetscReal A[],PetscInt ldx,const PetscScalar x[],PetscScalar y[])234d71ae5a4SJacob Faibussowitsch static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
235d71ae5a4SJacob Faibussowitsch {
236d0c080abSJoseph Pusztay   PetscScalar z[3];
2379371c9d4SSatish Balay   z[0] = x[0];
2389371c9d4SSatish Balay   z[1] = x[ldx];
2399371c9d4SSatish Balay   z[2] = x[ldx * 2];
240d0c080abSJoseph Pusztay   y[0] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
241d0c080abSJoseph Pusztay   y[ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
242d0c080abSJoseph Pusztay   y[ldx * 2] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
243d0c080abSJoseph Pusztay   (void)PetscLogFlops(15.0);
244d0c080abSJoseph Pusztay }
245d0c080abSJoseph Pusztay 
246d0c080abSJoseph Pusztay /*
247d0c080abSJoseph Pusztay   Gaussian - The Gaussian function G(x)
248d0c080abSJoseph Pusztay 
249d0c080abSJoseph Pusztay   Input Parameters:
250d0c080abSJoseph Pusztay +  dim   - The number of dimensions, or size of x
251d0c080abSJoseph Pusztay .  mu    - The mean, or center
252d0c080abSJoseph Pusztay .  sigma - The standard deviation, or width
253d0c080abSJoseph Pusztay -  x     - The evaluation point of the function
254d0c080abSJoseph Pusztay 
255d0c080abSJoseph Pusztay   Output Parameter:
256d0c080abSJoseph Pusztay . ret - The value G(x)
257d0c080abSJoseph Pusztay */
Gaussian(PetscInt dim,const PetscReal mu[],PetscReal sigma,const PetscReal x[])258d71ae5a4SJacob Faibussowitsch static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[])
259d71ae5a4SJacob Faibussowitsch {
260d0c080abSJoseph Pusztay   PetscReal arg = 0.0;
261d0c080abSJoseph Pusztay   PetscInt  d;
262d0c080abSJoseph Pusztay 
263d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]);
264d0c080abSJoseph Pusztay   return PetscPowReal(2.0 * PETSC_PI * sigma, -dim / 2.0) * PetscExpReal(-arg / (2.0 * sigma));
265d0c080abSJoseph Pusztay }
266d0c080abSJoseph Pusztay 
267d0c080abSJoseph Pusztay /*
268d0c080abSJoseph Pusztay   ComputeGradS - Compute grad_v dS_eps/df
269d0c080abSJoseph Pusztay 
270d0c080abSJoseph Pusztay   Input Parameters:
271d0c080abSJoseph Pusztay + dim      - The dimension
272d0c080abSJoseph Pusztay . Np       - The number of particles
273d0c080abSJoseph Pusztay . vp       - The velocity v_p of the particle at which we evaluate
274d0c080abSJoseph Pusztay . velocity - The velocity field for all particles
275d0c080abSJoseph Pusztay . epsilon  - The regularization strength
276d0c080abSJoseph Pusztay 
277d0c080abSJoseph Pusztay   Output Parameter:
278d0c080abSJoseph Pusztay . integral - The output grad_v dS_eps/df (v_p)
279d0c080abSJoseph Pusztay 
280d0c080abSJoseph Pusztay   Note:
281d0c080abSJoseph Pusztay   This comes from (3.6) in [1], and we are computing
282d0c080abSJoseph Pusztay $   \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q)
283d0c080abSJoseph Pusztay   which is discretized by using a one-point quadrature in each box l at its center v^c_l
284d0c080abSJoseph 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)
285d0c080abSJoseph Pusztay   where h^d is the volume of each box.
286d0c080abSJoseph Pusztay */
ComputeGradS(PetscInt dim,PetscInt Np,const PetscReal vp[],const PetscReal velocity[],PetscReal integral[],AppCtx * ctx)287d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx)
288d71ae5a4SJacob Faibussowitsch {
289d0c080abSJoseph Pusztay   PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5 * h - L;
290d0c080abSJoseph Pusztay   PetscInt  nx = roundf(2. * L / h);
291d0c080abSJoseph Pusztay   PetscInt  ny = dim > 1 ? nx : 1;
292d0c080abSJoseph Pusztay   PetscInt  nz = dim > 2 ? nx : 1;
293d0c080abSJoseph Pusztay   PetscInt  i, j, k, d, q, dbg = 0;
294d0c080abSJoseph Pusztay 
295d0c080abSJoseph Pusztay   PetscFunctionBeginHot;
296d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) integral[d] = 0.0;
297d0c080abSJoseph Pusztay   for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) {
298d0c080abSJoseph Pusztay     for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) {
299d0c080abSJoseph Pusztay       for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) {
300d0c080abSJoseph Pusztay         PetscReal sum = 0.0;
301d0c080abSJoseph Pusztay 
30263a3b9bcSJacob 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]));
303d0c080abSJoseph Pusztay         /* \log \sum_k \psi(v - v_k)  */
304d0c080abSJoseph Pusztay         for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q * dim], epsilon, vc_l);
305d0c080abSJoseph Pusztay         sum = PetscLogReal(sum);
30657508eceSPierre Jolivet         for (d = 0; d < dim; ++d) integral[d] += (-1. / epsilon) * PetscAbsReal(vp[d] - vc_l[d]) * Gaussian(dim, vp, epsilon, vc_l) * sum;
307d0c080abSJoseph Pusztay       }
308d0c080abSJoseph Pusztay     }
309d0c080abSJoseph Pusztay   }
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311d0c080abSJoseph Pusztay }
312d0c080abSJoseph Pusztay 
313d0c080abSJoseph Pusztay /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */
QCompute(PetscInt dim,const PetscReal vp[],const PetscReal vq[],PetscReal Q[])314d71ae5a4SJacob Faibussowitsch static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[])
315d71ae5a4SJacob Faibussowitsch {
316d0c080abSJoseph Pusztay   PetscReal xi[3], xi2, xi3, mag;
317d0c080abSJoseph Pusztay   PetscInt  d, e;
318d0c080abSJoseph Pusztay 
319d0c080abSJoseph Pusztay   PetscFunctionBeginHot;
320d0c080abSJoseph Pusztay   DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi);
321d0c080abSJoseph Pusztay   xi2 = DMPlex_DotD_Internal(dim, xi, xi);
322d0c080abSJoseph Pusztay   mag = PetscSqrtReal(xi2);
323d0c080abSJoseph Pusztay   xi3 = xi2 * mag;
324d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) {
325ad540459SPierre Jolivet     for (e = 0; e < dim; ++e) Q[d * dim + e] = -xi[d] * xi[e] / xi3;
326d0c080abSJoseph Pusztay     Q[d * dim + d] += 1. / mag;
327d0c080abSJoseph Pusztay   }
3283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
329d0c080abSJoseph Pusztay }
330d0c080abSJoseph Pusztay 
RHSFunctionParticles(TS ts,PetscReal t,Vec U,Vec R,PetscCtx ctx)331*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, PetscCtx ctx)
332d71ae5a4SJacob Faibussowitsch {
333d0c080abSJoseph Pusztay   AppCtx            *user = (AppCtx *)ctx;
334d0c080abSJoseph Pusztay   PetscInt           dbg  = 0;
335d0c080abSJoseph Pusztay   DM                 sw;  /* Particles */
336d0c080abSJoseph Pusztay   Vec                sol; /* Solution vector at current time */
337d0c080abSJoseph Pusztay   const PetscScalar *u;   /* input solution vector */
338d0c080abSJoseph Pusztay   PetscScalar       *r;
339d0c080abSJoseph Pusztay   PetscReal         *velocity;
340d0c080abSJoseph Pusztay   PetscInt           dim, Np, p, q;
341d0c080abSJoseph Pusztay 
342d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
3439566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(R));
3449566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
3459566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
3469566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
3479566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &sol));
3489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(sol, &velocity));
3499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
3509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
351d0c080abSJoseph Pusztay   Np /= dim;
3529566063dSJacob Faibussowitsch   if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Part  ppr     x        y\n"));
353d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
354d0c080abSJoseph Pusztay     PetscReal gradS_p[3] = {0., 0., 0.};
355d0c080abSJoseph Pusztay 
3569566063dSJacob Faibussowitsch     PetscCall(ComputeGradS(dim, Np, &velocity[p * dim], velocity, gradS_p, user));
357d0c080abSJoseph Pusztay     for (q = 0; q < Np; ++q) {
358d0c080abSJoseph Pusztay       PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9];
359d0c080abSJoseph Pusztay 
360d0c080abSJoseph Pusztay       if (q == p) continue;
3619566063dSJacob Faibussowitsch       PetscCall(ComputeGradS(dim, Np, &velocity[q * dim], velocity, gradS_q, user));
362d0c080abSJoseph Pusztay       DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS);
3639566063dSJacob Faibussowitsch       PetscCall(QCompute(dim, &u[p * dim], &u[q * dim], Q));
364d0c080abSJoseph Pusztay       switch (dim) {
365d71ae5a4SJacob Faibussowitsch       case 2:
366d71ae5a4SJacob Faibussowitsch         DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p * dim]);
367d71ae5a4SJacob Faibussowitsch         break;
368d71ae5a4SJacob Faibussowitsch       case 3:
369d71ae5a4SJacob Faibussowitsch         DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p * dim]);
370d71ae5a4SJacob Faibussowitsch         break;
371d71ae5a4SJacob Faibussowitsch       default:
372d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %" PetscInt_FMT, dim);
373d0c080abSJoseph Pusztay       }
374d0c080abSJoseph Pusztay     }
37563a3b9bcSJacob 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]));
376d0c080abSJoseph Pusztay   }
3779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
3789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
3799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(sol, &velocity));
3809566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(R, NULL, "-residual_view"));
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
382d0c080abSJoseph Pusztay }
383d0c080abSJoseph Pusztay 
384d0c080abSJoseph Pusztay /*
385d0c080abSJoseph Pusztay  TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform
386d0c080abSJoseph Pusztay  the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid
387d0c080abSJoseph Pusztay  to migrate between.
388d0c080abSJoseph Pusztay */
UpdateSwarm(TS ts)389d71ae5a4SJacob Faibussowitsch static PetscErrorCode UpdateSwarm(TS ts)
390d71ae5a4SJacob Faibussowitsch {
391d0c080abSJoseph Pusztay   PetscInt           idx, n;
392d0c080abSJoseph Pusztay   const PetscScalar *u;
393d0c080abSJoseph Pusztay   PetscScalar       *velocity;
394d0c080abSJoseph Pusztay   DM                 sw;
395d0c080abSJoseph Pusztay   Vec                sol;
396d0c080abSJoseph Pusztay 
397d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
3989566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sw));
3999566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity));
4009566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &sol));
4019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(sol, &u));
4029566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(sol, &n));
403d0c080abSJoseph Pusztay   for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx];
4049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(sol, &u));
4059566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
407d0c080abSJoseph Pusztay }
408d0c080abSJoseph Pusztay 
InitializeSolve(TS ts,Vec u)409d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitializeSolve(TS ts, Vec u)
410d71ae5a4SJacob Faibussowitsch {
411d0c080abSJoseph Pusztay   DM      dm;
412d0c080abSJoseph Pusztay   AppCtx *user;
413d0c080abSJoseph Pusztay 
414d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
4159566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
4169566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &user));
4179566063dSJacob Faibussowitsch   PetscCall(SetInitialCoordinates(dm));
4189566063dSJacob Faibussowitsch   PetscCall(SetInitialConditions(dm, u));
4193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
420d0c080abSJoseph Pusztay }
421d0c080abSJoseph Pusztay 
main(int argc,char ** argv)422d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
423d71ae5a4SJacob Faibussowitsch {
424d0c080abSJoseph Pusztay   TS       ts;     /* nonlinear solver */
425d0c080abSJoseph Pusztay   DM       dm, sw; /* Velocity space mesh and Particle Swarm */
426d0c080abSJoseph Pusztay   Vec      u, v;   /* problem vector */
427d0c080abSJoseph Pusztay   MPI_Comm comm;
428d0c080abSJoseph Pusztay   AppCtx   user;
429d0c080abSJoseph Pusztay 
430327415f7SBarry Smith   PetscFunctionBeginUser;
4319566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
432d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
4339566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
434d0c080abSJoseph Pusztay   /* Initialize objects and set initial conditions */
4359566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
4369566063dSJacob Faibussowitsch   PetscCall(CreateParticles(dm, &sw, &user));
4379566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(sw, &user));
4389566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(sw, "velocity"));
4399566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
4409566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, sw));
4419566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10.0));
4429566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.1));
4439566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1));
4449566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
4459566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user));
4469566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
4479566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
4489566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v));
4499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v, &u));
4509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v));
4519566063dSJacob Faibussowitsch   PetscCall(TSComputeInitialCondition(ts, u));
4529566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, UpdateSwarm));
4539566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
4549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
4559566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
4569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
4579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
459b122ec5aSJacob Faibussowitsch   return 0;
460d0c080abSJoseph Pusztay }
461d0c080abSJoseph Pusztay 
462d0c080abSJoseph Pusztay /*TEST
463d0c080abSJoseph Pusztay    build:
464d0c080abSJoseph Pusztay      requires: triangle !single !complex
465d0c080abSJoseph Pusztay    test:
466d0c080abSJoseph Pusztay      suffix: midpoint
46730602db0SMatthew 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 \
46841d17464SJames Wright            -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_dmswarm_monitor_moments_interval 1 -snes_fd
469d0c080abSJoseph Pusztay TEST*/
470