MOM6
|
Controls where open boundary conditions are applied.
This module implements some aspects of internal open boundary conditions in MOM.
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 | file_obc_cs |
Control structure for open boundaries that read from files. Probably lots to update here. More... | |
type | obc_registry_type |
Type to carry basic OBC information needed for updating values. More... | |
type | obc_segment_data_type |
Open boundary segment data from files (mostly). More... | |
type | obc_segment_tracer_type |
Tracer on OBC segment data structure, for putting into a segment tracer registry. More... | |
type | obc_segment_type |
Open boundary segment data structure. More... | |
type | obc_struct_type |
Type to carry something (what] for the OBC registry. More... | |
type | ocean_obc_type |
Open-boundary data. More... | |
type | segment_tracer_registry_type |
Registry type for tracers on segments. More... | |
Functions/Subroutines | |
subroutine, public | open_boundary_config (G, US, param_file, OBC) |
Enables OBC module and reads configuration parameters This routine is called from MOM_initialize_fixed which occurs before the initialization of the vertical coordinate and ALE_init. Therefore segment data are not fully initialized here. The remainder of the segment data are initialized in a later call to update_open_boundary_data. More... | |
subroutine | initialize_segment_data (G, OBC, PF) |
Allocate space for reading OBC data from files. It sets up the required vertical remapping. In the process, it does funky stuff with the MPI processes. More... | |
subroutine | setup_segment_indices (G, seg, Is_obc, Ie_obc, Js_obc, Je_obc) |
Define indices for segment and store in hor_index_type using global segment bounds corresponding to q-points. More... | |
subroutine | setup_u_point_obc (OBC, G, US, segment_str, l_seg, PF, reentrant_y) |
Parse an OBC_SEGMENT_%%% string starting with "I=" and configure placement and type of OBC accordingly. More... | |
subroutine | setup_v_point_obc (OBC, G, US, segment_str, l_seg, PF, reentrant_x) |
Parse an OBC_SEGMENT_%%% string starting with "J=" and configure placement and type of OBC accordingly. More... | |
subroutine | parse_segment_str (ni_global, nj_global, segment_str, l, m, n, action_str, reentrant) |
Parse an OBC_SEGMENT_%%% string. More... | |
subroutine | parse_segment_data_str (segment_str, var, value, filenam, fieldnam, fields, num_fields, debug) |
Parse an OBC_SEGMENT_%%_DATA string. More... | |
subroutine | parse_for_tracer_reservoirs (OBC, PF, use_temperature) |
Parse all the OBC_SEGMENT_%%_DATA strings again to see which need tracer reservoirs (all pes need to know). More... | |
subroutine | parse_segment_param_real (segment_str, var, param_value, debug) |
Parse an OBC_SEGMENT_%%_PARAMS string. More... | |
subroutine, public | open_boundary_init (G, GV, US, param_file, OBC, restart_CSp) |
Initialize open boundary control structure and do any necessary rescaling of OBC fields that have been read from a restart file. More... | |
logical function, public | open_boundary_query (OBC, apply_open_OBC, apply_specified_OBC, apply_Flather_OBC, apply_nudged_OBC, needs_ext_seg_data) |
subroutine | open_boundary_dealloc (OBC) |
Deallocate open boundary data. More... | |
subroutine, public | open_boundary_end (OBC) |
Close open boundary data. More... | |
subroutine, public | open_boundary_impose_normal_slope (OBC, G, depth) |
Sets the slope of bathymetry normal to an open bounndary to zero. More... | |
subroutine, public | open_boundary_impose_land_mask (OBC, G, areaCu, areaCv, US) |
Reconcile masks and open boundaries, deallocate OBC on PEs where it is not needed. Also adjust u- and v-point cell area on specified open boundaries and mask all points outside open boundaries. More... | |
subroutine | setup_obc_tracer_reservoirs (G, OBC) |
Make sure the OBC tracer reservoirs are initialized. More... | |
subroutine, public | radiation_open_bdry_conds (OBC, u_new, u_old, v_new, v_old, G, US, dt) |
Apply radiation conditions to 3D u,v at open boundaries. More... | |
subroutine, public | open_boundary_apply_normal_flow (OBC, G, u, v) |
Applies OBC values stored in segments to 3d u,v fields. More... | |
subroutine, public | open_boundary_zero_normal_flow (OBC, G, u, v) |
Applies zero values to 3d u,v fields on OBC segments. More... | |
subroutine | gradient_at_q_points (G, segment, uvel, vvel) |
Calculate the tangential gradient of the normal flow at the boundary q-points. More... | |
subroutine, public | set_tracer_data (OBC, tv, h, G, PF, tracer_Reg) |
Sets the initial values of the tracer open boundary conditions. Redoing this elsewhere. More... | |
integer function | lookup_seg_field (OBC_seg, field) |
Needs documentation. More... | |
subroutine | allocate_obc_segment_data (OBC, segment) |
Allocate segment data fields. More... | |
subroutine | deallocate_obc_segment_data (OBC, segment) |
Deallocate segment data fields. More... | |
subroutine, public | open_boundary_test_extern_uv (G, OBC, u, v) |
Set tangential velocities outside of open boundaries to silly values (used for checking the interior state is independent of values outside of the domain). More... | |
subroutine, public | open_boundary_test_extern_h (G, GV, OBC, h) |
Set thicknesses outside of open boundaries to silly values (used for checking the interior state is independent of values outside of the domain). More... | |
subroutine, public | update_obc_segment_data (G, GV, US, OBC, tv, h, Time) |
Update the OBC values on the segments. More... | |
subroutine, public | register_obc (name, param_file, Reg) |
register open boundary objects for boundary updates. More... | |
subroutine, public | obc_registry_init (param_file, Reg) |
This routine include declares and sets the variable "version". More... | |
logical function, public | register_file_obc (param_file, CS, OBC_Reg) |
Add file to OBC registry. More... | |
subroutine, public | file_obc_end (CS) |
Clean up the file OBC from registry. More... | |
subroutine, public | segment_tracer_registry_init (param_file, segment) |
Initialize the segment tracer registry. More... | |
subroutine, public | register_segment_tracer (tr_ptr, param_file, GV, segment, OBC_scalar, OBC_array) |
subroutine, public | segment_tracer_registry_end (Reg) |
Clean up the segment tracer registry. More... | |
subroutine, public | register_temp_salt_segments (GV, OBC, tr_Reg, param_file) |
subroutine, public | fill_temp_salt_segments (G, OBC, tv) |
subroutine | mask_outside_obcs (G, US, param_file, OBC) |
Find the region outside of all open boundary segments and make sure it is set to land mask. Gonna need to know global land mask as well to get it right... More... | |
subroutine | flood_fill (G, color, cin, cout, cland) |
flood the cin, cout values More... | |
subroutine | flood_fill2 (G, color, cin, cout, cland) |
flood the cin, cout values More... | |
subroutine, public | open_boundary_register_restarts (HI, GV, OBC, Reg, param_file, restart_CSp, use_temperature) |
Register OBC segment data for restarts. More... | |
subroutine, public | update_segment_tracer_reservoirs (G, GV, uhr, vhr, h, OBC, dt, Reg) |
Update the OBC tracer reservoirs after the tracers have been updated. More... | |
subroutine | adjustsegmentetatofitbathymetry (G, GV, US, segment, fld) |
Adjust interface heights to fit the bathymetry and diagnose layer thickness. More... | |
Variables | |
integer, parameter, public | obc_none = 0 |
Indicates the use of no open boundary. More... | |
integer, parameter, public | obc_simple = 1 |
Indicates the use of a simple inflow open boundary. More... | |
integer, parameter, public | obc_wall = 2 |
Indicates the use of a closed wall. More... | |
integer, parameter, public | obc_flather = 3 |
Indicates the use of a Flather open boundary. More... | |
integer, parameter, public | obc_radiation = 4 |
Indicates the use of a radiation open boundary. More... | |
integer, parameter, public | obc_direction_n = 100 |
Indicates the boundary is an effective northern boundary. More... | |
integer, parameter, public | obc_direction_s = 200 |
Indicates the boundary is an effective southern boundary. More... | |
integer, parameter, public | obc_direction_e = 300 |
Indicates the boundary is an effective eastern boundary. More... | |
integer, parameter, public | obc_direction_w = 400 |
Indicates the boundary is an effective western boundary. More... | |
integer, parameter | max_obc_fields = 100 |
Maximum number of data fields needed for OBC segments. More... | |
integer | id_clock_pass |
A CPU time clock. More... | |
character(len=40) | mdl = "MOM_open_boundary" |
This module's name. More... | |
|
private |
Adjust interface heights to fit the bathymetry and diagnose layer thickness.
If the bottom most interface is below the topography then the bottom-most layers are contracted to GVAngstrom_m. If the bottom most interface is above the topography then the entire column is dilated (expanded) to fill the void.
[in] | g | The ocean's grid structure |
[in] | gv | The ocean's vertical grid structure |
[in] | us | A dimensional unit scaling type |
[in,out] | segment | pointer to segment type |
[in] | fld | field index to adjust thickness |
Definition at line 4593 of file MOM_open_boundary.F90.
Referenced by update_obc_segment_data().
|
private |
Allocate segment data fields.
obc | Open boundary structure | |
[in,out] | segment | Open boundary segment |
Definition at line 3171 of file MOM_open_boundary.F90.
Referenced by setup_u_point_obc(), and setup_v_point_obc().
|
private |
Deallocate segment data fields.
obc | Open boundary structure | |
[in,out] | segment | Open boundary segment |
Definition at line 3278 of file MOM_open_boundary.F90.
References segment_tracer_registry_end().
Referenced by open_boundary_dealloc().
subroutine, public mom_open_boundary::file_obc_end | ( | type(file_obc_cs), pointer | CS | ) |
Clean up the file OBC from registry.
cs | OBC file control structure. |
Definition at line 3972 of file MOM_open_boundary.F90.
Referenced by mom_boundary_update::obc_register_end().
subroutine, public mom_open_boundary::fill_temp_salt_segments | ( | type(ocean_grid_type), intent(inout) | G, |
type(ocean_obc_type), pointer | OBC, | ||
type(thermo_var_ptrs), intent(inout) | tv | ||
) |
[in,out] | g | Ocean grid structure |
obc | Open boundary structure | |
[in,out] | tv | Thermodynamics structure |
Definition at line 4122 of file MOM_open_boundary.F90.
References obc_direction_s, obc_direction_w, and setup_obc_tracer_reservoirs().
|
private |
flood the cin, cout values
[in,out] | g | Ocean grid structure |
[in,out] | color | For sorting inside from outside |
[in] | cin | color for inside the domain |
[in] | cout | color for outside the domain |
[in] | cland | color for inside the land mask |
Definition at line 4292 of file MOM_open_boundary.F90.
Referenced by mask_outside_obcs().
|
private |
flood the cin, cout values
[in,out] | g | Ocean grid structure |
[in,out] | color | For sorting inside from outside |
[in] | cin | color for inside the domain |
[in] | cout | color for outside the domain |
[in] | cland | color for inside the land mask |
Definition at line 4352 of file MOM_open_boundary.F90.
Referenced by mask_outside_obcs().
|
private |
Calculate the tangential gradient of the normal flow at the boundary q-points.
[in] | g | Ocean grid structure |
segment | OBC segment structure | |
[in] | uvel | zonal velocity [L T-1 ~> m s-1] |
[in] | vvel | meridional velocity [L T-1 ~> m s-1] |
Definition at line 2966 of file MOM_open_boundary.F90.
References obc_direction_e, and obc_direction_n.
Referenced by radiation_open_bdry_conds().
|
private |
Allocate space for reading OBC data from files. It sets up the required vertical remapping. In the process, it does funky stuff with the MPI processes.
[in] | g | Ocean grid structure |
[in,out] | obc | Open boundary control structure |
[in] | pf | Parameter file handle |
Definition at line 541 of file MOM_open_boundary.F90.
References mom_remapping::initialize_remapping(), mdl, mom_error_handler::mom_error(), mom_error_handler::mom_mesg(), and parse_segment_data_str().
Referenced by open_boundary_config().
|
private |
Needs documentation.
obc_seg | OBC segment | |
[in] | field | The field name |
Definition at line 3152 of file MOM_open_boundary.F90.
|
private |
Find the region outside of all open boundary segments and make sure it is set to land mask. Gonna need to know global land mask as well to get it right...
[in,out] | g | Ocean grid structure |
[in] | param_file | Parameter file handle |
obc | Open boundary structure | |
[in] | us | A dimensional unit scaling type |
Definition at line 4183 of file MOM_open_boundary.F90.
References flood_fill(), flood_fill2(), mdl, obc_direction_e, obc_direction_n, obc_direction_s, and obc_direction_w.
Referenced by open_boundary_config().
subroutine, public mom_open_boundary::obc_registry_init | ( | type(param_file_type), intent(in) | param_file, |
type(obc_registry_type), pointer | Reg | ||
) |
This routine include declares and sets the variable "version".
[in] | param_file | open file to parse for model parameters |
reg | pointer to OBC registry |
Definition at line 3925 of file MOM_open_boundary.F90.
Referenced by register_obc().
subroutine, public mom_open_boundary::open_boundary_apply_normal_flow | ( | type(ocean_obc_type), pointer | OBC, |
type(ocean_grid_type), intent(inout) | G, | ||
real, dimension(szib_(g),szj_(g),szk_(g)), intent(inout) | u, | ||
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(inout) | v | ||
) |
Applies OBC values stored in segments to 3d u,v fields.
obc | Open boundary control structure | |
[in,out] | g | Ocean grid structure |
[in,out] | u | u field to update on open boundaries [L T-1 ~> m s-1] |
[in,out] | v | v field to update on open boundaries [L T-1 ~> m s-1] |
Definition at line 2898 of file MOM_open_boundary.F90.
Referenced by radiation_open_bdry_conds().
subroutine, public mom_open_boundary::open_boundary_config | ( | type(dyn_horgrid_type), intent(inout) | G, |
type(unit_scale_type), intent(in) | US, | ||
type(param_file_type), intent(in) | param_file, | ||
type(ocean_obc_type), pointer | OBC | ||
) |
Enables OBC module and reads configuration parameters This routine is called from MOM_initialize_fixed which occurs before the initialization of the vertical coordinate and ALE_init. Therefore segment data are not fully initialized here. The remainder of the segment data are initialized in a later call to update_open_boundary_data.
[in,out] | g | Ocean grid structure |
[in] | us | A dimensional unit scaling type |
[in] | param_file | Parameter file handle |
obc | Open boundary control structure |
Definition at line 317 of file MOM_open_boundary.F90.
References initialize_segment_data(), mask_outside_obcs(), mdl, mom_error_handler::mom_error(), obc_none, open_boundary_dealloc(), open_boundary_query(), mom_string_functions::remove_spaces(), setup_u_point_obc(), and setup_v_point_obc().
|
private |
Deallocate open boundary data.
obc | Open boundary control structure |
Definition at line 1563 of file MOM_open_boundary.F90.
References deallocate_obc_segment_data().
Referenced by open_boundary_config(), and open_boundary_end().
subroutine, public mom_open_boundary::open_boundary_end | ( | type(ocean_obc_type), pointer | OBC | ) |
Close open boundary data.
obc | Open boundary control structure |
Definition at line 1588 of file MOM_open_boundary.F90.
References open_boundary_dealloc().
subroutine, public mom_open_boundary::open_boundary_impose_land_mask | ( | type(ocean_obc_type), pointer | OBC, |
type(dyn_horgrid_type), intent(inout) | G, | ||
real, dimension(szib_(g),szj_(g)), intent(inout) | areaCu, | ||
real, dimension(szi_(g),szjb_(g)), intent(inout) | areaCv, | ||
type(unit_scale_type), intent(in) | US | ||
) |
Reconcile masks and open boundaries, deallocate OBC on PEs where it is not needed. Also adjust u- and v-point cell area on specified open boundaries and mask all points outside open boundaries.
obc | Open boundary control structure | |
[in,out] | g | Ocean grid structure |
[in] | us | A dimensional unit scaling type |
[in,out] | areacu | Area of a u-cell [L2 ~> m2] |
[in,out] | areacv | Area of a u-cell [L2 ~> m2] |
Definition at line 1639 of file MOM_open_boundary.F90.
References obc_direction_e, obc_direction_s, obc_direction_w, and obc_none.
Referenced by mom_fixed_initialization::mom_initialize_fixed().
subroutine, public mom_open_boundary::open_boundary_impose_normal_slope | ( | type(ocean_obc_type), pointer | OBC, |
type(dyn_horgrid_type), intent(in) | G, | ||
real, dimension(szi_(g),szj_(g)), intent(inout) | depth | ||
) |
Sets the slope of bathymetry normal to an open bounndary to zero.
obc | Open boundary control structure | |
[in] | g | Ocean grid structure |
[in,out] | depth | Bathymetry at h-points |
Definition at line 1594 of file MOM_open_boundary.F90.
References obc_direction_e, obc_direction_n, obc_direction_s, and obc_direction_w.
subroutine, public mom_open_boundary::open_boundary_init | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(unit_scale_type), intent(in) | US, | ||
type(param_file_type), intent(in) | param_file, | ||
type(ocean_obc_type), pointer | OBC, | ||
type(mom_restart_cs), pointer | restart_CSp | ||
) |
Initialize open boundary control structure and do any necessary rescaling of OBC fields that have been read from a restart file.
[in] | g | Ocean grid structure |
[in] | gv | Container for vertical grid information |
[in] | us | A dimensional unit scaling type |
[in] | param_file | Parameter file handle |
obc | Open boundary control structure | |
restart_csp | Restart structure, data intent(inout) |
Definition at line 1479 of file MOM_open_boundary.F90.
References id_clock_pass.
logical function, public mom_open_boundary::open_boundary_query | ( | type(ocean_obc_type), pointer | OBC, |
logical, intent(in), optional | apply_open_OBC, | ||
logical, intent(in), optional | apply_specified_OBC, | ||
logical, intent(in), optional | apply_Flather_OBC, | ||
logical, intent(in), optional | apply_nudged_OBC, | ||
logical, intent(in), optional | needs_ext_seg_data | ||
) |
obc | Open boundary control structure | |
[in] | apply_open_obc | Returns True if open_*_BCs_exist_globally is true |
[in] | apply_specified_obc | Returns True if specified_*_BCs_exist_globally is true |
[in] | apply_flather_obc | Returns True if Flather_*_BCs_exist_globally is true |
[in] | apply_nudged_obc | Returns True if nudged_*_BCs_exist_globally is true |
[in] | needs_ext_seg_data | Returns True if external segment data needed |
Definition at line 1541 of file MOM_open_boundary.F90.
Referenced by open_boundary_config().
subroutine, public mom_open_boundary::open_boundary_register_restarts | ( | type(hor_index_type), intent(in) | HI, |
type(verticalgrid_type), pointer | GV, | ||
type(ocean_obc_type), pointer | OBC, | ||
type(tracer_registry_type), pointer | Reg, | ||
type(param_file_type), intent(in) | param_file, | ||
type(mom_restart_cs), pointer | restart_CSp, | ||
logical, intent(in) | use_temperature | ||
) |
Register OBC segment data for restarts.
[in] | hi | Horizontal indices |
gv | Container for vertical grid information | |
obc | OBC data structure, data intent(inout) | |
reg | pointer to tracer registry | |
[in] | param_file | Parameter file handle |
restart_csp | Restart structure, data intent(inout) | |
[in] | use_temperature | If true, T and S are used |
Definition at line 4413 of file MOM_open_boundary.F90.
References parse_for_tracer_reservoirs().
subroutine, public mom_open_boundary::open_boundary_test_extern_h | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(ocean_obc_type), pointer | OBC, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout) | h | ||
) |
Set thicknesses outside of open boundaries to silly values (used for checking the interior state is independent of values outside of the domain).
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
obc | Open boundary structure | |
[in,out] | h | Layer thickness [H ~> m or kg m-2] |
Definition at line 3357 of file MOM_open_boundary.F90.
References obc_direction_e, and obc_direction_n.
Referenced by mom_dynamics_split_rk2::step_mom_dyn_split_rk2().
subroutine, public mom_open_boundary::open_boundary_test_extern_uv | ( | type(ocean_grid_type), intent(in) | G, |
type(ocean_obc_type), pointer | OBC, | ||
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(inout) | u, | ||
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(inout) | v | ||
) |
Set tangential velocities outside of open boundaries to silly values (used for checking the interior state is independent of values outside of the domain).
[in] | g | Ocean grid structure |
obc | Open boundary structure | |
[in,out] | u | Zonal velocity [L T-1 ~> m s-1] |
[in,out] | v | Meridional velocity [L T-1 ~> m s-1] |
Definition at line 3314 of file MOM_open_boundary.F90.
References obc_direction_e, and obc_direction_n.
subroutine, public mom_open_boundary::open_boundary_zero_normal_flow | ( | type(ocean_obc_type), pointer | OBC, |
type(ocean_grid_type), intent(inout) | G, | ||
real, dimension(szib_(g),szj_(g),szk_(g)), intent(inout) | u, | ||
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(inout) | v | ||
) |
Applies zero values to 3d u,v fields on OBC segments.
obc | Open boundary control structure | |
[in,out] | g | Ocean grid structure |
[in,out] | u | u field to update on open boundaries |
[in,out] | v | v field to update on open boundaries |
Definition at line 2934 of file MOM_open_boundary.F90.
Referenced by mom_dynamics_unsplit::step_mom_dyn_unsplit(), and mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2().
|
private |
Parse all the OBC_SEGMENT_%%_DATA strings again to see which need tracer reservoirs (all pes need to know).
[in,out] | obc | Open boundary control structure |
[in] | pf | Parameter file handle |
[in] | use_temperature | If true, T and S are used |
Definition at line 1332 of file MOM_open_boundary.F90.
References mdl, and parse_segment_data_str().
Referenced by open_boundary_register_restarts().
|
private |
Parse an OBC_SEGMENT_%%_DATA string.
[in] | segment_str | A string in form of "VAR1=file:foo1.nc(varnam1),VAR2=file:foo2.nc(varnam2),..." |
[in] | var | The name of the variable for which parameters are needed |
[out] | filenam | The name of the input file if using "file" method |
[out] | fieldnam | The name of the variable in the input file if using "file" method |
[out] | value | A constant value if using the "value" method |
[out] | fields | List of fieldnames for each segment |
[out] | num_fields | The number of fields in the segment data |
[in] | debug | If present and true, write verbose debugging messages |
Definition at line 1243 of file MOM_open_boundary.F90.
References mom_string_functions::extract_word(), and mom_error_handler::mom_error().
Referenced by initialize_segment_data(), and parse_for_tracer_reservoirs().
|
private |
Parse an OBC_SEGMENT_%%_PARAMS string.
[in] | segment_str | A string in form of "VAR1=file:foo1.nc(varnam1),VAR2=file:foo2.nc(varnam2),..." |
[in] | var | The name of the variable for which parameters are needed |
[out] | param_value | The value of the parameter |
[in] | debug | If present and true, write verbose debugging messages |
Definition at line 1396 of file MOM_open_boundary.F90.
|
private |
Parse an OBC_SEGMENT_%%% string.
[in] | ni_global | Number of h-points in zonal direction |
[in] | nj_global | Number of h-points in meridional direction |
[in] | segment_str | A string in form of "I=l,J=m:n,string" or "J=l,I=m,n,string" |
[out] | l | The value of I=l, if segment_str begins with I=l, or the value of J=l |
[out] | m | The value of J=m, if segment_str begins with I=, or the value of I=m |
[out] | n | The value of J=n, if segment_str begins with I=, or the value of I=n |
[out] | action_str | The "string" part of segment_str |
[in] | reentrant | is domain reentrant in relevant direction? |
Definition at line 1129 of file MOM_open_boundary.F90.
References mom_string_functions::extract_word(), interpret_int_expr(), and mom_error_handler::mom_error().
Referenced by setup_u_point_obc(), and setup_v_point_obc().
subroutine, public mom_open_boundary::radiation_open_bdry_conds | ( | type(ocean_obc_type), pointer | OBC, |
real, dimension(szib_(g),szj_(g),szk_(g)), intent(inout) | u_new, | ||
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in) | u_old, | ||
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(inout) | v_new, | ||
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in) | v_old, | ||
type(ocean_grid_type), intent(inout) | G, | ||
type(unit_scale_type), intent(in) | US, | ||
real, intent(in) | dt | ||
) |
Apply radiation conditions to 3D u,v at open boundaries.
[in,out] | g | Ocean grid structure |
obc | Open boundary control structure | |
[in,out] | u_new | On exit, new u values on open boundaries On entry, the old time-level v but including barotropic accelerations [L T-1 ~> m s-1]. |
[in] | u_old | Original unadjusted u [L T-1 ~> m s-1] |
[in,out] | v_new | On exit, new v values on open boundaries. On entry, the old time-level v but including barotropic accelerations [L T-1 ~> m s-1]. |
[in] | v_old | Original unadjusted v [L T-1 ~> m s-1] |
[in] | us | A dimensional unit scaling type |
[in] | dt | Appropriate timestep [T ~> s] |
Definition at line 1788 of file MOM_open_boundary.F90.
References gradient_at_q_points(), id_clock_pass, obc_direction_e, obc_direction_n, obc_direction_s, obc_direction_w, and open_boundary_apply_normal_flow().
logical function, public mom_open_boundary::register_file_obc | ( | type(param_file_type), intent(in) | param_file, |
type(file_obc_cs), pointer | CS, | ||
type(obc_registry_type), pointer | OBC_Reg | ||
) |
Add file to OBC registry.
[in] | param_file | parameter file. |
cs | file control structure. | |
obc_reg | OBC registry. |
Definition at line 3951 of file MOM_open_boundary.F90.
References register_obc().
Referenced by mom_boundary_update::call_obc_register().
subroutine, public mom_open_boundary::register_obc | ( | character(len=32), intent(in) | name, |
type(param_file_type), intent(in) | param_file, | ||
type(obc_registry_type), pointer | Reg | ||
) |
register open boundary objects for boundary updates.
[in] | name | OBC name used for error messages |
[in] | param_file | file to parse for model parameter values |
reg | pointer to the tracer registry |
Definition at line 3899 of file MOM_open_boundary.F90.
References obc_registry_init().
Referenced by dyed_channel_initialization::register_dyed_channel_obc(), and register_file_obc().
subroutine, public mom_open_boundary::register_segment_tracer | ( | type(tracer_type), target | tr_ptr, |
type(param_file_type), intent(in) | param_file, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(obc_segment_type), intent(inout) | segment, | ||
real, intent(in), optional | OBC_scalar, | ||
logical, intent(in), optional | OBC_array | ||
) |
[in] | gv | ocean vertical grid structure |
tr_ptr | A target that can be used to set a pointer to the stored value of tr. This target must be an enduring part of the control structure, because the tracer registry will use this memory, but it also means that any updates to this structure in the calling module will be available subsequently to the tracer registry. | |
[in] | param_file | file to parse for model parameter values |
[in,out] | segment | current segment data structure |
[in] | obc_scalar | If present, use scalar value for segment tracer inflow concentration. |
[in] | obc_array | If true, use array values for segment tracer inflow concentration. |
Definition at line 4013 of file MOM_open_boundary.F90.
References segment_tracer_registry_init().
Referenced by dome_initialization::dome_set_obc_data(), dyed_obcs_initialization::dyed_obcs_set_obc_data(), and register_temp_salt_segments().
subroutine, public mom_open_boundary::register_temp_salt_segments | ( | type(verticalgrid_type), intent(in) | GV, |
type(ocean_obc_type), pointer | OBC, | ||
type(tracer_registry_type), pointer | tr_Reg, | ||
type(param_file_type), intent(in) | param_file | ||
) |
[in] | gv | ocean vertical grid structure |
obc | Open boundary structure | |
tr_reg | Tracer registry | |
[in] | param_file | file to parse for model parameter values |
Definition at line 4088 of file MOM_open_boundary.F90.
References register_segment_tracer().
subroutine, public mom_open_boundary::segment_tracer_registry_end | ( | type(segment_tracer_registry_type), pointer | Reg | ) |
Clean up the segment tracer registry.
reg | pointer to tracer registry |
Definition at line 4074 of file MOM_open_boundary.F90.
Referenced by deallocate_obc_segment_data().
subroutine, public mom_open_boundary::segment_tracer_registry_init | ( | type(param_file_type), intent(in) | param_file, |
type(obc_segment_type), intent(inout) | segment | ||
) |
Initialize the segment tracer registry.
[in] | param_file | open file to parse for model parameters |
[in,out] | segment | the segment |
Definition at line 3981 of file MOM_open_boundary.F90.
References mdl.
Referenced by register_segment_tracer().
subroutine, public mom_open_boundary::set_tracer_data | ( | type(ocean_obc_type), pointer | OBC, |
type(thermo_var_ptrs), intent(inout) | tv, | ||
real, dimension(szi_(g),szj_(g), szk_(g)), intent(inout) | h, | ||
type(ocean_grid_type), intent(inout) | G, | ||
type(param_file_type), intent(in) | PF, | ||
type(tracer_registry_type), pointer | tracer_Reg | ||
) |
Sets the initial values of the tracer open boundary conditions. Redoing this elsewhere.
[in,out] | g | Ocean grid structure |
obc | Open boundary structure | |
[in,out] | tv | Thermodynamics structure |
[in,out] | h | Thickness |
[in] | pf | Parameter file handle |
tracer_reg | Tracer registry |
Definition at line 3090 of file MOM_open_boundary.F90.
References obc_direction_e, obc_direction_n, obc_direction_s, and obc_direction_w.
|
private |
Make sure the OBC tracer reservoirs are initialized.
[in,out] | g | Ocean grid structure |
obc | Open boundary control structure |
Definition at line 1749 of file MOM_open_boundary.F90.
Referenced by fill_temp_salt_segments().
subroutine mom_open_boundary::setup_segment_indices | ( | type(dyn_horgrid_type), intent(in) | G, |
type(obc_segment_type), intent(inout) | seg, | ||
integer, intent(in) | Is_obc, | ||
integer, intent(in) | Ie_obc, | ||
integer, intent(in) | Js_obc, | ||
integer, intent(in) | Je_obc | ||
) |
Define indices for segment and store in hor_index_type using global segment bounds corresponding to q-points.
[in] | g | grid type |
[in,out] | seg | Open boundary segment |
[in] | is_obc | Q-point global i-index of start of segment |
[in] | ie_obc | Q-point global i-index of end of segment |
[in] | js_obc | Q-point global j-index of start of segment |
[in] | je_obc | Q-point global j-index of end of segment |
Definition at line 799 of file MOM_open_boundary.F90.
Referenced by setup_u_point_obc(), and setup_v_point_obc().
|
private |
Parse an OBC_SEGMENT_%%% string starting with "I=" and configure placement and type of OBC accordingly.
obc | Open boundary control structure | |
[in] | g | Ocean grid structure |
[in] | us | A dimensional unit scaling type |
[in] | segment_str | A string in form of "I=%,J=%:%,string" |
[in] | l_seg | which segment is this? |
[in] | pf | Parameter file handle |
[in] | reentrant_y | is the domain reentrant in y? |
Definition at line 858 of file MOM_open_boundary.F90.
References allocate_obc_segment_data(), mdl, mom_error_handler::mom_error(), obc_direction_e, obc_direction_w, parse_segment_str(), and setup_segment_indices().
Referenced by open_boundary_config().
|
private |
Parse an OBC_SEGMENT_%%% string starting with "J=" and configure placement and type of OBC accordingly.
obc | Open boundary control structure | |
[in] | g | Ocean grid structure |
[in] | us | A dimensional unit scaling type |
[in] | segment_str | A string in form of "J=%,I=%:%,string" |
[in] | l_seg | which segment is this? |
[in] | pf | Parameter file handle |
[in] | reentrant_x | is the domain reentrant in x? |
Definition at line 994 of file MOM_open_boundary.F90.
References allocate_obc_segment_data(), mdl, mom_error_handler::mom_error(), obc_direction_n, obc_direction_s, parse_segment_str(), and setup_segment_indices().
Referenced by open_boundary_config().
subroutine, public mom_open_boundary::update_obc_segment_data | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(unit_scale_type), intent(in) | US, | ||
type(ocean_obc_type), pointer | OBC, | ||
type(thermo_var_ptrs), intent(in) | tv, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout) | h, | ||
type(time_type), intent(in) | Time | ||
) |
Update the OBC values on the segments.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | us | A dimensional unit scaling type |
obc | Open boundary structure | |
[in] | tv | Thermodynamics structure |
[in,out] | h | Thickness [m] |
[in] | time | Model time |
Definition at line 3401 of file MOM_open_boundary.F90.
References adjustsegmentetatofitbathymetry(), obc_direction_e, obc_direction_n, obc_direction_s, and obc_direction_w.
Referenced by mom_state_initialization::mom_initialize_state().
subroutine, public mom_open_boundary::update_segment_tracer_reservoirs | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
real, dimension(szib_(g),szj_(g),szk_(g)), intent(in) | uhr, | ||
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in) | vhr, | ||
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in) | h, | ||
type(ocean_obc_type), pointer | OBC, | ||
real, intent(in) | dt, | ||
type(tracer_registry_type), pointer | Reg | ||
) |
Update the OBC tracer reservoirs after the tracers have been updated.
[in] | g | The ocean's grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | uhr | accumulated volume/mass flux through the zonal face [H L2 ~> m3 or kg] |
[in] | vhr | accumulated volume/mass flux through the meridional face [H L2 ~> m3 or kg] |
[in] | h | layer thickness after advection [H ~> m or kg m-2] |
obc | Open boundary structure | |
[in] | dt | time increment [T ~> s] |
reg | pointer to tracer registry |
Definition at line 4513 of file MOM_open_boundary.F90.
References obc_direction_s, and obc_direction_w.
Referenced by mom::step_mom_tracer_dyn().
integer mom_open_boundary::id_clock_pass |
A CPU time clock.
Definition at line 301 of file MOM_open_boundary.F90.
Referenced by open_boundary_init(), and radiation_open_bdry_conds().
|
private |
Maximum number of data fields needed for OBC segments.
Definition at line 67 of file MOM_open_boundary.F90.
|
private |
This module's name.
Definition at line 303 of file MOM_open_boundary.F90.
Referenced by initialize_segment_data(), mask_outside_obcs(), open_boundary_config(), parse_for_tracer_reservoirs(), segment_tracer_registry_init(), setup_u_point_obc(), and setup_v_point_obc().
integer, parameter, public mom_open_boundary::obc_direction_e = 300 |
Indicates the boundary is an effective eastern boundary.
Definition at line 65 of file MOM_open_boundary.F90.
Referenced by gradient_at_q_points(), mask_outside_obcs(), open_boundary_impose_land_mask(), open_boundary_impose_normal_slope(), open_boundary_test_extern_h(), open_boundary_test_extern_uv(), mom_continuity_ppm::ppm_reconstruction_x(), radiation_open_bdry_conds(), set_tracer_data(), setup_u_point_obc(), update_obc_segment_data(), mom_sum_output::write_energy(), mom_continuity_ppm::zonal_face_thickness(), mom_continuity_ppm::zonal_flux_layer(), and mom_continuity_ppm::zonal_mass_flux().
integer, parameter, public mom_open_boundary::obc_direction_n = 100 |
Indicates the boundary is an effective northern boundary.
Definition at line 63 of file MOM_open_boundary.F90.
Referenced by mom_barotropic::apply_velocity_obcs(), mom_barotropic::btcalc(), mom_barotropic::btstep(), mom_coriolisadv::coradcalc(), mom_vert_friction::find_coupling_coef(), gradient_at_q_points(), mom_hor_visc::horizontal_viscosity(), mask_outside_obcs(), mom_continuity_ppm::merid_face_thickness(), mom_continuity_ppm::merid_flux_layer(), mom_continuity_ppm::meridional_mass_flux(), open_boundary_impose_normal_slope(), open_boundary_test_extern_h(), open_boundary_test_extern_uv(), mom_continuity_ppm::ppm_reconstruction_y(), radiation_open_bdry_conds(), mom_ale::remap_all_state_vars(), set_tracer_data(), mom_barotropic::set_up_bt_obc(), setup_v_point_obc(), update_obc_segment_data(), mom_vert_friction::vertvisc_coef(), and mom_sum_output::write_energy().
integer, parameter, public mom_open_boundary::obc_direction_s = 200 |
Indicates the boundary is an effective southern boundary.
Definition at line 64 of file MOM_open_boundary.F90.
Referenced by mom_barotropic::apply_velocity_obcs(), mom_barotropic::btcalc(), mom_barotropic::btstep(), fill_temp_salt_segments(), mom_vert_friction::find_coupling_coef(), mom_hor_visc::horizontal_viscosity(), mask_outside_obcs(), open_boundary_impose_land_mask(), open_boundary_impose_normal_slope(), mom_continuity_ppm::ppm_reconstruction_y(), radiation_open_bdry_conds(), set_tracer_data(), mom_barotropic::set_up_bt_obc(), setup_v_point_obc(), update_obc_segment_data(), update_segment_tracer_reservoirs(), mom_vert_friction::vertvisc_coef(), and mom_sum_output::write_energy().
integer, parameter, public mom_open_boundary::obc_direction_w = 400 |
Indicates the boundary is an effective western boundary.
Definition at line 66 of file MOM_open_boundary.F90.
Referenced by fill_temp_salt_segments(), mom_vert_friction::find_coupling_coef(), mask_outside_obcs(), open_boundary_impose_land_mask(), open_boundary_impose_normal_slope(), mom_continuity_ppm::ppm_reconstruction_x(), radiation_open_bdry_conds(), set_tracer_data(), setup_u_point_obc(), update_obc_segment_data(), update_segment_tracer_reservoirs(), mom_vert_friction::vertvisc_coef(), and mom_sum_output::write_energy().
integer, parameter, public mom_open_boundary::obc_flather = 3 |
Indicates the use of a Flather open boundary.
Definition at line 61 of file MOM_open_boundary.F90.
integer, parameter, public mom_open_boundary::obc_none = 0 |
Indicates the use of no open boundary.
Definition at line 58 of file MOM_open_boundary.F90.
Referenced by open_boundary_config(), and open_boundary_impose_land_mask().
integer, parameter, public mom_open_boundary::obc_radiation = 4 |
Indicates the use of a radiation open boundary.
Definition at line 62 of file MOM_open_boundary.F90.
integer, parameter, public mom_open_boundary::obc_simple = 1 |
Indicates the use of a simple inflow open boundary.
Definition at line 59 of file MOM_open_boundary.F90.
integer, parameter, public mom_open_boundary::obc_wall = 2 |
Indicates the use of a closed wall.
Definition at line 60 of file MOM_open_boundary.F90.