xref: /petsc/src/ts/tutorials/ex77.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1a4e36a32SMatthew G. Knepley static char help[] = "Time-dependent reactive low Mach Flow in 2d and 3d channels with finite elements.\n\
2a4e36a32SMatthew G. Knepley We solve the reactive low Mach flow problem in a rectangular domain\n\
3a4e36a32SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to discretize the flow\n\
4a4e36a32SMatthew G. Knepley and particles (DWSWARM) to discretize the chemical species.\n\n\n";
5a4e36a32SMatthew G. Knepley 
6a4e36a32SMatthew G. Knepley /*F
7a4e36a32SMatthew G. Knepley This low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the
8a4e36a32SMatthew G. Knepley finite element method on an unstructured mesh. The weak form equations are
9a4e36a32SMatthew G. Knepley 
10a4e36a32SMatthew G. Knepley \begin{align*}
11a4e36a32SMatthew G. Knepley     < q, \nabla\cdot u > = 0
12a4e36a32SMatthew G. Knepley     <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p >  - < v, f  >  = 0
13a4e36a32SMatthew G. Knepley     < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0
14a4e36a32SMatthew G. Knepley \end{align*}
15a4e36a32SMatthew G. Knepley 
16a4e36a32SMatthew G. Knepley where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity.
17a4e36a32SMatthew G. Knepley 
18a4e36a32SMatthew G. Knepley For visualization, use
19a4e36a32SMatthew G. Knepley 
20a4e36a32SMatthew G. Knepley   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
21a4e36a32SMatthew G. Knepley 
22a4e36a32SMatthew G. Knepley The particles can be visualized using
23a4e36a32SMatthew G. Knepley 
24a4e36a32SMatthew G. Knepley   -part_dm_view draw -part_dm_view_swarm_radius 0.03
25a4e36a32SMatthew G. Knepley 
26a4e36a32SMatthew G. Knepley F*/
27a4e36a32SMatthew G. Knepley 
28a4e36a32SMatthew G. Knepley #include <petscdmplex.h>
29a4e36a32SMatthew G. Knepley #include <petscdmswarm.h>
30a4e36a32SMatthew G. Knepley #include <petscts.h>
31a4e36a32SMatthew G. Knepley #include <petscds.h>
32a4e36a32SMatthew G. Knepley #include <petscbag.h>
33a4e36a32SMatthew G. Knepley 
349371c9d4SSatish Balay typedef enum {
359371c9d4SSatish Balay   SOL_TRIG_TRIG,
369371c9d4SSatish Balay   NUM_SOL_TYPES
379371c9d4SSatish Balay } SolType;
38a4e36a32SMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES + 1] = {"trig_trig", "unknown"};
39a4e36a32SMatthew G. Knepley 
409371c9d4SSatish Balay typedef enum {
419371c9d4SSatish Balay   PART_LAYOUT_CELL,
429371c9d4SSatish Balay   PART_LAYOUT_BOX,
439371c9d4SSatish Balay   NUM_PART_LAYOUT_TYPES
449371c9d4SSatish Balay } PartLayoutType;
45a4e36a32SMatthew G. Knepley const char *partLayoutTypes[NUM_PART_LAYOUT_TYPES + 1] = {"cell", "box", "unknown"};
46a4e36a32SMatthew G. Knepley 
47a4e36a32SMatthew G. Knepley typedef struct {
48a4e36a32SMatthew G. Knepley   PetscReal nu;    /* Kinematic viscosity */
49a4e36a32SMatthew G. Knepley   PetscReal alpha; /* Thermal diffusivity */
50a4e36a32SMatthew G. Knepley   PetscReal T_in;  /* Inlet temperature*/
51a4e36a32SMatthew G. Knepley   PetscReal omega; /* Rotation speed in MMS benchmark */
52a4e36a32SMatthew G. Knepley } Parameter;
53a4e36a32SMatthew G. Knepley 
54a4e36a32SMatthew G. Knepley typedef struct {
55a4e36a32SMatthew G. Knepley   /* Problem definition */
56a4e36a32SMatthew G. Knepley   PetscBag       bag;          /* Holds problem parameters */
57a4e36a32SMatthew G. Knepley   SolType        solType;      /* MMS solution type */
58a4e36a32SMatthew G. Knepley   PartLayoutType partLayout;   /* Type of particle distribution */
59a4e36a32SMatthew G. Knepley   PetscInt       Npc;          /* The initial number of particles per cell */
60a4e36a32SMatthew G. Knepley   PetscReal      partLower[3]; /* Lower left corner of particle box */
61a4e36a32SMatthew G. Knepley   PetscReal      partUpper[3]; /* Upper right corner of particle box */
62a4e36a32SMatthew G. Knepley   PetscInt       Npb;          /* The initial number of particles per box dimension */
63a4e36a32SMatthew G. Knepley } AppCtx;
64a4e36a32SMatthew G. Knepley 
65a4e36a32SMatthew G. Knepley typedef struct {
66a4e36a32SMatthew G. Knepley   PetscReal ti; /* The time for ui, at the beginning of the advection solve */
67a4e36a32SMatthew G. Knepley   PetscReal tf; /* The time for uf, at the end of the advection solve */
68a4e36a32SMatthew G. Knepley   Vec       ui; /* The PDE solution field at ti */
69a4e36a32SMatthew G. Knepley   Vec       uf; /* The PDE solution field at tf */
70a4e36a32SMatthew G. Knepley   Vec       x0; /* The initial particle positions at t = 0 */
71a4e36a32SMatthew G. Knepley   PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
72a4e36a32SMatthew G. Knepley   AppCtx *ctx; /* Context for exact solution */
73a4e36a32SMatthew G. Knepley } AdvCtx;
74a4e36a32SMatthew G. Knepley 
zero(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)75*2a8381b2SBarry Smith static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
76d71ae5a4SJacob Faibussowitsch {
77a4e36a32SMatthew G. Knepley   PetscInt d;
78a4e36a32SMatthew G. Knepley   for (d = 0; d < Nc; ++d) u[d] = 0.0;
793ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
80a4e36a32SMatthew G. Knepley }
81a4e36a32SMatthew G. Knepley 
constant(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)82*2a8381b2SBarry Smith static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
83d71ae5a4SJacob Faibussowitsch {
84a4e36a32SMatthew G. Knepley   PetscInt d;
85a4e36a32SMatthew G. Knepley   for (d = 0; d < Nc; ++d) u[d] = 1.0;
863ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
87a4e36a32SMatthew G. Knepley }
88a4e36a32SMatthew G. Knepley 
89a4e36a32SMatthew G. Knepley /*
90a4e36a32SMatthew G. Knepley   CASE: trigonometric-trigonometric
91a4e36a32SMatthew G. Knepley   In 2D we use exact solution:
92a4e36a32SMatthew G. Knepley 
93a4e36a32SMatthew G. Knepley     x = r0 cos(w t + theta0)  r0     = sqrt(x0^2 + y0^2)
94a4e36a32SMatthew G. Knepley     y = r0 sin(w t + theta0)  theta0 = arctan(y0/x0)
95a4e36a32SMatthew G. Knepley     u = -w r0 sin(theta0) = -w y
96a4e36a32SMatthew G. Knepley     v =  w r0 cos(theta0) =  w x
97a4e36a32SMatthew G. Knepley     p = x + y - 1
98a4e36a32SMatthew G. Knepley     T = t + x + y
99a4e36a32SMatthew G. Knepley     f = <1, 1>
100a4e36a32SMatthew G. Knepley     Q = 1 + w (x - y)/r
101a4e36a32SMatthew G. Knepley 
102a4e36a32SMatthew G. Knepley   so that
103a4e36a32SMatthew G. Knepley 
104a4e36a32SMatthew G. Knepley     \nabla \cdot u = 0 + 0 = 0
105a4e36a32SMatthew G. Knepley 
106a4e36a32SMatthew G. Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
107a4e36a32SMatthew G. Knepley     = <0, 0> + u_i d_i u_j - \nu 0 + <1, 1>
108a4e36a32SMatthew G. Knepley     = <1, 1> + w^2 <-y, x> . <<0, 1>, <-1, 0>>
109a4e36a32SMatthew G. Knepley     = <1, 1> + w^2 <-x, -y>
110a4e36a32SMatthew G. Knepley     = <1, 1> - w^2 <x, y>
111a4e36a32SMatthew G. Knepley 
112a4e36a32SMatthew G. Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
113a4e36a32SMatthew G. Knepley     = 1 + <u, v> . <1, 1> - \alpha 0
114a4e36a32SMatthew G. Knepley     = 1 + u + v
115a4e36a32SMatthew G. Knepley */
trig_trig_x(PetscInt dim,PetscReal time,const PetscReal X[],PetscInt Nf,PetscScalar * x,PetscCtx ctx)116*2a8381b2SBarry Smith static PetscErrorCode trig_trig_x(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *x, PetscCtx ctx)
117d71ae5a4SJacob Faibussowitsch {
118a4e36a32SMatthew G. Knepley   const PetscReal x0     = X[0];
119a4e36a32SMatthew G. Knepley   const PetscReal y0     = X[1];
120a4e36a32SMatthew G. Knepley   const PetscReal R0     = PetscSqrtReal(x0 * x0 + y0 * y0);
121a4e36a32SMatthew G. Knepley   const PetscReal theta0 = PetscAtan2Real(y0, x0);
122a4e36a32SMatthew G. Knepley   Parameter      *p      = (Parameter *)ctx;
123a4e36a32SMatthew G. Knepley 
124a4e36a32SMatthew G. Knepley   x[0] = R0 * PetscCosReal(p->omega * time + theta0);
125a4e36a32SMatthew G. Knepley   x[1] = R0 * PetscSinReal(p->omega * time + theta0);
1263ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
127a4e36a32SMatthew G. Knepley }
trig_trig_u(PetscInt dim,PetscReal time,const PetscReal X[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)128*2a8381b2SBarry Smith static PetscErrorCode trig_trig_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
129d71ae5a4SJacob Faibussowitsch {
130a4e36a32SMatthew G. Knepley   Parameter *p = (Parameter *)ctx;
131a4e36a32SMatthew G. Knepley 
132a4e36a32SMatthew G. Knepley   u[0] = -p->omega * X[1];
133a4e36a32SMatthew G. Knepley   u[1] = p->omega * X[0];
1343ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
135a4e36a32SMatthew G. Knepley }
trig_trig_u_t(PetscInt dim,PetscReal time,const PetscReal X[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)136*2a8381b2SBarry Smith static PetscErrorCode trig_trig_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
137d71ae5a4SJacob Faibussowitsch {
138a4e36a32SMatthew G. Knepley   u[0] = 0.0;
139a4e36a32SMatthew G. Knepley   u[1] = 0.0;
1403ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
141a4e36a32SMatthew G. Knepley }
142a4e36a32SMatthew G. Knepley 
trig_trig_p(PetscInt dim,PetscReal time,const PetscReal X[],PetscInt Nf,PetscScalar * p,PetscCtx ctx)143*2a8381b2SBarry Smith static PetscErrorCode trig_trig_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, PetscCtx ctx)
144d71ae5a4SJacob Faibussowitsch {
145a4e36a32SMatthew G. Knepley   p[0] = X[0] + X[1] - 1.0;
1463ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
147a4e36a32SMatthew G. Knepley }
148a4e36a32SMatthew G. Knepley 
trig_trig_T(PetscInt dim,PetscReal time,const PetscReal X[],PetscInt Nf,PetscScalar * T,PetscCtx ctx)149*2a8381b2SBarry Smith static PetscErrorCode trig_trig_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, PetscCtx ctx)
150d71ae5a4SJacob Faibussowitsch {
151a4e36a32SMatthew G. Knepley   T[0] = time + X[0] + X[1];
1523ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
153a4e36a32SMatthew G. Knepley }
trig_trig_T_t(PetscInt dim,PetscReal time,const PetscReal X[],PetscInt Nf,PetscScalar * T,PetscCtx ctx)154*2a8381b2SBarry Smith static PetscErrorCode trig_trig_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, PetscCtx ctx)
155d71ae5a4SJacob Faibussowitsch {
156a4e36a32SMatthew G. Knepley   T[0] = 1.0;
1573ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
158a4e36a32SMatthew G. Knepley }
159a4e36a32SMatthew G. Knepley 
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[])160d71ae5a4SJacob Faibussowitsch 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[])
161d71ae5a4SJacob Faibussowitsch {
162a4e36a32SMatthew G. Knepley   const PetscReal omega = PetscRealPart(constants[3]);
163a4e36a32SMatthew G. Knepley   PetscInt        Nc    = dim;
164a4e36a32SMatthew G. Knepley   PetscInt        c, d;
165a4e36a32SMatthew G. Knepley 
166a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0] + d];
167a4e36a32SMatthew G. Knepley 
168a4e36a32SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
169a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
170a4e36a32SMatthew G. Knepley   }
171a4e36a32SMatthew G. Knepley   f0[0] -= 1.0 - omega * omega * X[0];
172a4e36a32SMatthew G. Knepley   f0[1] -= 1.0 - omega * omega * X[1];
173a4e36a32SMatthew G. Knepley }
174a4e36a32SMatthew G. Knepley 
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[])175d71ae5a4SJacob Faibussowitsch 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[])
176d71ae5a4SJacob Faibussowitsch {
177a4e36a32SMatthew G. Knepley   const PetscReal omega = PetscRealPart(constants[3]);
178a4e36a32SMatthew G. Knepley   PetscInt        d;
179a4e36a32SMatthew G. Knepley 
180a4e36a32SMatthew G. Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0] + d] * u_x[uOff_x[2] + d];
181a4e36a32SMatthew G. Knepley   f0[0] += u_t[uOff[2]] - (1.0 + omega * (X[0] - X[1]));
182a4e36a32SMatthew G. Knepley }
183a4e36a32SMatthew G. Knepley 
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[])184d71ae5a4SJacob Faibussowitsch 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[])
185d71ae5a4SJacob Faibussowitsch {
186a4e36a32SMatthew G. Knepley   PetscInt d;
187a4e36a32SMatthew G. Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d];
188a4e36a32SMatthew G. Knepley }
189a4e36a32SMatthew G. Knepley 
190a4e36a32SMatthew G. Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */
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[])191d71ae5a4SJacob Faibussowitsch 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[])
192d71ae5a4SJacob Faibussowitsch {
193a4e36a32SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[0]);
194a4e36a32SMatthew G. Knepley   const PetscInt  Nc = dim;
195a4e36a32SMatthew G. Knepley   PetscInt        c, d;
196a4e36a32SMatthew G. Knepley 
197a4e36a32SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
198ad540459SPierre Jolivet     for (d = 0; d < dim; ++d) f1[c * dim + d] = nu * (u_x[c * dim + d] + u_x[d * dim + c]);
199a4e36a32SMatthew G. Knepley     f1[c * dim + c] -= u[uOff[1]];
200a4e36a32SMatthew G. Knepley   }
201a4e36a32SMatthew G. Knepley }
202a4e36a32SMatthew G. Knepley 
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[])203d71ae5a4SJacob Faibussowitsch 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[])
204d71ae5a4SJacob Faibussowitsch {
205a4e36a32SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
206a4e36a32SMatthew G. Knepley   PetscInt        d;
207a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha * u_x[uOff_x[2] + d];
208a4e36a32SMatthew G. Knepley }
209a4e36a32SMatthew G. Knepley 
210a4e36a32SMatthew G. Knepley /* Jacobians */
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[])211d71ae5a4SJacob Faibussowitsch 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[])
212d71ae5a4SJacob Faibussowitsch {
213a4e36a32SMatthew G. Knepley   PetscInt d;
214a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
215a4e36a32SMatthew G. Knepley }
216a4e36a32SMatthew G. Knepley 
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[])217d71ae5a4SJacob Faibussowitsch 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[])
218d71ae5a4SJacob Faibussowitsch {
219a4e36a32SMatthew G. Knepley   PetscInt       c, d;
220a4e36a32SMatthew G. Knepley   const PetscInt Nc = dim;
221a4e36a32SMatthew G. Knepley 
222a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift;
223a4e36a32SMatthew G. Knepley 
224a4e36a32SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
225ad540459SPierre Jolivet     for (d = 0; d < dim; ++d) g0[c * Nc + d] += u_x[c * Nc + d];
226a4e36a32SMatthew G. Knepley   }
227a4e36a32SMatthew G. Knepley }
228a4e36a32SMatthew G. Knepley 
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[])229d71ae5a4SJacob Faibussowitsch 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[])
230d71ae5a4SJacob Faibussowitsch {
231a4e36a32SMatthew G. Knepley   PetscInt NcI = dim;
232a4e36a32SMatthew G. Knepley   PetscInt NcJ = dim;
233a4e36a32SMatthew G. Knepley   PetscInt c, d, e;
234a4e36a32SMatthew G. Knepley 
235a4e36a32SMatthew G. Knepley   for (c = 0; c < NcI; ++c) {
236a4e36a32SMatthew G. Knepley     for (d = 0; d < NcJ; ++d) {
237a4e36a32SMatthew G. Knepley       for (e = 0; e < dim; ++e) {
238ad540459SPierre Jolivet         if (c == d) g1[(c * NcJ + d) * dim + e] += u[e];
239a4e36a32SMatthew G. Knepley       }
240a4e36a32SMatthew G. Knepley     }
241a4e36a32SMatthew G. Knepley   }
242a4e36a32SMatthew G. Knepley }
243a4e36a32SMatthew G. Knepley 
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[])244d71ae5a4SJacob Faibussowitsch 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[])
245d71ae5a4SJacob Faibussowitsch {
246a4e36a32SMatthew G. Knepley   PetscInt d;
247a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0;
248a4e36a32SMatthew G. Knepley }
249a4e36a32SMatthew G. Knepley 
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[])250d71ae5a4SJacob Faibussowitsch 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[])
251d71ae5a4SJacob Faibussowitsch {
252a4e36a32SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[0]);
253a4e36a32SMatthew G. Knepley   const PetscInt  Nc = dim;
254a4e36a32SMatthew G. Knepley   PetscInt        c, d;
255a4e36a32SMatthew G. Knepley 
256a4e36a32SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
257a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
258a4e36a32SMatthew G. Knepley       g3[((c * Nc + c) * dim + d) * dim + d] += nu; // gradU
259a4e36a32SMatthew G. Knepley       g3[((c * Nc + d) * dim + d) * dim + c] += nu; // gradU transpose
260a4e36a32SMatthew G. Knepley     }
261a4e36a32SMatthew G. Knepley   }
262a4e36a32SMatthew G. Knepley }
263a4e36a32SMatthew G. Knepley 
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[])264d71ae5a4SJacob Faibussowitsch 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[])
265d71ae5a4SJacob Faibussowitsch {
266a4e36a32SMatthew G. Knepley   PetscInt d;
267a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g0[d] = u_tShift;
268a4e36a32SMatthew G. Knepley }
269a4e36a32SMatthew G. Knepley 
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[])270d71ae5a4SJacob Faibussowitsch 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[])
271d71ae5a4SJacob Faibussowitsch {
272a4e36a32SMatthew G. Knepley   PetscInt d;
273a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2] + d];
274a4e36a32SMatthew G. Knepley }
275a4e36a32SMatthew G. Knepley 
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[])276d71ae5a4SJacob Faibussowitsch 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[])
277d71ae5a4SJacob Faibussowitsch {
278a4e36a32SMatthew G. Knepley   PetscInt d;
279a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0] + d];
280a4e36a32SMatthew G. Knepley }
281a4e36a32SMatthew G. Knepley 
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[])282d71ae5a4SJacob Faibussowitsch 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[])
283d71ae5a4SJacob Faibussowitsch {
284a4e36a32SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
285a4e36a32SMatthew G. Knepley   PetscInt        d;
286a4e36a32SMatthew G. Knepley 
287a4e36a32SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d * dim + d] = alpha;
288a4e36a32SMatthew G. Knepley }
289a4e36a32SMatthew G. Knepley 
ProcessOptions(MPI_Comm comm,AppCtx * options)290d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
291d71ae5a4SJacob Faibussowitsch {
292a4e36a32SMatthew G. Knepley   PetscInt sol, pl, n;
293a4e36a32SMatthew G. Knepley 
294a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
295a4e36a32SMatthew G. Knepley   options->solType      = SOL_TRIG_TRIG;
296a4e36a32SMatthew G. Knepley   options->partLayout   = PART_LAYOUT_CELL;
297a4e36a32SMatthew G. Knepley   options->Npc          = 1;
298a4e36a32SMatthew G. Knepley   options->Npb          = 1;
299a4e36a32SMatthew G. Knepley   options->partLower[0] = options->partLower[1] = options->partLower[2] = 0.;
300a4e36a32SMatthew G. Knepley   options->partUpper[0] = options->partUpper[1] = options->partUpper[2] = 1.;
301d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");
302a4e36a32SMatthew G. Knepley   sol = options->solType;
3039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL));
304a4e36a32SMatthew G. Knepley   options->solType = (SolType)sol;
305a4e36a32SMatthew G. Knepley   pl               = options->partLayout;
3069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL));
307a4e36a32SMatthew G. Knepley   options->partLayout = (PartLayoutType)pl;
3089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL));
309a4e36a32SMatthew G. Knepley   n = 3;
3109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL));
311a4e36a32SMatthew G. Knepley   n = 3;
3129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL));
3139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL));
314d0609cedSBarry Smith   PetscOptionsEnd();
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
316a4e36a32SMatthew G. Knepley }
317a4e36a32SMatthew G. Knepley 
SetupParameters(AppCtx * user)318d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(AppCtx *user)
319d71ae5a4SJacob Faibussowitsch {
320a4e36a32SMatthew G. Knepley   PetscBag   bag;
321a4e36a32SMatthew G. Knepley   Parameter *p;
322a4e36a32SMatthew G. Knepley 
323a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
324a4e36a32SMatthew G. Knepley   /* setup PETSc parameter bag */
325*2a8381b2SBarry Smith   PetscCall(PetscBagGetData(user->bag, &p));
3269566063dSJacob Faibussowitsch   PetscCall(PetscBagSetName(user->bag, "par", "Low Mach flow parameters"));
327a4e36a32SMatthew G. Knepley   bag = user->bag;
3289566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity"));
3299566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity"));
3309566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature"));
3319566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark"));
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
333a4e36a32SMatthew G. Knepley }
334a4e36a32SMatthew G. Knepley 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)335d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
336d71ae5a4SJacob Faibussowitsch {
337a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
3389566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
3399566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
3409566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
3419566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
343a4e36a32SMatthew G. Knepley }
344a4e36a32SMatthew G. Knepley 
SetupProblem(DM dm,AppCtx * user)345d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
346d71ae5a4SJacob Faibussowitsch {
347*2a8381b2SBarry Smith   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
348*2a8381b2SBarry Smith   PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
349a4e36a32SMatthew G. Knepley   PetscDS    prob;
35045480ffeSMatthew G. Knepley   DMLabel    label;
351a4e36a32SMatthew G. Knepley   Parameter *ctx;
352a4e36a32SMatthew G. Knepley   PetscInt   id;
353a4e36a32SMatthew G. Knepley 
354a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
3559566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
3569566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
357a4e36a32SMatthew G. Knepley   switch (user->solType) {
358a4e36a32SMatthew G. Knepley   case SOL_TRIG_TRIG:
3599566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v));
3609566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w));
361a4e36a32SMatthew G. Knepley 
362a4e36a32SMatthew G. Knepley     exactFuncs[0]   = trig_trig_u;
363a4e36a32SMatthew G. Knepley     exactFuncs[1]   = trig_trig_p;
364a4e36a32SMatthew G. Knepley     exactFuncs[2]   = trig_trig_T;
365a4e36a32SMatthew G. Knepley     exactFuncs_t[0] = trig_trig_u_t;
366a4e36a32SMatthew G. Knepley     exactFuncs_t[1] = NULL;
367a4e36a32SMatthew G. Knepley     exactFuncs_t[2] = trig_trig_T_t;
368a4e36a32SMatthew G. Knepley     break;
369d71ae5a4SJacob Faibussowitsch   default:
370d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
371a4e36a32SMatthew G. Knepley   }
372a4e36a32SMatthew G. Knepley 
3739566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL));
374a4e36a32SMatthew G. Knepley 
3759566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu));
3769566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL));
3779566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL));
3789566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL));
3799566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT));
380a4e36a32SMatthew G. Knepley   /* Setup constants */
381a4e36a32SMatthew G. Knepley   {
382a4e36a32SMatthew G. Knepley     Parameter  *param;
383a4e36a32SMatthew G. Knepley     PetscScalar constants[4];
384a4e36a32SMatthew G. Knepley 
385*2a8381b2SBarry Smith     PetscCall(PetscBagGetData(user->bag, &param));
386a4e36a32SMatthew G. Knepley 
387a4e36a32SMatthew G. Knepley     constants[0] = param->nu;
388a4e36a32SMatthew G. Knepley     constants[1] = param->alpha;
389a4e36a32SMatthew G. Knepley     constants[2] = param->T_in;
390a4e36a32SMatthew G. Knepley     constants[3] = param->omega;
3919566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(prob, 4, constants));
392a4e36a32SMatthew G. Knepley   }
393a4e36a32SMatthew G. Knepley   /* Setup Boundary Conditions */
394*2a8381b2SBarry Smith   PetscCall(PetscBagGetData(user->bag, &ctx));
395a4e36a32SMatthew G. Knepley   id = 3;
39657d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)exactFuncs[0], (PetscVoidFn *)exactFuncs_t[0], ctx, NULL));
397a4e36a32SMatthew G. Knepley   id = 1;
39857d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)exactFuncs[0], (PetscVoidFn *)exactFuncs_t[0], ctx, NULL));
399a4e36a32SMatthew G. Knepley   id = 2;
40057d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)exactFuncs[0], (PetscVoidFn *)exactFuncs_t[0], ctx, NULL));
401a4e36a32SMatthew G. Knepley   id = 4;
40257d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)exactFuncs[0], (PetscVoidFn *)exactFuncs_t[0], ctx, NULL));
403a4e36a32SMatthew G. Knepley   id = 3;
40457d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)exactFuncs[2], (PetscVoidFn *)exactFuncs_t[2], ctx, NULL));
405a4e36a32SMatthew G. Knepley   id = 1;
40657d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)exactFuncs[2], (PetscVoidFn *)exactFuncs_t[2], ctx, NULL));
407a4e36a32SMatthew G. Knepley   id = 2;
40857d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)exactFuncs[2], (PetscVoidFn *)exactFuncs_t[2], ctx, NULL));
409a4e36a32SMatthew G. Knepley   id = 4;
41057d50842SBarry Smith   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)exactFuncs[2], (PetscVoidFn *)exactFuncs_t[2], ctx, NULL));
411a4e36a32SMatthew G. Knepley 
412a4e36a32SMatthew G. Knepley   /*setup exact solution.*/
4139566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx));
4149566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx));
4159566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx));
4169566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx));
4179566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx));
4189566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx));
4193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
420a4e36a32SMatthew G. Knepley }
421a4e36a32SMatthew G. Knepley 
422a4e36a32SMatthew G. Knepley /* x_t = v
423a4e36a32SMatthew G. Knepley 
424a4e36a32SMatthew G. Knepley    Note that here we use the velocity field at t_{n+1} to advect the particles from
425a4e36a32SMatthew G. Knepley    t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or
426a4e36a32SMatthew G. Knepley    the method of characteristics.
427a4e36a32SMatthew G. Knepley */
FreeStreaming(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)428*2a8381b2SBarry Smith static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
429d71ae5a4SJacob Faibussowitsch {
430a4e36a32SMatthew G. Knepley   AdvCtx             *adv = (AdvCtx *)ctx;
431a4e36a32SMatthew G. Knepley   Vec                 u   = adv->ui;
432a4e36a32SMatthew G. Knepley   DM                  sdm, dm, vdm;
433a4e36a32SMatthew G. Knepley   Vec                 vel, locvel, pvel;
434a4e36a32SMatthew G. Knepley   IS                  vis;
435a4e36a32SMatthew G. Knepley   DMInterpolationInfo ictx;
436a4e36a32SMatthew G. Knepley   const PetscScalar  *coords, *v;
437a4e36a32SMatthew G. Knepley   PetscScalar        *f;
438a4e36a32SMatthew G. Knepley   PetscInt            vf[1] = {0};
439a4e36a32SMatthew G. Knepley   PetscInt            dim, Np;
440a4e36a32SMatthew G. Knepley 
441a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
4429566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sdm));
4439566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sdm, &dm));
4449566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(sdm, &pvel));
4459566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(sdm, &Np));
4469566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
447a4e36a32SMatthew G. Knepley   /* Get local velocity */
4489566063dSJacob Faibussowitsch   PetscCall(DMCreateSubDM(dm, 1, vf, &vis, &vdm));
4499566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(u, vis, &vel));
4509566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(vdm, &locvel));
4519566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL));
4529566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel));
4539566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel));
4549566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(u, vis, &vel));
4559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&vis));
456a4e36a32SMatthew G. Knepley   /* Interpolate velocity */
4579566063dSJacob Faibussowitsch   PetscCall(DMInterpolationCreate(PETSC_COMM_SELF, &ictx));
4589566063dSJacob Faibussowitsch   PetscCall(DMInterpolationSetDim(ictx, dim));
4599566063dSJacob Faibussowitsch   PetscCall(DMInterpolationSetDof(ictx, dim));
4609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &coords));
4619566063dSJacob Faibussowitsch   PetscCall(DMInterpolationAddPoints(ictx, Np, (PetscReal *)coords));
4629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &coords));
46352aa1562SMatthew G. Knepley   /* Particles that lie outside the domain should be dropped,
46452aa1562SMatthew G. Knepley      whereas particles that move to another partition should trigger a migration */
4659566063dSJacob Faibussowitsch   PetscCall(DMInterpolationSetUp(ictx, vdm, PETSC_FALSE, PETSC_TRUE));
4669566063dSJacob Faibussowitsch   PetscCall(VecSet(pvel, 0.));
4679566063dSJacob Faibussowitsch   PetscCall(DMInterpolationEvaluate(ictx, vdm, locvel, pvel));
4689566063dSJacob Faibussowitsch   PetscCall(DMInterpolationDestroy(&ictx));
4699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(vdm, &locvel));
4709566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&vdm));
471a4e36a32SMatthew G. Knepley 
4729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
4739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pvel, &v));
4749566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(f, v, Np * dim));
4759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pvel, &v));
4769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
4779566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(sdm, &pvel));
4783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
479a4e36a32SMatthew G. Knepley }
480a4e36a32SMatthew G. Knepley 
SetInitialParticleConditions(TS ts,Vec u)481d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u)
482d71ae5a4SJacob Faibussowitsch {
483a4e36a32SMatthew G. Knepley   AppCtx      *user;
484a4e36a32SMatthew G. Knepley   void        *ctx;
485a4e36a32SMatthew G. Knepley   DM           dm;
486a4e36a32SMatthew G. Knepley   PetscScalar *coords;
487a4e36a32SMatthew G. Knepley   PetscReal    x[3], dx[3];
488a4e36a32SMatthew G. Knepley   PetscInt     n[3];
4895f80ce2aSJacob Faibussowitsch   PetscInt     dim, d, i, j, k;
490a4e36a32SMatthew G. Knepley 
4917510d9b0SBarry Smith   PetscFunctionBeginUser;
4929566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &ctx));
493a4e36a32SMatthew G. Knepley   user = ((AdvCtx *)ctx)->ctx;
4949566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
4959566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
496a4e36a32SMatthew G. Knepley   switch (user->partLayout) {
497d71ae5a4SJacob Faibussowitsch   case PART_LAYOUT_CELL:
498d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmSetPointCoordinatesRandom(dm, user->Npc));
499d71ae5a4SJacob Faibussowitsch     break;
500a4e36a32SMatthew G. Knepley   case PART_LAYOUT_BOX:
501a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
502a4e36a32SMatthew G. Knepley       n[d]  = user->Npb;
503a4e36a32SMatthew G. Knepley       dx[d] = (user->partUpper[d] - user->partLower[d]) / PetscMax(1, n[d] - 1);
504a4e36a32SMatthew G. Knepley     }
5059566063dSJacob Faibussowitsch     PetscCall(VecGetArray(u, &coords));
506a4e36a32SMatthew G. Knepley     switch (dim) {
507a4e36a32SMatthew G. Knepley     case 2:
508a4e36a32SMatthew G. Knepley       x[0] = user->partLower[0];
509a4e36a32SMatthew G. Knepley       for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
510a4e36a32SMatthew G. Knepley         x[1] = user->partLower[1];
511a4e36a32SMatthew G. Knepley         for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
512a4e36a32SMatthew G. Knepley           const PetscInt p = j * n[0] + i;
513a4e36a32SMatthew G. Knepley           for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d];
514a4e36a32SMatthew G. Knepley         }
515a4e36a32SMatthew G. Knepley       }
516a4e36a32SMatthew G. Knepley       break;
517a4e36a32SMatthew G. Knepley     case 3:
518a4e36a32SMatthew G. Knepley       x[0] = user->partLower[0];
519a4e36a32SMatthew G. Knepley       for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
520a4e36a32SMatthew G. Knepley         x[1] = user->partLower[1];
521a4e36a32SMatthew G. Knepley         for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
522a4e36a32SMatthew G. Knepley           x[2] = user->partLower[2];
523a4e36a32SMatthew G. Knepley           for (k = 0; k < n[2]; ++k, x[2] += dx[2]) {
524a4e36a32SMatthew G. Knepley             const PetscInt p = (k * n[1] + j) * n[0] + i;
525a4e36a32SMatthew G. Knepley             for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d];
526a4e36a32SMatthew G. Knepley           }
527a4e36a32SMatthew G. Knepley         }
528a4e36a32SMatthew G. Knepley       }
529a4e36a32SMatthew G. Knepley       break;
530d71ae5a4SJacob Faibussowitsch     default:
531d71ae5a4SJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %" PetscInt_FMT, dim);
532a4e36a32SMatthew G. Knepley     }
5339566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(u, &coords));
534a4e36a32SMatthew G. Knepley     break;
535d71ae5a4SJacob Faibussowitsch   default:
536d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]);
537a4e36a32SMatthew G. Knepley   }
5383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
539a4e36a32SMatthew G. Knepley }
540a4e36a32SMatthew G. Knepley 
SetupDiscretization(DM dm,DM sdm,AppCtx * user)541d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user)
542d71ae5a4SJacob Faibussowitsch {
543a4e36a32SMatthew G. Knepley   DM             cdm = dm;
54419307e5cSMatthew G. Knepley   DMSwarmCellDM  celldm;
545a4e36a32SMatthew G. Knepley   PetscFE        fe[3];
546a4e36a32SMatthew G. Knepley   Parameter     *param;
54719307e5cSMatthew G. Knepley   PetscInt      *swarm_cellid, n[3];
548a4e36a32SMatthew G. Knepley   PetscReal      x[3], dx[3];
549a4e36a32SMatthew G. Knepley   PetscScalar   *coords;
550a4e36a32SMatthew G. Knepley   DMPolytopeType ct;
551a4e36a32SMatthew G. Knepley   PetscInt       dim, d, cStart, cEnd, c, Np, p, i, j, k;
552a4e36a32SMatthew G. Knepley   PetscBool      simplex;
553a4e36a32SMatthew G. Knepley   MPI_Comm       comm;
55419307e5cSMatthew G. Knepley   const char    *cellid;
555a4e36a32SMatthew G. Knepley 
556a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
5579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
5599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
560a4e36a32SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
561a4e36a32SMatthew G. Knepley   /* Create finite element */
5629566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5639566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
5649566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
565a4e36a32SMatthew G. Knepley 
5669566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
5679566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
5689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
569a4e36a32SMatthew G. Knepley 
5709566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]));
5719566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
5729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[2], "temperature"));
573a4e36a32SMatthew G. Knepley 
574a4e36a32SMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
5759566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
5769566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
5779566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe[2]));
5789566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
5799566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, user));
580*2a8381b2SBarry Smith   PetscCall(PetscBagGetData(user->bag, &param));
581a4e36a32SMatthew G. Knepley   while (cdm) {
5829566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
5839566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
584a4e36a32SMatthew G. Knepley   }
5859566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[0]));
5869566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[1]));
5879566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[2]));
588a4e36a32SMatthew G. Knepley 
589a4e36a32SMatthew G. Knepley   {
590a4e36a32SMatthew G. Knepley     PetscObject  pressure;
591a4e36a32SMatthew G. Knepley     MatNullSpace nullspacePres;
592a4e36a32SMatthew G. Knepley 
5939566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, 1, NULL, &pressure));
5949566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
5959566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nullspacePres));
5969566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullspacePres));
597a4e36a32SMatthew G. Knepley   }
598a4e36a32SMatthew G. Knepley 
599a4e36a32SMatthew G. Knepley   /* Setup particle information */
6009566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(sdm, DMSWARM_PIC));
6019566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL));
6029566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(sdm));
60319307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(sdm, &celldm));
60419307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
605a4e36a32SMatthew G. Knepley   switch (user->partLayout) {
606a4e36a32SMatthew G. Knepley   case PART_LAYOUT_CELL:
6079566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0));
6089566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(sdm));
60919307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(sdm, cellid, NULL, NULL, (void **)&swarm_cellid));
610a4e36a32SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
611a4e36a32SMatthew G. Knepley       for (p = 0; p < user->Npc; ++p) {
612a4e36a32SMatthew G. Knepley         const PetscInt n = c * user->Npc + p;
613a4e36a32SMatthew G. Knepley 
61419307e5cSMatthew G. Knepley         swarm_cellid[n] = c;
615a4e36a32SMatthew G. Knepley       }
616a4e36a32SMatthew G. Knepley     }
61719307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sdm, cellid, NULL, NULL, (void **)&swarm_cellid));
6189566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetPointCoordinatesRandom(sdm, user->Npc));
619a4e36a32SMatthew G. Knepley     break;
620a4e36a32SMatthew G. Knepley   case PART_LAYOUT_BOX:
621a4e36a32SMatthew G. Knepley     Np = 1;
622a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
623a4e36a32SMatthew G. Knepley       n[d]  = user->Npb;
624a4e36a32SMatthew G. Knepley       dx[d] = (user->partUpper[d] - user->partLower[d]) / PetscMax(1, n[d] - 1);
625a4e36a32SMatthew G. Knepley       Np *= n[d];
626a4e36a32SMatthew G. Knepley     }
6279566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetLocalSizes(sdm, Np, 0));
6289566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(sdm));
6299566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
630a4e36a32SMatthew G. Knepley     switch (dim) {
631a4e36a32SMatthew G. Knepley     case 2:
632a4e36a32SMatthew G. Knepley       x[0] = user->partLower[0];
633a4e36a32SMatthew G. Knepley       for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
634a4e36a32SMatthew G. Knepley         x[1] = user->partLower[1];
635a4e36a32SMatthew G. Knepley         for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
636a4e36a32SMatthew G. Knepley           const PetscInt p = j * n[0] + i;
637a4e36a32SMatthew G. Knepley           for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d];
638a4e36a32SMatthew G. Knepley         }
639a4e36a32SMatthew G. Knepley       }
640a4e36a32SMatthew G. Knepley       break;
641a4e36a32SMatthew G. Knepley     case 3:
642a4e36a32SMatthew G. Knepley       x[0] = user->partLower[0];
643a4e36a32SMatthew G. Knepley       for (i = 0; i < n[0]; ++i, x[0] += dx[0]) {
644a4e36a32SMatthew G. Knepley         x[1] = user->partLower[1];
645a4e36a32SMatthew G. Knepley         for (j = 0; j < n[1]; ++j, x[1] += dx[1]) {
646a4e36a32SMatthew G. Knepley           x[2] = user->partLower[2];
647a4e36a32SMatthew G. Knepley           for (k = 0; k < n[2]; ++k, x[2] += dx[2]) {
648a4e36a32SMatthew G. Knepley             const PetscInt p = (k * n[1] + j) * n[0] + i;
649a4e36a32SMatthew G. Knepley             for (d = 0; d < dim; ++d) coords[p * dim + d] = x[d];
650a4e36a32SMatthew G. Knepley           }
651a4e36a32SMatthew G. Knepley         }
652a4e36a32SMatthew G. Knepley       }
653a4e36a32SMatthew G. Knepley       break;
654d71ae5a4SJacob Faibussowitsch     default:
655d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %" PetscInt_FMT, dim);
656a4e36a32SMatthew G. Knepley     }
6579566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
65819307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(sdm, cellid, NULL, NULL, (void **)&swarm_cellid));
65919307e5cSMatthew G. Knepley     for (p = 0; p < Np; ++p) swarm_cellid[p] = 0;
66019307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(sdm, cellid, NULL, NULL, (void **)&swarm_cellid));
6619566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(sdm, PETSC_TRUE));
662a4e36a32SMatthew G. Knepley     break;
663d71ae5a4SJacob Faibussowitsch   default:
664d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]);
665a4e36a32SMatthew G. Knepley   }
6669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sdm, "Particles"));
6679566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view"));
6683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
669a4e36a32SMatthew G. Knepley }
670a4e36a32SMatthew G. Knepley 
CreatePressureNullSpace(DM dm,PetscInt ofield,PetscInt nfield,MatNullSpace * nullSpace)671d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
672d71ae5a4SJacob Faibussowitsch {
673a4e36a32SMatthew G. Knepley   Vec vec;
674a4e36a32SMatthew G. Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
675a4e36a32SMatthew G. Knepley 
676a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
67763a3b9bcSJacob Faibussowitsch   PetscCheck(ofield == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %" PetscInt_FMT, ofield);
678a4e36a32SMatthew G. Knepley   funcs[nfield] = constant;
6799566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &vec));
6809566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
6819566063dSJacob Faibussowitsch   PetscCall(VecNormalize(vec, NULL));
6829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)vec, "Pressure Null Space"));
6839566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view"));
6849566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace));
6859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vec));
6863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
687a4e36a32SMatthew G. Knepley }
688a4e36a32SMatthew G. Knepley 
RemoveDiscretePressureNullspace_Private(TS ts,Vec u)689d71ae5a4SJacob Faibussowitsch static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u)
690d71ae5a4SJacob Faibussowitsch {
691a4e36a32SMatthew G. Knepley   DM           dm;
692a4e36a32SMatthew G. Knepley   MatNullSpace nullsp;
693a4e36a32SMatthew G. Knepley 
6947510d9b0SBarry Smith   PetscFunctionBeginUser;
6959566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
6969566063dSJacob Faibussowitsch   PetscCall(CreatePressureNullSpace(dm, 1, 1, &nullsp));
6979566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceRemove(nullsp, u));
6989566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&nullsp));
6993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
700a4e36a32SMatthew G. Knepley }
701a4e36a32SMatthew G. Knepley 
702a4e36a32SMatthew G. Knepley /* Make the discrete pressure discretely divergence free */
RemoveDiscretePressureNullspace(TS ts)703d71ae5a4SJacob Faibussowitsch static PetscErrorCode RemoveDiscretePressureNullspace(TS ts)
704d71ae5a4SJacob Faibussowitsch {
705a4e36a32SMatthew G. Knepley   Vec u;
706a4e36a32SMatthew G. Knepley 
7077510d9b0SBarry Smith   PetscFunctionBeginUser;
7089566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &u));
7099566063dSJacob Faibussowitsch   PetscCall(RemoveDiscretePressureNullspace_Private(ts, u));
7103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
711a4e36a32SMatthew G. Knepley }
712a4e36a32SMatthew G. Knepley 
SetInitialConditions(TS ts,Vec u)713d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditions(TS ts, Vec u)
714d71ae5a4SJacob Faibussowitsch {
715a4e36a32SMatthew G. Knepley   DM        dm;
716a4e36a32SMatthew G. Knepley   PetscReal t;
717a4e36a32SMatthew G. Knepley 
7187510d9b0SBarry Smith   PetscFunctionBeginUser;
7199566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
7209566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
7219566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
7229566063dSJacob Faibussowitsch   PetscCall(RemoveDiscretePressureNullspace_Private(ts, u));
7233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
724a4e36a32SMatthew G. Knepley }
725a4e36a32SMatthew G. Knepley 
MonitorError(TS ts,PetscInt step,PetscReal crtime,Vec u,PetscCtx ctx)726*2a8381b2SBarry Smith static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, PetscCtx ctx)
727d71ae5a4SJacob Faibussowitsch {
728*2a8381b2SBarry Smith   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
729a4e36a32SMatthew G. Knepley   void     *ctxs[3];
730a4e36a32SMatthew G. Knepley   DM        dm;
731a4e36a32SMatthew G. Knepley   PetscDS   ds;
732a4e36a32SMatthew G. Knepley   Vec       v;
733a4e36a32SMatthew G. Knepley   PetscReal ferrors[3];
734b22b7e2eSMatthew G. Knepley   PetscInt  tl, l, f;
735a4e36a32SMatthew G. Knepley 
736a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
7379566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
7389566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
739a4e36a32SMatthew G. Knepley 
7409566063dSJacob Faibussowitsch   for (f = 0; f < 3; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
7419566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors));
7429566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tl));
7439566063dSJacob Faibussowitsch   for (l = 0; l < tl; ++l) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t"));
7449566063dSJacob Faibussowitsch   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]));
745a4e36a32SMatthew G. Knepley 
7469566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
7479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution"));
7489566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
7499566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
750a4e36a32SMatthew G. Knepley 
7519566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &v));
7529566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
7539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
7549566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
7559566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &v));
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757a4e36a32SMatthew G. Knepley }
758a4e36a32SMatthew G. Knepley 
759d2ae0b49SMatthew G. Knepley /* Note that adv->x0 will not be correct after migration */
ComputeParticleError(TS ts,Vec u,Vec e)760d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e)
761d71ae5a4SJacob Faibussowitsch {
762a4e36a32SMatthew G. Knepley   AdvCtx            *adv;
763a4e36a32SMatthew G. Knepley   DM                 sdm;
764a4e36a32SMatthew G. Knepley   Parameter         *param;
765a4e36a32SMatthew G. Knepley   const PetscScalar *xp0, *xp;
766a4e36a32SMatthew G. Knepley   PetscScalar       *ep;
767a4e36a32SMatthew G. Knepley   PetscReal          time;
768a4e36a32SMatthew G. Knepley   PetscInt           dim, Np, p;
769a4e36a32SMatthew G. Knepley   MPI_Comm           comm;
770a4e36a32SMatthew G. Knepley 
771a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
7729566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &time));
7739566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &adv));
774*2a8381b2SBarry Smith   PetscCall(PetscBagGetData(adv->ctx->bag, &param));
7759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
7769566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sdm));
7779566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sdm, &dim));
7789566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(sdm, &Np));
7799566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(adv->x0, &xp0));
7809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &xp));
7819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(e, &ep));
782a4e36a32SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
783a4e36a32SMatthew G. Knepley     PetscScalar x[3];
784a4e36a32SMatthew G. Knepley     PetscReal   x0[3];
785a4e36a32SMatthew G. Knepley     PetscInt    d;
786a4e36a32SMatthew G. Knepley 
787a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p * dim + d]);
7889566063dSJacob Faibussowitsch     PetscCall(adv->exact(dim, time, x0, 1, x, param));
789a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) ep[p * dim + d] += x[d] - xp[p * dim + d];
790a4e36a32SMatthew G. Knepley   }
7919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(adv->x0, &xp0));
7929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &xp));
7939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(e, &ep));
7943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
795a4e36a32SMatthew G. Knepley }
796a4e36a32SMatthew G. Knepley 
MonitorParticleError(TS ts,PetscInt step,PetscReal time,Vec u,PetscCtx ctx)797*2a8381b2SBarry Smith static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, PetscCtx ctx)
798d71ae5a4SJacob Faibussowitsch {
799a4e36a32SMatthew G. Knepley   AdvCtx            *adv = (AdvCtx *)ctx;
800a4e36a32SMatthew G. Knepley   DM                 sdm;
801a4e36a32SMatthew G. Knepley   Parameter         *param;
802a4e36a32SMatthew G. Knepley   const PetscScalar *xp0, *xp;
803a4e36a32SMatthew G. Knepley   PetscReal          error = 0.0;
804b22b7e2eSMatthew G. Knepley   PetscInt           dim, tl, l, Np, p;
805a4e36a32SMatthew G. Knepley   MPI_Comm           comm;
806a4e36a32SMatthew G. Knepley 
807a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
808*2a8381b2SBarry Smith   PetscCall(PetscBagGetData(adv->ctx->bag, &param));
8099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
8109566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &sdm));
8119566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sdm, &dim));
8129566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(sdm, &Np));
8139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(adv->x0, &xp0));
8149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &xp));
815a4e36a32SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
816a4e36a32SMatthew G. Knepley     PetscScalar x[3];
817a4e36a32SMatthew G. Knepley     PetscReal   x0[3];
818a4e36a32SMatthew G. Knepley     PetscReal   perror = 0.0;
819a4e36a32SMatthew G. Knepley     PetscInt    d;
820a4e36a32SMatthew G. Knepley 
821a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p * dim + d]);
8229566063dSJacob Faibussowitsch     PetscCall(adv->exact(dim, time, x0, 1, x, param));
823a4e36a32SMatthew G. Knepley     for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p * dim + d]));
824a4e36a32SMatthew G. Knepley     error += perror;
825a4e36a32SMatthew G. Knepley   }
8269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(adv->x0, &xp0));
8279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &xp));
8289566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tl));
8299566063dSJacob Faibussowitsch   for (l = 0; l < tl; ++l) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t"));
8309566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int)step, (double)time, (double)error));
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
832a4e36a32SMatthew G. Knepley }
833a4e36a32SMatthew G. Knepley 
AdvectParticles(TS ts)834d71ae5a4SJacob Faibussowitsch static PetscErrorCode AdvectParticles(TS ts)
835d71ae5a4SJacob Faibussowitsch {
836a4e36a32SMatthew G. Knepley   TS        sts;
837a4e36a32SMatthew G. Knepley   DM        sdm;
838ad95fc09SMatthew G. Knepley   Vec       coordinates;
839a4e36a32SMatthew G. Knepley   AdvCtx   *adv;
840a4e36a32SMatthew G. Knepley   PetscReal time;
841d2ae0b49SMatthew G. Knepley   PetscBool lreset, reset;
842d2ae0b49SMatthew G. Knepley   PetscInt  dim, n, N, newn, newN;
843a4e36a32SMatthew G. Knepley 
844a4e36a32SMatthew G. Knepley   PetscFunctionBeginUser;
8459566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ts, "_SwarmTS", (PetscObject *)&sts));
8469566063dSJacob Faibussowitsch   PetscCall(TSGetDM(sts, &sdm));
847*2a8381b2SBarry Smith   PetscCall(TSGetRHSFunction(sts, NULL, NULL, &adv));
8489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sdm, &dim));
8499566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(sdm, &N));
8509566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(sdm, &n));
8519566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates));
8529566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &time));
8539566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(sts, time));
854a4e36a32SMatthew G. Knepley   adv->tf = time;
8559566063dSJacob Faibussowitsch   PetscCall(TSSolve(sts, coordinates));
8569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates));
8579566063dSJacob Faibussowitsch   PetscCall(VecCopy(adv->uf, adv->ui));
858a4e36a32SMatthew G. Knepley   adv->ti = adv->tf;
859a4e36a32SMatthew G. Knepley 
8609566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate(sdm, PETSC_TRUE));
8619566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetSize(sdm, &newN));
8629566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(sdm, &newn));
863d2ae0b49SMatthew G. Knepley   lreset = (n != newn || N != newN) ? PETSC_TRUE : PETSC_FALSE;
8645440e5dcSBarry Smith   PetscCallMPI(MPIU_Allreduce(&lreset, &reset, 1, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sts)));
865ad95fc09SMatthew G. Knepley   if (reset) {
8669566063dSJacob Faibussowitsch     PetscCall(TSReset(sts));
8679566063dSJacob Faibussowitsch     PetscCall(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor));
868ad95fc09SMatthew G. Knepley   }
8699566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view"));
8703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
871a4e36a32SMatthew G. Knepley }
872a4e36a32SMatthew G. Knepley 
main(int argc,char ** argv)873d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
874d71ae5a4SJacob Faibussowitsch {
875a4e36a32SMatthew G. Knepley   DM        dm, sdm;
876a4e36a32SMatthew G. Knepley   TS        ts, sts;
877d2ae0b49SMatthew G. Knepley   Vec       u, xtmp;
878a4e36a32SMatthew G. Knepley   AppCtx    user;
879a4e36a32SMatthew G. Knepley   AdvCtx    adv;
880a4e36a32SMatthew G. Knepley   PetscReal t;
881a4e36a32SMatthew G. Knepley   PetscInt  dim;
882a4e36a32SMatthew G. Knepley 
883327415f7SBarry Smith   PetscFunctionBeginUser;
8849566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
8859566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
8869566063dSJacob Faibussowitsch   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
8879566063dSJacob Faibussowitsch   PetscCall(SetupParameters(&user));
8889566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
8899566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
8909566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
8919566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &user));
892a4e36a32SMatthew G. Knepley   /* Discretize chemical species */
8939566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &sdm));
8949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sdm, "part_"));
8959566063dSJacob Faibussowitsch   PetscCall(DMSetType(sdm, DMSWARM));
8969566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
8979566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(sdm, dim));
8989566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(sdm, dm));
899a4e36a32SMatthew G. Knepley   /* Setup problem */
9009566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, sdm, &user));
9019566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
902a4e36a32SMatthew G. Knepley 
9039566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
9049566063dSJacob Faibussowitsch   PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace));
905a4e36a32SMatthew G. Knepley 
9069566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user));
9079566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user));
9089566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user));
9099566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
9109566063dSJacob Faibussowitsch   PetscCall(TSSetPreStep(ts, RemoveDiscretePressureNullspace));
9119566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MonitorError, &user, NULL));
9129566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
913a4e36a32SMatthew G. Knepley 
9149566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions)); /* Must come after SetFromOptions() */
9159566063dSJacob Faibussowitsch   PetscCall(SetInitialConditions(ts, u));
9169566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
9179566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, 0, t));
9189566063dSJacob Faibussowitsch   PetscCall(DMTSCheckFromOptions(ts, u));
919a4e36a32SMatthew G. Knepley 
920a4e36a32SMatthew G. Knepley   /* Setup particle position integrator */
9219566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &sts));
9229566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)sts, "part_"));
9239566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)sts, (PetscObject)ts, 1));
9249566063dSJacob Faibussowitsch   PetscCall(TSSetDM(sts, sdm));
9259566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(sts, TS_NONLINEAR));
9269566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP));
9279566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(sts, MonitorParticleError, &adv, NULL));
9289566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(sts));
9299566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(sts, &adv));
9309566063dSJacob Faibussowitsch   PetscCall(TSSetComputeExactError(sts, ComputeParticleError));
9319566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(sts, SetInitialParticleConditions));
932a4e36a32SMatthew G. Knepley   adv.ti = t;
933a4e36a32SMatthew G. Knepley   adv.uf = u;
9349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(adv.uf, &adv.ui));
9359566063dSJacob Faibussowitsch   PetscCall(VecCopy(u, adv.ui));
9369566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(sts, NULL, FreeStreaming, &adv));
9379566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, AdvectParticles));
9389566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)ts, "_SwarmTS", (PetscObject)sts));
9399566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor));
9409566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(sdm, &adv.x0));
9419566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp));
9429566063dSJacob Faibussowitsch   PetscCall(VecCopy(xtmp, adv.x0));
9439566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp));
944a4e36a32SMatthew G. Knepley   switch (user.solType) {
945d71ae5a4SJacob Faibussowitsch   case SOL_TRIG_TRIG:
946d71ae5a4SJacob Faibussowitsch     adv.exact = trig_trig_x;
947d71ae5a4SJacob Faibussowitsch     break;
948d71ae5a4SJacob Faibussowitsch   default:
949d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType);
950a4e36a32SMatthew G. Knepley   }
951a4e36a32SMatthew G. Knepley   adv.ctx = &user;
952a4e36a32SMatthew G. Knepley 
9539566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
9549566063dSJacob Faibussowitsch   PetscCall(DMTSCheckFromOptions(ts, u));
9559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution"));
956a4e36a32SMatthew G. Knepley 
9579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
9589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&adv.x0));
9599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&adv.ui));
9609566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
9619566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sdm));
9629566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
9639566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&sts));
9649566063dSJacob Faibussowitsch   PetscCall(PetscBagDestroy(&user.bag));
9659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
966b122ec5aSJacob Faibussowitsch   return 0;
967a4e36a32SMatthew G. Knepley }
968a4e36a32SMatthew G. Knepley 
969a4e36a32SMatthew G. Knepley /*TEST
970a4e36a32SMatthew G. Knepley 
971a4e36a32SMatthew G. Knepley   # Swarm does not work with complex
972a4e36a32SMatthew G. Knepley   test:
973a4e36a32SMatthew G. Knepley     suffix: 2d_tri_p2_p1_p1_tconvp
974a4e36a32SMatthew G. Knepley     requires: triangle !single !complex
975a4e36a32SMatthew G. Knepley     args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \
976a4e36a32SMatthew G. Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
977188af4bfSBarry Smith       -dmts_check .001 -ts_max_steps 4 -ts_time_step 0.1 -ts_monitor_cancel \
978a4e36a32SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
979a4e36a32SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
980a4e36a32SMatthew G. Knepley         -fieldsplit_0_pc_type lu \
981a4e36a32SMatthew G. Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
982a4e36a32SMatthew G. Knepley       -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \
983188af4bfSBarry Smith       -part_ts_max_steps 2 -part_ts_time_step 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel
984b22b7e2eSMatthew G. Knepley   test:
985b22b7e2eSMatthew G. Knepley     suffix: 2d_tri_p2_p1_p1_exit
986b22b7e2eSMatthew G. Knepley     requires: triangle !single !complex
987b22b7e2eSMatthew G. Knepley     args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 1 \
988b22b7e2eSMatthew G. Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
989188af4bfSBarry Smith       -dmts_check .001 -ts_max_steps 10 -ts_time_step 0.1 \
990b22b7e2eSMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
991b22b7e2eSMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
992b22b7e2eSMatthew G. Knepley         -fieldsplit_0_pc_type lu \
993b22b7e2eSMatthew G. Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
994b22b7e2eSMatthew G. Knepley       -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \
995188af4bfSBarry Smith       -part_ts_max_steps 20 -part_ts_time_step 0.05
996a4e36a32SMatthew G. Knepley 
997a4e36a32SMatthew G. Knepley TEST*/
998