Page cover

PROB

The prob.h file creates the namespace PROB, that defines the probelm to be solved and the methods to do it. Check the file from one of the examples(LINK) as a starting point

Problem parameters

The code requires a ProbParm structure, which is passed to user deifned functions and allows to centralised one place where varibles are stored.

struct ProbParm {
  Real p_l = 1.0;
  Real p_r = 0.1;
  Real rho_l = 1.0;
  Real rho_r = 0.125;
  Real u_l = 0.0;
  Real u_r = 0.0;
};

The file uses thet AMReX type Real instead of floats, which can be single or double precision. Check GNU_Makefilefor set-up Makefile

Global Constants

In some cases in convenient to defined a a few global constants , such as

static constexpr Real Reynolds = 3000.0;  

These constants are defined as static constexpr, allowing them to be evaluated at compile time. This makes them efficient for use in other compile-time expressions. Global constants can be utilized within methods, the ProbParm structure, and as parameters to methods Unlike the ProbParm structure, global constants are directly accessible in all user-defined functions.

NOTE: By default Cerisse uses SI units (unlike PeleC) or no-units, include the file "Constants.h" to access conversion factors (CGS to SI and viceversa) as well as universal constants with appropiate precision.

Themodynamic and Transport Closures

The namespace PROB needs to define the closures_dt template, which creates the physical "closures" of the problems, that means which themrodynamics modesl are used as well as transport properties For example:

typedef closures_dt<indicies_t, visc_suth_t, cond_suth_t,
                    calorifically_perfect_gas_t<indicies_t>> ProbClosures;

