MOM6
BFB_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 !********+*********+*********+*********+*********+*********+*********+**
23 !* *
24 !* By Robert Hallberg, April 1994 - June 2002 *
25 !* *
26 !* This subroutine initializes the fields for the simulations. *
27 !* The one argument passed to initialize, Time, is set to the *
28 !* current time of the simulation. The fields which are initialized *
29 !* here are: *
30 !* G%g_prime - The reduced gravity at each interface, in m s-2. *
31 !* G%Rlay - Layer potential density (coordinate variable) in kg m-3.*
32 !* If SPONGE is defined: *
33 !* A series of subroutine calls are made to set up the damping *
34 !* rates and reference profiles for all variables that are damped *
35 !* in the sponge. *
36 !* *
37 !* These variables are all set in the set of subroutines (in this *
38 !* file) BFB_initialize_sponges_southonly and BFB_set_coord. *
39 !* *
40 !********+*********+*********+*********+*********+*********+*********+**
41 
44 use mom_get_input, only : directories
45 use mom_grid, only : ocean_grid_type
46 use mom_io, only : close_file, create_file, fieldtype, file_exists
47 use mom_io, only : open_file, read_data, read_axis_data, single_file
48 use mom_io, only : write_field, slasher
54 implicit none ; private
55 
56 #include <MOM_memory.h>
57 
58 public bfb_set_coord
60 
61 logical :: first_call = .true.
62 
63 contains
64 
65 subroutine bfb_set_coord(Rlay, g_prime, GV, param_file, eqn_of_state)
66 ! This subroutine specifies the vertical coordinate in terms of temperature at the surface and at the bottom. This case is set up in
67 ! such a way that the temperature of the topmost layer is equal to the SST at the southern edge of the domain. The temperatures are
68 ! then converted to densities of the top and bottom layers and linearly interpolated for the intermediate layers.
69  real, dimension(NKMEM_), intent(out) :: Rlay, g_prime
70  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
71  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
72  type(eos_type), pointer :: eqn_of_state
73  real :: drho_dt, SST_s, T_bot, rho_top, rho_bot
74  integer :: k, nz
75  character(len=40) :: mdl = "BFB_set_coord" ! This subroutine's name.
76 
77  call get_param(param_file, mdl, "DRHO_DT", drho_dt, &
78  "Rate of change of density with temperature.", &
79  units="kg m-3 K-1", default=-0.2)
80  call get_param(param_file, mdl, "SST_S", sst_s, &
81  "SST at the suothern edge of the domain.", units="C", default=20.0)
82  call get_param(param_file, mdl, "T_BOT", t_bot, &
83  "Bottom Temp", units="C", default=5.0)
84  rho_top = gv%rho0 + drho_dt*sst_s
85  rho_bot = gv%rho0 + drho_dt*t_bot
86  nz = gv%ke
87 
88  !call MOM_error(FATAL, &
89  ! "BFB_initialization.F90, BFB_set_coord: " // &
90  ! "Unmodified user routine called - you must edit the routine to use it")
91  do k = 1,nz
92  rlay(k) = (rho_bot - rho_top)/(nz-1)*real(k-1) + rho_top
93  if (k >1) then
94  g_prime(k) = (rlay(k) - rlay(k-1))*gv%g_earth/gv%rho0
95  else
96  g_prime(k) = gv%g_earth
97  end if
98  !Rlay(:) = 0.0
99  !g_prime(:) = 0.0
100  end do
101 
102  if (first_call) call write_bfb_log(param_file)
103 
104 end subroutine bfb_set_coord
105 
106 subroutine bfb_initialize_sponges_southonly(G, use_temperature, tv, param_file, CSp, h)
107 ! This subroutine sets up the sponges for the southern bouundary of the domain. Maximum damping occurs within 2 degrees lat of the
108 ! boundary. The damping linearly decreases northward over the next 2 degrees.
109  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
110  logical, intent(in) :: use_temperature
111  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
112  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
113  type(sponge_cs), pointer :: CSp
114  real, dimension(NIMEM_, NJMEM_, NKMEM_), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
115  !call MOM_error(FATAL, &
116  ! "BFB_initialization.F90, BFB_initialize_sponges: " // &
117  ! "Unmodified user routine called - you must edit the routine to use it")
118 
119  real :: eta(szi_(g),szj_(g),szk_(g)+1) ! A temporary array for eta.
120  real :: Idamp(szi_(g),szj_(g)) ! The inverse damping rate, in s-1.
121 
122  real :: H0(szk_(g))
123  real :: min_depth
124  real :: damp, e_dense, damp_new, slat, wlon, lenlat, lenlon, nlat
125  character(len=40) :: mdl = "BFB_initialize_sponges_southonly" ! This subroutine's name.
126  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
127 
128  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
129  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
130 
131  eta(:,:,:) = 0.0 ; idamp(:,:) = 0.0
132 
133 ! Here the inverse damping time, in s-1, is set. Set Idamp to 0 !
134 ! wherever there is no sponge, and the subroutines that are called !
135 ! will automatically set up the sponges only where Idamp is positive!
136 ! and mask2dT is 1. !
137 
138 ! Set up sponges for DOME configuration
139  call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, &
140  "The minimum depth of the ocean.", units="m", default=0.0)
141 
142  call get_param(param_file, mdl, "SOUTHLAT", slat, &
143  "The southern latitude of the domain.", units="degrees")
144  call get_param(param_file, mdl, "LENLAT", lenlat, &
145  "The latitudinal length of the domain.", units="degrees")
146  call get_param(param_file, mdl, "WESTLON", wlon, &
147  "The western longitude of the domain.", units="degrees", default=0.0)
148  call get_param(param_file, mdl, "LENLON", lenlon, &
149  "The longitudinal length of the domain.", units="degrees")
150  nlat = slat + lenlat
151  do k=1,nz ; h0(k) = -g%max_depth * real(k-1) / real(nz) ; enddo
152 ! do k=1,nz ; H0(k) = -G%max_depth * real(k-1) / real(nz-1) ; enddo ! Use for meridional thickness profile initialization
153  do i=is,ie; do j=js,je
154  if (g%geoLatT(i,j) < slat+2.0) then ; damp = 1.0
155  elseif (g%geoLatT(i,j) < slat+4.0) then
156  damp_new = 1.0*(slat+4.0-g%geoLatT(i,j))/2.0
157  else ; damp = 0.0
158  endif
159 
160  ! These will be streched inside of apply_sponge, so they can be in
161  ! depth space for Boussinesq or non-Boussinesq models.
162 
163  ! This section is used for uniform thickness initialization
164  do k = 1,nz; eta(i,j,k) = h0(k); enddo
165 
166  ! The below section is used for meridional temperature profile thickness initiation
167  ! do k = 1,nz; eta(i,j,k) = H0(k); enddo
168  ! if (G%geoLatT(i,j) > 40.0) then
169  ! do k = 1,nz
170  ! eta(i,j,k) = -G%Angstrom_z*(k-1)
171  ! enddo
172  ! elseif (G%geoLatT(i,j) > 20.0) then
173  ! do k = 1,nz
174  ! eta(i,j,k) = min(H0(k) + (G%geoLatT(i,j) - 20.0)*(G%max_depth - nz*G%Angstrom_z)/20.0, -(k-1)*G%angstrom_z)
175  ! enddo
176  ! endif
177  eta(i,j,nz+1) = -g%max_depth
178 
179  if (g%bathyT(i,j) > min_depth) then
180  idamp(i,j) = damp/86400.0
181  else ; idamp(i,j) = 0.0 ; endif
182  enddo ; enddo
183 
184 ! This call sets up the damping rates and interface heights.
185 ! This sets the inverse damping timescale fields in the sponges. !
186  call initialize_sponge(idamp, eta, g, param_file, csp)
187 
188 ! Now register all of the fields which are damped in the sponge. !
189 ! By default, momentum is advected vertically within the sponge, but !
190 ! momentum is typically not damped within the sponge. !
191 
192  if (first_call) call write_bfb_log(param_file)
193 
195 
196 !> Write output about the parameter values being used.
197 subroutine write_bfb_log(param_file)
198  type(param_file_type), intent(in) :: param_file !< A structure indicating the
199  !! open file to parse for model
200  !! parameter values.
201 
202 ! This include declares and sets the variable "version".
203 #include "version_variable.h"
204  character(len=40) :: mdl = "BFB_initialization" ! This module's name.
205 
206  call log_version(param_file, mdl, version)
207  first_call = .false.
208 
209 end subroutine write_bfb_log
210 
211 end module bfb_initialization
subroutine, public read_axis_data(filename, axis_name, var)
Definition: MOM_io.F90:424
subroutine write_bfb_log(param_file)
Write output about the parameter values being used.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:45
Provides the ocean grid type.
Definition: MOM_grid.F90:2
This module contains I/O framework code.
Definition: MOM_io.F90:2
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
subroutine, public calculate_density_derivs(T, S, pressure, drho_dT, drho_dS, start, npts, EOS)
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
Definition: MOM_EOS.F90:214
subroutine, public initialize_sponge(Iresttime, int_height, G, param_file, CS, Iresttime_i_mean, int_height_i_mean)
Definition: MOM_sponge.F90:142
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
Type to carry basic tracer information.
logical function, public is_root_pe()
subroutine, public add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
This subroutine adds open boundary condition concentrations for a tracer that has previously been reg...
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
subroutine, public bfb_set_coord(Rlay, g_prime, GV, param_file, eqn_of_state)
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
subroutine, public create_file(unit, filename, vars, novars, fields, threading, timeunit, G, dG, GV)
Routine creates a new NetCDF file. It also sets up structures that describe this file and variables t...
Definition: MOM_io.F90:82
subroutine, public mom_error(level, message, all_print)
subroutine, public bfb_initialize_sponges_southonly(G, use_temperature, tv, param_file, CSp, h)
A control structure for the equation of state.
Definition: MOM_EOS.F90:55