|
96 | 96 | # [Zhang, Shu (2011)](https://doi.org/10.1098/rspa.2011.0153). |
97 | 97 |
|
98 | 98 | # It works the following way. For every passed (scalar) variable and for every DG element we calculate |
99 | | -# the minimal value $value_{min}$. If this value falls below the given threshold $\varepsilon$, |
| 99 | +# the minimal value $value_\text{min}$. If this value falls below the given threshold $\varepsilon$, |
100 | 100 | # the approximation is slightly adapted such that the minimal value of the relevant variable lies |
101 | 101 | # now above the threshold. |
102 | 102 | # ```math |
103 | | -# \underline{u}^{new} = \theta * \underline{u} + (1-\theta) * u_{mean} |
| 103 | +# \underline{u}^\text{new} = \theta * \underline{u} + (1-\theta) * u_\text{mean} |
104 | 104 | # ``` |
105 | 105 | # where $\underline{u}$ are the collected pointwise evaluation coefficients in element $e$ and |
106 | | -# $u_{mean}$ the integral mean of the quantity in $e$. The new coefficients are a convex combination |
| 106 | +# $u_\text{mean}$ the integral mean of the quantity in $e$. The new coefficients are a convex combination |
107 | 107 | # of these two values with factor |
108 | 108 | # ```math |
109 | | -# \theta = \frac{value_{mean} - \varepsilon}{value_{mean} - value_{min}}, |
| 109 | +# \theta = \frac{value_\text{mean} - \varepsilon}{value_\text{mean} - value_\text{min}}, |
110 | 110 | # ``` |
111 | | -# where $value_{mean}$ is the relevant variable evaluated for the mean value $u_{mean}$. |
| 111 | +# where $value_\text{mean}$ is the relevant variable evaluated for the mean value $u_\text{mean}$. |
112 | 112 |
|
113 | 113 | # The adapted approximation keeps the exact same mean value, but the relevant variable is now greater |
114 | 114 | # or equal the threshold $\varepsilon$ at every node in every element. |
@@ -225,6 +225,137 @@ sol = solve(ode, CarpenterKennedy2N54(stage_limiter!, williamson_condition=false |
225 | 225 | using Plots |
226 | 226 | plot(sol) |
227 | 227 |
|
| 228 | +# # Entropy bounded limiter |
| 229 | + |
| 230 | +# As argued in the description of the positivity preserving limiter above it might sometimes be |
| 231 | +# necessary to apply advanced techniques to ensure a physically meaningful solution. |
| 232 | +# Apart from the positivity of pressure and density, the physical entropy of the system should increase |
| 233 | +# over the course of a simulation, see e.g. [this](https://doi.org/10.1016/0168-9274(86)90029-2) paper by Tadmor where this property is |
| 234 | +# shown for the compressible Euler equations. |
| 235 | +# As this is not necessarily the case for the numerical approximation (especially for the high-order, non-diffusive DG discretizations), |
| 236 | +# Lv and Ihme devised an a-posteriori limiter in [this paper](https://doi.org/10.1016/j.jcp.2015.04.026) which can be applied after each Runge-Kutta stage. |
| 237 | +# This limiter enforces a non-decrease in the physical, thermodynamic entropy $S$ |
| 238 | +# by bounding the entropy decrease (entropy increase is always tolerated) $\Delta S$ in each grid cell. |
| 239 | +# |
| 240 | +# This translates into a requirement that the entropy of the limited approximation $S\Big(\mathcal{L}\big[\boldsymbol u(\Delta t) \big] \Big)$ should be |
| 241 | +# greater or equal than the previous iterates' entropy $S\big(\boldsymbol u(0) \big)$, enforced at each quadrature point: |
| 242 | +# ```math |
| 243 | +# S\Big(\mathcal{L}\big[\boldsymbol u(\Delta t, \boldsymbol{x}_i) \big] \Big) \overset{!}{\geq} S\big(\boldsymbol u(0, \boldsymbol{x}_i) \big), \quad i = 1, \dots, (k+1)^d |
| 244 | +# ``` |
| 245 | +# where $k$ denotes the polynomial degree of the element-wise solution and $d$ is the spatial dimension. |
| 246 | +# For an ideal gas, the thermodynamic entropy $S(\boldsymbol u) = S(p, \rho)$ is given by |
| 247 | +# ```math |
| 248 | +# S(p, \rho) = \ln \left( \frac{p}{\rho^\gamma} \right) \: . |
| 249 | +# ``` |
| 250 | +# Thus, the non-decrease in entropy can be reformulated as |
| 251 | +# ```math |
| 252 | +# p(\boldsymbol{x}_i) - e^{ S\big(\boldsymbol u(0, \boldsymbol{x}_i) \big)} \cdot \rho(\boldsymbol{x}_i)^\gamma \overset{!}{\geq} 0, \quad i = 1, \dots, (k+1)^d \: . |
| 253 | +# ``` |
| 254 | +# In a practical simulation, we might tolerate a maximum (exponentiated) entropy decrease per element, i.e., |
| 255 | +# ```math |
| 256 | +# \Delta e^S \coloneqq \min_{i} \left\{ p(\boldsymbol{x}_i) - e^{ S\big(\boldsymbol u(0, \boldsymbol{x}_i) \big)} \cdot \rho(\boldsymbol{x}_i)^\gamma \right\} < c |
| 257 | +# ``` |
| 258 | +# with hyper-parameter $c$ which is to be specified by the user. |
| 259 | +# The default value for the corresponding parameter $c=$ `exp_entropy_decrease_max` is set to $-10^{-13}$, i.e., slightly less than zero to |
| 260 | +# avoid spurious limiter actions for cells in which the entropy remains effectively constant. |
| 261 | +# Other values can be specified by setting the `exp_entropy_decrease_max` keyword in the constructor of the limiter: |
| 262 | +# ```julia |
| 263 | +# stage_limiter! = EntropyBoundedLimiter(exp_entropy_decrease_max=-1e-9) |
| 264 | +# ``` |
| 265 | +# Smaller values (larger in absolute value) for `exp_entropy_decrease_max` relax the entropy increase requirement and are thus less diffusive. |
| 266 | +# On the other hand, for larger values (smaller in absolute value) of `exp_entropy_decrease_max` the limiter acts more often and the solution becomes more diffusive. |
| 267 | +# |
| 268 | +# In particular, we compute again a limiting parameter $\vartheta \in [0, 1]$ which is then used to blend the |
| 269 | +# unlimited nodal values $\boldsymbol u$ with the mean value $\boldsymbol u_{\text{mean}}$ of the element: |
| 270 | +# ```math |
| 271 | +# \mathcal{L} [\boldsymbol u](\vartheta) \coloneqq (1 - \vartheta) \boldsymbol u + \vartheta \cdot \boldsymbol u_{\text{mean}} |
| 272 | +# ``` |
| 273 | +# For the exact definition of $\vartheta$ the interested user is referred to section 4.4 of the paper by Lv and Ihme. |
| 274 | +# Note that therein the limiting parameter is denoted by $\epsilon$, which is not to be confused with the threshold $\varepsilon$ of the Zhang-Shu limiter. |
| 275 | + |
| 276 | +# As for the positivity preserving limiter, the entropy bounded limiter may be applied after every Runge-Kutta stage. |
| 277 | +# Both fixed timestep methods such as [`CarpenterKennedy2N54`](https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/) and |
| 278 | +# adaptive timestep methods such as [`SSPRK43`](https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/) are supported. |
| 279 | +# We would like to remark that of course every `stage_limiter!` can also be used as a `step_limiter!`, i.e., |
| 280 | +# acting only after the full time step has been taken. |
| 281 | + |
| 282 | +# As an example, we consider a variant of the [1D medium blast wave example](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/tree_1d_dgsem/elixir_euler_blast_wave.jl) |
| 283 | +# wherein the shock capturing method discussed above is employed to handle the shock. |
| 284 | +# Here, we use a standard DG solver with HLLC surface flux: |
| 285 | +using Trixi |
| 286 | + |
| 287 | +solver = DGSEM(polydeg = 3, surface_flux = flux_hllc) |
| 288 | + |
| 289 | +# The remaining setup is the same as in the standard example: |
| 290 | + |
| 291 | +using OrdinaryDiffEq |
| 292 | + |
| 293 | +############################################################################### |
| 294 | +# semidiscretization of the compressible Euler equations |
| 295 | + |
| 296 | +equations = CompressibleEulerEquations1D(1.4) |
| 297 | + |
| 298 | +function initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations1D) |
| 299 | + ## Modified From Hennemann & Gassner JCP paper 2020 (Sec. 6.3) -> "medium blast wave" |
| 300 | + ## Set up polar coordinates |
| 301 | + inicenter = SVector(0.0) |
| 302 | + x_norm = x[1] - inicenter[1] |
| 303 | + r = abs(x_norm) |
| 304 | + ## The following code is equivalent to |
| 305 | + ## phi = atan(0.0, x_norm) |
| 306 | + ## cos_phi = cos(phi) |
| 307 | + ## in 1D but faster |
| 308 | + cos_phi = x_norm > 0 ? one(x_norm) : -one(x_norm) |
| 309 | + |
| 310 | + ## Calculate primitive variables |
| 311 | + rho = r > 0.5 ? 1.0 : 1.1691 |
| 312 | + v1 = r > 0.5 ? 0.0 : 0.1882 * cos_phi |
| 313 | + p = r > 0.5 ? 1.0E-3 : 1.245 |
| 314 | + |
| 315 | + return prim2cons(SVector(rho, v1, p), equations) |
| 316 | +end |
| 317 | +initial_condition = initial_condition_blast_wave |
| 318 | + |
| 319 | +coordinates_min = (-2.0,) |
| 320 | +coordinates_max = (2.0,) |
| 321 | +mesh = TreeMesh(coordinates_min, coordinates_max, |
| 322 | + initial_refinement_level = 6, |
| 323 | + n_cells_max = 10_000) |
| 324 | + |
| 325 | +semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver) |
| 326 | + |
| 327 | +############################################################################### |
| 328 | +# ODE solvers, callbacks etc. |
| 329 | + |
| 330 | +tspan = (0.0, 12.5) |
| 331 | +ode = semidiscretize(semi, tspan) |
| 332 | + |
| 333 | +summary_callback = SummaryCallback() |
| 334 | + |
| 335 | +analysis_interval = 100 |
| 336 | + |
| 337 | +analysis_callback = AnalysisCallback(semi, interval = analysis_interval) |
| 338 | + |
| 339 | +alive_callback = AliveCallback(analysis_interval = analysis_interval) |
| 340 | + |
| 341 | +stepsize_callback = StepsizeCallback(cfl = 0.5) |
| 342 | + |
| 343 | +callbacks = CallbackSet(summary_callback, |
| 344 | + analysis_callback, alive_callback, |
| 345 | + stepsize_callback) |
| 346 | + |
| 347 | +# We specify the `stage_limiter!` supplied to the classic SSPRK33 integrator |
| 348 | +# with strict entropy increase enforcement |
| 349 | +# (recall that the tolerated exponentiated entropy decrease is set to -1e-13). |
| 350 | +stage_limiter! = EntropyBoundedLimiter() |
| 351 | + |
| 352 | +# We run the simulation with the SSPRK33 method and the entropy bounded limiter: |
| 353 | +sol = solve(ode, SSPRK33(stage_limiter!); |
| 354 | + dt = 1.0, |
| 355 | + callback = callbacks); |
| 356 | + |
| 357 | +using Plots |
| 358 | +plot(sol) |
228 | 359 |
|
229 | 360 | # ## Package versions |
230 | 361 |
|
|
0 commit comments