Reference

Model grid

ROMS.create_gridFunction
ROMS.create_grid(fname,h,f,lon_r,lat_r,mask_r,angle,pm,pn,dndx,dmde)

Create a NetCDF grid file fname using the bathymetry h, Coriolis parameter f and longitude, latitude, mask, angle and strechting factors are rho-points.

Note

This function currently only works for non-rotated grids (angle = 0) and the spherical grids.

source
ROMS.smoothgridFunction
hsmooth = smoothgrid(h,hmin,rmax)

Smooth the topography to get a maximum r factor equalt to rmax

The original code is from Pierrick Penven from Roms_tools. Copyright (c) 2002-2006 by Pierrick Penven, GPL

source
ROMS.reduce_resFunction
 field = reduce_res(field,red);

reduce resolution of field field by a factor "red"

source
ROMS.stretchingFunction
s,C = stretching(Vstretching, theta_s, theta_b, hc, N, kgrid, report)

Compute ROMS vertical coordinate stretching function.

Given vertical terrain-following vertical stretching parameters, this routine computes the vertical stretching function used in ROMS vertical coordinate transformation. Check the following (link)[https://www.myroms.org/wiki/index.php/Vertical_S-coordinate] for details:

On Input:

  • Vstretching: Vertical stretching function:
                   Vstretching = 1,  original (Song and Haidvogel, 1994)
                   Vstretching = 2,  A. Shchepetkin (UCLA-ROMS, 2005)
                   Vstretching = 3,  R. Geyer BBL refinement
                   Vstretching = 4,  A. Shchepetkin (UCLA-ROMS, 2010)
                   Vstretching = 5,  Quadractic (Souza et al., 2015)
  • theta_s: S-coordinate surface control parameter (scalar)
  • theta_b: S-coordinate bottom control parameter (scalar)
  • hc: Width (m) of surface or bottom boundary layer in which higher vertical resolution is required during stretching (scalar)
  • N: Number of vertical levels (scalar)
  • kgrid Depth grid type logical switch: kgrid = 0, function at vertical RHO-points kgrid = 1, function at vertical W-points
  • report Flag to report detailed information (OPTIONAL): report = false, do not report report = true, report information

On Output:

  • s: S-coordinate independent variable, [-1 <= s <= 0] at vertical RHO- or W-points (vector)
  • C: Nondimensional, monotonic, vertical stretching function, C(s), 1D array, [-1 <= C(s) <= 0]

The code is ported from the matlab code stretching.m from Hernan G. Arango. Copyright (c) 2002-2020 The ROMS/TOMS Group, Licensed under a MIT/X style license

source
ROMS.GridType
grid = ROMS.Grid(grid_fname,opt)

Loads the ROMS grid from a NetCDF file grid_fname. The grid structure contains the longitude, latitude, depth (z), angle and mask at rho, u, v and w points for a C-grid.

Note

grid.z corresponds to an elevation equal to zero.

Example

opt = (
    Tcline = 50,   # m
    theta_s = 5,   # surface refinement
    theta_b = 0.4, # bottom refinement
    nlevels = 32,  # number of vertical levels
    Vtransform  = 2,
    Vstretching = 4,
)
grid = ROMS.Grid(expanduser("~/Models/LS2v/LS2v.nc"),opt)
source
ROMS.generate_gridFunction
ROMS.generate_grid(grid_fname,bath_name,xr,yr,red,opt,hmin,rmax)

Generate the netCDF file grid_fname representing the model grid from the bathymetry file bath_name. The domain is bounded by the longitude range xr and the latitude range yr. The resolution of the bathymetry is reduced by the factors red[1] and red[2] in the longitude and latitude directions. opt is tuple with additional parameters describing the vertical axis. hmin and rmax is the minimum depth and rmax is the smoothness parameter of the bathymetry. The $r$ parameter is defined as:

$r = \max\left( \underset{ij}{\max} \frac{|h_{i,j} - h_{i+1,j}|}{h_{i,j} + h_{i+1,j}}, \underset{ij}{\max} \frac{|h_{i,j} - h_{i,j+1}|}{h_{i,j} + h_{i,j+1}} \right)$

The parameter opt contains for example:

# model specific parameters
opt = (
    Tcline = 50,   # m
    theta_s = 5,   # surface refinement
    theta_b = 0.4, # bottom refinement
    nlevels = 32,  # number of vertical levels
    Vtransform  = 2,
    Vstretching = 4,
)
source
ROMS.nudgecoefFunction
nudgecoef(domain,nudge_filename,alpha,Niter,halo,tscale; max_tscale = 5e5)

Generate trace nudging coefficients with a value of 1/tscale at the open boundaries over a halo grid cells. coefficients field smoothed spatially with a diffusion coefficient alpha over Niter iterations. For grid cells where the nudging time scale exceeds max_tscale, nudging is disabled (coefficient is set to zero).

source
ROMS.diffusion2!Function
 diffusion2!(f,alpha,Niter)

Two-dimensional diffusion of field f. f is the initial condition of field. alpha is a vector or tuple with two elements corresponding to the diffusion coefficient (multiplied by grid spacing) for the two dimensions. Niter is the number of iterations. On output, f is the field after Niter iterations.

source

Bathymetry

ROMS.gebco_loadFunction
x,y,b = gebco_load(bath_name,xr,yr)

Loads GEBCO bathymetry with lon and lat range xr, yr.

source

Boundary conditions and initial conditions

ROMS.CMEMSType
ds = ROMS.CMEMS(username,password,service_id,mapping,cachedir;
               motu_server = "http://nrt.cmems-du.eu/motu-web/Motu",
               # Put here the path of the script motuclient
               motu_program = "motuclient",
           )

Returns a structure ds to connect to a CMEMS Motu server using the python tools motuclient (which must be available in your PATH). The mapping parameter contains a dictorary linking the NetCDF CF standard namer to the underlying NetCDF variable names and the product identifers (more information is available in the product user manual). cachedir is a directory where the products are downloaded for caching.

Example

The values of service_id and mapping below are specific to the Mediterranean Sea and must be adapted for other domains.

cmems_username = "Alice"
cmems_password = "rabbit"
outdir = "/tmp"
service_id = "MEDSEA_ANALYSISFORECAST_PHY_006_013-TDS"
mapping = Dict(
    # var  product_id
    :sea_surface_height_above_geoid => ("zos","med-cmcc-ssh-an-fc-d"),
    :sea_water_potential_temperature => ("thetao", "med-cmcc-tem-an-fc-d"),
    :sea_water_salinity => ("so","med-cmcc-sal-an-fc-d"),
    :eastward_sea_water_velocity => ("uo", "med-cmcc-cur-an-fc-d"),
    :northward_sea_water_velocity => ("vo", "med-cmcc-cur-an-fc-d"),
)
dataset = ROMS.CMEMS(cmems_username,cmems_password,service_id,mapping,outdir)
source
ROMS.loadFunction
v,(x,y,z,t) = ROMS.load(ds::CMEMS,name::Symbol; kwargs...)

Loads a variable from a CMEMS remote resource. name is the NetCDF CF standard name.

source
ROMS.interp_climFunction
ROMS.interp_clim(domain,clim_filename,dataset,timerange;
                 padding = 0.5,
                 missing_value = -9999.)

Interpolate dataset on the the model grid domain and creating the climatology file clim_filename for all dates between timerange[1] and timerange[2].

source
ROMS.extract_icFunction
ROMS.extract_ic(domain,clim_filename,icfile,t0::DateTime;
                time_origin = DateTime(1858,11,17))

From the climatology clim_filename extract a single time instance at the time t0 (or the nearest) and save the result into icfile.

source
ROMS.extract_bcFunction
ROMS.extract_bc(domain,clim_filename,bc_filename; missing_value = 9999)

From the climatology clim_filename extract the boundary conditions at all open boundaries (using the mask in domain) and save the result into the netCDF nc_filename.

source

Atmospheric forcings

ROMS.prepare_ecmwfFunction
ROMS.prepare_ecmwf(
   atmo_fname,Vnames,filename_prefix,domain_name;
   time_origin = DateTime(1858,11,17),
   reset_accumulation = time -> Time(time) in (Time(0,0),Time(12,0)),
)

Generate ROMS forcing fields from the ECMWF data file atmo_fname. Note that some variables are accumulated. Per default, the accumulation is reset at 00:00 and 12:00 UTC. The accumulation period is determined from the time resolution (usually 3 hours).

Note ERA5 reanalysis (hourly data): Accumulations are performed over the hour.

Example

datadir = "..."
atmo_fname = joinpath(datadir,"ecmwf_sample_data.nc")
filename_prefix = joinpath(datadir,"liguriansea_")
domain_name = "Ligurian Sea Region"
Vnames = ["sustr","svstr","shflux","swflux","swrad","Uwind","Vwind","lwrad",
    "lwrad_down","latent","sensible","cloud","rain","Pair","Tair","Qair"]
ROMS.prepare_ecmwf(atmo_fname,Vnames,filename_prefix,domain_name)
)

