xref: /petsc/src/ts/tutorials/ex77.c (revision efa12513287cff49a2b9648ae83199dcbfaad71a)
1 static char help[] = "Time-dependent reactive low Mach Flow in 2d and 3d channels with finite elements.\n\
2 We solve the reactive low Mach flow problem in a rectangular domain\n\
3 using a parallel unstructured mesh (DMPLEX) to discretize the flow\n\
4 and particles (DWSWARM) to discretize the chemical species.\n\n\n";
5 
6 /*F
7 This low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the
8 finite element method on an unstructured mesh. The weak form equations are
9 
10 \begin{align*}
11     < q, \nabla\cdot u > = 0
12     <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p >  - < v, f  >  = 0
13     < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0
14 \end{align*}
15 
16 where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity.
17 
18 For visualization, use
19 
20   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
21 
22 The particles can be visualized using
23 
24   -part_dm_view draw -part_dm_view_swarm_radius 0.03
25 
26 F*/
27 
28 #include <petscdmplex.h>
29 #include <petscdmswarm.h>
30 #include <petscts.h>
31 #include <petscds.h>
32 #include <petscbag.h>
33 
34 typedef enum {SOL_TRIG_TRIG, NUM_SOL_TYPES} SolType;
35 const char *solTypes[NUM_SOL_TYPES+1] = {"trig_trig",  "unknown"};
36 
37 typedef enum {PART_LAYOUT_CELL, PART_LAYOUT_BOX, NUM_PART_LAYOUT_TYPES} PartLayoutType;
38 const char *partLayoutTypes[NUM_PART_LAYOUT_TYPES+1] = {"cell", "box",  "unknown"};
39 
40 typedef struct {
41   PetscReal nu;    /* Kinematic viscosity */
42   PetscReal alpha; /* Thermal diffusivity */
43   PetscReal T_in;  /* Inlet temperature*/
44   PetscReal omega; /* Rotation speed in MMS benchmark */
45 } Parameter;
46 
47 typedef struct {
48   /* Problem definition */
49   PetscBag       bag;          /* Holds problem parameters */
50   SolType        solType;      /* MMS solution type */
51   PartLayoutType partLayout;   /* Type of particle distribution */
52   PetscInt       Npc;          /* The initial number of particles per cell */
53   PetscReal      partLower[3]; /* Lower left corner of particle box */
54   PetscReal      partUpper[3]; /* Upper right corner of particle box */
55   PetscInt       Npb;          /* The initial number of particles per box dimension */
56 } AppCtx;
57 
58 typedef struct {
59   PetscReal ti; /* The time for ui, at the beginning of the advection solve */
60   PetscReal tf; /* The time for uf, at the end of the advection solve */
61   Vec       ui; /* The PDE solution field at ti */
62   Vec       uf; /* The PDE solution field at tf */
63   Vec       x0; /* The initial particle positions at t = 0 */
64   PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
65   AppCtx   *ctx; /* Context for exact solution */
66 } AdvCtx;
67 
68 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
69 {
70   PetscInt d;
71   for (d = 0; d < Nc; ++d) u[d] = 0.0;
72   return 0;
73 }
74 
75 static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
76 {
77   PetscInt d;
78   for (d = 0; d < Nc; ++d) u[d] = 1.0;
79   return 0;
80 }
81 
82 /*
83   CASE: trigonometric-trigonometric
84   In 2D we use exact solution:
85 
86     x = r0 cos(w t + theta0)  r0     = sqrt(x0^2 + y0^2)
87     y = r0 sin(w t + theta0)  theta0 = arctan(y0/x0)
88     u = -w r0 sin(theta0) = -w y
89     v =  w r0 cos(theta0) =  w x
90     p = x + y - 1
91     T = t + x + y
92     f = <1, 1>
93     Q = 1 + w (x - y)/r
94 
95   so that
96 
97     \nabla \cdot u = 0 + 0 = 0
98 
99   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
100     = <0, 0> + u_i d_i u_j - \nu 0 + <1, 1>
101     = <1, 1> + w^2 <-y, x> . <<0, 1>, <-1, 0>>
102     = <1, 1> + w^2 <-x, -y>
103     = <1, 1> - w^2 <x, y>
104 
105   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
106     = 1 + <u, v> . <1, 1> - \alpha 0
107     = 1 + u + v
108 */
109 static PetscErrorCode trig_trig_x(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *x, void *ctx)
110 {
111   const PetscReal x0     = X[0];
112   const PetscReal y0     = X[1];
113   const PetscReal R0     = PetscSqrtReal(x0*x0 + y0*y0);
114   const PetscReal theta0 = PetscAtan2Real(y0, x0);
115   Parameter      *p      = (Parameter *) ctx;
116 
117   x[0] = R0*PetscCosReal(p->omega*time + theta0);
118   x[1] = R0*PetscSinReal(p->omega*time + theta0);
119   return 0;
120 }
121 static PetscErrorCode trig_trig_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
122 {
123   Parameter *p = (Parameter *) ctx;
124 
125   u[0] = -p->omega*X[1];
126   u[1] =  p->omega*X[0];
127   return 0;
128 }
129 static PetscErrorCode trig_trig_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
130 {
131   u[0] = 0.0;
132   u[1] = 0.0;
133   return 0;
134 }
135 
136 static PetscErrorCode trig_trig_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
137 {
138   p[0] = X[0] + X[1] - 1.0;
139   return 0;
140 }
141 
142 static PetscErrorCode trig_trig_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
143 {
144   T[0] = time + X[0] + X[1];
145   return 0;
146 }
147 static PetscErrorCode trig_trig_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
148 {
149   T[0] = 1.0;
150   return 0;
151 }
152 
153 static void f0_trig_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
154                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
155                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
156                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
157 {
158   const PetscReal omega = PetscRealPart(constants[3]);
159   PetscInt        Nc    = dim;
160   PetscInt        c, d;
161 
162   for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d];
163 
164   for (c = 0; c < Nc; ++c) {
165     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
166   }
167   f0[0] -= 1.0 - omega*omega*X[0];
168   f0[1] -= 1.0 - omega*omega*X[1];
169 }
170 
171 static void f0_trig_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
172                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
173                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
174                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
175 {
176   const PetscReal omega = PetscRealPart(constants[3]);
177   PetscInt        d;
178 
179   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
180   f0[0] += u_t[uOff[2]] - (1.0 + omega*(X[0] - X[1]));
181 }
182 
183 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
184                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
185                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
186                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
187 {
188   PetscInt d;
189   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
190 }
191 
192 /*f1_v = \nu[grad(u) + grad(u)^T] - pI */
193 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
194                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
195                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
196                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
197 {
198   const PetscReal nu = PetscRealPart(constants[0]);
199   const PetscInt    Nc = dim;
200   PetscInt        c, d;
201 
202   for (c = 0; c < Nc; ++c) {
203     for (d = 0; d < dim; ++d) {
204       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
205       //f1[c*dim+d] = nu*u_x[c*dim+d];
206     }
207     f1[c*dim+c] -= u[uOff[1]];
208   }
209 }
210 
211 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
212                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
213                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
214                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
215 {
216   const PetscReal alpha = PetscRealPart(constants[1]);
217   PetscInt d;
218   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
219 }
220 
221 /*Jacobians*/
222 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
223                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
224                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
225                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
226 {
227   PetscInt d;
228   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
229 }
230 
231 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
232                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
233                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
234                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
235 {
236   PetscInt c, d;
237   const PetscInt  Nc = dim;
238 
239   for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift;
240 
241   for (c = 0; c < Nc; ++c) {
242     for (d = 0; d < dim; ++d) {
243       g0[c*Nc+d] += u_x[c*Nc+d];
244     }
245   }
246 }
247 
248 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
249                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
250                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
251                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
252 {
253   PetscInt NcI = dim;
254   PetscInt NcJ = dim;
255   PetscInt c, d, e;
256 
257   for (c = 0; c < NcI; ++c) {
258     for (d = 0; d < NcJ; ++d) {
259       for (e = 0; e < dim; ++e) {
260         if (c == d) {
261           g1[(c*NcJ+d)*dim+e] += u[e];
262         }
263       }
264     }
265   }
266 }
267 
268 
269 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
270                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
271                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
272                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
273 {
274   PetscInt d;
275   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
276 }
277 
278 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
279                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
280                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
281                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
282 {
283    const PetscReal nu = PetscRealPart(constants[0]);
284    const PetscInt  Nc = dim;
285    PetscInt        c, d;
286 
287   for (c = 0; c < Nc; ++c) {
288     for (d = 0; d < dim; ++d) {
289       g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU
290       g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose
291     }
292   }
293 }
294 
295 static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
296                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
297                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
298                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
299 {
300   PetscInt d;
301   for (d = 0; d < dim; ++d) g0[d] = u_tShift;
302 }
303 
304 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
305                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
306                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
307                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
308 {
309   PetscInt d;
310   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
311 }
312 
313 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
314                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
315                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
316                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
317 {
318   PetscInt d;
319   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
320 }
321 
322 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
323                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
324                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
325                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
326 {
327   const PetscReal alpha = PetscRealPart(constants[1]);
328   PetscInt               d;
329 
330   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
331 }
332 
333 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
334 {
335   PetscInt       sol, pl, n;
336   PetscErrorCode ierr;
337 
338   PetscFunctionBeginUser;
339   options->solType    = SOL_TRIG_TRIG;
340   options->partLayout = PART_LAYOUT_CELL;
341   options->Npc        = 1;
342   options->Npb        = 1;
343 
344   options->partLower[0] = options->partLower[1] = options->partLower[2] = 0.;
345   options->partUpper[0] = options->partUpper[1] = options->partUpper[2] = 1.;
346   ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr);
347   sol  = options->solType;
348   ierr = PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
349   options->solType = (SolType) sol;
350   pl   = options->partLayout;
351   ierr = PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL);CHKERRQ(ierr);
352   options->partLayout = (PartLayoutType) pl;
353   ierr = PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL);CHKERRQ(ierr);
354   n    = 3;
355   ierr = PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL);CHKERRQ(ierr);
356   n    = 3;
357   ierr = PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL);CHKERRQ(ierr);
358   ierr = PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL);CHKERRQ(ierr);
359   ierr = PetscOptionsEnd();
360   PetscFunctionReturn(0);
361 }
362 
363 static PetscErrorCode SetupParameters(AppCtx *user)
364 {
365   PetscBag       bag;
366   Parameter     *p;
367   PetscErrorCode ierr;
368 
369   PetscFunctionBeginUser;
370   /* setup PETSc parameter bag */
371   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
372   ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr);
373   bag  = user->bag;
374   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0, "nu",    "Kinematic viscosity");CHKERRQ(ierr);
375   ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr);
376   ierr = PetscBagRegisterReal(bag, &p->T_in,  1.0, "T_in",  "Inlet temperature");CHKERRQ(ierr);
377   ierr = PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark");CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
382 {
383   PetscErrorCode ierr;
384 
385   PetscFunctionBeginUser;
386   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
387   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
388   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
393 {
394   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
395   PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
396   PetscDS          prob;
397   Parameter       *ctx;
398   PetscInt         id;
399   PetscErrorCode   ierr;
400 
401   PetscFunctionBeginUser;
402   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
403   switch(user->solType){
404   case SOL_TRIG_TRIG:
405     ierr = PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v);CHKERRQ(ierr);
406     ierr = PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w);CHKERRQ(ierr);
407 
408     exactFuncs[0]   = trig_trig_u;
409     exactFuncs[1]   = trig_trig_p;
410     exactFuncs[2]   = trig_trig_T;
411     exactFuncs_t[0] = trig_trig_u_t;
412     exactFuncs_t[1] = NULL;
413     exactFuncs_t[2] = trig_trig_T_t;
414     break;
415    default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
416   }
417 
418   ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr);
419 
420   ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu,  NULL,  g3_vu);CHKERRQ(ierr);
421   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_vp, NULL);CHKERRQ(ierr);
422   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL,  NULL);CHKERRQ(ierr);
423   ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL,  NULL);CHKERRQ(ierr);
424   ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL,  g3_wT);CHKERRQ(ierr);
425   /* Setup constants */
426   {
427     Parameter  *param;
428     PetscScalar constants[4];
429 
430     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
431 
432     constants[0] = param->nu;
433     constants[1] = param->alpha;
434     constants[2] = param->T_in;
435     constants[3] = param->omega;
436     ierr = PetscDSSetConstants(prob, 4, constants);CHKERRQ(ierr);
437   }
438   /* Setup Boundary Conditions */
439   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
440   id   = 3;
441   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity",    "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
442   id   = 1;
443   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
444   id   = 2;
445   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity",  "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
446   id   = 4;
447   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity",   "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
448   id   = 3;
449   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp",    "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
450   id   = 1;
451   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
452   id   = 2;
453   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp",  "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
454   id   = 4;
455   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp",   "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
456 
457   /*setup exact solution.*/
458   ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr);
459   ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr);
460   ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr);
461   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr);
462   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr);
463   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 /* x_t = v
468 
469    Note that here we use the velocity field at t_{n+1} to advect the particles from
470    t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or
471    the method of characteristics.
472 */
473 static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, void *ctx)
474 {
475   AdvCtx             *adv = (AdvCtx *) ctx;
476   Vec                 u   = adv->ui;
477   DM                  sdm, dm, vdm;
478   Vec                 vel, locvel, pvel;
479   IS                  vis;
480   DMInterpolationInfo ictx;
481   const PetscScalar  *coords, *v;
482   PetscScalar        *f;
483   PetscInt            vf[1] = {0};
484   PetscInt            dim, Np;
485   PetscErrorCode      ierr;
486 
487   PetscFunctionBeginUser;
488   ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr);
489   ierr = DMSwarmGetCellDM(sdm, &dm);CHKERRQ(ierr);
490   ierr = DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor);CHKERRQ(ierr);
491   ierr = DMGetGlobalVector(sdm, &pvel);CHKERRQ(ierr);
492   ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr);
493   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
494   /* Get local velocity */
495   ierr = DMCreateSubDM(dm, 1, vf, &vis, &vdm);CHKERRQ(ierr);
496   ierr = VecGetSubVector(u, vis, &vel);CHKERRQ(ierr);
497   ierr = DMGetLocalVector(vdm, &locvel);CHKERRQ(ierr);
498   ierr = DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL);CHKERRQ(ierr);
499   ierr = DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel);CHKERRQ(ierr);
500   ierr = DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel);CHKERRQ(ierr);
501   ierr = VecRestoreSubVector(u, vis, &vel);CHKERRQ(ierr);
502   ierr = ISDestroy(&vis);CHKERRQ(ierr);
503   /* Interpolate velocity */
504   ierr = DMInterpolationCreate(PETSC_COMM_SELF, &ictx);CHKERRQ(ierr);
505   ierr = DMInterpolationSetDim(ictx, dim);CHKERRQ(ierr);
506   ierr = DMInterpolationSetDof(ictx, dim);CHKERRQ(ierr);
507   ierr = VecGetArrayRead(X, &coords);CHKERRQ(ierr);
508   ierr = DMInterpolationAddPoints(ictx, Np, (PetscReal *) coords);CHKERRQ(ierr);
509   ierr = VecRestoreArrayRead(X, &coords);CHKERRQ(ierr);
510   ierr = DMInterpolationSetUp(ictx, vdm, PETSC_FALSE);CHKERRQ(ierr);
511   ierr = DMInterpolationEvaluate(ictx, vdm, locvel, pvel);CHKERRQ(ierr);
512   ierr = DMInterpolationDestroy(&ictx);CHKERRQ(ierr);
513   ierr = DMRestoreLocalVector(vdm, &locvel);CHKERRQ(ierr);
514   ierr = DMDestroy(&vdm);CHKERRQ(ierr);
515 
516   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
517   ierr = VecGetArrayRead(pvel, &v);CHKERRQ(ierr);
518   ierr = PetscArraycpy(f, v, Np*dim);CHKERRQ(ierr);
519   ierr = VecRestoreArrayRead(pvel, &v);CHKERRQ(ierr);
520   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
521   ierr = DMRestoreGlobalVector(sdm, &pvel);CHKERRQ(ierr);
522   PetscFunctionReturn(0);
523 }
524 
525 static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u)
526 {
527   AppCtx        *user;
528   void          *ctx;
529   DM             dm;
530   PetscScalar   *coords;
531   PetscReal      x[3], dx[3];
532   PetscInt       n[3];
533   PetscInt       Np, dim, d, i, j, k;
534   PetscErrorCode ierr;
535 
536   PetscFunctionBegin;
537   ierr = TSGetApplicationContext(ts, &ctx);CHKERRQ(ierr);
538   user = ((AdvCtx *) ctx)->ctx;
539   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
540   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
541   switch (user->partLayout) {
542     case PART_LAYOUT_CELL:
543       ierr = DMSwarmSetPointCoordinatesRandom(dm, user->Npc);CHKERRQ(ierr);
544       break;
545     case PART_LAYOUT_BOX:
546       Np = 1;
547       for (d = 0; d < dim; ++d) {
548         n[d]  = user->Npb;
549         dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1);
550         Np   *= n[d];
551       }
552       ierr = VecGetArray(u, &coords);CHKERRQ(ierr);
553       switch (dim) {
554         case 2:
555           x[0] = user->partLower[0];
556           for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
557             x[1] = user->partLower[1];
558             for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
559               const PetscInt p = j*n[0] + i;
560               for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d];
561             }
562           }
563           break;
564         case 3:
565           x[0] = user->partLower[0];
566           for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
567             x[1] = user->partLower[1];
568             for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
569               x[2] = user->partLower[2];
570               for (k = 0; k < n[2]; ++k, x[2] += dx[2]) {
571                 const PetscInt p = (k*n[1] + j)*n[0] + i;
572                 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d];
573               }
574             }
575           }
576           break;
577         default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim);
578       }
579       ierr = VecRestoreArray(u, &coords);CHKERRQ(ierr);
580       break;
581     default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]);
582   }
583   PetscFunctionReturn(0);
584 }
585 
586 static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user)
587 {
588   DM              cdm = dm;
589   PetscFE         fe[3];
590   Parameter      *param;
591   PetscInt       *cellid, n[3];
592   PetscReal       x[3], dx[3];
593   PetscScalar    *coords;
594   DMPolytopeType  ct;
595   PetscInt        dim, d, cStart, cEnd, c, Np, p, i, j, k;
596   PetscBool       simplex;
597   MPI_Comm        comm;
598   PetscErrorCode  ierr;
599 
600   PetscFunctionBeginUser;
601   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
602   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
603   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
604   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
605   /* Create finite element */
606   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
607   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
608   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
609 
610   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
611   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
612   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
613 
614   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr);
615   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
616   ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr);
617 
618   /* Set discretization and boundary conditions for each mesh */
619   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
620   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
621   ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr);
622   ierr = DMCreateDS(dm);CHKERRQ(ierr);
623   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
624   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
625   while (cdm) {
626     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
627     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
628   }
629   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
630   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
631   ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr);
632 
633   {
634     PetscObject  pressure;
635     MatNullSpace nullspacePres;
636 
637     ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr);
638     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr);
639     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr);
640     ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr);
641   }
642 
643   /* Setup particle information */
644   ierr = DMSwarmSetType(sdm, DMSWARM_PIC);CHKERRQ(ierr);
645   ierr = DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL);CHKERRQ(ierr);
646   ierr = DMSwarmFinalizeFieldRegister(sdm);CHKERRQ(ierr);
647   switch (user->partLayout) {
648     case PART_LAYOUT_CELL:
649       ierr = DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0);CHKERRQ(ierr);
650       ierr = DMSetFromOptions(sdm);CHKERRQ(ierr);
651       ierr = DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
652       for (c = cStart; c < cEnd; ++c) {
653         for (p = 0; p < user->Npc; ++p) {
654           const PetscInt n = c*user->Npc + p;
655 
656           cellid[n] = c;
657         }
658       }
659       ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
660       ierr = DMSwarmSetPointCoordinatesRandom(sdm, user->Npc);CHKERRQ(ierr);
661       break;
662     case PART_LAYOUT_BOX:
663       Np = 1;
664       for (d = 0; d < dim; ++d) {
665         n[d]  = user->Npb;
666         dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1);
667         Np   *= n[d];
668       }
669       ierr = DMSwarmSetLocalSizes(sdm, Np, 0);CHKERRQ(ierr);
670       ierr = DMSetFromOptions(sdm);CHKERRQ(ierr);
671       ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
672       switch (dim) {
673         case 2:
674           x[0] = user->partLower[0];
675           for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
676             x[1] = user->partLower[1];
677             for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
678               const PetscInt p = j*n[0] + i;
679               for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d];
680             }
681           }
682           break;
683         case 3:
684           x[0] = user->partLower[0];
685           for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
686             x[1] = user->partLower[1];
687             for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
688               x[2] = user->partLower[2];
689               for (k = 0; k < n[2]; ++k, x[2] += dx[2]) {
690                 const PetscInt p = (k*n[1] + j)*n[0] + i;
691                 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d];
692               }
693             }
694           }
695           break;
696         default: SETERRQ1(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim);
697       }
698       ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
699       ierr = DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
700       for (p = 0; p < Np; ++p) cellid[p] = 0;
701       ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
702       ierr = DMSwarmMigrate(sdm, PETSC_TRUE);CHKERRQ(ierr);
703       break;
704     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]);
705   }
706   ierr = PetscObjectSetName((PetscObject) sdm, "Particles");CHKERRQ(ierr);
707   ierr = DMViewFromOptions(sdm, NULL, "-dm_view");CHKERRQ(ierr);
708   PetscFunctionReturn(0);
709 }
710 
711 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
712 {
713   Vec              vec;
714   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
715   PetscErrorCode   ierr;
716 
717   PetscFunctionBeginUser;
718   if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield);
719   funcs[nfield] = constant;
720   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
721   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
722   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
723   ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr);
724   ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr);
725   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr);
726   ierr = VecDestroy(&vec);CHKERRQ(ierr);
727   PetscFunctionReturn(0);
728 }
729 
730 static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u)
731 {
732   DM             dm;
733   MatNullSpace   nullsp;
734   PetscErrorCode ierr;
735 
736   PetscFunctionBegin;
737   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
738   ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr);
739   ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr);
740   ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
741   PetscFunctionReturn(0);
742 }
743 
744 /* Make the discrete pressure discretely divergence free */
745 static PetscErrorCode RemoveDiscretePressureNullspace(TS ts)
746 {
747   Vec            u;
748   PetscErrorCode ierr;
749 
750   PetscFunctionBegin;
751   ierr = TSGetSolution(ts, &u);CHKERRQ(ierr);
752   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
753   PetscFunctionReturn(0);
754 }
755 
756 static PetscErrorCode SetInitialConditions(TS ts, Vec u)
757 {
758   DM             dm;
759   PetscReal      t;
760   PetscErrorCode ierr;
761 
762   PetscFunctionBegin;
763   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
764   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
765   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
766   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
767   PetscFunctionReturn(0);
768 }
769 
770 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
771 {
772   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
773   void            *ctxs[3];
774   DM               dm;
775   PetscDS          ds;
776   Vec              v;
777   PetscReal        ferrors[3];
778   PetscInt         f;
779   PetscErrorCode   ierr;
780 
781   PetscFunctionBeginUser;
782   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
783   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
784 
785   for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);}
786   ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr);
787   ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g]\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1], (double) ferrors[2]);CHKERRQ(ierr);
788 
789   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
790   //ierr = TSGetSolution(ts, &u);CHKERRQ(ierr);
791   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
792   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
793   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
794 
795   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
796   // ierr = VecSet(v, 0.0);CHKERRQ(ierr);
797   ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr);
798   ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr);
799   ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr);
800   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
801 
802   PetscFunctionReturn(0);
803 }
804 
805 static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e)
806 {
807   AdvCtx            *adv;
808   DM                 sdm;
809   Parameter         *param;
810   const PetscScalar *xp0, *xp;
811   PetscScalar       *ep;
812   PetscReal          time;
813   PetscInt           dim, Np, p;
814   MPI_Comm           comm;
815   PetscErrorCode     ierr;
816 
817   PetscFunctionBeginUser;
818   ierr = TSGetTime(ts, &time);CHKERRQ(ierr);
819   ierr = TSGetApplicationContext(ts, (void **) &adv);CHKERRQ(ierr);
820   ierr = PetscBagGetData(adv->ctx->bag, (void **) &param);CHKERRQ(ierr);
821   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
822   ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr);
823   ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr);
824   ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr);
825   ierr = VecGetArrayRead(adv->x0, &xp0);CHKERRQ(ierr);
826   ierr = VecGetArrayRead(u, &xp);CHKERRQ(ierr);
827   ierr = VecGetArrayWrite(e, &ep);CHKERRQ(ierr);
828   for (p = 0; p < Np; ++p) {
829     PetscScalar x[3];
830     PetscReal   x0[3];
831     PetscInt    d;
832 
833     for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]);
834     ierr = adv->exact(dim, time, x0, 1, x, param);CHKERRQ(ierr);
835     for (d = 0; d < dim; ++d) ep[p*dim+d] += x[d] - xp[p*dim+d];
836   }
837   ierr = VecRestoreArrayRead(adv->x0, &xp0);CHKERRQ(ierr);
838   ierr = VecRestoreArrayRead(u, &xp);CHKERRQ(ierr);
839   ierr = VecRestoreArrayWrite(e, &ep);CHKERRQ(ierr);
840   PetscFunctionReturn(0);
841 }
842 
843 static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx)
844 {
845   AdvCtx            *adv = (AdvCtx *) ctx;
846   DM                 sdm;
847   Parameter         *param;
848   const PetscScalar *xp0, *xp;
849   PetscReal          error = 0.0;
850   PetscInt           dim, Np, p;
851   MPI_Comm           comm;
852   PetscErrorCode     ierr;
853 
854   PetscFunctionBeginUser;
855   ierr = PetscBagGetData(adv->ctx->bag, (void **) &param);CHKERRQ(ierr);
856   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
857   ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr);
858   ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr);
859   ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr);
860   ierr = VecGetArrayRead(adv->x0, &xp0);CHKERRQ(ierr);
861   ierr = VecGetArrayRead(u, &xp);CHKERRQ(ierr);
862   for (p = 0; p < Np; ++p) {
863     PetscScalar x[3];
864     PetscReal   x0[3];
865     PetscReal   perror = 0.0;
866     PetscInt    d;
867 
868     for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]);
869     ierr = adv->exact(dim, time, x0, 1, x, param);CHKERRQ(ierr);
870     for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p*dim+d]));
871     error += perror;
872   }
873   ierr = VecRestoreArrayRead(adv->x0, &xp0);CHKERRQ(ierr);
874   ierr = VecRestoreArrayRead(u, &xp);CHKERRQ(ierr);
875   ierr = PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int) step, (double) time, (double) error);CHKERRQ(ierr);
876   PetscFunctionReturn(0);
877 }
878 
879 static PetscErrorCode AdvectParticles(TS ts)
880 {
881   TS                 sts;
882   DM                 sdm;
883   Vec                p;
884   AdvCtx            *adv;
885   PetscScalar       *coord, *a;
886   const PetscScalar *ca;
887   PetscReal          time;
888   PetscInt           n;
889   PetscErrorCode     ierr;
890 
891   PetscFunctionBeginUser;
892   ierr = PetscObjectQuery((PetscObject) ts, "_SwarmTS",  (PetscObject *) &sts);CHKERRQ(ierr);
893   ierr = PetscObjectQuery((PetscObject) ts, "_SwarmSol", (PetscObject *) &p);CHKERRQ(ierr);
894   ierr = TSGetDM(sts, &sdm);CHKERRQ(ierr);
895   ierr = TSGetRHSFunction(sts, NULL, NULL, (void **) &adv);CHKERRQ(ierr);
896   ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coord);CHKERRQ(ierr);
897   ierr = VecGetLocalSize(p, &n);CHKERRQ(ierr);
898   ierr = VecGetArray(p, &a);CHKERRQ(ierr);
899   ierr = PetscArraycpy(a, coord, n);CHKERRQ(ierr);
900   ierr = VecRestoreArray(p, &a);CHKERRQ(ierr);
901   ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coord);CHKERRQ(ierr);
902   ierr = TSGetTime(ts, &time);CHKERRQ(ierr);
903   ierr = TSSetMaxTime(sts, time);CHKERRQ(ierr);
904   adv->tf = time;
905   ierr = TSSolve(sts, p);CHKERRQ(ierr);
906   ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coord);CHKERRQ(ierr);
907   ierr = VecGetLocalSize(p, &n);CHKERRQ(ierr);
908   ierr = VecGetArrayRead(p, &ca);CHKERRQ(ierr);
909   ierr = PetscArraycpy(coord, ca, n);CHKERRQ(ierr);
910   ierr = VecRestoreArrayRead(p, &ca);CHKERRQ(ierr);
911   ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coord);CHKERRQ(ierr);
912 
913   ierr = VecCopy(adv->uf, adv->ui);CHKERRQ(ierr);
914   adv->ti = adv->tf;
915 
916   ierr = DMSwarmMigrate(sdm, PETSC_TRUE);CHKERRQ(ierr);
917   ierr = DMViewFromOptions(sdm, NULL, "-dm_view");CHKERRQ(ierr);
918   PetscFunctionReturn(0);
919 }
920 
921 int main(int argc, char **argv)
922 {
923   DM              dm, sdm;
924   TS              ts, sts;
925   Vec             u, p, xtmp;
926   AppCtx          user;
927   AdvCtx          adv;
928   PetscReal       t;
929   PetscInt        dim;
930   PetscErrorCode  ierr;
931 
932   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
933   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
934   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
935   ierr = SetupParameters(&user);CHKERRQ(ierr);
936   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
937   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
938   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
939   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
940   /* Discretize chemical species */
941   ierr = DMCreate(PETSC_COMM_WORLD, &sdm);CHKERRQ(ierr);
942   ierr = PetscObjectSetOptionsPrefix((PetscObject) sdm, "part_");CHKERRQ(ierr);
943   ierr = DMSetType(sdm, DMSWARM);CHKERRQ(ierr);
944   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
945   ierr = DMSetDimension(sdm, dim);CHKERRQ(ierr);
946   ierr = DMSwarmSetCellDM(sdm, dm);CHKERRQ(ierr);
947   /* Setup problem */
948   ierr = SetupDiscretization(dm, sdm, &user);CHKERRQ(ierr);
949   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
950 
951   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
952   ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
953 
954   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr);
955   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr);
956   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr);
957   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
958   ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr);
959   ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
960   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
961 
962   ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */
963   ierr = SetInitialConditions(ts, u);CHKERRQ(ierr);
964   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
965   ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr);
966   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
967 
968   /* Setup particle position integrator */
969   ierr = TSCreate(PETSC_COMM_WORLD, &sts);CHKERRQ(ierr);
970   ierr = PetscObjectSetOptionsPrefix((PetscObject) sts, "part_");CHKERRQ(ierr);
971   ierr = TSSetDM(sts, sdm);CHKERRQ(ierr);
972   ierr = TSSetProblemType(sts, TS_NONLINEAR);CHKERRQ(ierr);
973   ierr = TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
974   ierr = TSMonitorSet(sts, MonitorParticleError, &adv, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
975   ierr = TSSetFromOptions(sts);CHKERRQ(ierr);
976   ierr = TSSetApplicationContext(sts, &adv);CHKERRQ(ierr);
977   ierr = TSSetComputeExactError(sts, ComputeParticleError);CHKERRQ(ierr);
978   ierr = TSSetComputeInitialCondition(sts, SetInitialParticleConditions);CHKERRQ(ierr);
979   adv.ti = t;
980   adv.uf = u;
981   ierr = VecDuplicate(adv.uf, &adv.ui);
982   ierr = VecCopy(u, adv.ui);CHKERRQ(ierr);
983   ierr = TSSetRHSFunction(sts, NULL, FreeStreaming, &adv);CHKERRQ(ierr);
984   ierr = TSSetPostStep(ts, AdvectParticles);CHKERRQ(ierr);
985   ierr = PetscObjectCompose((PetscObject) ts, "_SwarmTS", (PetscObject) sts);CHKERRQ(ierr);
986   ierr = DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor);CHKERRQ(ierr);
987   ierr = DMCreateGlobalVector(sdm, &p);CHKERRQ(ierr);
988   ierr = PetscObjectCompose((PetscObject) ts, "_SwarmSol", (PetscObject) p);CHKERRQ(ierr);
989   ierr = DMCreateGlobalVector(sdm, &adv.x0);CHKERRQ(ierr);
990   ierr = DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp);CHKERRQ(ierr);
991   ierr = VecCopy(xtmp, adv.x0);CHKERRQ(ierr);
992   ierr = DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp);CHKERRQ(ierr);
993   switch(user.solType){
994     case SOL_TRIG_TRIG: adv.exact = trig_trig_x;break;
995     default: SETERRQ2(PetscObjectComm((PetscObject) sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType);
996   }
997   adv.ctx = &user;
998 
999   ierr = TSSolve(ts, u);CHKERRQ(ierr);
1000   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
1001   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
1002 
1003   ierr = VecDestroy(&u);CHKERRQ(ierr);
1004   ierr = VecDestroy(&p);CHKERRQ(ierr);
1005   ierr = VecDestroy(&adv.x0);CHKERRQ(ierr);
1006   ierr = VecDestroy(&adv.ui);CHKERRQ(ierr);
1007   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1008   ierr = DMDestroy(&sdm);CHKERRQ(ierr);
1009   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1010   ierr = TSDestroy(&sts);CHKERRQ(ierr);
1011   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
1012   ierr = PetscFinalize();
1013   return ierr;
1014 }
1015 
1016 /*TEST
1017 
1018   # Swarm does not work with complex
1019   test:
1020     suffix: 2d_tri_p2_p1_p1_tconvp
1021     requires: triangle !single !complex
1022     args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \
1023       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1024       -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 -ts_monitor_cancel \
1025       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
1026       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
1027         -fieldsplit_0_pc_type lu \
1028         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
1029       -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \
1030       -part_ts_max_steps 2 -part_ts_dt 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel
1031 
1032 TEST*/
1033