Reference
Model grid
ROMS.create_grid
— FunctionROMS.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.
This function currently only works for non-rotated grids (angle = 0) and the spherical grids.
ROMS.smoothgrid
— Functionhsmooth = 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
ROMS.stiffness_ratios
— Functionrx0,rx1 = ROMS.stiffness_ratios(mask_u,mask_v,z_w)
rx0,rx1 = ROMS.stiffness_ratios(domain)
Compute maximum grid stiffness ratios following Beckmann and Haidvogel (rx0
) and Haney (rx1
).
ROMS.reduce_res
— Function field = reduce_res(field,red);
reduce resolution of field field
by a factor "red"
ROMS.stretching
— Functions,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-pointsreport
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
ROMS.Grid
— Typegrid = 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.
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)
ROMS.generate_grid
— FunctionROMS.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,
)
ROMS.nudgecoef
— Functionnudgecoef(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).
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.
Bathymetry
ROMS.gebco_load
— Functionx,y,b = gebco_load(bath_name,xr,yr)
Loads GEBCO bathymetry with lon and lat range xr, yr.
Boundary conditions and initial conditions
ROMS.CMEMS_zarr
— Functionds = ROMS.CMEMS_zarr(product_id,mapping,cachedir;
chunks = 60,
time_shift = 0,
kwargs...
)
Returns a structure ds
to connect to a CMEMS zarr server. 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.
While for most datasets (and CMEMS in the past) the time represents the central time the time axis. However since 2024, the time in the CMEMS Zarr data represents now the beginning of the time interval. Therefore time_shift
has to be added to the time variable to account for this difference. For example, if for a daily dataset, the first time instance is the average from 2000-01-01:00:00:00 to 2000-01-02:00:00:00, then the Zarr file records 2000-01-01:00:00:00
(the beginning for the averaging interval) rather than 2000-01-01:12:00:00
(the center for the averaging interval). In this case, time_shift
should be 12*60*60
(12 hours in seconds).
Example
The values of product_id
and mapping
(with dataset_id
) below are specific to the Mediterranean Sea and must be adapted for other domains.
outdir = "/tmp"
product_id = "MEDSEA_MULTIYEAR_PHY_006_004"
mapping = Dict(
# var dataset_id
:sea_surface_height_above_geoid => "med-cmcc-ssh-an-fc-d",
:sea_water_potential_temperature => "med-cmcc-tem-an-fc-d",
:sea_water_salinity => "med-cmcc-sal-an-fc-d",
:eastward_sea_water_velocity => "med-cmcc-cur-an-fc-d",
:northward_sea_water_velocity => "med-cmcc-cur-an-fc-d",
)
dataset = ROMS.CMEMS_zarr(product_id,mapping,outdir,time_shift = 12*60*60)
ROMS.load
— Functionv,(x,y,z,t) = ROMS.load(ds,name::Symbol; kwargs...)
Loads a variable from a remote resource ds
. name
is the NetCDF CF standard name.
ROMS.interp_clim
— FunctionROMS.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]
.
ROMS.extract_ic
— FunctionROMS.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
.
ROMS.extract_bc
— FunctionROMS.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
.
ROMS.whenopen
— FunctionLBC = ROMS.whenopen(domain::Grid,BC::AbstractString)
ROMS.whenopen
returns a list of the four lateral boundary condition (corresponding to "west","south","east" and "north"). The corresponding elements of LBC will be "Clo" for closed boundaries for the provided BC
for open boundary conditions (e.g. Cha, Fla, Rad, RadNud, as defined in the roms.in glossary).
Atmospheric forcings
ROMS.prepare_ecmwf
— FunctionROMS.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
ROMS.gfs_url
— Functionurl = 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
.
ROMS.download_gfs
— Functionatmo_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.
ROMS.prepare_gfs
— FunctionROMS.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)
)
Thermodynamics
ROMS.vapor_pressure
— Functione = 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
ROMS.vapor_pressure_Buck
— Functione = 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)$
ROMS.relative_humidity
— Functionrh = 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
ROMS.latent_heat_of_vaporization
— Functionλ = ROMS.latent_heat_of_vaporization(T)
Compute the latent heat of vaporization λ
(J/kg) of water at temperature T
(in degree Celsius).
The function implements equation 2.55 (page 38) of Foken, T, 2008: Micrometeorology. Springer, Berlin, Germany.
ROMS.build
— FunctionROMS.build(roms_application,modeldir;
stdout = stdout,
use_mpi = false,
use_openmp = false,
use_mpif90 = use_mpi,
use_netcdf4 = true,
my_header_dir = modeldir,
my_analytical_dir = modeldir,
fortran_compiler = "ifort",
clean = true,
jobs = 1,
bindir = modeldir,
extra_env = Dict(),
)
Compile the ROMS source code with the application name roms_application
and the ROMS project directory is the directory modeldir
which will contain the produced binaries.
See build_roms.sh for more information.
ROMS.run_model
— FunctionROMS.run_model(modeldir::AbstractString,romsin::AbstractString;
use_mpi = false,
use_openmp = false,
stdout = stdout,
mpiexec = "mpiexec",
np = 1)
Executes ROMS with the model directory modeldir
(containing the ROMS binary) and the input file romsin
using np
processes (or threads for OpenMP).