Based on forcing/d_ecmwf2roms.m (svn revision 1102):

Copyright (c) 2002-2017 The ROMS/TOMS Group      Hernan G. Arango
Licensed under a MIT/X style license             John Wilkin
See License_ROMS.txt
source
ROMS.gfs_urlFunction
url = ROMS.gfs_url(time,tau;
         modelname = "gfs",
         resolution = 0.25,
         baseurl = "https://rda.ucar.edu/thredds/dodsC/files/g/ds084.1/")

Returns the OPeNDAP url for the GFS data at time time (DateTime) and the forecast time tau (hours) from the archive specified at baseurl.

source
ROMS.download_gfsFunction
atmo_src = ROMS.download_gfs(xr,yr,tr,cachedir)

Downloads GFS 0.25° model results from the NCAR Research Data Archive within the longitude range xr, latitude range yr and time range tr. Ranges are list of two elements with the start and end value. Results are saved in cachedir.

See ROMS.prepage_gfs for an example.

source
ROMS.prepare_gfsFunction
ROMS.prepare_gfs(
   atmo_src,Vnames,filename_prefix,domain_name;
   time_origin = DateTime(1858,11,17),
)

Generate ROMS forcing fields from GFS data atmo_src (a generated by ROMS.download_gfs). The other arguments are the same as for ROMS.prepage_ecmwf. The example below shows all currently supported values for Vnames.

