MOM6
tidal_bay_initialization Module Reference

Detailed Description

Configures the model for the "tidal_bay" experiment. tidal_bay = Tidally resonant bay from Zygmunt Kowalik's class on tides.

Data Types

type  tidal_bay_obc_cs
 Control structure for tidal bay open boundaries. More...
 

Functions/Subroutines

logical function, public register_tidal_bay_obc (param_file, CS, OBC_Reg)
 Add tidal bay to OBC registry. More...
 
subroutine, public tidal_bay_obc_end (CS)
 Clean up the tidal bay OBC from registry. More...
 
subroutine, public tidal_bay_set_obc_data (OBC, CS, G, h, Time)
 This subroutine sets the properties of flow at open boundary conditions. More...
 

Function/Subroutine Documentation

◆ register_tidal_bay_obc()

logical function, public tidal_bay_initialization::register_tidal_bay_obc ( type(param_file_type), intent(in)  param_file,
type(tidal_bay_obc_cs), pointer  CS,
type(obc_registry_type), pointer  OBC_Reg 
)

Add tidal bay to OBC registry.

Parameters
[in]param_fileparameter file.
cstidal bay control structure.
obc_regOBC registry.

Definition at line 34 of file tidal_bay_initialization.F90.

34  type(param_file_type), intent(in) :: param_file !< parameter file.
35  type(tidal_bay_OBC_CS), pointer :: CS !< tidal bay control structure.
36  type(OBC_registry_type), pointer :: OBC_Reg !< OBC registry.
37  logical :: register_tidal_bay_OBC
38  character(len=32) :: casename = "tidal bay" !< This case's name.
39 
40  if (associated(cs)) then
41  call mom_error(warning, "register_tidal_bay_OBC called with an "// &
42  "associated control structure.")
43  return
44  endif
45  allocate(cs)
46 
47  ! Register the open boundaries.
48  call register_obc(casename, param_file, obc_reg)
49  register_tidal_bay_obc = .true.
50 

References mom_error_handler::mom_error().

Here is the call graph for this function:

◆ tidal_bay_obc_end()

subroutine, public tidal_bay_initialization::tidal_bay_obc_end ( type(tidal_bay_obc_cs), pointer  CS)

Clean up the tidal bay OBC from registry.

Parameters
cstidal bay control structure.

Definition at line 55 of file tidal_bay_initialization.F90.

55  type(tidal_bay_OBC_CS), pointer :: CS !< tidal bay control structure.
56 
57  if (associated(cs)) then
58  deallocate(cs)
59  endif

Referenced by mom_boundary_update::obc_register_end().

Here is the caller graph for this function:

◆ tidal_bay_set_obc_data()

subroutine, public tidal_bay_initialization::tidal_bay_set_obc_data ( type(ocean_obc_type), pointer  OBC,
type(tidal_bay_obc_cs), pointer  CS,
type(ocean_grid_type), intent(in)  G,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(time_type), intent(in)  Time 
)

This subroutine sets the properties of flow at open boundary conditions.

Parameters
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used.
cstidal bay control structure.
[in]gThe ocean's grid structure.
[in]hlayer thickness.
[in]timemodel time.

Definition at line 64 of file tidal_bay_initialization.F90.

64  type(ocean_OBC_type), pointer :: OBC !< This open boundary condition type specifies
65  !! whether, where, and what open boundary
66  !! conditions are used.
67  type(tidal_bay_OBC_CS), pointer :: CS !< tidal bay control structure.
68  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
69  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< layer thickness.
70  type(time_type), intent(in) :: Time !< model time.
71 
72  ! The following variables are used to set up the transport in the tidal_bay example.
73  real :: time_sec, cff
74  real :: my_flux, total_area
75  real :: PI
76  real, allocatable :: my_area(:,:)
77  character(len=40) :: mdl = "tidal_bay_set_OBC_data" ! This subroutine's name.
78  integer :: i, j, k, itt, is, ie, js, je, isd, ied, jsd, jed, nz, n
79  integer :: IsdB, IedB, JsdB, JedB
80  type(OBC_segment_type), pointer :: segment => null()
81 
82  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
83  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
84  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
85 
86  pi = 4.0*atan(1.0)
87 
88  if (.not.associated(obc)) return
89 
90  allocate(my_area(1:1,js:je))
91 
92  time_sec = time_type_to_real(time)
93  cff = 0.1*sin(2.0*pi*time_sec/(12.0*3600.0))
94  my_area=0.0
95  my_flux=0.0
96  segment => obc%segment(1)
97 
98  do j=segment%HI%jsc,segment%HI%jec ; do i=segment%HI%IscB,segment%HI%IecB
99  if (obc%segnum_u(i,j) /= obc_none) then
100  do k=1,nz
101  my_area(1,j) = my_area(1,j) + h(i,j,k)*g%US%L_to_m*g%dyCu(i,j)
102  enddo
103  endif
104  enddo ; enddo
105  total_area = reproducing_sum(my_area)
106  my_flux = - cs%tide_flow*sin(2.0*pi*time_sec/(12.0*3600.0))
107 
108  do n = 1, obc%number_of_segments
109  segment => obc%segment(n)
110 
111  if (.not. segment%on_pe) cycle
112 
113  segment%normal_vel_bt(:,:) = g%US%m_s_to_L_T*my_flux/total_area
114  segment%eta(:,:) = cff
115 
116  enddo ! end segment loop
117 
mom_coms::reproducing_sum
Find an accurate and order-invariant sum of distributed 2d or 3d fields.
Definition: MOM_coms.F90:54