Preparation of the input files for ROMS

The code here is also available as a notebook 02_prep_roms.ipynb.

ROMS needs several input files in the NetCDF format, in paricular:

  • the model grid
  • the initial conditions
  • the boundary conditions
  • the atmospheric forcing fields

Optionally

  • the climatology file
  • the field defining the nudging strength

This script can use multiple threads if julia was started with multi-threading (option -t/--threads or the environement variable JULIA_NUM_THREADS)

using Dates
using ROMS
using ROMS: whenopen
using Downloads: download

The model bathymetry

While the full GEBCO bathymetry is relatively large, where use here a subset of the global bathymetry to reduce the downloading time. (longitude from 5°E to 15°E and latitude from 40°N to 45°N)

bath_name = expanduser("~/Data/Bathymetry/gebco_30sec_1_ligurian_sea.nc")

if !isfile(bath_name)
    mkpath(dirname(bath_name))
    download("https://dox.ulg.ac.be/index.php/s/piwSaFP3nhM8jSD/download",bath_name)
end;

The time range for the simulation:

  • t0 start time
  • t1 end time
t0 = DateTime(2023,1,1);
t1 = DateTime(2023,1,4);

Define the bounding box the of the grid

# range of longitude
xr = [7.6, 12.2];

# range of latitude
yr = [42, 44.5];

# reduce bathymetry in x and y direction
red = (4, 4)

# maximum normalized topographic variations
rmax = 0.4;

# minimal depth
hmin = 2; # m

# name of folders and files
modeldir = expanduser("~/ROMS-implementation-test")

# The model grid (`GRDNAME` in roms.in)
grd_name = joinpath(modeldir,"roms_grd_liguriansea.nc");

Some parameters specific to the vertical coordinate system

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,
)
(Tcline = 50, theta_s = 5, theta_b = 0.4, nlevels = 32, Vtransform = 2, Vstretching = 4)

Create the model directory and generate the model grid

mkpath(modeldir);

domain = ROMS.generate_grid(grd_name,bath_name,xr,yr,red,opt,hmin,rmax);

@info "domain size $(size(domain.mask))"
 Smooth the topography...
  36 iterations - rmax = 0.39828464368061617
 Smooth the topography a bit more...
  hmin = 2.0
Vtransform  = 2   ROMS-UCLA
   igrid     = 1 at horizontal RHO-points

Vstretching = 4   Shchepetkin (2010)
   kgrid    = 0   at vertical RHO-points
   theta_s  = 5
   theta_b  = 0.4
   hc       = 50

 S-coordinate curves: k, s(k), C(k)

     32    -1.562500000000e-02    -5.060164946710e-05
     31    -4.687500000000e-02    -4.572401105559e-04
     30    -7.812500000000e-02    -1.280298055894e-03
     29    -1.093750000000e-01    -2.539560500182e-03
     28    -1.406250000000e-01    -4.265266261424e-03
     27    -1.718750000000e-01    -6.498792311199e-03
     26    -2.031250000000e-01    -9.293584054181e-03
     25    -2.343750000000e-01    -1.271634639800e-02
     24    -2.656250000000e-01    -1.684851356914e-02
     23    -2.968750000000e-01    -2.178801812466e-02
     22    -3.281250000000e-01    -2.765138117566e-02
     21    -3.593750000000e-01    -3.457614600596e-02
     20    -3.906250000000e-01    -4.272367537770e-02
     19    -4.218750000000e-01    -5.228232795041e-02
     18    -4.531250000000e-01    -6.347102015626e-02
     17    -4.843750000000e-01    -7.654316490419e-02
     16    -5.156250000000e-01    -9.179095542844e-02
     15    -5.468750000000e-01    -1.095499286006e-01
     14    -5.781250000000e-01    -1.302036934689e-01
     13    -6.093750000000e-01    -1.541886431904e-01
     12    -6.406250000000e-01    -1.819983765227e-01
     11    -6.718750000000e-01    -2.141874325067e-01
     10    -7.031250000000e-01    -2.513737824031e-01
      9    -7.343750000000e-01    -2.942393202666e-01
      8    -7.656250000000e-01    -3.435273436662e-01
      7    -7.968750000000e-01    -4.000357194613e-01
      6    -8.281250000000e-01    -4.646040954176e-01
      5    -8.593750000000e-01    -5.380931709292e-01
      4    -8.906250000000e-01    -6.213537270252e-01
      3    -9.218750000000e-01    -7.151829201676e-01
      2    -9.531250000000e-01    -8.202653975901e-01
      1    -9.843750000000e-01    -9.370972868551e-01

