xref: /petsc/src/ts/tests/ex27.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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 
32d0c080abSJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33d0c080abSJoseph Pusztay {
34d0c080abSJoseph Pusztay   PetscErrorCode ierr;
35d0c080abSJoseph Pusztay 
36d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
37d0c080abSJoseph Pusztay   options->N         = 1;
38d0c080abSJoseph Pusztay   options->momentTol = 100.0*PETSC_MACHINE_EPSILON;
39d0c080abSJoseph Pusztay   options->L         = 1.0;
40d0c080abSJoseph Pusztay   options->h         = -1.0;
41d0c080abSJoseph Pusztay   options->epsilon   = -1.0;
42d0c080abSJoseph Pusztay 
43d0c080abSJoseph Pusztay   ierr = PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");CHKERRQ(ierr);
44d0c080abSJoseph Pusztay   ierr = PetscOptionsInt("-N", "Number of particles per spatial cell", "ex27.c", options->N, &options->N, NULL);CHKERRQ(ierr);
45d0c080abSJoseph Pusztay   ierr = PetscOptionsReal("-L", "Velocity-space extent", "ex27.c", options->L, &options->L, NULL);CHKERRQ(ierr);
46d0c080abSJoseph Pusztay   ierr = PetscOptionsReal("-h", "Velocity-space resolution", "ex27.c", options->h, &options->h, NULL);CHKERRQ(ierr);
47d0c080abSJoseph Pusztay   ierr = PetscOptionsReal("-epsilon", "Mollifier regularization parameter", "ex27.c", options->epsilon, &options->epsilon, NULL);CHKERRQ(ierr);
48d0c080abSJoseph Pusztay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
49d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
50d0c080abSJoseph Pusztay }
51d0c080abSJoseph Pusztay 
52d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
53d0c080abSJoseph Pusztay {
54d0c080abSJoseph Pusztay   PetscErrorCode ierr;
55d0c080abSJoseph Pusztay 
56d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
5730602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
5830602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
59d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
60d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
61d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
62d0c080abSJoseph Pusztay }
63d0c080abSJoseph Pusztay 
64d0c080abSJoseph Pusztay static PetscErrorCode SetInitialCoordinates(DM sw)
65d0c080abSJoseph Pusztay {
66d0c080abSJoseph Pusztay   AppCtx        *user;
67d0c080abSJoseph Pusztay   PetscRandom    rnd, rndv;
68d0c080abSJoseph Pusztay   DM             dm;
69d0c080abSJoseph Pusztay   DMPolytopeType ct;
70d0c080abSJoseph Pusztay   PetscBool      simplex;
71d0c080abSJoseph Pusztay   PetscReal     *centroid, *coords, *velocity, *xi0, *v0, *J, *invJ, detJ, *vals;
72d0c080abSJoseph Pusztay   PetscInt       dim, d, cStart, cEnd, c, Np, p;
73d0c080abSJoseph Pusztay   PetscErrorCode ierr;
74d0c080abSJoseph Pusztay 
75d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
76d0c080abSJoseph Pusztay   /* Randomization for coordinates */
77d0c080abSJoseph Pusztay   ierr = PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rnd);CHKERRQ(ierr);
78d0c080abSJoseph Pusztay   ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr);
79d0c080abSJoseph Pusztay   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
80d0c080abSJoseph Pusztay   ierr = PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rndv);CHKERRQ(ierr);
81d0c080abSJoseph Pusztay   ierr = PetscRandomSetInterval(rndv, -1., 1.);CHKERRQ(ierr);
82d0c080abSJoseph Pusztay   ierr = PetscRandomSetFromOptions(rndv);CHKERRQ(ierr);
83d0c080abSJoseph Pusztay   ierr = DMGetApplicationContext(sw, (void **) &user);CHKERRQ(ierr);
84d0c080abSJoseph Pusztay   Np   = user->N;
85d0c080abSJoseph Pusztay   ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr);
86d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellDM(sw, &dm);CHKERRQ(ierr);
87d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
88d0c080abSJoseph Pusztay   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
89d0c080abSJoseph Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
90d0c080abSJoseph Pusztay   ierr = PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr);
91d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
92d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
93d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &velocity);CHKERRQ(ierr);
94d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr);
95d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
96d0c080abSJoseph Pusztay     if (Np == 1) {
97d0c080abSJoseph Pusztay       ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
98d0c080abSJoseph Pusztay       for (d = 0; d < dim; ++d) {
99d0c080abSJoseph Pusztay         coords[c*dim+d] = centroid[d];
100d0c080abSJoseph Pusztay       }
101d0c080abSJoseph Pusztay       vals[c] = 1.0;
102d0c080abSJoseph Pusztay     } else {
103d0c080abSJoseph Pusztay       ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */
104d0c080abSJoseph Pusztay       for (p = 0; p < Np; ++p) {
105d0c080abSJoseph Pusztay         const PetscInt n   = c*Np + p;
106d0c080abSJoseph Pusztay         PetscReal      sum = 0.0, refcoords[3];
107d0c080abSJoseph Pusztay 
108d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) {
109d0c080abSJoseph Pusztay           ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr);
110d0c080abSJoseph Pusztay           sum += refcoords[d];
111d0c080abSJoseph Pusztay         }
112d0c080abSJoseph Pusztay         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
113d0c080abSJoseph Pusztay         vals[n] = 1.0;
114d0c080abSJoseph Pusztay         ierr = DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n*dim]);CHKERRQ(ierr);
115d0c080abSJoseph Pusztay       }
116d0c080abSJoseph Pusztay     }
117d0c080abSJoseph Pusztay   }
118d0c080abSJoseph Pusztay   /* Random velocity IC */
119d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
120d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
121d0c080abSJoseph Pusztay       for (d = 0; d < dim; ++d) {
122d0c080abSJoseph Pusztay         PetscReal v_val;
123d0c080abSJoseph Pusztay 
124*1e1ea65dSPierre Jolivet         ierr = PetscRandomGetValueReal(rndv, &v_val);CHKERRQ(ierr);
125d0c080abSJoseph Pusztay         velocity[p*dim+d] = v_val;
126d0c080abSJoseph Pusztay       }
127d0c080abSJoseph Pusztay     }
128d0c080abSJoseph Pusztay   }
129d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
130d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &velocity);CHKERRQ(ierr);
131d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr);
132d0c080abSJoseph Pusztay   ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr);
133d0c080abSJoseph Pusztay   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
134d0c080abSJoseph Pusztay   ierr = PetscRandomDestroy(&rndv);CHKERRQ(ierr);
135d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
136d0c080abSJoseph Pusztay }
137d0c080abSJoseph Pusztay 
138d0c080abSJoseph Pusztay /* Get velocities from swarm and place in solution vector */
139d0c080abSJoseph Pusztay static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
140d0c080abSJoseph Pusztay {
141d0c080abSJoseph Pusztay   DM             dm;
142d0c080abSJoseph Pusztay   AppCtx        *user;
143d0c080abSJoseph Pusztay   PetscReal     *velocity;
144d0c080abSJoseph Pusztay   PetscScalar   *initialConditions;
145d0c080abSJoseph Pusztay   PetscInt       dim, d, cStart, cEnd, c, Np, p, n;
146d0c080abSJoseph Pusztay   PetscErrorCode ierr;
147d0c080abSJoseph Pusztay 
148d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
149d0c080abSJoseph Pusztay   ierr = VecGetLocalSize(u, &n);CHKERRQ(ierr);
150d0c080abSJoseph Pusztay   ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr);
151d0c080abSJoseph Pusztay   Np   = user->N;
152d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
153d0c080abSJoseph Pusztay   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
154d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
155d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **) &velocity);CHKERRQ(ierr);
156d0c080abSJoseph Pusztay   ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr);
157d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
158d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
159d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
160d0c080abSJoseph Pusztay       for (d = 0; d < dim; d++) {
161d0c080abSJoseph Pusztay         initialConditions[n*dim+d] = velocity[n*dim+d];
162d0c080abSJoseph Pusztay       }
163d0c080abSJoseph Pusztay     }
164d0c080abSJoseph Pusztay   }
165d0c080abSJoseph Pusztay   ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr);
166d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **) &velocity);CHKERRQ(ierr);
167d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
168d0c080abSJoseph Pusztay }
169d0c080abSJoseph Pusztay 
170d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
171d0c080abSJoseph Pusztay {
172d0c080abSJoseph Pusztay   PetscInt      *cellid;
173d0c080abSJoseph Pusztay   PetscInt       dim, cStart, cEnd, c, Np = user->N, p;
174d0c080abSJoseph Pusztay   PetscBool      view = PETSC_FALSE;
175d0c080abSJoseph Pusztay   PetscErrorCode ierr;
176d0c080abSJoseph Pusztay 
177d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
178d0c080abSJoseph Pusztay   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
179d0c080abSJoseph Pusztay   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
180d0c080abSJoseph Pusztay   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
181d0c080abSJoseph Pusztay   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
182d0c080abSJoseph Pusztay   /* h = 2L/n and N = n^d */
183d0c080abSJoseph Pusztay   if (user->h < 0.) user->h = 2.*user->L / PetscPowReal(user->N, 1./dim);
184d0c080abSJoseph Pusztay   /* From Section 4 in [1], \epsilon = 0.64 h^.98 */
185d0c080abSJoseph Pusztay   if (user->epsilon < 0.) user->epsilon = 0.64*pow(user->h, 1.98);
186d0c080abSJoseph Pusztay   ierr = PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL);CHKERRQ(ierr);
187d0c080abSJoseph Pusztay   if (view) {
188d0c080abSJoseph Pusztay     ierr = PetscPrintf(PETSC_COMM_SELF, "N: %D L: %g h: %g eps: %g\n", user->N, user->L, user->h, user->epsilon);CHKERRQ(ierr);
189d0c080abSJoseph Pusztay   }
190d0c080abSJoseph Pusztay   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
191d0c080abSJoseph Pusztay   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
192d0c080abSJoseph Pusztay   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL);CHKERRQ(ierr);
193d0c080abSJoseph Pusztay   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR);CHKERRQ(ierr);
194d0c080abSJoseph Pusztay   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
195d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
196d0c080abSJoseph Pusztay   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
197d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
198d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
199d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
200d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
201d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
202d0c080abSJoseph Pusztay       cellid[n] = c;
203d0c080abSJoseph Pusztay     }
204d0c080abSJoseph Pusztay   }
205d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
206d0c080abSJoseph Pusztay   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
207d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
208d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
209d0c080abSJoseph Pusztay }
210d0c080abSJoseph Pusztay 
211d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
212d0c080abSJoseph Pusztay static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
213d0c080abSJoseph Pusztay {
214d0c080abSJoseph Pusztay   PetscInt d;
215d0c080abSJoseph Pusztay 
216d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d];
217d0c080abSJoseph Pusztay }
218d0c080abSJoseph Pusztay 
219d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
220d0c080abSJoseph Pusztay static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
221d0c080abSJoseph Pusztay {
222d0c080abSJoseph Pusztay   PetscReal sum = 0.0;
223d0c080abSJoseph Pusztay   PetscInt d;
224d0c080abSJoseph Pusztay 
225d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d])*y[d];
226d0c080abSJoseph Pusztay   return sum;
227d0c080abSJoseph Pusztay }
228d0c080abSJoseph Pusztay 
229d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */
230d0c080abSJoseph Pusztay static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
231d0c080abSJoseph Pusztay {
232d0c080abSJoseph Pusztay   PetscScalar z[2];
233d0c080abSJoseph Pusztay   z[0] = x[0]; z[1] = x[ldx];
234d0c080abSJoseph Pusztay   y[0]   += A[0]*z[0] + A[1]*z[1];
235d0c080abSJoseph Pusztay   y[ldx] += A[2]*z[0] + A[3]*z[1];
236d0c080abSJoseph Pusztay   (void)PetscLogFlops(6.0);
237d0c080abSJoseph Pusztay }
238d0c080abSJoseph Pusztay 
239d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */
240d0c080abSJoseph Pusztay static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
241d0c080abSJoseph Pusztay {
242d0c080abSJoseph Pusztay   PetscScalar z[3];
243d0c080abSJoseph Pusztay   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
244d0c080abSJoseph Pusztay   y[0]     += A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
245d0c080abSJoseph Pusztay   y[ldx]   += A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
246d0c080abSJoseph Pusztay   y[ldx*2] += A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
247d0c080abSJoseph Pusztay   (void)PetscLogFlops(15.0);
248d0c080abSJoseph Pusztay }
249d0c080abSJoseph Pusztay 
250d0c080abSJoseph Pusztay /*
251d0c080abSJoseph Pusztay   Gaussian - The Gaussian function G(x)
252d0c080abSJoseph Pusztay 
253d0c080abSJoseph Pusztay   Input Parameters:
254d0c080abSJoseph Pusztay +  dim   - The number of dimensions, or size of x
255d0c080abSJoseph Pusztay .  mu    - The mean, or center
256d0c080abSJoseph Pusztay .  sigma - The standard deviation, or width
257d0c080abSJoseph Pusztay -  x     - The evaluation point of the function
258d0c080abSJoseph Pusztay 
259d0c080abSJoseph Pusztay   Output Parameter:
260d0c080abSJoseph Pusztay . ret - The value G(x)
261d0c080abSJoseph Pusztay */
262d0c080abSJoseph Pusztay static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[])
263d0c080abSJoseph Pusztay {
264d0c080abSJoseph Pusztay   PetscReal arg = 0.0;
265d0c080abSJoseph Pusztay   PetscInt  d;
266d0c080abSJoseph Pusztay 
267d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]);
268d0c080abSJoseph Pusztay   return PetscPowReal(2.0*PETSC_PI*sigma, -dim/2.0) * PetscExpReal(-arg/(2.0*sigma));
269d0c080abSJoseph Pusztay }
270d0c080abSJoseph Pusztay 
271d0c080abSJoseph Pusztay /*
272d0c080abSJoseph Pusztay   ComputeGradS - Compute grad_v dS_eps/df
273d0c080abSJoseph Pusztay 
274d0c080abSJoseph Pusztay   Input Parameters:
275d0c080abSJoseph Pusztay + dim      - The dimension
276d0c080abSJoseph Pusztay . Np       - The number of particles
277d0c080abSJoseph Pusztay . vp       - The velocity v_p of the particle at which we evaluate
278d0c080abSJoseph Pusztay . velocity - The velocity field for all particles
279d0c080abSJoseph Pusztay . epsilon  - The regularization strength
280d0c080abSJoseph Pusztay 
281d0c080abSJoseph Pusztay   Output Parameter:
282d0c080abSJoseph Pusztay . integral - The output grad_v dS_eps/df (v_p)
283d0c080abSJoseph Pusztay 
284d0c080abSJoseph Pusztay   Note:
285d0c080abSJoseph Pusztay   This comes from (3.6) in [1], and we are computing
286d0c080abSJoseph Pusztay $   \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q)
287d0c080abSJoseph Pusztay   which is discretized by using a one-point quadrature in each box l at its center v^c_l
288d0c080abSJoseph 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)
289d0c080abSJoseph Pusztay   where h^d is the volume of each box.
290d0c080abSJoseph Pusztay */
291d0c080abSJoseph Pusztay static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx)
292d0c080abSJoseph Pusztay {
293d0c080abSJoseph Pusztay   PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5*h - L;
294d0c080abSJoseph Pusztay   PetscInt  nx = roundf(2.*L / h);
295d0c080abSJoseph Pusztay   PetscInt  ny = dim > 1 ? nx : 1;
296d0c080abSJoseph Pusztay   PetscInt  nz = dim > 2 ? nx : 1;
297d0c080abSJoseph Pusztay   PetscInt  i, j, k, d, q, dbg = 0;
298d0c080abSJoseph Pusztay   PetscErrorCode ierr;
299d0c080abSJoseph Pusztay 
300d0c080abSJoseph Pusztay   PetscFunctionBeginHot;
301d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) integral[d] = 0.0;
302d0c080abSJoseph Pusztay   for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) {
303d0c080abSJoseph Pusztay     for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) {
304d0c080abSJoseph Pusztay       for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) {
305d0c080abSJoseph Pusztay         PetscReal sum = 0.0;
306d0c080abSJoseph Pusztay 
307d0c080abSJoseph Pusztay         if (dbg) {ierr = PetscPrintf(PETSC_COMM_SELF, "(%D %D) vc_l: %g %g\n", i, j, vc_l[0], vc_l[1]);CHKERRQ(ierr);}
308d0c080abSJoseph Pusztay         /* \log \sum_k \psi(v - v_k)  */
309d0c080abSJoseph Pusztay         for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q*dim], epsilon, vc_l);
310d0c080abSJoseph Pusztay         sum = PetscLogReal(sum);
311d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) integral[d] += (-1./(epsilon))*PetscAbsReal(vp[d] - vc_l[d])*(Gaussian(dim, vp, epsilon, vc_l)) * sum;
312d0c080abSJoseph Pusztay       }
313d0c080abSJoseph Pusztay     }
314d0c080abSJoseph Pusztay   }
315d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
316d0c080abSJoseph Pusztay }
317d0c080abSJoseph Pusztay 
318d0c080abSJoseph Pusztay /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */
319d0c080abSJoseph Pusztay static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[])
320d0c080abSJoseph Pusztay {
321d0c080abSJoseph Pusztay   PetscReal xi[3], xi2, xi3, mag;
322d0c080abSJoseph Pusztay   PetscInt  d, e;
323d0c080abSJoseph Pusztay 
324d0c080abSJoseph Pusztay   PetscFunctionBeginHot;
325d0c080abSJoseph Pusztay   DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi);
326d0c080abSJoseph Pusztay   xi2 = DMPlex_DotD_Internal(dim, xi, xi);
327d0c080abSJoseph Pusztay   mag = PetscSqrtReal(xi2);
328d0c080abSJoseph Pusztay   xi3 = xi2 * mag;
329d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) {
330d0c080abSJoseph Pusztay     for (e = 0; e < dim; ++e) {
331d0c080abSJoseph Pusztay       Q[d*dim+e] = -xi[d]*xi[e] / xi3;
332d0c080abSJoseph Pusztay     }
333d0c080abSJoseph Pusztay     Q[d*dim+d] += 1. / mag;
334d0c080abSJoseph Pusztay   }
335d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
336d0c080abSJoseph Pusztay }
337d0c080abSJoseph Pusztay 
338d0c080abSJoseph Pusztay static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
339d0c080abSJoseph Pusztay {
340d0c080abSJoseph Pusztay   AppCtx            *user = (AppCtx *) ctx;
341d0c080abSJoseph Pusztay   PetscInt           dbg  = 0;
342d0c080abSJoseph Pusztay   DM                 sw;                  /* Particles */
343d0c080abSJoseph Pusztay   Vec                sol;                 /* Solution vector at current time */
344d0c080abSJoseph Pusztay   const PetscScalar *u;                   /* input solution vector */
345d0c080abSJoseph Pusztay   PetscScalar       *r;
346d0c080abSJoseph Pusztay   PetscReal         *velocity;
347d0c080abSJoseph Pusztay   PetscInt           dim, Np, p, q;
348d0c080abSJoseph Pusztay   PetscErrorCode     ierr;
349d0c080abSJoseph Pusztay 
350d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
351d0c080abSJoseph Pusztay   ierr = VecZeroEntries(R);CHKERRQ(ierr);
352d0c080abSJoseph Pusztay   ierr = TSGetDM(ts, &sw);CHKERRQ(ierr);CHKERRQ(ierr);
353*1e1ea65dSPierre Jolivet   ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr);
354d0c080abSJoseph Pusztay   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
355d0c080abSJoseph Pusztay   ierr = TSGetSolution(ts, &sol);CHKERRQ(ierr);
356*1e1ea65dSPierre Jolivet   ierr = VecGetArray(sol, &velocity);CHKERRQ(ierr);
357*1e1ea65dSPierre Jolivet   ierr = VecGetArray(R, &r);CHKERRQ(ierr);
358*1e1ea65dSPierre Jolivet   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
359d0c080abSJoseph Pusztay   Np  /= dim;
360d0c080abSJoseph Pusztay   if (dbg) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Part  ppr     x        y\n");CHKERRQ(ierr);}
361d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
362d0c080abSJoseph Pusztay     PetscReal gradS_p[3] = {0., 0., 0.};
363d0c080abSJoseph Pusztay 
364d0c080abSJoseph Pusztay     ierr = ComputeGradS(dim, Np, &velocity[p*dim], velocity, gradS_p, user);CHKERRQ(ierr);
365d0c080abSJoseph Pusztay     for (q = 0; q < Np; ++q) {
366d0c080abSJoseph Pusztay       PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9];
367d0c080abSJoseph Pusztay 
368d0c080abSJoseph Pusztay       if (q == p) continue;
369d0c080abSJoseph Pusztay       ierr = ComputeGradS(dim, Np, &velocity[q*dim], velocity, gradS_q, user);CHKERRQ(ierr);
370d0c080abSJoseph Pusztay       DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS);
371d0c080abSJoseph Pusztay       ierr = QCompute(dim, &u[p*dim], &u[q*dim], Q);CHKERRQ(ierr);
372d0c080abSJoseph Pusztay       switch (dim) {
373d0c080abSJoseph Pusztay         case 2: DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p*dim]);break;
374d0c080abSJoseph Pusztay         case 3: DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p*dim]);break;
375d0c080abSJoseph Pusztay         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %D", dim);
376d0c080abSJoseph Pusztay       }
377d0c080abSJoseph Pusztay     }
378d0c080abSJoseph Pusztay     if (dbg) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final %4D %10.8lf %10.8lf\n", p, r[p*dim+0], r[p*dim+1]);CHKERRQ(ierr);}
379d0c080abSJoseph Pusztay   }
380d0c080abSJoseph Pusztay   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
381d0c080abSJoseph Pusztay   ierr = VecRestoreArray(R, &r);CHKERRQ(ierr);
382d0c080abSJoseph Pusztay   ierr = VecRestoreArray(sol, &velocity);CHKERRQ(ierr);
383d0c080abSJoseph Pusztay   ierr = VecViewFromOptions(R, NULL, "-residual_view");CHKERRQ(ierr);
384d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
385d0c080abSJoseph Pusztay }
386d0c080abSJoseph Pusztay 
387d0c080abSJoseph Pusztay /*
388d0c080abSJoseph Pusztay  TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform
389d0c080abSJoseph Pusztay  the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid
390d0c080abSJoseph Pusztay  to migrate between.
391d0c080abSJoseph Pusztay */
392d0c080abSJoseph Pusztay static PetscErrorCode UpdateSwarm(TS ts)
393d0c080abSJoseph Pusztay {
394d0c080abSJoseph Pusztay   PetscInt idx, n;
395d0c080abSJoseph Pusztay   const PetscScalar *u;
396d0c080abSJoseph Pusztay   PetscScalar *velocity;
397d0c080abSJoseph Pusztay   DM sw;
398d0c080abSJoseph Pusztay   Vec sol;
399d0c080abSJoseph Pusztay   PetscErrorCode ierr;
400d0c080abSJoseph Pusztay 
401d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
402d0c080abSJoseph Pusztay   ierr = TSGetDM(ts, &sw);CHKERRQ(ierr);
403d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &velocity);CHKERRQ(ierr);
404d0c080abSJoseph Pusztay   ierr = TSGetSolution(ts, &sol);CHKERRQ(ierr);
405d0c080abSJoseph Pusztay   ierr = VecGetArrayRead(sol, &u);CHKERRQ(ierr);
406*1e1ea65dSPierre Jolivet   ierr = VecGetLocalSize(sol, &n);CHKERRQ(ierr);
407d0c080abSJoseph Pusztay   for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx];
408d0c080abSJoseph Pusztay   ierr = VecRestoreArrayRead(sol, &u);CHKERRQ(ierr);
409d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &velocity);CHKERRQ(ierr);
410d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
411d0c080abSJoseph Pusztay }
412d0c080abSJoseph Pusztay 
413d0c080abSJoseph Pusztay static PetscErrorCode InitializeSolve(TS ts, Vec u)
414d0c080abSJoseph Pusztay {
415d0c080abSJoseph Pusztay   DM             dm;
416d0c080abSJoseph Pusztay   AppCtx        *user;
417d0c080abSJoseph Pusztay   PetscErrorCode ierr;
418d0c080abSJoseph Pusztay 
419d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
420d0c080abSJoseph Pusztay   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
421d0c080abSJoseph Pusztay   ierr = DMGetApplicationContext(dm, (void **) &user);CHKERRQ(ierr);
422d0c080abSJoseph Pusztay   ierr = SetInitialCoordinates(dm);CHKERRQ(ierr);
423d0c080abSJoseph Pusztay   ierr = SetInitialConditions(dm, u);CHKERRQ(ierr);
424d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
425d0c080abSJoseph Pusztay }
426d0c080abSJoseph Pusztay 
427d0c080abSJoseph Pusztay int main(int argc,char **argv)
428d0c080abSJoseph Pusztay {
429d0c080abSJoseph Pusztay   TS             ts;     /* nonlinear solver */
430d0c080abSJoseph Pusztay   DM             dm, sw; /* Velocity space mesh and Particle Swarm */
431d0c080abSJoseph Pusztay   Vec            u, v;   /* problem vector */
432d0c080abSJoseph Pusztay   MPI_Comm       comm;
433d0c080abSJoseph Pusztay   AppCtx         user;
434d0c080abSJoseph Pusztay   PetscErrorCode ierr;
435d0c080abSJoseph Pusztay 
436d0c080abSJoseph Pusztay   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
437d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
438d0c080abSJoseph Pusztay   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
439d0c080abSJoseph Pusztay   /* Initialize objects and set initial conditions */
440d0c080abSJoseph Pusztay   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
441d0c080abSJoseph Pusztay   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
442d0c080abSJoseph Pusztay   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
443d0c080abSJoseph Pusztay   ierr = DMSwarmVectorDefineField(sw, "velocity");CHKERRQ(ierr);
444d0c080abSJoseph Pusztay   ierr = TSCreate(comm, &ts);CHKERRQ(ierr);
445d0c080abSJoseph Pusztay   ierr = TSSetDM(ts, sw);CHKERRQ(ierr);
446d0c080abSJoseph Pusztay   ierr = TSSetMaxTime(ts, 10.0);CHKERRQ(ierr);
447d0c080abSJoseph Pusztay   ierr = TSSetTimeStep(ts, 0.1);CHKERRQ(ierr);
448d0c080abSJoseph Pusztay   ierr = TSSetMaxSteps(ts, 1);CHKERRQ(ierr);
449d0c080abSJoseph Pusztay   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
450d0c080abSJoseph Pusztay   ierr = TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user);CHKERRQ(ierr);
451d0c080abSJoseph Pusztay   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
452d0c080abSJoseph Pusztay   ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr);
453d0c080abSJoseph Pusztay   ierr = DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v);CHKERRQ(ierr);
454d0c080abSJoseph Pusztay   ierr = VecDuplicate(v, &u);CHKERRQ(ierr);
455d0c080abSJoseph Pusztay   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v);CHKERRQ(ierr);
456d0c080abSJoseph Pusztay   ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr);
457d0c080abSJoseph Pusztay   ierr = TSSetPostStep(ts, UpdateSwarm);CHKERRQ(ierr);
458d0c080abSJoseph Pusztay   ierr = TSSolve(ts, u);CHKERRQ(ierr);
459d0c080abSJoseph Pusztay   ierr = VecDestroy(&u);CHKERRQ(ierr);
460d0c080abSJoseph Pusztay   ierr = TSDestroy(&ts);CHKERRQ(ierr);
461d0c080abSJoseph Pusztay   ierr = DMDestroy(&sw);CHKERRQ(ierr);
462d0c080abSJoseph Pusztay   ierr = DMDestroy(&dm);CHKERRQ(ierr);
463d0c080abSJoseph Pusztay   ierr = PetscFinalize();
464d0c080abSJoseph Pusztay   return ierr;
465d0c080abSJoseph Pusztay }
466d0c080abSJoseph Pusztay 
467d0c080abSJoseph Pusztay /*TEST
468d0c080abSJoseph Pusztay    build:
469d0c080abSJoseph Pusztay      requires: triangle !single !complex
470d0c080abSJoseph Pusztay    test:
471d0c080abSJoseph Pusztay      suffix: midpoint
47230602db0SMatthew 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 \
473d0c080abSJoseph Pusztay            -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_monitor_frequency 1 -snes_fd
474d0c080abSJoseph Pusztay TEST*/
475