Skip to content

Commit 5831468

Browse files
committed
test: add 2D advection WENO7/TENO7 convergence cases
The existing 2D vortex test can't measure 7th-order rate because the primitive->conserved IC has an O(eps^3 h^2) covariance floor that dominates WENO7's O(eps^2 h^7) scheme error at testable resolutions (h > eps^(1/5) needed; eps=0.01 gives N < 25). Add a 2D smooth-advection case with rho = 1 + 0.2 sin(2pi(x+y)) and v=(1,1) on the unit square with periodic BCs; linear primitives mean no covariance floor. Cases registered: Convergence -> 2D Adv -> {WENO7, TENO7} at N=[80, 96], CFL=0.005, 4 ranks, T=0.5 (one diagonal period). Both measure rate 6.86 (threshold 6.5) in ~5-8 min wall time per case.
1 parent 4334fa5 commit 5831468

2 files changed

Lines changed: 153 additions & 0 deletions

File tree

Lines changed: 116 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,116 @@
1+
#!/usr/bin/env python3
2+
"""
3+
2D single-fluid Euler convergence case (smooth diagonal advection).
4+
5+
Density wave rho = 1 + 0.2*sin(2*pi*(x+y)) with uniform velocity (u, v) = (1, 1)
6+
on the unit square with periodic BCs. Pressure p = 1 throughout. For this IC
7+
the Euler equations reduce to pure advection of all variables along the
8+
diagonal at speed sqrt(2); the wave phase x+y - 2t returns to x+y after
9+
T = 0.5. L2(rho(T) - rho(0)) measures accumulated scheme spatial truncation
10+
error.
11+
12+
The point of this case is to exercise WENO7 / TENO7 in 2D without the
13+
primitive-to-conserved covariance floor that dominates the 2D isentropic
14+
vortex (run_convergence.py) at testable resolutions. Here all primitives are
15+
linear in the conserved variables, so there is no floor.
16+
"""
17+
18+
import argparse
19+
import json
20+
import math
21+
22+
parser = argparse.ArgumentParser(description="2D Euler diagonal-advection convergence case")
23+
parser.add_argument("--mfc", type=json.loads, default="{}", metavar="DICT")
24+
parser.add_argument("-N", type=int, default=64, help="Grid points per dim (default: 64)")
25+
parser.add_argument("--order", type=int, default=5, help="WENO order: 1, 3, 5, or 7")
26+
parser.add_argument("--muscl", action="store_true", help="Use MUSCL-2 instead of WENO")
27+
parser.add_argument("--teno", action="store_true", help="Use TENO instead of WENO")
28+
parser.add_argument("--teno-ct", type=float, default=1e-6, help="TENO CT threshold (default: 1e-6)")
29+
parser.add_argument("--cfl", type=float, default=0.4, help="CFL number (default: 0.4; use 0.005 for WENO7/TENO7)")
30+
parser.add_argument("--no-mapped", action="store_true", help="Disable mapped WENO")
31+
parser.add_argument("--muscl-lim", type=int, default=0, help="MUSCL limiter: 0=unlimited 1=minmod ... (default: 0)")
32+
parser.add_argument("--time-stepper", type=int, default=3, help="Time stepper: 1=Euler 2=RK2 3=RK3 (default: 3)")
33+
args = parser.parse_args()
34+
35+
gamma = 1.4
36+
N = args.N
37+
m = N - 1
38+
L = 1.0
39+
dx = L / N
40+
41+
# Max wave speed: |u| + c with c = sqrt(gamma * p / rho_min). rho_min = 1 - 0.2.
42+
c_max = math.sqrt(gamma / 0.8) + 1.0
43+
dt = args.cfl * dx / c_max
44+
T_end = 0.5 # one full diagonal period: phase shift = 2*T_end = 1.0
45+
Nt = max(4, math.ceil(T_end / dt))
46+
dt = T_end / Nt
47+
48+
if args.muscl:
49+
scheme_params = {
50+
"recon_type": 2,
51+
"muscl_order": 2,
52+
"muscl_lim": args.muscl_lim,
53+
}
54+
else:
55+
scheme_params = {
56+
"recon_type": 1,
57+
"weno_order": args.order,
58+
"weno_eps": 1.0e-40,
59+
"weno_Re_flux": "F",
60+
"weno_avg": "F",
61+
"mapped_weno": "F" if (args.order == 1 or args.no_mapped or args.teno) else "T",
62+
"null_weights": "F",
63+
"mp_weno": "F",
64+
"teno": "T" if args.teno else "F",
65+
**({"teno_CT": args.teno_ct} if args.teno else {}),
66+
}
67+
68+
print(
69+
json.dumps(
70+
{
71+
"run_time_info": "F",
72+
"x_domain%beg": 0.0,
73+
"x_domain%end": L,
74+
"y_domain%beg": 0.0,
75+
"y_domain%end": L,
76+
"m": m,
77+
"n": m,
78+
"p": 0,
79+
"dt": dt,
80+
"t_step_start": 0,
81+
"t_step_stop": Nt,
82+
"t_step_save": Nt,
83+
"num_patches": 1,
84+
"model_eqns": 2,
85+
"alt_soundspeed": "F",
86+
"num_fluids": 1,
87+
"mpp_lim": "F",
88+
"mixture_err": "F",
89+
"time_stepper": args.time_stepper,
90+
"riemann_solver": 2,
91+
"wave_speeds": 1,
92+
"avg_state": 2,
93+
"bc_x%beg": -1,
94+
"bc_x%end": -1,
95+
"bc_y%beg": -1,
96+
"bc_y%end": -1,
97+
"format": 1,
98+
"precision": 2,
99+
"prim_vars_wrt": "T",
100+
"parallel_io": "F",
101+
"patch_icpp(1)%geometry": 3,
102+
"patch_icpp(1)%x_centroid": 0.5,
103+
"patch_icpp(1)%y_centroid": 0.5,
104+
"patch_icpp(1)%length_x": L,
105+
"patch_icpp(1)%length_y": L,
106+
"patch_icpp(1)%vel(1)": 1.0,
107+
"patch_icpp(1)%vel(2)": 1.0,
108+
"patch_icpp(1)%pres": 1.0,
109+
"patch_icpp(1)%alpha_rho(1)": "1.0 + 0.2 * sin(2.0 * pi * (x / lx + y / ly))",
110+
"patch_icpp(1)%alpha(1)": 1.0,
111+
"fluid_pp(1)%gamma": 1.0 / (gamma - 1.0),
112+
"fluid_pp(1)%pi_inf": 0.0,
113+
**scheme_params,
114+
}
115+
)
116+
)

