MOM6
MOM_open_boundary.F90 File Reference
#include <MOM_memory.h>
Include dependency graph for MOM_open_boundary.F90:

Go to the source code of this file.

Data Types

type  mom_open_boundary::obc_segment_data_type
 Open boundary segment data from files (mostly). More...
 
type  mom_open_boundary::obc_segment_type
 Open boundary segment data structure. More...
 
type  mom_open_boundary::ocean_obc_type
 Open-boundary data. More...
 
type  mom_open_boundary::file_obc_cs
 Control structure for open boundaries that read from files. Probably lots to update here. More...
 
type  mom_open_boundary::obc_struct_type
 Type to carry something (what] for the OBC registry. More...
 
type  mom_open_boundary::obc_registry_type
 Type to carry basic OBC information needed for updating values. More...
 

Modules

module  mom_open_boundary
 Controls where open boundary conditions are applied.
 

Functions/Subroutines

subroutine, public mom_open_boundary::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 mom_open_boundary::initialize_segment_data (G, OBC, PF)
 
subroutine mom_open_boundary::setup_segment_indices (G, seg, Is_obc, Ie_obc, Js_obc, Je_obc)
 
subroutine mom_open_boundary::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 mom_open_boundary::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 mom_open_boundary::parse_segment_str (ni_global, nj_global, segment_str, l, m, n, action_str)
 Parse an OBC_SEGMENT_%%% string. More...
 
integer function interpret_int_expr (string, imax)
 
subroutine mom_open_boundary::parse_segment_data_str (segment_str, var, value, filenam, fieldnam, fields, num_fields, debug)
 Parse an OBC_SEGMENT_%%_DATA string. More...
 
subroutine, public mom_open_boundary::open_boundary_init (G, param_file, OBC)
 Initialize open boundary control structure. More...
 
logical function, public mom_open_boundary::open_boundary_query (OBC, apply_open_OBC, apply_specified_OBC, apply_Flather_OBC, apply_nudged_OBC, needs_ext_seg_data)
 
subroutine mom_open_boundary::open_boundary_dealloc (OBC)
 Deallocate open boundary data. More...
 
subroutine, public mom_open_boundary::open_boundary_end (OBC)
 Close open boundary data. More...
 
subroutine, public mom_open_boundary::open_boundary_impose_normal_slope (OBC, G, depth)
 Sets the slope of bathymetry normal to an open bounndary to zero. More...
 
subroutine, public mom_open_boundary::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 mom_open_boundary::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 mom_open_boundary::open_boundary_apply_normal_flow (OBC, G, u, v)
 Applies OBC values stored in segments to 3d u,v fields. More...
 
subroutine, public mom_open_boundary::open_boundary_zero_normal_flow (OBC, G, u, v)
 Applies zero values to 3d u,v fields on OBC segments. More...
 
subroutine mom_open_boundary::gradient_at_q_points (G, segment, uvel, vvel)
 Calculate the tangential gradient of the normal flow at the boundary q-points. More...
 
subroutine, public mom_open_boundary::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 mom_open_boundary::lookup_seg_field (OBC_seg, field)
 
subroutine mom_open_boundary::allocate_obc_segment_data (OBC, segment)
 Allocate segment data fields. More...
 
subroutine, public mom_open_boundary::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 mom_open_boundary::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 mom_open_boundary::update_obc_segment_data (G, GV, OBC, tv, h, Time)
 Update the OBC values on the segments. More...
 
subroutine, public mom_open_boundary::register_obc (name, param_file, Reg)
 register open boundary objects for boundary updates. More...
 
subroutine, public mom_open_boundary::obc_registry_init (param_file, Reg)
 This routine include declares and sets the variable "version". More...
 
logical function, public mom_open_boundary::register_file_obc (param_file, CS, OBC_Reg)
 Add file to OBC registry. More...
 
subroutine, public mom_open_boundary::file_obc_end (CS)
 Clean up the file OBC from registry. More...
 

Variables

integer, parameter, public mom_open_boundary::obc_none = 0
 
integer, parameter, public mom_open_boundary::obc_simple = 1
 
integer, parameter, public mom_open_boundary::obc_wall = 2
 
integer, parameter, public mom_open_boundary::obc_flather = 3
 
integer, parameter, public mom_open_boundary::obc_radiation = 4
 
integer, parameter, public mom_open_boundary::obc_direction_n = 100
 Indicates the boundary is an effective northern boundary. More...
 
integer, parameter, public mom_open_boundary::obc_direction_s = 200
 Indicates the boundary is an effective southern boundary. More...
 
integer, parameter, public mom_open_boundary::obc_direction_e = 300
 Indicates the boundary is an effective eastern boundary. More...
 
integer, parameter, public mom_open_boundary::obc_direction_w = 400
 Indicates the boundary is an effective western boundary. More...
 
integer, parameter mom_open_boundary::max_obc_fields = 100
 Maximum number of data fields needed for OBC segments. More...
 
integer mom_open_boundary::id_clock_pass
 
character(len=40) mom_open_boundary::mdl = "MOM_open_boundary"
 

Function/Subroutine Documentation

◆ interpret_int_expr()

integer function parse_segment_str::interpret_int_expr ( character(len=*), intent(in)  string,
integer, intent(in)  imax 
)
private
Parameters
[in]stringInteger in form or I, N or N-I
[in]imaxValue to replace 'N' with

Definition at line 839 of file MOM_open_boundary.F90.

References mom_error_handler::mom_error().

Referenced by mom_open_boundary::parse_segment_str().

839  character(len=*), intent(in) :: string !< Integer in form or %I, N or N-%I
840  integer, intent(in) :: imax !< Value to replace 'N' with
841  ! Local variables
842  integer slen
843 
844  slen = len_trim(string)
845  if (slen==0) call mom_error(fatal, "MOM_open_boundary.F90, parse_segment_str"//&
846  "Parsed string was empty!")
847  if (len_trim(string)==1 .and. string(1:1)=='N') then
848  interpret_int_expr = imax
849  elseif (string(1:1)=='N') then
850  read(string(2:slen),*,err=911) interpret_int_expr
852  else
853  read(string(1:slen),*,err=911) interpret_int_expr
854  endif
855  return
856  911 call mom_error(fatal, "MOM_open_boundary.F90, parse_segment_str"//&
857  "Problem reading value from string '"//trim(string)//"'.")
integer function interpret_int_expr(string, imax)
Here is the call graph for this function:
Here is the caller graph for this function: