MOM6
BFB_surface_forcing.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 !* Rewritten by Robert Hallberg, June 2009 *
25 !* *
26 !* This file contains subroutines for specifying surface buoyancy *
27 !* forcing for the buoyancy-forced basin (BFB) case. *
28 !* BFB_buoyancy_forcing is used to restore the surface buoayncy to *
29 !* a linear meridional ramp of temperature. The extent of the ramp *
30 !* can be specified by LFR_SLAT (linear forcing ramp southern *
31 !* latitude) and LFR_NLAT. The temperatures at these edges of the *
32 !* ramp can be specified by SST_S and SST_N. *
33 !* *
34 !********+*********+*********+*********+*********+*********+*********+**
37 use mom_domains, only : pass_var, pass_vector, agrid
38 use mom_error_handler, only : mom_error, fatal, warning, is_root_pe
41 use mom_grid, only : ocean_grid_type
42 use mom_io, only : file_exists, read_data
43 use mom_time_manager, only : time_type, operator(+), operator(/), get_time
46 use mom_variables, only : surface
47 
48 implicit none ; private
49 
51 
52 type, public :: bfb_surface_forcing_cs ; private
53  ! This control structure should be used to store any run-time variables
54  ! associated with the user-specified forcing. It can be readily modified
55  ! for a specific case, and because it is private there will be no changes
56  ! needed in other code (although they will have to be recompiled).
57  ! The variables in the cannonical example are used for some common
58  ! cases, but do not need to be used.
59 
60  logical :: use_temperature ! If true, temperature and salinity are used as
61  ! state variables.
62  logical :: restorebuoy ! If true, use restoring surface buoyancy forcing.
63  real :: rho0 ! The density used in the Boussinesq
64  ! approximation, in kg m-3.
65  real :: g_earth ! The gravitational acceleration in m s-2.
66  real :: flux_const ! The restoring rate at the surface, in m s-1.
67  real :: gust_const ! A constant unresolved background gustiness
68  ! that contributes to ustar, in Pa.
69  real :: sst_s ! SST at the southern edge of the linear
70  ! forcing ramp
71  real :: sst_n ! SST at the northern edge of the linear
72  ! forcing ramp
73  real :: lfrslat ! Southern latitude where the linear forcing ramp
74  ! begins
75  real :: lfrnlat ! Northern latitude where the linear forcing ramp
76  ! ends
77  real :: drho_dt ! Rate of change of density with temperature.
78  ! Note that temperature is being used as a dummy
79  ! variable here. All temperatures are converted
80  ! into density.
81 
82  type(diag_ctrl), pointer :: diag ! A structure that is used to regulate the
83  ! timing of diagnostic output.
85 
86 contains
87 
88 subroutine bfb_buoyancy_forcing(state, fluxes, day, dt, G, CS)
89  type(surface), intent(inout) :: state
90  type(forcing), intent(inout) :: fluxes
91  type(time_type), intent(in) :: day
92  real, intent(in) :: dt !< The amount of time over which
93  !! the fluxes apply, in s
94  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
95  type(bfb_surface_forcing_cs), pointer :: CS
96 
97 ! This subroutine specifies the current surface fluxes of buoyancy or
98 ! temperature and fresh water. It may also be modified to add
99 ! surface fluxes of user provided tracers.
100 
101 ! When temperature is used, there are long list of fluxes that need to be
102 ! set - essentially the same as for a full coupled model, but most of these
103 ! can be simply set to zero. The net fresh water flux should probably be
104 ! set in fluxes%evap and fluxes%lprec, with any salinity restoring
105 ! appearing in fluxes%vprec, and the other water flux components
106 ! (fprec, lrunoff and frunoff) left as arrays full of zeros.
107 ! Evap is usually negative and precip is usually positive. All heat fluxes
108 ! are in W m-2 and positive for heat going into the ocean. All fresh water
109 ! fluxes are in kg m-2 s-1 and positive for water moving into the ocean.
110 
111 ! Arguments: state - A structure containing fields that describe the
112 ! surface state of the ocean.
113 ! (out) fluxes - A structure containing pointers to any possible
114 ! forcing fields. Unused fields have NULL ptrs.
115 ! (in) day_start - Start time of the fluxes.
116 ! (in) day_interval - Length of time over which these fluxes
117 ! will be applied.
118 ! (in) G - The ocean's grid structure.
119 ! (in) CS - A pointer to the control structure returned by a previous
120 ! call to user_surface_forcing_init
121 
122  real :: Temp_restore ! The temperature that is being restored toward, in C.
123  real :: Salin_restore ! The salinity that is being restored toward, in PSU.
124  real :: density_restore ! The potential density that is being restored
125  ! toward, in kg m-3.
126  real :: rhoXcp ! The mean density times the heat capacity, in J m-3 K-1.
127  real :: buoy_rest_const ! A constant relating density anomalies to the
128  ! restoring buoyancy flux, in m5 s-3 kg-1.
129  integer :: i, j, is, ie, js, je
130  integer :: isd, ied, jsd, jed
131 
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
134 
135  ! When modifying the code, comment out this error message. It is here
136  ! so that the original (unmodified) version is not accidentally used.
137  ! call MOM_error(FATAL, "User_buoyancy_surface_forcing: " // &
138  ! "User forcing routine called without modification." )
139 
140  ! Allocate and zero out the forcing arrays, as necessary. This portion is
141  ! usually not changed.
142  if (cs%use_temperature) then
143  call alloc_if_needed(fluxes%evap, isd, ied, jsd, jed)
144  call alloc_if_needed(fluxes%lprec, isd, ied, jsd, jed)
145  call alloc_if_needed(fluxes%fprec, isd, ied, jsd, jed)
146  call alloc_if_needed(fluxes%lrunoff, isd, ied, jsd, jed)
147  call alloc_if_needed(fluxes%frunoff, isd, ied, jsd, jed)
148  call alloc_if_needed(fluxes%vprec, isd, ied, jsd, jed)
149 
150  call alloc_if_needed(fluxes%sw, isd, ied, jsd, jed)
151  call alloc_if_needed(fluxes%lw, isd, ied, jsd, jed)
152  call alloc_if_needed(fluxes%latent, isd, ied, jsd, jed)
153  call alloc_if_needed(fluxes%sens, isd, ied, jsd, jed)
154  else ! This is the buoyancy only mode.
155  call alloc_if_needed(fluxes%buoy, isd, ied, jsd, jed)
156  endif
157 
158 
159  ! MODIFY THE CODE IN THE FOLLOWING LOOPS TO SET THE BUOYANCY FORCING TERMS.
160 
161  if ( cs%use_temperature ) then
162  ! Set whichever fluxes are to be used here. Any fluxes that
163  ! are always zero do not need to be changed here.
164  do j=js,je ; do i=is,ie
165  ! Fluxes of fresh water through the surface are in units of kg m-2 s-1
166  ! and are positive downward - i.e. evaporation should be negative.
167  fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
168  fluxes%lprec(i,j) = 0.0 * g%mask2dT(i,j)
169 
170  ! vprec will be set later, if it is needed for salinity restoring.
171  fluxes%vprec(i,j) = 0.0
172 
173  ! Heat fluxes are in units of W m-2 and are positive into the ocean.
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)
178  enddo ; enddo
179  else ! This is the buoyancy only mode.
180  do j=js,je ; do i=is,ie
181  ! fluxes%buoy is the buoyancy flux into the ocean in m2 s-3. A positive
182  ! buoyancy flux is of the same sign as heating the ocean.
183  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
184  enddo ; enddo
185  endif
186 
187  if (cs%restorebuoy) then
188  if (cs%use_temperature) then
189  call alloc_if_needed(fluxes%heat_added, isd, ied, jsd, jed)
190  ! When modifying the code, comment out this error message. It is here
191  ! so that the original (unmodified) version is not accidentally used.
192  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
193  "Temperature and salinity restoring used without modification." )
194 
195  rhoxcp = cs%Rho0 * fluxes%C_p
196  do j=js,je ; do i=is,ie
197  ! Set Temp_restore and Salin_restore to the temperature (in C) and
198  ! salinity (in PSU) that are being restored toward.
199  temp_restore = 0.0
200  salin_restore = 0.0
201 
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))))
207  enddo ; enddo
208  else
209  ! When modifying the code, comment out this error message. It is here
210  ! so that the original (unmodified) version is not accidentally used.
211  ! call MOM_error(FATAL, "User_buoyancy_surface_forcing: " // &
212  ! "Buoyancy restoring used without modification." )
213 
214  ! The -1 is because density has the opposite sign to buoyancy.
215  buoy_rest_const = -1.0 * (cs%G_Earth * cs%Flux_const) / cs%Rho0
216  temp_restore = 0.0
217  do j=js,je ; do i=is,ie
218  ! Set density_restore to an expression for the surface potential
219  ! density in kg m-3 that is being restored toward.
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
224  else
225  temp_restore = (cs%SST_s - cs%SST_n)/(cs%lfrslat - cs%lfrnlat) * &
226  (g%geoLatT(i,j) - cs%lfrslat) + cs%SST_s
227  end if
228 
229  density_restore = temp_restore*cs%drho_dt + cs%Rho0
230 
231  fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
232  (density_restore - state%sfc_density(i,j))
233  enddo ; enddo
234  endif
235  endif ! end RESTOREBUOY
236 
237 end subroutine bfb_buoyancy_forcing
238 
239 subroutine alloc_if_needed(ptr, isd, ied, jsd, jed)
240  ! If ptr is not associated, this routine allocates it with the given size
241  ! and zeros out its contents. This is equivalent to safe_alloc_ptr in
242  ! MOM_diag_mediator, but is here so as to be completely transparent.
243  real, pointer :: ptr(:,:)
244  integer :: isd, ied, jsd, jed
245  if (.not.ASSOCIATED(ptr)) then
246  allocate(ptr(isd:ied,jsd:jed))
247  ptr(:,:) = 0.0
248  endif
249 end subroutine alloc_if_needed
250 
251 subroutine bfb_surface_forcing_init(Time, G, param_file, diag, CS)
252  type(time_type), intent(in) :: Time
253  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
254  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
255  type(diag_ctrl), target, intent(in) :: diag
256  type(bfb_surface_forcing_cs), pointer :: CS
257 ! Arguments: Time - The current model time.
258 ! (in) G - The ocean's grid structure.
259 ! (in) param_file - A structure indicating the open file to parse for
260 ! model parameter values.
261 ! (in) diag - A structure that is used to regulate diagnostic output.
262 ! (in/out) CS - A pointer that is set to point to the control structure
263 ! for this module
264 
265 ! This include declares and sets the variable "version".
266 #include "version_variable.h"
267  character(len=40) :: mdl = "BFB_surface_forcing" ! This module's name.
268 
269  if (associated(cs)) then
270  call mom_error(warning, "BFB_surface_forcing_init called with an associated "// &
271  "control structure.")
272  return
273  endif
274  allocate(cs)
275  cs%diag => diag
276 
277  ! Read all relevant parameters and write them to the model log.
278  call log_version(param_file, mdl, version, "")
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.)
282 
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", &
309  default=0.02)
310 
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.)
321  ! Convert CS%Flux_const from m day-1 to m s-1.
322  cs%Flux_const = cs%Flux_const / 86400.0
323  endif
324 
325 end subroutine bfb_surface_forcing_init
326 
327 end module bfb_surface_forcing
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.
Definition: MOM_grid.F90:19
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Provides the ocean grid type.
Definition: MOM_grid.F90:2
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.
Definition: MOM_io.F90:2
subroutine, public bfb_surface_forcing_init(Time, G, param_file, diag, CS)
logical function, public query_averaging_enabled(diag_cs, time_int, time_end)
subroutine, public call_tracer_set_forcing(state, fluxes, day_start, day_interval, G, CS)
This subroutine calls the individual tracer modules&#39; 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 ...
integer function, public register_diag_field(module_name, field_name, axes, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name, cell_methods, x_cell_method, y_cell_method, v_cell_method, conversion, v_extensive)
Returns the "diag_mediator" handle for a group (native, CMOR, z-coord, ...) of diagnostics derived fr...
subroutine, public mom_error(level, message, all_print)