DALES interface

This is the OMUSE interface to DALES, the Dutch Atmospheric Large-Eddy Simulation.

Example:

from omuse.community.dales.interface import Dales
from omuse.units import units

d = Dales(workdir='dales-workdir', channel_type='sockets', number_of_workers=1)

# set parameters
d.parameters_DYNAMICS.iadv_mom = 6  # 6th order advection for momentum
d.parameters_DYNAMICS.iadv_thl = 5  # 5th order advection for scalars, less overshoots than 6th order
d.parameters_DYNAMICS.iadv_qt = 5
d.parameters_DYNAMICS.iadv_tke = 5

# access model variables
d.fields[:, :, :].U = 0 | units.m / units.s

# evolve the model
d.evolve_model(120 | units.s, exactEnd=True)

OMUSE models use variables with units, see Units in OMUSE.

The parameters are doumented in the DALES documentation: options. Setting parameters is possible before the model is time stepped or the model variables have been accessed, but not after.

Variable grids

various variables of the DALES model can be accessed from the OMUSE interface. The variables are grouped in grids, according to their dimensionality. For example, the U-velocity component in the model component in the model can be accessed as d.fields[i,j,k].U.

The following grids are defined:

Grid Description
fields 3D grid
profiles horizontal averages of the 3D grid (vertical profiles)
nudging_profiles vertical profiles of U, V, THL, QT to nudge the model towards.
forcing_profiles external forcings of U, V, THL, QT - vertical profiles
scalars surface fluxes and roughness lengths, scalars.
surface_fields 2D surface fields, e.g. liquid water path.

The grids contain the following variables:

Grid Variables
fields U, V, W, THL, QT, QL, QL_ice, QR, E12, T, pi, rswd, rswdir, rswdif,
rswu, rlwd, rlwu, rswdcs, rswucs, rlwdcs, rlwucs
profiles U, V, W, THL, QT, QL, QL_ice, QR, E12, T, P, rho
nudging_profiles U, V, THL, QT
forcing_profiles U, V, THL, QT
scalars wt, wq, z0m, z0h, Ps, QR
surface_fields LWP, RWP, TWP, ustar, z0m, z0h, tskin, qskin, LE, H, obl

Description of the variables:

Variable Unit Description
U,V,W m/s velocity components in x, y, and z directions
THL K liquid water potential temperature
QT kg/kg specific total humidity (i.e. vapor and cloud condensate but excluding precipitation)
QL kg/kg specific cloud condensate
QL_ice kg/kg specific cloud condensate in the form of ice
QR kg/kg precipitation
E12 m/s sqrt(subgrid-scale turbulent kinetic energy)
T K temperature
P Pa pressure
Ps Pa surface pressure
rho kg/m^3 density
pi Pa modified air pressure
rswu,rswd W/m^2 up-/down-welling short-wave radiative flux
rlwu,rlwd W/m^2 up-/down-welling long-wave radiative flux
r{s,l}w{u,d}cs W/m^2 up/down-welling long/short-wave radiative flux for clear sky
rswdir W/m^2 downwelling shortwave direct radiative flux
rswdif W/m^2 downwelling shortwave diffusive radiative flux
wt K m/s surface flux of heat
wq m/s surface flux of humidity
z0m, z0h m surface roughness length
ustar m/s friction velocity
LWP, RWP, TWP kg/m^2 liquid-, rain-, total water path
tskin K skin temperature
qskin kg/kg skin humidity
H, LE W/m^2 sensible and latent heat fluxes
obl m Obukhov length

The forcing and nudging profiles, and the 3D grids for prognostic variables can be read and written. Diagnosed quantities, horizontal averages, and water paths can only be read.

Note: the indexing brackets can also be placed on the variable name instead of on the grid name, e.g. d.fields.U[:, :, :] vs d.fields[:, :, :].U. Avoid using the first form, since assigning new values to the grid this way does not work.

Example scripts

Some examples of using Dales with OMUSE are bundled included in omuse/src/omuse/community/dales/example/ .

bubble.py

A Dales experiment with a warm air bubble. Shows model initialization, setting model parameters, setting initial conditions, evolving the model, and retreiving the state.

async.py

Shows how to use OMUSE’s asynchronous function calls, to let several model instances run concurrently. See also Asynchronous calls.

class omuse.community.dales.interface.Dales(**options)

OMUSE Dales Interface.

Parameters:
  • workdir (str, optional) – Working directory for DALES. Output files are placed here. If workdir doesn’t exist, and either inputdir or case is given, input files are copied from the provided directory to workdir. If workdir doesn’t exist, and no case or input files are provided, built-in defaults are used.
  • exp (int, optional) – Experiment number, used to number input files, e.g. prof.inp.001 .
  • inputdir (str, optional) – Directory for input files, copied to workdir.
  • case (str, optional) – Specify one of the cases bundled with DALES. Valid names include ‘bomex’, ‘rico’, ‘atex’, ‘fog’. Input files are copied to workdir.
  • z (numpy array, optional) – Override the vertical discretization with an array of increasing heights. The array is assumed to have a length unit attached to it.
  • interpolator (function handle, optional) – In case of a user-specified vertical discretization, this function determines the interpolation method used to obtain the initial profiles at the desired resolution. Should be a function f(str,z_new,z_old,y) where the first argument denotes the variable, the second and third resp. the new and old z-axes and the latter the profile values on the opriginal axis, returning an array of values on the new axis. Default = None, meaning linear interpolation.
  • number_of_workers (int, optional) – Number of MPI tasks to use. General OMUSE option. Default = 1.
  • channel_type (str, optional) – Communication channel between Python and DALES worker processes. Options ‘mpi’ or ‘sockets’. General OMUSE option.
  • redirection (str, optional) – Options for re-directing stdout and stderr of the DALES process. General OMUSE option. ‘none’ for no redirection, ‘file’ for redirection to files.
  • redirect_stdout_file (str, optional.) – File name for redirection of stdout, see above. General OMUSE option.
  • redirect_stderr_file (str, optional.) – File name for redirection of stderr, see above. General OMUSE option.
