MOM6
tidal_bay_initialization Module Reference

Detailed Description

The module 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 49 of file tidal_bay_initialization.F90.

References mom_error_handler::mom_error().

49  type(param_file_type), intent(in) :: param_file !< parameter file.
50  type(tidal_bay_obc_cs), pointer :: cs !< tidal bay control structure.
51  type(obc_registry_type), pointer :: obc_reg !< OBC registry.
52  logical :: register_tidal_bay_obc
53  character(len=32) :: casename = "tidal bay" !< This case's name.
54 
55  if (associated(cs)) then
56  call mom_error(warning, "register_tidal_bay_OBC called with an "// &
57  "associated control structure.")
58  return
59  endif
60  allocate(cs)
61 
62  ! Register the tracer for horizontal advection & diffusion.
63  call register_obc(casename, param_file, obc_reg)
64  register_tidal_bay_obc = .true.
65 
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 70 of file tidal_bay_initialization.F90.

Referenced by mom_boundary_update::obc_register_end().

70  type(tidal_bay_obc_cs), pointer :: cs !< tidal bay control structure.
71 
72  if (associated(cs)) then
73  deallocate(cs)
74  endif
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( g %isd: g %ied, g %jsd: g %jed, g %ke), 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 79 of file tidal_bay_initialization.F90.

79  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
80  !! whether, where, and what open boundary
81  !! conditions are used.
82  type(tidal_bay_obc_cs), pointer :: cs !< tidal bay control structure.
83  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
84  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< layer thickness.
85  type(time_type), intent(in) :: time !< model time.
86 
87  ! The following variables are used to set up the transport in the tidal_bay example.
88  real :: time_sec, cff
89  real :: my_flux, total_area
90  real :: pi
91  real, allocatable :: my_area(:,:)
92  character(len=40) :: mdl = "tidal_bay_set_OBC_data" ! This subroutine's name.
93  integer :: i, j, k, itt, is, ie, js, je, isd, ied, jsd, jed, nz, n
94  integer :: isdb, iedb, jsdb, jedb
95  type(obc_segment_type), pointer :: segment
96 
97  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
98  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
99  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
100 
101  pi = 4.0*atan(1.0) ;
102 
103  if (.not.associated(obc)) return
104 
105  allocate(my_area(1:1,js:je))
106 
107  time_sec = time_type_to_real(time)
108  cff = 0.1*sin(2.0*pi*time_sec/(12.0*3600.0))
109  my_area=0.0
110  my_flux=0.0
111  segment => obc%segment(1)
112 
113  do j=segment%HI%jsc,segment%HI%jec ; do i=segment%HI%IscB,segment%HI%IecB
114  if (obc%segnum_u(i,j) /= obc_none) then
115  do k=1,nz
116  my_area(1,j) = my_area(1,j) + h(i,j,k)*g%dyCu(i,j)
117  enddo
118  endif
119  enddo ; enddo
120  total_area = reproducing_sum(my_area)
121  my_flux = - cs%tide_flow*sin(2.0*pi*time_sec/(12.0*3600.0))
122 
123  do n = 1, obc%number_of_segments
124  segment => obc%segment(n)
125 
126  if (.not. segment%on_pe) cycle
127 
128  segment%normal_vel_bt(:,:) = my_flux/total_area
129  segment%eta(:,:) = cff
130 
131  enddo ! end segment loop
132 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19