The line above selects a problem that uses the structure indicies_t, which variables to solve(and how are they stored) how many to solve. The above also selects the Sutherland viscosity model visc_suth_t and condictivity model cond_suth_t. As well as use a perdect ideal gas as thermodynamics model. All this is wrapped in the class ```ProbClosures``, which will then pass to the equations. Some of the tranpsort model required input from the user, which can use default values (or probelem specific). Options in closures_dt, all units in SI

Closure
Options
defaults
Description

visc_suth_t

no

μ=1.458106\mu =1.458\cdot 10^{-6} at T=110.4T=110.4

Sutherland Viscosity

cond_suth_t

no

λ=2.495103\lambda =2.495 \cdot 10^{-3} at T=194T=194

Sutherland Conductivity

visc_const_t

viscosity

μ=1.85105\mu= 1.85 \cdot 10^{-5}

Constant Viscosity

cond_const_t

conductivity

λ=0.0262\lambda = 0.0262

Constant Conductivity

transport_const_t

viscosity

μ=1.85105\mu = 1.85 \cdot 10^{-5} and λ=0.0262\lambda = 0.0262

Constant Viscosity and

calorifically_perfect_gas_t

no

γ=1.4\gamma =1.4 M=28.96103M=28.96 \cdot 10^{-3}

perfect

multispecies_perfect_gas_t

no

PelePhysics

Used for PelePhysics options

Smagorinsky_t

yes

none

Used for LES

Passing Arguments

Oprions and parameters can be passed to the closures by creating a small structure methodparm_t before the assembly of closures_dt, so all parameters are changed directly in prob.h. For example:

struct methodparm_t {

  public:

  static constexpr Real viscosity   = 2.85e-5; 
  
};
typedef closures_dt<indicies_t, visc_const_t<methodparm_t>, cond_const_t<defaultparm_t>,
                    calorifically_perfect_gas_t<indicies_t>> ProbClosures;

The problem will use a viscosity of 2.85e-5 and a conductivity of 0.0262 (default values using the structure defaultparm_t ). Check Options for detailed description on arguments.

To use the default values, the line #include <NumParam.h> has to be included in the headers

Numerical Method and Equations

The namespace PROB needs to define as well the rhs_dt template, which defines which tre,s in the equatison to solve and how. In general

where the RHS includes the inviscid (Euler) terms viscous (for Navier-Stokes) and source terms. The C++ template follows:

template <typename euler, typename diffusive, typename source>

Splitting the RHS into Euler terms, diffusive (or viscous) terms and source terms, which need to be defined For example, the line below

typedef rhs_dt<skew_t<methodparm_t, ProbClosures>, no_diffusive_t, no_source_t > ProbRHS;

will solve

Ut=F(U)\frac{\partial U}{\partial t} = F(U)

corresponding to the Euler equations, with skew-symmetric numerical scheme (with order defined in methodparm_t, similar to closures_dt (see Options). Available options in rhs_dt are:

euler
Options
IBM
defaults
Description

skew_t

yes

yes

4th order, no dissipation

keep_euler_t

yes

no

no default

weno_t

yes

no

no default

WENO or TENO 5th order scheme

no_euler_t

no

yes

-

0 (not solving Euler)

diffusive options in rhs_dt (October 2024)

diffusive
Options
IBM
defaults
Description

diffusiveheat_t

yes

no

4th order

2/4/6 central scheme only heat

skew_t

yes

yes

4th order

2/4/6 central scheme

no_diffusive_t

no

yes

0

0 (not diffusive part)

source options in rhs_dt (October 2024)

diffusive
Options
IBM
defaults
Description

reactor_t

yes

yes

-

PelePhysics chemcail raection

user_source_t

no

yes

-

user-given source term

no_source_t

no

yes

0

0 (no source term part)

Initial Conditions

The prob_initdata function is called for every i,j,k cell

// initial condition
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void
prob_initdata(int i, int j, int k, Array4<Real> const &state,
              GeometryData const &geomdata, ProbClosures const &cls,
              ProbParm const &prob_parm) {

In most examples, a few auxiliar definitions follow that extract the size the domain and allow to define the spatial coordinate x (and y and z) of the cell

  const Real *prob_lo = geomdata.ProbLo();
  const Real *prob_hi = geomdata.ProbHi();
  const Real *dx = geomdata.CellSize();
  Real x = prob_lo[0] + (i + Real(0.5)) * dx[0];

This allow to define, for example, different condition depending on the x coordinate. The structure prob_parmof type ProbParm is used to recall generic problem parameters

  Real Pt, rhot, uxt;
  if (x < prob_hi[0] / 2) {
    Pt = prob_parm.p_l;
    rhot = prob_parm.rho_l;
    uxt = prob_parm.u_l;
  } else {
    Pt = prob_parm.p_r;
    rhot = prob_parm.rho_r;
    uxt = prob_parm.u_r;
  }
  Real et = Pt / (cls.gamma - Real(1.0));

The function needs to fille the state array, where the conservative variables exist. The index URHO, UMX are defined in the template indicies_t and are accessed through the class as cls.URHO and so on.

  state(i, j, k, cls.URHO) = rhot;
  state(i, j, k, cls.UMX)  = rhot * uxt;
  state(i, j, k, cls.UMY)  = Real(0.0);
  state(i, j, k, cls.UMZ)  = Real(0.0);
  state(i, j, k, cls.UET)  = et + Real(0.5) * rhot * uxt * uxt;
}

Since this is pure C++ code, it allows for the construction of complex initial conditions. Additionally the Utility class (see below) can be used to incorporate chemistry profiles, turbulent inflows, etc.

Source

A user-specific source term can be created by defining a template

template <typename cls_t > class user_source_t;
template <typename cls_t>
class user_source_t {
  public:
  void inline src(const amrex::MFIter &mfi,
                  const amrex::Array4<const amrex::Real> &prims,
                  const amrex::Array4<amrex::Real> &rhs, const cls_t *cls_d,
                  amrex::Real dt){
    const Box bx = mfi.tilebox();
    ProbParm const prob_parm;
    amrex::ParallelFor(bx,
      [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
      {
      const auto& cls = *cls_d;
      //  Pressure Source           
      rhs(i,j,k,cls.UMX) += prob_parm.dpdx;  
      });
  };
};

Boundary Condition

The prob.h can be used to implement user-specific boundary conditions (transient, turbulent, etc). To do that, a local template is defined (see boundaries tab for detail).

Utility

The utility class can be used to incorporate chemistry profiles, turbulent inflows etc, in the initial and boundary conditions. It inspired in the Utility options in PelePhysics and is a placeholder for complex interactions not specified in this fie, such as new data read from files. Utility can access data from input and external files. Conditions and parameters defined in PROB are known at compile time. To use the Utility requires an additional line in GNU_Makefile and to change the prob_initdata to allow an additional argument

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void prob_initdata(
    int i, int j, int k, Array4<Real> const &state,
    GeometryData const &geomdata, ProbClosures const &cls,
    ProbParm const &prob_parm, Utility* util = nullptr) {

This class allows access to multiples functions, for example PMF in PelePhysics.

Other (optional)

The names follow in the cons_var_names array

inline Vector<std::string> cons_vars_names={"Xmom","Ymom","Zmom","Energy","Energy"};

The type of variables, keep as it is, scalars are 0 and vectors are given by their components

inline Vector<int> cons_vars_type={1,2,3,0,0};

Data missed from input file

void inline inputs() {
  ParmParse pp;

  pp.add("cns.order_rk", 3);   // -2, 1, 2 or 3"
  pp.add("cns.stages_rk", 3);  // 1, 2 or 3
}

Last updated