MOM6
bfb_initialization Module Reference

Functions/Subroutines

subroutine, public bfb_set_coord (Rlay, g_prime, GV, param_file, eqn_of_state)
 
subroutine, public bfb_initialize_sponges_southonly (G, use_temperature, tv, param_file, CSp, h)
 
subroutine write_bfb_log (param_file)
 Write output about the parameter values being used. More...
 

Variables

logical first_call = .true.
 

Function/Subroutine Documentation

◆ bfb_initialize_sponges_southonly()

subroutine, public bfb_initialization::bfb_initialize_sponges_southonly ( type(ocean_grid_type), intent(in)  G,
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(nimem_, njmem_, nkmem_), intent(in)  h 
)
Parameters
[in]gThe ocean's grid structure
[in]tvA structure pointing to various thermodynamic variables
[in]param_fileA structure to parse for run-time parameters
[in]hLayer thicknesses, in H (usually m or kg m-2)

Definition at line 107 of file BFB_initialization.F90.

References first_call, mom_sponge::initialize_sponge(), and write_bfb_log().

Referenced by mom_state_initialization::mom_initialize_state().

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 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ bfb_set_coord()

subroutine, public bfb_initialization::bfb_set_coord ( real, dimension( : ), intent(out)  Rlay,
real, dimension( : ), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state 
)
Parameters
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters

Definition at line 66 of file BFB_initialization.F90.

References first_call, and write_bfb_log().

Referenced by mom_coord_initialization::mom_initialize_coord().

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 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 198 of file BFB_initialization.F90.

References first_call.

Referenced by bfb_initialize_sponges_southonly(), and bfb_set_coord().

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 
Here is the caller graph for this function:

Variable Documentation

◆ first_call

logical bfb_initialization::first_call = .true.
private

Definition at line 61 of file BFB_initialization.F90.

Referenced by bfb_initialize_sponges_southonly(), bfb_set_coord(), and write_bfb_log().

61 logical :: first_call = .true.