MOM6
bfb_initialization Module Reference

Detailed Description

Initialization of the boundary-forced-basing configuration.

Functions/Subroutines

subroutine, public bfb_set_coord (Rlay, g_prime, GV, US, param_file, eqn_of_state)
 This subroutine specifies the vertical coordinate in terms of temperature at the surface and at the bottom. This case is set up in 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 then converted to densities of the top and bottom layers and linearly interpolated for the intermediate layers. More...
 
subroutine, public bfb_initialize_sponges_southonly (G, GV, US, use_temperature, tv, param_file, CSp, h)
 This subroutine sets up the sponges for the southern bouundary of the domain. Maximum damping occurs within 2 degrees lat of the boundary. The damping linearly decreases northward over the next 2 degrees. More...
 
subroutine write_bfb_log (param_file)
 Write output about the parameter values being used. More...
 

Variables

logical first_call = .true.
 Unsafe model variable. More...
 

Function/Subroutine Documentation

◆ bfb_initialize_sponges_southonly()

subroutine, public bfb_initialization::bfb_initialize_sponges_southonly ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
logical, intent(in)  use_temperature,
type(thermo_var_ptrs), intent(in)  tv,
type(param_file_type), intent(in)  param_file,
type(sponge_cs), pointer  CSp,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h 
)

This subroutine sets up the sponges for the southern bouundary of the domain. Maximum damping occurs within 2 degrees lat of the boundary. The damping linearly decreases northward over the next 2 degrees.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[in]use_temperatureIf true, temperature and salinity are used as state variables.
[in]tvA structure pointing to various thermodynamic variables
[in]param_fileA structure to parse for run-time parameters
cspA pointer to the sponge control structure
[in]hLayer thicknesses [H ~> m or kg m-2]

Definition at line 79 of file BFB_initialization.F90.

80  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
81  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
82  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
83  logical, intent(in) :: use_temperature !< If true, temperature and salinity are used as
84  !! state variables.
85  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
86  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
87  type(sponge_CS), pointer :: CSp !< A pointer to the sponge control structure
88  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
89  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
90 
91  ! Local variables
92  real :: eta(SZI_(G),SZJ_(G),SZK_(GV)+1) ! A temporary array for eta, in depth units [Z ~> m].
93  real :: Idamp(SZI_(G),SZJ_(G)) ! The inverse damping rate [T-1 ~> s-1].
94  real :: H0(SZK_(GV)) ! Resting layer thicknesses in depth units [Z ~> m].
95  real :: min_depth ! The minimum ocean depth in depth units [Z ~> m].
96  real :: slat, wlon, lenlat, lenlon, nlat
97  real :: max_damping ! The maximum damping rate [T-1 ~> s-1]
98  character(len=40) :: mdl = "BFB_initialize_sponges_southonly" ! This subroutine's name.
99  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
100 
101  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
102  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
103 
104  eta(:,:,:) = 0.0 ; idamp(:,:) = 0.0
105 
106 ! Here the inverse damping time [T-1 ~> s-1], is set. Set Idamp to 0
107 ! wherever there is no sponge, and the subroutines that are called
108 ! will automatically set up the sponges only where Idamp is positive
109 ! and mask2dT is 1.
110 
111 ! Set up sponges for DOME configuration
112  call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, &
113  "The minimum depth of the ocean.", units="m", default=0.0, scale=us%m_to_Z)
114 
115  call get_param(param_file, mdl, "SOUTHLAT", slat, &
116  "The southern latitude of the domain.", units="degrees")
117  call get_param(param_file, mdl, "LENLAT", lenlat, &
118  "The latitudinal length of the domain.", units="degrees")
119  call get_param(param_file, mdl, "WESTLON", wlon, &
120  "The western longitude of the domain.", units="degrees", default=0.0)
121  call get_param(param_file, mdl, "LENLON", lenlon, &
122  "The longitudinal length of the domain.", units="degrees")
123  nlat = slat + lenlat
124  do k=1,nz ; h0(k) = -g%max_depth * real(k-1) / real(nz) ; enddo
125 
126  ! Use for meridional thickness profile initialization
127 ! do k=1,nz ; H0(k) = -G%max_depth * real(k-1) / real(nz-1) ; enddo
128 
129  max_damping = 1.0 / (86400.0*us%s_to_T)
130 
131  do i=is,ie; do j=js,je
132  if (g%bathyT(i,j) <= min_depth) then ; idamp(i,j) = 0.0
133  elseif (g%geoLatT(i,j) < slat+2.0) then ; idamp(i,j) = max_damping
134  elseif (g%geoLatT(i,j) < slat+4.0) then
135  idamp(i,j) = max_damping * (slat+4.0-g%geoLatT(i,j))/2.0
136  else ; idamp(i,j) = 0.0
137  endif
138 
139  ! These will be streched inside of apply_sponge, so they can be in
140  ! depth space for Boussinesq or non-Boussinesq models.
141 
142  ! This section is used for uniform thickness initialization
143  do k = 1,nz; eta(i,j,k) = h0(k); enddo
144 
145  ! The below section is used for meridional temperature profile thickness initiation
146  ! do k = 1,nz; eta(i,j,k) = H0(k); enddo
147  ! if (G%geoLatT(i,j) > 40.0) then
148  ! do k = 1,nz
149  ! eta(i,j,k) = -G%Angstrom_Z*(k-1)
150  ! enddo
151  ! elseif (G%geoLatT(i,j) > 20.0) then
152  ! do k = 1,nz
153  ! eta(i,j,k) = min(H0(k) + (G%geoLatT(i,j) - 20.0)*(G%max_depth - nz*G%Angstrom_Z)/20.0, &
154  ! -(k-1)*G%Angstrom_Z)
155  ! enddo
156  ! endif
157  eta(i,j,nz+1) = -g%max_depth
158 
159  enddo ; enddo
160 
161 ! This call sets up the damping rates and interface heights.
162 ! This sets the inverse damping timescale fields in the sponges. !
163  call initialize_sponge(idamp, eta, g, param_file, csp, gv)
164 
165 ! Now register all of the fields which are damped in the sponge. !
166 ! By default, momentum is advected vertically within the sponge, but !
167 ! momentum is typically not damped within the sponge. !
168 
169  if (first_call) call write_bfb_log(param_file)
170 

