MOM6
bfb_surface_forcing Module Reference

Detailed Description

Surface forcing for the boundary-forced-basin (BFB) configuration.

Data Types

type  bfb_surface_forcing_cs
 Control structure for BFB_surface_forcing. More...
 

Functions/Subroutines

subroutine, public bfb_buoyancy_forcing (sfc_state, fluxes, day, dt, G, US, CS)
 Bouyancy forcing for the boundary-forced-basin (BFB) configuration. More...
 
subroutine, public bfb_surface_forcing_init (Time, G, US, param_file, diag, CS)
 Initialization for forcing the boundary-forced-basin (BFB) configuration. More...
 

Function/Subroutine Documentation

◆ bfb_buoyancy_forcing()

subroutine, public bfb_surface_forcing::bfb_buoyancy_forcing ( type(surface), intent(inout)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(bfb_surface_forcing_cs), pointer  CS 
)

Bouyancy forcing for the boundary-forced-basin (BFB) configuration.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesA structure containing pointers to any possible forcing fields. Unused fields have NULL ptrs.
[in]dayTime of the fluxes.
[in]dtThe amount of time over which the fluxes apply [s]
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
csA pointer to the control structure returned by a previous call to BFB_surface_forcing_init.

Definition at line 51 of file BFB_surface_forcing.F90.

51  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
52  !! describe the surface state of the ocean.
53  type(forcing), intent(inout) :: fluxes !< A structure containing pointers to any
54  !! possible forcing fields. Unused fields
55  !! have NULL ptrs.
56  type(time_type), intent(in) :: day !< Time of the fluxes.
57  real, intent(in) :: dt !< The amount of time over which
58  !! the fluxes apply [s]
59  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
60  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
61  type(bfb_surface_forcing_cs), pointer :: cs !< A pointer to the control structure
62  !! returned by a previous call to
63  !! BFB_surface_forcing_init.
64  ! Local variables
65  real :: temp_restore ! The temperature that is being restored toward [degC].
66  real :: salin_restore ! The salinity that is being restored toward [ppt].
67  real :: density_restore ! The potential density that is being restored
68  ! toward [R ~> kg m-3].
69  real :: rhoxcp ! Reference density times heat capacity times unit scaling
70  ! factors [Q R degC-1 ~> J m-3 degC-1]
71  real :: buoy_rest_const ! A constant relating density anomalies to the
72  ! restoring buoyancy flux [L2 T-3 R-1 ~> m5 s-3 kg-1].
73  integer :: i, j, is, ie, js, je
74  integer :: isd, ied, jsd, jed
75 
76  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
77  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
78 
79  ! Allocate and zero out the forcing arrays, as necessary. This portion is
80  ! usually not changed.
81  if (cs%use_temperature) then
82  call safe_alloc_ptr(fluxes%evap, isd, ied, jsd, jed)
83  call safe_alloc_ptr(fluxes%lprec, isd, ied, jsd, jed)
84  call safe_alloc_ptr(fluxes%fprec, isd, ied, jsd, jed)
85  call safe_alloc_ptr(fluxes%lrunoff, isd, ied, jsd, jed)
86  call safe_alloc_ptr(fluxes%frunoff, isd, ied, jsd, jed)
87  call safe_alloc_ptr(fluxes%vprec, isd, ied, jsd, jed)
88 
89  call safe_alloc_ptr(fluxes%sw, isd, ied, jsd, jed)
90  call safe_alloc_ptr(fluxes%lw, isd, ied, jsd, jed)
91  call safe_alloc_ptr(fluxes%latent, isd, ied, jsd, jed)
92  call safe_alloc_ptr(fluxes%sens, isd, ied, jsd, jed)
93  else ! This is the buoyancy only mode.
94  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
95  endif
96 
97  if ( cs%use_temperature ) then
98  ! Set whichever fluxes are to be used here. Any fluxes that
99  ! are always zero do not need to be changed here.
100  do j=js,je ; do i=is,ie
101  ! Fluxes of fresh water through the surface are in units of [R Z T-1 ~> kg m-2 s-1]
102  ! and are positive downward - i.e. evaporation should be negative.
103  fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
104  fluxes%lprec(i,j) = 0.0 * g%mask2dT(i,j)
105 
106  ! vprec will be set later, if it is needed for salinity restoring.
107  fluxes%vprec(i,j) = 0.0
108 
109  ! Heat fluxes are in units of [Q R Z T-1 ~> W m-2] and are positive into the ocean.
110  fluxes%lw(i,j) = 0.0 * g%mask2dT(i,j)
111  fluxes%latent(i,j) = 0.0 * g%mask2dT(i,j)
112  fluxes%sens(i,j) = 0.0 * g%mask2dT(i,j)
113  fluxes%sw(i,j) = 0.0 * g%mask2dT(i,j)
114  enddo ; enddo
115  else ! This is the buoyancy only mode.
116  do j=js,je ; do i=is,ie
117  ! fluxes%buoy is the buoyancy flux into the ocean [L2 T-3 ~> m2 s-3]. A positive
118  ! buoyancy flux is of the same sign as heating the ocean.
119  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
120  enddo ; enddo
121  endif
122 
123  if (cs%restorebuoy) then
124  if (cs%use_temperature) then
125  call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
126  ! When modifying the code, comment out this error message. It is here
127  ! so that the original (unmodified) version is not accidentally used.
128  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
129  "Temperature and salinity restoring used without modification." )
130 
131  rhoxcp = cs%Rho0 * fluxes%C_p
132  do j=js,je ; do i=is,ie
133  ! Set Temp_restore and Salin_restore to the temperature (in degC) and
134  ! salinity (in ppt) that are being restored toward.
135  temp_restore = 0.0
136  salin_restore = 0.0
137 
138  fluxes%heat_added(i,j) = (g%mask2dT(i,j) * (rhoxcp * cs%Flux_const)) * &
139  (temp_restore - sfc_state%SST(i,j))
140  fluxes%vprec(i,j) = - (g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)) * &
141  ((salin_restore - sfc_state%SSS(i,j)) / (0.5 * (salin_restore + sfc_state%SSS(i,j))))
142  enddo ; enddo
143  else
144  ! When modifying the code, comment out this error message. It is here
145  ! so that the original (unmodified) version is not accidentally used.
146  ! call MOM_error(FATAL, "User_buoyancy_surface_forcing: " // &
147  ! "Buoyancy restoring used without modification." )
148 
149  ! The -1 is because density has the opposite sign to buoyancy.
150  buoy_rest_const = -1.0 * (cs%G_Earth * cs%Flux_const) / cs%Rho0
151  temp_restore = 0.0
152  do j=js,je ; do i=is,ie
153  ! Set density_restore to an expression for the surface potential
154  ! density [R ~> kg m-3] that is being restored toward.
155  if (g%geoLatT(i,j) < cs%lfrslat) then
156  temp_restore = cs%SST_s
157  elseif (g%geoLatT(i,j) > cs%lfrnlat) then
158  temp_restore = cs%SST_n
159  else
160  temp_restore = (cs%SST_s - cs%SST_n)/(cs%lfrslat - cs%lfrnlat) * &
161  (g%geoLatT(i,j) - cs%lfrslat) + cs%SST_s
162  endif
163 
164  density_restore = temp_restore*cs%drho_dt + cs%Rho0
165 
166  fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
167  (density_restore - sfc_state%sfc_density(i,j))
168  enddo ; enddo
169  endif
170  endif ! end RESTOREBUOY
171 

