xref: /petsc/src/ts/tutorials/hamiltonian/ex2.c (revision 84467f862f2de26368b63ea79dccd665bcda30cd)
1 #include "petscsys.h"
2 static char help[] = "Landau Damping/Two Stream instability test using Vlasov-Poisson equations\n";
3 
4 /*
5   To run the code with particles sinusoidally perturbed in x space use the test "pp_poisson_bsi_1d_4" or "pp_poisson_bsi_2d_4"
6   According to Lukas, good damping results come at ~16k particles per cell
7 
8   To visualize the efield use
9 
10     -monitor_efield
11 
12   To visualize the swarm distribution use
13 
14     -ts_monitor_hg_swarm
15 
16   To visualize the particles, we can use
17 
18     -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain 0 -ts_monitor_sp_swarm_phase 1 -draw_size 500,500
19 
20 For a Landau Damping verification run, we use
21 
22     -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 10,1 \
23       -dm_plex_box_lower 0.,-0.5 -dm_plex_box_upper 12.5664,0.5 -dm_plex_box_bd periodic,none \
24     -vdm_plex_dim 1 -vdm_plex_simplex 0 -vdm_plex_box_faces 2000 -vdm_plex_box_lower -10 -vdm_plex_box_upper 10 \
25     -dm_swarm_num_species 1 -charges -1.,1. \
26     -cosine_coefficients 0.01,0.5 -perturbed_weights -total_weight 1. \
27     -ts_type basicsymplectic -ts_basicsymplectic_type 1 -ts_dt 0.03 -ts_max_time 500 -ts_max_steps 500 \
28     -em_type primal -petscspace_degree 1 -em_snes_atol 1.e-12 -em_snes_error_if_not_converged -em_ksp_error_if_not_converged -em_pc_type svd \
29     -output_step 100 -check_vel_res -monitor_efield -ts_monitor -log_view
30 
31 */
32 #include <petscts.h>
33 #include <petscdmplex.h>
34 #include <petscdmswarm.h>
35 #include <petscfe.h>
36 #include <petscds.h>
37 #include <petscbag.h>
38 #include <petscdraw.h>
39 #include <petsc/private/dmpleximpl.h>  /* For norm and dot */
40 #include <petsc/private/petscfeimpl.h> /* For interpolation */
41 #include <petsc/private/dmswarmimpl.h> /* For swarm debugging */
42 #include "petscdm.h"
43 #include "petscdmlabel.h"
44 
45 PETSC_EXTERN PetscErrorCode stream(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
46 PETSC_EXTERN PetscErrorCode line(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
47 
48 const char *EMTypes[] = {"primal", "mixed", "coulomb", "none", "EMType", "EM_", NULL};
49 typedef enum {
50   EM_PRIMAL,
51   EM_MIXED,
52   EM_COULOMB,
53   EM_NONE
54 } EMType;
55 
56 typedef enum {
57   V0,
58   X0,
59   T0,
60   M0,
61   Q0,
62   PHI0,
63   POISSON,
64   VLASOV,
65   SIGMA,
66   NUM_CONSTANTS
67 } ConstantType;
68 typedef struct {
69   PetscScalar v0; /* Velocity scale, often the thermal velocity */
70   PetscScalar t0; /* Time scale */
71   PetscScalar x0; /* Space scale */
72   PetscScalar m0; /* Mass scale */
73   PetscScalar q0; /* Charge scale */
74   PetscScalar kb;
75   PetscScalar epsi0;
76   PetscScalar phi0;          /* Potential scale */
77   PetscScalar poissonNumber; /* Non-Dimensional Poisson Number */
78   PetscScalar vlasovNumber;  /* Non-Dimensional Vlasov Number */
79   PetscReal   sigma;         /* Nondimensional charge per length in x */
80 } Parameter;
81 
82 typedef struct {
83   PetscBag    bag;            /* Problem parameters */
84   PetscBool   error;          /* Flag for printing the error */
85   PetscBool   efield_monitor; /* Flag to show electric field monitor */
86   PetscBool   initial_monitor;
87   PetscBool   fake_1D;           /* Run simulation in 2D but zeroing second dimension */
88   PetscBool   perturbed_weights; /* Uniformly sample x,v space with gaussian weights */
89   PetscBool   poisson_monitor;
90   PetscInt    ostep; /* print the energy at each ostep time steps */
91   PetscInt    numParticles;
92   PetscReal   timeScale;              /* Nondimensionalizing time scale */
93   PetscReal   charges[2];             /* The charges of each species */
94   PetscReal   masses[2];              /* The masses of each species */
95   PetscReal   thermal_energy[2];      /* Thermal Energy (used to get other constants)*/
96   PetscReal   cosine_coefficients[2]; /*(alpha, k)*/
97   PetscReal   totalWeight;
98   PetscReal   stepSize;
99   PetscInt    steps;
100   PetscReal   initVel;
101   EMType      em; /* Type of electrostatic model */
102   SNES        snes;
103   PetscDraw   drawef;
104   PetscDrawLG drawlg_ef;
105   PetscDraw   drawic_x;
106   PetscDraw   drawic_v;
107   PetscDraw   drawic_w;
108   PetscDrawHG drawhgic_x;
109   PetscDrawHG drawhgic_v;
110   PetscDrawHG drawhgic_w;
111   PetscDraw   EDraw;
112   PetscDraw   RhoDraw;
113   PetscDraw   PotDraw;
114   PetscDrawSP EDrawSP;
115   PetscDrawSP RhoDrawSP;
116   PetscDrawSP PotDrawSP;
117   PetscBool   monitor_positions; /* Flag to show particle positins at each time step */
118   PetscDraw   positionDraw;
119   PetscDrawSP positionDrawSP;
120   DM          swarm;
121   PetscRandom random;
122   PetscBool   twostream;
123   PetscBool   checkweights;
124   PetscInt    checkVRes; /* Flag to check/output velocity residuals for nightly tests */
125 
126   PetscLogEvent RhsXEvent, RhsVEvent, ESolveEvent, ETabEvent;
127 } AppCtx;
128 
129 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
130 {
131   PetscFunctionBeginUser;
132   PetscInt d                      = 2;
133   PetscInt maxSpecies             = 2;
134   options->error                  = PETSC_FALSE;
135   options->efield_monitor         = PETSC_FALSE;
136   options->initial_monitor        = PETSC_FALSE;
137   options->fake_1D                = PETSC_FALSE;
138   options->perturbed_weights      = PETSC_FALSE;
139   options->poisson_monitor        = PETSC_FALSE;
140   options->ostep                  = 100;
141   options->timeScale              = 2.0e-14;
142   options->charges[0]             = -1.0;
143   options->charges[1]             = 1.0;
144   options->masses[0]              = 1.0;
145   options->masses[1]              = 1000.0;
146   options->thermal_energy[0]      = 1.0;
147   options->thermal_energy[1]      = 1.0;
148   options->cosine_coefficients[0] = 0.01;
149   options->cosine_coefficients[1] = 0.5;
150   options->initVel                = 1;
151   options->totalWeight            = 1.0;
152   options->drawef                 = NULL;
153   options->drawlg_ef              = NULL;
154   options->drawic_x               = NULL;
155   options->drawic_v               = NULL;
156   options->drawic_w               = NULL;
157   options->drawhgic_x             = NULL;
158   options->drawhgic_v             = NULL;
159   options->drawhgic_w             = NULL;
160   options->EDraw                  = NULL;
161   options->RhoDraw                = NULL;
162   options->PotDraw                = NULL;
163   options->EDrawSP                = NULL;
164   options->RhoDrawSP              = NULL;
165   options->PotDrawSP              = NULL;
166   options->em                     = EM_COULOMB;
167   options->numParticles           = 32768;
168   options->monitor_positions      = PETSC_FALSE;
169   options->positionDraw           = NULL;
170   options->positionDrawSP         = NULL;
171   options->twostream              = PETSC_FALSE;
172   options->checkweights           = PETSC_FALSE;
173   options->checkVRes              = 0;
174 
175   PetscOptionsBegin(comm, "", "Landau Damping and Two Stream options", "DMSWARM");
176   PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex2.c", options->error, &options->error, NULL));
177   PetscCall(PetscOptionsBool("-monitor_efield", "Flag to show efield plot", "ex2.c", options->efield_monitor, &options->efield_monitor, NULL));
178   PetscCall(PetscOptionsBool("-monitor_ics", "Flag to show initial condition histograms", "ex2.c", options->initial_monitor, &options->initial_monitor, NULL));
179   PetscCall(PetscOptionsBool("-monitor_positions", "The flag to show particle positions", "ex2.c", options->monitor_positions, &options->monitor_positions, NULL));
180   PetscCall(PetscOptionsBool("-monitor_poisson", "The flag to show charges, Efield and potential solve", "ex2.c", options->poisson_monitor, &options->poisson_monitor, NULL));
181   PetscCall(PetscOptionsBool("-fake_1D", "Flag to run a 1D simulation (but really in 2D)", "ex2.c", options->fake_1D, &options->fake_1D, NULL));
182   PetscCall(PetscOptionsBool("-twostream", "Run two stream instability", "ex2.c", options->twostream, &options->twostream, NULL));
183   PetscCall(PetscOptionsBool("-perturbed_weights", "Flag to run uniform sampling with perturbed weights", "ex2.c", options->perturbed_weights, &options->perturbed_weights, NULL));
184   PetscCall(PetscOptionsBool("-check_weights", "Ensure all particle weights are positive", "ex2.c", options->checkweights, &options->checkweights, NULL));
185   PetscCall(PetscOptionsInt("-check_vel_res", "Check particle velocity residuals for nightly tests", "ex2.c", options->checkVRes, &options->checkVRes, NULL));
186   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex2.c", options->ostep, &options->ostep, NULL));
187   PetscCall(PetscOptionsReal("-timeScale", "Nondimensionalizing time scale", "ex2.c", options->timeScale, &options->timeScale, NULL));
188   PetscCall(PetscOptionsReal("-initial_velocity", "Initial velocity of perturbed particle", "ex2.c", options->initVel, &options->initVel, NULL));
189   PetscCall(PetscOptionsReal("-total_weight", "Total weight of all particles", "ex2.c", options->totalWeight, &options->totalWeight, NULL));
190   PetscCall(PetscOptionsRealArray("-cosine_coefficients", "Amplitude and frequency of cosine equation used in initialization", "ex2.c", options->cosine_coefficients, &d, NULL));
191   PetscCall(PetscOptionsRealArray("-charges", "Species charges", "ex2.c", options->charges, &maxSpecies, NULL));
192   PetscCall(PetscOptionsEnum("-em_type", "Type of electrostatic solver", "ex2.c", EMTypes, (PetscEnum)options->em, (PetscEnum *)&options->em, NULL));
193   PetscOptionsEnd();
194 
195   PetscCall(PetscLogEventRegister("RhsX", TS_CLASSID, &options->RhsXEvent));
196   PetscCall(PetscLogEventRegister("RhsV", TS_CLASSID, &options->RhsVEvent));
197   PetscCall(PetscLogEventRegister("ESolve", TS_CLASSID, &options->ESolveEvent));
198   PetscCall(PetscLogEventRegister("ETab", TS_CLASSID, &options->ETabEvent));
199   PetscFunctionReturn(PETSC_SUCCESS);
200 }
201 
202 static PetscErrorCode SetupContext(DM dm, DM sw, AppCtx *user)
203 {
204   PetscFunctionBeginUser;
205   if (user->efield_monitor) {
206     PetscDrawAxis axis_ef;
207     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_efield", 0, 300, 400, 300, &user->drawef));
208     PetscCall(PetscDrawSetSave(user->drawef, "ex9_Efield.png"));
209     PetscCall(PetscDrawSetFromOptions(user->drawef));
210     PetscCall(PetscDrawLGCreate(user->drawef, 1, &user->drawlg_ef));
211     PetscCall(PetscDrawLGGetAxis(user->drawlg_ef, &axis_ef));
212     PetscCall(PetscDrawAxisSetLabels(axis_ef, "Electron Electric Field", "time", "E_max"));
213     PetscCall(PetscDrawLGSetLimits(user->drawlg_ef, 0., user->steps * user->stepSize, -10., 0.));
214     PetscCall(PetscDrawAxisSetLimits(axis_ef, 0., user->steps * user->stepSize, -10., 0.));
215   }
216   if (user->initial_monitor) {
217     PetscDrawAxis axis1, axis2, axis3;
218     PetscReal     dmboxlower[2], dmboxupper[2];
219     PetscInt      dim, cStart, cEnd;
220     PetscCall(DMGetDimension(sw, &dim));
221     PetscCall(DMGetBoundingBox(dm, dmboxlower, dmboxupper));
222     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
223 
224     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_initial_conditions_x", 0, 300, 400, 300, &user->drawic_x));
225     PetscCall(PetscDrawSetSave(user->drawic_x, "ex9_ic_x.png"));
226     PetscCall(PetscDrawSetFromOptions(user->drawic_x));
227     PetscCall(PetscDrawHGCreate(user->drawic_x, (int)dim, &user->drawhgic_x));
228     PetscCall(PetscDrawHGGetAxis(user->drawhgic_x, &axis1));
229     PetscCall(PetscDrawHGSetNumberBins(user->drawhgic_x, (int)(cEnd - cStart)));
230     PetscCall(PetscDrawAxisSetLabels(axis1, "Initial X Distribution", "X", "counts"));
231     PetscCall(PetscDrawAxisSetLimits(axis1, dmboxlower[0], dmboxupper[0], 0, 1500));
232 
233     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_initial_conditions_v", 400, 300, 400, 300, &user->drawic_v));
234     PetscCall(PetscDrawSetSave(user->drawic_v, "ex9_ic_v.png"));
235     PetscCall(PetscDrawSetFromOptions(user->drawic_v));
236     PetscCall(PetscDrawHGCreate(user->drawic_v, (int)dim, &user->drawhgic_v));
237     PetscCall(PetscDrawHGGetAxis(user->drawhgic_v, &axis2));
238     PetscCall(PetscDrawHGSetNumberBins(user->drawhgic_v, 1000));
239     PetscCall(PetscDrawAxisSetLabels(axis2, "Initial V_x Distribution", "V", "counts"));
240     PetscCall(PetscDrawAxisSetLimits(axis2, -1, 1, 0, 1500));
241 
242     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_initial_conditions_w", 800, 300, 400, 300, &user->drawic_w));
243     PetscCall(PetscDrawSetSave(user->drawic_w, "ex9_ic_w.png"));
244     PetscCall(PetscDrawSetFromOptions(user->drawic_w));
245     PetscCall(PetscDrawHGCreate(user->drawic_w, (int)dim, &user->drawhgic_w));
246     PetscCall(PetscDrawHGGetAxis(user->drawhgic_w, &axis3));
247     PetscCall(PetscDrawHGSetNumberBins(user->drawhgic_w, 10));
248     PetscCall(PetscDrawAxisSetLabels(axis3, "Initial W Distribution", "weight", "counts"));
249     PetscCall(PetscDrawAxisSetLimits(axis3, 0, 0.01, 0, 5000));
250   }
251   if (user->monitor_positions) {
252     PetscDrawAxis axis;
253 
254     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "position_monitor_species1", 0, 0, 400, 300, &user->positionDraw));
255     PetscCall(PetscDrawSetFromOptions(user->positionDraw));
256     PetscCall(PetscDrawSPCreate(user->positionDraw, 10, &user->positionDrawSP));
257     PetscCall(PetscDrawSPSetDimension(user->positionDrawSP, 1));
258     PetscCall(PetscDrawSPGetAxis(user->positionDrawSP, &axis));
259     PetscCall(PetscDrawSPReset(user->positionDrawSP));
260     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "x", "v"));
261     PetscCall(PetscDrawSetSave(user->positionDraw, "ex9_pos.png"));
262   }
263   if (user->poisson_monitor) {
264     PetscDrawAxis axis_E, axis_Rho, axis_Pot;
265 
266     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "Efield_monitor", 0, 0, 400, 300, &user->EDraw));
267     PetscCall(PetscDrawSetFromOptions(user->EDraw));
268     PetscCall(PetscDrawSPCreate(user->EDraw, 10, &user->EDrawSP));
269     PetscCall(PetscDrawSPSetDimension(user->EDrawSP, 1));
270     PetscCall(PetscDrawSPGetAxis(user->EDrawSP, &axis_E));
271     PetscCall(PetscDrawSPReset(user->EDrawSP));
272     PetscCall(PetscDrawAxisSetLabels(axis_E, "Particles", "x", "E"));
273     PetscCall(PetscDrawSetSave(user->EDraw, "ex9_E_spatial.png"));
274 
275     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "rho_monitor", 0, 0, 400, 300, &user->RhoDraw));
276     PetscCall(PetscDrawSetFromOptions(user->RhoDraw));
277     PetscCall(PetscDrawSPCreate(user->RhoDraw, 10, &user->RhoDrawSP));
278     PetscCall(PetscDrawSPSetDimension(user->RhoDrawSP, 1));
279     PetscCall(PetscDrawSPGetAxis(user->RhoDrawSP, &axis_Rho));
280     PetscCall(PetscDrawSPReset(user->RhoDrawSP));
281     PetscCall(PetscDrawAxisSetLabels(axis_Rho, "Particles", "x", "rho"));
282     PetscCall(PetscDrawSetSave(user->RhoDraw, "ex9_rho_spatial.png"));
283 
284     PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "potential_monitor", 0, 0, 400, 300, &user->PotDraw));
285     PetscCall(PetscDrawSetFromOptions(user->PotDraw));
286     PetscCall(PetscDrawSPCreate(user->PotDraw, 10, &user->PotDrawSP));
287     PetscCall(PetscDrawSPSetDimension(user->PotDrawSP, 1));
288     PetscCall(PetscDrawSPGetAxis(user->PotDrawSP, &axis_Pot));
289     PetscCall(PetscDrawSPReset(user->PotDrawSP));
290     PetscCall(PetscDrawAxisSetLabels(axis_Pot, "Particles", "x", "potential"));
291     PetscCall(PetscDrawSetSave(user->PotDraw, "ex9_phi_spatial.png"));
292   }
293   PetscFunctionReturn(PETSC_SUCCESS);
294 }
295 
296 static PetscErrorCode DestroyContext(AppCtx *user)
297 {
298   PetscFunctionBeginUser;
299   PetscCall(PetscDrawLGDestroy(&user->drawlg_ef));
300   PetscCall(PetscDrawDestroy(&user->drawef));
301   PetscCall(PetscDrawHGDestroy(&user->drawhgic_x));
302   PetscCall(PetscDrawDestroy(&user->drawic_x));
303   PetscCall(PetscDrawHGDestroy(&user->drawhgic_v));
304   PetscCall(PetscDrawDestroy(&user->drawic_v));
305   PetscCall(PetscDrawHGDestroy(&user->drawhgic_w));
306   PetscCall(PetscDrawDestroy(&user->drawic_w));
307   PetscCall(PetscDrawSPDestroy(&user->positionDrawSP));
308   PetscCall(PetscDrawDestroy(&user->positionDraw));
309 
310   PetscCall(PetscDrawSPDestroy(&user->EDrawSP));
311   PetscCall(PetscDrawDestroy(&user->EDraw));
312   PetscCall(PetscDrawSPDestroy(&user->RhoDrawSP));
313   PetscCall(PetscDrawDestroy(&user->RhoDraw));
314   PetscCall(PetscDrawSPDestroy(&user->PotDrawSP));
315   PetscCall(PetscDrawDestroy(&user->PotDraw));
316 
317   PetscCall(PetscBagDestroy(&user->bag));
318   PetscFunctionReturn(PETSC_SUCCESS);
319 }
320 
321 static PetscErrorCode CheckNonNegativeWeights(DM sw, AppCtx *user)
322 {
323   const PetscScalar *w;
324   PetscInt           Np;
325 
326   PetscFunctionBeginUser;
327   if (!user->checkweights) PetscFunctionReturn(PETSC_SUCCESS);
328   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w));
329   PetscCall(DMSwarmGetLocalSize(sw, &Np));
330   for (PetscInt p = 0; p < Np; ++p) PetscCheck(w[p] >= 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle %" PetscInt_FMT " has negative weight %g", p, w[p]);
331   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
332   PetscFunctionReturn(PETSC_SUCCESS);
333 }
334 
335 static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user)
336 {
337   DM                 dm;
338   const PetscReal   *coords;
339   const PetscScalar *w;
340   PetscReal          mom[3] = {0.0, 0.0, 0.0};
341   PetscInt           cell, cStart, cEnd, dim;
342 
343   PetscFunctionBeginUser;
344   PetscCall(DMGetDimension(sw, &dim));
345   PetscCall(DMSwarmGetCellDM(sw, &dm));
346   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
347   PetscCall(DMSwarmSortGetAccess(sw));
348   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&coords));
349   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w));
350   for (cell = cStart; cell < cEnd; ++cell) {
351     PetscInt *pidx;
352     PetscInt  Np, p, d;
353 
354     PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
355     for (p = 0; p < Np; ++p) {
356       const PetscInt   idx = pidx[p];
357       const PetscReal *c   = &coords[idx * dim];
358 
359       mom[0] += PetscRealPart(w[idx]);
360       mom[1] += PetscRealPart(w[idx]) * c[0];
361       for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d];
362       //if (w[idx] < 0. ) PetscPrintf(PETSC_COMM_WORLD, "error, negative weight %" PetscInt_FMT " \n", idx);
363     }
364     PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Np, &pidx));
365   }
366   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&coords));
367   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
368   PetscCall(DMSwarmSortRestoreAccess(sw));
369   PetscCallMPI(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
370   PetscFunctionReturn(PETSC_SUCCESS);
371 }
372 
373 static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
374 {
375   f0[0] = u[0];
376 }
377 
378 static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
379 {
380   f0[0] = x[0] * u[0];
381 }
382 
383 static void f0_r2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
384 {
385   PetscInt d;
386 
387   f0[0] = 0.0;
388   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0];
389 }
390 
391 static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
392 {
393   PetscDS     prob;
394   PetscScalar mom;
395   PetscInt    field = 0;
396 
397   PetscFunctionBeginUser;
398   PetscCall(DMGetDS(dm, &prob));
399   PetscCall(PetscDSSetObjective(prob, field, &f0_1));
400   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
401   moments[0] = PetscRealPart(mom);
402   PetscCall(PetscDSSetObjective(prob, field, &f0_x));
403   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
404   moments[1] = PetscRealPart(mom);
405   PetscCall(PetscDSSetObjective(prob, field, &f0_r2));
406   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
407   moments[2] = PetscRealPart(mom);
408   PetscFunctionReturn(PETSC_SUCCESS);
409 }
410 
411 static PetscErrorCode MonitorEField(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
412 {
413   AppCtx    *user = (AppCtx *)ctx;
414   DM         dm, sw;
415   PetscReal *E;
416   PetscReal  Enorm = 0., lgEnorm, lgEmax, sum = 0., Emax = 0., temp = 0., *weight, chargesum = 0.;
417   PetscReal *x, *v;
418   PetscInt  *species, dim, p, d, Np, cStart, cEnd;
419   PetscReal  pmoments[3]; /* \int f, \int x f, \int r^2 f */
420   PetscReal  fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */
421   Vec        rho;
422 
423   PetscFunctionBeginUser;
424   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
425   PetscCall(TSGetDM(ts, &sw));
426   PetscCall(DMSwarmGetCellDM(sw, &dm));
427   PetscCall(DMGetDimension(sw, &dim));
428   PetscCall(DMSwarmGetLocalSize(sw, &Np));
429   PetscCall(DMSwarmSortGetAccess(sw));
430   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
431   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v));
432   PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
433   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
434   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
435   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
436 
437   for (p = 0; p < Np; ++p) {
438     for (d = 0; d < 1; ++d) {
439       temp = PetscAbsReal(E[p * dim + d]);
440       if (temp > Emax) Emax = temp;
441     }
442     Enorm += PetscSqrtReal(E[p * dim] * E[p * dim]);
443     sum += E[p * dim];
444     chargesum += user->charges[0] * weight[p];
445   }
446   lgEnorm = Enorm != 0 ? PetscLog10Real(Enorm) : -16.;
447   lgEmax  = Emax != 0 ? PetscLog10Real(Emax) : -16.;
448 
449   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
450   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v));
451   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
452   PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
453   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
454 
455   Parameter *param;
456   PetscCall(PetscBagGetData(user->bag, (void **)&param));
457   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "charges", &rho));
458   if (user->em == EM_PRIMAL) {
459     PetscCall(computeParticleMoments(sw, pmoments, user));
460     PetscCall(computeFEMMoments(dm, rho, fmoments, user));
461   } else if (user->em == EM_MIXED) {
462     DM       potential_dm;
463     IS       potential_IS;
464     PetscInt fields = 1;
465     PetscCall(DMCreateSubDM(dm, 1, &fields, &potential_IS, &potential_dm));
466 
467     PetscCall(computeParticleMoments(sw, pmoments, user));
468     PetscCall(computeFEMMoments(potential_dm, rho, fmoments, user));
469     PetscCall(DMDestroy(&potential_dm));
470     PetscCall(ISDestroy(&potential_IS));
471   }
472   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "charges", &rho));
473 
474   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%f\t%+e\t%e\t%f\t%f\t%f\t%f\t%f\t%f\t%f\t%f\t%f\t%f\n", (double)t, (double)sum, (double)Enorm, (double)lgEnorm, (double)Emax, (double)lgEmax, (double)chargesum, (double)pmoments[0], (double)pmoments[1], (double)pmoments[2], (double)fmoments[0], (double)fmoments[1], (double)fmoments[2]));
475   PetscCall(PetscDrawLGAddPoint(user->drawlg_ef, &t, &lgEmax));
476   PetscCall(PetscDrawLGDraw(user->drawlg_ef));
477   PetscCall(PetscDrawSave(user->drawef));
478   PetscFunctionReturn(PETSC_SUCCESS);
479 }
480 
481 PetscErrorCode MonitorInitialConditions(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
482 {
483   AppCtx            *user = (AppCtx *)ctx;
484   DM                 dm, sw;
485   const PetscScalar *u;
486   PetscReal         *weight, *pos, *vel;
487   PetscInt           dim, p, Np, cStart, cEnd;
488 
489   PetscFunctionBegin;
490   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */
491   PetscCall(TSGetDM(ts, &sw));
492   PetscCall(DMSwarmGetCellDM(sw, &dm));
493   PetscCall(DMGetDimension(sw, &dim));
494   PetscCall(DMSwarmGetLocalSize(sw, &Np));
495   PetscCall(DMSwarmSortGetAccess(sw));
496   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
497 
498   if (step == 0) {
499     PetscCall(PetscDrawHGReset(user->drawhgic_x));
500     PetscCall(PetscDrawHGGetDraw(user->drawhgic_x, &user->drawic_x));
501     PetscCall(PetscDrawClear(user->drawic_x));
502     PetscCall(PetscDrawFlush(user->drawic_x));
503 
504     PetscCall(PetscDrawHGReset(user->drawhgic_v));
505     PetscCall(PetscDrawHGGetDraw(user->drawhgic_v, &user->drawic_v));
506     PetscCall(PetscDrawClear(user->drawic_v));
507     PetscCall(PetscDrawFlush(user->drawic_v));
508 
509     PetscCall(PetscDrawHGReset(user->drawhgic_w));
510     PetscCall(PetscDrawHGGetDraw(user->drawhgic_w, &user->drawic_w));
511     PetscCall(PetscDrawClear(user->drawic_w));
512     PetscCall(PetscDrawFlush(user->drawic_w));
513 
514     PetscCall(VecGetArrayRead(U, &u));
515     PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&vel));
516     PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
517     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&pos));
518 
519     PetscCall(VecGetLocalSize(U, &Np));
520     Np /= dim * 2;
521     for (p = 0; p < Np; ++p) {
522       PetscCall(PetscDrawHGAddValue(user->drawhgic_x, pos[p * dim]));
523       PetscCall(PetscDrawHGAddValue(user->drawhgic_v, vel[p * dim]));
524       PetscCall(PetscDrawHGAddValue(user->drawhgic_w, weight[p]));
525     }
526 
527     PetscCall(VecRestoreArrayRead(U, &u));
528     PetscCall(PetscDrawHGDraw(user->drawhgic_x));
529     PetscCall(PetscDrawHGSave(user->drawhgic_x));
530 
531     PetscCall(PetscDrawHGDraw(user->drawhgic_v));
532     PetscCall(PetscDrawHGSave(user->drawhgic_v));
533 
534     PetscCall(PetscDrawHGDraw(user->drawhgic_w));
535     PetscCall(PetscDrawHGSave(user->drawhgic_w));
536 
537     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&pos));
538     PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&vel));
539     PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
540   }
541   PetscFunctionReturn(PETSC_SUCCESS);
542 }
543 
544 static PetscErrorCode MonitorPositions_2D(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
545 {
546   AppCtx         *user = (AppCtx *)ctx;
547   DM              dm, sw;
548   PetscScalar    *x, *v, *weight;
549   PetscReal       lower[3], upper[3], speed;
550   const PetscInt *s;
551   PetscInt        dim, cStart, cEnd, c;
552 
553   PetscFunctionBeginUser;
554   if (step > 0 && step % user->ostep == 0) {
555     PetscCall(TSGetDM(ts, &sw));
556     PetscCall(DMSwarmGetCellDM(sw, &dm));
557     PetscCall(DMGetDimension(dm, &dim));
558     PetscCall(DMGetBoundingBox(dm, lower, upper));
559     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
560     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
561     PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v));
562     PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
563     PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&s));
564     PetscCall(DMSwarmSortGetAccess(sw));
565     PetscCall(PetscDrawSPReset(user->positionDrawSP));
566     PetscCall(PetscDrawSPSetLimits(user->positionDrawSP, lower[0], upper[0], lower[1], upper[1]));
567     PetscCall(PetscDrawSPSetLimits(user->positionDrawSP, lower[0], upper[0], -12, 12));
568     for (c = 0; c < cEnd - cStart; ++c) {
569       PetscInt *pidx, Npc, q;
570       PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx));
571       for (q = 0; q < Npc; ++q) {
572         const PetscInt p = pidx[q];
573         if (s[p] == 0) {
574           speed = PetscSqrtReal(PetscSqr(v[p * dim]) + PetscSqr(v[p * dim + 1]));
575           if (dim == 1 || user->fake_1D) {
576             PetscCall(PetscDrawSPAddPointColorized(user->positionDrawSP, &x[p * dim], &v[p * dim], &speed));
577           } else {
578             PetscCall(PetscDrawSPAddPointColorized(user->positionDrawSP, &x[p * dim], &x[p * dim + 1], &speed));
579           }
580         } else if (s[p] == 1) {
581           PetscCall(PetscDrawSPAddPoint(user->positionDrawSP, &x[p * dim], &v[p * dim]));
582         }
583       }
584       PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx));
585     }
586     PetscCall(PetscDrawSPDraw(user->positionDrawSP, PETSC_TRUE));
587     PetscCall(PetscDrawSave(user->positionDraw));
588     PetscCall(DMSwarmSortRestoreAccess(sw));
589     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
590     PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
591     PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v));
592     PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&s));
593   }
594   PetscFunctionReturn(PETSC_SUCCESS);
595 }
596 
597 static PetscErrorCode MonitorPoisson(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
598 {
599   AppCtx      *user = (AppCtx *)ctx;
600   DM           dm, sw;
601   PetscScalar *x, *E, *weight, *pot, *charges;
602   PetscReal    lower[3], upper[3], xval;
603   PetscInt     dim, cStart, cEnd, c;
604 
605   PetscFunctionBeginUser;
606   if (step > 0 && step % user->ostep == 0) {
607     PetscCall(TSGetDM(ts, &sw));
608     PetscCall(DMSwarmGetCellDM(sw, &dm));
609     PetscCall(DMGetDimension(dm, &dim));
610     PetscCall(DMGetBoundingBox(dm, lower, upper));
611     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
612 
613     PetscCall(PetscDrawSPReset(user->RhoDrawSP));
614     PetscCall(PetscDrawSPReset(user->EDrawSP));
615     PetscCall(PetscDrawSPReset(user->PotDrawSP));
616     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
617     PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
618     PetscCall(DMSwarmGetField(sw, "potential", NULL, NULL, (void **)&pot));
619     PetscCall(DMSwarmGetField(sw, "charges", NULL, NULL, (void **)&charges));
620     PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
621 
622     PetscCall(DMSwarmSortGetAccess(sw));
623     for (c = 0; c < cEnd - cStart; ++c) {
624       PetscReal Esum = 0.0;
625       PetscInt *pidx, Npc, q;
626       PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx));
627       for (q = 0; q < Npc; ++q) {
628         const PetscInt p = pidx[q];
629         Esum += E[p * dim];
630       }
631       xval = (c + 0.5) * ((upper[0] - lower[0]) / (cEnd - cStart));
632       PetscCall(PetscDrawSPAddPoint(user->EDrawSP, &xval, &Esum));
633       PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx));
634     }
635     for (c = 0; c < (cEnd - cStart); ++c) {
636       xval = (c + 0.5) * ((upper[0] - lower[0]) / (cEnd - cStart));
637       PetscCall(PetscDrawSPAddPoint(user->RhoDrawSP, &xval, &charges[c]));
638       PetscCall(PetscDrawSPAddPoint(user->PotDrawSP, &xval, &pot[c]));
639     }
640     PetscCall(PetscDrawSPDraw(user->RhoDrawSP, PETSC_TRUE));
641     PetscCall(PetscDrawSave(user->RhoDraw));
642     PetscCall(PetscDrawSPDraw(user->EDrawSP, PETSC_TRUE));
643     PetscCall(PetscDrawSave(user->EDraw));
644     PetscCall(PetscDrawSPDraw(user->PotDrawSP, PETSC_TRUE));
645     PetscCall(PetscDrawSave(user->PotDraw));
646     PetscCall(DMSwarmSortRestoreAccess(sw));
647     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
648     PetscCall(DMSwarmRestoreField(sw, "potential", NULL, NULL, (void **)&pot));
649     PetscCall(DMSwarmRestoreField(sw, "charges", NULL, NULL, (void **)&charges));
650     PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
651     PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
652   }
653   PetscFunctionReturn(PETSC_SUCCESS);
654 }
655 
656 static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
657 {
658   PetscBag   bag;
659   Parameter *p;
660 
661   PetscFunctionBeginUser;
662   /* setup PETSc parameter bag */
663   PetscCall(PetscBagGetData(ctx->bag, (void **)&p));
664   PetscCall(PetscBagSetName(ctx->bag, "par", "Vlasov-Poisson Parameters"));
665   bag = ctx->bag;
666   PetscCall(PetscBagRegisterScalar(bag, &p->v0, 1.0, "v0", "Velocity scale, m/s"));
667   PetscCall(PetscBagRegisterScalar(bag, &p->t0, 1.0, "t0", "Time scale, s"));
668   PetscCall(PetscBagRegisterScalar(bag, &p->x0, 1.0, "x0", "Space scale, m"));
669   PetscCall(PetscBagRegisterScalar(bag, &p->v0, 1.0, "phi0", "Potential scale, kg*m^2/A*s^3"));
670   PetscCall(PetscBagRegisterScalar(bag, &p->q0, 1.0, "q0", "Charge Scale, A*s"));
671   PetscCall(PetscBagRegisterScalar(bag, &p->m0, 1.0, "m0", "Mass Scale, kg"));
672   PetscCall(PetscBagRegisterScalar(bag, &p->epsi0, 1.0, "epsi0", "Permittivity of Free Space, kg"));
673   PetscCall(PetscBagRegisterScalar(bag, &p->kb, 1.0, "kb", "Boltzmann Constant, m^2 kg/s^2 K^1"));
674 
675   PetscCall(PetscBagRegisterScalar(bag, &p->sigma, 1.0, "sigma", "Charge per unit area, C/m^3"));
676   PetscCall(PetscBagRegisterScalar(bag, &p->poissonNumber, 1.0, "poissonNumber", "Non-Dimensional Poisson Number"));
677   PetscCall(PetscBagRegisterScalar(bag, &p->vlasovNumber, 1.0, "vlasovNumber", "Non-Dimensional Vlasov Number"));
678   PetscCall(PetscBagSetFromOptions(bag));
679   {
680     PetscViewer       viewer;
681     PetscViewerFormat format;
682     PetscBool         flg;
683 
684     PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
685     if (flg) {
686       PetscCall(PetscViewerPushFormat(viewer, format));
687       PetscCall(PetscBagView(bag, viewer));
688       PetscCall(PetscViewerFlush(viewer));
689       PetscCall(PetscViewerPopFormat(viewer));
690       PetscCall(PetscViewerDestroy(&viewer));
691     }
692   }
693   PetscFunctionReturn(PETSC_SUCCESS);
694 }
695 
696 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
697 {
698   PetscFunctionBeginUser;
699   PetscCall(DMCreate(comm, dm));
700   PetscCall(DMSetType(*dm, DMPLEX));
701   PetscCall(DMSetFromOptions(*dm));
702   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
703   PetscFunctionReturn(PETSC_SUCCESS);
704 }
705 
706 static void ion_f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
707 {
708   f0[0] = -constants[SIGMA];
709 }
710 
711 static void laplacian_f1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
712 {
713   PetscInt d;
714   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
715 }
716 
717 static void laplacian_g3(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
718 {
719   PetscInt d;
720   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
721 }
722 
723 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
724 {
725   *u = 0.0;
726   return PETSC_SUCCESS;
727 }
728 
729 /*
730    /  I   -grad\ / q \ = /0\
731    \-div    0  / \phi/   \f/
732 */
733 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
734 {
735   for (PetscInt d = 0; d < dim; ++d) f0[d] += u[uOff[0] + d];
736 }
737 
738 static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
739 {
740   for (PetscInt d = 0; d < dim; ++d) f1[d * dim + d] = u[uOff[1]];
741 }
742 
743 static void f0_phi_backgroundCharge(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
744 {
745   f0[0] += constants[SIGMA];
746   for (PetscInt d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d];
747 }
748 
749 /* Boundary residual. Dirichlet boundary for u means u_bdy=p*n */
750 static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
751 {
752   for (PetscInt d = 0; d < dim; ++d) g0[d * dim + d] = 1.0;
753 }
754 
755 static void g2_qphi(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
756 {
757   for (PetscInt d = 0; d < dim; ++d) g2[d * dim + d] = 1.0;
758 }
759 
760 static void g1_phiq(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
761 {
762   for (PetscInt d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
763 }
764 
765 static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
766 {
767   PetscFE   fephi, feq;
768   PetscDS   ds;
769   PetscBool simplex;
770   PetscInt  dim;
771 
772   PetscFunctionBeginUser;
773   PetscCall(DMGetDimension(dm, &dim));
774   PetscCall(DMPlexIsSimplex(dm, &simplex));
775   if (user->em == EM_MIXED) {
776     DMLabel        label;
777     const PetscInt id = 1;
778 
779     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", PETSC_DETERMINE, &feq));
780     PetscCall(PetscObjectSetName((PetscObject)feq, "field"));
781     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", PETSC_DETERMINE, &fephi));
782     PetscCall(PetscObjectSetName((PetscObject)fephi, "potential"));
783     PetscCall(PetscFECopyQuadrature(feq, fephi));
784     PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq));
785     PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fephi));
786     PetscCall(DMCreateDS(dm));
787     PetscCall(PetscFEDestroy(&fephi));
788     PetscCall(PetscFEDestroy(&feq));
789 
790     PetscCall(DMGetLabel(dm, "marker", &label));
791     PetscCall(DMGetDS(dm, &ds));
792 
793     PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q));
794     PetscCall(PetscDSSetResidual(ds, 1, f0_phi_backgroundCharge, NULL));
795     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL));
796     PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qphi, NULL));
797     PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_phiq, NULL, NULL));
798 
799     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, NULL, NULL));
800 
801   } else if (user->em == EM_PRIMAL) {
802     MatNullSpace nullsp;
803     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fephi));
804     PetscCall(PetscObjectSetName((PetscObject)fephi, "potential"));
805     PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fephi));
806     PetscCall(DMCreateDS(dm));
807     PetscCall(DMGetDS(dm, &ds));
808     PetscCall(PetscDSSetResidual(ds, 0, ion_f0, laplacian_f1));
809     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, laplacian_g3));
810     PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullsp));
811     PetscCall(PetscObjectCompose((PetscObject)fephi, "nullspace", (PetscObject)nullsp));
812     PetscCall(MatNullSpaceDestroy(&nullsp));
813     PetscCall(PetscFEDestroy(&fephi));
814   }
815   PetscFunctionReturn(PETSC_SUCCESS);
816 }
817 
818 static PetscErrorCode CreatePoisson(DM dm, AppCtx *user)
819 {
820   SNES         snes;
821   Mat          J;
822   MatNullSpace nullSpace;
823 
824   PetscFunctionBeginUser;
825   PetscCall(CreateFEM(dm, user));
826   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes));
827   PetscCall(SNESSetOptionsPrefix(snes, "em_"));
828   PetscCall(SNESSetDM(snes, dm));
829   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, user));
830   PetscCall(SNESSetFromOptions(snes));
831 
832   PetscCall(DMCreateMatrix(dm, &J));
833   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullSpace));
834   PetscCall(MatSetNullSpace(J, nullSpace));
835   PetscCall(MatNullSpaceDestroy(&nullSpace));
836   PetscCall(SNESSetJacobian(snes, J, J, NULL, NULL));
837   PetscCall(MatDestroy(&J));
838   user->snes = snes;
839   PetscFunctionReturn(PETSC_SUCCESS);
840 }
841 
842 PetscErrorCode PetscPDFPertubedConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
843 {
844   p[0] = (1 + 0.01 * PetscCosReal(0.5 * x[0])) / (2 * PETSC_PI);
845   p[1] = (1 + 0.01 * PetscCosReal(0.5 * x[1])) / (2 * PETSC_PI);
846   return PETSC_SUCCESS;
847 }
848 PetscErrorCode PetscPDFPertubedConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
849 {
850   p[0] = (1. + 0.01 * PetscCosReal(0.5 * x[0])) / (2 * PETSC_PI);
851   return PETSC_SUCCESS;
852 }
853 
854 PetscErrorCode PetscPDFCosine1D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
855 {
856   const PetscReal alpha = scale ? scale[0] : 0.0;
857   const PetscReal k     = scale ? scale[1] : 1.;
858   p[0]                  = (1 + alpha * PetscCosReal(k * x[0]));
859   return PETSC_SUCCESS;
860 }
861 
862 PetscErrorCode PetscPDFCosine2D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
863 {
864   const PetscReal alpha = scale ? scale[0] : 0.;
865   const PetscReal k     = scale ? scale[0] : 1.;
866   p[0]                  = (1 + alpha * PetscCosReal(k * (x[0] + x[1])));
867   return PETSC_SUCCESS;
868 }
869 
870 static PetscErrorCode CreateVelocityDM(DM sw, DM *vdm)
871 {
872   PetscFunctionBegin;
873   PetscCall(DMCreate(PETSC_COMM_SELF, vdm));
874   PetscCall(DMSetType(*vdm, DMPLEX));
875   PetscCall(DMPlexSetOptionsPrefix(*vdm, "v"));
876   PetscCall(DMSetFromOptions(*vdm));
877   PetscCall(DMViewFromOptions(*vdm, NULL, "-dm_view"));
878   PetscFunctionReturn(PETSC_SUCCESS);
879 }
880 
881 static PetscErrorCode InitializeParticles_Centroid(DM sw, PetscBool force1D)
882 {
883   DM_Swarm  *swarm = (DM_Swarm *)sw->data;
884   DM         xdm, vdm;
885   PetscReal  vmin[3], vmax[3];
886   PetscReal *x, *v;
887   PetscInt  *species, *cellid;
888   PetscInt   dim, xcStart, xcEnd, vcStart, vcEnd, Ns, Np, Npc, debug;
889   PetscBool  flg;
890   MPI_Comm   comm;
891 
892   PetscFunctionBegin;
893   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
894 
895   PetscOptionsBegin(comm, "", "DMSwarm Options", "DMSWARM");
896   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
897   PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg));
898   if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns));
899   PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_coords", "Debug output level for particle coordinate computations", "InitializeParticles", 0, &swarm->printCoords, NULL, 0));
900   PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_weights", "Debug output level for particle weight computations", "InitializeWeights", 0, &swarm->printWeights, NULL, 0));
901   PetscOptionsEnd();
902   debug = swarm->printCoords;
903 
904   PetscCall(DMGetDimension(sw, &dim));
905   PetscCall(DMSwarmGetCellDM(sw, &xdm));
906   PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd));
907 
908   PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm));
909   if (!vdm) {
910     PetscCall(CreateVelocityDM(sw, &vdm));
911     PetscCall(PetscObjectCompose((PetscObject)sw, "__vdm__", (PetscObject)vdm));
912     PetscCall(DMDestroy(&vdm));
913     PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm));
914   }
915   PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd));
916 
917   // One particle per centroid on the tensor product grid
918   Npc = (vcEnd - vcStart) * Ns;
919   Np  = (xcEnd - xcStart) * Npc;
920   PetscCall(DMSwarmSetLocalSizes(sw, Np, 0));
921   if (debug) {
922     PetscInt gNp;
923     PetscCallMPI(MPIU_Allreduce(&Np, &gNp, 1, MPIU_INT, MPIU_SUM, comm));
924     PetscCall(PetscPrintf(comm, "Global Np = %" PetscInt_FMT "\n", gNp));
925   }
926 
927   // Set species and cellid
928   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
929   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
930   for (PetscInt c = 0, p = 0; c < xcEnd - xcStart; ++c) {
931     for (PetscInt s = 0; s < Ns; ++s) {
932       for (PetscInt q = 0; q < Npc / Ns; ++q, ++p) {
933         species[p] = s;
934         cellid[p]  = c;
935       }
936     }
937   }
938   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
939   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
940 
941   // Set particle coordinates
942   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
943   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v));
944   PetscCall(DMSwarmSortGetAccess(sw));
945   PetscCall(DMGetBoundingBox(vdm, vmin, vmax));
946   PetscCall(DMGetCoordinatesLocalSetUp(xdm));
947   for (PetscInt c = 0; c < xcEnd - xcStart; ++c) {
948     const PetscInt cell = c + xcStart;
949     PetscInt      *pidx, Npc;
950     PetscReal      centroid[3], volume;
951 
952     PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx));
953     PetscCall(DMPlexComputeCellGeometryFVM(xdm, cell, &volume, centroid, NULL));
954     for (PetscInt s = 0; s < Ns; ++s) {
955       for (PetscInt q = 0; q < Npc / Ns; ++q) {
956         const PetscInt p = pidx[q * Ns + s];
957 
958         for (PetscInt d = 0; d < dim; ++d) {
959           x[p * dim + d] = centroid[d];
960           v[p * dim + d] = vmin[0] + (q + 0.5) * ((vmax[0] - vmin[0]) / (Npc / Ns));
961           if (force1D && d > 0) v[p * dim + d] = 0.;
962         }
963       }
964     }
965     PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx));
966   }
967   PetscCall(DMSwarmSortRestoreAccess(sw));
968   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
969   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v));
970   PetscFunctionReturn(PETSC_SUCCESS);
971 }
972 
973 /*
974   InitializeWeights - Compute weight for each local particle
975 
976   Input Parameters:
977 + sw          - The `DMSwarm`
978 . totalWeight - The sum of all particle weights
979 . force1D     - Flag to treat the problem as 1D
980 . func        - The PDF for the particle spatial distribution
981 - param       - The PDF parameters
982 
983   Notes:
984   The PDF for velocity is assumed to be a Gaussian
985 
986   The particle weights are returned in the `w_q` field of `sw`.
987 */
988 static PetscErrorCode InitializeWeights(DM sw, PetscReal totalWeight, PetscBool force1D, PetscProbFunc func, const PetscReal param[])
989 {
990   DM               xdm, vdm;
991   PetscScalar     *weight;
992   PetscQuadrature  xquad;
993   const PetscReal *xq, *xwq;
994   const PetscInt   order = 5;
995   PetscReal       *xqd   = NULL, xi0[3];
996   PetscReal        xwtot = 0., pwtot = 0.;
997   PetscInt         xNq;
998   PetscInt         dim, Ns, xcStart, xcEnd, vcStart, vcEnd, debug = ((DM_Swarm *)sw->data)->printWeights;
999   MPI_Comm         comm;
1000   PetscMPIInt      rank;
1001 
1002   PetscFunctionBegin;
1003   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
1004   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1005   PetscCall(DMGetDimension(sw, &dim));
1006   PetscCall(DMSwarmGetCellDM(sw, &xdm));
1007   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
1008   PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd));
1009   PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm));
1010   PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd));
1011 
1012   // Setup Quadrature for spatial and velocity weight calculations
1013   if (force1D) PetscCall(PetscDTGaussTensorQuadrature(1, 1, order, -1.0, 1.0, &xquad));
1014   else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, order, -1.0, 1.0, &xquad));
1015   PetscCall(PetscQuadratureGetData(xquad, NULL, NULL, &xNq, &xq, &xwq));
1016   if (force1D) {
1017     PetscCall(PetscCalloc1(xNq * dim, &xqd));
1018     for (PetscInt q = 0; q < xNq; ++q) xqd[q * dim] = xq[q];
1019   }
1020   for (PetscInt d = 0; d < dim; ++d) xi0[d] = -1.0;
1021 
1022   // Integrate the density function to get the weights of particles in each cell
1023   PetscCall(DMGetCoordinatesLocalSetUp(vdm));
1024   PetscCall(DMSwarmSortGetAccess(sw));
1025   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
1026   for (PetscInt c = xcStart; c < xcEnd; ++c) {
1027     PetscReal          xv0[3], xJ[9], xinvJ[9], xdetJ, xqr[3], xden, xw = 0.;
1028     PetscInt          *pidx, Npc;
1029     PetscInt           xNc;
1030     const PetscScalar *xarray;
1031     PetscScalar       *xcoords = NULL;
1032     PetscBool          xisDG;
1033 
1034     PetscCall(DMPlexGetCellCoordinates(xdm, c, &xisDG, &xNc, &xarray, &xcoords));
1035     PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx));
1036     PetscCheck(Npc == (vcEnd - vcStart) * Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of particles %" PetscInt_FMT " in cell (rank %d) != %" PetscInt_FMT " number of velocity vertices", Npc, rank, (vcEnd - vcStart) * Ns);
1037     PetscCall(DMPlexComputeCellGeometryFEM(xdm, c, NULL, xv0, xJ, xinvJ, &xdetJ));
1038     for (PetscInt q = 0; q < xNq; ++q) {
1039       // Transform quadrature points from ref space to real space
1040       if (force1D) CoordinatesRefToReal(dim, dim, xi0, xv0, xJ, &xqd[q * dim], xqr);
1041       else CoordinatesRefToReal(dim, dim, xi0, xv0, xJ, &xq[q * dim], xqr);
1042       // Get probability density at quad point
1043       //   No need to scale xqr since PDF will be periodic
1044       PetscCall((*func)(xqr, param, &xden));
1045       if (force1D) xdetJ = xJ[0]; // Only want x contribution
1046       xw += xden * (xwq[q] * xdetJ);
1047     }
1048     xwtot += xw;
1049     if (debug) {
1050       IS              globalOrdering;
1051       const PetscInt *ordering;
1052 
1053       PetscCall(DMPlexGetCellNumbering(xdm, &globalOrdering));
1054       PetscCall(ISGetIndices(globalOrdering, &ordering));
1055       PetscCall(PetscSynchronizedPrintf(comm, "c:%" PetscInt_FMT " [x_a,x_b] = %1.15f,%1.15f -> cell weight = %1.15f\n", ordering[c], (double)PetscRealPart(xcoords[0]), (double)PetscRealPart(xcoords[2]), (double)xw));
1056       PetscCall(ISRestoreIndices(globalOrdering, &ordering));
1057     }
1058     // Set weights to be Gaussian in velocity cells
1059     for (PetscInt vc = vcStart; vc < vcEnd; ++vc) {
1060       const PetscInt     p  = pidx[vc * Ns + 0];
1061       PetscReal          vw = 0.;
1062       PetscInt           vNc;
1063       const PetscScalar *varray;
1064       PetscScalar       *vcoords = NULL;
1065       PetscBool          visDG;
1066 
1067       PetscCall(DMPlexGetCellCoordinates(vdm, vc, &visDG, &vNc, &varray, &vcoords));
1068       // TODO: Fix 2 stream Ask Joe
1069       //   Two stream function from 1/2pi v^2 e^(-v^2/2)
1070       //   vw = 1. / (PetscSqrtReal(2 * PETSC_PI)) * (((coords_v[0] * PetscExpReal(-PetscSqr(coords_v[0]) / 2.)) - (coords_v[1] * PetscExpReal(-PetscSqr(coords_v[1]) / 2.)))) - 0.5 * PetscErfReal(coords_v[0] / PetscSqrtReal(2.)) + 0.5 * (PetscErfReal(coords_v[1] / PetscSqrtReal(2.)));
1071       vw = 0.5 * (PetscErfReal(vcoords[1] / PetscSqrtReal(2.)) - PetscErfReal(vcoords[0] / PetscSqrtReal(2.)));
1072 
1073       weight[p] = totalWeight * vw * xw;
1074       pwtot += weight[p];
1075       PetscCheck(weight[p] <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle %" PetscInt_FMT " weight exceeeded 1: %g, %g, %g", p, xw, vw, totalWeight);
1076       PetscCall(DMPlexRestoreCellCoordinates(vdm, vc, &visDG, &vNc, &varray, &vcoords));
1077       if (debug > 1) PetscCall(PetscPrintf(comm, "particle %" PetscInt_FMT ": %g, vw: %g xw: %g\n", p, weight[p], vw, xw));
1078     }
1079     PetscCall(DMPlexRestoreCellCoordinates(xdm, c, &xisDG, &xNc, &xarray, &xcoords));
1080     PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx));
1081   }
1082   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
1083   PetscCall(DMSwarmSortRestoreAccess(sw));
1084   PetscCall(PetscQuadratureDestroy(&xquad));
1085   if (force1D) PetscCall(PetscFree(xqd));
1086 
1087   if (debug) {
1088     PetscReal wtot[2] = {pwtot, xwtot}, gwtot[2];
1089 
1090     PetscCall(PetscSynchronizedFlush(comm, NULL));
1091     PetscCallMPI(MPIU_Allreduce(wtot, gwtot, 2, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD));
1092     PetscCall(PetscPrintf(comm, "particle weight sum = %1.10f cell weight sum = %1.10f\n", (double)gwtot[0], (double)gwtot[1]));
1093   }
1094   PetscFunctionReturn(PETSC_SUCCESS);
1095 }
1096 
1097 static PetscErrorCode InitializeParticles_PerturbedWeights(DM sw, AppCtx *user)
1098 {
1099   PetscReal scale[2] = {user->cosine_coefficients[0], user->cosine_coefficients[1]};
1100 
1101   PetscFunctionBegin;
1102   PetscCall(InitializeParticles_Centroid(sw, user->fake_1D));
1103   PetscCall(InitializeWeights(sw, user->totalWeight, user->fake_1D, user->fake_1D ? PetscPDFCosine1D : PetscPDFCosine2D, scale));
1104   PetscFunctionReturn(PETSC_SUCCESS);
1105 }
1106 
1107 static PetscErrorCode InitializeConstants(DM sw, AppCtx *user)
1108 {
1109   DM         dm;
1110   PetscInt  *species;
1111   PetscReal *weight, totalCharge = 0., totalWeight = 0., gmin[3], gmax[3], global_charge, global_weight;
1112   PetscInt   Np, dim;
1113 
1114   PetscFunctionBegin;
1115   PetscCall(DMSwarmGetCellDM(sw, &dm));
1116   PetscCall(DMGetDimension(sw, &dim));
1117   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1118   PetscCall(DMGetBoundingBox(dm, gmin, gmax));
1119   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
1120   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
1121   for (PetscInt p = 0; p < Np; ++p) {
1122     totalWeight += weight[p];
1123     totalCharge += user->charges[species[p]] * weight[p];
1124   }
1125   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
1126   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1127   {
1128     Parameter *param;
1129     PetscReal  Area;
1130 
1131     PetscCall(PetscBagGetData(user->bag, (void **)&param));
1132     switch (dim) {
1133     case 1:
1134       Area = (gmax[0] - gmin[0]);
1135       break;
1136     case 2:
1137       if (user->fake_1D) {
1138         Area = (gmax[0] - gmin[0]);
1139       } else {
1140         Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]);
1141       }
1142       break;
1143     case 3:
1144       Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]) * (gmax[2] - gmin[2]);
1145       break;
1146     default:
1147       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension %" PetscInt_FMT " not supported", dim);
1148     }
1149     PetscCallMPI(MPIU_Allreduce(&totalWeight, &global_weight, 1, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD));
1150     PetscCallMPI(MPIU_Allreduce(&totalCharge, &global_charge, 1, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD));
1151     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "dim = %" PetscInt_FMT "\ttotalWeight = %f, user->charges[species[0]] = %f\ttotalCharge = %f, Total Area = %f\n", dim, (double)global_weight, (double)user->charges[0], (double)global_charge, (double)Area));
1152     param->sigma = PetscAbsReal(global_charge / (Area));
1153 
1154     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sigma: %g\n", (double)param->sigma));
1155     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(x0,v0,t0,m0,q0,phi0): (%e, %e, %e, %e, %e, %e) - (P, V) = (%e, %e)\n", (double)param->x0, (double)param->v0, (double)param->t0, (double)param->m0, (double)param->q0, (double)param->phi0, (double)param->poissonNumber,
1156                           (double)param->vlasovNumber));
1157   }
1158   /* Setup Constants */
1159   {
1160     PetscDS    ds;
1161     Parameter *param;
1162     PetscCall(PetscBagGetData(user->bag, (void **)&param));
1163     PetscScalar constants[NUM_CONSTANTS];
1164     constants[SIGMA]   = param->sigma;
1165     constants[V0]      = param->v0;
1166     constants[T0]      = param->t0;
1167     constants[X0]      = param->x0;
1168     constants[M0]      = param->m0;
1169     constants[Q0]      = param->q0;
1170     constants[PHI0]    = param->phi0;
1171     constants[POISSON] = param->poissonNumber;
1172     constants[VLASOV]  = param->vlasovNumber;
1173     PetscCall(DMGetDS(dm, &ds));
1174     PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants));
1175   }
1176   PetscFunctionReturn(PETSC_SUCCESS);
1177 }
1178 
1179 static PetscErrorCode InitializeVelocities_Fake1D(DM sw, AppCtx *user)
1180 {
1181   DM         dm;
1182   PetscReal *v;
1183   PetscInt  *species, cStart, cEnd;
1184   PetscInt   dim, Np;
1185 
1186   PetscFunctionBegin;
1187   PetscCall(DMGetDimension(sw, &dim));
1188   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1189   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v));
1190   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
1191   PetscCall(DMSwarmGetCellDM(sw, &dm));
1192   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1193   PetscRandom rnd;
1194   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
1195   PetscCall(PetscRandomSetInterval(rnd, 0, 1.));
1196   PetscCall(PetscRandomSetFromOptions(rnd));
1197 
1198   for (PetscInt p = 0; p < Np; ++p) {
1199     PetscReal a[3] = {0., 0., 0.}, vel[3] = {0., 0., 0.};
1200 
1201     PetscCall(PetscRandomGetValueReal(rnd, &a[0]));
1202     if (user->perturbed_weights) {
1203       PetscCall(PetscPDFSampleConstant1D(a, NULL, vel));
1204     } else {
1205       PetscCall(PetscPDFSampleGaussian1D(a, NULL, vel));
1206     }
1207     v[p * dim] = vel[0];
1208   }
1209   PetscCall(PetscRandomDestroy(&rnd));
1210   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v));
1211   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1212   PetscFunctionReturn(PETSC_SUCCESS);
1213 }
1214 
1215 static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
1216 {
1217   PetscReal v0[2] = {1., 0.};
1218   PetscInt  dim;
1219 
1220   PetscFunctionBeginUser;
1221   PetscCall(DMGetDimension(dm, &dim));
1222   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1223   PetscCall(DMSetType(*sw, DMSWARM));
1224   PetscCall(DMSetDimension(*sw, dim));
1225   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1226   PetscCall(DMSwarmSetCellDM(*sw, dm));
1227   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1228   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
1229   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
1230   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL));
1231   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initVelocity", dim, PETSC_REAL));
1232   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "E_field", dim, PETSC_REAL));
1233   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "potential", dim, PETSC_REAL));
1234   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "charges", dim, PETSC_REAL));
1235   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
1236   PetscCall(DMSetApplicationContext(*sw, user));
1237   PetscCall(DMSetFromOptions(*sw));
1238   user->swarm = *sw;
1239   if (user->perturbed_weights) {
1240     PetscCall(InitializeParticles_PerturbedWeights(*sw, user));
1241   } else {
1242     PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
1243     PetscCall(DMSwarmInitializeCoordinates(*sw));
1244     if (user->fake_1D) {
1245       PetscCall(InitializeVelocities_Fake1D(*sw, user));
1246     } else {
1247       PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0));
1248     }
1249   }
1250   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
1251   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
1252   {
1253     Vec gc, gc0, gv, gv0;
1254 
1255     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
1256     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0));
1257     PetscCall(VecCopy(gc, gc0));
1258     PetscCall(VecViewFromOptions(gc, NULL, "-ic_x_view"));
1259     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
1260     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0));
1261     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "velocity", &gv));
1262     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initVelocity", &gv0));
1263     PetscCall(VecCopy(gv, gv0));
1264     PetscCall(VecViewFromOptions(gv, NULL, "-ic_v_view"));
1265     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "velocity", &gv));
1266     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initVelocity", &gv0));
1267   }
1268   {
1269     const char *fieldnames[2] = {DMSwarmPICField_coor, "velocity"};
1270 
1271     PetscCall(DMSwarmVectorDefineField(*sw, 2, fieldnames));
1272   }
1273   PetscFunctionReturn(PETSC_SUCCESS);
1274 }
1275 
1276 static PetscErrorCode ComputeFieldAtParticles_Coulomb(SNES snes, DM sw, PetscReal E[])
1277 {
1278   AppCtx     *user;
1279   PetscReal  *coords;
1280   PetscInt   *species, dim, Np, Ns;
1281   PetscMPIInt size;
1282 
1283   PetscFunctionBegin;
1284   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)snes), &size));
1285   PetscCheck(size == 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Coulomb code only works in serial");
1286   PetscCall(DMGetDimension(sw, &dim));
1287   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1288   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
1289   PetscCall(DMGetApplicationContext(sw, (void *)&user));
1290 
1291   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1292   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
1293   for (PetscInt p = 0; p < Np; ++p) {
1294     PetscReal *pcoord = &coords[p * dim];
1295     PetscReal  pE[3]  = {0., 0., 0.};
1296 
1297     /* Calculate field at particle p due to particle q */
1298     for (PetscInt q = 0; q < Np; ++q) {
1299       PetscReal *qcoord = &coords[q * dim];
1300       PetscReal  rpq[3], r, r3, q_q;
1301 
1302       if (p == q) continue;
1303       q_q = user->charges[species[q]] * 1.;
1304       for (PetscInt d = 0; d < dim; ++d) rpq[d] = pcoord[d] - qcoord[d];
1305       r = DMPlex_NormD_Internal(dim, rpq);
1306       if (r < PETSC_SQRT_MACHINE_EPSILON) continue;
1307       r3 = PetscPowRealInt(r, 3);
1308       for (PetscInt d = 0; d < dim; ++d) pE[d] += q_q * rpq[d] / r3;
1309     }
1310     for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = pE[d];
1311   }
1312   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1313   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1314   PetscFunctionReturn(PETSC_SUCCESS);
1315 }
1316 
1317 static PetscErrorCode ComputeFieldAtParticles_Primal(SNES snes, DM sw, PetscReal E[])
1318 {
1319   DM         dm;
1320   AppCtx    *user;
1321   PetscDS    ds;
1322   PetscFE    fe;
1323   Mat        M_p, M;
1324   Vec        phi, locPhi, rho, f;
1325   PetscReal *coords;
1326   PetscInt   dim, cStart, cEnd, Np;
1327 
1328   PetscFunctionBegin;
1329   PetscCall(DMGetApplicationContext(sw, (void *)&user));
1330   PetscCall(PetscLogEventBegin(user->ESolveEvent, snes, sw, 0, 0));
1331   PetscCall(DMGetDimension(sw, &dim));
1332   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1333 
1334   KSP          ksp;
1335   Vec          rho0;
1336   const char **oldFields;
1337   const char  *fields[1] = {"w_q"};
1338   PetscInt     Nf;
1339   const char **tmp;
1340 
1341   /* Create the charges rho */
1342   PetscCall(SNESGetDM(snes, &dm));
1343   PetscCall(DMSwarmVectorGetField(sw, &Nf, &tmp));
1344   PetscCall(PetscMalloc1(Nf, &oldFields));
1345   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscStrallocpy(tmp[f], (char **)&oldFields[f]));
1346   PetscCall(DMSwarmVectorDefineField(sw, 1, fields));
1347   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
1348   PetscCall(DMSwarmVectorDefineField(sw, Nf, oldFields));
1349   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(oldFields[f]));
1350   PetscCall(PetscFree(oldFields));
1351 
1352   PetscCall(DMCreateMassMatrix(dm, dm, &M));
1353   PetscCall(DMGetGlobalVector(dm, &rho0));
1354   PetscCall(PetscObjectSetName((PetscObject)rho0, "Charge density (rho0) from Primal Compute"));
1355   PetscCall(DMGetGlobalVector(dm, &rho));
1356   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
1357   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
1358 
1359   PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
1360   PetscCall(MatMultTranspose(M_p, f, rho));
1361   PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
1362   PetscCall(MatViewFromOptions(M, NULL, "-m_view"));
1363   PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
1364   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
1365 
1366   PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &ksp));
1367   PetscCall(KSPSetOptionsPrefix(ksp, "em_proj_"));
1368   PetscCall(KSPSetOperators(ksp, M, M));
1369   PetscCall(KSPSetFromOptions(ksp));
1370   PetscCall(KSPSolve(ksp, rho, rho0));
1371   PetscCall(VecViewFromOptions(rho0, NULL, "-rho0_view"));
1372 
1373   PetscInt           rhosize;
1374   PetscReal         *charges;
1375   const PetscScalar *rho_vals;
1376   PetscCall(DMSwarmGetField(sw, "charges", NULL, NULL, (void **)&charges));
1377   PetscCall(VecGetLocalSize(rho0, &rhosize));
1378   PetscCall(VecGetArrayRead(rho0, &rho_vals));
1379   for (PetscInt c = 0; c < rhosize; ++c) charges[c] = rho_vals[c];
1380   PetscCall(VecRestoreArrayRead(rho0, &rho_vals));
1381   PetscCall(DMSwarmRestoreField(sw, "charges", NULL, NULL, (void **)&charges));
1382 
1383   PetscCall(VecScale(rho, -1.0));
1384 
1385   PetscCall(VecViewFromOptions(rho0, NULL, "-rho0_view"));
1386   PetscCall(VecViewFromOptions(rho, NULL, "-rho_view"));
1387   PetscCall(DMRestoreGlobalVector(dm, &rho0));
1388   PetscCall(KSPDestroy(&ksp));
1389   PetscCall(MatDestroy(&M_p));
1390   PetscCall(MatDestroy(&M));
1391 
1392   PetscCall(DMGetGlobalVector(dm, &phi));
1393   PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
1394   PetscCall(VecSet(phi, 0.0));
1395   PetscCall(SNESSolve(snes, rho, phi));
1396   PetscCall(DMRestoreGlobalVector(dm, &rho));
1397   PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
1398 
1399   PetscInt           phisize;
1400   PetscReal         *pot;
1401   const PetscScalar *phi_vals;
1402   PetscCall(DMSwarmGetField(sw, "potential", NULL, NULL, (void **)&pot));
1403   PetscCall(VecGetLocalSize(phi, &phisize));
1404   PetscCall(VecGetArrayRead(phi, &phi_vals));
1405   for (PetscInt c = 0; c < phisize; ++c) pot[c] = phi_vals[c];
1406   PetscCall(VecRestoreArrayRead(phi, &phi_vals));
1407   PetscCall(DMSwarmRestoreField(sw, "potential", NULL, NULL, (void **)&pot));
1408 
1409   PetscCall(DMGetLocalVector(dm, &locPhi));
1410   PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi));
1411   PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi));
1412   PetscCall(DMRestoreGlobalVector(dm, &phi));
1413   PetscCall(PetscLogEventEnd(user->ESolveEvent, snes, sw, 0, 0));
1414 
1415   PetscCall(DMGetDS(dm, &ds));
1416   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
1417   PetscCall(DMSwarmSortGetAccess(sw));
1418   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1419   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1420 
1421   PetscCall(PetscLogEventBegin(user->ETabEvent, snes, sw, 0, 0));
1422   for (PetscInt c = cStart; c < cEnd; ++c) {
1423     PetscTabulation tab;
1424     PetscScalar    *clPhi = NULL;
1425     PetscReal      *pcoord, *refcoord;
1426     PetscReal       v[3], J[9], invJ[9], detJ;
1427     PetscInt       *points;
1428     PetscInt        Ncp;
1429 
1430     PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
1431     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
1432     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
1433     for (PetscInt cp = 0; cp < Ncp; ++cp)
1434       for (PetscInt d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d];
1435     PetscCall(DMPlexCoordinatesToReference(dm, PETSC_TRUE, c, Ncp, pcoord, refcoord));
1436     PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
1437     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v, J, invJ, &detJ));
1438     PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi));
1439     for (PetscInt cp = 0; cp < Ncp; ++cp) {
1440       const PetscReal *basisDer = tab->T[1];
1441       const PetscInt   p        = points[cp];
1442 
1443       for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = 0.;
1444       PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, clPhi, dim, invJ, NULL, cp, &E[p * dim]));
1445       for (PetscInt d = 0; d < dim; ++d) {
1446         E[p * dim + d] *= -1.0;
1447         if (user->fake_1D && d > 0) E[p * dim + d] = 0;
1448       }
1449     }
1450     PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi));
1451     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
1452     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
1453     PetscCall(PetscTabulationDestroy(&tab));
1454     PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points));
1455   }
1456   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1457   PetscCall(DMSwarmSortRestoreAccess(sw));
1458   PetscCall(DMRestoreLocalVector(dm, &locPhi));
1459   PetscCall(PetscLogEventEnd(user->ETabEvent, snes, sw, 0, 0));
1460   PetscFunctionReturn(PETSC_SUCCESS);
1461 }
1462 
1463 static PetscErrorCode ComputeFieldAtParticles_Mixed(SNES snes, DM sw, PetscReal E[])
1464 {
1465   AppCtx         *user;
1466   DM              dm, potential_dm;
1467   KSP             ksp;
1468   IS              potential_IS;
1469   PetscDS         ds;
1470   PetscFE         fe;
1471   PetscFEGeom     feGeometry;
1472   Mat             M_p, M;
1473   Vec             phi, locPhi, rho, f, temp_rho, rho0;
1474   PetscQuadrature q;
1475   PetscReal      *coords, *pot;
1476   PetscInt        dim, cStart, cEnd, Np, pot_field = 1;
1477   const char    **oldFields;
1478   const char     *fields[1] = {"w_q"};
1479   PetscInt        Nf;
1480   const char    **tmp;
1481 
1482   PetscFunctionBegin;
1483   PetscCall(DMGetApplicationContext(sw, &user));
1484   PetscCall(PetscLogEventBegin(user->ESolveEvent, snes, sw, 0, 0));
1485   PetscCall(DMGetDimension(sw, &dim));
1486   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1487 
1488   /* Create the charges rho */
1489   PetscCall(SNESGetDM(snes, &dm));
1490   PetscCall(DMGetGlobalVector(dm, &rho));
1491   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
1492 
1493   PetscCall(DMCreateSubDM(dm, 1, &pot_field, &potential_IS, &potential_dm));
1494 
1495   PetscCall(DMSwarmVectorGetField(sw, &Nf, &tmp));
1496   PetscCall(PetscMalloc1(Nf, &oldFields));
1497   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscStrallocpy(tmp[f], (char **)&oldFields[f]));
1498   PetscCall(DMSwarmVectorDefineField(sw, 1, fields));
1499   PetscCall(DMCreateMassMatrix(sw, potential_dm, &M_p));
1500   PetscCall(DMSwarmVectorDefineField(sw, Nf, oldFields));
1501   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(oldFields[f]));
1502   PetscCall(PetscFree(oldFields));
1503 
1504   PetscCall(DMCreateMassMatrix(potential_dm, potential_dm, &M));
1505   PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
1506   PetscCall(MatViewFromOptions(M, NULL, "-m_view"));
1507   PetscCall(DMGetGlobalVector(potential_dm, &temp_rho));
1508   PetscCall(PetscObjectSetName((PetscObject)temp_rho, "Mf"));
1509   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
1510   PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
1511   PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
1512   PetscCall(MatMultTranspose(M_p, f, temp_rho));
1513   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
1514   PetscCall(DMGetGlobalVector(potential_dm, &rho0));
1515   PetscCall(PetscObjectSetName((PetscObject)rho0, "Charge density (rho0) from Mixed Compute"));
1516 
1517   PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &ksp));
1518   PetscCall(KSPSetOptionsPrefix(ksp, "em_proj"));
1519   PetscCall(KSPSetOperators(ksp, M, M));
1520   PetscCall(KSPSetFromOptions(ksp));
1521   PetscCall(KSPSolve(ksp, temp_rho, rho0));
1522   PetscCall(VecViewFromOptions(rho0, NULL, "-rho0_view"));
1523 
1524   PetscInt           rhosize;
1525   PetscReal         *charges;
1526   const PetscScalar *rho_vals;
1527   Parameter         *param;
1528   PetscCall(PetscBagGetData(user->bag, (void **)&param));
1529   PetscCall(DMSwarmGetField(sw, "charges", NULL, NULL, (void **)&charges));
1530   PetscCall(VecGetLocalSize(rho0, &rhosize));
1531 
1532   /* Integral over reference element is size 1.  Reference element area is 4.  Scale rho0 by 1/4 because the basis function is 1/4 */
1533   PetscCall(VecScale(rho0, 0.25));
1534   PetscCall(VecGetArrayRead(rho0, &rho_vals));
1535   for (PetscInt c = 0; c < rhosize; ++c) charges[c] = rho_vals[c];
1536   PetscCall(VecRestoreArrayRead(rho0, &rho_vals));
1537   PetscCall(DMSwarmRestoreField(sw, "charges", NULL, NULL, (void **)&charges));
1538 
1539   PetscCall(VecISCopy(rho, potential_IS, SCATTER_FORWARD, temp_rho));
1540   PetscCall(VecScale(rho, 0.25));
1541   PetscCall(VecViewFromOptions(rho0, NULL, "-rho0_view"));
1542   PetscCall(VecViewFromOptions(temp_rho, NULL, "-temprho_view"));
1543   PetscCall(VecViewFromOptions(rho, NULL, "-rho_view"));
1544   PetscCall(DMRestoreGlobalVector(potential_dm, &temp_rho));
1545   PetscCall(DMRestoreGlobalVector(potential_dm, &rho0));
1546 
1547   PetscCall(MatDestroy(&M_p));
1548   PetscCall(MatDestroy(&M));
1549   PetscCall(KSPDestroy(&ksp));
1550   PetscCall(DMDestroy(&potential_dm));
1551   PetscCall(ISDestroy(&potential_IS));
1552 
1553   PetscCall(DMGetGlobalVector(dm, &phi));
1554   PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
1555   PetscCall(VecSet(phi, 0.0));
1556   PetscCall(SNESSolve(snes, rho, phi));
1557   PetscCall(DMRestoreGlobalVector(dm, &rho));
1558 
1559   PetscInt           phisize;
1560   const PetscScalar *phi_vals;
1561   PetscCall(DMSwarmGetField(sw, "potential", NULL, NULL, (void **)&pot));
1562   PetscCall(VecGetLocalSize(phi, &phisize));
1563   PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
1564   PetscCall(VecGetArrayRead(phi, &phi_vals));
1565   for (PetscInt c = 0; c < phisize; ++c) pot[c] = phi_vals[c];
1566   PetscCall(VecRestoreArrayRead(phi, &phi_vals));
1567   PetscCall(DMSwarmRestoreField(sw, "potential", NULL, NULL, (void **)&pot));
1568 
1569   PetscCall(DMGetLocalVector(dm, &locPhi));
1570   PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi));
1571   PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi));
1572   PetscCall(DMRestoreGlobalVector(dm, &phi));
1573   PetscCall(PetscLogEventEnd(user->ESolveEvent, snes, sw, 0, 0));
1574 
1575   PetscCall(PetscLogEventBegin(user->ETabEvent, snes, sw, 0, 0));
1576   PetscCall(DMGetDS(dm, &ds));
1577   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
1578   PetscCall(DMSwarmSortGetAccess(sw));
1579   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1580   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1581   PetscCall(PetscFEGetQuadrature(fe, &q));
1582   PetscCall(PetscFECreateCellGeometry(fe, q, &feGeometry));
1583   for (PetscInt c = cStart; c < cEnd; ++c) {
1584     PetscTabulation tab;
1585     PetscScalar    *clPhi = NULL;
1586     PetscReal      *pcoord, *refcoord;
1587     PetscInt       *points;
1588     PetscInt        Ncp;
1589 
1590     PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
1591     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
1592     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
1593     for (PetscInt cp = 0; cp < Ncp; ++cp)
1594       for (PetscInt d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d];
1595     PetscCall(DMPlexCoordinatesToReference(dm, PETSC_TRUE, c, Ncp, pcoord, refcoord));
1596     PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
1597     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, q, feGeometry.v, feGeometry.J, feGeometry.invJ, feGeometry.detJ));
1598     PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi));
1599 
1600     for (PetscInt cp = 0; cp < Ncp; ++cp) {
1601       const PetscInt p = points[cp];
1602 
1603       for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = 0.;
1604       PetscCall(PetscFEInterpolateAtPoints_Static(fe, tab, clPhi, &feGeometry, cp, &E[p * dim]));
1605       PetscCall(PetscFEPushforward(fe, &feGeometry, 1, &E[p * dim]));
1606       for (PetscInt d = 0; d < dim; ++d) {
1607         E[p * dim + d] *= -2.0;
1608         if (user->fake_1D && d > 0) E[p * dim + d] = 0;
1609       }
1610     }
1611     PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi));
1612     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
1613     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
1614     PetscCall(PetscTabulationDestroy(&tab));
1615     PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points));
1616   }
1617   PetscCall(PetscFEDestroyCellGeometry(fe, &feGeometry));
1618   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1619   PetscCall(DMSwarmSortRestoreAccess(sw));
1620   PetscCall(DMRestoreLocalVector(dm, &locPhi));
1621   PetscCall(PetscLogEventEnd(user->ETabEvent, snes, sw, 0, 0));
1622   PetscFunctionReturn(PETSC_SUCCESS);
1623 }
1624 
1625 static PetscErrorCode ComputeFieldAtParticles(SNES snes, DM sw, PetscReal E[])
1626 {
1627   AppCtx  *ctx;
1628   PetscInt dim, Np;
1629 
1630   PetscFunctionBegin;
1631   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1632   PetscValidHeaderSpecific(sw, DM_CLASSID, 2);
1633   PetscAssertPointer(E, 3);
1634   PetscCall(DMGetDimension(sw, &dim));
1635   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1636   PetscCall(DMGetApplicationContext(sw, &ctx));
1637   PetscCall(PetscArrayzero(E, Np * dim));
1638 
1639   switch (ctx->em) {
1640   case EM_PRIMAL:
1641     PetscCall(ComputeFieldAtParticles_Primal(snes, sw, E));
1642     break;
1643   case EM_COULOMB:
1644     PetscCall(ComputeFieldAtParticles_Coulomb(snes, sw, E));
1645     break;
1646   case EM_MIXED:
1647     PetscCall(ComputeFieldAtParticles_Mixed(snes, sw, E));
1648     break;
1649   case EM_NONE:
1650     break;
1651   default:
1652     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No solver for electrostatic model %s", EMTypes[ctx->em]);
1653   }
1654   PetscFunctionReturn(PETSC_SUCCESS);
1655 }
1656 
1657 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
1658 {
1659   DM                 sw;
1660   SNES               snes = ((AppCtx *)ctx)->snes;
1661   const PetscReal   *coords, *vel;
1662   const PetscScalar *u;
1663   PetscScalar       *g;
1664   PetscReal         *E, m_p = 1., q_p = -1.;
1665   PetscInt           dim, d, Np, p;
1666 
1667   PetscFunctionBeginUser;
1668   PetscCall(TSGetDM(ts, &sw));
1669   PetscCall(DMGetDimension(sw, &dim));
1670   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1671   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1672   PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
1673   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1674   PetscCall(VecGetArrayRead(U, &u));
1675   PetscCall(VecGetArray(G, &g));
1676 
1677   PetscCall(ComputeFieldAtParticles(snes, sw, E));
1678 
1679   Np /= 2 * dim;
1680   for (p = 0; p < Np; ++p) {
1681     for (d = 0; d < dim; ++d) {
1682       g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d];
1683       g[(p * 2 + 1) * dim + d] = q_p * E[p * dim + d] / m_p;
1684     }
1685   }
1686   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1687   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1688   PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
1689   PetscCall(VecRestoreArrayRead(U, &u));
1690   PetscCall(VecRestoreArray(G, &g));
1691   PetscFunctionReturn(PETSC_SUCCESS);
1692 }
1693 
1694 /* J_{ij} = dF_i/dx_j
1695    J_p = (  0   1)
1696          (-w^2  0)
1697    TODO Now there is another term with w^2 from the electric field. I think we will need to invert the operator.
1698         Perhaps we can approximate the Jacobian using only the cellwise P-P gradient from Coulomb
1699 */
1700 static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, void *ctx)
1701 {
1702   DM               sw;
1703   const PetscReal *coords, *vel;
1704   PetscInt         dim, d, Np, p, rStart;
1705 
1706   PetscFunctionBeginUser;
1707   PetscCall(TSGetDM(ts, &sw));
1708   PetscCall(DMGetDimension(sw, &dim));
1709   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1710   PetscCall(MatGetOwnershipRange(J, &rStart, NULL));
1711   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1712   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1713   Np /= 2 * dim;
1714   for (p = 0; p < Np; ++p) {
1715     const PetscReal x0      = coords[p * dim + 0];
1716     const PetscReal vy0     = vel[p * dim + 1];
1717     const PetscReal omega   = vy0 / x0;
1718     PetscScalar     vals[4] = {0., 1., -PetscSqr(omega), 0.};
1719 
1720     for (d = 0; d < dim; ++d) {
1721       const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
1722       PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
1723     }
1724   }
1725   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1726   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1727   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1728   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
1729   PetscFunctionReturn(PETSC_SUCCESS);
1730 }
1731 
1732 static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
1733 {
1734   AppCtx            *user = (AppCtx *)ctx;
1735   DM                 sw;
1736   const PetscScalar *v;
1737   PetscScalar       *xres;
1738   PetscInt           Np, p, d, dim;
1739 
1740   PetscFunctionBeginUser;
1741   PetscCall(PetscLogEventBegin(user->RhsXEvent, ts, 0, 0, 0));
1742   PetscCall(TSGetDM(ts, &sw));
1743   PetscCall(DMGetDimension(sw, &dim));
1744   PetscCall(VecGetLocalSize(Xres, &Np));
1745   PetscCall(VecGetArrayRead(V, &v));
1746   PetscCall(VecGetArray(Xres, &xres));
1747   Np /= dim;
1748   for (p = 0; p < Np; ++p) {
1749     for (d = 0; d < dim; ++d) {
1750       xres[p * dim + d] = v[p * dim + d];
1751       if (user->fake_1D && d > 0) xres[p * dim + d] = 0;
1752     }
1753   }
1754   PetscCall(VecRestoreArrayRead(V, &v));
1755   PetscCall(VecRestoreArray(Xres, &xres));
1756   PetscCall(PetscLogEventEnd(user->RhsXEvent, ts, 0, 0, 0));
1757   PetscFunctionReturn(PETSC_SUCCESS);
1758 }
1759 
1760 static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
1761 {
1762   DM                 sw;
1763   AppCtx            *user = (AppCtx *)ctx;
1764   SNES               snes = ((AppCtx *)ctx)->snes;
1765   const PetscScalar *x;
1766   const PetscReal   *coords, *vel;
1767   PetscReal         *E, m_p, q_p;
1768   PetscScalar       *vres;
1769   PetscInt           Np, p, dim, d;
1770   Parameter         *param;
1771 
1772   PetscFunctionBeginUser;
1773   PetscCall(PetscLogEventBegin(user->RhsVEvent, ts, 0, 0, 0));
1774   PetscCall(TSGetDM(ts, &sw));
1775   PetscCall(DMGetDimension(sw, &dim));
1776   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1777   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1778   PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
1779   PetscCall(PetscBagGetData(user->bag, (void **)&param));
1780   m_p = user->masses[0] * param->m0;
1781   q_p = user->charges[0] * param->q0;
1782   PetscCall(VecGetLocalSize(Vres, &Np));
1783   PetscCall(VecGetArrayRead(X, &x));
1784   PetscCall(VecGetArray(Vres, &vres));
1785   PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension must be 2");
1786   PetscCall(ComputeFieldAtParticles(snes, sw, E));
1787 
1788   Np /= dim;
1789   for (p = 0; p < Np; ++p) {
1790     for (d = 0; d < dim; ++d) {
1791       vres[p * dim + d] = q_p * E[p * dim + d] / m_p;
1792       if (user->fake_1D && d > 0) vres[p * dim + d] = 0.;
1793     }
1794   }
1795   PetscCall(VecRestoreArrayRead(X, &x));
1796   /*
1797     Synchronized, ordered output for parallel/sequential test cases.
1798     In the 1D (on the 2D mesh) case, every y component should be zero.
1799   */
1800   if (user->checkVRes) {
1801     PetscBool pr = user->checkVRes > 1 ? PETSC_TRUE : PETSC_FALSE;
1802     PetscInt  step;
1803 
1804     PetscCall(TSGetStepNumber(ts, &step));
1805     if (pr) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "step: %" PetscInt_FMT "\n", step));
1806     for (PetscInt p = 0; p < Np; ++p) {
1807       if (pr) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Residual: %.12g %.12g\n", (double)PetscRealPart(vres[p * dim + 0]), (double)PetscRealPart(vres[p * dim + 1])));
1808       PetscCheck(PetscAbsScalar(vres[p * dim + 1]) < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Y velocity should be 0., not %g", (double)PetscRealPart(vres[p * dim + 1]));
1809     }
1810     if (pr) PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
1811   }
1812   PetscCall(VecRestoreArray(Vres, &vres));
1813   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1814   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1815   PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
1816   PetscCall(PetscLogEventEnd(user->RhsVEvent, ts, 0, 0, 0));
1817   PetscFunctionReturn(PETSC_SUCCESS);
1818 }
1819 
1820 static PetscErrorCode CreateSolution(TS ts)
1821 {
1822   DM       sw;
1823   Vec      u;
1824   PetscInt dim, Np;
1825 
1826   PetscFunctionBegin;
1827   PetscCall(TSGetDM(ts, &sw));
1828   PetscCall(DMGetDimension(sw, &dim));
1829   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1830   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
1831   PetscCall(VecSetBlockSize(u, dim));
1832   PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE));
1833   PetscCall(VecSetUp(u));
1834   PetscCall(TSSetSolution(ts, u));
1835   PetscCall(VecDestroy(&u));
1836   PetscFunctionReturn(PETSC_SUCCESS);
1837 }
1838 
1839 static PetscErrorCode SetProblem(TS ts)
1840 {
1841   AppCtx *user;
1842   DM      sw;
1843 
1844   PetscFunctionBegin;
1845   PetscCall(TSGetDM(ts, &sw));
1846   PetscCall(DMGetApplicationContext(sw, (void **)&user));
1847   // Define unified system for (X, V)
1848   {
1849     Mat      J;
1850     PetscInt dim, Np;
1851 
1852     PetscCall(DMGetDimension(sw, &dim));
1853     PetscCall(DMSwarmGetLocalSize(sw, &Np));
1854     PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
1855     PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE));
1856     PetscCall(MatSetBlockSize(J, 2 * dim));
1857     PetscCall(MatSetFromOptions(J));
1858     PetscCall(MatSetUp(J));
1859     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user));
1860     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user));
1861     PetscCall(MatDestroy(&J));
1862   }
1863   /* Define split system for X and V */
1864   {
1865     Vec             u;
1866     IS              isx, isv, istmp;
1867     const PetscInt *idx;
1868     PetscInt        dim, Np, rstart;
1869 
1870     PetscCall(TSGetSolution(ts, &u));
1871     PetscCall(DMGetDimension(sw, &dim));
1872     PetscCall(DMSwarmGetLocalSize(sw, &Np));
1873     PetscCall(VecGetOwnershipRange(u, &rstart, NULL));
1874     PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 0, 2, &istmp));
1875     PetscCall(ISGetIndices(istmp, &idx));
1876     PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx));
1877     PetscCall(ISRestoreIndices(istmp, &idx));
1878     PetscCall(ISDestroy(&istmp));
1879     PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 1, 2, &istmp));
1880     PetscCall(ISGetIndices(istmp, &idx));
1881     PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv));
1882     PetscCall(ISRestoreIndices(istmp, &idx));
1883     PetscCall(ISDestroy(&istmp));
1884     PetscCall(TSRHSSplitSetIS(ts, "position", isx));
1885     PetscCall(TSRHSSplitSetIS(ts, "momentum", isv));
1886     PetscCall(ISDestroy(&isx));
1887     PetscCall(ISDestroy(&isv));
1888     PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user));
1889     PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user));
1890   }
1891   PetscFunctionReturn(PETSC_SUCCESS);
1892 }
1893 
1894 static PetscErrorCode DMSwarmTSRedistribute(TS ts)
1895 {
1896   DM        sw;
1897   Vec       u;
1898   PetscReal t, maxt, dt;
1899   PetscInt  n, maxn;
1900 
1901   PetscFunctionBegin;
1902   PetscCall(TSGetDM(ts, &sw));
1903   PetscCall(TSGetTime(ts, &t));
1904   PetscCall(TSGetMaxTime(ts, &maxt));
1905   PetscCall(TSGetTimeStep(ts, &dt));
1906   PetscCall(TSGetStepNumber(ts, &n));
1907   PetscCall(TSGetMaxSteps(ts, &maxn));
1908 
1909   PetscCall(TSReset(ts));
1910   PetscCall(TSSetDM(ts, sw));
1911   PetscCall(TSSetFromOptions(ts));
1912   PetscCall(TSSetTime(ts, t));
1913   PetscCall(TSSetMaxTime(ts, maxt));
1914   PetscCall(TSSetTimeStep(ts, dt));
1915   PetscCall(TSSetStepNumber(ts, n));
1916   PetscCall(TSSetMaxSteps(ts, maxn));
1917 
1918   PetscCall(CreateSolution(ts));
1919   PetscCall(SetProblem(ts));
1920   PetscCall(TSGetSolution(ts, &u));
1921   PetscFunctionReturn(PETSC_SUCCESS);
1922 }
1923 
1924 PetscErrorCode line(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx)
1925 {
1926   DM        sw, cdm;
1927   PetscInt  Np;
1928   PetscReal low[2], high[2];
1929   AppCtx   *user = (AppCtx *)ctx;
1930 
1931   sw = user->swarm;
1932   PetscCall(DMSwarmGetCellDM(sw, &cdm));
1933   // Get the bounding box so we can equally space the particles
1934   PetscCall(DMGetLocalBoundingBox(cdm, low, high));
1935   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1936   // shift it by h/2 so nothing is initialized directly on a boundary
1937   x[0] = ((high[0] - low[0]) / Np) * (p + 0.5);
1938   x[1] = 0.;
1939   return PETSC_SUCCESS;
1940 }
1941 
1942 /*
1943   InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values.
1944 
1945   Input Parameters:
1946 + ts         - The TS
1947 - useInitial - Flag to also set the initial conditions to the current coordinates and velocities and setup the problem
1948 
1949   Output Parameters:
1950 . u - The initialized solution vector
1951 
1952   Level: advanced
1953 
1954 .seealso: InitializeSolve()
1955 */
1956 static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial)
1957 {
1958   DM       sw;
1959   Vec      u, gc, gv, gc0, gv0;
1960   IS       isx, isv;
1961   PetscInt dim;
1962   AppCtx  *user;
1963 
1964   PetscFunctionBeginUser;
1965   PetscCall(TSGetDM(ts, &sw));
1966   PetscCall(DMGetApplicationContext(sw, &user));
1967   PetscCall(DMGetDimension(sw, &dim));
1968   if (useInitial) {
1969     PetscReal v0[2] = {1., 0.};
1970     if (user->perturbed_weights) {
1971       PetscCall(InitializeParticles_PerturbedWeights(sw, user));
1972     } else {
1973       PetscCall(DMSwarmComputeLocalSizeFromOptions(sw));
1974       PetscCall(DMSwarmInitializeCoordinates(sw));
1975       if (user->fake_1D) {
1976         PetscCall(InitializeVelocities_Fake1D(sw, user));
1977       } else {
1978         PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0));
1979       }
1980     }
1981     PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
1982     PetscCall(DMSwarmTSRedistribute(ts));
1983   }
1984   PetscCall(TSGetSolution(ts, &u));
1985   PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
1986   PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
1987   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
1988   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0));
1989   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
1990   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initVelocity", &gv0));
1991   if (useInitial) {
1992     PetscCall(VecCopy(gc, gc0));
1993     PetscCall(VecCopy(gv, gv0));
1994   }
1995   PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc));
1996   PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv));
1997   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
1998   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0));
1999   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
2000   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initVelocity", &gv0));
2001   PetscFunctionReturn(PETSC_SUCCESS);
2002 }
2003 
2004 static PetscErrorCode InitializeSolve(TS ts, Vec u)
2005 {
2006   PetscFunctionBegin;
2007   PetscCall(TSSetSolution(ts, u));
2008   PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE));
2009   PetscFunctionReturn(PETSC_SUCCESS);
2010 }
2011 
2012 static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
2013 {
2014   MPI_Comm           comm;
2015   DM                 sw;
2016   AppCtx            *user;
2017   const PetscScalar *u;
2018   const PetscReal   *coords, *vel;
2019   PetscScalar       *e;
2020   PetscReal          t;
2021   PetscInt           dim, Np, p;
2022 
2023   PetscFunctionBeginUser;
2024   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
2025   PetscCall(TSGetDM(ts, &sw));
2026   PetscCall(DMGetApplicationContext(sw, &user));
2027   PetscCall(DMGetDimension(sw, &dim));
2028   PetscCall(TSGetSolveTime(ts, &t));
2029   PetscCall(VecGetArray(E, &e));
2030   PetscCall(VecGetArrayRead(U, &u));
2031   PetscCall(DMSwarmGetLocalSize(sw, &Np));
2032   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
2033   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
2034   Np /= 2 * dim;
2035   for (p = 0; p < Np; ++p) {
2036     /* TODO generalize initial conditions and project into plane instead of assuming x-y */
2037     const PetscReal    r0    = DMPlex_NormD_Internal(dim, &coords[p * dim]);
2038     const PetscReal    th0   = PetscAtan2Real(coords[p * dim + 1], coords[p * dim + 0]);
2039     const PetscReal    v0    = DMPlex_NormD_Internal(dim, &vel[p * dim]);
2040     const PetscReal    omega = v0 / r0;
2041     const PetscReal    ct    = PetscCosReal(omega * t + th0);
2042     const PetscReal    st    = PetscSinReal(omega * t + th0);
2043     const PetscScalar *x     = &u[(p * 2 + 0) * dim];
2044     const PetscScalar *v     = &u[(p * 2 + 1) * dim];
2045     const PetscReal    xe[3] = {r0 * ct, r0 * st, 0.0};
2046     const PetscReal    ve[3] = {-v0 * st, v0 * ct, 0.0};
2047     PetscInt           d;
2048 
2049     for (d = 0; d < dim; ++d) {
2050       e[(p * 2 + 0) * dim + d] = x[d] - xe[d];
2051       e[(p * 2 + 1) * dim + d] = v[d] - ve[d];
2052     }
2053     if (user->error) {
2054       const PetscReal en   = 0.5 * DMPlex_DotRealD_Internal(dim, v, v);
2055       const PetscReal exen = 0.5 * PetscSqr(v0);
2056       PetscCall(PetscPrintf(comm, "t %.4g: p%" PetscInt_FMT " error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g (%.10lf%%)\n", (double)t, p, (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 0) * dim]), (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 1) * dim]), (double)x[0], (double)x[1], (double)v[0], (double)v[1], (double)xe[0], (double)xe[1], (double)ve[0], (double)ve[1], (double)en, (double)exen, (double)(PetscAbsReal(exen - en) * 100. / exen)));
2057     }
2058   }
2059   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
2060   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
2061   PetscCall(VecRestoreArrayRead(U, &u));
2062   PetscCall(VecRestoreArray(E, &e));
2063   PetscFunctionReturn(PETSC_SUCCESS);
2064 }
2065 
2066 static PetscErrorCode MigrateParticles(TS ts)
2067 {
2068   DM               sw, cdm;
2069   const PetscReal *L;
2070 
2071   PetscFunctionBeginUser;
2072   PetscCall(TSGetDM(ts, &sw));
2073   PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre"));
2074   {
2075     Vec        u, gc, gv, position, momentum;
2076     IS         isx, isv;
2077     PetscReal *pos, *mom;
2078 
2079     PetscCall(TSGetSolution(ts, &u));
2080     PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
2081     PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
2082     PetscCall(VecGetSubVector(u, isx, &position));
2083     PetscCall(VecGetSubVector(u, isv, &momentum));
2084     PetscCall(VecGetArray(position, &pos));
2085     PetscCall(VecGetArray(momentum, &mom));
2086     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
2087     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
2088     PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc));
2089     PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv));
2090 
2091     PetscCall(DMSwarmGetCellDM(sw, &cdm));
2092     PetscCall(DMGetPeriodicity(cdm, NULL, NULL, &L));
2093     if ((L[0] || L[1]) >= 0.) {
2094       PetscReal *x, *v, upper[3], lower[3];
2095       PetscInt   Np, dim;
2096 
2097       PetscCall(DMSwarmGetLocalSize(sw, &Np));
2098       PetscCall(DMGetDimension(cdm, &dim));
2099       PetscCall(DMGetBoundingBox(cdm, lower, upper));
2100       PetscCall(VecGetArray(gc, &x));
2101       PetscCall(VecGetArray(gv, &v));
2102       for (PetscInt p = 0; p < Np; ++p) {
2103         for (PetscInt d = 0; d < dim; ++d) {
2104           if (pos[p * dim + d] < lower[d]) {
2105             x[p * dim + d] = pos[p * dim + d] + (upper[d] - lower[d]);
2106           } else if (pos[p * dim + d] > upper[d]) {
2107             x[p * dim + d] = pos[p * dim + d] - (upper[d] - lower[d]);
2108           } else {
2109             x[p * dim + d] = pos[p * dim + d];
2110           }
2111           PetscCheck(x[p * dim + d] >= lower[d] && x[p * dim + d] <= upper[d], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "p: %" PetscInt_FMT "x[%" PetscInt_FMT "] %g", p, d, (double)x[p * dim + d]);
2112           v[p * dim + d] = mom[p * dim + d];
2113         }
2114       }
2115       PetscCall(VecRestoreArray(gc, &x));
2116       PetscCall(VecRestoreArray(gv, &v));
2117     }
2118     PetscCall(VecRestoreArray(position, &pos));
2119     PetscCall(VecRestoreArray(momentum, &mom));
2120     PetscCall(VecRestoreSubVector(u, isx, &position));
2121     PetscCall(VecRestoreSubVector(u, isv, &momentum));
2122     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
2123     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
2124   }
2125   PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
2126   PetscCall(DMSwarmTSRedistribute(ts));
2127   PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE));
2128   PetscFunctionReturn(PETSC_SUCCESS);
2129 }
2130 
2131 int main(int argc, char **argv)
2132 {
2133   DM        dm, sw;
2134   TS        ts;
2135   Vec       u;
2136   PetscReal dt;
2137   PetscInt  maxn;
2138   AppCtx    user;
2139 
2140   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2141   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2142   PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag));
2143   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
2144   PetscCall(CreatePoisson(dm, &user));
2145   PetscCall(CreateSwarm(dm, &user, &sw));
2146   PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
2147   PetscCall(InitializeConstants(sw, &user));
2148   PetscCall(DMSetApplicationContext(sw, &user));
2149 
2150   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2151   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2152   PetscCall(TSSetDM(ts, sw));
2153   PetscCall(TSSetMaxTime(ts, 0.1));
2154   PetscCall(TSSetTimeStep(ts, 0.00001));
2155   PetscCall(TSSetMaxSteps(ts, 100));
2156   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2157 
2158   if (user.efield_monitor) PetscCall(TSMonitorSet(ts, MonitorEField, &user, NULL));
2159   if (user.initial_monitor) PetscCall(TSMonitorSet(ts, MonitorInitialConditions, &user, NULL));
2160   if (user.monitor_positions) PetscCall(TSMonitorSet(ts, MonitorPositions_2D, &user, NULL));
2161   if (user.poisson_monitor) PetscCall(TSMonitorSet(ts, MonitorPoisson, &user, NULL));
2162 
2163   PetscCall(TSSetFromOptions(ts));
2164   PetscCall(TSGetTimeStep(ts, &dt));
2165   PetscCall(TSGetMaxSteps(ts, &maxn));
2166   user.steps    = maxn;
2167   user.stepSize = dt;
2168   PetscCall(SetupContext(dm, sw, &user));
2169   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
2170   PetscCall(TSSetComputeExactError(ts, ComputeError));
2171   PetscCall(TSSetPostStep(ts, MigrateParticles));
2172   PetscCall(CreateSolution(ts));
2173   PetscCall(TSGetSolution(ts, &u));
2174   PetscCall(TSComputeInitialCondition(ts, u));
2175   PetscCall(CheckNonNegativeWeights(sw, &user));
2176   PetscCall(TSSolve(ts, NULL));
2177 
2178   PetscCall(SNESDestroy(&user.snes));
2179   PetscCall(TSDestroy(&ts));
2180   PetscCall(DMDestroy(&sw));
2181   PetscCall(DMDestroy(&dm));
2182   PetscCall(DestroyContext(&user));
2183   PetscCall(PetscFinalize());
2184   return 0;
2185 }
2186 
2187 /*TEST
2188 
2189    build:
2190     requires: !complex double
2191 
2192   # This tests that we can put particles in a box and compute the Coulomb force
2193   # Recommend -draw_size 500,500
2194    testset:
2195      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
2196      args: -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 20,1 \
2197              -dm_plex_box_lower 0,-1 -dm_plex_box_upper 12.5664,1 \
2198              -dm_plex_box_bd periodic,none \
2199            -dm_swarm_coordinate_density constant -dm_swarm_num_particles 100 \
2200            -sigma 1.0e-8 -timeScale 2.0e-14 \
2201            -ts_type basicsymplectic -ts_basicsymplectic_type 1 \
2202              -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain 0 -ts_monitor_sp_swarm_phase 0 \
2203            -output_step 50 -check_vel_res
2204      test:
2205        suffix: none_1d
2206        args: -em_type none -error
2207      test:
2208        suffix: coulomb_1d
2209        args: -em_type coulomb
2210 
2211    # for viewers
2212    #-ts_monitor_sp_swarm_phase -ts_monitor_sp_swarm -em_snes_monitor -ts_monitor_sp_swarm_multi_species 0 -ts_monitor_sp_swarm_retain 0
2213    testset:
2214      nsize: {{1 2}}
2215      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
2216      args: -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 36,1 \
2217              -dm_plex_box_lower 0.,-0.5 -dm_plex_box_upper 12.5664,0.5 \
2218              -dm_plex_box_bd periodic,none \
2219            -vdm_plex_dim 1 -vdm_plex_simplex 0 -vdm_plex_box_faces 10 \
2220              -vdm_plex_box_lower -3 -vdm_plex_box_upper 3 \
2221            -dm_swarm_num_species 1 -twostream -charges -1.,1. -sigma 1.0e-8 \
2222              -cosine_coefficients 0.01,0.5 -perturbed_weights -total_weight 1. \
2223            -ts_type basicsymplectic -ts_basicsymplectic_type 2 \
2224              -ts_dt 0.01 -ts_max_time 5 -ts_max_steps 10 \
2225            -em_snes_atol 1.e-15 -em_snes_error_if_not_converged -em_ksp_error_if_not_converged \
2226            -output_step 1 -check_vel_res -dm_swarm_print_coords 1 -dm_swarm_print_weights 1
2227      test:
2228        suffix: two_stream_c0
2229        args: -em_type primal -petscfe_default_quadrature_order 2 -petscspace_degree 2 -em_pc_type svd
2230      test:
2231        suffix: two_stream_rt
2232        requires: superlu_dist
2233        args: -em_type mixed \
2234                -potential_petscspace_degree 0 \
2235                -potential_petscdualspace_lagrange_use_moments \
2236                -potential_petscdualspace_lagrange_moment_order 2 \
2237                -field_petscspace_degree 2 -field_petscfe_default_quadrature_order 1 \
2238                -field_petscspace_type sum \
2239                  -field_petscspace_variables 2 \
2240                  -field_petscspace_components 2 \
2241                  -field_petscspace_sum_spaces 2 \
2242                  -field_petscspace_sum_concatenate true \
2243                  -field_sumcomp_0_petscspace_variables 2 \
2244                  -field_sumcomp_0_petscspace_type tensor \
2245                  -field_sumcomp_0_petscspace_tensor_spaces 2 \
2246                  -field_sumcomp_0_petscspace_tensor_uniform false \
2247                  -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
2248                  -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
2249                  -field_sumcomp_1_petscspace_variables 2 \
2250                  -field_sumcomp_1_petscspace_type tensor \
2251                  -field_sumcomp_1_petscspace_tensor_spaces 2 \
2252                  -field_sumcomp_1_petscspace_tensor_uniform false \
2253                  -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
2254                  -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
2255                -field_petscdualspace_form_degree -1 \
2256                -field_petscdualspace_order 1 \
2257                -field_petscdualspace_lagrange_trimmed true \
2258              -em_snes_error_if_not_converged \
2259              -em_ksp_type preonly -em_ksp_error_if_not_converged \
2260              -em_pc_type fieldsplit -em_pc_fieldsplit_type schur \
2261                -em_pc_fieldsplit_schur_fact_type full -em_pc_fieldsplit_schur_precondition full \
2262                -em_fieldsplit_field_pc_type lu \
2263                  -em_fieldsplit_field_pc_factor_mat_solver_type superlu_dist \
2264                -em_fieldsplit_potential_pc_type svd
2265 
2266    # For an eyeball check, we use
2267    # -ts_max_steps 1000 -dm_plex_box_faces 10,1 -vdm_plex_box_faces 2000 -monitor_efield
2268    # For verification, we use
2269    # -ts_max_steps 1000 -dm_plex_box_faces 100,1 -vdm_plex_box_faces 8000 -dm_swarm_num_particles 800000 -monitor_efield
2270    # -ts_monitor_sp_swarm_multi_species 0 -ts_monitor_sp_swarm_retain 0 -ts_monitor_sp_swarm_phase 1 -draw_size 500,500
2271    testset:
2272      nsize: {{1 2}}
2273      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
2274      args: -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 10,1 \
2275              -dm_plex_box_lower 0.,-0.5 -dm_plex_box_upper 12.5664,0.5 \
2276              -dm_plex_box_bd periodic,none \
2277            -vdm_plex_dim 1 -vdm_plex_simplex 0 -vdm_plex_box_faces 10 \
2278              -vdm_plex_box_lower -10 -vdm_plex_box_upper 10 \
2279            -dm_swarm_num_species 1 -charges -1.,1. \
2280              -cosine_coefficients 0.01,0.5 -perturbed_weights -total_weight 1. \
2281            -ts_type basicsymplectic -ts_basicsymplectic_type 1 \
2282              -ts_dt 0.03 -ts_max_time 500 -ts_max_steps 1 \
2283            -em_snes_atol 1.e-12 -em_snes_error_if_not_converged -em_ksp_error_if_not_converged \
2284            -output_step 1 -check_vel_res -dm_swarm_print_coords 1 -dm_swarm_print_weights 1
2285 
2286      test:
2287        suffix: uniform_equilibrium_1d
2288        args: -cosine_coefficients 0.0,0.5 -em_type primal -petscspace_degree 1 -em_pc_type svd
2289      test:
2290        suffix: uniform_primal_1d
2291        args: -em_type primal -petscspace_degree 1 -em_pc_type svd
2292      test:
2293        requires: superlu_dist
2294        suffix: uniform_mixed_1d
2295        args: -em_type mixed \
2296                -potential_petscspace_degree 0 \
2297                -potential_petscdualspace_lagrange_use_moments \
2298                -potential_petscdualspace_lagrange_moment_order 2 \
2299                -field_petscspace_degree 2 -field_petscfe_default_quadrature_order 1 \
2300                -field_petscspace_type sum \
2301                  -field_petscspace_variables 2 \
2302                  -field_petscspace_components 2 \
2303                  -field_petscspace_sum_spaces 2 \
2304                  -field_petscspace_sum_concatenate true \
2305                  -field_sumcomp_0_petscspace_variables 2 \
2306                  -field_sumcomp_0_petscspace_type tensor \
2307                  -field_sumcomp_0_petscspace_tensor_spaces 2 \
2308                  -field_sumcomp_0_petscspace_tensor_uniform false \
2309                  -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
2310                  -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
2311                  -field_sumcomp_1_petscspace_variables 2 \
2312                  -field_sumcomp_1_petscspace_type tensor \
2313                  -field_sumcomp_1_petscspace_tensor_spaces 2 \
2314                  -field_sumcomp_1_petscspace_tensor_uniform false \
2315                  -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
2316                  -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
2317                -field_petscdualspace_form_degree -1 \
2318                -field_petscdualspace_order 1 \
2319                -field_petscdualspace_lagrange_trimmed true \
2320              -em_snes_error_if_not_converged \
2321              -em_ksp_type preonly -em_ksp_error_if_not_converged \
2322              -em_pc_type fieldsplit -em_pc_fieldsplit_type schur \
2323                -em_pc_fieldsplit_schur_fact_type full -em_pc_fieldsplit_schur_precondition full \
2324                -em_fieldsplit_field_pc_type lu \
2325                  -em_fieldsplit_field_pc_factor_mat_solver_type superlu_dist \
2326                -em_fieldsplit_potential_pc_type svd
2327 
2328 TEST*/
2329