Vtransform  = 2   ROMS-UCLA
   igrid     = 5 at horizontal RHO-points

Vstretching = 4   Shchepetkin (2010)
   kgrid    = 1   at vertical W-points
   theta_s  = 5
   theta_b  = 0.4
   hc       = 50

 S-coordinate curves: k, s(k), C(k)

     32     0.000000000000e+00     0.000000000000e+00
     31    -3.125000000000e-02    -2.027105200195e-04
     30    -6.250000000000e-02    -8.157187016713e-04
     29    -9.375000000000e-02    -1.853765651128e-03
     28    -1.250000000000e-01    -3.341792624089e-03
     27    -1.562500000000e-01    -5.315506910025e-03
     26    -1.875000000000e-01    -7.822187499588e-03
     25    -2.187500000000e-01    -1.092174369724e-02
     24    -2.500000000000e-01    -1.468804314646e-02
     23    -2.812500000000e-01    -1.921052856298e-02
     22    -3.125000000000e-01    -2.459614455050e-02
     21    -3.437500000000e-01    -3.097159680944e-02
     20    -3.750000000000e-01    -3.848596528424e-02
     19    -4.062500000000e-01    -4.731368954908e-02
     18    -4.375000000000e-01    -5.765793793687e-02
     17    -4.687500000000e-01    -6.975436012900e-02
     16    -5.000000000000e-01    -8.387520422218e-02
     15    -5.312500000000e-01    -1.003337511642e-01
     14    -5.625000000000e-01    -1.194889786791e-01
     13    -5.937500000000e-01    -1.417503093415e-01
     12    -6.250000000000e-01    -1.675822183846e-01
     11    -6.562500000000e-01    -1.975083703627e-01
     10    -6.875000000000e-01    -2.321148135475e-01
      9    -7.187500000000e-01    -2.720515804961e-01
      8    -7.500000000000e-01    -3.180318172552e-01
      7    -7.812500000000e-01    -3.708272899397e-01
      6    -8.125000000000e-01    -4.312588001567e-01
      5    -8.437500000000e-01    -5.001796956603e-01
      4    -8.750000000000e-01    -5.784503244456e-01
      3    -9.062500000000e-01    -6.669010130264e-01
      2    -9.375000000000e-01    -7.662810584266e-01
      1    -9.687500000000e-01    -8.771914691872e-01
      0    -1.000000000000e+00    -1.000000000000e+00

