48 implicit none ;
private 60 logical :: use_temperature
62 logical :: restorebuoy
90 type(
forcing),
intent(inout) :: fluxes
91 type(time_type),
intent(in) :: day
92 real,
intent(in) :: dt
123 real :: Salin_restore
124 real :: density_restore
127 real :: buoy_rest_const
129 integer :: i, j, is, ie, js, je
130 integer :: isd, ied, jsd, jed
132 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
133 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
142 if (cs%use_temperature)
then 161 if ( cs%use_temperature )
then 164 do j=js,je ;
do i=is,ie
167 fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
168 fluxes%lprec(i,j) = 0.0 * g%mask2dT(i,j)
171 fluxes%vprec(i,j) = 0.0
174 fluxes%lw(i,j) = 0.0 * g%mask2dT(i,j)
175 fluxes%latent(i,j) = 0.0 * g%mask2dT(i,j)
176 fluxes%sens(i,j) = 0.0 * g%mask2dT(i,j)
177 fluxes%sw(i,j) = 0.0 * g%mask2dT(i,j)
180 do j=js,je ;
do i=is,ie
183 fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
187 if (cs%restorebuoy)
then 188 if (cs%use_temperature)
then 192 call mom_error(fatal,
"User_buoyancy_surface_forcing: " // &
193 "Temperature and salinity restoring used without modification." )
195 rhoxcp = cs%Rho0 * fluxes%C_p
196 do j=js,je ;
do i=is,ie
202 fluxes%heat_added(i,j) = (g%mask2dT(i,j) * (rhoxcp * cs%Flux_const)) * &
203 (temp_restore - state%SST(i,j))
204 fluxes%vprec(i,j) = - (g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)) * &
205 ((salin_restore - state%SSS(i,j)) / &
206 (0.5 * (salin_restore + state%SSS(i,j))))
215 buoy_rest_const = -1.0 * (cs%G_Earth * cs%Flux_const) / cs%Rho0
217 do j=js,je ;
do i=is,ie
220 if (g%geoLatT(i,j) < cs%lfrslat)
then 221 temp_restore = cs%SST_s
222 else if (g%geoLatT(i,j) > cs%lfrnlat)
then 223 temp_restore = cs%SST_n
225 temp_restore = (cs%SST_s - cs%SST_n)/(cs%lfrslat - cs%lfrnlat) * &
226 (g%geoLatT(i,j) - cs%lfrslat) + cs%SST_s
229 density_restore = temp_restore*cs%drho_dt + cs%Rho0
231 fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
232 (density_restore - state%sfc_density(i,j))
243 real,
pointer :: ptr(:,:)
244 integer :: isd, ied, jsd, jed
245 if (.not.
ASSOCIATED(ptr))
then 246 allocate(ptr(isd:ied,jsd:jed))
252 type(time_type),
intent(in) :: Time
255 type(
diag_ctrl),
target,
intent(in) :: diag
266 #include "version_variable.h" 267 character(len=40) :: mdl =
"BFB_surface_forcing" 269 if (
associated(cs))
then 270 call mom_error(warning,
"BFB_surface_forcing_init called with an associated "// &
271 "control structure.")
279 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", cs%use_temperature, &
280 "If true, Temperature and salinity are used as state \n"//&
281 "variables.", default=.true.)
283 call get_param(param_file, mdl,
"G_EARTH", cs%G_Earth, &
284 "The gravitational acceleration of the Earth.", &
285 units=
"m s-2", default = 9.80)
286 call get_param(param_file, mdl,
"RHO_0", cs%Rho0, &
287 "The mean ocean density used with BOUSSINESQ true to \n"//&
288 "calculate accelerations and the mass for conservation \n"//&
289 "properties, or with BOUSSINSEQ false to convert some \n"//&
290 "parameters from vertical units of m to kg m-2.", &
291 units=
"kg m-3", default=1035.0)
292 call get_param(param_file, mdl,
"LFR_SLAT", cs%lfrslat, &
293 "Southern latitude where the linear forcing ramp begins.", &
294 units=
"degrees", default = 20.0)
295 call get_param(param_file, mdl,
"LFR_NLAT", cs%lfrnlat, &
296 "Northern latitude where the linear forcing ramp ends.", &
297 units=
"degrees", default = 40.0)
298 call get_param(param_file, mdl,
"SST_S", cs%SST_s, &
299 "SST at the southern edge of the linear forcing ramp.", &
300 units=
"C", default = 20.0)
301 call get_param(param_file, mdl,
"SST_N", cs%SST_n, &
302 "SST at the northern edge of the linear forcing ramp.", &
303 units=
"C", default = 10.0)
304 call get_param(param_file, mdl,
"DRHO_DT", cs%drho_dt, &
305 "The rate of change of density with temperature.", &
306 units=
"kg m-3 K-1", default = -0.2)
307 call get_param(param_file, mdl,
"GUST_CONST", cs%gust_const, &
308 "The background gustiness in the winds.", units=
"Pa", &
311 call get_param(param_file, mdl,
"RESTOREBUOY", cs%restorebuoy, &
312 "If true, the buoyancy fluxes drive the model back \n"//&
313 "toward some specified surface state with a rate \n"//&
314 "given by FLUXCONST.", default= .false.)
315 if (cs%restorebuoy)
then 316 call get_param(param_file, mdl,
"FLUXCONST", cs%Flux_const, &
317 "The constant that relates the restoring surface fluxes \n"//&
318 "to the relative surface anomalies (akin to a piston \n"//&
319 "velocity). Note the non-MKS units.", units=
"m day-1", &
320 fail_if_missing=.true.)
322 cs%Flux_const = cs%Flux_const / 86400.0
The following structure contains pointers to various fields which may be used describe the surface st...
subroutine alloc_if_needed(ptr, isd, ied, jsd, jed)
subroutine, public bfb_buoyancy_forcing(state, fluxes, day, dt, G, CS)
This module implements boundary forcing for MOM6.
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
subroutine, public allocate_forcing_type(G, fluxes, stress, ustar, water, heat, shelf, press, iceberg)
Conditionally allocate fields within the forcing type.
This module contains I/O framework code.
subroutine, public bfb_surface_forcing_init(Time, G, param_file, diag, CS)
subroutine, public call_tracer_set_forcing(state, fluxes, day_start, day_interval, G, CS)
This subroutine calls the individual tracer modules' subroutines to specify or read quantities relate...
logical function, public is_root_pe()
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
subroutine, public mom_error(level, message, all_print)