xref: /libCEED/examples/fluids/qfunctions/advection.h (revision 5ebd836c59d60a2e5e1cb67f6731404c7da26f85)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Advection initial condition and operator for Navier-Stokes example using PETSc
10 #include <ceed.h>
11 #include <math.h>
12 
13 #include "advection_types.h"
14 #include "newtonian_state.h"
15 #include "newtonian_types.h"
16 #include "stabilization_types.h"
17 #include "utils.h"
18 
19 // *****************************************************************************
20 // This QFunction sets the initial conditions and the boundary conditions
21 //   for two test cases: ROTATION and TRANSLATION
22 //
23 // -- ROTATION (default)
24 //      Initial Conditions:
25 //        Mass Density:
26 //          Constant mass density of 1.0
27 //        Momentum Density:
28 //          Rotational field in x,y
29 //        Energy Density:
30 //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
31 //            increases to (1.-r/rc), then 0. everywhere else
32 //
33 //      Boundary Conditions:
34 //        Mass Density:
35 //          0.0 flux
36 //        Momentum Density:
37 //          0.0
38 //        Energy Density:
39 //          0.0 flux
40 //
41 // -- TRANSLATION
42 //      Initial Conditions:
43 //        Mass Density:
44 //          Constant mass density of 1.0
45 //        Momentum Density:
46 //           Constant rectilinear field in x,y
47 //        Energy Density:
48 //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
49 //            increases to (1.-r/rc), then 0. everywhere else
50 //
51 //      Boundary Conditions:
52 //        Mass Density:
53 //          0.0 flux
54 //        Momentum Density:
55 //          0.0
56 //        Energy Density:
57 //          Inflow BCs:
58 //            E = E_wind
59 //          Outflow BCs:
60 //            E = E(boundary)
61 //          Both In/Outflow BCs for E are applied weakly in the
62 //            QFunction "Advection2d_Sur"
63 //
64 // *****************************************************************************
65 
66 // *****************************************************************************
67 // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
68 // *****************************************************************************
69 CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
70   const SetupContextAdv context = (SetupContextAdv)ctx;
71   const CeedScalar      rc      = context->rc;
72   const CeedScalar      lx      = context->lx;
73   const CeedScalar      ly      = context->ly;
74   const CeedScalar      lz      = dim == 2 ? 0. : context->lz;
75   const CeedScalar     *wind    = context->wind;
76 
77   const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
78   const CeedScalar theta     = dim == 2 ? M_PI / 3 : M_PI;
79   const CeedScalar x0[3]     = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
80 
81   const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
82 
83   CeedScalar r = 0.;
84   switch (context->initial_condition_type) {
85     case ADVECTIONIC_BUBBLE_SPHERE:
86     case ADVECTIONIC_BUBBLE_CYLINDER:
87       r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
88       break;
89     case ADVECTIONIC_COSINE_HILL:
90       r = sqrt(Square(x - center[0]) + Square(y - center[1]));
91       break;
92     case ADVECTIONIC_SKEW:
93       break;
94   }
95 
96   switch (context->wind_type) {
97     case WIND_ROTATION:
98       q[0] = 1.;
99       q[1] = -(y - center[1]);
100       q[2] = (x - center[0]);
101       q[3] = 0;
102       break;
103     case WIND_TRANSLATION:
104       q[0] = 1.;
105       q[1] = wind[0];
106       q[2] = wind[1];
107       q[3] = dim == 2 ? 0. : wind[2];
108       break;
109     default:
110       return 1;
111   }
112 
113   switch (context->initial_condition_type) {
114     case ADVECTIONIC_BUBBLE_SPHERE:
115     case ADVECTIONIC_BUBBLE_CYLINDER:
116       switch (context->bubble_continuity_type) {
117         // original continuous, smooth shape
118         case BUBBLE_CONTINUITY_SMOOTH:
119           q[4] = r <= rc ? (1. - r / rc) : 0.;
120           break;
121         // discontinuous, sharp back half shape
122         case BUBBLE_CONTINUITY_BACK_SHARP:
123           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
124           break;
125         // attempt to define a finite thickness that will get resolved under grid refinement
126         case BUBBLE_CONTINUITY_THICK:
127           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
128           break;
129         case BUBBLE_CONTINUITY_COSINE:
130           q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
131           break;
132       }
133       break;
134     case ADVECTIONIC_COSINE_HILL: {
135       CeedScalar half_width = context->lx / 2;
136       q[4]                  = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
137     } break;
138     case ADVECTIONIC_SKEW: {
139       CeedScalar       skewed_barrier[3]  = {wind[0], wind[1], 0};
140       CeedScalar       inflow_to_point[3] = {x - context->lx / 2, y, 0};
141       CeedScalar       cross_product[3]   = {0};
142       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
143       Cross3(skewed_barrier, inflow_to_point, cross_product);
144 
145       q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
146       if ((x < boundary_threshold && wind[0] < boundary_threshold) ||                // outflow at -x boundary
147           (y < boundary_threshold && wind[1] < boundary_threshold) ||                // outflow at -y boundary
148           (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) ||  // outflow at +x boundary
149           (y > context->ly - boundary_threshold && wind[1] > boundary_threshold)     // outflow at +y boundary
150       ) {
151         q[4] = 0;
152       }
153     } break;
154   }
155   return 0;
156 }
157 
158 // *****************************************************************************
159 // This QFunction sets the initial conditions for 3D advection
160 // *****************************************************************************
161 CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
162   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
163   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
164 
165   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
166     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
167     CeedScalar       q[5] = {0.};
168 
169     Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
170     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
171   }
172   return 0;
173 }
174 
175 // *****************************************************************************
176 // This QFunction sets the initial conditions for 2D advection
177 // *****************************************************************************
178 CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
179   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
180   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
181   const SetupContextAdv context    = (SetupContextAdv)ctx;
182 
183   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
184     const CeedScalar x[]  = {X[0][i], X[1][i]};
185     CeedScalar       q[5] = {0.};
186 
187     Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
188     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
189   }
190   return 0;
191 }
192 
193 CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) {
194   switch (N) {
195     case 2:
196       QdataUnpack_2D(Q, i, q_data, wdetJ, (CeedScalar(*)[2])dXdx);
197       break;
198     case 3:
199       QdataUnpack_3D(Q, i, q_data, wdetJ, (CeedScalar(*)[3])dXdx);
200       break;
201   }
202 }
203 
204 CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx,
205                                                  CeedScalar *normal) {
206   switch (N) {
207     case 2:
208       QdataBoundaryUnpack_2D(Q, i, q_data, wdetJ, normal);
209       break;
210     case 3:
211       QdataBoundaryUnpack_3D(Q, i, q_data, wdetJ, (CeedScalar(*)[3])dXdx, normal);
212       break;
213   }
214   return CEED_ERROR_SUCCESS;
215 }
216 
217 CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
218                                                                  StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
219                                                                  State *grad_s) {
220   switch (N) {
221     case 2: {
222       for (CeedInt k = 0; k < 2; k++) {
223         CeedScalar dqi[5];
224         for (CeedInt j = 0; j < 5; j++) {
225           dqi[j] = grad_q[(Q * 5) * 0 + Q * j + i] * dXdx[0 * N + k] + grad_q[(Q * 5) * 1 + Q * j + i] * dXdx[1 * N + k];
226         }
227         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
228       }
229       CeedScalar U[5] = {0.};
230       grad_s[2]       = StateFromU(gas, U);
231     } break;
232     case 3:
233       StatePhysicalGradientFromReference(Q, i, gas, s, state_var, grad_q, (CeedScalar(*)[3])dXdx, grad_s);
234       break;
235   }
236 }
237 
238 // *****************************************************************************
239 // This QFunction implements Advection for implicit time stepping method
240 // *****************************************************************************
241 CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
242   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
243   const CeedScalar(*grad_q)            = in[1];
244   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
245   const CeedScalar(*q_data)            = in[3];
246 
247   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
248   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
249   CeedScalar *jac_data               = out[2];
250 
251   AdvectionContext                 context   = (AdvectionContext)ctx;
252   const CeedScalar                 CtauS     = context->CtauS;
253   const CeedScalar                 zeros[14] = {0.};
254   NewtonianIdealGasContext         gas;
255   struct NewtonianIdealGasContext_ gas_struct = {0};
256   gas                                         = &gas_struct;
257 
258   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
259     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
260     const State      s     = StateFromU(gas, qi);
261 
262     CeedScalar wdetJ, dXdx[9];
263     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
264     State grad_s[3];
265     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
266 
267     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
268 
269     for (CeedInt f = 0; f < 4; f++) {
270       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
271       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
272     }
273 
274     CeedScalar div_u = 0;
275     for (CeedInt j = 0; j < dim; j++) {
276       for (CeedInt k = 0; k < dim; k++) {
277         div_u += grad_s[k].Y.velocity[j];
278       }
279     }
280     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
281     CeedScalar strong_res  = q_dot[4][i] + strong_conv;
282 
283     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
284 
285     CeedScalar uX[3] = {0.};
286     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
287 
288     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
289       v[4][i] += wdetJ * strong_conv;
290     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
291       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
292     }
293 
294     CeedScalar TauS = 0;
295     switch (context->stabilization_tau) {
296       case STAB_TAU_CTAU:
297         TauS = CtauS / sqrt(Dot3(uX, uX));
298         break;
299       case STAB_TAU_ADVDIFF_SHAKIB: {
300         CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
301         MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
302 
303         MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
304         TauS = 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a);
305       } break;
306     }
307 
308     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
309         case STAB_NONE:
310           break;
311         case STAB_SU:
312           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
313           break;
314         case STAB_SUPG:
315           grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j];
316           break;
317       }
318     StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
319   }
320 }
321 
322 CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
323   IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
324   return 0;
325 }
326 
327 CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
328   IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
329   return 0;
330 }
331 
332 // *****************************************************************************
333 // This QFunction implements Advection for explicit time stepping method
334 // *****************************************************************************
335 CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
336   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
337   const CeedScalar(*grad_q)        = in[1];
338   const CeedScalar(*q_data)        = in[2];
339 
340   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
341   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
342 
343   AdvectionContext                 context = (AdvectionContext)ctx;
344   const CeedScalar                 CtauS   = context->CtauS;
345   NewtonianIdealGasContext         gas;
346   struct NewtonianIdealGasContext_ gas_struct = {0};
347   gas                                         = &gas_struct;
348 
349   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
350     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
351     const State      s     = StateFromU(gas, qi);
352 
353     CeedScalar wdetJ, dXdx[9];
354     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
355     State grad_s[3];
356     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
357 
358     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
359 
360     for (CeedInt f = 0; f < 4; f++) {
361       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
362       v[f][i] = 0.;
363     }
364 
365     CeedScalar div_u = 0;
366     for (CeedInt j = 0; j < dim; j++) {
367       for (CeedInt k = 0; k < dim; k++) {
368         div_u += grad_s[k].Y.velocity[j];
369       }
370     }
371     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
372 
373     CeedScalar uX[3] = {0.};
374     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
375 
376     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
377       v[4][i] = -wdetJ * strong_conv;
378       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
379     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
380       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
381       v[4][i] = 0.;
382     }
383 
384     const CeedScalar TauS = CtauS / sqrt(Dot3(uX, uX));
385     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
386         case STAB_NONE:
387           break;
388         case STAB_SU:
389         case STAB_SUPG:
390           grad_v[j][4][i] -= wdetJ * TauS * strong_conv * uX[j];
391           break;
392       }
393   }
394 }
395 
396 CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
397   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
398   return 0;
399 }
400 
401 CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
402   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
403   return 0;
404 }
405 
406 // *****************************************************************************
407 // This QFunction implements consistent outflow and inflow BCs
408 //      for advection
409 //
410 //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
411 //    sign(dot(wind, normal)) > 0 : outflow BCs
412 //    sign(dot(wind, normal)) < 0 : inflow BCs
413 //
414 //  Outflow BCs:
415 //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
416 //
417 //  Inflow BCs:
418 //    A prescribed Total Energy (E_wind) is applied weakly.
419 // *****************************************************************************
420 CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
421   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
422   const CeedScalar(*q_data_sur)    = in[2];
423 
424   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
425   AdvectionContext context     = (AdvectionContext)ctx;
426   const CeedScalar E_wind      = context->E_wind;
427   const CeedScalar strong_form = context->strong_form;
428   const bool       is_implicit = context->implicit;
429 
430   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
431     const CeedScalar rho  = q[0][i];
432     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
433     const CeedScalar E    = q[4][i];
434 
435     CeedScalar wdetJb, norm[3];
436     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, norm);
437     wdetJb *= is_implicit ? -1. : 1.;
438 
439     const CeedScalar u_normal = DotN(norm, u, dim);
440 
441     // No Change in density or momentum
442     for (CeedInt j = 0; j < 4; j++) {
443       v[j][i] = 0;
444     }
445     // Implementing in/outflow BCs
446     if (u_normal > 0) {  // outflow
447       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
448     } else {  // inflow
449       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
450     }
451   }
452   return 0;
453 }
454 
455 CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
456   Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
457   return 0;
458 }
459 
460 CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
461   Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
462   return 0;
463 }
464