[ Info: domain size (138, 75)

The boundary and initial conditions

# GCM interpolated on model grid (`CLMNAME` in roms.in)
clm_name =  joinpath(modeldir,"roms_clm_2023.nc")

# initial conditions (`ININAME` in roms.in)
ini_name =  joinpath(modeldir,"roms_ini_2023.nc")

# boundary conditions (`BRYNAME` in roms.in)
bry_name =  joinpath(modeldir,"roms_bry_2023.nc")

# temporary directory of the OGCM data
outdir = joinpath(modeldir,"OGCM")
mkpath(outdir)
"/home/runner/ROMS-implementation-test/OGCM"
  • For CMEMS boundary conditions https://marine.copernicus.eu/:
    • You may need to adapt the CMEMS product_id and mapping (if the model domain is outside of the Mediterranean Sea)
    • Data will be downloaded and saved in NetCDF by "chunks" of 60 days in the folder OGCM under the content of the variable basedir
    • You need to remove the files in this directory if you rerun the script with a different time range.

Here we use the following dataset: https://doi.org/10.25423/CMCC/MEDSEAMULTIYEARPHY006004_E3R1

product_id = "MEDSEA_MULTIYEAR_PHY_006_004"
"MEDSEA_MULTIYEAR_PHY_006_004"

mapping the variable (CF names) with the CMEMS dataset_id

mapping = Dict(
    :sea_surface_height_above_geoid => "med-cmcc-ssh-rean-d",
    :sea_water_potential_temperature => "med-cmcc-tem-rean-d",
    :sea_water_salinity => "med-cmcc-sal-rean-d",
    :eastward_sea_water_velocity => "med-cmcc-cur-rean-d",
    :northward_sea_water_velocity => "med-cmcc-cur-rean-d",
)

dataset = ROMS.CMEMS_zarr(product_id,mapping,outdir, time_shift = 12*60*60)
ROMS.CDMDataset{ZarrDatasets.ZarrDataset}(DataStructures.DefaultDict(:sea_water_salinity => "https://s3.waw3-1.cloudferro.com/mdl-arco-time-035/arco/MEDSEA_MULTIYEAR_PHY_006_004/med-cmcc-sal-rean-d_202012/timeChunked.zarr", :eastward_sea_water_velocity => "https://s3.waw3-1.cloudferro.com/mdl-arco-time-035/arco/MEDSEA_MULTIYEAR_PHY_006_004/med-cmcc-cur-rean-d_202012/timeChunked.zarr", :sea_surface_height_above_geoid => "https://s3.waw3-1.cloudferro.com/mdl-arco-time-035/arco/MEDSEA_MULTIYEAR_PHY_006_004/med-cmcc-ssh-rean-d_202012/timeChunked.zarr", :sea_water_potential_temperature => "https://s3.waw3-1.cloudferro.com/mdl-arco-time-036/arco/MEDSEA_MULTIYEAR_PHY_006_004/med-cmcc-tem-rean-d_202012/timeChunked.zarr", :northward_sea_water_velocity => "https://s3.waw3-1.cloudferro.com/mdl-arco-time-035/arco/MEDSEA_MULTIYEAR_PHY_006_004/med-cmcc-cur-rean-d_202012/timeChunked.zarr"), "/home/runner/ROMS-implementation-test/OGCM", Dict{Symbol, String}(), 60, Dict{Symbol, Any}(:_omitcode => [404, 403]), 43200)

Extent the time range by one extra day

ROMS.interp_clim(domain,clm_name,dataset,[t0-Dates.Day(1), t1+Dates.Day(1)])

ROMS.extract_ic(domain,clm_name,ini_name, t0);
ROMS.extract_bc(domain,clm_name,bry_name)
[ Info: download sea_surface_height_above_geoid in /home/runner/ROMS-implementation-test/OGCM/bc19f9716470e5766f561323c8885cbed77a46447b452f176394ca77e0c65331.nc (13149:13153)
[ Info: download sea_water_potential_temperature in /home/runner/ROMS-implementation-test/OGCM/278703240e4f8873236513a37b3b39c6ccfdbd87ee127f6dddf7188cedc630b2.nc (13149:13153)
[ Info: download sea_water_salinity in /home/runner/ROMS-implementation-test/OGCM/d35e5ad1eb56014e42eb24ad006e38931bab46134e7c4939caeeb9baa73d2a12.nc (13149:13153)
[ Info: download eastward_sea_water_velocity in /home/runner/ROMS-implementation-test/OGCM/ab64aa3534519da56071be2caf1ecd6c4df88e5298dc1494bae6e2014f77881e.nc (13149:13153)
[ Info: download northward_sea_water_velocity in /home/runner/ROMS-implementation-test/OGCM/0d19493ccd67b7f829dceaa66c2300299a80275ef965f5268e7e03aac0ab3b1a.nc (13149:13153)
[ Info: load 2022-12-31T12:00:00
[ Info: interpolate 2022-12-31T12:00:00
[ Info: load 2023-01-01T12:00:00
[ Info: interpolate 2023-01-01T12:00:00
[ Info: load 2023-01-02T12:00:00
[ Info: interpolate 2023-01-02T12:00:00
[ Info: load 2023-01-03T12:00:00
[ Info: interpolate 2023-01-03T12:00:00
[ Info: load 2023-01-04T12:00:00
[ Info: interpolate 2023-01-04T12:00:00
[ Info: Request to initalize from 2023-01-01T00:00:00
[ Info: Try to initalize from 2022-12-31T12:00:00
copy zeta
copy ubar
copy vbar
copy u
copy v
copy temp
copy salt
[ Info: DSTART = 59945.0

Nudging coefficients (NUDNAME)

tscale = 7; # days
alpha = 0.3;
halo = 2;
Niter = 50
max_tscale = 5e5

nud_name = joinpath(modeldir,"roms_nud_$(tscale)_$(Niter).nc")
tracer_NudgeCoef = ROMS.nudgecoef(domain,nud_name,alpha,Niter,
          halo,tscale; max_tscale = max_tscale);
[ Info: Number of sea grid cell:              5497
[ Info: Number sea grid cell with nudging:    2807
[ Info: Number sea grid cell without nudging: 2690

The atmospheric forcings

Prepare atmospheric forcings (FRCNAME)

ecmwf_fname = expanduser("~/Data/Atmosphere/ecmwf_operational_archive_2022-12-01_2024-02-01.nc")

if !isfile(ecmwf_fname)
    mkpath(dirname(ecmwf_fname))
    download("https://data-assimilation.net/upload/OCEA0036/ecmwf_operational_archive_2022-12-01_2024-02-01.nc",ecmwf_fname)
end

frc_name_prefix = joinpath(modeldir,"roms_frc_2023_")
domain_name = "Ligurian Sea Region"
Vnames = ["sustr","svstr","shflux","swflux","swrad","Uwind","Vwind",
    "lwrad","lwrad_down","latent","sensible","cloud","rain","Pair","Tair","Qair"]

# forcing_filenames corresponds to `FRCNAME` in roms.in
forcing_filenames = ROMS.prepare_ecmwf(ecmwf_fname,Vnames,frc_name_prefix,domain_name)
16-element Vector{Tuple{String, String}}:
 ("sustr", "/home/runner/ROMS-implementation-test/roms_frc_2023_sms.nc")
 ("svstr", "/home/runner/ROMS-implementation-test/roms_frc_2023_sms.nc")
 ("shflux", "/home/runner/ROMS-implementation-test/roms_frc_2023_shflux.nc")
 ("swflux", "/home/runner/ROMS-implementation-test/roms_frc_2023_swflux.nc")
 ("swrad", "/home/runner/ROMS-implementation-test/roms_frc_2023_swrad.nc")
 ("Uwind", "/home/runner/ROMS-implementation-test/roms_frc_2023_wind.nc")
 ("Vwind", "/home/runner/ROMS-implementation-test/roms_frc_2023_wind.nc")
 ("lwrad", "/home/runner/ROMS-implementation-test/roms_frc_2023_lwrad.nc")
 ("lwrad_down", "/home/runner/ROMS-implementation-test/roms_frc_2023_lwrad.nc")
 ("latent", "/home/runner/ROMS-implementation-test/roms_frc_2023_latent.nc")
 ("sensible", "/home/runner/ROMS-implementation-test/roms_frc_2023_sensible.nc")
 ("cloud", "/home/runner/ROMS-implementation-test/roms_frc_2023_cloud.nc")
 ("rain", "/home/runner/ROMS-implementation-test/roms_frc_2023_rain.nc")
 ("Pair", "/home/runner/ROMS-implementation-test/roms_frc_2023_Pair.nc")
 ("Tair", "/home/runner/ROMS-implementation-test/roms_frc_2023_Tair.nc")
 ("Qair", "/home/runner/ROMS-implementation-test/roms_frc_2023_Qair.nc")

We print a list of all generated files.

fn(name) = basename(name) # use relative file path
# fn(name) = name         # use absolute file path

println()
println("The created netCDF files are in $modeldir.");
println("The following information has to be added to roms.in. A template of this file is")
println("provided in the directory User/External of your ROMS source code")
println("You can also use relative or absolute file names.")
println()
println("! grid file ")
println("     GRDNAME == $(fn(grd_name))")
println()
println("! initial conditions")
println("     ININAME == $(fn(ini_name))")
println()
println("! boundary conditions")
println("     NBCFILES == 1")
println("     BRYNAME == $(fn(bry_name))")
println()
println("! climatology or large-scale circulatio model")
println("     NCLMFILES == 1")
println("     CLMNAME == $(fn(clm_name))")
println()
println("! nudging coefficients file (optional)")
println("     NUDNAME == $(fn(nud_name))")
println()
println("! forcing files")
println("     NFFILES == $(length(Vnames))")

for i in 1:length(Vnames)
    if i == 1
        print("     FRCNAME == ")
    else
        print("                ")
    end
    print("$(fn(frc_name_prefix))$(Vnames[i]).nc")
    if i < length(Vnames)
        print(" \\")
    end
    println()
end

The created netCDF files are in /home/runner/ROMS-implementation-test.
The following information has to be added to roms.in. A template of this file is
provided in the directory User/External of your ROMS source code
You can also use relative or absolute file names.

! grid file 
     GRDNAME == roms_grd_liguriansea.nc

! initial conditions
     ININAME == roms_ini_2023.nc

! boundary conditions
     NBCFILES == 1
     BRYNAME == roms_bry_2023.nc

! climatology or large-scale circulatio model
     NCLMFILES == 1
     CLMNAME == roms_clm_2023.nc

! nudging coefficients file (optional)
     NUDNAME == roms_nud_7_50.nc

! forcing files
     NFFILES == 16
     FRCNAME == roms_frc_2023_sustr.nc \
                roms_frc_2023_svstr.nc \
                roms_frc_2023_shflux.nc \
                roms_frc_2023_swflux.nc \
                roms_frc_2023_swrad.nc \
                roms_frc_2023_Uwind.nc \
                roms_frc_2023_Vwind.nc \
                roms_frc_2023_lwrad.nc \
                roms_frc_2023_lwrad_down.nc \
                roms_frc_2023_latent.nc \
                roms_frc_2023_sensible.nc \
                roms_frc_2023_cloud.nc \
                roms_frc_2023_rain.nc \
                roms_frc_2023_Pair.nc \
                roms_frc_2023_Tair.nc \
                roms_frc_2023_Qair.nc

Check the resulting files such as bathymetry, initial conditions, boundary conditions, interpolated model (clm_name file) and visualizing them.

Configuration files

Beside the created NetCDF files, ROMS needs two configuration files (roms.in and varinfo.yaml)

romsdir = expanduser("~/src/roms")
modeldir = expanduser("~/ROMS-implementation-test")
simulationdir = joinpath(modeldir,"Simulation1")
mkpath(simulationdir)

frc_name = joinpath.(modeldir,sort(filter(startswith("roms_frc"),readdir(modeldir))));

Copy varinfo.yaml from ~/src/roms/ROMS/External/varinfo.yaml in your directory for your simulation (e.g. ROMS-implementation-test). This file does not need to be changed.

var_name_template = joinpath(romsdir,"ROMS","External","varinfo.yaml")
var_name = joinpath(simulationdir,"varinfo.yaml")
cp(var_name_template,var_name; force=true);

Load the ROMS grid

domain = ROMS.Grid(grd_name);

We use roms.in from ~/src/roms/User/External/roms.in as a template

intemplate = joinpath(romsdir,"User","External","roms.in")
infile = joinpath(simulationdir,"roms.in");

This file is typicall edited with a text editor (when editing this file, do not use "tabs".). Check the glossary at the end of this file for the meaning of the keys that we will change.

Here we edit the file programmatically. These are the changes that are done in the following:

  • adapt MyAppCPP and change it to LIGURIANSEA

  • adapt file names VARNAME, GRDNAME, ININAME, BRYNAME, CLMNAME, FRCNAME and NFFILES (varinfo.yaml, LS2v.nc, ic2019.nc, bc2019.nc, clim2019.nc, liguriansea2019_*.nc, * means the different variables). NFFILES is the number of forcing files.

  • change Lm, Mm and N based on the dimensions of your grid (make sure to read the glossary for these variable in roms.in)

  • read the desciption about "lateral boundary conditions" and adapt boundaries LBC:

    • use closed (Clo) for boundaries without sea-point
    • for open boundaries:
      • free-surface: Chapman implicit (Cha)
      • 2D U/V-momentum: Flather (Fla)
      • 3D U/V-momentum, temperature, salinity: Radiation with nudging (RadNud)
      • mixing TKE: Radiation (Rad)
  • set the starting time and time reference

DSTART = ...
TIME_REF =  18581117

where DSTART is here the number of days since 1858-11-17 or November 17, 1858 (see also modified Julia day) of the start of the model simulation (t0 in the julia script). For instance the number of days since 2014-01-01 (year-month-day) can be computed by of following commands in Julia:

using Dates
Date(2020,1,1) - Date(1858,11,17)

The inverse operation can be done with:

using Dates
Date(1858,11,17) + Day(58849)

You can use DateTime if you want to specify hour, minutes or seconds.

  • Adapt the length of a time step DT (in seconds) and number of time steps NTIMES
  • DT can be 300 seconds
  • Initially we choose:
    • NTIMES -> number of time step corresponding to 2 days (e.g. 2*24*60*60/DT where DT is the time steps in seconds)
    • NHIS, NAVG-> number of time steps corresponding to 1 hour
    • NRST -> number of time steps correspond to 1 hour
# time step (seconds)
DT = 300.
# output frequency of ROMS in time steps
NHIS = round(Int,24*60*60 / DT)
NRST = NAVG = NHIS

# number of time steps
t0 = DateTime(2023,1,1);
t1 = DateTime(2023,1,4);
NTIMES = floor(Int,Dates.value(t1-t0) / (DT * 1000))
864

How many CPU cores does your machine have? You can use the command top in a shell terminal followed by 1. The number of CPU cores should be NtileI * NtileJ. The parameters NtileI and NtileJ are defined in roms.in.

NtileI = 1
NtileJ = 1

substitutions = Dict(
    "TITLE" => "My test",
    "NtileI" => NtileI,
    "NtileJ" => NtileJ,
    "TIME_REF" => "18581117",
    "VARNAME" => var_name,
    "GRDNAME" => grd_name,
    "ININAME" => ini_name,
    "BRYNAME" => bry_name,
    "CLMNAME" => clm_name,
    "NFFILES" => length(frc_name),
    "FRCNAME" => join(frc_name,"  \\\n       "),
    "Vtransform" => domain.Vtransform,
    "Vstretching" => domain.Vstretching,
    "THETA_S" => domain.theta_s,
    "THETA_B" => domain.theta_b,
    "TCLINE" => domain.Tcline,
    "Lm" => size(domain.h,1)-2,
    "Mm" => size(domain.h,2)-2,
    "N" => domain.nlevels,
    "LBC(isFsur)" => whenopen(domain,"Cha"),
    "LBC(isUbar)" => whenopen(domain,"Fla"),
    "LBC(isVbar)" => whenopen(domain,"Fla"),
    "LBC(isUvel)" => whenopen(domain,"RadNud"),
    "LBC(isVvel)" => whenopen(domain,"RadNud"),
    "LBC(isMtke)" => whenopen(domain,"Rad"),
    "LBC(isTvar)" => whenopen(domain,"RadNud") * " \\\n" * whenopen(domain,"RadNud"),
    "DT" => DT,
    "NHIS" => NHIS,
    "NAVG" => NAVG,
    "NRST" => NRST,
    "NTIMES" => NTIMES,
    "NUDNAME" => nud_name,
    "TNUDG" => "10.0d0 10.0d0",
    "LtracerCLM" => "T T",
    "LnudgeTCLM" => "T T",
    "OBCFAC" => 10.0,
)

ROMS.infilereplace(intemplate,infile,substitutions)
[ Info: setting TITLE to: My test
[ Info: setting VARNAME to: /home/runner/ROMS-implementation-test/Simulation1/varinfo.yaml
[ Info: setting Lm to: 136
[ Info: setting Mm to: 73
[ Info: setting N to: 32
[ Info: setting NtileI to: 1
[ Info: setting NtileJ to: 1
[ Info: setting LBC(isFsur) to: Cha Cha Clo Clo
[ Info: setting LBC(isUbar) to: Fla Fla Clo Clo
[ Info: setting LBC(isVbar) to: Fla Fla Clo Clo
[ Info: setting LBC(isUvel) to: RadNud RadNud Clo Clo
[ Info: setting LBC(isVvel) to: RadNud RadNud Clo Clo
[ Info: setting LBC(isMtke) to: Rad Rad Clo Clo
┌ Info: setting LBC(isTvar) to: RadNud RadNud Clo Clo \
└ RadNud RadNud Clo Clo
[ Info: setting NTIMES to: 864
[ Info: setting DT to: 300.0
[ Info: setting NRST to: 288
[ Info: setting NHIS to: 288
[ Info: setting NAVG to: 288
[ Info: setting Vtransform to: 2
[ Info: setting Vstretching to: 4
[ Info: setting THETA_S to: 5.0
[ Info: setting THETA_B to: 0.4
[ Info: setting TCLINE to: 50.0
[ Info: setting TIME_REF to: 18581117
[ Info: setting TNUDG to: 10.0d0 10.0d0
[ Info: setting OBCFAC to: 10.0
[ Info: setting LtracerCLM to: T T
[ Info: setting LnudgeTCLM to: T T
[ Info: setting GRDNAME to: /home/runner/ROMS-implementation-test/roms_grd_liguriansea.nc
[ Info: setting ININAME to: /home/runner/ROMS-implementation-test/roms_ini_2023.nc
[ Info: setting BRYNAME to: /home/runner/ROMS-implementation-test/roms_bry_2023.nc
[ Info: setting CLMNAME to: /home/runner/ROMS-implementation-test/roms_clm_2023.nc
[ Info: setting NUDNAME to: /home/runner/ROMS-implementation-test/roms_nud_7_50.nc
[ Info: setting NFFILES to: 13
┌ Info: setting FRCNAME to: /home/runner/ROMS-implementation-test/roms_frc_2023_Pair.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_Qair.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_Tair.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_cloud.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_latent.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_lwrad.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_rain.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_sensible.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_shflux.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_sms.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_swflux.nc  \
│        /home/runner/ROMS-implementation-test/roms_frc_2023_swrad.nc  \
└        /home/runner/ROMS-implementation-test/roms_frc_2023_wind.nc
68c68
<        TITLE = Wind-Driven Upwelling/Downwelling over a Periodic Channel
---
>        TITLE = My test
77c77
<      VARNAME = ROMS/External/varinfo.yaml
---
>      VARNAME = /home/runner/ROMS-implementation-test/Simulation1/varinfo.yaml
95,97c95,97
<           Lm == 41            ! Number of I-direction INTERIOR RHO-points
<           Mm == 80            ! Number of J-direction INTERIOR RHO-points
<            N == 16            ! Number of vertical levels
---
>           Lm == 136            ! Number of I-direction INTERIOR RHO-points
>           Mm == 73            ! Number of J-direction INTERIOR RHO-points
>            N == 32            ! Number of vertical levels
185,190c185,190
<    LBC(isFsur) ==   Per     Clo     Per     Clo         ! free-surface
<    LBC(isUbar) ==   Per     Clo     Per     Clo         ! 2D U-momentum
<    LBC(isVbar) ==   Per     Clo     Per     Clo         ! 2D V-momentum
<    LBC(isUvel) ==   Per     Clo     Per     Clo         ! 3D U-momentum
<    LBC(isVvel) ==   Per     Clo     Per     Clo         ! 3D V-momentum
<    LBC(isMtke) ==   Per     Clo     Per     Clo         ! mixing TKE
---
>    LBC(isFsur) == Cha Cha Clo Clo         ! free-surface
>    LBC(isUbar) == Fla Fla Clo Clo         ! 2D U-momentum
>    LBC(isVbar) == Fla Fla Clo Clo         ! 2D V-momentum
>    LBC(isUvel) == RadNud RadNud Clo Clo         ! 3D U-momentum
>    LBC(isVvel) == RadNud RadNud Clo Clo         ! 3D V-momentum
>    LBC(isMtke) == Rad Rad Clo Clo         ! mixing TKE
192c192,193
<    LBC(isTvar) ==   Per     Clo     Per     Clo \       ! temperature
---
>    LBC(isTvar) == RadNud RadNud Clo Clo \
> RadNud RadNud Clo Clo       ! temperature
225,226c226,227
<       NTIMES == 1440
<           DT == 300.0d0
---
>       NTIMES == 864
>           DT == 300.0
263c264
<         NHIS == 72
---
>         NHIS == 288
268c269
<         NAVG == 72
---
>         NAVG == 288
405,407c406,408
<      THETA_S == 3.0d0                      ! surface stretching parameter
<      THETA_B == 0.0d0                      ! bottom  stretching parameter
<       TCLINE == 25.0d0                     ! critical depth (m)
---
>      THETA_S == 5.0                      ! surface stretching parameter
>      THETA_B == 0.4                      ! bottom  stretching parameter
>       TCLINE == 50.0                     ! critical depth (m)
425c426
<     TIME_REF =  0.0d0                      ! yyyymmdd.dd
---
>     TIME_REF = 18581117                      ! yyyymmdd.dd
430c431
<        TNUDG == 2*0.0d0                    ! days
---
>        TNUDG == 10.0d0 10.0d0                    ! days
439c440
<       OBCFAC == 0.0d0                      ! nondimensional
---
>       OBCFAC == 10.0                      ! nondimensional
473c474
<   LtracerCLM == F F                        ! temperature, salinity, inert
---
>   LtracerCLM == T T                        ! temperature, salinity, inert
483c484
<   LnudgeTCLM == F F                        ! temperature, salinity, inert
---
>   LnudgeTCLM == T T                        ! temperature, salinity, inert
964,965c965,966
<      GRDNAME == roms_grd.nc
<      ININAME == roms_ini.nc
---
>      GRDNAME == /home/runner/ROMS-implementation-test/roms_grd_liguriansea.nc
>      ININAME == /home/runner/ROMS-implementation-test/roms_ini_2023.nc
1007c1008
<      BRYNAME == roms_bry.nc
---
>      BRYNAME == /home/runner/ROMS-implementation-test/roms_bry_2023.nc
1019c1020
<      CLMNAME == roms_clm.nc
---
>      CLMNAME == /home/runner/ROMS-implementation-test/roms_clm_2023.nc
1023c1024
<      NUDNAME == roms_nud.nc
---
>      NUDNAME == /home/runner/ROMS-implementation-test/roms_nud_7_50.nc
1051c1052
<      NFFILES == 1                          ! number of unique forcing files
---
>      NFFILES == 13                          ! number of unique forcing files
1053c1054,1066
<      FRCNAME == roms_frc.nc                ! forcing file 1, grid 1
---
>      FRCNAME == /home/runner/ROMS-implementation-test/roms_frc_2023_Pair.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_Qair.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_Tair.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_cloud.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_latent.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_lwrad.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_rain.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_sensible.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_shflux.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_sms.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_swflux.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_swrad.nc  \
>        /home/runner/ROMS-implementation-test/roms_frc_2023_wind.nc                ! forcing file 1, grid 1

Always make make sure that THETA_S, THETA_B, TCLINE, Vtransform and Vstretching match the values in your julia script. We can review the changes with the shell command:

diff -u --color ~/src/roms/User/External/roms.in roms.in