FJORD TIDAL CASE
Note: In SVN revision 933 (January 26, 2019) all ocean_*.in files were renamed to roms_*.in and all ocean* ROMS executables were renamed to roms* in order to facilitate and clarify model coupling efforts. More information can be found in the ROMS repository Trac ticket #794. If you are working with a ROMS release prior to revision 933 you will need to replace roms_upwelling.in with ocean_upwelling.in and romsS, romsM, or romsO with oceanS, oceanM, or oceanO in all commands below.
There are a few ways to setup tidal forcing in ROMS (...learn more about tides in ROMS). This tutorial will explain the simplest way, which is to prescribe (analytically) tidal variations in sea-surface height at a open boundary (this uses a FLATHER/CHAPMAN combination).
I learned this trick from hetland in this forum post (at the bottom).
RESTRICTIONS: This technique is ONLY applicable to embayments with a relatively small open boundary (say < 10 km). In this cases, the tidal height along the open boundary is essentially uniform and can be prescribed analytically.
PREREQUISITES: This tutorial assumes that you have (1) downloaded ROMS, (2) installed it in your computer (or cluster), and (3) tested it by compiling and running one of the included test cases. If you haven't done all the above, check the "Getting Started" and "Tutorials" WikiROMS sections.
WATCH HERE a YouTube Video with a simulation product of this tutorial. The plotted colormap is current velocity in the u direction.
Geographical Preamble
This application is for Ship Harbour, an estuarine fjord in Nova Scotia, Canada. Click here to see the location in Google. Tides are semidiurnal and tidal range is 1.4 m on average and 2 m on spring tides. Ship Harbour is a long embayment (~7 km) with an open (EAST) boundary of ~1 km, therefore it is an ideal candidate for the type of analytical tidal forcing taught in this tutorial. For now I will only include tidal forcing, however, there is a river at the uppermost end of the estuary, which discharges freshwater at an annual average rate of 18 m3 s-1. I plan to write another tutorial on how to add a river, but for now is only tides.
Grid Generation
The first step to set up a realistic application is to set up a realistic grid. There are several software packages to generate ROMS grids; SEAGRID and GRIDGEN being the most popular ones. I used the much-less-fancy EASYGRID, which is a bit easier to get up and running.
DOWNLOAD HERE the grid and initialization files required for this tutorial. Alternatively, you can create your own grid and initialization files for this tutorial using the grid-generation software of your choice. You can download the bathymetry and coastline data for Ship Harbour Fjord HERE.
NOTE: I also wrote a tutorial for grid generation using EASYGRID. The output of that tutorial are the grid and initialization files used in this tutorial (see above). So if you want to start from scratch, you may want to begin with that tutorial.
Compiling ROMS with tides
Before you compile ROMS, you need to create a header (.h) file with the appropriate cpp definitions that turn ON the analytical tides. Also, you need to modify some analytical Fortran files, where you are going to specify how to create the analytical tide.
Creating header (.h) file
Create a file named fjord.h and copy-paste the cpp definitions below:
/* ** ** Options for Tidal Fjord. ** ** Application flag: FJORD ** Input script: roms_fjord.in */ #define UV_ADV /* use to turn ON or OFF advection terms */ #define UV_COR /* use to turn ON or OFF Coriolis term */ #define UV_QDRAG /* use to turn ON or OFF quadratic bottom friction */ #define UV_VIS4 /* use to turn ON or OFF harmonic horizontal mixing */ #define MIX_S_UV /* momentum mixing on s-surfaces */ #define DJ_GRADPS /* use if splines density Jacobian (Shchepetkin, 2000) */ #define TS_U3HADVECTION /* use if 3rd-order upstream horiz. advection */ #define TS_C4VADVECTION /* use if 4th-order centered vertical advection */ #define SOLVE3D /* use if solving 3D primitive equations */ #define SPLINES /* use to activate parabolic splines reconstruction */ /* */ #define ANA_SMFLUX /* use if analytical surface momentum stress */ #define ANA_STFLUX /* use if analytical surface temperature flux */ #define ANA_SSFLUX /* use if analytical surface salinity flux */ #define ANA_BSFLUX /* use if analytical bottom salinity flux */ #define ANA_BTFLUX /* use if analytical bottom temperature flux */ #define MASKING /* use if analytical masking is enabled */ #define EAST_FSCHAPMAN /* use if free-surface Chapman condition*/ #define EAST_M2FLATHER /* use if 2D momentum Flather condition*/ #define EAST_M3RADIATION /* use if 3D momentum radiation condition*/ #define EAST_TRADIATION /* use if tracers radiation condition*/ #define ANA_FSOBC /* use if analytical free-surface boundary conditions*/ #define ANA_M2OBC /* use if analytical 2D momentum boundary conditions*/
The last 6 cpp definitions are the responsible for the analytical tidal forcing. I used the basin.h header file as a starting template... the first 10 cpp definitions are from it.
If your open boundary is not EAST... change EAST_FSCHAPMAN, EAST_M2FLATHER, EAST_M3RADIATION and EAST_TRADIATION to represent your open boundary (e.g. WEST__FSCHAPMAN for a west open boundary... etc.).
Modify ana_fsobc.h
You will have to edit the file ana_fsobc.h (located in trunk/ROMS/Functionals) to tell ROMS to estimate surface height analytically (at the open boundary) when running the FJORD case. Below is a snippet of the end of the file... you have to ADD the green code.
#elif defined TEST_CHAN
IF (WESTERN_EDGE) THEN
cff=0.0_r8
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_west(j)=cff
END DO
END IF
IF (EASTERN_EDGE) THEN
cff=-0.4040_r8*MIN(time(ng)/150000.0_r8,1.0_r8)
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_east(j)=cff
END DO
END IF
#elif defined WEDDELL
IF (WESTERN_EDGE) THEN
fac=TANH((tdays(ng)-dstart)/1.0_r8)
omega=2.0_r8*pi*time(ng)/(12.42_r8*3600.0_r8) ! M2 Tide period
val=0.53_r8+(0.53_r8-0.48_r8)/REAL(Iend+1,r8)
phase=(277.0_r8+(277.0_r8-240.0_r8)/REAL(Iend+1,r8))*deg2rad
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_west(j)=fac*val*COS(omega-phase)
END DO
END IF
IF (EASTERN_EDGE) THEN
fac=TANH((tdays(ng)-dstart)/1.0_r8)
omega=2.0_r8*pi*time(ng)/(12.42_r8*3600.0_r8) ! M2 Tide period
val=0.53_r8+(0.53_r8-0.48_r8)
phase=(277.0_r8+(277.0_r8-240.0_r8))*deg2rad
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_east(j)=fac*val*COS(omega-phase)
END DO
END IF
#elif defined FJORD
IF (EASTERN_EDGE) THEN
fac=TANH((tdays(ng)-dstart)/1.0_r8)
omega=2.0_r8*pi*time(ng)/(12.42_r8*3600.0_r8) ! M2 Tide period
val=0.53_r8+(0.53_r8-0.48_r8)
phase=(277.0_r8+(277.0_r8-240.0_r8))*deg2rad
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_east(j)=fac*val*COS(omega-phase)
END DO
END IF
#else
IF (EASTERN_EDGE) THEN
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_east(j)=0.0_r8
END DO
END IF
IF (WESTERN_EDGE) THEN
DO j=JstrR,JendR
BOUNDARY(ng)%zeta_west(j)=0.0_r8
END DO
END IF
IF (SOUTHERN_EDGE) THEN
DO i=IstrR,IendR
BOUNDARY(ng)%zeta_south(i)=0.0_r8
END DO
END IF
IF (NORTHERN_EDGE) THEN
DO i=IstrR,IendR
BOUNDARY(ng)%zeta_north(i)=0.0_r8
END DO
END IF
#endif
RETURN
END SUBROUTINE ana_fsobc_tile
Modify ana_m2obc.h
You will have to edit the file ana_m2obc.h (located in trunk/ROMS/Functionals) to tell ROMS to estimate currents analytically (at the open boundary) when running the FJORD case. Below is a snippet of the end of the file... you have to ADD the green code.
#elif defined WEDDELL
IF (WESTERN_EDGE) THEN
fac=TANH((tdays(ng)-dstart)/1.0_r8)
omega=2.0_r8*pi*time(ng)/(12.42_r8*3600.0_r8) ! M2 Tide period
minor=0.0143_r8+(0.0143_r8+0.010_r8)/REAL(Iend+1,r8)
major=0.1144_r8+(0.1144_r8-0.013_r8)/REAL(Iend+1,r8)
phase=(318.0_r8+(318.0_r8-355.0_r8)/REAL(Iend+1,r8))*deg2rad
angle=(125.0_r8+(125.0_r8- 25.0_r8)/REAL(Iend+1,r8))*deg2rad
DO j=JstrR,JendR
val=0.5_r8*(angler(Istr-1,j)+angler(Istr,j))
BOUNDARY(ng)%ubar_west(j)=fac*(major*COS(angle-val)* &
& COS(omega-phase)- &
& minor*SIN(angle-val)* &
& SIN(omega-phase))
END DO
DO j=Jstr,JendR
val=0.5_r8*(angler(Istr-1,j-1)+angler(Istr-1,j))
BOUNDARY(ng)%vbar_west(j)=fac*(major*SIN(angle-val)* &
& COS(omega-phase)- &
& minor*SIN(angle-val)* &
& COS(omega-phase))
END DO
END IF
IF (EASTERN_EDGE) THEN
fac=TANH((tdays(ng)-dstart)/1.0_r8)
omega=2.0_r8*pi*time(ng)/(12.42_r8*3600.0_r8) ! M2 Tide period
minor=0.0143_r8+(0.0143_r8+0.010_r8)
major=0.1144_r8+(0.1144_r8-0.013_r8)
phase=(318.0_r8+(318.0_r8-355.0_r8))*deg2rad
angle=(125.0_r8+(125.0_r8- 25.0_r8))*deg2rad
DO j=JstrR,JendR
val=0.5_r8*(angler(Iend,j)+angler(Iend+1,j))
BOUNDARY(ng)%ubar_east(j)=fac*(major*COS(angle-val)* &
& COS(omega-phase)- &
& minor*SIN(angle-val)* &
& SIN(omega-phase))
END DO
DO j=Jstr,JendR
val=0.5_r8*(angler(Iend+1,j-1)+angler(Iend+1,j))
BOUNDARY(ng)%vbar_east(j)=fac*(major*SIN(angle-val)* &
& COS(omega-phase)- &
& minor*SIN(angle-val)* &
& COS(omega-phase))
END DO
END IF
#elif defined FJORD
IF (EASTERN_EDGE) THEN
fac=TANH((tdays(ng)-dstart)/1.0_r8)
omega=2.0_r8*pi*time(ng)/(12.42_r8*3600.0_r8) ! M2 Tide period
minor=0.0143_r8+(0.0143_r8+0.010_r8)
major=0.1144_r8+(0.1144_r8-0.013_r8)
phase=(318.0_r8+(318.0_r8-355.0_r8))*deg2rad
angle=(125.0_r8+(125.0_r8- 25.0_r8))*deg2rad
DO j=JstrR,JendR
val=0.5_r8*(angler(Iend,j)+angler(Iend+1,j))
BOUNDARY(ng)%ubar_east(j)=fac*(major*COS(angle-val)* &
& COS(omega-phase)- &
& minor*SIN(angle-val)* &
& SIN(omega-phase))
END DO
DO j=Jstr,JendR
val=0.5_r8*(angler(Iend+1,j-1)+angler(Iend+1,j))
BOUNDARY(ng)%vbar_east(j)=fac*(major*SIN(angle-val)* &
& COS(omega-phase)- &
& minor*SIN(angle-val)* &
& COS(omega-phase))
END DO
END IF
#else
IF (EASTERN_EDGE) THEN
DO j=JstrR,JendR
BOUNDARY(ng)%ubar_east(j)=0.0_r8
END DO
DO j=Jstr,JendR
BOUNDARY(ng)%vbar_east(j)=0.0_r8
END DO
END IF
IF (WESTERN_EDGE) THEN
DO j=JstrR,JendR
BOUNDARY(ng)%ubar_west(j)=0.0_r8
END DO
DO j=Jstr,JendR
BOUNDARY(ng)%vbar_west(j)=0.0_r8
END DO
END IF
IF (SOUTHERN_EDGE) THEN
DO i=Istr,IendR
BOUNDARY(ng)%ubar_south(i)=0.0_r8
END DO
DO i=IstrR,IendR
BOUNDARY(ng)%vbar_south(i)=0.0_r8
END DO
END IF
IF (NORTHERN_EDGE) THEN
DO i=Istr,IendR
BOUNDARY(ng)%ubar_north(i)=0.0_r8
END DO
DO i=IstrR,IendR
BOUNDARY(ng)%vbar_north(i)=0.0_r8
END DO
END IF
#endif
RETURN
END SUBROUTINE ana_m2obc_tile
Now you are ready to compile ROMS!
Same as you did with the test cases... now you have to compile ROMS using the newly created fjord.h header file. If you are unsure of how to compile ROMS, you may want to take a look to the build.sh / build.bash page and to the example within the "contributions" section of that page.
Creating roms_fjord.in
Once ROMS has been compiled and you have your oceanS (or oceanM) executable file, you need to created an input file to run ROMS. You may want to make a copy of roms_basin.in (trunk/ROMS/External) and use it as a starting template.
Don't forget to rename the copy as roms_fjord.in
Below are snippets of roms_basin.in. In GREEN are the parts that you need to change in the copy (roms_fjord.in). Parts in RED also need to be changed in the copy, but the actual replaced text will depend on your own configuration.
Application Title & RHO-grid dimensions
! Application title.
!
TITLE = Tidal Fjord
!
! C-preprocessing Flag.
!
MyAppCPP = FJORD
!
! Input variable information file name. This file needs to be processed
! first so all information arrays can be initialized properly.
!
VARNAME = ../PATH2VARINFO/varinfo.dat
!
! Grid dimension parameters. See notes below in the Glossary for how to set
! these parameters correctly.
!
Lm == 90 ! Number of I-direction INTERIOR RHO-points
Mm == 24 ! Number of J-direction INTERIOR RHO-points
N == 10 ! Number of vertical levels
NOTE: The values of Lm, Mm and N are output to the screen when creating the grid file using EASYGRID.
Parallelizations parameters
If Running ROMS in parallel, use the configuration below... otherwise you may leave NtileI == 1 and NtileJ == 1
! Domain decomposition parameters for serial, distributed-memory or
! shared-memory configurations used to determine tile horizontal range
! indices (Istr,Iend) and (Jstr,Jend), [1:Ngrids].
!
NtileI == 4 ! I-direction partition
NtileJ == 1 ! J-direction partition
Time steps
! Time-Stepping parameters.
!
NTIMES == 1866241 ! 3 days
DT == 10
NDTFAST == 20
NOTE: The values of DT is also output to the screen when creating the grid file using EASYGRID.
Other parameters
These changes are just to speed up simulations.
! Number of eigenvalues (NEV) and eigenvectors (NCV) to compute for the
! Lanczos/Arnoldi problem in the Generalized Stability Theory (GST)
! analysis. NCV must be greater than NEV (see documentation below).
!
NEV = 2 ! Number of eigenvalues
NCV = 10 ! Number of eigenvectors
!
! Input/Output parameters.
!
NRREC == 0
LcycleRST == T
NRST == 360 ! Every 1 hour
NSTA == 360 ! Every 1 hour
NFLT == 360 ! Every 1 hour
NINFO == 360 ! Every 1 hour
!
! Output history, average, diagnostic files parameters.
!
LDEFOUT == T
NHIS == 60 ! Every 10 minutes
NDEFHIS == 0
NTSAVG == 1
NAVG == 60 ! Every 10 minutes
NDEFAVG == 0
NTSDIA == 1
NDIA == 60 ! Every 10 minutes
NDEFDIA == 0
Diffusion coefficients
I'm not sure why, but the values below work ok, while other values make ROMS to blow up
! Harmonic/biharmonic horizontal diffusion of tracer: [1:NAT+NPT,Ngrids].
!
TNU2 == 20.0d0 20.0d0 ! m2/s
TNU4 == 2*0.0d0 ! m4/s
!
! Harmononic/biharmonic, horizontal viscosity coefficient: [Ngrids].
!
VISC2 == 100.0d0 ! m2/s
VISC4 == 0.0d0 ! m4/s
!
! Vertical mixing coefficients for active tracers: [1:NAT+NPT,Ngrids]
!
AKT_BAK == 1.0d-6 1.0d-6 ! m2/s
Vertical S-coordinates parameters
Adjust parameters for the depth of the our fjord
! Vertical S-coordinates parameters, [1:Ngrids].
!
THETA_S == 5.0d0 ! 0 < THETA_S < 20
THETA_B == 0.4d0 ! 0 < THETA_B < 1
TCLINE == 30.0d0 ! m
Time-stamp at start and reference time
! Time-stamp assigned for model initialization, reference time
! origin for tidal forcing, and model reference time for output
! NetCDF units attribute.
!
DSTART = 52791.0d0 ! days
TIDE_START = 52791.0d0 ! days
TIME_REF = 18581117.0 ! yyyymmdd.dd
Logical switches
! Logical switches (TRUE/FALSE) to specify the state surface forcing ! variable whose stochastic optimals is required. ! SOstate(isUstr) == F ! surface u-stress SOstate(isVstr) == F ! surface v-stress
To speed up simulations
Hout(inert) == F ! inert passive tracers
Hout(idBott) == F F F F F F F F F F F F F F F F
Input files
Write path to the GRID and INITIALIZATION input files. !Comment-out the ones that we don't need.
! Input NetCDF file names, [1:Ngrids].
!
GRDNAME == ../PATH2yourFILE/Fjord_grd.nc
ININAME == ../PATH2yourFILE/Fjord_ini.nc
! ITLNAME == ocean_itl.nc
! IRPNAME == ocean_irp.nc
! IADNAME == ocean_iad.nc
! CLMNAME == ocean_clm.nc
! BRYNAME == ocean_bry.nc
! FWDNAME == ocean_fwd.nc
! ADSNAME == ocean_ads.nc
Forcing files
No forcing. !Comment-out the what we don't need.
! NFFILES == 1 ! number of forcing files ! ! FRCNAME == ocean_frc.nc ! forcing file 1, grid 1
Output files
Finally, the output file names (may want to change all, eventhough only his, avg and rst will be used)
! Output NetCDF file names, [1:Ngrids].
!
GSTNAME == Fjord_ocean_gst.nc
RSTNAME == Fjord_ocean_rst.nc
HISNAME == Fjord_ocean_his.nc
TLMNAME == Fjord_ocean_tlm.nc
TLFNAME == Fjord_ocean_tlf.nc
ADJNAME == Fjord_ocean_adj.nc
AVGNAME == Fjord_ocean_avg.nc
DIANAME == Fjord_ocean_dia.nc
STANAME == Fjord_ocean_sta.nc
FLTNAME == Fjord_ocean_flt.nc
Running ROMS
Similar to the "Getting Started" tutorial...
- To run ROMS in serial, just type:
oceanS < ROMS/External/roms_fjord.in > & log &
- or to run in parallel (distributed-memory) on four processors:
mpirun -np 4 oceanM ROMS/External/roms_fjord.in > & log &