33 implicit none ;
private 35 #include <MOM_memory.h> 42 real :: tide_flow = 3.0e6
52 logical :: register_tidal_bay_OBC
53 character(len=32) :: casename =
"tidal bay" 55 if (
associated(cs))
then 56 call mom_error(warning,
"register_tidal_bay_OBC called with an "// &
57 "associated control structure.")
63 call register_obc(casename, param_file, obc_reg)
64 register_tidal_bay_obc = .true.
72 if (
associated(cs))
then 84 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
85 type(time_type),
intent(in) :: Time
89 real :: my_flux, total_area
91 real,
allocatable :: my_area(:,:)
92 character(len=40) :: mdl =
"tidal_bay_set_OBC_data" 93 integer :: i, j, k, itt, is, ie, js, je, isd, ied, jsd, jed, nz, n
94 integer :: IsdB, IedB, JsdB, JedB
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
103 if (.not.
associated(obc))
return 105 allocate(my_area(1:1,js:je))
107 time_sec = time_type_to_real(time)
108 cff = 0.1*sin(2.0*pi*time_sec/(12.0*3600.0))
111 segment => obc%segment(1)
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 116 my_area(1,j) = my_area(1,j) + h(i,j,k)*g%dyCu(i,j)
121 my_flux = - cs%tide_flow*sin(2.0*pi*time_sec/(12.0*3600.0))
123 do n = 1, obc%number_of_segments
124 segment => obc%segment(n)
126 if (.not. segment%on_pe) cycle
128 segment%normal_vel_bt(:,:) = my_flux/total_area
129 segment%eta(:,:) = cff
subroutine, public tidal_bay_set_obc_data(OBC, CS, G, h, Time)
This subroutine sets the properties of flow at open boundary conditions.
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
Open boundary segment data structure.
subroutine, public register_obc(name, param_file, Reg)
register open boundary objects for boundary updates.
The module configures the model for the "tidal_bay" experiment. tidal_bay = Tidally resonant bay from...
Type to carry basic OBC information needed for updating values.
logical function, public is_root_pe()
subroutine, public mom_mesg(message, verb, all_print)
Control structure for tidal bay open boundaries.
subroutine, public tidal_bay_obc_end(CS)
Clean up the tidal bay OBC from registry.
integer, parameter, public obc_none
Controls where open boundary conditions are applied.
subroutine, public mom_error(level, message, all_print)
logical function, public register_tidal_bay_obc(param_file, CS, OBC_Reg)
Add tidal bay to OBC registry.