◆ bfb_set_coord()

subroutine, public bfb_initialization::bfb_set_coord ( real, dimension(gv%ke), intent(out)  Rlay,
real, dimension(gv%ke+1), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state 
)

This subroutine specifies the vertical coordinate in terms of temperature at the surface and at the bottom. This case is set up in 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 then converted to densities of the top and bottom layers and linearly interpolated for the intermediate layers.

Parameters
[in]gvThe ocean's vertical grid structure
[out]rlayLayer potential density [R ~> kg m-3].
[out]g_primeThe reduced gravity at each interface [L2 Z-1 T-2 ~> m s-2].
[in]usA dimensional unit scaling type
[in]param_fileA structure to parse for run-time parameters
eqn_of_stateEquation of state structure

Definition at line 38 of file BFB_initialization.F90.

39  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
40  real, dimension(GV%ke), intent(out) :: Rlay !< Layer potential density [R ~> kg m-3].
41  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity at each
42  !! interface [L2 Z-1 T-2 ~> m s-2].
43  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
44  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
45  type(EOS_type), pointer :: eqn_of_state !< Equation of state structure
46  ! Local variables
47  real :: drho_dt, SST_s, T_bot, rho_top, rho_bot
48  integer :: k, nz
49  character(len=40) :: mdl = "BFB_set_coord" ! This subroutine's name.
50 
51  call get_param(param_file, mdl, "DRHO_DT", drho_dt, &
52  "Rate of change of density with temperature.", &
53  units="kg m-3 K-1", default=-0.2, scale=us%kg_m3_to_R)
54  call get_param(param_file, mdl, "SST_S", sst_s, &
55  "SST at the suothern edge of the domain.", units="C", default=20.0)
56  call get_param(param_file, mdl, "T_BOT", t_bot, &
57  "Bottom Temp", units="C", default=5.0)
58  rho_top = gv%Rho0 + drho_dt*sst_s
59  rho_bot = gv%Rho0 + drho_dt*t_bot
60  nz = gv%ke
61 
62  do k = 1,nz
63  rlay(k) = (rho_bot - rho_top)/(nz-1)*real(k-1) + rho_top
64  if (k >1) then
65  g_prime(k) = (rlay(k) - rlay(k-1)) * gv%g_Earth / (gv%Rho0)
66  else
67  g_prime(k) = gv%g_Earth
68  endif
69  !Rlay(:) = 0.0
70  !g_prime(:) = 0.0
71  enddo
72 
73  if (first_call) call write_bfb_log(param_file)
74 

◆ write_bfb_log()

subroutine bfb_initialization::write_bfb_log ( type(param_file_type), intent(in)  param_file)
private

Write output about the parameter values being used.

Parameters
[in]param_fileA structure indicating the open file to parse for model parameter values.

Definition at line 174 of file BFB_initialization.F90.

175  type(param_file_type), intent(in) :: param_file !< A structure indicating the
176  !! open file to parse for model
177  !! parameter values.
178 
179 ! This include declares and sets the variable "version".
180 #include "version_variable.h"
181  character(len=40) :: mdl = "BFB_initialization" ! This module's name.
182 
183  call log_version(param_file, mdl, version)
184  first_call = .false.
185 

Variable Documentation

◆ first_call

logical bfb_initialization::first_call = .true.
private

Unsafe model variable.

Todo:
Remove this module variable

Definition at line 30 of file BFB_initialization.F90.

30 logical :: first_call = .true.