MOM6
BFB_surface_forcing.F90
Go to the documentation of this file.
1 !> Surface forcing for the boundary-forced-basin (BFB) configuration
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
8 use mom_domains, only : pass_var, pass_vector, agrid
9 use mom_error_handler, only : mom_error, fatal, warning, is_root_pe
12 use mom_grid, only : ocean_grid_type
13 use mom_io, only : file_exists, read_data
15 use mom_time_manager, only : time_type, operator(+), operator(/)
19 use mom_variables, only : surface
20 
21 implicit none ; private
22 
24 
25 !> Control structure for BFB_surface_forcing
26 type, public :: bfb_surface_forcing_cs ; private
27 
28  logical :: use_temperature !< If true, temperature and salinity are used as state variables.
29  logical :: restorebuoy !< If true, use restoring surface buoyancy forcing.
30  real :: rho0 !< The density used in the Boussinesq approximation [R ~> kg m-3].
31  real :: g_earth !< The gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
32  real :: flux_const !< The restoring rate at the surface [Z T-1 ~> m s-1].
33  real :: gust_const !< A constant unresolved background gustiness
34  !! that contributes to ustar [Pa].
35  real :: sst_s !< SST at the southern edge of the linear forcing ramp [degC]
36  real :: sst_n !< SST at the northern edge of the linear forcing ramp [degC]
37  real :: lfrslat !< Southern latitude where the linear forcing ramp begins [degLat]
38  real :: lfrnlat !< Northern latitude where the linear forcing ramp ends [degLat]
39  real :: drho_dt !< Rate of change of density with temperature [R degC-1 ~> kg m-3 degC-1].
40  !! Note that temperature is being used as a dummy variable here.
41  !! All temperatures are converted into density.
42 
43  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
44  !! regulate the timing of diagnostic output.
46 
47 contains
48 
49 !> Bouyancy forcing for the boundary-forced-basin (BFB) configuration
50 subroutine bfb_buoyancy_forcing(state, fluxes, day, dt, G, US, CS)
51  type(surface), intent(inout) :: 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 [J T s-1 Z-1 m-2 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 [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 [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 = us%R_to_kg_m3*us%Z_to_m*us%s_to_T * 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 - state%SST(i,j))
140  fluxes%vprec(i,j) = - (g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)) * &
141  ((salin_restore - state%SSS(i,j)) / (0.5 * (salin_restore + 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 [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 - us%kg_m3_to_R*state%sfc_density(i,j))
168  enddo ; enddo
169  endif
170  endif ! end RESTOREBUOY
171 
172 end subroutine bfb_buoyancy_forcing
173 
174 !> Initialization for forcing the boundary-forced-basin (BFB) configuration
175 subroutine bfb_surface_forcing_init(Time, G, US, param_file, diag, CS)
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 "//&
199  "variables.", default=.true.)
200 
201  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
202  "The gravitational acceleration of the Earth.", &
203  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
204  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
205  "The mean ocean density used with BOUSSINESQ true to "//&
206  "calculate accelerations and the mass for conservation "//&
207  "properties, or with BOUSSINSEQ false to convert some "//&
208  "parameters from vertical units of m to kg m-2.", &
209  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
210  call get_param(param_file, mdl, "LFR_SLAT", cs%lfrslat, &
211  "Southern latitude where the linear forcing ramp begins.", &
212  units="degrees", default=20.0)
213  call get_param(param_file, mdl, "LFR_NLAT", cs%lfrnlat, &
214  "Northern latitude where the linear forcing ramp ends.", &
215  units="degrees", default=40.0)
216  call get_param(param_file, mdl, "SST_S", cs%SST_s, &
217  "SST at the southern edge of the linear forcing ramp.", &
218  units="C", default=20.0)
219  call get_param(param_file, mdl, "SST_N", cs%SST_n, &
220  "SST at the northern edge of the linear forcing ramp.", &
221  units="C", default=10.0)
222  call get_param(param_file, mdl, "DRHO_DT", cs%drho_dt, &
223  "The rate of change of density with temperature.", &
224  units="kg m-3 K-1", default=-0.2, scale=us%kg_m3_to_R)
225  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
226  "The background gustiness in the winds.", units="Pa", &
227  default=0.02)
228 
229  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
230  "If true, the buoyancy fluxes drive the model back "//&
231  "toward some specified surface state with a rate "//&
232  "given by FLUXCONST.", default= .false.)
233  if (cs%restorebuoy) then
234  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
235  "The constant that relates the restoring surface fluxes "//&
236  "to the relative surface anomalies (akin to a piston "//&
237  "velocity). Note the non-MKS units.", &
238  units="m day-1", scale=us%m_to_Z*us%T_to_s, fail_if_missing=.true.)
239  ! Convert CS%Flux_const from m day-1 to m s-1.
240  cs%Flux_const = cs%Flux_const / 86400.0
241  endif
242 
243 end subroutine bfb_surface_forcing_init
244 
245 end module bfb_surface_forcing
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_forcing_type::allocate_forcing_type
subroutine, public allocate_forcing_type(G, fluxes, water, heat, ustar, press, shelf, iceberg, salt)
Conditionally allocate fields within the forcing type.
Definition: MOM_forcing_type.F90:2811
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_safe_alloc
Convenience functions for safely allocating memory without accidentally reallocating pointer and caus...
Definition: MOM_safe_alloc.F90:3
mom_diag_mediator::query_averaging_enabled
logical function, public query_averaging_enabled(diag_cs, time_int, time_end)
Call this subroutine to determine whether the averaging is currently enabled. .true....
Definition: MOM_diag_mediator.F90:1850
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_diag_mediator::register_diag_field
integer function, public register_diag_field(module_name, field_name, axes_in, 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...
Definition: MOM_diag_mediator.F90:1878
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
bfb_surface_forcing::bfb_surface_forcing_init
subroutine, public bfb_surface_forcing_init(Time, G, US, param_file, diag, CS)
Initialization for forcing the boundary-forced-basin (BFB) configuration.
Definition: BFB_surface_forcing.F90:176
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
bfb_surface_forcing::bfb_buoyancy_forcing
subroutine, public bfb_buoyancy_forcing(state, fluxes, day, dt, G, US, CS)
Bouyancy forcing for the boundary-forced-basin (BFB) configuration.
Definition: BFB_surface_forcing.F90:51
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_tracer_flow_control
Orchestrates the registration and calling of tracer packages.
Definition: MOM_tracer_flow_control.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_tracer_flow_control::tracer_flow_control_cs
The control structure for orchestrating the calling of tracer packages.
Definition: MOM_tracer_flow_control.F90:75
mom_error_handler::is_root_pe
logical function, public is_root_pe()
This returns .true. if the current PE is the root PE.
Definition: MOM_error_handler.F90:44
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:50
bfb_surface_forcing
Surface forcing for the boundary-forced-basin (BFB) configuration.
Definition: BFB_surface_forcing.F90:2
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_error_handler::mom_error
subroutine, public mom_error(level, message, all_print)
This provides a convenient interface for writing an mpp_error message with run-time filter based on a...
Definition: MOM_error_handler.F90:72
mom_safe_alloc::safe_alloc_ptr
Allocate a pointer to a 1-d, 2-d or 3-d array.
Definition: MOM_safe_alloc.F90:12
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
mom_tracer_flow_control::call_tracer_set_forcing
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...
Definition: MOM_tracer_flow_control.F90:382
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239
bfb_surface_forcing::bfb_surface_forcing_cs
Control structure for BFB_surface_forcing.
Definition: BFB_surface_forcing.F90:26