MOM6
mom_boundary_update Module Reference

Detailed Description

Controls where open boundary conditions are applied.

This module updates the open boundary arrays when time-varying. It caused a circular dependency with the tidal_bay setup when MOM_open_boundary.

A small fragment of the grid is shown below:

j+1 x ^ x ^ x At x: q, CoriolisBu j+1 > o > o > At ^: v, tauy j x ^ x ^ x At >: u, taux j > o > o > At o: h, bathyT, buoy, tr, T, S, Rml, ustar j-1 x ^ x ^ x i-1 i i+1 At x & ^: i i+1 At > & o:

The boundaries always run through q grid points (x).

Data Types

type  update_obc_cs
 

Functions/Subroutines

subroutine, public call_obc_register (param_file, CS, OBC)
 The following subroutines and associated definitions provide the machinery to register and call the subroutines that initialize open boundary conditions. More...
 
subroutine, public update_obc_data (OBC, G, GV, tv, h, CS, Time)
 Calls appropriate routine to update the open boundary conditions. More...
 
subroutine, public obc_register_end (CS)
 Clean up the OBC registry. More...
 

Variables

integer id_clock_pass
 
character(len=40) mdl = "MOM_boundary_update"
 

Function/Subroutine Documentation

◆ call_obc_register()

subroutine, public mom_boundary_update::call_obc_register ( type(param_file_type), intent(in)  param_file,
type(update_obc_cs), pointer  CS,
type(ocean_obc_type), pointer  OBC 
)

The following subroutines and associated definitions provide the machinery to register and call the subroutines that initialize open boundary conditions.

Parameters
[in]param_fileParameter file to parse
csControl structure for OBCs
obcOpen boundary structure

Definition at line 56 of file MOM_boundary_update.F90.

References mom_error_handler::mom_error(), and mom_open_boundary::register_file_obc().

Referenced by mom::initialize_mom().

56  type(param_file_type), intent(in) :: param_file !< Parameter file to parse
57  type(update_obc_cs), pointer :: cs !< Control structure for OBCs
58  type(ocean_obc_type), pointer :: obc !< Open boundary structure
59  character(len=40) :: mdl = "MOM_boundary_update" ! This module's name.
60 
61  if (associated(cs)) then
62  call mom_error(warning, "call_OBC_register called with an associated "// &
63  "control structure.")
64  return
65  else ; allocate(cs) ; endif
66 
67  call log_version(param_file, mdl, version, "")
68 
69  call get_param(param_file, mdl, "USE_FILE_OBC", cs%use_files, &
70  "If true, use external files for the open boundary.", &
71  default=.false.)
72  call get_param(param_file, mdl, "USE_TIDAL_BAY_OBC", cs%use_tidal_bay, &
73  "If true, use the tidal_bay open boundary.", &
74  default=.false.)
75  call get_param(param_file, mdl, "USE_KELVIN_WAVE_OBC", cs%use_Kelvin, &
76  "If true, use the Kelvin wave open boundary.", &
77  default=.false.)
78  call get_param(param_file, mdl, "USE_SHELFWAVE_OBC", cs%use_shelfwave, &
79  "If true, use the shelfwave open boundary.", &
80  default=.false.)
81 
82  if (cs%use_files) cs%use_files = &
83  register_file_obc(param_file, cs%file_OBC_CSp, &
84  obc%OBC_Reg)
85  if (cs%use_tidal_bay) cs%use_tidal_bay = &
86  register_tidal_bay_obc(param_file, cs%tidal_bay_OBC_CSp, &
87  obc%OBC_Reg)
88  if (cs%use_Kelvin) cs%use_Kelvin = &
89  register_kelvin_obc(param_file, cs%Kelvin_OBC_CSp, &
90  obc%OBC_Reg)
91  if (cs%use_shelfwave) cs%use_shelfwave = &
92  register_shelfwave_obc(param_file, cs%shelfwave_OBC_CSp, &
93  obc%OBC_Reg)
94 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ obc_register_end()