◆ bfb_surface_forcing_init()

subroutine, public bfb_surface_forcing::bfb_surface_forcing_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(in), target  diag,
type(bfb_surface_forcing_cs), pointer  CS 
)

Initialization for forcing the boundary-forced-basin (BFB) configuration.

Parameters
[in]timeThe current model time.
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
[in]param_fileA structure to parse for run-time parameters
[in]diagA structure that is used to regulate diagnostic output.
csA pointer to the control structure for this module

Definition at line 176 of file BFB_surface_forcing.F90.

176  type(time_type), intent(in) :: time !< The current model time.
177  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
178  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
179  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
180  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to
181  !! regulate diagnostic output.
182  type(bfb_surface_forcing_cs), pointer :: cs !< A pointer to the control structure for this module
183 ! This include declares and sets the variable "version".
184 #include "version_variable.h"
185  character(len=40) :: mdl = "BFB_surface_forcing" ! This module's name.
186 
187  if (associated(cs)) then
188  call mom_error(warning, "BFB_surface_forcing_init called with an associated "// &
189  "control structure.")
190  return
191  endif
192  allocate(cs)
193  cs%diag => diag
194 
195  ! Read all relevant parameters and write them to the model log.
196  call log_version(param_file, mdl, version, "")
197  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
198  "If true, Temperature and salinity are used as state variables.", default=.true.)
199 
200  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
201  "The gravitational acceleration of the Earth.", &
202  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
203  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
204  "The mean ocean density used with BOUSSINESQ true to "//&
205  "calculate accelerations and the mass for conservation "//&
206  "properties, or with BOUSSINSEQ false to convert some "//&
207  "parameters from vertical units of m to kg m-2.", &
208  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
209  call get_param(param_file, mdl, "LFR_SLAT", cs%lfrslat, &
210  "Southern latitude where the linear forcing ramp begins.", &
211  units="degrees", default=20.0)
212  call get_param(param_file, mdl, "LFR_NLAT", cs%lfrnlat, &
213  "Northern latitude where the linear forcing ramp ends.", &
214  units="degrees", default=40.0)
215  call get_param(param_file, mdl, "SST_S", cs%SST_s, &
216  "SST at the southern edge of the linear forcing ramp.", &
217  units="C", default=20.0)
218  call get_param(param_file, mdl, "SST_N", cs%SST_n, &
219  "SST at the northern edge of the linear forcing ramp.", &
220  units="C", default=10.0)
221  call get_param(param_file, mdl, "DRHO_DT", cs%drho_dt, &
222  "The rate of change of density with temperature.", &
223  units="kg m-3 K-1", default=-0.2, scale=us%kg_m3_to_R)
224  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
225  "The background gustiness in the winds.", units="Pa", default=0.0)
226 
227  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
228  "If true, the buoyancy fluxes drive the model back "//&
229  "toward some specified surface state with a rate "//&
230  "given by FLUXCONST.", default= .false.)
231  if (cs%restorebuoy) then
232  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
233  "The constant that relates the restoring surface fluxes to the relative "//&
234  "surface anomalies (akin to a piston velocity). Note the non-MKS units.", &
235  default=0.0, units="m day-1", scale=us%m_to_Z*us%T_to_s)
236  ! Convert CS%Flux_const from m day-1 to m s-1.
237  cs%Flux_const = cs%Flux_const / 86400.0
238  endif
239