toolchain/mfc/test/cases.py

Lines changed: 37 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -65,6 +65,18 @@
6565
("RK3", ["--order", "5", "--time-stepper", "3"], 3, 0.3, [0.50, 0.25]),
6666
]
6767

68+
# 2D diagonal advection (smooth density wave; no covariance floor): exercises
69+
# WENO7 / TENO7 in 2D. CFL=0.005 keeps RK3 temporal error below O(h^7).
70+
# With 4 ranks (2x2 decomp) the per-rank stencil constraint n+1 >= 5*weno_order=35
71+
# sets min_N=70 (rounded up to 80). N capped at 96 to keep CI cost reasonable:
72+
# 7th-order error scaling between N=80 and N=96 is ~3.6x, well above the
73+
# machine-precision floor and clean for rate measurement.
74+
_RES_2D_ADV_DEFAULT = [80, 96]
75+
_CONVERGENCE_2D_ADV_SCHEMES = [
76+
("WENO7", ["--order", "7", "--cfl", "0.005"], 7, 0.5, 80, 96),
77+
("TENO7", ["--order", "7", "--teno", "--teno-ct", "1e-9", "--cfl", "0.005"], 7, 0.5, 80, 96),
78+
]
79+
6880

6981
def add_convergence_cases(cases):
7082
num_ranks = 4
@@ -153,6 +165,30 @@ def add_convergence_cases(cases):
153165
)
154166
)
155167

168+
# 2D diagonal advection — same 4-rank decomposition as the rest of the suite.
169+
for label, extra_args, expected, tol, min_N, max_N in _CONVERGENCE_2D_ADV_SCHEMES:
170+
cases.append(
171+
define_convergence_case(
172+
f"Convergence -> 2D Adv -> {label}",
173+
spec={
174+
"runner": "2d_vortex",
175+
"case_path": "examples/2D_advection_convergence/case.py",
176+
"extra_args": extra_args,
177+
"expected_order": expected,
178+
"tol": tol,
179+
"resolutions": _RES_2D_ADV_DEFAULT,
180+
"min_N": min_N,
181+
"max_N": max_N,
182+
"ndim": 2,
183+
"domain_len": 1.0,
184+
"cons_vars": _CONS_VARS_2D,
185+
"primary_idx": 1,
186+
"num_ranks": num_ranks,
187+
},
188+
ppn=num_ranks,
189+
)
190+
)
191+
156192

157193
def get_bc_mods(bc: int, dimInfo):
158194
params = {}
@@ -1731,6 +1767,7 @@ def foreach_example():
17311767
"1D_euler_convergence",
17321768
"1D_advection_convergence",
17331769
"1D_sod_convergence",
1770+
"2D_advection_convergence",
17341771
"2D_zero_circ_vortex_analytical",
17351772
"3D_TaylorGreenVortex_analytical",
17361773
"3D_IGR_TaylorGreenVortex_nvidia",

0 commit comments

Comments
 (0)