MOM6
tidal_bay_initialization.F90
Go to the documentation of this file.
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !***********************************************************************
21 
22 use mom_coms, only : reproducing_sum
24 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
26 use mom_grid, only : ocean_grid_type
31 use mom_time_manager, only : time_type, set_time, time_type_to_real
32 
33 implicit none ; private
34 
35 #include <MOM_memory.h>
36 
39 
40 !> Control structure for tidal bay open boundaries.
41 type, public :: tidal_bay_obc_cs ; private
42  real :: tide_flow = 3.0e6 !< Maximum tidal flux.
43 end type tidal_bay_obc_cs
44 
45 contains
46 
47 !> Add tidal bay to OBC registry.
48 function register_tidal_bay_obc(param_file, CS, OBC_Reg)
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 
66 end function register_tidal_bay_obc
67 
68 !> Clean up the tidal bay OBC from registry.
69 subroutine tidal_bay_obc_end(CS)
70  type(tidal_bay_obc_cs), pointer :: CS !< tidal bay control structure.
71 
72  if (associated(cs)) then
73  deallocate(cs)
74  endif
75 end subroutine tidal_bay_obc_end
76 
77 !> This subroutine sets the properties of flow at open boundary conditions.
78 subroutine tidal_bay_set_obc_data(OBC, CS, G, h, Time)
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 
133 end subroutine tidal_bay_set_obc_data
134 
135 !> \namespace tidal_bay_initialization
136 !!
137 !! The module configures the model for the "tidal_bay" experiment.
138 !! tidal_bay = Tidally resonant bay from Zygmunt Kowalik's class on tides.
139 end module tidal_bay_initialization
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.
Definition: MOM_grid.F90:19
Provides the ocean grid type.
Definition: MOM_grid.F90:2
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.