Example

tr = (DateTime(2019,1,1),DateTime(2019,1,7))
xr = (7.5, 12.375)
yr = (41.875, 44.625)

cachedir = expanduser("~/tmp/GFS")
atmo_src = ROMS.download_gfs(xr,yr,tr,cachedir)

filename_prefix = "liguriansea_"
domain_name = "Ligurian Sea Region"
Vnames = ["sustr","svstr","swflux","swrad","Uwind","Vwind",
          "sensible","cloud","rain","Pair","Tair","Qair"]
ROMS.prepare_gfs(atmo_src,Vnames,filename_prefix,domain_name)
)
source

Thermodynamics

ROMS.vapor_pressureFunction
e = ROMS.vapor_pressure(T)

actual vapor pressure in hPa (millibars) from dewpoint temperature T in degree Celsius using using [1]. If T is the air temperature, then e is the saturated vapor pressure over liquid water is given by:

$e(T) = 6.11 \cdot 10 ^ {\left( \frac{7.5 T}{237.7 + T} \right)}$

[1] https://web.archive.org/web/20200926200733/https://www.weather.gov/media/epz/wxcalc/vaporPressure.pdf

source
ROMS.vapor_pressure_BuckFunction
e = ROMS.vapor_pressure_Buck(T)

actual vapor pressure in hPa (millibars) from dewpoint temperature T in degree Celsius using using Buck (1996). If T is the air temperature, then e is the saturated vapor pressure over liquid water is given by:

$e(T) = 6.1121 \exp \left(\left( 18.678 - \frac{T} {234.5}\right)\left( \frac{T} {257.14 + T} \right)\right)$

source
ROMS.relative_humidityFunction
rh = ROMS.relative_humidity(temperature_2m_C,dew_temperature_2m_C)

Compute the relative humidity (between 0 and 100) from temperature at 2 m, and dew_temperature at 2 m) both in degree Celsius)

[1] https://web.archive.org/web/20200926200733/https://www.weather.gov/media/epz/wxcalc/vaporPressure.pdf

source