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_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... | |
Functions/Subroutines | |
subroutine, public | open_boundary_config (G, 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) |
subroutine | setup_segment_indices (G, seg, Is_obc, Ie_obc, Js_obc, Je_obc) |
subroutine | setup_u_point_obc (OBC, G, segment_str, l_seg) |
Parse an OBC_SEGMENT_%%% string starting with "I=" and configure placement and type of OBC accordingly. More... | |
subroutine | setup_v_point_obc (OBC, G, segment_str, l_seg) |
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) |
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, public | open_boundary_init (G, param_file, OBC) |
Initialize open boundary control structure. 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) |
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. More... | |
subroutine, public | radiation_open_bdry_conds (OBC, u_new, u_old, v_new, v_old, G, 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 and h open boundary conditions. Also allocates and fills the segmentT and segmentS arrays, but they are not yet used anywhere. More... | |
integer function | lookup_seg_field (OBC_seg, field) |
subroutine | allocate_obc_segment_data (OBC, segment) |
Allocate 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, 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, 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... | |
Variables | |
integer, parameter, public | obc_none = 0 |
integer, parameter, public | obc_simple = 1 |
integer, parameter, public | obc_wall = 2 |
integer, parameter, public | obc_flather = 3 |
integer, parameter, public | obc_radiation = 4 |
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 |
character(len=40) | mdl = "MOM_open_boundary" |
|
private |
Allocate segment data fields.
obc | Open boundary structure | |
[in,out] | segment | Open boundary segment |
Definition at line 1670 of file MOM_open_boundary.F90.
Referenced by setup_u_point_obc(), and setup_v_point_obc().
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 2079 of file MOM_open_boundary.F90.
Referenced by mom_boundary_update::obc_register_end().
|
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 |
[in] | vvel | meridional velocity |
Definition at line 1455 of file MOM_open_boundary.F90.
References obc_direction_e, and obc_direction_n.
Referenced by radiation_open_bdry_conds().
|
private |
[in] | g | Ocean grid structure |
[in,out] | obc | Open boundary control structure |
[in] | pf | Parameter file handle |
Definition at line 331 of file MOM_open_boundary.F90.
References mdl, mom_error_handler::mom_error(), and parse_segment_data_str().
Referenced by open_boundary_config().
|
private |
Definition at line 1649 of file MOM_open_boundary.F90.
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 2032 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( 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 | ||
) |
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 |
[in,out] | v | v field to update on open boundaries |
Definition at line 1389 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(in) | G, |
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] | g | Ocean grid structure |
[in] | param_file | Parameter file handle |
obc | Open boundary control structure |
Definition at line 216 of file MOM_open_boundary.F90.
References initialize_segment_data(), mdl, mom_error_handler::mom_error(), obc_none, open_boundary_dealloc(), 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 1003 of file MOM_open_boundary.F90.
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 1013 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(in) | G, | ||
real, dimension( g %isdb: g %iedb, g %jsd: g %jed), intent(inout) | areaCu, | ||
real, dimension( g %isd: g %ied, g %jsdb: g %jedb), intent(inout) | areaCv | ||
) |
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.
obc | Open boundary control structure | |
[in] | g | Ocean grid structure |
[in,out] | areacu | Area of a u-cell (m2) |
[in,out] | areacv | Area of a u-cell (m2) |
Definition at line 1063 of file MOM_open_boundary.F90.
References obc_direction_e, obc_direction_s, 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( g %isd: g %ied, g %jsd: g %jed), 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 1019 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(param_file_type), intent(in) | param_file, | ||
type(ocean_obc_type), pointer | OBC | ||
) |
Initialize open boundary control structure.
[in] | g | Ocean grid structure |
[in] | param_file | Parameter file handle |
obc | Open boundary control structure |
Definition at line 947 of file MOM_open_boundary.F90.
References id_clock_pass, and mdl.
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 | If present, returns True if specified_*_BCs_exist_globally is true |
[in] | apply_specified_obc | If present, returns True if specified_*_BCs_exist_globally is true |
[in] | apply_flather_obc | If present, returns True if Flather_*_BCs_exist_globally is true |
[in] | apply_nudged_obc | If present, returns True if nudged_*_BCs_exist_globally is true |
[in] | needs_ext_seg_data | If present, returns True if external segment data needed |
Definition at line 981 of file MOM_open_boundary.F90.
Referenced by mom_state_initialization::mom_initialize_state().
subroutine, public mom_open_boundary::open_boundary_test_extern_h | ( | type(ocean_grid_type), intent(in) | G, |
type(ocean_obc_type), pointer | OBC, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %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 |
obc | Open boundary structure | |
[in,out] | h | Layer thickness (m or kg/m2) |
Definition at line 1766 of file MOM_open_boundary.F90.
References obc_direction_e, and obc_direction_n.
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 (m/s) |
[in,out] | v | Meridional velocity (m/s) |
Definition at line 1722 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( 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 | ||
) |
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 1423 of file MOM_open_boundary.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().
|
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 |
Definition at line 863 of file MOM_open_boundary.F90.
References mom_string_functions::extract_word(), and mom_error_handler::mom_error().
Referenced by initialize_segment_data().
|
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 |
Definition at line 771 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( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(inout) | u_new, | ||
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in) | u_old, | ||
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(inout) | v_new, | ||
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in) | v_old, | ||
type(ocean_grid_type), intent(inout) | G, | ||
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 | New u values on open boundaries |
[in] | u_old | Original unadjusted u |
[in,out] | v_new | New v values on open boundaries |
[in] | v_old | Original unadjusted v |
[in] | dt | Appropriate timestep |
Definition at line 1151 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 2058 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 2006 of file MOM_open_boundary.F90.
References obc_registry_init().
Referenced by register_file_obc().
subroutine, public mom_open_boundary::set_tracer_data | ( | type(ocean_obc_type), pointer | OBC, |
type(thermo_var_ptrs), intent(inout) | tv, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), 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 and h open boundary conditions. Also allocates and fills the segmentT and segmentS arrays, but they are not yet used anywhere.
[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 1518 of file MOM_open_boundary.F90.
References obc_direction_e, obc_direction_n, obc_direction_s, and obc_direction_w.
Referenced by mom_state_initialization::mom_initialize_state().
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 | ||
) |
[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 503 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] | segment_str | A string in form of "I=%,J=%:%,string" |
[in] | l_seg | which segment is this? |
Definition at line 565 of file MOM_open_boundary.F90.
References allocate_obc_segment_data(), 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] | segment_str | A string in form of "J=%,I=%:%,string" |
[in] | l_seg | which segment is this? |
Definition at line 668 of file MOM_open_boundary.F90.
References allocate_obc_segment_data(), 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(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 |
obc | Open boundary structure | |
[in] | tv | Thermodynamics structure |
[in,out] | h | Thickness |
Definition at line 1807 of file MOM_open_boundary.F90.
References obc_direction_s, and obc_direction_w.
integer mom_open_boundary::id_clock_pass |
Definition at line 200 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 53 of file MOM_open_boundary.F90.
|
private |
Definition at line 202 of file MOM_open_boundary.F90.
Referenced by initialize_segment_data(), open_boundary_config(), and open_boundary_init().
integer, parameter, public mom_open_boundary::obc_direction_e = 300 |
Indicates the boundary is an effective eastern boundary.
Definition at line 51 of file MOM_open_boundary.F90.
Referenced by mom_continuity_ppm::continuity_ppm(), gradient_at_q_points(), open_boundary_impose_land_mask(), open_boundary_impose_normal_slope(), open_boundary_test_extern_h(), open_boundary_test_extern_uv(), radiation_open_bdry_conds(), set_tracer_data(), setup_u_point_obc(), 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 49 of file MOM_open_boundary.F90.
Referenced by mom_barotropic::apply_eta_obcs(), mom_barotropic::apply_velocity_obcs(), mom_barotropic::btcalc(), mom_barotropic::btstep(), mom_continuity_ppm::continuity_ppm(), mom_coriolisadv::coradcalc(), mom_vert_friction::find_coupling_coef(), gradient_at_q_points(), mom_hor_visc::horizontal_viscosity(), open_boundary_impose_normal_slope(), open_boundary_test_extern_h(), open_boundary_test_extern_uv(), radiation_open_bdry_conds(), set_tracer_data(), mom_barotropic::set_up_bt_obc(), setup_v_point_obc(), and mom_vert_friction::vertvisc_coef().
integer, parameter, public mom_open_boundary::obc_direction_s = 200 |
Indicates the boundary is an effective southern boundary.
Definition at line 50 of file MOM_open_boundary.F90.
Referenced by mom_barotropic::apply_eta_obcs(), mom_barotropic::apply_velocity_obcs(), mom_barotropic::btcalc(), mom_barotropic::btstep(), mom_continuity_ppm::continuity_ppm(), mom_vert_friction::find_coupling_coef(), mom_hor_visc::horizontal_viscosity(), open_boundary_impose_land_mask(), open_boundary_impose_normal_slope(), radiation_open_bdry_conds(), set_tracer_data(), mom_barotropic::set_up_bt_obc(), setup_v_point_obc(), update_obc_segment_data(), and mom_vert_friction::vertvisc_coef().
integer, parameter, public mom_open_boundary::obc_direction_w = 400 |
Indicates the boundary is an effective western boundary.
Definition at line 52 of file MOM_open_boundary.F90.
Referenced by mom_continuity_ppm::continuity_ppm(), mom_vert_friction::find_coupling_coef(), open_boundary_impose_normal_slope(), radiation_open_bdry_conds(), set_tracer_data(), setup_u_point_obc(), update_obc_segment_data(), mom_vert_friction::vertvisc_coef(), and mom_continuity_ppm::zonal_mass_flux().
integer, parameter, public mom_open_boundary::obc_flather = 3 |
Definition at line 47 of file MOM_open_boundary.F90.
integer, parameter, public mom_open_boundary::obc_none = 0 |
Definition at line 46 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 |
Definition at line 48 of file MOM_open_boundary.F90.
integer, parameter, public mom_open_boundary::obc_simple = 1 |
Definition at line 46 of file MOM_open_boundary.F90.
integer, parameter, public mom_open_boundary::obc_wall = 2 |
Definition at line 46 of file MOM_open_boundary.F90.