Peano
Loading...
Searching...
No Matches
coupling.py
Go to the documentation of this file.
1import peano4, exahype2
2import os, sys
3import argparse
4import subprocess
5import mpmath as mp
6
7sys.path.insert(0, os.path.abspath("../aderdg/"))
8sys.path.insert(0, os.path.abspath("../aderdg/scenarios"))
9import scenarios
10
11modes = {
12 "release": peano4.output.CompileMode.Release,
13 "trace": peano4.output.CompileMode.Trace,
14 "assert": peano4.output.CompileMode.Asserts,
15 "debug": peano4.output.CompileMode.Debug,
16}
17
18available_prec = {
19 "bf16", "fp16", "fp32", "fp64"
20}
21
22available_scenarios = {
23 "AcousticPlanarWaves": scenarios.AcousticPlanarWaves(dimensions=2),
24 "AdvectionLinear": scenarios.AdvectionLinear(),
25 "ElasticPlanarWaves": scenarios.ElasticPlanarWaves(dimensions=2),
26 "EulerGaussianBell": scenarios.EulerGaussianBell(),
27 "EulerIsotropicVortex": scenarios.EulerIsotropicVortex(),
28 "SWERadialDamBreak": scenarios.SWERadialDamBreak(),
29 "SWERestingLake": scenarios.SWERestingLake(),
30}
31
32parser = argparse.ArgumentParser(description="ExaHyPE 2 - ADER testing script")
33
34parser.add_argument(
35 "-md",
36 "--mesh-depth",
37 dest="md",
38 type=int,
39 default=3,
40 help="Depth of coarsest mesh level, i.e if 2 is specified there will be 9 cells per dimension",
41)
42parser.add_argument(
43 "-amr",
44 "--adaptive-levels",
45 dest="adaptivity_levels",
46 type=int,
47 default=0,
48 help="Number of AMR grid levels on top of hmax (0 by default)",
49)
50parser.add_argument("-o", "--order", dest="order", type=int, default=3, help="DG Order")
51parser.add_argument(
52 "-p",
53 "--p",
54 dest="polynomials",
55 type=int,
56 default=0,
57 help="Polynomial type, 0 is Gauss-Legendre, 1 is Gauss-Lobatto",
58)
59parser.add_argument(
60 "-m", "--mode", dest="mode", default="release", help="|".join(modes.keys())
61)
62parser.add_argument(
63 "-pr", "--precision", dest="precision", default="fp64", help="|".join(available_prec)
64)
65parser.add_argument(
66 "-s",
67 "--scenario",
68 dest="s",
69 default="SWERadialDamBreak",
70 help="|".join(available_scenarios.keys()),
71)
72
73args = parser.parse_args()
74
75if args.s is None:
76 while True:
77 try:
78 s = input(
79 "Which of the following scenarios would you like to try out?\n"
80 + " - ".join(available_scenarios.keys())
81 + "\n"
82 )
83 scenario = available_scenarios[s]
84 except KeyError:
85 continue
86 else:
87 # User has specified a valid scenario
88 break
89else:
90 scenario = available_scenarios[args.s]
91
92order = args.order
93max_h = 1.1 * scenario._domain_size / (3.0**args.md)
94min_h = max_h * 3.0 ** (-args.adaptivity_levels)
95
96polynomials = (
97 exahype2.solvers.aderdg.Polynomials.Gauss_Legendre
98 if args.polynomials == 0
99 else exahype2.solvers.aderdg.Polynomials.Gauss_Lobatto
100)
101
102project = exahype2.Project(
103 ["tests", "exahype2", "aderdg"],
104 ".",
105 executable=scenario.__class__.__name__,
106)
107
108solver = exahype2.solvers.aderdg.GlobalAdaptiveTimeStep(
109 name=scenario.__class__.__name__,
110 order=order,
111 min_cell_h=min_h,
112 max_cell_h=max_h,
113 time_step_relaxation=0.9,
114 unknowns=scenario._equation.num_unknowns,
115 auxiliary_variables=scenario._equation.num_auxiliary_variables
116)
117
118solver.add_kernel_optimisations(
119 polynomials=polynomials, is_linear=scenario._equation.is_linear,
120 riemann_solver_implementation = scenario._equation.riemann_solver(),
121 precision = args.precision
122)
123
124solver.set_implementation(
125 initial_conditions=scenario.initial_conditions(),
126 boundary_conditions=scenario.boundary_conditions(),
127 eigenvalues=scenario._equation.eigenvalues(),
128 flux=scenario._equation.flux(),
129 ncp=scenario._equation.ncp(),
130)
131
132if scenario.analytical_solution() != exahype2.solvers.PDETerms.None_Implementation:
133 exahype2.solvers.aderdg.ErrorMeasurement(
134 solver,
135 error_measurement_implementation=scenario.analytical_solution(),
136 output_file_name="Error_" + scenario.__class__.__name__,
137 )
138
139project.add_solver(solver)
140
141
142
143solver2 = exahype2.solvers.aderdg.GlobalAdaptiveTimeStep(
144 name=scenario.__class__.__name__+"2",
145 order=order,
146 min_cell_h=min_h,
147 max_cell_h=max_h,
148 time_step_relaxation=0.9,
149 unknowns=scenario._equation.num_unknowns,
150 auxiliary_variables=scenario._equation.num_auxiliary_variables
151)
152
153solver2.add_kernel_optimisations(
154 polynomials=polynomials, is_linear=scenario._equation.is_linear,
155 riemann_solver_implementation = scenario._equation.riemann_solver(),
156 precision = "fp32"#precision = args.precision
157)
158
159solver2.set_implementation(
160 initial_conditions=scenario.initial_conditions(),
161 boundary_conditions=scenario.boundary_conditions(),
162 eigenvalues=scenario._equation.eigenvalues(),
163 flux=scenario._equation.flux(),
164 ncp=scenario._equation.ncp(),
165)
166
167
168
169if scenario.analytical_solution() != exahype2.solvers.PDETerms.None_Implementation:
170 exahype2.solvers.aderdg.ErrorMeasurement(
171 solver2,
172 error_measurement_implementation=scenario.analytical_solution(),
173 output_file_name="Error_" + scenario.__class__.__name__,
174 )
175
176project.add_solver(solver2)
177
178# from fuseADERSolvers import fuseADERSolvers
179# fuseADERSolvers(solver, solver2, "marker.getX()[0]<=0.0")
180
181sys.path.insert(0, os.path.abspath("./coupling"))
182from coupling import StaticCoupling
183
185 name = "limitingSolver",
186 regularSolver = solver,
187 limitingSolver = solver2,
188 limiting_criterion_implementation = """
189 return x[0]<=0.0;
190"""
191)
192project.add_solver(limiter)
193
194
195
196
197project.set_output_path("solutions")
198scenario.set_global_simulation_parameters(project)
199
200project.set_load_balancing(
201 "toolbox::loadbalancing::strategies::SpreadOutOnceGridStagnates",
202 "new ::exahype2::LoadBalancingConfiguration(0.98)",
203)
204project.set_build_mode(modes[args.mode])
205project = project.generate_Peano4_project("False")
206project.build(make_clean_first=True)
Scenario reproduced from Dumbser & Käser, https://doi.org/10.1111/j.1365-246X.2006....
Very simple scenario in which the initial value of x is shifted in each spatial dimension.
Scenario reproduced from Dumbser & Käser, https://doi.org/10.1111/j.1365-246X.2006....
Scenario reproduced from Ioratti, Dumbser & Loubère, https://doi.org/10.1007/s10915-020-01209-w (p.
Classic radial dam break SWE equations, with constant initial water height but a bump in the bathymet...
Resting lake scenario for the shallow water equations.