Branch data Line data Source code
1 : : // *****************************************************************************
2 : : /*!
3 : : \file src/Control/Inciter/InputDeck/InputDeck.hpp
4 : : \copyright 2012-2015 J. Bakosi,
5 : : 2016-2018 Los Alamos National Security, LLC.,
6 : : 2019-2021 Triad National Security, LLC.
7 : : All rights reserved. See the LICENSE file for details.
8 : : \brief Inciter's new input deck definition
9 : : \details This file defines the heterogeneous struct that is used for storing
10 : : the data from user input during the control file parsing of the
11 : : computational shock hydrodynamics tool, Inciter.
12 : : */
13 : : // *****************************************************************************
14 : : #ifndef InputDeck_h
15 : : #define InputDeck_h
16 : :
17 : : #include <getopt.h>
18 : : #include "Types.hpp"
19 : : #include "TaggedTuple.hpp"
20 : : #include "Inciter/CmdLine/CmdLine.hpp"
21 : : #include "Transfer.hpp"
22 : : #include "Inciter/OutVar.hpp"
23 : : #include "Inciter/Options/PDE.hpp"
24 : : #include "Inciter/Options/Problem.hpp"
25 : : #include "Inciter/Options/Scheme.hpp"
26 : : #include "Inciter/Options/Limiter.hpp"
27 : : #include "Inciter/Options/Flux.hpp"
28 : : #include "Inciter/Options/Initiate.hpp"
29 : : #include "Inciter/Options/AMRInitial.hpp"
30 : : #include "Inciter/Options/AMRError.hpp"
31 : : #include "Inciter/Options/PrefIndicator.hpp"
32 : : #include "Inciter/Options/MeshVelocity.hpp"
33 : : #include "Inciter/Options/MeshVelocitySmoother.hpp"
34 : : #include "Inciter/Options/Material.hpp"
35 : : #include "Options/PartitioningAlgorithm.hpp"
36 : : #include "Options/TxtFloatFormat.hpp"
37 : : #include "Options/FieldFile.hpp"
38 : : #include "Options/Error.hpp"
39 : : #include "Options/UserTable.hpp"
40 : :
41 : : namespace inciter {
42 : :
43 : : namespace ctr {
44 : :
45 : : using ncomp_t = std::size_t;
46 : :
47 : : using bclist = tk::TaggedTuple< brigand::list<
48 : : tag::dirichlet, std::vector< std::size_t >,
49 : : tag::symmetry, std::vector< std::size_t >,
50 : : tag::outlet, std::vector< std::size_t >,
51 : : tag::farfield, std::vector< std::size_t >,
52 : : tag::extrapolate, std::vector< std::size_t >,
53 : : tag::noslipwall, std::vector< std::size_t >,
54 : : tag::slipwall, std::vector< std::size_t >
55 : : > >;
56 : :
57 : : // Transport
58 : : using transportList = tk::TaggedTuple< brigand::list<
59 : : tag::physics, PhysicsType,
60 : : tag::ncomp, std::size_t,
61 : : tag::intsharp, int,
62 : : tag::intsharp_param, tk::real,
63 : : tag::problem, ProblemType,
64 : : tag::diffusivity, std::vector< tk::real >,
65 : : tag::lambda, std::vector< tk::real >,
66 : : tag::u0, std::vector< tk::real >
67 : : > >;
68 : :
69 : : // CompFlow
70 : : using compflowList = tk::TaggedTuple< brigand::list<
71 : : tag::physics, PhysicsType,
72 : : tag::problem, ProblemType,
73 : : tag::alpha, tk::real,
74 : : tag::beta, tk::real,
75 : : tag::betax, tk::real,
76 : : tag::betay, tk::real,
77 : : tag::betaz, tk::real,
78 : : tag::r0, tk::real,
79 : : tag::p0, tk::real,
80 : : tag::ce, tk::real,
81 : : tag::kappa, tk::real
82 : : > >;
83 : :
84 : : // MultiMat
85 : : using multimatList = tk::TaggedTuple< brigand::list<
86 : : tag::physics, PhysicsType,
87 : : tag::nmat, std::size_t,
88 : : tag::prelax, uint64_t,
89 : : tag::prelax_timescale, tk::real,
90 : : tag::intsharp, int,
91 : : tag::intsharp_param, tk::real,
92 : : tag::rho0constraint, uint64_t,
93 : : tag::dt_sos_massavg, int,
94 : : tag::problem, ProblemType,
95 : : tag::viscous, bool
96 : : > >;
97 : :
98 : : // MultiSpecies
99 : : using multispeciesList = tk::TaggedTuple< brigand::list<
100 : : tag::physics, PhysicsType,
101 : : tag::nspec, std::size_t,
102 : : tag::problem, ProblemType,
103 : : tag::viscous, bool
104 : : > >;
105 : :
106 : : // Material/EOS object
107 : : using materialList = tk::TaggedTuple< brigand::list<
108 : : tag::eos, MaterialType,
109 : : tag::id, std::vector< uint64_t >,
110 : : tag::gamma, std::vector< tk::real >,
111 : : tag::pstiff, std::vector< tk::real >,
112 : : tag::w_gru, std::vector< tk::real >,
113 : : tag::A_jwl, std::vector< tk::real >,
114 : : tag::B_jwl, std::vector< tk::real >,
115 : : tag::C_jwl, std::vector< tk::real >,
116 : : tag::R1_jwl, std::vector< tk::real >,
117 : : tag::R2_jwl, std::vector< tk::real >,
118 : : tag::rho0_jwl, std::vector< tk::real >,
119 : : tag::de_jwl, std::vector< tk::real >,
120 : : tag::rhor_jwl, std::vector< tk::real >,
121 : : tag::Tr_jwl, std::vector< tk::real >,
122 : : tag::Pr_jwl, std::vector< tk::real >,
123 : : tag::mu, std::vector< tk::real >,
124 : : tag::yield_stress, std::vector< tk::real >,
125 : : tag::alpha, std::vector< tk::real >,
126 : : tag::K0, std::vector< tk::real >,
127 : : tag::cv, std::vector< tk::real >,
128 : : tag::k, std::vector< tk::real >
129 : : > >;
130 : :
131 : : // Species/EOS object
132 : : using speciesList = tk::TaggedTuple< brigand::list<
133 : : tag::id, std::vector< uint64_t >,
134 : : tag::gamma, std::vector< tk::real >,
135 : : tag::R, std::vector< tk::real >,
136 : : tag::cp_coeff, std::vector< std::vector< std::vector< tk::real > > >,
137 : : tag::t_range, std::vector< std::vector< tk::real > >,
138 : : tag::dH_ref, std::vector< tk::real >
139 : : > >;
140 : :
141 : : // Boundary conditions block
142 : : using bcList = tk::TaggedTuple< brigand::list<
143 : : tag::mesh, std::vector< std::size_t >,
144 : : tag::dirichlet, std::vector< std::size_t >,
145 : : tag::symmetry, std::vector< std::size_t >,
146 : : tag::outlet, std::vector< std::size_t >,
147 : : tag::farfield, std::vector< std::size_t >,
148 : : tag::extrapolate, std::vector< std::size_t >,
149 : : tag::noslipwall, std::vector< std::size_t >,
150 : : tag::slipwall, std::vector< std::size_t >,
151 : : tag::velocity, std::vector< tk::real >,
152 : : tag::pressure, tk::real,
153 : : tag::density, tk::real,
154 : : tag::temperature, tk::real,
155 : : tag::mass_fractions, std::vector< tk::real >,
156 : : tag::materialid, std::size_t,
157 : : tag::inlet, std::vector<
158 : : tk::TaggedTuple< brigand::list<
159 : : tag::sideset, std::vector< uint64_t >,
160 : : tag::velocity, std::vector< tk::real >,
161 : : tag::pressure, tk::real,
162 : : tag::temperature, tk::real,
163 : : tag::materialid, std::size_t
164 : : > >
165 : : >,
166 : : tag::timedep, std::vector<
167 : : tk::TaggedTuple< brigand::list<
168 : : tag::sideset, std::vector< uint64_t >,
169 : : tag::fn, std::vector< tk::real >
170 : : > >
171 : : >,
172 : : tag::back_pressure, tk::TaggedTuple< brigand::list<
173 : : tag::sideset, std::vector< std::size_t >,
174 : : tag::pressure, tk::real
175 : : > >
176 : : > >;
177 : :
178 : : // IC box
179 : : using boxList = tk::TaggedTuple< brigand::list<
180 : : tag::materialid, std::size_t,
181 : : tag::volume, tk::real,
182 : : tag::mass, tk::real,
183 : : tag::density, tk::real,
184 : : tag::velocity, std::vector< tk::real >,
185 : : tag::pressure, tk::real,
186 : : tag::energy, tk::real,
187 : : tag::energy_content, tk::real,
188 : : tag::temperature, tk::real,
189 : : tag::mass_fractions, std::vector< tk::real >,
190 : : tag::xmin, tk::real,
191 : : tag::xmax, tk::real,
192 : : tag::ymin, tk::real,
193 : : tag::ymax, tk::real,
194 : : tag::zmin, tk::real,
195 : : tag::zmax, tk::real,
196 : : tag::orientation, std::vector< tk::real >,
197 : : tag::initiate, inciter::ctr::InitiateType,
198 : : tag::point, std::vector< tk::real >,
199 : : tag::init_time, tk::real,
200 : : tag::front_width, tk::real,
201 : : tag::front_speed, tk::real
202 : : > >;
203 : :
204 : : // IC meshblock
205 : : using meshblockList = tk::TaggedTuple< brigand::list<
206 : : tag::blockid, std::uint64_t,
207 : : tag::materialid, std::size_t,
208 : : tag::volume, tk::real,
209 : : tag::mass, tk::real,
210 : : tag::density, tk::real,
211 : : tag::velocity, std::vector< tk::real >,
212 : : tag::pressure, tk::real,
213 : : tag::energy, tk::real,
214 : : tag::energy_content, tk::real,
215 : : tag::temperature, tk::real,
216 : : tag::mass_fractions, std::vector< tk::real >,
217 : : tag::initiate, inciter::ctr::InitiateType,
218 : : tag::point, std::vector< tk::real >,
219 : : tag::init_time, tk::real,
220 : : tag::front_width, tk::real,
221 : : tag::front_speed, tk::real
222 : : > >;
223 : :
224 : : // Initial conditions (ic) block
225 : : using icList = tk::TaggedTuple< brigand::list<
226 : : tag::materialid, std::size_t,
227 : : tag::pressure, tk::real,
228 : : tag::temperature, tk::real,
229 : : tag::mass_fractions, std::vector< tk::real >,
230 : : tag::density, tk::real,
231 : : tag::energy, tk::real,
232 : : tag::velocity, std::vector< tk::real >,
233 : : tag::box, std::vector< boxList >,
234 : : tag::meshblock, std::vector< meshblockList >
235 : : > >;
236 : :
237 : : // Overset mesh block
238 : : using meshList = tk::TaggedTuple< brigand::list<
239 : : tag::filename, std::string,
240 : : tag::location, std::vector< tk::real >,
241 : : tag::orientation, std::vector< tk::real >,
242 : : tag::mass, tk::real,
243 : : tag::moment_of_inertia, tk::real,
244 : : tag::center_of_mass, std::vector< tk::real >
245 : : > >;
246 : :
247 : : // Field output block
248 : : using fieldOutputList = tk::TaggedTuple< brigand::list<
249 : : tag::interval, uint32_t,
250 : : tag::time_interval, tk::real,
251 : : tag::time_range, std::vector< tk::real >,
252 : : tag::refined, bool,
253 : : tag::filetype, tk::ctr::FieldFileType,
254 : : tag::sideset, std::vector< uint64_t >,
255 : : tag::outvar, std::vector< OutVar >,
256 : : tag::elemalias, std::vector< std::string >, // only for error checking
257 : : tag::elemvar, std::vector< std::string >, // only for error checking
258 : : tag::nodealias, std::vector< std::string >, // only for error checking
259 : : tag::nodevar, std::vector< std::string > // only for error checking
260 : : > >;
261 : :
262 : : // Diagnostics block
263 : : using diagnosticsList = tk::TaggedTuple< brigand::list<
264 : : tag::interval, uint32_t,
265 : : tag::error, tk::ctr::ErrorType,
266 : : tag::format, tk::ctr::TxtFloatFormatType,
267 : : tag::precision, std::streamsize
268 : : > >;
269 : :
270 : : // History output block
271 : : using historyOutputList = tk::TaggedTuple< brigand::list<
272 : : tag::interval, uint32_t,
273 : : tag::time_interval, tk::real,
274 : : tag::time_range, std::vector< tk::real >,
275 : : tag::format, tk::ctr::TxtFloatFormatType,
276 : : tag::precision, std::streamsize,
277 : : tag::point, std::vector<
278 : : tk::TaggedTuple< brigand::list<
279 : : tag::id, std::string,
280 : : tag::coord, std::vector< tk::real >
281 : : > >
282 : : >
283 : : > >;
284 : :
285 : : using ConfigMembers = brigand::list<
286 : :
287 : : tag::title, std::string,
288 : :
289 : : // Command line parameters
290 : : tag::cmd, CmdLine,
291 : :
292 : : // time stepping options
293 : : tag::nstep, uint64_t,
294 : : tag::term, tk::real,
295 : : tag::t0, tk::real,
296 : : tag::dt, tk::real,
297 : : tag::cfl, tk::real,
298 : : tag::ttyi, uint32_t,
299 : : tag::imex_runge_kutta, uint32_t,
300 : : tag::imex_maxiter, uint32_t,
301 : : tag::imex_reltol, tk::real,
302 : : tag::imex_abstol, tk::real,
303 : :
304 : : // steady-state solver options
305 : : tag::steady_state, bool,
306 : : tag::residual, tk::real,
307 : : tag::rescomp, uint32_t,
308 : :
309 : : // mesh partitioning and reordering/sorting choices
310 : : tag::partitioning, tk::ctr::PartitioningAlgorithmType,
311 : : tag::pelocal_reorder, bool,
312 : : tag::operator_reorder, bool,
313 : :
314 : : // discretization scheme choices
315 : : tag::scheme, SchemeType,
316 : : tag::ndof, std::size_t,
317 : : tag::rdof, std::size_t,
318 : : tag::flux, FluxType,
319 : : tag::lowspeed_kp, tk::real,
320 : : tag::lowspeed_ku, tk::real,
321 : :
322 : : // limiter options
323 : : tag::limiter, LimiterType,
324 : : tag::cweight, tk::real,
325 : : tag::shock_detector_coeff, tk::real,
326 : : tag::accuracy_test, bool,
327 : : tag::limsol_projection, bool,
328 : :
329 : : // PDE options
330 : : tag::ncomp, std::size_t,
331 : : tag::pde, PDEType,
332 : : tag::transport, transportList,
333 : : tag::compflow, compflowList,
334 : : tag::multimat, multimatList,
335 : : tag::multispecies, multispeciesList,
336 : :
337 : : // Dependent variable name
338 : : tag::depvar, std::vector< char >,
339 : :
340 : : tag::sys, std::map< std::size_t, std::size_t >,
341 : :
342 : : tag::material, std::vector< materialList >,
343 : :
344 : : tag::species, std::vector< speciesList >,
345 : :
346 : : tag::matidxmap, tk::TaggedTuple< brigand::list<
347 : : tag::eosidx, std::vector< std::size_t >,
348 : : tag::matidx, std::vector< std::size_t >,
349 : : tag::solidx, std::vector< std::size_t >
350 : : > >,
351 : :
352 : : // Conditions
353 : : tag::bc, std::vector< bcList >,
354 : : tag::ic, icList,
355 : : tag::mesh, std::vector< meshList >,
356 : : tag::transfer, std::vector< Transfer >,
357 : :
358 : : // Rigid-body motion solver
359 : : tag::rigid_body_motion, tk::TaggedTuple< brigand::list<
360 : : tag::rigid_body_movt, bool,
361 : : tag::rigid_body_dof, std::size_t,
362 : : tag::symmetry_plane, std::size_t
363 : : > >,
364 : :
365 : : // ALE block
366 : : // ---------------------------------------------------------------------------
367 : : tag::ale, tk::TaggedTuple< brigand::list<
368 : : tag::ale, bool,
369 : : tag::smoother, MeshVelocitySmootherType,
370 : : tag::mesh_velocity, MeshVelocityType,
371 : : tag::mesh_motion, std::vector< std::size_t >,
372 : : tag::meshforce, std::vector< tk::real >,
373 : : tag::dvcfl, tk::real,
374 : : tag::vortmult, tk::real,
375 : : tag::maxit, std::size_t,
376 : : tag::tolerance, tk::real,
377 : : tag::dirichlet, std::vector< std::size_t >,
378 : : tag::symmetry, std::vector< std::size_t >,
379 : : tag::move, std::vector<
380 : : tk::TaggedTuple< brigand::list<
381 : : tag::sideset, std::vector< uint64_t >,
382 : : tag::fntype, tk::ctr::UserTableType,
383 : : tag::fn, std::vector< tk::real >
384 : : > >
385 : : >
386 : : > >,
387 : :
388 : : // p-refinement block
389 : : // ---------------------------------------------------------------------------
390 : : tag::pref, tk::TaggedTuple< brigand::list<
391 : : tag::pref, bool,
392 : : tag::indicator, PrefIndicatorType,
393 : : tag::ndofmax, std::size_t,
394 : : tag::tolref, tk::real
395 : : > >,
396 : :
397 : : // AMR block
398 : : // ---------------------------------------------------------------------------
399 : : tag::amr, tk::TaggedTuple< brigand::list<
400 : : tag::amr, bool,
401 : : tag::t0ref, bool,
402 : : tag::dtref, bool,
403 : : tag::dtref_uniform, bool,
404 : : tag::dtfreq, std::size_t,
405 : : tag::maxlevels, std::size_t,
406 : : tag::initial, std::vector< AMRInitialType >,
407 : : tag::edgelist, std::vector< std::size_t >,
408 : : tag::coords, tk::TaggedTuple< brigand::list<
409 : : tag::xminus, tk::real,
410 : : tag::xplus, tk::real,
411 : : tag::yminus, tk::real,
412 : : tag::yplus, tk::real,
413 : : tag::zminus, tk::real,
414 : : tag::zplus, tk::real
415 : : > >,
416 : : tag::error, AMRErrorType,
417 : : tag::refvar, std::vector< char >,
418 : : tag::tol_refine, tk::real,
419 : : tag::tol_derefine, tk::real
420 : : > >,
421 : :
422 : : // Output options
423 : : tag::field_output, fieldOutputList,
424 : : tag::diagnostics, diagnosticsList,
425 : : tag::history_output, historyOutputList
426 : : >;
427 : :
428 : : // Class storing the Config params
429 : : class InputDeck : public tk::TaggedTuple< ConfigMembers > {
430 : :
431 : : public:
432 : : //! Set of tags to ignore when printing this InputDeck
433 : : using ignore = CmdLine::ignore;
434 : :
435 : : //! \brief Constructor: set defaults
436 : : //! \param[in] cl Previously parsed and store command line
437 : : //! \details Anything not set here is initialized by the compiler using the
438 : : //! default constructor for the corresponding type.
439 : 1348 : explicit InputDeck( const CmdLine& cl = {} ) {
440 : : // Set previously parsed command line
441 [ + - ]: 1348 : get< tag::cmd >() = cl;
442 : : // Default time stepping params
443 : 1348 : get< tag::dt >() = 0.0;
444 : 1348 : get< tag::cfl >() = 0.0;
445 : : // Default AMR settings
446 : 1348 : auto rmax =
447 : : std::numeric_limits< tk::real >::max() / 100;
448 : 1348 : get< tag::amr, tag::coords, tag::xminus >() = rmax;
449 : 1348 : get< tag::amr, tag::coords, tag::xplus >() = -rmax;
450 : 1348 : get< tag::amr, tag::coords, tag::yminus >() = rmax;
451 : 1348 : get< tag::amr, tag::coords, tag::yplus >() = -rmax;
452 : 1348 : get< tag::amr, tag::coords, tag::zminus >() = rmax;
453 : 1348 : get< tag::amr, tag::coords, tag::zplus >() = -rmax;
454 : :
455 : : // -----------------------------------------------------------------------
456 : : /*
457 : : Keyword vector
458 : : The following code generates a vector of keywords, for the sole purpose
459 : : of documentation and user-help (accessible via the --helpctr and
460 : : --helpkw cmdline arguments. If entries for a keyword are not added to
461 : : this vector, help will not be output for it. The entries follow the
462 : : function signature of the tk::entry_t constructor (defined in
463 : : Base/Types.hpp).
464 : : */
465 : : // -----------------------------------------------------------------------
466 : 2696 : std::set< tk::entry_t > keywords;
467 : :
468 : : keywords.insert({"inciter",
469 : : "Start configuration block for inciter",
470 : : R"(This keyword is used to select inciter. Inciter, is a continuum-realm
471 : : shock hydrodynamics tool, solving a system of PDEs. The entire control
472 [ + - ][ + - ]: 1348 : file must be enclosed within the inciter block)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
473 : :
474 : : keywords.insert({"title", "Title", R"(The title may be specified in
475 [ + - ][ + - ]: 1348 : the input file. It is optional.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
476 : :
477 : : // -----------------------------------------------------------------------
478 : : // time stepping options
479 : : // -----------------------------------------------------------------------
480 : :
481 : : keywords.insert({"nstep", "Set number of time steps to take",
482 : : R"(This keyword is used to specify the number of time steps to take in a
483 : : simulation. The number of time steps are used in conjunction with the
484 : : maximmum time specified by keyword 'term': the simulation stops whichever
485 : : is reached first. Both 'nstep' and 'term' can be left unspecified, in
486 [ + - ][ + - ]: 1348 : which case their default values are used. See also 'term'.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
487 : :
488 : : keywords.insert({"term", "Set maximum physical time to simulate",
489 : : R"(This keyword is used to specify the termination time in a simulation.
490 : : The termination time and number of time steps, specified by 'nstep', are
491 : : used in conjunction to determine when to stop a simulation: whichever is
492 : : reached first. Both 'nstep' and 'term' can be left unspecified, in which
493 [ + - ][ + - ]: 1348 : case their default values are used. See also 'nstep'.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
494 : :
495 : : keywords.insert({"t0", "Set starting non-dimensional time",
496 : : R"(This keyword is used to specify the starting time in a simulation.)",
497 [ + - ][ + - ]: 1348 : "real"});
[ + - ][ + - ]
[ + - ][ + - ]
498 : :
499 : : keywords.insert({"dt", "Select constant time step size",
500 : : R"(This keyword is used to specify the time step size that used as a
501 : : constant during simulation. Setting 'cfl' and 'dt' are mutually
502 [ + - ][ + - ]: 1348 : exclusive. If both 'cfl' and 'dt' are set, 'dt' wins.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
503 : :
504 : : keywords.insert({"cfl",
505 : : "Set the Courant-Friedrichs-Lewy (CFL) coefficient",
506 : : R"(This keyword is used to specify the CFL coefficient for
507 : : variable-time-step-size simulations. Setting 'cfl' and 'dt' are mutually
508 [ + - ][ + - ]: 1348 : exclusive. If both 'cfl' and 'dt' are set, 'dt' wins.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
509 : :
510 : : keywords.insert({"ttyi", "Set screen output interval",
511 : : R"(This keyword is used to specify the interval in time steps for screen
512 [ + - ][ + - ]: 1348 : output during a simulation.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
513 : :
514 : : keywords.insert({"imex_runge_kutta",
515 : : "Toggle use of IMplicit-EXplicit Runge-Kutta scheme",
516 : : R"(This keywords is used to turn IMEX integrator on/off for solid materials
517 : : in a multimat run. Plastic terms are integrated implicitly in time. This
518 : : flag will activate an Implicit-Explicit Runge-Kutta scheme to replace the
519 : : explicit one that is usually used. Scheme taken from Cavaglieri, D., &
520 : : Bewley, T. (2015). Low-storage implicit/explicit Runge–Kutta schemes for
521 : : the simulation of stiff high-dimensional ODE systems. Journal of
522 [ + - ][ + - ]: 1348 : Computational Physics, 286, 172-193.)", "uint 0/1"});
[ + - ][ + - ]
[ + - ][ + - ]
523 : :
524 : : keywords.insert({"imex_maxiter",
525 : : "Set maximum number of iterations for non-linear solver with IMEX-RK scheme",
526 : : R"(This keywords is used to specify the maximum number of iterations that
527 : : the non-linear solver uses to obtain the implicit unknowns within the
528 [ + - ][ + - ]: 1348 : Implicit-Explicit Runge-Kutta scheme.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
529 : :
530 : : keywords.insert({"imex_reltol",
531 : : "Set relative tolerance for non-linear solver with IMEX-RK scheme",
532 : : R"(This keywords is used to specify the relative tolerance that
533 : : the non-linear solver uses to obtain the implicit unknowns within the
534 [ + - ][ + - ]: 1348 : Implicit-Explicit Runge-Kutta scheme.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
535 : :
536 : : keywords.insert({"imex_abstol",
537 : : "Set absolute tolerance for non-linear solver with IMEX-RK scheme",
538 : : R"(This keywords is used to specify the absolute tolerance that
539 : : the non-linear solver uses to obtain the implicit unknowns within the
540 [ + - ][ + - ]: 1348 : Implicit-Explicit Runge-Kutta scheme.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
541 : :
542 : : // -----------------------------------------------------------------------
543 : : // steady-state solver options
544 : : // -----------------------------------------------------------------------
545 : :
546 : : keywords.insert({"steady_state", "March to steady state",
547 : : R"(This keyword is used indicate that local time stepping should be used
548 [ + - ][ + - ]: 1348 : to march towards a stationary solution.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
549 : :
550 : : keywords.insert({"residual",
551 : : "Set the convergence criterion for the residual to reach",
552 : : R"(This keyword is used to specify a convergence criterion for the local
553 : : time stepping marching to steady state, below which the simulation is
554 [ + - ][ + - ]: 1348 : considered converged.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
555 : :
556 : : keywords.insert({"rescomp",
557 : : "Equation system component index for convergence",
558 : : R"(This keyword is used to specify a single integer that is used to denote
559 : : the equation component index in the complete system of equations
560 : : configured, to use for the convergence criterion for local
561 [ + - ][ + - ]: 1348 : time stepping marching towards steady state.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
562 : :
563 : : // -----------------------------------------------------------------------
564 : : // mesh partitioning and reordering/sorting choices
565 : : // -----------------------------------------------------------------------
566 : :
567 : : keywords.insert({"partitioning",
568 : : "Select mesh partitioning algorithm",
569 : : R"(This keyword is used to select a mesh partitioning algorithm. See
570 : : Control/Options/PartitioningAlgorithm.hpp for valid options.)",
571 [ + - ][ + - ]: 1348 : "string"});
[ + - ][ + - ]
[ + - ][ + - ]
572 : :
573 : : keywords.insert({"rcb",
574 : : "Select recursive coordinate bisection mesh partitioner",
575 : : R"(This keyword is used to select the recursive coordinate bisection (RCB)
576 : : mesh partitioner. RCB is a geometry-based partitioner used to distribute
577 : : an input mesh among processing elements. See
578 [ + - ][ + - ]: 1348 : Control/Options/PartitioningAlgorithm.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
579 : :
580 : : keywords.insert({"rib",
581 : : "Select recursive inertial bisection mesh partitioner",
582 : : R"(This keyword is used to select the recursive inertial bisection (RIB)
583 : : mesh partitioner. RIB is a geometry-based partitioner used to distribute
584 : : an input mesh among processing elements. See
585 [ + - ][ + - ]: 1348 : Control/Options/PartitioningAlgorithm.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
586 : :
587 : : keywords.insert({"hsfc",
588 : : "Select Hilbert Space Filling Curve (HSFC) mesh partitioner",
589 : : R"(This keyword is used to select the Hilbert Space Filling Curve (HSFC)
590 : : mesh partitioner. HSFC is a geometry-based partitioner used to distribute
591 : : an input mesh among processing elements. See
592 [ + - ][ + - ]: 1348 : Control/Options/PartitioningAlgorithm.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
593 : :
594 : : keywords.insert({"phg",
595 : : "Select parallel hypergraph mesh partitioner",
596 : : R"(This keyword is used to select the parallel hypergraph (PHG)
597 : : mesh partitioner. PHG is a graph-based partitioner used to distribute an
598 : : input mesh among processing elements. See
599 [ + - ][ + - ]: 1348 : Control/Options/PartitioningAlgorithm.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
600 : :
601 : : keywords.insert({"mj",
602 : : "Select multi-jagged (MJ) mesh partitioner",
603 : : R"(This keyword is used to select the multi-jagged (MJ) mesh partitioner.
604 : : MJ is a geometry-based partitioner used to distribute an input mesh among
605 : : processing elements. See
606 [ + - ][ + - ]: 1348 : Control/Options/PartitioningAlgorithm.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
607 : :
608 : : keywords.insert({"pelocal_reorder",
609 : : "PE-local reorder",
610 : : R"(This keyword is used in inciter as a keyword in the inciter...end block
611 : : as "pelocal_reorder true" (or false) to do (or not do) a global
612 : : distributed mesh reordering across all PEs that yields an approximately
613 : : continuous mesh node ID order as mesh partitions are assigned to PEs after
614 [ + - ][ + - ]: 1348 : mesh partitioning. This reordering is optional.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
615 : :
616 : : keywords.insert({"operator_reorder",
617 : : "Operator-access reorder",
618 : : R"(This keyword is used in inciter as a keyword in the inciter...end block
619 : : as "operator_reorder on" (or off) to do (or not do) a local mesh node
620 : : reordering based on the PDE operator access pattern. This reordering is
621 [ + - ][ + - ]: 1348 : optional.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
622 : :
623 : : // -----------------------------------------------------------------------
624 : : // discretization scheme choices
625 : : // -----------------------------------------------------------------------
626 : :
627 : : keywords.insert({"scheme", "Select discretization scheme",
628 : : R"(This keyword is used to select a spatial discretization scheme,
629 : : necessarily connected to the temporal discretization scheme. See
630 [ + - ][ + - ]: 1348 : Control/Inciter/Options/Scheme.hpp for valid options.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
631 : :
632 : : keywords.insert({"alecg",
633 : : "Select continuous Galerkin with ALE + Runge-Kutta",
634 : : R"(This keyword is used to select the continuous Galerkin finite element
635 : : scheme in the arbitrary Lagrangian-Eulerian (ALE) reference frame combined
636 : : with Runge-Kutta (RK) time stepping.
637 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Scheme.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
638 : :
639 : : keywords.insert({"oversetfe",
640 : : "Select continuous Galerkin finite element with overset meshes + "
641 : : "Runge-Kutta",
642 : : R"(This keyword is used to select the continuous Galerkin finite element
643 : : scheme with Runge-Kutta (RK) time stepping, combined with overset grids.
644 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Scheme.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
645 : :
646 : : keywords.insert({"dg",
647 : : "Select 1st-order discontinuous Galerkin discretization + Runge-Kutta",
648 : : R"(This keyword is used to select the first-order accurate discontinuous
649 : : Galerkin, DG(P0), spatial discretiztaion used in Inciter. As this is first
650 : : order accurate, it is intended for testing and debugging purposes only.
651 : : Selecting this spatial discretization also selects the Runge-Kutta scheme
652 : : for time discretization. See Control/Inciter/Options/Scheme.hpp for other
653 [ + - ][ + - ]: 1348 : valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
654 : :
655 : : keywords.insert({"p0p1",
656 : : "Select 2nd-order finite volume discretization + Runge-Kutta",
657 : : R"(This keyword is used to select the second-order accurate finite volume,
658 : : P0P1, spatial discretiztaion used in Inciter. This method uses a
659 : : least-squares procedure to reconstruct the second-order solution from the
660 : : first-order one. Selecting this spatial discretization also selects the
661 : : Runge-Kutta scheme for time discretization.
662 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Scheme.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
663 : :
664 : : keywords.insert({"dgp1",
665 : : "Select 2nd-order discontinuous Galerkin discretization + Runge-Kutta",
666 : : R"(This keyword is used to select the second-order accurate discontinuous
667 : : Galerkin, DG(P1), spatial discretiztaion used in Inciter. Selecting this
668 : : spatial discretization also selects the Runge-Kutta scheme for time
669 : : discretization. See Control/Inciter/Options/Scheme.hpp for other
670 [ + - ][ + - ]: 1348 : valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
671 : :
672 : : keywords.insert({"dgp2",
673 : : "Select 3nd-order discontinuous Galerkin discretization + Runge-Kutta",
674 : : R"(This keyword is used to select the third-order accurate discontinuous
675 : : Galerkin, DG(P2), spatial discretiztaion used in Inciter. Selecting this
676 : : spatial discretization also selects the Runge-Kutta scheme for time
677 : : discretization. See Control/Inciter/Options/Scheme.hpp for other
678 [ + - ][ + - ]: 1348 : valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
679 : :
680 : : keywords.insert({"pdg",
681 : : "Select p-adaptive discontinuous Galerkin discretization + Runge-Kutta",
682 : : R"(This keyword is used to select the polynomial adaptive discontinuous
683 : : Galerkin spatial discretizaion used in Inciter. Selecting this spatial
684 : : discretization also selects the Runge-Kutta scheme for time
685 : : discretization. See Control/Inciter/Options/Scheme.hpp for other valid
686 [ + - ][ + - ]: 1348 : options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
687 : :
688 : : keywords.insert({"fv",
689 : : "Select 2nd-order finite volume discretization + Runge-Kutta",
690 : : R"(This keyword is used to select the second-order accurate finite volume,
691 : : P0P1, spatial discretiztaion used in Inciter. This method uses a
692 : : least-squares procedure to reconstruct the second-order solution from the
693 : : first-order one. See Control/Inciter/Options/Scheme.hpp for other valid
694 [ + - ][ + - ]: 1348 : options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
695 : :
696 : : keywords.insert({"ndof", "Number of evolved solution DOFs",
697 [ + - ][ + - ]: 1348 : R"(The number of solution DOFs that are evolved.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
698 : :
699 : : keywords.insert({"rdof", "Total number of solution DOFs",
700 : : R"(The total number of solution DOFs, including the reconstructed and the
701 [ + - ][ + - ]: 1348 : evolved ones.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
702 : :
703 : : // -----------------------------------------------------------------------
704 : : // limiter options
705 : : // -----------------------------------------------------------------------
706 : :
707 : : keywords.insert({"limiter", "Select limiter function",
708 : : R"(This keyword is used to select a limiter function, used for
709 : : discontinuous Galerkin (DG) spatial discretization used in inciter. See
710 [ + - ][ + - ]: 1348 : Control/Inciter/Options/Limiter.hpp for valid options.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
711 : :
712 : : keywords.insert({"nolimiter", "No limiter used",
713 : : R"(This keyword is used for discontinuous Galerkin (DG) spatial
714 : : discretization without any limiter in inciter. See
715 [ + - ][ + - ]: 1348 : Control/Inciter/Options/Limiter.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
716 : :
717 : : keywords.insert({"wenop1",
718 : : "Select the Weighted Essentially Non-Oscillatory (WENO) limiter for DGP1",
719 : : R"(This keyword is used to select the Weighted Essentially Non-Oscillatory
720 : : limiter used for discontinuous Galerkin (DG) P1 spatial discretization
721 : : used in inciter. See Control/Inciter/Options/Limiter.hpp for other valid
722 [ + - ][ + - ]: 1348 : options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
723 : :
724 : : keywords.insert({"cweight",
725 : : "Set value for central linear weight used by WENO, cweight",
726 : : R"(This keyword is used to set the central linear weight used for the
727 : : central stencil in the Weighted Essentially Non-Oscillatory (WENO) limiter
728 [ + - ][ + - ]: 1348 : for discontinuous Galerkin (DG) methods.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
729 : :
730 : : keywords.insert({"superbeep1",
731 : : "Select the Superbee limiter for DGP1",
732 : : R"(This keyword is used to select the Superbee limiter used for
733 : : discontinuous Galerkin (DG) P1 spatial discretization used in inciter.
734 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Limiter.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
735 : :
736 : : keywords.insert({"shock_detector_coeff",
737 : : "Configure the coefficient used in shock indicator",
738 : : R"(This keyword can be used to configure the coefficient used in the
739 [ + - ][ + - ]: 1348 : threshold calculation for the shock indicator.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
740 : :
741 : : keywords.insert({"vertexbasedp1",
742 : : "Select the vertex-based limiter for DGP1",
743 : : R"(This keyword is used to select the vertex-based limiter used for
744 : : discontinuous Galerkin (DG) P1 spatial discretization used in inciter.
745 : : Ref. Kuzmin, D. (2010). A vertex-based hierarchical slope limiter for
746 : : p-adaptive discontinuous Galerkin methods. Journal of computational and
747 : : applied mathematics, 233(12), 3077-3085.
748 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Limiter.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
749 : :
750 : : keywords.insert({"accuracy_test", "Toggle accuracy test setup",
751 : : R"(This keyword is used to specify if the current setup is for an
752 : : order-of-accuracy testing, used for discontinuous Galerkin (DG) spatial
753 : : discretization in inciter. This deactivates certain robustness corrections
754 : : which might impact order-of-accuracy. Only intended for simple test
755 [ + - ][ + - ]: 1348 : problems and not for real problems.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
756 : :
757 : : keywords.insert({"limsol_projection",
758 : : "Toggle limited solution projection",
759 : : R"(This keyword is used to specify limited solution projection.
760 : : This is used for discontinuous Galerkin (DG) spatial discretization in
761 : : inciter, for multi-material hydrodynamics. This uses a projection to
762 : : obtain bulk momentum and material energies from the limited primitive
763 : : quantities. This step is essential to obtain closure-law obeying limited
764 : : quantities. See Pandare et al. (2023). On the Design of Stable,
765 : : Consistent, and Conservative High-Order Methods for Multi-Material
766 [ + - ][ + - ]: 1348 : Hydrodynamics. J Comp Phys (490).)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
767 : :
768 : : // -----------------------------------------------------------------------
769 : : // flux options
770 : : // -----------------------------------------------------------------------
771 : :
772 : : keywords.insert({"flux", "Select flux function",
773 : : R"(This keyword is used to select a flux function, used for
774 : : discontinuous Galerkin (DG) spatial discretization used in inciter. See
775 [ + - ][ + - ]: 1348 : Control/Inciter/Options/Flux.hpp for valid options.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
776 : :
777 : : keywords.insert({"laxfriedrichs",
778 : : "Select Lax-Friedrichs flux function",
779 : : R"(This keyword is used to select the Lax-Friedrichs flux function used
780 : : for discontinuous Galerkin (DG) spatial discretization used in inciter.
781 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Flux.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
782 : :
783 : : keywords.insert({"hllc",
784 : : "Select the Harten-Lax-van Leer-Contact (HLLC) flux function",
785 : : R"(This keyword is used to select the Harten-Lax-van Leer-Contact flux
786 : : function used for discontinuous Galerkin (DG) spatial discretization
787 : : used in inciter. See Control/Inciter/Options/Flux.hpp for other valid
788 [ + - ][ + - ]: 1348 : options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
789 : :
790 : : keywords.insert({"upwind", "Select the upwind flux function",
791 : : R"(This keyword is used to select the upwind flux
792 : : function used for discontinuous Galerkin (DG) spatial discretization
793 : : used in inciter. It is only usable for scalar transport.
794 [ + - ][ + - ]: 1348 : See Control/Inciter/Options/Flux.hpp for other valid options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
795 : :
796 : : keywords.insert({"ausm",
797 : : "Select the Advection Upstream Splitting Method (AUSM) flux function",
798 : : R"(This keyword is used to select the AUSM flux
799 : : function used for discontinuous Galerkin (DG) spatial discretization
800 : : used in inciter. It is only set up for for multi-material hydro, and
801 [ + - ][ + - ]: 1348 : not selectable for anything else.)"});
[ + - ][ + - ]
[ + - ][ + - ]
802 : :
803 : : keywords.insert({"ldfss",
804 : : "Select the Low Diffusion Flux Splitting Scheme (LDFSS)",
805 : : R"(This keyword is used to select the LDFSS flux
806 : : function used for discontinuous Galerkin (DG) spatial discretization
807 : : used in inciter. It is only set up for for multi-material hydro, and
808 [ + - ][ + - ]: 1348 : not selectable for anything else.)"});
[ + - ][ + - ]
[ + - ][ + - ]
809 : :
810 : : keywords.insert({"lowspeed_kp",
811 : : "Select the low-speed coefficient K_p in the AUSM+up flux function",
812 : : R"(This keyword is used to select the low-speed coefficient K_p in the
813 : : AUSM+up flux function used for the DG or FV spatial discretization for
814 : : multi-material hydro, and not used for anything else. The default
815 : : value is 0, and recommended value for low speed flows (Mach < 0.1) is
816 [ + - ][ + - ]: 1348 : 1.)"});
[ + - ][ + - ]
[ + - ][ + - ]
817 : :
818 : : keywords.insert({"lowspeed_ku",
819 : : "Select the low-speed coefficient K_u in the AUSM+up flux function",
820 : : R"(This keyword is used to select the low-speed coefficient K_u in the
821 : : AUSM+up flux function used for the DG or FV spatial discretization for
822 : : multi-material hydro, and not used for anything else. The default
823 : : value is 1, and recommended value for low speed flows (Mach < 0.1) is
824 [ + - ][ + - ]: 1348 : 1.)"});
[ + - ][ + - ]
[ + - ][ + - ]
825 : :
826 : : keywords.insert({"hll",
827 : : "Select the Harten-Lax-vanLeer (HLL) flux function",
828 : : R"(This keyword is used to select the HLL flux
829 : : function used for discontinuous Galerkin (DG) spatial discretization
830 : : used in inciter. It is only set up for for multi-material hydro, and
831 [ + - ][ + - ]: 1348 : not selectable for anything else.)"});
[ + - ][ + - ]
[ + - ][ + - ]
832 : :
833 : : keywords.insert({"hlld",
834 : : "Select the Harten-Lax-vanLeer-Discontinuities (HLLD) flux function",
835 : : R"(This keyword is used to select the HLLD flux
836 : : function used for discontinuous Galerkin (DG) spatial discretization
837 : : used in inciter. It is only set up for for multi-material runs. This
838 : : flux is designed to handle normal and shear waves within solid
839 [ + - ][ + - ]: 1348 : materials)"});
[ + - ][ + - ]
[ + - ][ + - ]
840 : :
841 : : // -----------------------------------------------------------------------
842 : : // PDE keywords
843 : : // -----------------------------------------------------------------------
844 : :
845 : : keywords.insert({"problem",
846 : : "Specify problem configuration for partial differential equation solver",
847 : : R"(This keyword is used to specify the problem configuration for the
848 [ + - ][ + - ]: 1348 : partial differential equation solver.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
849 : :
850 : : keywords.insert({"transport",
851 : : "Start configuration block for an transport equation",
852 : : R"(This keyword is used to introduce a transport block, used to
853 : : specify the configuration for a transport equation type.)",
854 [ + - ][ + - ]: 1348 : "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
855 : :
856 : : keywords.insert({"ncomp",
857 : : "Set number of scalar components for a system of transport equations",
858 : : R"(This keyword is used to specify the number of scalar
859 [ + - ][ + - ]: 1348 : components of transport (linear advection) equations.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
860 : :
861 : : keywords.insert({"compflow",
862 : : "Start configuration block for the compressible flow equations",
863 : : R"(This keyword is used to introduce the compflow block, used to
864 : : specify the configuration for a system of partial differential equations,
865 [ + - ][ + - ]: 1348 : governing single material compressible fluid flow.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
866 : :
867 : : keywords.insert({"multimat",
868 : : "Start configuration block for the compressible multi-material equations",
869 : : R"(This keyword is used to introduce the multimat block,
870 : : used to specify the configuration for a system of partial differential
871 : : equations, governing compressible multi-material hydrodynamics assuming
872 [ + - ][ + - ]: 1348 : velocity equilibrium (single velocity).)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
873 : :
874 : : keywords.insert({"multispecies",
875 : : "Start configuration block for the compressible multi-species equations",
876 : : R"(This keyword is used to introduce the multispecies block,
877 : : used to specify the configuration for a system of partial differential
878 : : equations, governing compressible multi-species fluid dynamics.)",
879 [ + - ][ + - ]: 1348 : "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
880 : :
881 : : keywords.insert({"nmat",
882 : : "Set number of materials for the multi-material system",
883 : : R"(This keyword is used to specify the number of materials for
884 [ + - ][ + - ]: 1348 : multi-material flow, see also the keyword 'multimat'.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
885 : :
886 : : keywords.insert({"nspec",
887 : : "Set number of species for the multi-species system",
888 : : R"(This keyword is used to specify the number of species for
889 [ + - ][ + - ]: 1348 : multi-species flow, see also the keyword 'multispecies'.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
890 : :
891 : : keywords.insert({"prelax",
892 : : "Toggle multi-material finite pressure relaxation",
893 : : R"(This keyword is used to turn finite pressure relaxation between
894 : : multiple materials on/off. It is used only for the multi-material solver,
895 [ + - ][ + - ]: 1348 : and has no effect when used for the other PDE types.)", "uint 0/1"});
[ + - ][ + - ]
[ + - ][ + - ]
896 : :
897 : : keywords.insert({"prelax_timescale",
898 : : "Time-scale for multi-material finite pressure relaxation",
899 : : R"(This keyword is used to specify the time-scale at which finite pressure
900 : : relaxation between multiple materials occurs. The default value of 0.25
901 : : corresponds to a relaxation time that is 4 times the time required for a
902 : : sound wave to pass through a computational element. It is used only for
903 [ + - ][ + - ]: 1348 : multimat, and has no effect for the other PDE types.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
904 : :
905 : : keywords.insert({"intsharp",
906 : : "Toggle multi-material interface sharpening",
907 : : R"(This keyword is used to turn interface sharpening on/off. It uses the
908 : : multi-material THINC interface reconstruction.
909 : : Ref. Pandare A. K., Waltz J., & Bakosi J. (2021) Multi-Material
910 : : Hydrodynamics with Algebraic Sharp Interface Capturing. Computers &
911 : : Fluids, doi: https://doi.org/10.1016/j.compfluid.2020.104804. It is used
912 : : for the multi-material and the transport solver, and has no effect when
913 [ + - ][ + - ]: 1348 : used for the other PDE types.)", "uint 0/1"});
[ + - ][ + - ]
[ + - ][ + - ]
914 : :
915 : : keywords.insert({"intsharp_param",
916 : : "Parameter for multi-material interface sharpening",
917 : : R"(This keyword is used to specify the parameter for the interface
918 : : sharpening. This parameter affects how many cells the material interfaces
919 : : span, after the use of sharpening. It is used for multimat and transport,
920 [ + - ][ + - ]: 1348 : and has no effect for the other PDE types.)", "real" });
[ + - ][ + - ]
[ + - ][ + - ]
921 : :
922 : : keywords.insert({"rho0constraint",
923 : : "Toggle the density constraint correction",
924 : : R"(This keyword is used to toggle the density constraint in solid
925 : : dynamics on/off. It is used only for the multi-material solver in the
926 [ + - ][ + - ]: 1348 : presence of solids. The default is 1 (on).)", "uint 0/1"});
[ + - ][ + - ]
[ + - ][ + - ]
927 : :
928 : : keywords.insert({"dt_sos_massavg",
929 : : "Toggle method for calculating speed of sound used for time step in a cell",
930 : : R"(This keyword is used to specify the method to calculate the speed of
931 : : sound in a cell used for the time step. If set to 1, the speed of sound
932 : : will be calculated using the mass average, rather than the maximum value
933 : : across materials. It is used for multimat, and has no effect for the
934 [ + - ][ + - ]: 1348 : other PDE types.)", "uint 0/1" });
[ + - ][ + - ]
[ + - ][ + - ]
935 : :
936 : : // Dependent variable name
937 : : keywords.insert({"depvar",
938 : : "Select dependent variable name for PDE.",
939 [ + - ][ + - ]: 1348 : R"(Select dependent variable name for PDE.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
940 : :
941 : : // -----------------------------------------------------------------------
942 : : // physics choices
943 : : // -----------------------------------------------------------------------
944 : :
945 : : keywords.insert({"physics",
946 : : "Specify the physics configuration for a system of PDEs",
947 : : R"(This keyword is used to select the physics configuration for a
948 : : particular PDE system. Valid options depend on the system of PDEs in
949 [ + - ][ + - ]: 1348 : which the keyword is used.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
950 : :
951 : : keywords.insert({"advection",
952 : : "Specify the advection physics",
953 : : R"(This keyword is used to select the advection physics for the transport
954 [ + - ][ + - ]: 1348 : PDE system. Only usable for 'transport'.)"});
[ + - ][ + - ]
[ + - ][ + - ]
955 : :
956 : : keywords.insert({"advdiff",
957 : : "Specify the advection + diffusion physics",
958 : : R"(This keyword is used to select the advection + diffusion physics
959 [ + - ][ + - ]: 1348 : for transport PDEs. Only usable for 'transport'.)"});
[ + - ][ + - ]
[ + - ][ + - ]
960 : :
961 : : keywords.insert({"euler",
962 : : "Specify the Euler (inviscid) compressible flow physics",
963 : : R"(This keyword is used to select the Euler (inviscid) compressible
964 [ + - ][ + - ]: 1348 : flow physics configuration. Usable for 'compflow' and 'multimat')"});
[ + - ][ + - ]
[ + - ][ + - ]
965 : :
966 : : keywords.insert({"energy_pill",
967 : : "Specify the energy pill physics",
968 : : R"(This keyword is used to select an energy pill initialization as physics
969 : : configuration for multiple material compressible flow. Parameters for the
970 : : linearly traveling front are required to be specified when energy_pill is
971 : : selected. See 'linear' for more details. Currently setup only for
972 [ + - ][ + - ]: 1348 : 'multimat')"});
[ + - ][ + - ]
[ + - ][ + - ]
973 : :
974 : : // -----------------------------------------------------------------------
975 : : // material/eos object
976 : : // -----------------------------------------------------------------------
977 : :
978 : : keywords.insert({"material",
979 : : "Start configuration block for material (eos) properties",
980 : : R"(This keyword is used to introduce a material block, used to
981 [ + - ][ + - ]: 1348 : specify material properties.)", "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
982 : :
983 : : keywords.insert({"species",
984 : : "Start configuration block for species (eos) properties",
985 : : R"(This keyword is used to introduce a species block, used to
986 : : specify species properties in the multi species solver.)",
987 [ + - ][ + - ]: 1348 : "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
988 : :
989 : : keywords.insert({"id", "ID",
990 : : R"(This keyword is used to specify an ID, a positive integer. Usage is
991 : : context specific, i.e. what block it is specified in. E.g. Inside the
992 : : material block, it is used to specify a block consisting of IDs of
993 [ + - ][ + - ]: 1348 : materials of that EOS type)", "vector of uints"});
[ + - ][ + - ]
[ + - ][ + - ]
994 : :
995 : : keywords.insert({"eos", "Select equation of state (type)",
996 : : R"(This keyword is used to select an equation of state for a material.)",
997 [ + - ][ + - ]: 1348 : "string"});
[ + - ][ + - ]
[ + - ][ + - ]
998 : :
999 : : keywords.insert({"gamma", "ratio of specific heats",
1000 : : R"(This keyword is used to specify the material property, ratio of
1001 [ + - ][ + - ]: 1348 : specific heats.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1002 : :
1003 : : keywords.insert({"pstiff", "EoS stiffness parameter",
1004 : : R"(This keyword is used to specify the material property, stiffness
1005 [ + - ][ + - ]: 1348 : parameter in the stiffened gas equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1006 : :
1007 : : keywords.insert({"w_gru", "Grueneisen coefficient",
1008 : : R"(This keyword is used to specify the material property, Gruneisen
1009 : : coefficient for the Jones-Wilkins-Lee equation of state.)",
1010 [ + - ][ + - ]: 1348 : "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1011 : :
1012 : : keywords.insert({"A_jwl", "JWL EoS A parameter",
1013 : : R"(This keyword is used to specify the material property A (units: Pa)
1014 [ + - ][ + - ]: 1348 : for the Jones-Wilkins-Lee equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1015 : :
1016 : : keywords.insert({"B_jwl", "JWL EoS B parameter",
1017 : : R"(This keyword is used to specify the material property B (units: Pa)
1018 [ + - ][ + - ]: 1348 : for the Jones-Wilkins-Lee equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1019 : :
1020 : : keywords.insert({"C_jwl", "JWL EoS C parameter",
1021 : : R"(This keyword is used to specify the material property C (units: Pa)
1022 [ + - ][ + - ]: 1348 : for the Jones-Wilkins-Lee equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1023 : :
1024 : : keywords.insert({"R1_jwl", "JWL EoS R1 parameter",
1025 : : R"(This keyword is used to specify the material property R1 for the
1026 [ + - ][ + - ]: 1348 : Jones-Wilkins-Lee equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1027 : :
1028 : : keywords.insert({"R2_jwl", "JWL EoS R2 parameter",
1029 : : R"(This keyword is used to specify the material property R2 for the
1030 [ + - ][ + - ]: 1348 : Jones-Wilkins-Lee equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1031 : :
1032 : : keywords.insert({"rho0_jwl", "JWL EoS rho0 parameter",
1033 : : R"(This keyword is used to specify the material property rho0, which is
1034 : : the density of initial state (units: kg/m3) for the Jones-Wilkins-Lee
1035 [ + - ][ + - ]: 1348 : equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1036 : :
1037 : : keywords.insert({"de_jwl", "JWL EoS de parameter",
1038 : : R"(This keyword is used to specify the material property de, which is the
1039 : : heat of detonation for products; and for reactants, it is chosen such that
1040 : : the ambient internal energy (e0) is 0 (units: J/kg). Used for the
1041 [ + - ][ + - ]: 1348 : Jones-Wilkins-Lee equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1042 : :
1043 : : keywords.insert({"rhor_jwl", "JWL EoS rhor parameter",
1044 : : R"(This keyword is used to specify the material property rhor, which is
1045 : : the density of reference state (units: kg/m3) for the Jones-Wilkins-Lee
1046 [ + - ][ + - ]: 1348 : equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1047 : :
1048 : : keywords.insert({"Tr_jwl", "JWL EoS Tr parameter",
1049 : : R"(This keyword is used to specify the material property Tr, which is the
1050 : : temperature of reference state (units: K) for the Jones-Wilkins-Lee
1051 [ + - ][ + - ]: 1348 : equation of state.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1052 : :
1053 : : keywords.insert({"Pr_jwl", "JWL EoS er parameter",
1054 : : R"(This keyword is used to specify the material property Pr, which is the
1055 : : pressure at the reference state (units: Pa) for the Jones-Wilkins-Lee
1056 : : equation of state. It is used to calculate the reference temperature for
1057 [ + - ][ + - ]: 1348 : the EoS.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1058 : :
1059 : : keywords.insert({"mu", "shear modulus/dynamic viscosity",
1060 : : R"(This keyword is used to specify the material property, shear modulus
1061 [ + - ][ + - ]: 1348 : for solids, or dynamic viscosity for fluids.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1062 : :
1063 : : keywords.insert({"yield_stress", "Yield stress of solid material",
1064 : : R"(This keyword is used to specify the material property yield stress,
1065 : : which indicates the stress (units: Pa) after which the material begins
1066 [ + - ][ + - ]: 1348 : plastic flow.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1067 : :
1068 : : keywords.insert({"alpha", "alpha parameter for Godunov-Romenski EOS",
1069 : : R"(This keyword is used to specify the alpha parameter for
1070 [ + - ][ + - ]: 1348 : Godunov-Romenski EOS for solids.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1071 : :
1072 : : keywords.insert({"K0", "K0 parameter for Godunov-Romenski EOS",
1073 : : R"(This keyword is used to specify the K0 parameter for
1074 [ + - ][ + - ]: 1348 : Godunov-Romenski EOS for solids.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1075 : :
1076 : : keywords.insert({"cv", "specific heat at constant volume",
1077 : : R"(This keyword is used to specify the material property, specific heat at
1078 [ + - ][ + - ]: 1348 : constant volume.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1079 : :
1080 : : keywords.insert({"k", "heat conductivity",
1081 : : R"(This keyword is used to specify the material property, heat
1082 [ + - ][ + - ]: 1348 : conductivity.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1083 : :
1084 : : keywords.insert({"cp_coeff", "specific heat coefficients for TPG",
1085 : : R"(This keyword is used to specify species' coefficients in the
1086 : : thermally perfect gas polynomial fit (per temperature range)
1087 : : for specific heat at constant volume. The outer vector is per
1088 : : species, the middle vector is per temperature range, and the
1089 : : inner vector contains 8 coefficients to describe the polynomial.)",
1090 [ + - ][ + - ]: 1348 : "vector of vector of vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1091 : :
1092 : : keywords.insert({"t_range", "temperature range for TPG specific heat polynomials",
1093 : : R"(This keyword is used to specify the temperature range for each specific
1094 : : heat polynomial given in cp_coeff. The outer vector is per species,
1095 : : and the inner vector gives the temperature ranges. The inner vector must
1096 : : be sized 1 larger than the number of temperature ranges, as the ranges
1097 [ + - ][ + - ]: 1348 : are read as [T0, T1],[T1, T2],...)", "vector of vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1098 : :
1099 : : keywords.insert({"dH_ref", "reference enthalpy for TPG specific heat polynomials",
1100 : : R"(This keyword is used to specify the reference enthalpy at 298.15 K so that
1101 : : the reference temperature in the enthalpy calculations is 0 K. This number
1102 : : is taken from the NASA Glenn 2002 report, and is the heat of formation
1103 [ + - ][ + - ]: 1348 : divided by the species molar mass.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1104 : :
1105 : : keywords.insert({"R", "Specific gas constant",
1106 : : R"(This keyword is used to specify the species property, specific gas
1107 [ + - ][ + - ]: 1348 : constant, in units J/kg.K.)", "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1108 : :
1109 : : keywords.insert({"stiffenedgas",
1110 : : "Select the stiffened gas equation of state",
1111 [ + - ][ + - ]: 1348 : R"(This keyword is used to select the stiffened gas equation of state.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1112 : :
1113 : : keywords.insert({"jwl", "Select the JWL equation of state",
1114 : : R"(This keyword is used to select the Jones, Wilkins, Lee equation of
1115 [ + - ][ + - ]: 1348 : state.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1116 : :
1117 : : keywords.insert({"smallshearsolid",
1118 : : "Select the SMALLSHEARSOLID equation of state",
1119 : : R"(This keyword is used to select the small shear strain equation of state
1120 : : for solids. This EOS uses a small-shear approximation for the elastic
1121 : : contribution, and a stiffened gas EOS for the hydrodynamic contribution of
1122 : : the internal energy See Plohr, J. N., & Plohr, B. J. (2005). Linearized
1123 : : analysis of Richtmyer–Meshkov flow for elastic materials. Journal of Fluid
1124 [ + - ][ + - ]: 1348 : Mechanics, 537, 55-89 for further details.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1125 : :
1126 : : keywords.insert({"wilkins_aluminum",
1127 : : "Select Wilkins' equation of state for aluminum",
1128 : : R"(This keyword is used to select Wilkin's equation of state for solids
1129 : : and a hydro EoS for aluminum. These functions were taken from Example 4
1130 : : of Barton, Philip T. "An interface-capturing Godunov method
1131 : : for the simulation of compressible solid-fluid problems." Journal
1132 [ + - ][ + - ]: 1348 : of Computational Physics 390 (2019): 25-50.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1133 : :
1134 : : keywords.insert({"godunovromenski",
1135 : : "Select godunovromenski equation of state for solids",
1136 : : R"(This keyword is used to select Godunov-Romenski equation of state
1137 : : for solids. These functions were taken from Example 1
1138 : : of Barton, Philip T. "An interface-capturing Godunov method
1139 : : for the simulation of compressible solid-fluid problems." Journal
1140 [ + - ][ + - ]: 1348 : of Computational Physics 390 (2019): 25-50.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1141 : :
1142 : : keywords.insert({"matidxmap",
1143 : : "AUTO-GENERATED Material index map for EOS",
1144 : : R"(The following AUTO-GENERATED data structure is used to index into the
1145 : : correct material vector entry. This is done using the following three maps:
1146 : : 1. eosidx: This vector provides the eos-index (value) in the
1147 : : vector<tag::material> for the given user-spec material id (index).
1148 : : 2. matidx: This vector provides the material-index (value) inside the
1149 : : vector<tag::material>[eosidx] block for the given user-specified
1150 : : material id (index).
1151 : : 3. solidx: This vector provides the solid-index (value) assigned to
1152 [ + - ][ + - ]: 1348 : the given user-specified material id (index). It is 0 for fluids.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1153 : :
1154 : : // -----------------------------------------------------------------------
1155 : : // output object
1156 : : // -----------------------------------------------------------------------
1157 : :
1158 : : keywords.insert({"field_output",
1159 : : "Start of field_output input block",
1160 : : R"(This keyword is used to start a block in the input file containing the
1161 [ + - ][ + - ]: 1348 : list and settings of requested field output.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1162 : :
1163 : : keywords.insert({"interval",
1164 : : "Set interval (in units of iteration count)",
1165 : : R"(This keyword is used to specify an interval in units of iteration count
1166 : : (i.e., number of time steps). This must be used within a relevant
1167 [ + - ][ + - ]: 1348 : block.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1168 : :
1169 : : keywords.insert({"time_interval",
1170 : : "Set interval (in units of physics time)",
1171 : : R"(This keyword is used to specify an interval in units of physics time.
1172 [ + - ][ + - ]: 1348 : This must be used within a relevant block.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1173 : :
1174 : : keywords.insert({"time_range",
1175 : : "Configure physics time range for output (in units of physics time)",
1176 : : R"(This keyword is used to configure field-, or history-output, specifying
1177 : : a start time, a stop time, and an output frequency in physics time units.
1178 : : Example: 'time_range = {0.2, 0.3, 0.001}', which specifies that from t=0.2 to
1179 : : t=0.3 output should happen at physics time units of dt=0.001. This must be
1180 [ + - ][ + - ]: 1348 : used within a relevant block.)", "vector of 3 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1181 : :
1182 : : keywords.insert({"refined", "Toggle refined field output on/off",
1183 : : R"(This keyword can be used to turn on/off refined field output, which
1184 : : refines the mesh and evaluates the solution on the refined mesh for saving
1185 [ + - ][ + - ]: 1348 : the solution.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
1186 : :
1187 : : keywords.insert({"filetype", "Select output file type",
1188 : : R"(This keyword is used to specify the output file type of
1189 [ + - ][ + - ]: 1348 : mesh-based field output in a field_output block.)", "string" });
[ + - ][ + - ]
[ + - ][ + - ]
1190 : :
1191 : : keywords.insert({"elemvar",
1192 : : "Specify list of elem-centered variables for output",
1193 : : R"(This keyword is used to specify elem-centered variables for output to
1194 [ + - ][ + - ]: 1348 : file. It is used in field_output blocks.)", "vector of string"});
[ + - ][ + - ]
[ + - ][ + - ]
1195 : :
1196 : : keywords.insert({"nodevar",
1197 : : "Specify list of node-centered variables for output",
1198 : : R"(This keyword is used to specify node-centered variables for output to
1199 [ + - ][ + - ]: 1348 : file. It is used in field_output blocks.)", "vector of string"});
[ + - ][ + - ]
[ + - ][ + - ]
1200 : :
1201 : : keywords.insert({"sideset",
1202 : : "Specify list of side sets",
1203 : : R"(This keyword is used to specify side sets. Usage is context specific,
1204 : : i.e. depends on what block it is specified in. Eg. in the field_output
1205 : : block it specifies the sidesets on which field output is desired.)",
1206 [ + - ][ + - ]: 1348 : "vector of uints"});
[ + - ][ + - ]
[ + - ][ + - ]
1207 : :
1208 : : keywords.insert({"diagnostics",
1209 : : "Specify the diagnostics block",
1210 : : R"(This keyword is used to introduce the dagnostics block, used to
1211 [ + - ][ + - ]: 1348 : configure diagnostics output.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1212 : :
1213 : : keywords.insert({"interval",
1214 : : "Set interval (in units of iteration count)",
1215 : : R"(This keyword is used to specify an interval in units of iteration count
1216 : : (i.e., number of time steps). Usage is context specific, and this must
1217 [ + - ][ + - ]: 1348 : be used within a relevant block.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1218 : :
1219 : : keywords.insert({"error", "Select an error type",
1220 : : R"(This keyword is used to select the error type. Used either in the
1221 : : diagnostics block to specify error norm, or in the AMR block to specify
1222 [ + - ][ + - ]: 1348 : the error for solution-adaptive mesh refinement.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1223 : :
1224 : : keywords.insert({"format",
1225 : : "Specify the ASCII floating-point output format",
1226 : : R"(This keyword is used to select the
1227 : : floating-point output format for ASCII floating-point number output.
1228 : : Valid options are 'default', 'fixed', and 'scientific'. For more info on
1229 : : these various formats, see
1230 [ + - ][ + - ]: 1348 : http://en.cppreference.com/w/cpp/io/manip/fixed.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1231 : :
1232 : : keywords.insert({"precision",
1233 : : "Precision in digits for ASCII floating-point output",
1234 : : R"(This keyword is used to select
1235 : : the precision in digits for ASCII floating-point real number output.
1236 : : Example: "precision=10", which selects ten digits for floating-point
1237 : : output, e.g., 3.141592654. The number of digits must be larger than zero
1238 : : and lower than the maximum representable digits for the given
1239 : : floating-point type. For more info on setting the precision in C++, see
1240 : : http://en.cppreference.com/w/cpp/io/manip/setprecision, and
1241 [ + - ][ + - ]: 1348 : http://en.cppreference.com/w/cpp/types/numeric_limits/digits10)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1242 : :
1243 : : keywords.insert({"history_output",
1244 : : "Start of history_output input block",
1245 : : R"(This keyword is used to start a block in the input file containing the
1246 [ + - ][ + - ]: 1348 : descriptions and settings of requested history output.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1247 : :
1248 : : keywords.insert({"point",
1249 : : "Start configuration block for history point, or a single point in IC blocks",
1250 : : R"(This keyword is used to either introduce a vector block used to
1251 : : specify probes for history output, or in the IC/BC block to specify
1252 : : a single point. When used in history output, it takes sub-entries of
1253 : : 'id' and 'coord'. When used in IC/BC, directly takes three reals as
1254 [ + - ][ + - ]: 1348 : coordinates.)", "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1255 : :
1256 : : keywords.insert({"coord", "Specify point coordinates",
1257 : : R"(This keyword is used to specify coordinates of the history-point.)",
1258 [ + - ][ + - ]: 1348 : "3 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1259 : :
1260 : : keywords.insert({"exodusii", "Select ExodusII output",
1261 : : R"(This keyword is used to select the
1262 : : ExodusII output file type readable by, e.g., ParaView of either a requested
1263 : : probability density function (PDF) within a pdfs ... end block or for
1264 : : mesh-based field output in a field_output ... end block. Example:
1265 : : "filetype exodusii", which selects ExodusII file output. For more info on
1266 [ + - ][ + - ]: 1348 : ExodusII, see http://sourceforge.net/projects/exodusii.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1267 : :
1268 : : keywords.insert({"density", "Request/specify density",
1269 : : R"(This keyword is used to request/specify the density. Usage is
1270 : : context specific. When specifed as 'elemvar' or 'nodevar' inside
1271 : : field_output, requests density as an output quantity. Otherwise,
1272 [ + - ][ + - ]: 1348 : specifies density at IC/BC.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1273 : :
1274 : : keywords.insert({"x-momentum",
1275 : : "Request x-momentum",
1276 : : R"(This keyword is used to request the fluid x-momentum as an output
1277 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1278 : :
1279 : : keywords.insert({"y-momentum",
1280 : : "Request y-momentum",
1281 : : R"(This keyword is used to request the fluid y-momentum as an output
1282 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1283 : :
1284 : : keywords.insert({"z-momentum",
1285 : : "Request z-momentum",
1286 : : R"(This keyword is used to request the fluid z-momentum as an output
1287 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1288 : :
1289 : : keywords.insert({
1290 : : "specific_total_energy", "Request specific total energy",
1291 : : R"(This keyword is used to request the specific total energy as an output
1292 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1293 : :
1294 : : keywords.insert({
1295 : : "volumetric_total_energy", "Request total volumetric energy",
1296 : : R"(This keyword is used to request the volumetric total energy as an output
1297 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1298 : :
1299 : : keywords.insert({"x-velocity",
1300 : : "Request x-velocity",
1301 : : R"(This keyword is used to request the fluid x-velocity as an output
1302 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1303 : :
1304 : : keywords.insert({"y-velocity",
1305 : : "Request y-velocity",
1306 : : R"(This keyword is used to request the fluid y-velocity as an output
1307 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1308 : :
1309 : : keywords.insert({"z-velocity",
1310 : : "Request z-velocity",
1311 : : R"(This keyword is used to request the fluid z-velocity as an output
1312 [ + - ][ + - ]: 1348 : variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1313 : :
1314 : : keywords.insert({"pressure", "Request/specify pressure",
1315 : : R"(This keyword is used to request/specify the pressure. Usage is
1316 : : context specific. When specifed as 'elemvar' or 'nodevar' inside
1317 : : field_output, requests pressure as an output quantity. Otherwise,
1318 [ + - ][ + - ]: 1348 : specifies pressure at IC/BC.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1319 : :
1320 : : keywords.insert({"material_indicator",
1321 : : "Request material_indicator",
1322 : : R"(This keyword is used to request the material indicator function as an
1323 [ + - ][ + - ]: 1348 : output variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1324 : :
1325 : : keywords.insert({"analytic",
1326 : : "Request analytic solution",
1327 : : R"(This keyword is used to request the analytic solution (if exist) as an
1328 [ + - ][ + - ]: 1348 : output variable.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1329 : :
1330 : : keywords.insert({"l2", "Select the L2 norm",
1331 [ + - ][ + - ]: 1348 : R"(This keyword is used to enable computing the L2 norm.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1332 : :
1333 : : keywords.insert({"linf", "Select the L_{infinity} norm",
1334 [ + - ][ + - ]: 1348 : R"(This keyword is used to enable computing the L-infinity norm.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1335 : :
1336 : : keywords.insert({"default",
1337 : : "Select the default ASCII floating-point output",
1338 : : R"(This keyword is used to select the
1339 : : 'default' floating-point output format for ASCII floating-point real
1340 : : number output. For more info on these various formats, see
1341 [ + - ][ + - ]: 1348 : http://en.cppreference.com/w/cpp/io/manip/fixed.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1342 : :
1343 : : keywords.insert({"fixed",
1344 : : "Select the fixed ASCII floating-point output",
1345 : : R"(This keyword is used to select the
1346 : : 'fixed' floating-point output format for ASCII floating-point real
1347 : : number output. For more info on these various formats, see
1348 [ + - ][ + - ]: 1348 : http://en.cppreference.com/w/cpp/io/manip/fixed.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1349 : :
1350 : : keywords.insert({"scientific",
1351 : : "Select the scientific ASCII floating-point output",
1352 : : R"(This keyword is used to select the
1353 : : 'scientific' floating-point output format for ASCII floating-point real
1354 : : number output. For more info on these various formats, see
1355 [ + - ][ + - ]: 1348 : http://en.cppreference.com/w/cpp/io/manip/fixed.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1356 : :
1357 : : // -----------------------------------------------------------------------
1358 : : // ALE options
1359 : : // -----------------------------------------------------------------------
1360 : :
1361 : : keywords.insert({"ale", "Start block configuring ALE",
1362 : : R"(This keyword is used to introduce the ale block, used to
1363 : : configure arbitrary Lagrangian-Eulerian (ALE) mesh movement.)",
1364 [ + - ][ + - ]: 1348 : "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1365 : :
1366 : : keywords.insert({"smoother", "Select mesh velocity smoother",
1367 : : R"(This keyword is used to select a mesh velocity smoother option, used
1368 : : for Arbitrary-Lagrangian-Eulerian (ALE) mesh motion. Valid options are
1369 [ + - ][ + - ]: 1348 : 'laplace', 'helmholtz', and 'none')", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1370 : :
1371 : : keywords.insert({"mesh_velocity", "Select mesh velocity",
1372 : : R"(This keyword is used to select a mesh velocity option, used for
1373 : : Arbitrary-Lagrangian-Eulerian (ALE) mesh motion. Valid options are
1374 [ + - ][ + - ]: 1348 : 'sine', 'fluid', and 'user_defined".)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1375 : :
1376 : : keywords.insert({"mesh_motion",
1377 : : "List of dimension indices that are allowed to move in ALE calculations",
1378 : : R"(This keyword is used to specify a list of integers (0, 1, or 2) whose
1379 : : coordinate directions corresponding to x, y, or z are allowed to move with
1380 : : the mesh velocity in ALE calculations. Useful for 1D/2D problems.)",
1381 [ + - ][ + - ]: 1348 : "vector of uints"});
[ + - ][ + - ]
[ + - ][ + - ]
1382 : :
1383 : : keywords.insert({"meshforce", "Set ALE mesh force model parameter(s)",
1384 : : R"(This keyword is used to specify a vector of real numbers used to
1385 : : parameterize a mesh force model for ALE. The length of the vector must
1386 [ + - ][ + - ]: 1348 : exactly be 4.)", "vector of 4 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1387 : :
1388 : : keywords.insert({"dvcfl",
1389 : : "Set the volume-change Courant-Friedrichs-Lewy (CFL) coefficient",
1390 : : R"(This keyword is used to specify the volume-change (dV/dt) CFL coefficient
1391 : : for variable-time-step-size simulations due to volume change in time in
1392 : : arbitrary-Lagrangian-Eulerian (ALE) calculations. Setting 'dvcfl' only has
1393 : : effect in ALE calculations and used together with 'cfl'. See also J. Waltz,
1394 : : N.R. Morgan, T.R. Canfield, M.R.J. Charest, L.D. Risinger, J.G. Wohlbier, A
1395 : : three-dimensional finite element arbitrary Lagrangian–Eulerian method for
1396 : : shock hydrodynamics on unstructured grids, Computers & Fluids, 92: 172-187,
1397 [ + - ][ + - ]: 1348 : 2014.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1398 : :
1399 : : keywords.insert({"vortmult",
1400 : : "Configure vorticity multiplier for ALE mesh velocity",
1401 : : R"(This keyword is used to configure the multiplier for the vorticity term
1402 : : in the mesh velocity smoother (mesh_velocity=fluid) or for the potential
1403 : : gradient for the Helmholtz mesh velocity (mesh_velocity=helmholtz) for ALE
1404 : : mesh motion. For 'fluid' this is coefficient c2 in Eq.(36) of Waltz,
1405 : : Morgan, Canfield, Charest, Risinger, Wohlbier, A three-dimensional finite
1406 : : element arbitrary Lagrangian–Eulerian method for shock hydrodynamics on
1407 : : unstructured grids, Computers & Fluids, 2014, and for 'helmholtz', this
1408 : : is coefficient a1 in Eq.(23) of Bakosi, Waltz, Morgan, Improved ALE mesh
1409 : : velocities for complex flows, International Journal for Numerical Methods
1410 [ + - ][ + - ]: 1348 : in Fluids, 2017. )", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1411 : :
1412 : : keywords.insert({"maxit",
1413 : : "Set the max number of iterations for the ALE mesh velocity linear solve",
1414 : : R"(This keyword is used to specify the maximum number of linear solver
1415 : : iterations taken to converge the mesh velocity linear solve in
1416 : : arbitrary-Lagrangian-Eulerian (ALE) calculations. See also J. Waltz,
1417 : : N.R. Morgan, T.R. Canfield, M.R.J. Charest, L.D. Risinger, J.G. Wohlbier, A
1418 : : three-dimensional finite element arbitrary Lagrangian–Eulerian method for
1419 : : shock hydrodynamics on unstructured grids, Computers & Fluids, 92: 172-187,
1420 [ + - ][ + - ]: 1348 : 2014.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1421 : :
1422 : : keywords.insert({"tolerance",
1423 : : "Set the tolerance for the ALE mesh velocity linear solve",
1424 : : R"(This keyword is used to specify the tolerance to converge the mesh
1425 : : velocity linear solve for in
1426 : : arbitrary-Lagrangian-Eulerian (ALE) calculations. See also J. Waltz,
1427 : : N.R. Morgan, T.R. Canfield, M.R.J. Charest, L.D. Risinger, J.G. Wohlbier, A
1428 : : three-dimensional finite element arbitrary Lagrangian–Eulerian method for
1429 : : shock hydrodynamics on unstructured grids, Computers & Fluids, 92: 172-187,
1430 [ + - ][ + - ]: 1348 : 2014.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1431 : :
1432 : : keywords.insert({"move",
1433 : : "Start configuration block configuring surface movement in ALE",
1434 : : R"(This keyword is used to introduce a move block, used to
1435 [ + - ][ + - ]: 1348 : configure surface movement for ALE simulations.)", "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1436 : :
1437 : : keywords.insert({"fntype",
1438 : : "Select how a user-defined function is interpreted",
1439 : : R"(This keyword is used to select how a user-defined function should be
1440 [ + - ][ + - ]: 1348 : interpreted.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1441 : :
1442 : : keywords.insert({"fn", "Specify a discrete user-defined function",
1443 : : R"(This keyword is used to specify a user-defined function block with
1444 : : discrete points, listed inside the fn block. Used in ale mesh motion and
1445 [ + - ][ + - ]: 1348 : time-dependent BC specification)", "reals" });
[ + - ][ + - ]
[ + - ][ + - ]
1446 : :
1447 : : keywords.insert({"laplace",
1448 : : "Select the Laplace mesh velocity smoother for ALE",
1449 : : R"(This keyword is used to select the 'Laplace' mesh velocity smoother
1450 [ + - ][ + - ]: 1348 : for Arbitrary-Lagrangian-Eulerian (ALE) mesh motion.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1451 : :
1452 : : keywords.insert({"helmholtz",
1453 : : "Select the Helmholtz velocity for ALE",
1454 : : R"(This keyword is used to select the a velocity, computed from the
1455 : : Helmholtz-decomposition as the mesh velocity for
1456 : : Arbitrary-Lagrangian-Eulerian (ALE) mesh motion. See J. Bakosi, J. Waltz,
1457 : : N. Morgan, Improved ALE mesh velocities for complex flows, Int. J. Numer.
1458 [ + - ][ + - ]: 1348 : Meth. Fl., 1-10, 2017, https://doi.org/10.1002/fld.4403.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1459 : :
1460 : : keywords.insert({"none", "Select none option",
1461 : : R"(This keyword is used to select the 'none' option from a list of
1462 [ + - ][ + - ]: 1348 : configuration options.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1463 : :
1464 : : keywords.insert({"sine",
1465 : : "Prescribe sinusoidal mesh velocity for ALE",
1466 : : R"(This keyword is used to prescribe a sinusoidal mesh velocity
1467 [ + - ][ + - ]: 1348 : for Arbitrary-Lagrangian-Eulerian (ALE) mesh motion.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1468 : :
1469 : : keywords.insert({"fluid", "Select the fluid velocity for ALE",
1470 : : R"(This keyword is used to select the 'fluid' velocity as the mesh velocity
1471 [ + - ][ + - ]: 1348 : for Arbitrary-Lagrangian-Eulerian (ALE) mesh motion.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1472 : :
1473 : : // -----------------------------------------------------------------------
1474 : : // h/p adaptation objects
1475 : : // -----------------------------------------------------------------------
1476 : :
1477 : : keywords.insert({"amr",
1478 : : "Start configuration block configuring adaptive mesh refinement",
1479 : : R"(This keyword is used to introduce the amr block, used to
1480 [ + - ][ + - ]: 1348 : configure adaptive mesh refinement.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1481 : :
1482 : : keywords.insert({"t0ref", "Enable mesh refinement at t<0",
1483 : : R"(This keyword is used to enable initial mesh refinement, which can be
1484 : : configured to perform multiple levels of mesh refinement based on various
1485 [ + - ][ + - ]: 1348 : refinement criteria and configuration settings.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
1486 : :
1487 : : keywords.insert({"dtref", "Enable mesh refinement at t>0",
1488 : : R"(This keyword is used to enable solution-adaptive mesh refinement during
1489 [ + - ][ + - ]: 1348 : time stepping.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
1490 : :
1491 : : keywords.insert({"dtref_uniform",
1492 : : "Enable mesh refinement at t>0 but only perform uniform refinement",
1493 : : R"(This keyword is used to force uniform-only solution-adaptive mesh
1494 [ + - ][ + - ]: 1348 : refinement during time stepping.)", "bool"});
[ + - ][ + - ]
[ + - ][ + - ]
1495 : :
1496 : : keywords.insert({"dtfreq",
1497 : : "Set mesh refinement frequency during time stepping",
1498 : : R"(This keyword is used to configure the frequency of mesh refinement
1499 [ + - ][ + - ]: 1348 : during time stepping.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1500 : :
1501 : : keywords.insert({"maxlevels",
1502 : : "Set maximum allowed mesh refinement levels",
1503 : : R"(This keyword is used to configure the maximum allowed mesh refinement
1504 [ + - ][ + - ]: 1348 : levels.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1505 : :
1506 : : keywords.insert({"initial",
1507 : : "Configure initial mesh refinement (before time stepping)",
1508 : : R"(This keyword is used to add to a list of initial mesh refinement types
1509 : : that happens before t = 0. Allowed options are 'uniform',
1510 : : 'uniform_derefine', 'initial_conditions', 'coords', 'edgelist')",
1511 [ + - ][ + - ]: 1348 : "vector of strings"});
[ + - ][ + - ]
[ + - ][ + - ]
1512 : :
1513 : : keywords.insert({"coords",
1514 : : "Configure initial refinement using coordinate planes",
1515 : : R"(This keyword can be used to configure entire volumes on a given side of
1516 : : a plane in 3D space. The keyword introduces an coords block within
1517 : : an amr block. All edges of the input mesh will be tagged for refinement
1518 : : whose end-points lie within the given ranges.
1519 : : Example: 'xminus 0.5' refines all edges whose end-point coordinates are
1520 : : less than 0.5. Multiple specifications are understood by combining with
1521 : : a logical AND. That is: 'xminus 0.5 yplus 0.3' refines all edges whose
1522 : : end-point x coordinates are less than 0.5 AND y coordinates are larger than
1523 [ + - ][ + - ]: 1348 : 0.3.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1524 : :
1525 : : keywords.insert({"xminus",
1526 : : "Configure initial refinement for coordinates lower than an x-normal plane",
1527 : : R"(This keyword can be used to configure a mesh refinement volume for edges
1528 : : whose end-points are less than the x coordinate of a plane perpendicular
1529 : : to coordinate x in 3D space. The keyword must be used in a coords-block
1530 : : within an amr-block with syntax 'xminus <real>'. All edges of the
1531 : : input mesh will be tagged for refinement whose end-points lie less than (-)
1532 : : the real number given. Example: 'xminus 0.5' refines all edges whose end-point
1533 [ + - ][ + - ]: 1348 : x-coordinates are less than 0.5.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1534 : :
1535 : : keywords.insert({"xplus",
1536 : : "Configure initial refinement for coordinates larger than an x-normal plane",
1537 : : R"(This keyword can be used to configure a mesh refinement volume for edges
1538 : : whose end-points are larger than the x coordinate of a plane perpendicular
1539 : : to coordinate x in 3D space. The keyword must be used in a coords-block
1540 : : within an amr-block with syntax 'xplus <real>'. All edges of the
1541 : : input mesh will be tagged for refinement whose end-points lie larger than
1542 : : (+) the real number given. Example: 'xplus 0.5' refines all edges whose
1543 [ + - ][ + - ]: 1348 : end-point coordinates are larger than 0.5.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1544 : :
1545 : : keywords.insert({"yminus",
1546 : : "Configure initial refinement for coordinates lower than an y-normal plane",
1547 : : R"(This keyword can be used to configure a mesh refinement volume for edges
1548 : : whose end-points are less than the y coordinate of a plane perpendicular
1549 : : to coordinate y in 3D space. The keyword must be used in a coords-block
1550 : : within an amr-block with syntax 'yminus <real>'. All edges of the
1551 : : input mesh will be tagged for refinement whose end-points lie less than (-)
1552 : : the real number given. Example: 'yminus 0.5' refines all edges whose end-point
1553 [ + - ][ + - ]: 1348 : coordinates are less than 0.5.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1554 : :
1555 : : keywords.insert({"yplus",
1556 : : "Configure initial refinement for coordinates larger than an y-normal plane",
1557 : : R"(This keyword can be used to configure a mesh refinement volume for edges
1558 : : whose end-points are larger than the y coordinate of a plane perpendicular
1559 : : to coordinate y in 3D space. The keyword must be used in a coords-block
1560 : : within an amr-block with syntax 'yplus <real>'. All edges of the
1561 : : input mesh will be tagged for refinement whose end-points lie larger than
1562 : : (+) the real number given. Example: 'yplus 0.5' refines all edges whose
1563 [ + - ][ + - ]: 1348 : end-point coordinates are larger than 0.5.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1564 : :
1565 : : keywords.insert({"zminus",
1566 : : "Configure initial refinement for coordinates lower than an z-normal plane",
1567 : : R"(This keyword can be used to configure a mesh refinement volume for edges
1568 : : whose end-points are less than the z coordinate of a plane perpendicular
1569 : : to coordinate z in 3D space. The keyword must be used in a coords-block
1570 : : within an amr-block with syntax 'zminus <real>'. All edges of the
1571 : : input mesh will be tagged for refinement whose end-points lie less than (-)
1572 : : the real number given. Example: 'zminus 0.5' refines all edges whose end-point
1573 [ + - ][ + - ]: 1348 : coordinates are less than 0.5.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1574 : :
1575 : : keywords.insert({"zplus",
1576 : : "Configure initial refinement for coordinates larger than an z-normal plane",
1577 : : R"(This keyword can be used to configure a mesh refinement volume for edges
1578 : : whose end-points are larger than the z coordinate of a plane perpendicular
1579 : : to coordinate z in 3D space. The keyword must be used in a coords-block
1580 : : within an amr-block with syntax 'zplus <real>'. All edges of the
1581 : : input mesh will be tagged for refinement whose end-points lie larger than
1582 : : (+) the real number given. Example: 'zplus 0.5' refines all edges whose
1583 [ + - ][ + - ]: 1348 : end-point coordinates are larger than 0.5.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1584 : :
1585 : : keywords.insert({"edgelist",
1586 : : "Configure edge-node pairs for initial refinement",
1587 : : R"(This keyword can be used to configure a list of edges that are explicitly
1588 : : tagged for initial refinement during setup in inciter. The keyword
1589 : : introduces an edgelist block within an amr block and must
1590 : : contain a list of integer pairs, i.e., the number of ids must be even,
1591 : : denoting the end-points of the nodes (=edge) which should be tagged for
1592 [ + - ][ + - ]: 1348 : refinement.)", "vector of uints"});
[ + - ][ + - ]
[ + - ][ + - ]
1593 : :
1594 : : keywords.insert({"error",
1595 : : "Configure the error type for solution-adaptive mesh refinement",
1596 : : R"(This keyword is used to select the algorithm used to estimate the error
1597 : : for solution-adaptive mesh refinement. Available options are 'jump' and
1598 [ + - ][ + - ]: 1348 : 'hessian')", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1599 : :
1600 : : keywords.insert({"refvar",
1601 : : "Configure dependent variables used for adaptive mesh refinement",
1602 : : R"(This keyword is used to configured a list of dependent variables that
1603 : : trigger adaptive mesh refinement based on estimating their numerical error.
1604 : : These refinement variables are used for both initial (i.e., before time
1605 : : stepping) mesh refinement as well as during time stepping. Only previously
1606 : : (i.e., earlier in the input file) selected dependent variables can be
1607 : : configured as refinement variables. Dependent variables are required to be
1608 : : defined in all equation system configuration blocks, e.g., transport ...
1609 : : end, by using the 'depvar' keyword. Example: transport depvar c end amr
1610 : : refvar c end end. Selecting a particular scalar component in a system is
1611 : : done by appending the equation number to the refvar: Example: transport
1612 : : depvar q ncomp 3 end amr refvar q1 q2 end end, which configures two
1613 : : refinement variables: the first and third scalar component of the previously
1614 [ + - ][ + - ]: 1348 : configured transport equation system.)", "vector of char"});
[ + - ][ + - ]
[ + - ][ + - ]
1615 : :
1616 : : keywords.insert({"tol_refine", "Configure refine tolerance",
1617 : : R"(This keyword is used to set the tolerance used to tag an edge for
1618 [ + - ][ + - ]: 1348 : refinement if the relative error exceeds this value.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1619 : :
1620 : : keywords.insert({"tol_derefine",
1621 : : "Configure derefine tolerance",
1622 : : R"(This keyword is used to set the tolerance used to tag an edge for
1623 [ + - ][ + - ]: 1348 : derefinement if the relative error is below this value.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1624 : :
1625 : : keywords.insert({"uniform",
1626 : : "Select uniform initial mesh refinement",
1627 : : R"(This keyword is used to select uniform initial mesh refinement.)",
1628 [ + - ][ + - ]: 1348 : "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1629 : :
1630 : : keywords.insert({"uniform_derefine",
1631 : : "Select uniform initial mesh de-refinement",
1632 : : R"(This keyword is used to select uniform initial mesh de-refinement.)",
1633 [ + - ][ + - ]: 1348 : "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1634 : :
1635 : : keywords.insert({"initial_conditions",
1636 : : "Select initial-conditions-based initial mesh refinement",
1637 : : R"(This keyword is used to select initial-conditions-based initial mesh
1638 [ + - ][ + - ]: 1348 : refinement.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1639 : :
1640 : : keywords.insert({"jump",
1641 : : "Error estimation based on the solution jump normalized by solution value",
1642 : : R"(This keyword is used to select the jump-based error indicator for
1643 : : solution-adaptive mesh refinement. The error is estimated by computing the
1644 : : magnitude of the jump in the solution value normalized by the solution
1645 [ + - ][ + - ]: 1348 : value.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1646 : :
1647 : : keywords.insert({"hessian",
1648 : : "Error estimation based on the Hessian normalized by solution value",
1649 : : R"(This keyword is used to select the Hessian-based error indicator for
1650 : : solution-adaptive mesh refinement. The error is estimated by computing the
1651 : : Hessian (2nd derivative matrix) of the solution normalized by sum of the
1652 [ + - ][ + - ]: 1348 : absolute values of the gradients at edges-end points.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1653 : :
1654 : : keywords.insert({"pref",
1655 : : "Start configuration block configuring p-adaptive refinement",
1656 : : R"(This keyword is used to introduce the pref block, to
1657 [ + - ][ + - ]: 1348 : configure p-adaptive refinement)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1658 : :
1659 : : keywords.insert({"indicator",
1660 : : "Configure the specific adaptive indicator for p-adaptive DG scheme",
1661 : : R"(This keyword can be used to configure a specific type of adaptive
1662 : : indicator for p-adaptive refinement of the DG scheme. The keyword must
1663 : : be used in a pref block. Available options are 'pref_spectral_decay' and
1664 [ + - ][ + - ]: 1348 : 'pref_non_conformity'.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1665 : :
1666 : : keywords.insert({"ndofmax",
1667 : : "Configure the maximum number of degree of freedom for p-adaptive DG",
1668 : : R"(This keyword can be used to configure a maximum number of degree of
1669 : : freedom for p-adaptive refinement of the DG scheme. The keyword must
1670 [ + - ][ + - ]: 1348 : be used in a pref block.)", "uint either 4 or 10"});
[ + - ][ + - ]
[ + - ][ + - ]
1671 : :
1672 : : keywords.insert({"tolref",
1673 : : "Configure the tolerance for p-refinement for p-adaptive DG",
1674 : : R"(This keyword can be used to configure a tolerance for p-adaptive
1675 : : refinement for the DG scheme. The keyword must be used in a pref
1676 : : block. All elements with a refinement indicator larger than this
1677 [ + - ][ + - ]: 1348 : tolerance will be p-refined.)", "real between 0 and 1"});
[ + - ][ + - ]
[ + - ][ + - ]
1678 : :
1679 : : keywords.insert({"spectral_decay",
1680 : : "Select the spectral-decay indicator for p-adaptive DG scheme",
1681 : : R"(This keyword is used to select the spectral-decay indicator used for
1682 : : p-adaptive discontinuous Galerkin (DG) discretization used in inciter.
1683 : : See Control/Inciter/Options/PrefIndicator.hpp for other valid options.)",
1684 [ + - ][ + - ]: 1348 : "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1685 : :
1686 : : keywords.insert({"non_conformity",
1687 : : "Select the non-conformity indicator for p-adaptive DG scheme",
1688 : : R"(This keyword is used to select the non-conformity indicator used for
1689 : : p-adaptive discontinuous Galerkin (DG) discretization used in inciter.
1690 : : See Control/Inciter/Options/PrefIndicator.hpp for other valid options.)",
1691 [ + - ][ + - ]: 1348 : "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1692 : :
1693 : : // -----------------------------------------------------------------------
1694 : : // boundary condition options
1695 : : // -----------------------------------------------------------------------
1696 : :
1697 : : keywords.insert({"bc",
1698 : : "Start configuration block for boundary conditions",
1699 : : R"(This keyword is used to introduce the bc block, used for
1700 : : boundary conditions. This is a vector block, where each vector entry
1701 [ + - ][ + - ]: 1348 : specifies BCs for a particular mesh)", "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1702 : :
1703 : : keywords.insert({"mesh",
1704 : : "List meshes on which the following BCs apply",
1705 : : R"(This keyword is used to list multiple meshes on which the boundary
1706 : : conditions listed in this particular bc-block apply.)",
1707 [ + - ][ + - ]: 1348 : "vector of uints"});
[ + - ][ + - ]
[ + - ][ + - ]
1708 : :
1709 : : keywords.insert({"dirichlet",
1710 : : "List sidesets with Dirichlet boundary conditions",
1711 : : R"(This keyword is used to list Dirichlet sidesets.
1712 : : This keyword is used to list multiple sidesets on
1713 : : which a prescribed Dirichlet BC is then applied. Such prescribed BCs
1714 : : at each point in space and time are evaluated using a built-in function,
1715 [ + - ][ + - ]: 1348 : e.g., using the method of manufactured solutions.)", "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1716 : :
1717 : : keywords.insert({"symmetry",
1718 : : "List sidesets with symmetry boundary conditions",
1719 : : R"(This keyword is used to list (multiple) symmetry BC sidesets.)",
1720 [ + - ][ + - ]: 1348 : "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1721 : :
1722 : : keywords.insert({"inlet",
1723 : : "List sidesets with inlet boundary conditions",
1724 : : R"(This keyword is used to list (multiple) inlet BC sidesets.)",
1725 [ + - ][ + - ]: 1348 : "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1726 : :
1727 : : keywords.insert({"outlet",
1728 : : "List sidesets with outlet boundary conditions",
1729 : : R"(This keyword is used to list (multiple) outlet BC sidesets.)",
1730 [ + - ][ + - ]: 1348 : "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1731 : :
1732 : : keywords.insert({"farfield",
1733 : : "List sidesets with farfield boundary conditions",
1734 : : R"(This keyword is used to list (multiple) farfield BC sidesets.
1735 : : Keywords allowed in a bc_farfield block are 'density', 'velocity',
1736 [ + - ][ + - ]: 1348 : 'pressure')", "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1737 : :
1738 : : keywords.insert({"extrapolate",
1739 : : "List sidesets with Extrapolation boundary conditions",
1740 : : R"(This keyword is used to list (multiple) extrapolate BC sidesets.)",
1741 [ + - ][ + - ]: 1348 : "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1742 : :
1743 : : keywords.insert({"noslipwall",
1744 : : "List sidesets with no-slip wall boundary conditions",
1745 : : R"(This keyword is used to list (multiple) no-slip wall BC sidesets.)",
1746 [ + - ][ + - ]: 1348 : "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1747 : :
1748 : : keywords.insert({"slipwall",
1749 : : "List sidesets with slip wall boundary conditions",
1750 : : R"(This keyword is used to list (multiple) slip wall BC sidesets.)",
1751 [ + - ][ + - ]: 1348 : "vector of uint(s)"});
[ + - ][ + - ]
[ + - ][ + - ]
1752 : :
1753 : : keywords.insert({"timedep",
1754 : : "Start configuration block describing time dependent boundary conditions",
1755 : : R"(This keyword is used to introduce a bc_timedep block, used to
1756 : : specify the configuration of time dependent boundary conditions for a
1757 : : partial differential equation. A discrete function in time t in the form
1758 : : of a table with 6 columns (t, pressure(t), density(t), vx(t), vy(t), vz(t))
1759 : : is expected inside a fn ... end block, specified within the bc_timedep
1760 : : block. Multiple such bc_timedep blocks can be specified for different
1761 [ + - ][ + - ]: 1348 : time dependent BCs on different groups of side sets.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1762 : :
1763 : : keywords.insert({"back_pressure",
1764 : : "Start configuration block describing back pressure boundary conditions",
1765 : : R"(This keyword is used to introduce a back pressure BC block. This
1766 : : block requires a 'sideset' vector and 'pressure' to be specified within
1767 [ + - ][ + - ]: 1348 : it.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1768 : :
1769 : : keywords.insert({"velocity", "Specify velocity",
1770 : : R"(This keyword is used to configure a velocity vector used in a
1771 : : context-specific way, e.g., for boundary or initial conditions.)",
1772 [ + - ][ + - ]: 1348 : "vector of 3 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1773 : :
1774 : : // -----------------------------------------------------------------------
1775 : : // Rigid-body motion solver
1776 : : // -----------------------------------------------------------------------
1777 : :
1778 : : keywords.insert({"rigid_body_motion", "Specify a rigid body motion block",
1779 : : R"(This keyword is used to specify a rigid body motion block, to move an
1780 : : overset mesh as a rigid body. Number of degrees of freedom and the
1781 : : symmetry plane (if any) are specified within this block.)",
1782 [ + - ][ + - ]: 1348 : "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1783 : :
1784 : : keywords.insert({"rigid_body_dof",
1785 : : "Number of rigid body degrees of freedom", R"(This keyword is used to
1786 : : specify the number of degrees of freedom the rigid body has. Valid
1787 : : options are 3 and 6. 3 DOFs indicate translation in two dimensions and
1788 : : rotation about symmetry plane axis; 6 DOFs indication translation in
1789 : : three dimensions and rotation about three axes. If 3 DOFs is specified,
1790 : : symmetry plane is required; if 6 DOFs is specified symmetry plane is not
1791 [ + - ][ + - ]: 1348 : used.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1792 : :
1793 : : keywords.insert({"symmetry_plane", "Symmetry plane for rigid body motion",
1794 : : R"(This keyword is used to specify the symmetry plane for a 3 DOF rigid
1795 [ + - ][ + - ]: 1348 : body motion solver. 1: x-plane, 2: y-plane, 3: z-plane.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1796 : :
1797 : : // -----------------------------------------------------------------------
1798 : : // IC object
1799 : : // -----------------------------------------------------------------------
1800 : :
1801 : : keywords.insert({"ic",
1802 : : "Introduce an ic block used to configure initial conditions",
1803 : : R"(This keyword is used to introduce an ic block used to set initial
1804 [ + - ][ + - ]: 1348 : conditions.)", "block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1805 : :
1806 : : keywords.insert({"materialid", "Specify material id",
1807 : : R"(This keyword is used to configure the material id within an IC box,
1808 : : IC mesh-block, farfield BC, or in the background as a part of the
1809 [ + - ][ + - ]: 1348 : initialization.)", "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1810 : :
1811 : : keywords.insert({"temperature", "Specify temperature",
1812 : : R"(This keyword is used to configure temperature, used for, e.g.,
1813 [ + - ][ + - ]: 1348 : boundary or initial conditions.)" , "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1814 : :
1815 : : keywords.insert({"mass_fractions", "Specify species mass fractions",
1816 : : R"(This keyword is used to configure species mass fractions, used for,
1817 [ + - ][ + - ]: 1348 : e.g., boundary or initial conditions.)" , "vector of reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1818 : :
1819 : : keywords.insert({"box",
1820 : : "Introduce a box block used to assign initial conditions",
1821 : : R"(This keyword is used to introduce a IC box block used to assign
1822 : : initial conditions within a box given by spatial coordinates.)",
1823 [ + - ][ + - ]: 1348 : "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1824 : :
1825 : : keywords.insert({"meshblock",
1826 : : "Introduce a meshblock block used to assign initial conditions",
1827 : : R"(This keyword is used to introduce a IC meshblock block used to
1828 : : assign initial conditions within a mesh block specified in the mesh file.)",
1829 [ + - ][ + - ]: 1348 : "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1830 : :
1831 : : keywords.insert({"blockid", "Specify mesh block id",
1832 : : R"(This keyword is used to configure the mesh block id within the
1833 : : meshblock-block as a part of the initialization. It is strongly
1834 : : recommended to use contiguous block ids in mesh file starting from 1.)",
1835 [ + - ][ + - ]: 1348 : "uint"});
[ + - ][ + - ]
[ + - ][ + - ]
1836 : :
1837 : : keywords.insert({"volume", "Specify volume",
1838 : : R"(This keyword is used to configure the volume of a meshblock.)",
1839 [ + - ][ + - ]: 1348 : "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1840 : :
1841 : : keywords.insert({"mass", "Specify mass",
1842 : : R"(This keyword is used to configure the mass within a box/meshblock,
1843 : : or mass of the rigid body which is conformally meshed using the overset
1844 [ + - ][ + - ]: 1348 : mesh.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1845 : :
1846 : : keywords.insert({"energy", "Specify energy per unit mass",
1847 : : R"(This keyword is used to configure energy per unit mass, used for, e.g.,
1848 [ + - ][ + - ]: 1348 : boundary or initial conditions.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1849 : :
1850 : : keywords.insert({"energy_content", "Specify energy per unit volume",
1851 : : R"(This keyword is used to configure energy per unit volume, used for
1852 [ + - ][ + - ]: 1348 : initial conditions.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1853 : :
1854 : : keywords.insert({"xmin", "Minimum x coordinate",
1855 : : R"(This keyword used to configure a minimum x coordinate to specify
1856 [ + - ][ + - ]: 1348 : a box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1857 : :
1858 : : keywords.insert({"xmax", "Maximum x coordinate",
1859 : : R"(This keyword used to configure a maximum x coordinate to specify
1860 [ + - ][ + - ]: 1348 : a box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1861 : :
1862 : : keywords.insert({"ymin", "Minimum y coordinate",
1863 : : R"(This keyword used to configure a minimum y coordinate to specify
1864 [ + - ][ + - ]: 1348 : a box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1865 : :
1866 : : keywords.insert({"ymax", "Maximum y coordinate",
1867 : : R"(This keyword used to configure a maximum y coordinate to specify
1868 [ + - ][ + - ]: 1348 : a box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1869 : :
1870 : : keywords.insert({"zmin", "Minimum z coordinate",
1871 : : R"(This keyword used to configure a minimum z coordinate to specify
1872 [ + - ][ + - ]: 1348 : a box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1873 : :
1874 : : keywords.insert({"zmax", "Maximum z coordinate",
1875 : : R"(This keyword used to configure a maximum z coordinate to specify
1876 [ + - ][ + - ]: 1348 : a box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1877 : :
1878 : : keywords.insert({"orientation", "Configure orientation",
1879 : : R"(Configure orientation of an IC box for rotation about centroid of
1880 : : box; or configure orientation of a mesh relative to another.)",
1881 [ + - ][ + - ]: 1348 : "vector of 3 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1882 : :
1883 : : keywords.insert({"initiate", "Initiation type",
1884 : : R"(This keyword is used to select an initiation type to configure how
1885 : : values are assigned for a box/meshblock initialization. This can be used
1886 : : to specify, how the values are assigned to mesh nodes within a box. Uses:
1887 : : (1) impulse: assign the full values at t=0 for all points in a box,
1888 : : (2) linear: use a linear function in time and space, configured with an
1889 : : initiation point in space, a constant velocity of the growing spherical
1890 : : front in time (and space) linearly, and width of the front and assigns
1891 : : values to mesh points falling within the growing spherical shell inside
1892 [ + - ][ + - ]: 1348 : a configured box.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1893 : :
1894 : : keywords.insert({"init_time","Specify the initialization time",
1895 : : R"(This keyword is used to specify the time at which the propagating front
1896 : : is initialized for a mesh block or box IC, with 'initiate linear' type.
1897 : : Delays in initializing separate mesh blocks or boxes can be achieved using
1898 [ + - ][ + - ]: 1348 : different initialization times.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1899 : :
1900 : : keywords.insert({"front_width", "Specify a front width",
1901 : : R"(This keyword is used to specify the width of the propagating front for
1902 : : a mesh block or box IC, with 'initiate linear' type. The suggested value
1903 : : of the front width is about 4-5 times the mesh size inside the mesh block
1904 [ + - ][ + - ]: 1348 : or box.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1905 : :
1906 : : keywords.insert({"front_speed", "Specify a front speed",
1907 : : R"(This keyword is used to specify the speed at which a front propagates
1908 [ + - ][ + - ]: 1348 : for a mesh block or box IC, with 'initiate linear' type.)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1909 : :
1910 : : keywords.insert({"impulse",
1911 : : "Select the impulse initiation type, for a box/meshblock IC",
1912 : : R"(This keyword can be used to select the 'impulse' initiation/assignment
1913 : : type for box initial conditions. It simply assigns the prescribed values
1914 [ + - ][ + - ]: 1348 : to all mesh points within a configured box at t=0.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1915 : :
1916 : : keywords.insert({"linear",
1917 : : "Select the linear initiation type, for a box/meshblock IC",
1918 : : R"(This keyword is be used to specify the 'linear' initiation parameters
1919 : : for a particular box or meshblock, as a part of the 'energy_pill'
1920 : : initialization. Linear initiation uses a linear function in time and space,
1921 : : configured with an initiation point in space, a constant velocity of the
1922 : : growing spherical front in time (and space) linearly, and width of the front
1923 : : and assigns values to mesh points falling within the growing spherical shell
1924 : : inside a configured box or meshblock. The following keywords are required
1925 : : in the box/meshblock block if 'linear' is used: 'init_time',
1926 [ + - ][ + - ]: 1348 : 'front_width', 'front_speed')"});
[ + - ][ + - ]
[ + - ][ + - ]
1927 : :
1928 : : // -----------------------------------------------------------------------
1929 : : // Overset mesh object
1930 : : // -----------------------------------------------------------------------
1931 : :
1932 : : keywords.insert({"mesh",
1933 : : "Start configuration block assigning a mesh to a solver",
1934 : : R"(This keyword is used to introduce a mesh block, used to
1935 [ + - ][ + - ]: 1348 : assign and configure a mesh to a solver.)", "vector block-title"});
[ + - ][ + - ]
[ + - ][ + - ]
1936 : :
1937 : : keywords.insert({"filename", "Set filename",
1938 [ + - ][ + - ]: 1348 : R"(Set filename, e.g., mesh filename for solver coupling.)", "string"});
[ + - ][ + - ]
[ + - ][ + - ]
1939 : :
1940 : : keywords.insert({"location", "Configure location",
1941 : : R"(Configure location of a mesh relative to another.)",
1942 [ + - ][ + - ]: 1348 : "vector of 3 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1943 : :
1944 : : keywords.insert({"moment_of_inertia", "Moment of inertia of rigid body",
1945 [ + - ][ + - ]: 1348 : R"(Moment of inertia of rigid body for rotational motion)", "real"});
[ + - ][ + - ]
[ + - ][ + - ]
1946 : :
1947 : : keywords.insert({"center_of_mass", "Center of mass of rigid body",
1948 : : R"(Center of mass of rigid body used to compute torque for rotational
1949 [ + - ][ + - ]: 1348 : motion)", "vector of 3 reals"});
[ + - ][ + - ]
[ + - ][ + - ]
1950 : :
1951 : : // -----------------------------------------------------------------------
1952 : : // pre-configured problems
1953 : : // -----------------------------------------------------------------------
1954 : :
1955 : : keywords.insert({"user_defined",
1956 : : "Select user-defined specification for a problem",
1957 : : R"(This keyword is used to select the user-defined specification for an
1958 : : option. This could be a 'problem' to be solved by a partial differential
1959 : : equation, but can also be a 'user-defined' mesh velocity specification for
1960 [ + - ][ + - ]: 1348 : ALE mesh motion.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1961 : :
1962 : : keywords.insert({"shear_diff",
1963 : : "Select the shear + diffusion test problem ",
1964 : : R"(This keyword is used to select the shear diffusion test problem. The
1965 : : initial and boundary conditions are specified to set up the test problem
1966 : : suitable to exercise and test the advection and diffusion terms of the
1967 [ + - ][ + - ]: 1348 : scalar transport equation.)" });
[ + - ][ + - ]
[ + - ][ + - ]
1968 : :
1969 : : keywords.insert({"slot_cyl",
1970 : : "Select Zalesak's slotted cylinder test problem",
1971 : : R"(This keyword is used to select Zalesak's slotted cylinder test
1972 : : problem. The initial and boundary conditions are specified to set up the
1973 : : test problem suitable to exercise and test the advection and diffusion
1974 [ + - ][ + - ]: 1348 : terms of the scalar transport equation.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1975 : :
1976 : : keywords.insert({"gauss_hump",
1977 : : "Select advection of 2D Gaussian hump test problem",
1978 : : R"(This keyword is used to select the advection of 2D Gaussian hump test
1979 : : problem. The initial and boundary conditions are specified to set up the
1980 : : test problem suitable to exercise and test the advection
1981 [ + - ][ + - ]: 1348 : terms of the scalar transport equation.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1982 : :
1983 : : keywords.insert({"cyl_advect",
1984 : : "Select advection of cylinder test problem",
1985 : : R"(This keyword is used to select the advection of cylinder test
1986 : : problem. The initial and boundary conditions are specified to set up the
1987 : : test problem suitable to exercise and test the advection
1988 [ + - ][ + - ]: 1348 : terms of the scalar transport equation.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1989 : :
1990 : : keywords.insert({"cyl_vortex",
1991 : : "Select deformation of cylinder in a vortex test problem",
1992 : : R"(This keyword is used to select the test problem which deforms a cylinder
1993 : : in a vortical velocity field. The initial and boundary conditions are
1994 : : specified to set up the test problem suitable to exercise and test the
1995 [ + - ][ + - ]: 1348 : advection terms of the scalar transport equation.)"});
[ + - ][ + - ]
[ + - ][ + - ]
1996 : :
1997 : : keywords.insert({"vortical_flow",
1998 : : "Select the vortical flow test problem ",
1999 : : R"(This keyword is used to select the vortical flow test problem. The
2000 : : purpose of this test problem is to test velocity errors generated by spatial
2001 : : operators in the presence of 3D vorticity and in particluar the
2002 : : superposition of planar and vortical flows, analogous to voritcity
2003 : : stretching. For more details, see Waltz,
2004 : : et. al, "Manufactured solutions for the three-dimensional Euler equations
2005 : : with relevance to Inertial Confinement Fusion", Journal of Computational
2006 [ + - ][ + - ]: 1348 : Physics 267 (2014) 196-209.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2007 : :
2008 : : keywords.insert({"nl_energy_growth",
2009 : : "Select the nonlinear energy growth test problem",
2010 : : R"(This keyword is used to select the nonlinear energy growth test problem.
2011 : : The purpose of this test problem is to test nonlinear, time dependent energy
2012 : : growth and the subsequent development of pressure gradients due to coupling
2013 : : between the internal energy and the equation of state. For more details,
2014 : : see Waltz, et. al, "Manufactured
2015 : : solutions for the three-dimensional Euler equations with relevance to
2016 : : Inertial Confinement Fusion", Journal of Computational Physics 267 (2014)
2017 [ + - ][ + - ]: 1348 : 196-209.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2018 : :
2019 : : keywords.insert({"rayleigh_taylor",
2020 : : "Select the Rayleigh-Taylor test problem ",
2021 : : R"(This keyword is used to select the Rayleigh-Taylor unstable configuration
2022 : : test problem. The purpose of this test problem is to assess time dependent
2023 : : fluid motion in the presence of Rayleigh-Taylor unstable conditions, i.e.
2024 : : opposing density and pressure gradients.
2025 : : For more details, see Waltz, et. al, "Manufactured solutions for the
2026 : : three-dimensional Euler equations with relevance to Inertial Confinement
2027 [ + - ][ + - ]: 1348 : Fusion", Journal of Computational Physics 267 (2014) 196-209.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2028 : :
2029 : : keywords.insert({"taylor_green",
2030 : : "Select the Taylor-Green test problem ",
2031 : : R"(This keyword is used to select the Taylor-Green vortex test problem. The
2032 : : purpose of this problem is to test time accuracy and the correctness of the
2033 : : discretization of the viscous term in the Navier-Stokes equation. For more
2034 : : details on the flow, see G.I. Taylor, A.E.
2035 : : Green, "Mechanism of the Production of Small Eddies from Large Ones", Proc.
2036 : : R. Soc. Lond. A 1937 158 499-521; DOI: 10.1098/rspa.1937.0036. Published 3
2037 [ + - ][ + - ]: 1348 : February 1937.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2038 : :
2039 : : keywords.insert({"sod_shocktube",
2040 : : "Select the Sod shock-tube test problem ",
2041 : : R"(This keyword is used to select the Sod shock-tube test problem. The
2042 : : purpose of this test problem is to test the correctness of the
2043 : : approximate Riemann solver and its shock and interface capturing
2044 : : capabilities. For more details, see
2045 : : G. A. Sod, "A Survey of Several Finite Difference Methods for Systems of
2046 : : Nonlinear Hyperbolic Conservation Laws", J. Comput. Phys., 27 (1978)
2047 [ + - ][ + - ]: 1348 : 1–31.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2048 : :
2049 : : keywords.insert({"rotated_sod_shocktube",
2050 : : "Select the rotated Sod shock-tube test problem ",
2051 : : R"(This keyword is used to select the rotated Sod shock-tube test problem.
2052 : : This the same as Sod shocktube but the geometry is rotated about X, Y, Z
2053 : : each by 45 degrees (in that order) so that none of the domain boundary align
2054 : : with any of the coordinate directions. The purpose of this test problem is
2055 : : to test the correctness of the approximate Riemann solver and its shock and
2056 : : interface capturing capabilities in an arbitrarily oriented geometry.
2057 : : For more details on the Sod
2058 : : problem, see G. A. Sod, "A Survey of Several Finite Difference Methods for
2059 : : Systems of Nonlinear Hyperbolic Conservation Laws", J. Comput. Phys., 27
2060 [ + - ][ + - ]: 1348 : (1978) 1–31.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2061 : :
2062 : : keywords.insert({"shedding_flow",
2063 : : "Select the Shedding flow test problem ",
2064 : : R"(This keyword is used to select the Shedding flow test problem. It
2065 : : describe a quasi-2D inviscid flow over a triangular wedge in tetrahedron
2066 : : grid. The purpose of this test problem is to test the capability of DG
2067 : : scheme for retaining the shape of vortices and also different error
2068 : : indicator behavior for this external flow problem when p-adaptive DG scheme
2069 [ + - ][ + - ]: 1348 : is applied.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2070 : :
2071 : : keywords.insert({"sedov_blastwave",
2072 : : "Select the Sedov blast-wave test problem ",
2073 : : R"(This keyword is used to select the Sedov blast-wave test problem. The
2074 : : purpose of this test problem is to test the correctness of the
2075 : : approximate Riemann solver and its strong shock and interface capturing
2076 [ + - ][ + - ]: 1348 : capabilities.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2077 : :
2078 : : keywords.insert({"interface_advection",
2079 : : "Select the interface advection test problem ",
2080 : : R"(This keyword is used to select the interface advection test problem.
2081 : : The purpose of this test problem is to test the well-balancedness of the
2082 : : multi-material discretization and its interface capturing
2083 [ + - ][ + - ]: 1348 : capabilities.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2084 : :
2085 : : keywords.insert({"gauss_hump_compflow",
2086 : : "Select advection of 2D Gaussian hump test problem",
2087 : : R"(This keyword is used to select the advection of 2D Gaussian hump test
2088 : : problem. The initial and boundary conditions are specified to set up the
2089 : : test problem suitable to exercise and test the advection terms of the
2090 : : Euler equations. The baseline of the density distribution in this testcase
2091 : : is 1 instead of 0 in gauss_hump_transport which enables it to be the
2092 [ + - ][ + - ]: 1348 : regression testcase for p-adaptive DG scheme.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2093 : :
2094 : : keywords.insert({"waterair_shocktube",
2095 : : "Select the water-air shock-tube test problem ",
2096 : : R"(This keyword is used to select the Water-air shock-tube test problem.
2097 : : The purpose of this test problem is to test the correctness of the
2098 : : multi-material pressure relaxation procedure and its interface capturing
2099 : : capabilities. For more details, see
2100 : : Chiapolino, A., Saurel, R., & Nkonga, B. (2017). Sharpening diffuse
2101 : : interfaces with compressible fluids on unstructured meshes. Journal of
2102 [ + - ][ + - ]: 1348 : Computational Physics, 340, 389-417.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2103 : :
2104 : : keywords.insert({"shock_hebubble",
2105 : : "Select the shock He-bubble test problem ",
2106 : : R"(This keyword is used to select the shock He-bubble test problem. The
2107 : : purpose of this test problem is to test the correctness of the
2108 : : multi-material algorithm and its shock-interface interaction
2109 : : capabilities. For more details, see
2110 : : Quirk, J. J., & Karni, S. (1996). On the dynamics of a shock–bubble
2111 [ + - ][ + - ]: 1348 : interaction. Journal of Fluid Mechanics, 318, 129-163.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2112 : :
2113 : : keywords.insert({"underwater_ex",
2114 : : "Select the underwater explosion test problem ",
2115 : : R"(This keyword is used to select the underwater explosion test problem.
2116 : : The purpose of this test problem is to test the correctness of the
2117 : : multi-material algorithm and its interface capturing capabilities in the
2118 : : presence of strong shocks and large deformations.
2119 : : For more details, see
2120 : : Chiapolino, A., Saurel, R., & Nkonga, B. (2017). Sharpening diffuse
2121 : : interfaces with compressible fluids on unstructured meshes. Journal of
2122 [ + - ][ + - ]: 1348 : Computational Physics, 340, 389-417.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2123 : :
2124 : : keywords.insert({"shockdensity_wave",
2125 : : "Select the shock-density wave test problem ",
2126 : : R"(This keyword is used to select the shock-density wave test problem.
2127 : : THe purpose of this test problem is to assess the accuracy of high order
2128 : : method in predicting the interaction of a density wave with a shock front.
2129 : : For more details, see Yu, L., Matthias
2130 : : I. (2014). Discontinuous Galerkin method for multicomponent chemically
2131 : : reacting flows and combustion. Journal of Computational Physics, 270,
2132 [ + - ][ + - ]: 1348 : 105-137.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2133 : :
2134 : : keywords.insert({"equilinterface_advect",
2135 : : "Select the advection of equilibrium interface problem ",
2136 : : R"(This keyword is used to select the advection of equilibrium interface
2137 : : problem. This is a manufactured problem with source terms with nonlinear
2138 : : solutions near the material interface. Source terms are used to ensure
2139 [ + - ][ + - ]: 1348 : that the conservation laws are satisfied by the manufactured solution.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2140 : :
2141 : : keywords.insert({"sinewave_packet",
2142 : : "Select the advection of sinewave packet problem ",
2143 : : R"(This keyword is used to select the advection of sinewave packet
2144 [ + - ][ + - ]: 1348 : problem.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2145 : :
2146 : : keywords.insert({"richtmyer_meshkov",
2147 : : "Select the Richtmyer-Meshkov instability problem ",
2148 : : R"(This keyword is used to select the Richtmyer-Meshkov instability
2149 [ + - ][ + - ]: 1348 : problem. In this problem, a shock hits a perturbed material interface.)"});
[ + - ][ + - ]
[ + - ][ + - ]
2150 : :
2151 : : // -----------------------------------------------------------------------
2152 : :
2153 : : // Initialize help: fill own keywords
2154 : 1348 : tk::ctr::Info ctrinfoFill(get< tag::cmd, tag::ctrinfo >());
2155 [ + + ]: 342392 : for (const auto& i : keywords) {
2156 [ + - ]: 341044 : ctrinfoFill.fill(i);
2157 : : }
2158 : 1348 : }
2159 : :
2160 : : //! Query scheme centering
2161 : : //! \return Scheme centering
2162 : 192 : tk::Centering centering() const
2163 [ + - ]: 192 : { return ctr::Scheme().centering( get< tag::scheme >() ); }
2164 : :
2165 : : /** @name Pack/Unpack: Serialize InputDeck object for Charm++ */
2166 : : ///@{
2167 : : //! \brief Pack/Unpack serialize member function
2168 : : //! \param[in,out] p Charm++'s PUP::er serializer object reference
2169 : 1790 : void pup( PUP::er& p ) { tk::TaggedTuple< ConfigMembers >::pup(p); }
2170 : : //! \brief Pack/Unpack serialize operator|
2171 : : //! \param[in,out] p Charm++'s PUP::er serializer object reference
2172 : : //! \param[in,out] c InputDeck object reference
2173 : 1790 : friend void operator|( PUP::er& p, InputDeck& c ) { c.pup(p); }
2174 : : //@}
2175 : :
2176 : : };
2177 : :
2178 : : } // ctr::
2179 : : } // inciter::
2180 : :
2181 : : #endif // InputDeck_h
|