subroutine, public mom_boundary_update::obc_register_end ( type(update_obc_cs), pointer  CS)

Clean up the OBC registry.

Parameters
csControl structure for OBCs

Definition at line 136 of file MOM_boundary_update.F90.

References mom_open_boundary::file_obc_end(), kelvin_initialization::kelvin_obc_end(), and tidal_bay_initialization::tidal_bay_obc_end().

136  type(update_obc_cs), pointer :: cs !< Control structure for OBCs
137 
138  if (cs%use_files) call file_obc_end(cs%file_OBC_CSp)
139  if (cs%use_tidal_bay) call tidal_bay_obc_end(cs%tidal_bay_OBC_CSp)
140  if (cs%use_Kelvin) call kelvin_obc_end(cs%Kelvin_OBC_CSp)
141 
142  if (associated(cs)) deallocate(cs)
Here is the call graph for this function:

◆ update_obc_data()

subroutine, public mom_boundary_update::update_obc_data ( type(ocean_obc_type), pointer  OBC,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  h,
type(update_obc_cs), pointer  CS,
type(time_type), intent(in)  Time 
)

Calls appropriate routine to update the open boundary conditions.

Parameters
[in]gOcean grid structure
[in]gvOcean vertical grid structure
[in]tvThermodynamics structure
[in,out]hlayer thickness
obcOpen boundary structure
csControl structure for OBCs
[in]timeModel time

Definition at line 99 of file MOM_boundary_update.F90.

Referenced by mom_dynamics_legacy_split::step_mom_dyn_legacy_split(), mom_dynamics_split_rk2::step_mom_dyn_split_rk2(), mom_dynamics_unsplit::step_mom_dyn_unsplit(), and mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2().

99  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
100  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
101  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
102  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< layer thickness
103  type(ocean_obc_type), pointer :: obc !< Open boundary structure
104  type(update_obc_cs), pointer :: cs !< Control structure for OBCs
105  type(time_type), intent(in) :: time !< Model time
106  ! Local variables
107  logical :: read_obc_eta = .false.
108  logical :: read_obc_uv = .false.
109  logical :: read_obc_ts = .false.
110  integer :: i, j, k, itt, is, ie, js, je, isd, ied, jsd, jed, nz
111  integer :: isd_off, jsd_off
112  integer :: isdb, iedb, jsdb, jedb
113  character(len=40) :: mdl = "update_OBC_data" ! This subroutine's name.
114  character(len=200) :: filename, obc_file, inputdir ! Strings for file/path
115 
116  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
117  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
118  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
119 
120 ! Something here... with CS%file_OBC_CSp?
121 ! if (CS%use_files) &
122 ! call update_OBC_segment_data(G, GV, OBC, tv, h, Time)
123  if (cs%use_tidal_bay) &
124  call tidal_bay_set_obc_data(obc, cs%tidal_bay_OBC_CSp, g, h, time)
125  if (cs%use_Kelvin) &
126  call kelvin_set_obc_data(obc, cs%Kelvin_OBC_CSp, g, h, time)
127  if (cs%use_shelfwave) &
128  call shelfwave_set_obc_data(obc, cs%shelfwave_OBC_CSp, g, h, time)
129  if (obc%needs_IO_for_data) &
130  call update_obc_segment_data(g, gv, obc, tv, h, time)
131 
Here is the caller graph for this function:

Variable Documentation

◆ id_clock_pass

integer mom_boundary_update::id_clock_pass

Definition at line 44 of file MOM_boundary_update.F90.

44 integer :: id_clock_pass

◆ mdl

character(len=40) mom_boundary_update::mdl = "MOM_boundary_update"
private

Definition at line 46 of file MOM_boundary_update.F90.

46 character(len=40) :: mdl = "MOM_boundary_update" ! This module's name.