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