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