get_dx()

Dales grid cell size in x-direction

Returns:Grid resolution (in m) along U-direction
Return type:float
get_dy()

Dales grid cell size in y-direction

Returns:Grid resolution (in m) along V-direction
Return type:float
get_field(field, imin=1, imax=None, jmin=1, jmax=None, kmin=1, kmax=None, **kwargs)

Dales volume field retrieval method

Parameters:
  • field (str) – Variable shortname: either LWP, RWP, TWP, U, V, W, THL, T, QT, QL, Qsat
  • imin (integer, optional) – Lower one-based x-bound of data block
  • imax (integer, optional) – Upper one-based x-bound of data block
  • jmin (integer, optional) – Lower one-based y-bound of data block
  • jmax (integer, optional) – Upper one-based y-bound of data block
  • kmin (integer, optional) – Lower one-based z-bound of data block
  • kmax (integer, optional) – Upper one-based z-bound of data block
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

3D block containing variable values

Return type:

numpy.array

get_itot()

Dales number of grid cells in x-direction

Returns:Number of grid cells along U-direction
Return type:int
get_jtot()

Dales number of grid cells in y-direction

Returns:Number of grid cells along V-direction
Return type:int
get_ktot()

Dales number of grid cells in z-direction

Returns:Number of vertical layers
Return type:int
get_presf(k=None, **kwargs)

Dales full level mean pressure retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Full level mean pressure profile

Return type:

numpy.array

get_presh(k=None, **kwargs)

Dales half level mean pressure retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Half level mean pressure profile

Return type:

numpy.array

get_profile(field, **kwargs)

Dales generic profile retrieval method

Parameters:
  • field (str) – Variable shortname: either U, V, W, THL, T, QT, QL, E12
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

1D array containing mean vertical profile values

Return type:

numpy.array

get_profile_E12(k=None, **kwargs)

Dales turbulence kinetic energy profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Turbulence kinetic energy vertical profile

Return type:

numpy.array

get_profile_QL(k=None, **kwargs)

Dales liquid water content profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Liquid water content vertical profile

Return type:

numpy.array

get_profile_QL_ice(k=None, **kwargs)

Dales ice water content profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Ice water content vertical profile

Return type:

numpy.array

get_profile_QR(k=None, **kwargs)

Dales rain water content profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Rain water content vertical profile

Return type:

numpy.array

get_profile_QT(k=None, **kwargs)

Dales total humidity profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Total humidity vertical profile

Return type:

numpy.array

get_profile_T(k=None, **kwargs)

Dales temperature profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Temperature vertical profile

Return type:

numpy.array

get_profile_THL(k=None, **kwargs)

Dales liquid water virtual temperature profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Liquid water virtual temperature vertical profile

Return type:

numpy.array

get_profile_U(k=None, **kwargs)

Dales eastward wind profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Eastward vertical wind profile

Return type:

numpy.array

get_profile_V(k=None, **kwargs)

Dales northward wind profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Northward vertical wind profile

Return type:

numpy.array

get_profile_W(k=None, **kwargs)

Dales upward wind profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Upward vertical wind profile

Return type:

numpy.array

get_rhobf(k=None, **kwargs)

Dales base density profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Full level density base profile

Return type:

numpy.array

get_rhof(k=None, **kwargs)

Dales mean density profile retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Full level mean density profile

Return type:

numpy.array

get_xsize()

Dales domain extent in x-direction

Returns:Domain length (in m) along U-direction
Return type:float
get_ysize()

Dales domain extent in y-direction

Returns:Domain length (in m) along V-direction
Return type:float
get_zf(k=None, **kwargs)

Dales full level heights retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Full level heights

Return type:

numpy.array

get_zh(k=None, **kwargs)

Dales half level heights retrieval method

Parameters:
  • k (integer array, optional) – Restrict profile to this set of vertical indices.
  • async (boolean, optional) – Execute function asynchronously, return request object
Returns:

Half level heights

Return type:

numpy.array

set_field(field, a, imin=1, jmin=1, kmin=1, **kwargs)

Dales volume field insertion method

Parameters:
  • field (str) – prognostic variable shortname: either U, V, W, THL, QT
  • a (numpy.array) – Block of values for substitution in state
  • imin (integer, optional) – Lower one-based x-bound of data block
  • jmin (integer, optional) – Lower one-based y-bound of data block
  • kmin (integer, optional) – Lower one-based z-bound of data block
  • async (boolean, optional) – Execute function asynchronously, return request object