MOM6
MOM_surface_forcing.F90
1 !> Functions that calculate the surface wind stresses and fluxes of buoyancy
2 !! or temperature/salinity andfresh water, in ocean-only (solo) mode.
3 !!
4 !! These functions are called every time step, even if the wind stresses
5 !! or buoyancy fluxes are constant in time - in that case these routines
6 !! return quickly without doing anything. In addition, any I/O of forcing
7 !! fields is controlled by surface_forcing_init, located in this file.
9 
10 ! This file is part of MOM6. See LICENSE.md for the license.
11 
12 use mom_constants, only : hlv, hlf
13 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
14 use mom_cpu_clock, only : clock_module
15 use mom_diag_mediator, only : post_data, query_averaging_enabled
16 use mom_diag_mediator, only : diag_ctrl, safe_alloc_ptr
17 use mom_domains, only : pass_var, pass_vector, agrid, to_south, to_west, to_all
18 use mom_domains, only : fill_symmetric_edges, cgrid_ne
19 use mom_error_handler, only : mom_error, fatal, warning, mom_mesg, is_root_pe
20 use mom_error_handler, only : calltree_enter, calltree_leave
22 use mom_string_functions, only : uppercase
24 use mom_forcing_type, only : set_net_mass_forcing, copy_common_forcing_fields
25 use mom_forcing_type, only : set_derived_forcing_fields
26 use mom_forcing_type, only : forcing_diags, mech_forcing_diags, register_forcing_type_diags
27 use mom_forcing_type, only : allocate_forcing_type, deallocate_forcing_type
28 use mom_forcing_type, only : allocate_mech_forcing, deallocate_mech_forcing
29 use mom_grid, only : ocean_grid_type
30 use mom_get_input, only : get_mom_input, directories
31 use mom_io, only : file_exists, mom_read_data, mom_read_vector, slasher
32 use mom_io, only : east_face, north_face, num_timelevels
33 use mom_restart, only : register_restart_field, restart_init, mom_restart_cs
34 use mom_restart, only : restart_init_end, save_restart, restore_state
35 use mom_time_manager, only : time_type, operator(+), operator(/), get_time, time_type_to_real
36 use mom_tracer_flow_control, only : call_tracer_set_forcing
39 use mom_variables, only : surface
40 use meso_surface_forcing, only : meso_buoyancy_forcing
41 use meso_surface_forcing, only : meso_surface_forcing_init, meso_surface_forcing_cs
42 use user_surface_forcing, only : user_wind_forcing, user_buoyancy_forcing
43 use user_surface_forcing, only : user_surface_forcing_init, user_surface_forcing_cs
44 use user_revise_forcing, only : user_alter_forcing, user_revise_forcing_init
46 use idealized_hurricane, only : idealized_hurricane_wind_init
47 use idealized_hurricane, only : idealized_hurricane_wind_forcing, scm_idealized_hurricane_wind_forcing
49 use scm_cvmix_tests, only : scm_cvmix_tests_surface_forcing_init
50 use scm_cvmix_tests, only : scm_cvmix_tests_wind_forcing
51 use scm_cvmix_tests, only : scm_cvmix_tests_buoyancy_forcing
53 use bfb_surface_forcing, only : bfb_buoyancy_forcing
54 use bfb_surface_forcing, only : bfb_surface_forcing_init, bfb_surface_forcing_cs
55 use dumbbell_surface_forcing, only : dumbbell_surface_forcing_init, dumbbell_surface_forcing_cs
56 use dumbbell_surface_forcing, only : dumbbell_buoyancy_forcing
57 use data_override_mod, only : data_override_init, data_override
58 
59 implicit none ; private
60 
61 #include <MOM_memory.h>
62 
63 public set_forcing
64 public surface_forcing_init
65 public forcing_save_restart
66 
67 !> Structure containing pointers to the forcing fields that may be used to drive MOM.
68 !! All fluxes are positive into the ocean.
69 type, public :: surface_forcing_cs ; private
70 
71  logical :: use_temperature !< if true, temp & salinity used as state variables
72  logical :: restorebuoy !< if true, use restoring surface buoyancy forcing
73  logical :: adiabatic !< if true, no diapycnal mass fluxes or surface buoyancy forcing
74  logical :: variable_winds !< if true, wind stresses vary with time
75  logical :: variable_buoyforce !< if true, buoyancy forcing varies with time.
76  real :: south_lat !< southern latitude of the domain
77  real :: len_lat !< domain length in latitude
78 
79  real :: rho0 !< Boussinesq reference density [R ~> kg m-3]
80  real :: g_earth !< gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
81  real :: flux_const !< piston velocity for surface restoring [Z T-1 ~> m s-1]
82  real :: flux_const_t !< piston velocity for surface temperature restoring [m s-1]
83  real :: flux_const_s !< piston velocity for surface salinity restoring [Z T-1 ~> m s-1]
84  real :: latent_heat_fusion !< latent heat of fusion times [Q ~> J kg-1]
85  real :: latent_heat_vapor !< latent heat of vaporization [Q ~> J kg-1]
86  real :: tau_x0 !< Constant zonal wind stress used in the WIND_CONFIG="const" forcing
87  real :: tau_y0 !< Constant meridional wind stress used in the WIND_CONFIG="const" forcing
88 
89  real :: gust_const !< constant unresolved background gustiness for ustar [R L Z T-1 ~> Pa]
90  logical :: read_gust_2d !< if true, use 2-dimensional gustiness supplied from a file
91  real, pointer :: gust(:,:) => null() !< spatially varying unresolved background gustiness [R L Z T-1 ~> Pa]
92  !! gust is used when read_gust_2d is true.
93 
94  real, pointer :: t_restore(:,:) => null() !< temperature to damp (restore) the SST to [degC]
95  real, pointer :: s_restore(:,:) => null() !< salinity to damp (restore) the SSS [ppt]
96  real, pointer :: dens_restore(:,:) => null() !< density to damp (restore) surface density [R ~> kg m-3]
97 
98  integer :: buoy_last_lev_read = -1 !< The last time level read from buoyancy input files
99 
100  ! if WIND_CONFIG=='gyres' then use the following as = A, B, C and n respectively for
101  ! taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L)
102  real :: gyres_taux_const !< A constant wind stress [Pa].
103  real :: gyres_taux_sin_amp !< The amplitude of cosine wind stress gyres [Pa], if WIND_CONFIG=='gyres'.
104  real :: gyres_taux_cos_amp !< The amplitude of cosine wind stress gyres [Pa], if WIND_CONFIG=='gyres'.
105  real :: gyres_taux_n_pis !< The number of sine lobes in the basin if if WIND_CONFIG=='gyres'
106  logical :: answers_2018 !< If true, use the order of arithmetic and expressions that recover
107  !! the answers from the end of 2018. Otherwise, use a form of the gyre
108  !! wind stresses that are rotationally invariant and more likely to be
109  !! the same between compilers.
110  logical :: fix_ustar_gustless_bug !< If true correct a bug in the time-averaging of the
111  !! gustless wind friction velocity.
112  ! if WIND_CONFIG=='scurves' then use the following to define a piecwise scurve profile
113  real :: scurves_ydata(20) = 90. !< Latitudes of scurve nodes [degreesN]
114  real :: scurves_taux(20) = 0. !< Zonal wind stress values at scurve nodes [Pa]
115 
116  real :: t_north !< target temperatures at north used in buoyancy_forcing_linear
117  real :: t_south !< target temperatures at south used in buoyancy_forcing_linear
118  real :: s_north !< target salinity at north used in buoyancy_forcing_linear
119  real :: s_south !< target salinity at south used in buoyancy_forcing_linear
120 
121  logical :: first_call_set_forcing = .true. !< True until after the first call to set_forcing
122  logical :: archaic_omip_file = .true. !< If true use the variable names and data fields from
123  !! a very old version of the OMIP forcing
124  logical :: dataoverrideisinitialized = .false. !< If true, data override has been initialized
125 
126  real :: wind_scale !< value by which wind-stresses are scaled, ND.
127  real :: constantheatforcing !< value used for sensible heat flux when buoy_config="const" [Q R Z T-1 ~> W m-2]
128 
129  character(len=8) :: wind_stagger !< A character indicating how the wind stress components
130  !! are staggered in WIND_FILE. Valid values are A or C for now.
131  type(tracer_flow_control_cs), pointer :: tracer_flow_csp => null() !< A pointer to the structure
132  !! that is used to orchestrate the calling of tracer packages
133 !#CTRL# type(ctrl_forcing_CS), pointer :: ctrl_forcing_CSp => NULL()
134  type(mom_restart_cs), pointer :: restart_csp => null() !< A pointer to the restart control structure
135 
136  type(diag_ctrl), pointer :: diag !< structure used to regulate timing of diagnostic output
137 
138  character(len=200) :: inputdir !< directory where NetCDF input files are.
139  character(len=200) :: wind_config !< indicator for wind forcing type (2gyre, USER, FILE..)
140  character(len=200) :: wind_file !< if wind_config is "file", file to use
141  character(len=200) :: buoy_config !< indicator for buoyancy forcing type
142 
143  character(len=200) :: longwave_file = '' !< The file from which the longwave heat flux is read
144  character(len=200) :: shortwave_file = '' !< The file from which the shortwave heat flux is read
145  character(len=200) :: evaporation_file = '' !< The file from which the evaporation is read
146  character(len=200) :: sensibleheat_file = '' !< The file from which the sensible heat flux is read
147  character(len=200) :: latentheat_file = '' !< The file from which the latent heat flux is read
148 
149  character(len=200) :: rain_file = '' !< The file from which the rainfall is read
150  character(len=200) :: snow_file = '' !< The file from which the snowfall is read
151  character(len=200) :: runoff_file = '' !< The file from which the runoff is read
152 
153  character(len=200) :: longwaveup_file = '' !< The file from which the upward longwave heat flux is read
154  character(len=200) :: shortwaveup_file = '' !< The file from which the upward shorwave heat flux is read
155 
156  character(len=200) :: sstrestore_file = '' !< The file from which to read the sea surface
157  !! temperature to restore toward
158  character(len=200) :: salinityrestore_file = '' !< The file from which to read the sea surface
159  !! salinity to restore toward
160 
161  character(len=80) :: stress_x_var = '' !< X-windstress variable name in the input file
162  character(len=80) :: stress_y_var = '' !< Y-windstress variable name in the input file
163  character(len=80) :: ustar_var = '' !< ustar variable name in the input file
164  character(len=80) :: lw_var = '' !< lonngwave heat flux variable name in the input file
165  character(len=80) :: sw_var = '' !< shortwave heat flux variable name in the input file
166  character(len=80) :: latent_var = '' !< latent heat flux variable name in the input file
167  character(len=80) :: sens_var = '' !< sensible heat flux variable name in the input file
168  character(len=80) :: evap_var = '' !< evaporation variable name in the input file
169  character(len=80) :: rain_var = '' !< rainfall variable name in the input file
170  character(len=80) :: snow_var = '' !< snowfall variable name in the input file
171  character(len=80) :: lrunoff_var = '' !< liquid runoff variable name in the input file
172  character(len=80) :: frunoff_var = '' !< frozen runoff variable name in the input file
173  character(len=80) :: sst_restore_var = '' !< target sea surface temeperature variable name in the input file
174  character(len=80) :: sss_restore_var = '' !< target sea surface salinity variable name in the input file
175 
176  ! These variables give the number of time levels in the various forcing files.
177  integer :: wind_nlev = -1 !< The number of time levels in the file of wind stress
178  integer :: sw_nlev = -1 !< The number of time levels in the file of shortwave heat flux
179  integer :: lw_nlev = -1 !< The number of time levels in the file of longwave heat flux
180  integer :: latent_nlev = -1 !< The number of time levels in the file of latent heat flux
181  integer :: sens_nlev = -1 !< The number of time levels in the file of sensible heat flux
182  integer :: evap_nlev = -1 !< The number of time levels in the file of evaporation
183  integer :: precip_nlev = -1 !< The number of time levels in the file of precipitation
184  integer :: runoff_nlev = -1 !< The number of time levels in the file of runoff
185  integer :: sst_nlev = -1 !< The number of time levels in the file of target SST
186  integer :: sss_nlev = -1 !< The number of time levels in the file of target SSS
187 
188  ! These variables give the last time level read for the various forcing files.
189  integer :: wind_last_lev = -1 !< The last time level read of wind stress
190  integer :: sw_last_lev = -1 !< The last time level read of shortwave heat flux
191  integer :: lw_last_lev = -1 !< The last time level read of longwave heat flux
192  integer :: latent_last_lev = -1 !< The last time level read of latent heat flux
193  integer :: sens_last_lev = -1 !< The last time level read of sensible heat flux
194  integer :: evap_last_lev = -1 !< The last time level read of evaporation
195  integer :: precip_last_lev = -1 !< The last time level read of precipitation
196  integer :: runoff_last_lev = -1 !< The last time level read of runoff
197  integer :: sst_last_lev = -1 !< The last time level read of target SST
198  integer :: sss_last_lev = -1 !< The last time level read of target SSS
199 
200  type(forcing_diags), public :: handles !< A structure with diagnostics handles
201 
202  !>@{ Control structures for named forcing packages
203  type(user_revise_forcing_cs), pointer :: urf_cs => null()
204  type(user_surface_forcing_cs), pointer :: user_forcing_csp => null()
205  type(bfb_surface_forcing_cs), pointer :: bfb_forcing_csp => null()
206  type(dumbbell_surface_forcing_cs), pointer :: dumbbell_forcing_csp => null()
207  type(meso_surface_forcing_cs), pointer :: meso_forcing_csp => null()
208  type(idealized_hurricane_cs), pointer :: idealized_hurricane_csp => null()
209  type(scm_cvmix_tests_cs), pointer :: scm_cvmix_tests_csp => null()
210  !>@}
211 
212 end type surface_forcing_cs
213 
214 integer :: id_clock_forcing !< A CPU time clock
215 
216 contains
217 
218 !> Calls subroutines in this file to get surface forcing fields.
219 !!
220 !! It also allocates and initializes the fields in the forcing and mech_forcing types
221 !! the first time it is called.
222 subroutine set_forcing(sfc_state, forces, fluxes, day_start, day_interval, G, US, CS)
223  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
224  !! describe the surface state of the ocean.
225  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
226  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
227  type(time_type), intent(in) :: day_start !< The start time of the fluxes
228  type(time_type), intent(in) :: day_interval !< Length of time over which these fluxes applied
229  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
230  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
231  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
232  !! a previous surface_forcing_init call
233  ! Local variables
234  real :: dt ! length of time over which fluxes applied [s]
235  type(time_type) :: day_center ! central time of the fluxes.
236  integer :: isd, ied, jsd, jed
237  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
238 
239  call cpu_clock_begin(id_clock_forcing)
240  call calltree_enter("set_forcing, MOM_surface_forcing.F90")
241 
242  day_center = day_start + day_interval/2
243  dt = time_type_to_real(day_interval)
244 
245  if (cs%first_call_set_forcing) then
246  ! Allocate memory for the mechanical and thermodyanmic forcing fields.
247  call allocate_mech_forcing(g, forces, stress=.true., ustar=.true., press=.true.)
248 
249  call allocate_forcing_type(g, fluxes, ustar=.true., fix_accum_bug=cs%fix_ustar_gustless_bug)
250  if (trim(cs%buoy_config) /= "NONE") then
251  if ( cs%use_temperature ) then
252  call allocate_forcing_type(g, fluxes, water=.true., heat=.true., press=.true.)
253  if (cs%restorebuoy) then
254  call safe_alloc_ptr(cs%T_Restore,isd, ied, jsd, jed)
255  call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
256  call safe_alloc_ptr(cs%S_Restore, isd, ied, jsd, jed)
257  endif
258  else ! CS%use_temperature false.
259  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
260 
261  if (cs%restorebuoy) call safe_alloc_ptr(cs%Dens_Restore, isd, ied, jsd, jed)
262  endif ! endif for CS%use_temperature
263  endif
264  endif
265 
266  ! calls to various wind options
267  if (cs%variable_winds .or. cs%first_call_set_forcing) then
268 
269  if (trim(cs%wind_config) == "file") then
270  call wind_forcing_from_file(sfc_state, forces, day_center, g, us, cs)
271  elseif (trim(cs%wind_config) == "data_override") then
272  call wind_forcing_by_data_override(sfc_state, forces, day_center, g, us, cs)
273  elseif (trim(cs%wind_config) == "2gyre") then
274  call wind_forcing_2gyre(sfc_state, forces, day_center, g, us, cs)
275  elseif (trim(cs%wind_config) == "1gyre") then
276  call wind_forcing_1gyre(sfc_state, forces, day_center, g, us, cs)
277  elseif (trim(cs%wind_config) == "gyres") then
278  call wind_forcing_gyres(sfc_state, forces, day_center, g, us, cs)
279  elseif (trim(cs%wind_config) == "zero") then
280  call wind_forcing_const(sfc_state, forces, 0., 0., day_center, g, us, cs)
281  elseif (trim(cs%wind_config) == "const") then
282  call wind_forcing_const(sfc_state, forces, cs%tau_x0, cs%tau_y0, day_center, g, us, cs)
283  elseif (trim(cs%wind_config) == "Neverworld" .or. trim(cs%wind_config) == "Neverland") then
284  call neverworld_wind_forcing(sfc_state, forces, day_center, g, us, cs)
285  elseif (trim(cs%wind_config) == "scurves") then
286  call scurve_wind_forcing(sfc_state, forces, day_center, g, us, cs)
287  elseif (trim(cs%wind_config) == "ideal_hurr") then
288  call idealized_hurricane_wind_forcing(sfc_state, forces, day_center, g, us, cs%idealized_hurricane_CSp)
289  elseif (trim(cs%wind_config) == "SCM_ideal_hurr") then
290  call scm_idealized_hurricane_wind_forcing(sfc_state, forces, day_center, g, us, cs%idealized_hurricane_CSp)
291  elseif (trim(cs%wind_config) == "SCM_CVmix_tests") then
292  call scm_cvmix_tests_wind_forcing(sfc_state, forces, day_center, g, us, cs%SCM_CVmix_tests_CSp)
293  elseif (trim(cs%wind_config) == "USER") then
294  call user_wind_forcing(sfc_state, forces, day_center, g, us, cs%user_forcing_CSp)
295  elseif (cs%variable_winds .and. .not.cs%first_call_set_forcing) then
296  call mom_error(fatal, &
297  "MOM_surface_forcing: Variable winds defined with no wind config")
298  else
299  call mom_error(fatal, &
300  "MOM_surface_forcing:Unrecognized wind config "//trim(cs%wind_config))
301  endif
302  endif
303 
304  ! calls to various buoyancy forcing options
305  if (cs%restorebuoy .and. .not.cs%variable_buoyforce) then
306  call mom_error(fatal, "With RESTOREBUOY = True, VARIABLE_BUOYFORCE = True should be used. "//&
307  "Otherwise, this can lead to diverging solutions when a simulation "//&
308  "is continued using a restart file.")
309  endif
310 
311  if ((cs%variable_buoyforce .or. cs%first_call_set_forcing) .and. &
312  (.not.cs%adiabatic)) then
313  if (trim(cs%buoy_config) == "file") then
314  call buoyancy_forcing_from_files(sfc_state, fluxes, day_center, dt, g, us, cs)
315  elseif (trim(cs%buoy_config) == "data_override") then
316  call buoyancy_forcing_from_data_override(sfc_state, fluxes, day_center, dt, g, us, cs)
317  elseif (trim(cs%buoy_config) == "zero") then
318  call buoyancy_forcing_zero(sfc_state, fluxes, day_center, dt, g, cs)
319  elseif (trim(cs%buoy_config) == "const") then
320  call buoyancy_forcing_const(sfc_state, fluxes, day_center, dt, g, us, cs)
321  elseif (trim(cs%buoy_config) == "linear") then
322  call buoyancy_forcing_linear(sfc_state, fluxes, day_center, dt, g, us, cs)
323  elseif (trim(cs%buoy_config) == "MESO") then
324  call meso_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%MESO_forcing_CSp)
325  elseif (trim(cs%buoy_config) == "SCM_CVmix_tests") then
326  call scm_cvmix_tests_buoyancy_forcing(sfc_state, fluxes, day_center, g, us, cs%SCM_CVmix_tests_CSp)
327  elseif (trim(cs%buoy_config) == "USER") then
328  call user_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%user_forcing_CSp)
329  elseif (trim(cs%buoy_config) == "BFB") then
330  call bfb_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%BFB_forcing_CSp)
331  elseif (trim(cs%buoy_config) == "dumbbell") then
332  call dumbbell_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%dumbbell_forcing_CSp)
333  elseif (trim(cs%buoy_config) == "NONE") then
334  call mom_mesg("MOM_surface_forcing: buoyancy forcing has been set to omitted.")
335  elseif (cs%variable_buoyforce .and. .not.cs%first_call_set_forcing) then
336  call mom_error(fatal, &
337  "MOM_surface_forcing: Variable buoy defined with no buoy config.")
338  else
339  call mom_error(fatal, &
340  "MOM_surface_forcing: Unrecognized buoy config "//trim(cs%buoy_config))
341  endif
342  endif
343 
344  if (associated(cs%tracer_flow_CSp)) then
345  call call_tracer_set_forcing(sfc_state, fluxes, day_start, day_interval, g, cs%tracer_flow_CSp)
346  endif
347 
348  ! Allow for user-written code to alter the fluxes after all the above
349  call user_alter_forcing(sfc_state, fluxes, day_center, g, cs%urf_CS)
350 
351  ! Fields that exist in both the forcing and mech_forcing types must be copied.
352  if (cs%variable_winds .or. cs%first_call_set_forcing) then
353  call copy_common_forcing_fields(forces, fluxes, g)
354  call set_derived_forcing_fields(forces, fluxes, g, us, cs%Rho0)
355  endif
356 
357  if ((cs%variable_buoyforce .or. cs%first_call_set_forcing) .and. &
358  (.not.cs%adiabatic)) then
359  call set_net_mass_forcing(fluxes, forces, g, us)
360  endif
361 
362  cs%first_call_set_forcing = .false.
363 
364  call cpu_clock_end(id_clock_forcing)
365  call calltree_leave("set_forcing")
366 
367 end subroutine set_forcing
368 
369 !> Sets the surface wind stresses to constant values
370 subroutine wind_forcing_const(sfc_state, forces, tau_x0, tau_y0, day, G, US, CS)
371  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
372  !! describe the surface state of the ocean.
373  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
374  real, intent(in) :: tau_x0 !< The zonal wind stress [Pa]
375  real, intent(in) :: tau_y0 !< The meridional wind stress [Pa]
376  type(time_type), intent(in) :: day !< The time of the fluxes
377  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
378  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
379  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
380  !! a previous surface_forcing_init call
381  ! Local variables
382  real :: Pa_conversion ! A unit conversion factor from Pa to the internal units [R Z L T-2 Pa-1 ~> 1]
383  real :: mag_tau
384  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
385 
386  call calltree_enter("wind_forcing_const, MOM_surface_forcing.F90")
387  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
388  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
389  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
390 
391  !set steady surface wind stresses, in units of Pa.
392  mag_tau = pa_conversion * sqrt( tau_x0**2 + tau_y0**2)
393 
394  do j=js,je ; do i=is-1,ieq
395  forces%taux(i,j) = tau_x0 * pa_conversion
396  enddo ; enddo
397 
398  do j=js-1,jeq ; do i=is,ie
399  forces%tauy(i,j) = tau_y0 * pa_conversion
400  enddo ; enddo
401 
402  if (cs%read_gust_2d) then
403  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
404  forces%ustar(i,j) = sqrt( us%L_to_Z * ( mag_tau + cs%gust(i,j) ) / cs%Rho0 )
405  enddo ; enddo ; endif
406  else
407  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
408  forces%ustar(i,j) = sqrt( us%L_to_Z * ( mag_tau + cs%gust_const ) / cs%Rho0 )
409  enddo ; enddo ; endif
410  endif
411 
412  call calltree_leave("wind_forcing_const")
413 end subroutine wind_forcing_const
414 
415 
416 !> Sets the surface wind stresses to set up two idealized gyres.
417 subroutine wind_forcing_2gyre(sfc_state, forces, day, G, US, CS)
418  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
419  !! describe the surface state of the ocean.
420  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
421  type(time_type), intent(in) :: day !< The time of the fluxes
422  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
423  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
424  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
425  !! a previous surface_forcing_init call
426  ! Local variables
427  real :: PI
428  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
429 
430  call calltree_enter("wind_forcing_2gyre, MOM_surface_forcing.F90")
431  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
432  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
433 
434  !set the steady surface wind stresses, in units of Pa.
435  pi = 4.0*atan(1.0)
436 
437  do j=js,je ; do i=is-1,ieq
438  forces%taux(i,j) = 0.1*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
439  (1.0 - cos(2.0*pi*(g%geoLatCu(i,j)-cs%South_lat) / cs%len_lat))
440  enddo ; enddo
441 
442  do j=js-1,jeq ; do i=is,ie
443  forces%tauy(i,j) = 0.0
444  enddo ; enddo
445 
446  call calltree_leave("wind_forcing_2gyre")
447 end subroutine wind_forcing_2gyre
448 
449 
450 !> Sets the surface wind stresses to set up a single idealized gyre.
451 subroutine wind_forcing_1gyre(sfc_state, forces, day, G, US, CS)
452  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
453  !! describe the surface state of the ocean.
454  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
455  type(time_type), intent(in) :: day !< The time of the fluxes
456  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
457  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
458  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
459  !! a previous surface_forcing_init call
460  ! Local variables
461  real :: PI
462  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
463 
464  call calltree_enter("wind_forcing_1gyre, MOM_surface_forcing.F90")
465  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
466  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
467 
468  ! set the steady surface wind stresses, in units of Pa.
469  pi = 4.0*atan(1.0)
470 
471  do j=js,je ; do i=is-1,ieq
472  forces%taux(i,j) = -0.2*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
473  cos(pi*(g%geoLatCu(i,j)-cs%South_lat)/cs%len_lat)
474  enddo ; enddo
475 
476  do j=js-1,jeq ; do i=is,ie
477  forces%tauy(i,j) = 0.0
478  enddo ; enddo
479 
480  call calltree_leave("wind_forcing_1gyre")
481 end subroutine wind_forcing_1gyre
482 
483 !> Sets the surface wind stresses to set up idealized gyres.
484 subroutine wind_forcing_gyres(sfc_state, forces, day, G, US, CS)
485  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
486  !! describe the surface state of the ocean.
487  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
488  type(time_type), intent(in) :: day !< The time of the fluxes
489  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
490  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
491  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
492  !! a previous surface_forcing_init call
493  ! Local variables
494  real :: PI, y, I_rho
495  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
496 
497  call calltree_enter("wind_forcing_gyres, MOM_surface_forcing.F90")
498  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
499  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
500 
501  ! steady surface wind stresses [Pa]
502  pi = 4.0*atan(1.0)
503 
504  do j=js-1,je+1 ; do i=is-1,ieq
505  y = (g%geoLatCu(i,j)-cs%South_lat) / cs%len_lat
506  forces%taux(i,j) = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
507  (cs%gyres_taux_const + &
508  ( cs%gyres_taux_sin_amp*sin(cs%gyres_taux_n_pis*pi*y) &
509  + cs%gyres_taux_cos_amp*cos(cs%gyres_taux_n_pis*pi*y) ))
510  enddo ; enddo
511 
512  do j=js-1,jeq ; do i=is-1,ie+1
513  forces%tauy(i,j) = 0.0
514  enddo ; enddo
515 
516  ! set the friction velocity
517  if (cs%answers_2018) then
518  do j=js,je ; do i=is,ie
519  forces%ustar(i,j) = sqrt(us%L_to_Z * ((cs%gust_const/cs%Rho0) + &
520  sqrt(0.5*(forces%tauy(i,j-1)*forces%tauy(i,j-1) + forces%tauy(i,j)*forces%tauy(i,j) + &
521  forces%taux(i-1,j)*forces%taux(i-1,j) + forces%taux(i,j)*forces%taux(i,j)))/cs%Rho0) )
522  enddo ; enddo
523  else
524  i_rho = us%L_to_Z / cs%Rho0
525  do j=js,je ; do i=is,ie
526  forces%ustar(i,j) = sqrt( (cs%gust_const + &
527  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
528  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
529  enddo ; enddo
530  endif
531 
532  call calltree_leave("wind_forcing_gyres")
533 end subroutine wind_forcing_gyres
534 
535 !> Sets the surface wind stresses, forces%taux and forces%tauy for the
536 !! Neverworld forcing configuration.
537 subroutine neverworld_wind_forcing(sfc_state, forces, day, G, US, CS)
538  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
539  !! describe the surface state of the ocean.
540  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
541  type(time_type), intent(in) :: day !< Time used for determining the fluxes.
542  type(ocean_grid_type), intent(inout) :: G !< Grid structure.
543  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
544  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
545  !! a previous surface_forcing_init call
546  ! Local variables
547  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
548  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
549  real :: PI, I_rho, y
550  real :: tau_max ! The magnitude of the wind stress [R Z L T-2 ~> Pa]
551  real :: off
552 
553  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
554  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
555  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
556  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
557 
558  ! Allocate the forcing arrays, if necessary.
559  call allocate_mech_forcing(g, forces, stress=.true.)
560 
561  ! Set the surface wind stresses, in units of Pa. A positive taux
562  ! accelerates the ocean to the (pseudo-)east.
563 
564  ! The i-loop extends to is-1 so that taux can be used later in the
565  ! calculation of ustar - otherwise the lower bound would be Isq.
566  pi = 4.0*atan(1.0)
567  forces%taux(:,:) = 0.0
568  tau_max = 0.2 * us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
569  off = 0.02
570  do j=js,je ; do i=is-1,ieq
571  y = (g%geoLatT(i,j)-g%south_lat)/g%len_lat
572 
573  if (y <= 0.29) then
574  forces%taux(i,j) = forces%taux(i,j) + tau_max * ( (1/0.29)*y - ( 1/(2*pi) )*sin( (2*pi*y) / 0.29 ) )
575  endif
576  if ((y > 0.29) .and. (y <= (0.8-off))) then
577  forces%taux(i,j) = forces%taux(i,j) + tau_max *(0.35+0.65*cos(pi*(y-0.29)/(0.51-off)) )
578  endif
579  if ((y > (0.8-off)) .and. (y <= (1-off))) then
580  forces%taux(i,j) = forces%taux(i,j) + tau_max *( 1.5*( (y-1+off) - (0.1/pi)*sin(10.0*pi*(y-0.8+off)) ) )
581  endif
582  forces%taux(i,j) = g%mask2dCu(i,j) * forces%taux(i,j)
583  enddo ; enddo
584 
585  do j=js-1,jeq ; do i=is,ie
586  forces%tauy(i,j) = g%mask2dCv(i,j) * 0.0
587  enddo ; enddo
588 
589  ! Set the surface friction velocity, in units of m s-1. ustar is always positive.
590  if (associated(forces%ustar)) then
591  i_rho = us%L_to_Z / cs%Rho0
592  do j=js,je ; do i=is,ie
593  forces%ustar(i,j) = sqrt( (cs%gust_const + &
594  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
595  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
596  enddo ; enddo
597  endif
598 
599 end subroutine neverworld_wind_forcing
600 
601 !> Sets the zonal wind stresses to a piecewise series of s-curves.
602 subroutine scurve_wind_forcing(sfc_state, forces, day, G, US, CS)
603  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
604  !! describe the surface state of the ocean.
605  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
606  type(time_type), intent(in) :: day !< Time used for determining the fluxes.
607  type(ocean_grid_type), intent(inout) :: G !< Grid structure.
608  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
609  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
610  !! a previous surface_forcing_init call
611  ! Local variables
612  integer :: i, j, kseg
613  real :: lon, lat, I_rho, y, L
614 ! real :: ydata(7) = (/ -70., -45., -15., 0., 15., 45., 70. /)
615 ! real :: taudt(7) = (/ 0., 0.2, -0.1, -0.02, -0.1, 0.1, 0. /)
616 
617  ! Allocate the forcing arrays, if necessary.
618  call allocate_mech_forcing(g, forces, stress=.true.)
619 
620  kseg = 1
621  do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB
622  lon = g%geoLonCu(i,j)
623  lat = g%geoLatCu(i,j)
624 
625  ! Find segment k s.t. ydata(k)<= lat < ydata(k+1)
626  do while (lat>=cs%scurves_ydata(kseg+1) .and. kseg<6)
627  kseg = kseg+1
628  enddo
629  do while (lat<cs%scurves_ydata(kseg) .and. kseg>1)
630  kseg = kseg-1
631  enddo
632 
633  y = lat - cs%scurves_ydata(kseg)
634  l = cs%scurves_ydata(kseg+1) - cs%scurves_ydata(kseg)
635  forces%taux(i,j) = cs%scurves_taux(kseg) + &
636  ( cs%scurves_taux(kseg+1) - cs%scurves_taux(kseg) ) * scurve(y, l)
637  forces%taux(i,j) = g%mask2dCu(i,j) * forces%taux(i,j)
638  enddo ; enddo
639 
640  do j=g%JsdB,g%JedB ; do i=g%isd,g%ied
641  forces%tauy(i,j) = g%mask2dCv(i,j) * 0.0
642  enddo ; enddo
643 
644  ! Set the surface friction velocity, in units of m s-1. ustar is always positive.
645  if (associated(forces%ustar)) then
646  i_rho = us%L_to_Z / cs%Rho0
647  do j=g%jsc,g%jec ; do i=g%isc,g%iec
648  forces%ustar(i,j) = sqrt( (cs%gust_const + &
649  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
650  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
651  enddo ; enddo
652  endif
653 
654 end subroutine scurve_wind_forcing
655 
656 !> Returns the value of a cosine-bell function evaluated at x/L
657 real function scurve(x,L)
658  real , intent(in) :: x !< non-dimensional position
659  real , intent(in) :: L !< non-dimensional width
660  real :: s
661 
662  s = x/l
663  scurve = (3. - 2.*s) * (s*s)
664 end function scurve
665 
666 ! Sets the surface wind stresses from input files.
667 subroutine wind_forcing_from_file(sfc_state, forces, day, G, US, CS)
668  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
669  !! describe the surface state of the ocean.
670  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
671  type(time_type), intent(in) :: day !< The time of the fluxes
672  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
673  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
674  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
675  !! a previous surface_forcing_init call
676  ! Local variables
677  character(len=200) :: filename ! The name of the input file.
678  real :: temp_x(szi_(g),szj_(g)) ! Pseudo-zonal and psuedo-meridional
679  real :: temp_y(szi_(g),szj_(g)) ! wind stresses at h-points [R L Z T-1 ~> Pa].
680  real :: Pa_conversion ! A unit conversion factor from Pa to the internal wind stress
681  ! units [R Z L T-2 Pa-1 ~> 1]
682  integer :: time_lev_daily ! The time levels to read for fields with
683  integer :: time_lev_monthly ! daily and montly cycles.
684  integer :: time_lev ! The time level that is used for a field.
685  integer :: days, seconds
686  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
687  logical :: read_Ustar
688 
689  call calltree_enter("wind_forcing_from_file, MOM_surface_forcing.F90")
690  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
691  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
692  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
693 
694  call get_time(day, seconds, days)
695  time_lev_daily = days - 365*floor(real(days) / 365.0)
696 
697  if (time_lev_daily < 31) then ; time_lev_monthly = 0
698  elseif (time_lev_daily < 59) then ; time_lev_monthly = 1
699  elseif (time_lev_daily < 90) then ; time_lev_monthly = 2
700  elseif (time_lev_daily < 120) then ; time_lev_monthly = 3
701  elseif (time_lev_daily < 151) then ; time_lev_monthly = 4
702  elseif (time_lev_daily < 181) then ; time_lev_monthly = 5
703  elseif (time_lev_daily < 212) then ; time_lev_monthly = 6
704  elseif (time_lev_daily < 243) then ; time_lev_monthly = 7
705  elseif (time_lev_daily < 273) then ; time_lev_monthly = 8
706  elseif (time_lev_daily < 304) then ; time_lev_monthly = 9
707  elseif (time_lev_daily < 334) then ; time_lev_monthly = 10
708  else ; time_lev_monthly = 11
709  endif
710 
711  time_lev_daily = time_lev_daily+1
712  time_lev_monthly = time_lev_monthly+1
713 
714  select case (cs%wind_nlev)
715  case (12) ; time_lev = time_lev_monthly
716  case (365) ; time_lev = time_lev_daily
717  case default ; time_lev = 1
718  end select
719 
720  if (time_lev /= cs%wind_last_lev) then
721  filename = trim(cs%wind_file)
722  read_ustar = (len_trim(cs%ustar_var) > 0)
723 ! if (is_root_pe()) &
724 ! write(*,'("Wind_forcing Reading time level ",I," last was ",I,".")')&
725 ! time_lev-1,CS%wind_last_lev-1
726  select case ( uppercase(cs%wind_stagger(1:1)) )
727  case ("A")
728  temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
729  call mom_read_vector(filename, cs%stress_x_var, cs%stress_y_var, &
730  temp_x(:,:), temp_y(:,:), g%Domain, stagger=agrid, &
731  timelevel=time_lev, scale=pa_conversion)
732 
733  call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
734  do j=js,je ; do i=is-1,ieq
735  forces%taux(i,j) = 0.5 * cs%wind_scale * (temp_x(i,j) + temp_x(i+1,j))
736  enddo ; enddo
737  do j=js-1,jeq ; do i=is,ie
738  forces%tauy(i,j) = 0.5 * cs%wind_scale * (temp_y(i,j) + temp_y(i,j+1))
739  enddo ; enddo
740 
741  if (.not.read_ustar) then
742  if (cs%read_gust_2d) then
743  do j=js,je ; do i=is,ie
744  forces%ustar(i,j) = sqrt((cs%gust(i,j) + &
745  sqrt(temp_x(i,j)*temp_x(i,j) + temp_y(i,j)*temp_y(i,j))) * us%L_to_Z / cs%Rho0)
746  enddo ; enddo
747  else
748  do j=js,je ; do i=is,ie
749  forces%ustar(i,j) = sqrt(us%L_to_Z * (cs%gust_const/cs%Rho0 + &
750  sqrt(temp_x(i,j)*temp_x(i,j) + temp_y(i,j)*temp_y(i,j)) / cs%Rho0) )
751  enddo ; enddo
752  endif
753  endif
754  case ("C")
755  if (g%symmetric) then
756  if (.not.associated(g%Domain_aux)) call mom_error(fatal, &
757  " wind_forcing_from_file with C-grid input and symmetric memory "//&
758  " called with a non-associated auxiliary domain in the grid type.")
759  ! Read the data as though symmetric memory were not being used, and
760  ! then translate it appropriately.
761  temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
762  call mom_read_vector(filename, cs%stress_x_var, cs%stress_y_var, &
763  temp_x(:,:), temp_y(:,:), &
764  g%Domain_aux, stagger=cgrid_ne, timelevel=time_lev, &
765  scale=pa_conversion)
766  do j=js,je ; do i=is,ie
767  forces%taux(i,j) = cs%wind_scale * temp_x(i,j)
768  forces%tauy(i,j) = cs%wind_scale * temp_y(i,j)
769  enddo ; enddo
770  call fill_symmetric_edges(forces%taux, forces%tauy, g%Domain, stagger=cgrid_ne)
771  else
772  call mom_read_vector(filename, cs%stress_x_var, cs%stress_y_var, &
773  forces%taux(:,:), forces%tauy(:,:), &
774  g%Domain, stagger=cgrid_ne, timelevel=time_lev, &
775  scale=pa_conversion)
776 
777  if (cs%wind_scale /= 1.0) then
778  do j=js,je ; do i=isq,ieq
779  forces%taux(i,j) = cs%wind_scale * forces%taux(i,j)
780  enddo ; enddo
781  do j=jsq,jeq ; do i=is,ie
782  forces%tauy(i,j) = cs%wind_scale * forces%tauy(i,j)
783  enddo ; enddo
784  endif
785  endif
786 
787  call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
788  if (.not.read_ustar) then
789  if (cs%read_gust_2d) then
790  do j=js, je ; do i=is, ie
791  forces%ustar(i,j) = sqrt((cs%gust(i,j) + &
792  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
793  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * us%L_to_Z / cs%Rho0 )
794  enddo ; enddo
795  else
796  do j=js, je ; do i=is, ie
797  forces%ustar(i,j) = sqrt(us%L_to_Z * ( (cs%gust_const/cs%Rho0) + &
798  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
799  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2)))/cs%Rho0))
800  enddo ; enddo
801  endif
802  endif
803  case default
804  call mom_error(fatal, "wind_forcing_from_file: Unrecognized stagger "//&
805  trim(cs%wind_stagger)//" is not 'A' or 'C'.")
806  end select
807 
808  if (read_ustar) then
809  call mom_read_data(filename, cs%Ustar_var, forces%ustar(:,:), &
810  g%Domain, timelevel=time_lev, scale=us%m_to_Z*us%T_to_s)
811  endif
812 
813  cs%wind_last_lev = time_lev
814 
815  endif ! time_lev /= CS%wind_last_lev
816 
817  call calltree_leave("wind_forcing_from_file")
818 end subroutine wind_forcing_from_file
819 
820 
821 ! Sets the surface wind stresses via the data override facility.
822 subroutine wind_forcing_by_data_override(sfc_state, forces, day, G, US, CS)
823  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
824  !! describe the surface state of the ocean.
825  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
826  type(time_type), intent(in) :: day !< The time of the fluxes
827  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
828  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
829  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
830  !! a previous surface_forcing_init call
831  ! Local variables
832  real :: temp_x(szi_(g),szj_(g)) ! Pseudo-zonal and psuedo-meridional
833  real :: temp_y(szi_(g),szj_(g)) ! wind stresses at h-points [Pa].
834  real :: temp_ustar(szi_(g),szj_(g)) ! ustar [m s-1] (not rescaled).
835  real :: Pa_conversion ! A unit conversion factor from Pa to the internal units [R Z L T-2 Pa-1 ~> 1]
836  integer :: i, j, is_in, ie_in, js_in, je_in
837  logical :: read_uStar
838 
839  call calltree_enter("wind_forcing_by_data_override, MOM_surface_forcing.F90")
840 
841  if (.not.cs%dataOverrideIsInitialized) then
842  call allocate_mech_forcing(g, forces, stress=.true., ustar=.true., press=.true.)
843  call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
844  cs%dataOverrideIsInitialized = .true.
845  endif
846 
847  is_in = g%isc - g%isd + 1 ; ie_in = g%iec - g%isd + 1
848  js_in = g%jsc - g%jsd + 1 ; je_in = g%jec - g%jsd + 1
849  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
850 
851  temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
852  call data_override('OCN', 'taux', temp_x, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
853  call data_override('OCN', 'tauy', temp_y, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
854  call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
855  ! Ignore CS%wind_scale when using data_override ?????
856  do j=g%jsc,g%jec ; do i=g%isc-1,g%IecB
857  forces%taux(i,j) = pa_conversion * 0.5 * (temp_x(i,j) + temp_x(i+1,j))
858  enddo ; enddo
859  do j=g%jsc-1,g%JecB ; do i=g%isc,g%iec
860  forces%tauy(i,j) = pa_conversion * 0.5 * (temp_y(i,j) + temp_y(i,j+1))
861  enddo ; enddo
862 
863  read_ustar = (len_trim(cs%ustar_var) > 0) ! Need better control higher up ????
864  if (read_ustar) then
865  do j=g%jsc,g%jec ; do i=g%isc,g%iec ; temp_ustar(i,j) = us%Z_to_m*us%s_to_T*forces%ustar(i,j) ; enddo ; enddo
866  call data_override('OCN', 'ustar', temp_ustar, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
867  do j=g%jsc,g%jec ; do i=g%isc,g%iec ; forces%ustar(i,j) = us%m_to_Z*us%T_to_s*temp_ustar(i,j) ; enddo ; enddo
868  else
869  if (cs%read_gust_2d) then
870  call data_override('OCN', 'gust', cs%gust, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
871  do j=g%jsc,g%jec ; do i=g%isc,g%iec
872  forces%ustar(i,j) = sqrt((pa_conversion * sqrt(temp_x(i,j)*temp_x(i,j) + &
873  temp_y(i,j)*temp_y(i,j)) + cs%gust(i,j)) * us%L_to_Z / cs%Rho0)
874  enddo ; enddo
875  else
876  do j=g%jsc,g%jec ; do i=g%isc,g%iec
877  forces%ustar(i,j) = sqrt(us%L_to_Z * (pa_conversion*sqrt(temp_x(i,j)*temp_x(i,j) + &
878  temp_y(i,j)*temp_y(i,j))/cs%Rho0 + cs%gust_const/cs%Rho0 ))
879  enddo ; enddo
880  endif
881  endif
882 
883  call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
884 ! call pass_var(forces%ustar, G%Domain, To_All) Not needed ?????
885 
886  call calltree_leave("wind_forcing_by_data_override")
887 end subroutine wind_forcing_by_data_override
888 
889 
890 !> Specifies zero surface bouyancy fluxes from input files.
891 subroutine buoyancy_forcing_from_files(sfc_state, fluxes, day, dt, G, US, CS)
892  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
893  !! describe the surface state of the ocean.
894  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
895  type(time_type), intent(in) :: day !< The time of the fluxes
896  real, intent(in) :: dt !< The amount of time over which
897  !! the fluxes apply [s]
898  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
899  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
900  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
901  !! a previous surface_forcing_init call
902  ! Local variables
903  real, dimension(SZI_(G),SZJ_(G)) :: &
904  temp, & ! A 2-d temporary work array with various units.
905  SST_anom, & ! Instantaneous sea surface temperature anomalies from a
906  ! target (observed) value [degC].
907  sss_anom, & ! Instantaneous sea surface salinity anomalies from a target
908  ! (observed) value [ppt].
909  sss_mean ! A (mean?) salinity about which to normalize local salinity
910  ! anomalies when calculating restorative precipitation
911  ! anomalies [ppt].
912 
913  real :: kg_m2_s_conversion ! A combination of unit conversion factors for rescaling
914  ! mass fluxes [R Z s m2 kg-1 T-1 ~> 1].
915  real :: rhoXcp ! reference density times heat capacity [Q R degC-1 ~> J m-3 degC-1]
916 
917  integer :: time_lev_daily ! time levels to read for fields with daily cycle
918  integer :: time_lev_monthly ! time levels to read for fields with monthly cycle
919  integer :: time_lev ! time level that for a field
920 
921  integer :: days, seconds
922  integer :: i, j, is, ie, js, je
923 
924  call calltree_enter("buoyancy_forcing_from_files, MOM_surface_forcing.F90")
925 
926  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
927  kg_m2_s_conversion = us%kg_m2s_to_RZ_T
928 
929  if (cs%use_temperature) rhoxcp = cs%Rho0 * fluxes%C_p
930 
931  ! Read the buoyancy forcing file
932  call get_time(day, seconds, days)
933 
934  time_lev_daily = days - 365*floor(real(days) / 365.0)
935 
936  if (time_lev_daily < 31) then ; time_lev_monthly = 0
937  elseif (time_lev_daily < 59) then ; time_lev_monthly = 1
938  elseif (time_lev_daily < 90) then ; time_lev_monthly = 2
939  elseif (time_lev_daily < 120) then ; time_lev_monthly = 3
940  elseif (time_lev_daily < 151) then ; time_lev_monthly = 4
941  elseif (time_lev_daily < 181) then ; time_lev_monthly = 5
942  elseif (time_lev_daily < 212) then ; time_lev_monthly = 6
943  elseif (time_lev_daily < 243) then ; time_lev_monthly = 7
944  elseif (time_lev_daily < 273) then ; time_lev_monthly = 8
945  elseif (time_lev_daily < 304) then ; time_lev_monthly = 9
946  elseif (time_lev_daily < 334) then ; time_lev_monthly = 10
947  else ; time_lev_monthly = 11
948  endif
949 
950  time_lev_daily = time_lev_daily +1
951  time_lev_monthly = time_lev_monthly+1
952 
953  if (time_lev_daily /= cs%buoy_last_lev_read) then
954 
955  ! longwave
956  select case (cs%LW_nlev)
957  case (12) ; time_lev = time_lev_monthly
958  case (365) ; time_lev = time_lev_daily
959  case default ; time_lev = 1
960  end select
961  call mom_read_data(cs%longwave_file, cs%LW_var, fluxes%lw(:,:), &
962  g%Domain, timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
963  if (cs%archaic_OMIP_file) then
964  call mom_read_data(cs%longwaveup_file, "lwup_sfc", temp(:,:), g%Domain, &
965  timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
966  do j=js,je ; do i=is,ie ; fluxes%LW(i,j) = fluxes%LW(i,j) - temp(i,j) ; enddo ; enddo
967  endif
968  cs%LW_last_lev = time_lev
969 
970  ! evaporation
971  select case (cs%evap_nlev)
972  case (12) ; time_lev = time_lev_monthly
973  case (365) ; time_lev = time_lev_daily
974  case default ; time_lev = 1
975  end select
976  if (cs%archaic_OMIP_file) then
977  call mom_read_data(cs%evaporation_file, cs%evap_var, fluxes%evap(:,:), &
978  g%Domain, timelevel=time_lev, scale=-kg_m2_s_conversion)
979  do j=js,je ; do i=is,ie
980  fluxes%latent(i,j) = cs%latent_heat_vapor*fluxes%evap(i,j)
981  fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
982  enddo ; enddo
983  else
984  call mom_read_data(cs%evaporation_file, cs%evap_var, fluxes%evap(:,:), &
985  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
986  endif
987  cs%evap_last_lev = time_lev
988 
989  select case (cs%latent_nlev)
990  case (12) ; time_lev = time_lev_monthly
991  case (365) ; time_lev = time_lev_daily
992  case default ; time_lev = 1
993  end select
994  if (.not.cs%archaic_OMIP_file) then
995  call mom_read_data(cs%latentheat_file, cs%latent_var, fluxes%latent(:,:), &
996  g%Domain, timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
997  do j=js,je ; do i=is,ie
998  fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
999  enddo ; enddo
1000  endif
1001  cs%latent_last_lev = time_lev
1002 
1003  select case (cs%sens_nlev)
1004  case (12) ; time_lev = time_lev_monthly
1005  case (365) ; time_lev = time_lev_daily
1006  case default ; time_lev = 1
1007  end select
1008  if (cs%archaic_OMIP_file) then
1009  call mom_read_data(cs%sensibleheat_file, cs%sens_var, fluxes%sens(:,:), &
1010  g%Domain, timelevel=time_lev, scale=-us%W_m2_to_QRZ_T)
1011  else
1012  call mom_read_data(cs%sensibleheat_file, cs%sens_var, fluxes%sens(:,:), &
1013  g%Domain, timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
1014  endif
1015  cs%sens_last_lev = time_lev
1016 
1017  select case (cs%SW_nlev)
1018  case (12) ; time_lev = time_lev_monthly
1019  case (365) ; time_lev = time_lev_daily
1020  case default ; time_lev = 1
1021  end select
1022  call mom_read_data(cs%shortwave_file, cs%SW_var, fluxes%sw(:,:), g%Domain, &
1023  timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
1024  if (cs%archaic_OMIP_file) then
1025  call mom_read_data(cs%shortwaveup_file, "swup_sfc", temp(:,:), g%Domain, &
1026  timelevel=time_lev, scale=us%W_m2_to_QRZ_T)
1027  do j=js,je ; do i=is,ie
1028  fluxes%sw(i,j) = fluxes%sw(i,j) - temp(i,j)
1029  enddo ; enddo
1030  endif
1031  cs%SW_last_lev = time_lev
1032 
1033  select case (cs%precip_nlev)
1034  case (12) ; time_lev = time_lev_monthly
1035  case (365) ; time_lev = time_lev_daily
1036  case default ; time_lev = 1
1037  end select
1038  call mom_read_data(cs%snow_file, cs%snow_var, &
1039  fluxes%fprec(:,:), g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1040  call mom_read_data(cs%rain_file, cs%rain_var, &
1041  fluxes%lprec(:,:), g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1042  if (cs%archaic_OMIP_file) then
1043  do j=js,je ; do i=is,ie
1044  fluxes%lprec(i,j) = fluxes%lprec(i,j) - fluxes%fprec(i,j)
1045  enddo ; enddo
1046  endif
1047  cs%precip_last_lev = time_lev
1048 
1049  select case (cs%runoff_nlev)
1050  case (12) ; time_lev = time_lev_monthly
1051  case (365) ; time_lev = time_lev_daily
1052  case default ; time_lev = 1
1053  end select
1054  if (cs%archaic_OMIP_file) then
1055  call mom_read_data(cs%runoff_file, cs%lrunoff_var, temp(:,:), &
1056  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1057  do j=js,je ; do i=is,ie
1058  fluxes%lrunoff(i,j) = temp(i,j)*us%m_to_L**2*g%IareaT(i,j)
1059  enddo ; enddo
1060  call mom_read_data(cs%runoff_file, cs%frunoff_var, temp(:,:), &
1061  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1062  do j=js,je ; do i=is,ie
1063  fluxes%frunoff(i,j) = temp(i,j)*us%m_to_L**2*g%IareaT(i,j)
1064  enddo ; enddo
1065  else
1066  call mom_read_data(cs%runoff_file, cs%lrunoff_var, fluxes%lrunoff(:,:), &
1067  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1068  call mom_read_data(cs%runoff_file, cs%frunoff_var, fluxes%frunoff(:,:), &
1069  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
1070  endif
1071  cs%runoff_last_lev = time_lev
1072 
1073 ! Read the SST and SSS fields for damping.
1074  if (cs%restorebuoy) then !#CTRL# .or. associated(CS%ctrl_forcing_CSp)) then
1075  select case (cs%SST_nlev)
1076  case (12) ; time_lev = time_lev_monthly
1077  case (365) ; time_lev = time_lev_daily
1078  case default ; time_lev = 1
1079  end select
1080  call mom_read_data(cs%SSTrestore_file, cs%SST_restore_var, &
1081  cs%T_Restore(:,:), g%Domain, timelevel=time_lev)
1082  cs%SST_last_lev = time_lev
1083 
1084  select case (cs%SSS_nlev)
1085  case (12) ; time_lev = time_lev_monthly
1086  case (365) ; time_lev = time_lev_daily
1087  case default ; time_lev = 1
1088  end select
1089  call mom_read_data(cs%salinityrestore_file, cs%SSS_restore_var, &
1090  cs%S_Restore(:,:), g%Domain, timelevel=time_lev)
1091  cs%SSS_last_lev = time_lev
1092  endif
1093  cs%buoy_last_lev_read = time_lev_daily
1094 
1095  ! mask out land points and compute heat content of water fluxes
1096  ! assume liquid precip enters ocean at SST
1097  ! assume frozen precip enters ocean at 0degC
1098  ! assume liquid runoff enters ocean at SST
1099  ! assume solid runoff (calving) enters ocean at 0degC
1100  ! mass leaving the ocean has heat_content determined in MOM_diabatic_driver.F90
1101  do j=js,je ; do i=is,ie
1102  fluxes%evap(i,j) = fluxes%evap(i,j) * g%mask2dT(i,j)
1103  fluxes%lprec(i,j) = fluxes%lprec(i,j) * g%mask2dT(i,j)
1104  fluxes%fprec(i,j) = fluxes%fprec(i,j) * g%mask2dT(i,j)
1105  fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * g%mask2dT(i,j)
1106  fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * g%mask2dT(i,j)
1107  fluxes%lw(i,j) = fluxes%lw(i,j) * g%mask2dT(i,j)
1108  fluxes%sens(i,j) = fluxes%sens(i,j) * g%mask2dT(i,j)
1109  fluxes%sw(i,j) = fluxes%sw(i,j) * g%mask2dT(i,j)
1110  fluxes%latent(i,j) = fluxes%latent(i,j) * g%mask2dT(i,j)
1111 
1112  fluxes%latent_evap_diag(i,j) = fluxes%latent_evap_diag(i,j) * g%mask2dT(i,j)
1113  fluxes%latent_fprec_diag(i,j) = -fluxes%fprec(i,j)*cs%latent_heat_fusion
1114  fluxes%latent_frunoff_diag(i,j) = -fluxes%frunoff(i,j)*cs%latent_heat_fusion
1115  enddo ; enddo
1116 
1117  endif ! time_lev /= CS%buoy_last_lev_read
1118 
1119 
1120  ! restoring surface boundary fluxes
1121  if (cs%restorebuoy) then
1122 
1123  if (cs%use_temperature) then
1124  do j=js,je ; do i=is,ie
1125  if (g%mask2dT(i,j) > 0) then
1126  fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1127  ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const_T)
1128  fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const_S) * &
1129  (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
1130  (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
1131  else
1132  fluxes%heat_added(i,j) = 0.0
1133  fluxes%vprec(i,j) = 0.0
1134  endif
1135  enddo ; enddo
1136  else
1137  do j=js,je ; do i=is,ie
1138  if (g%mask2dT(i,j) > 0) then
1139  fluxes%buoy(i,j) = (cs%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1140  (cs%G_Earth * cs%Flux_const / cs%Rho0)
1141  else
1142  fluxes%buoy(i,j) = 0.0
1143  endif
1144  enddo ; enddo
1145  endif
1146 
1147  else ! not RESTOREBUOY
1148  if (.not.cs%use_temperature) then
1149  call mom_error(fatal, "buoyancy_forcing in MOM_surface_forcing: "// &
1150  "The fluxes need to be defined without RESTOREBUOY.")
1151  endif
1152 
1153  endif ! end RESTOREBUOY
1154 
1155 !#CTRL# if (associated(CS%ctrl_forcing_CSp)) then
1156 !#CTRL# do j=js,je ; do i=is,ie
1157 !#CTRL# SST_anom(i,j) = sfc_state%SST(i,j) - CS%T_Restore(i,j)
1158 !#CTRL# SSS_anom(i,j) = sfc_state%SSS(i,j) - CS%S_Restore(i,j)
1159 !#CTRL# SSS_mean(i,j) = 0.5*(sfc_state%SSS(i,j) + CS%S_Restore(i,j))
1160 !#CTRL# enddo ; enddo
1161 !#CTRL# call apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, fluxes%heat_added, &
1162 !#CTRL# fluxes%vprec, day, dt, G, CS%ctrl_forcing_CSp)
1163 !#CTRL# endif
1164 
1165  call calltree_leave("buoyancy_forcing_from_files")
1166 end subroutine buoyancy_forcing_from_files
1167 
1168 !> Specifies zero surface bouyancy fluxes from data over-ride.
1169 subroutine buoyancy_forcing_from_data_override(sfc_state, fluxes, day, dt, G, US, CS)
1170  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1171  !! describe the surface state of the ocean.
1172  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1173  type(time_type), intent(in) :: day !< The time of the fluxes
1174  real, intent(in) :: dt !< The amount of time over which
1175  !! the fluxes apply [s]
1176  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1177  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1178  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1179  !! a previous surface_forcing_init call
1180  ! Local variables
1181  real, dimension(SZI_(G),SZJ_(G)) :: &
1182  temp, & ! A 2-d temporary work array with various units.
1183  SST_anom, & ! Instantaneous sea surface temperature anomalies from a
1184  ! target (observed) value [degC].
1185  sss_anom, & ! Instantaneous sea surface salinity anomalies from a target
1186  ! (observed) value [ppt].
1187  sss_mean ! A (mean?) salinity about which to normalize local salinity
1188  ! anomalies when calculating restorative precipitation
1189  ! anomalies [ppt].
1190  real :: kg_m2_s_conversion ! A combination of unit conversion factors for rescaling
1191  ! mass fluxes [R Z s m2 kg-1 T-1 ~> 1].
1192  real :: rhoXcp ! The mean density times the heat capacity [Q R degC-1 ~> J m-3 degC-1].
1193 
1194  integer :: time_lev_daily ! The time levels to read for fields with
1195  integer :: time_lev_monthly ! daily and montly cycles.
1196  integer :: itime_lev ! The time level that is used for a field.
1197 
1198  integer :: days, seconds
1199  integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
1200  integer :: is_in, ie_in, js_in, je_in
1201 
1202  call calltree_enter("buoyancy_forcing_from_data_override, MOM_surface_forcing.F90")
1203 
1204  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1205  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1206  kg_m2_s_conversion = us%kg_m2s_to_RZ_T
1207 
1208  if (cs%use_temperature) rhoxcp = cs%Rho0 * fluxes%C_p
1209 
1210  if (.not.cs%dataOverrideIsInitialized) then
1211  call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
1212  cs%dataOverrideIsInitialized = .true.
1213  endif
1214 
1215  is_in = g%isc - g%isd + 1
1216  ie_in = g%iec - g%isd + 1
1217  js_in = g%jsc - g%jsd + 1
1218  je_in = g%jec - g%jsd + 1
1219 
1220  call data_override('OCN', 'lw', fluxes%lw(:,:), day, &
1221  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=US%W_m2_to_QRZ_T
1222  if (us%QRZ_T_to_W_m2 /= 1.0) then ; do j=js,je ; do i=is,ie
1223  fluxes%lw(i,j) = fluxes%lw(i,j) * us%W_m2_to_QRZ_T
1224  enddo ; enddo ; endif
1225  call data_override('OCN', 'evap', fluxes%evap(:,:), day, &
1226  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1227 
1228  ! note the sign convention
1229  do j=js,je ; do i=is,ie
1230  ! The normal convention is that fluxes%evap positive into the ocean
1231  ! but evap is normally a positive quantity in the files
1232  ! This conversion is dangerous because it is not clear whether the data files have been read!
1233  fluxes%evap(i,j) = -kg_m2_s_conversion*fluxes%evap(i,j)
1234  fluxes%latent(i,j) = cs%latent_heat_vapor*fluxes%evap(i,j)
1235  fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
1236  enddo ; enddo
1237 
1238  call data_override('OCN', 'sens', fluxes%sens(:,:), day, &
1239  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1240 
1241  ! note the sign convention
1242  do j=js,je ; do i=is,ie
1243  fluxes%sens(i,j) = -us%W_m2_to_QRZ_T * fluxes%sens(i,j) ! Normal convention is positive into the ocean
1244  ! but sensible is normally a positive quantity in the files
1245  enddo ; enddo
1246 
1247  call data_override('OCN', 'sw', fluxes%sw(:,:), day, &
1248  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=US%W_m2_to_QRZ_T
1249  if (us%QRZ_T_to_W_m2 /= 1.0) then ; do j=js,je ; do i=is,ie
1250  fluxes%sw(i,j) = fluxes%sw(i,j) * us%W_m2_to_QRZ_T
1251  enddo ; enddo ; endif
1252 
1253  call data_override('OCN', 'snow', fluxes%fprec(:,:), day, &
1254  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1255 
1256  call data_override('OCN', 'rain', fluxes%lprec(:,:), day, &
1257  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1258 
1259  call data_override('OCN', 'runoff', fluxes%lrunoff(:,:), day, &
1260  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1261 
1262  call data_override('OCN', 'calving', fluxes%frunoff(:,:), day, &
1263  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1264 
1265  if (kg_m2_s_conversion /= 1.0) then ; do j=js,je ; do i=is,ie
1266  fluxes%lprec(i,j) = fluxes%lprec(i,j) * kg_m2_s_conversion
1267  fluxes%fprec(i,j) = fluxes%fprec(i,j) * kg_m2_s_conversion
1268  fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * kg_m2_s_conversion
1269  fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * kg_m2_s_conversion
1270  enddo ; enddo ; endif
1271 
1272 ! Read the SST and SSS fields for damping.
1273  if (cs%restorebuoy) then !#CTRL# .or. associated(CS%ctrl_forcing_CSp)) then
1274  call data_override('OCN', 'SST_restore', cs%T_restore(:,:), day, &
1275  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1276 
1277  call data_override('OCN', 'SSS_restore', cs%S_restore(:,:), day, &
1278  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1279 
1280  endif
1281 
1282  ! restoring boundary fluxes
1283  if (cs%restorebuoy) then
1284  if (cs%use_temperature) then
1285  do j=js,je ; do i=is,ie
1286  if (g%mask2dT(i,j) > 0) then
1287  fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1288  ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const_T)
1289  fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const_S) * &
1290  (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
1291  (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
1292  else
1293  fluxes%heat_added(i,j) = 0.0
1294  fluxes%vprec(i,j) = 0.0
1295  endif
1296  enddo ; enddo
1297  else
1298  do j=js,je ; do i=is,ie
1299  if (g%mask2dT(i,j) > 0) then
1300  fluxes%buoy(i,j) = (cs%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1301  (cs%G_Earth * cs%Flux_const / cs%Rho0)
1302  else
1303  fluxes%buoy(i,j) = 0.0
1304  endif
1305  enddo ; enddo
1306  endif
1307  else ! not RESTOREBUOY
1308  if (.not.cs%use_temperature) then
1309  call mom_error(fatal, "buoyancy_forcing in MOM_surface_forcing: "// &
1310  "The fluxes need to be defined without RESTOREBUOY.")
1311  endif
1312  endif ! end RESTOREBUOY
1313 
1314 
1315  ! mask out land points and compute heat content of water fluxes
1316  ! assume liquid precip enters ocean at SST
1317  ! assume frozen precip enters ocean at 0degC
1318  ! assume liquid runoff enters ocean at SST
1319  ! assume solid runoff (calving) enters ocean at 0degC
1320  ! mass leaving ocean has heat_content determined in MOM_diabatic_driver.F90
1321  do j=js,je ; do i=is,ie
1322  fluxes%evap(i,j) = fluxes%evap(i,j) * g%mask2dT(i,j)
1323  fluxes%lprec(i,j) = fluxes%lprec(i,j) * g%mask2dT(i,j)
1324  fluxes%fprec(i,j) = fluxes%fprec(i,j) * g%mask2dT(i,j)
1325  fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * g%mask2dT(i,j)
1326  fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * g%mask2dT(i,j)
1327  fluxes%lw(i,j) = fluxes%lw(i,j) * g%mask2dT(i,j)
1328  fluxes%latent(i,j) = fluxes%latent(i,j) * g%mask2dT(i,j)
1329  fluxes%sens(i,j) = fluxes%sens(i,j) * g%mask2dT(i,j)
1330  fluxes%sw(i,j) = fluxes%sw(i,j) * g%mask2dT(i,j)
1331 
1332  fluxes%latent_evap_diag(i,j) = fluxes%latent_evap_diag(i,j) * g%mask2dT(i,j)
1333  fluxes%latent_fprec_diag(i,j) = -fluxes%fprec(i,j)*cs%latent_heat_fusion
1334  fluxes%latent_frunoff_diag(i,j) = -fluxes%frunoff(i,j)*cs%latent_heat_fusion
1335  enddo ; enddo
1336 
1337 
1338 !#CTRL# if (associated(CS%ctrl_forcing_CSp)) then
1339 !#CTRL# do j=js,je ; do i=is,ie
1340 !#CTRL# SST_anom(i,j) = sfc_state%SST(i,j) - CS%T_Restore(i,j)
1341 !#CTRL# SSS_anom(i,j) = sfc_state%SSS(i,j) - CS%S_Restore(i,j)
1342 !#CTRL# SSS_mean(i,j) = 0.5*(sfc_state%SSS(i,j) + CS%S_Restore(i,j))
1343 !#CTRL# enddo ; enddo
1344 !#CTRL# call apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, fluxes%heat_added, &
1345 !#CTRL# fluxes%vprec, day, dt, G, CS%ctrl_forcing_CSp)
1346 !#CTRL# endif
1347 
1348  call calltree_leave("buoyancy_forcing_from_data_override")
1349 end subroutine buoyancy_forcing_from_data_override
1350 
1351 !> This subroutine specifies zero surface bouyancy fluxes
1352 subroutine buoyancy_forcing_zero(sfc_state, fluxes, day, dt, G, CS)
1353  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1354  !! describe the surface state of the ocean.
1355  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1356  type(time_type), intent(in) :: day !< The time of the fluxes
1357  real, intent(in) :: dt !< The amount of time over which
1358  !! the fluxes apply [s]
1359  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1360  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1361  !! a previous surface_forcing_init call
1362  ! Local variables
1363  integer :: i, j, is, ie, js, je
1364 
1365  call calltree_enter("buoyancy_forcing_zero, MOM_surface_forcing.F90")
1366  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1367 
1368  if (cs%use_temperature) then
1369  do j=js,je ; do i=is,ie
1370  fluxes%evap(i,j) = 0.0
1371  fluxes%lprec(i,j) = 0.0
1372  fluxes%fprec(i,j) = 0.0
1373  fluxes%vprec(i,j) = 0.0
1374  fluxes%lrunoff(i,j) = 0.0
1375  fluxes%frunoff(i,j) = 0.0
1376  fluxes%lw(i,j) = 0.0
1377  fluxes%latent(i,j) = 0.0
1378  fluxes%sens(i,j) = 0.0
1379  fluxes%sw(i,j) = 0.0
1380  fluxes%latent_evap_diag(i,j) = 0.0
1381  fluxes%latent_fprec_diag(i,j) = 0.0
1382  fluxes%latent_frunoff_diag(i,j) = 0.0
1383  enddo ; enddo
1384  else
1385  do j=js,je ; do i=is,ie
1386  fluxes%buoy(i,j) = 0.0
1387  enddo ; enddo
1388  endif
1389 
1390  call calltree_leave("buoyancy_forcing_zero")
1391 end subroutine buoyancy_forcing_zero
1392 
1393 
1394 !> Sets up spatially and temporally constant surface heat fluxes.
1395 subroutine buoyancy_forcing_const(sfc_state, fluxes, day, dt, G, US, CS)
1396  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1397  !! describe the surface state of the ocean.
1398  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1399  type(time_type), intent(in) :: day !< The time of the fluxes
1400  real, intent(in) :: dt !< The amount of time over which
1401  !! the fluxes apply [s]
1402  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1403  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1404  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1405  !! a previous surface_forcing_init call
1406  ! Local variables
1407  integer :: i, j, is, ie, js, je
1408  call calltree_enter("buoyancy_forcing_const, MOM_surface_forcing.F90")
1409  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1410 
1411  if (cs%use_temperature) then
1412  do j=js,je ; do i=is,ie
1413  fluxes%evap(i,j) = 0.0
1414  fluxes%lprec(i,j) = 0.0
1415  fluxes%fprec(i,j) = 0.0
1416  fluxes%vprec(i,j) = 0.0
1417  fluxes%lrunoff(i,j) = 0.0
1418  fluxes%frunoff(i,j) = 0.0
1419  fluxes%lw(i,j) = 0.0
1420  fluxes%latent(i,j) = 0.0
1421  fluxes%sens(i,j) = cs%constantHeatForcing * g%mask2dT(i,j)
1422  fluxes%sw(i,j) = 0.0
1423  fluxes%latent_evap_diag(i,j) = 0.0
1424  fluxes%latent_fprec_diag(i,j) = 0.0
1425  fluxes%latent_frunoff_diag(i,j) = 0.0
1426  enddo ; enddo
1427  else
1428  do j=js,je ; do i=is,ie
1429  fluxes%buoy(i,j) = 0.0
1430  enddo ; enddo
1431  endif
1432 
1433  call calltree_leave("buoyancy_forcing_const")
1434 end subroutine buoyancy_forcing_const
1435 
1436 !> Sets surface fluxes of heat and salinity by restoring to temperature and
1437 !! salinity profiles that vary linearly with latitude.
1438 subroutine buoyancy_forcing_linear(sfc_state, fluxes, day, dt, G, US, CS)
1439  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1440  !! describe the surface state of the ocean.
1441  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1442  type(time_type), intent(in) :: day !< The time of the fluxes
1443  real, intent(in) :: dt !< The amount of time over which
1444  !! the fluxes apply [s]
1445  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1446  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1447  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1448  !! a previous surface_forcing_init call
1449  ! Local variables
1450  real :: y, T_restore, S_restore
1451  integer :: i, j, is, ie, js, je
1452 
1453  call calltree_enter("buoyancy_forcing_linear, MOM_surface_forcing.F90")
1454  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1455 
1456  ! This case has no surface buoyancy forcing.
1457  if (cs%use_temperature) then
1458  do j=js,je ; do i=is,ie
1459  fluxes%evap(i,j) = 0.0
1460  fluxes%lprec(i,j) = 0.0
1461  fluxes%fprec(i,j) = 0.0
1462  fluxes%vprec(i,j) = 0.0
1463  fluxes%lrunoff(i,j) = 0.0
1464  fluxes%frunoff(i,j) = 0.0
1465  fluxes%lw(i,j) = 0.0
1466  fluxes%latent(i,j) = 0.0
1467  fluxes%sens(i,j) = 0.0
1468  fluxes%sw(i,j) = 0.0
1469  fluxes%latent_evap_diag(i,j) = 0.0
1470  fluxes%latent_fprec_diag(i,j) = 0.0
1471  fluxes%latent_frunoff_diag(i,j) = 0.0
1472  enddo ; enddo
1473  else
1474  do j=js,je ; do i=is,ie
1475  fluxes%buoy(i,j) = 0.0
1476  enddo ; enddo
1477  endif
1478 
1479  if (cs%restorebuoy) then
1480  if (cs%use_temperature) then
1481  do j=js,je ; do i=is,ie
1482  y = (g%geoLatCu(i,j)-cs%South_lat)/cs%len_lat
1483  t_restore = cs%T_south + (cs%T_north-cs%T_south)*y
1484  s_restore = cs%S_south + (cs%S_north-cs%S_south)*y
1485  if (g%mask2dT(i,j) > 0) then
1486  fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1487  ((t_restore - sfc_state%SST(i,j)) * ((cs%Rho0 * fluxes%C_p) * cs%Flux_const))
1488  fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const) * &
1489  (s_restore - sfc_state%SSS(i,j)) / &
1490  (0.5*(sfc_state%SSS(i,j) + s_restore))
1491  else
1492  fluxes%heat_added(i,j) = 0.0
1493  fluxes%vprec(i,j) = 0.0
1494  endif
1495  enddo ; enddo
1496  else
1497  call mom_error(fatal, "buoyancy_forcing_linear in MOM_surface_forcing: "// &
1498  "RESTOREBUOY to linear not written yet.")
1499  !do j=js,je ; do i=is,ie
1500  ! if (G%mask2dT(i,j) > 0) then
1501  ! fluxes%buoy(i,j) = (CS%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1502  ! (CS%G_Earth * CS%Flux_const / CS%Rho0)
1503  ! else
1504  ! fluxes%buoy(i,j) = 0.0
1505  ! endif
1506  !enddo ; enddo
1507  endif
1508  else ! not RESTOREBUOY
1509  if (.not.cs%use_temperature) then
1510  call mom_error(fatal, "buoyancy_forcing_linear in MOM_surface_forcing: "// &
1511  "The fluxes need to be defined without RESTOREBUOY.")
1512  endif
1513  endif ! end RESTOREBUOY
1514 
1515  call calltree_leave("buoyancy_forcing_linear")
1516 end subroutine buoyancy_forcing_linear
1517 
1518 !> Save a restart file for the forcing fields
1519 subroutine forcing_save_restart(CS, G, Time, directory, time_stamped, &
1520  filename_suffix)
1521  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1522  !! a previous surface_forcing_init call
1523  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1524  type(time_type), intent(in) :: Time !< model time at this call; needed for mpp_write calls
1525  character(len=*), intent(in) :: directory !< directory into which to write these restart files
1526  logical, optional, intent(in) :: time_stamped !< If true, the restart file names
1527  !! include a unique time stamp; the default is false.
1528  character(len=*), optional, intent(in) :: filename_suffix !< optional suffix (e.g., a time-stamp)
1529  !! to append to the restart fname
1530 
1531  if (.not.associated(cs)) return
1532  if (.not.associated(cs%restart_CSp)) return
1533 
1534  call save_restart(directory, time, g, cs%restart_CSp, time_stamped)
1535 
1536 end subroutine forcing_save_restart
1537 
1538 !> Initialize the surface forcing module
1539 subroutine surface_forcing_init(Time, G, US, param_file, diag, CS, tracer_flow_CSp)
1540  type(time_type), intent(in) :: Time !< The current model time
1541  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1542  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1543  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1544  type(diag_ctrl), target, intent(inout) :: diag !< structure used to regulate diagnostic output
1545  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1546  !! a previous surface_forcing_init call
1547  type(tracer_flow_control_cs), pointer :: tracer_flow_CSp !< Forcing for tracers?
1548 
1549  ! Local variables
1550  type(directories) :: dirs
1551  logical :: new_sim
1552  type(time_type) :: Time_frc
1553  ! This include declares and sets the variable "version".
1554 # include "version_variable.h"
1555  real :: flux_const_default ! The unscaled value of FLUXCONST [m day-1]
1556  logical :: default_2018_answers
1557  character(len=40) :: mdl = "MOM_surface_forcing" ! This module's name.
1558  character(len=200) :: filename, gust_file ! The name of the gustiness input file.
1559 
1560  if (associated(cs)) then
1561  call mom_error(warning, "surface_forcing_init called with an associated "// &
1562  "control structure.")
1563  return
1564  endif
1565  allocate(cs)
1566 
1567  id_clock_forcing=cpu_clock_id('(Ocean surface forcing)', grain=clock_module)
1568  call cpu_clock_begin(id_clock_forcing)
1569 
1570  cs%diag => diag
1571  if (associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
1572 
1573  ! Read all relevant parameters and write them to the model log.
1574  call log_version(param_file, mdl, version, '')
1575  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
1576  "If true, Temperature and salinity are used as state "//&
1577  "variables.", default=.true.)
1578  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, &
1579  "The directory in which all input files are found.", &
1580  default=".")
1581  cs%inputdir = slasher(cs%inputdir)
1582 
1583  call get_param(param_file, mdl, "ADIABATIC", cs%adiabatic, &
1584  "There are no diapycnal mass fluxes if ADIABATIC is "//&
1585  "true. This assumes that KD = KDML = 0.0 and that "//&
1586  "there is no buoyancy forcing, but makes the model "//&
1587  "faster by eliminating subroutine calls.", default=.false.)
1588  call get_param(param_file, mdl, "VARIABLE_WINDS", cs%variable_winds, &
1589  "If true, the winds vary in time after the initialization.", &
1590  default=.true.)
1591  call get_param(param_file, mdl, "VARIABLE_BUOYFORCE", cs%variable_buoyforce, &
1592  "If true, the buoyancy forcing varies in time after the "//&
1593  "initialization of the model.", default=.true.)
1594 
1595  call get_param(param_file, mdl, "BUOY_CONFIG", cs%buoy_config, &
1596  "The character string that indicates how buoyancy forcing "//&
1597  "is specified. Valid options include (file), (zero), "//&
1598  "(linear), (USER), (BFB) and (NONE).", default="zero")
1599  if (trim(cs%buoy_config) == "file") then
1600  call get_param(param_file, mdl, "ARCHAIC_OMIP_FORCING_FILE", cs%archaic_OMIP_file, &
1601  "If true, use the forcing variable decomposition from "//&
1602  "the old German OMIP prescription that predated CORE. If "//&
1603  "false, use the variable groupings available from MOM "//&
1604  "output diagnostics of forcing variables.", default=.true.)
1605  if (cs%archaic_OMIP_file) then
1606  call get_param(param_file, mdl, "LONGWAVEDOWN_FILE", cs%longwave_file, &
1607  "The file with the downward longwave heat flux, in "//&
1608  "variable lwdn_sfc.", fail_if_missing=.true.)
1609  call get_param(param_file, mdl, "LONGWAVEUP_FILE", cs%longwaveup_file, &
1610  "The file with the upward longwave heat flux, in "//&
1611  "variable lwup_sfc.", fail_if_missing=.true.)
1612  call get_param(param_file, mdl, "EVAPORATION_FILE", cs%evaporation_file, &
1613  "The file with the evaporative moisture flux, in "//&
1614  "variable evap.", fail_if_missing=.true.)
1615  call get_param(param_file, mdl, "SENSIBLEHEAT_FILE", cs%sensibleheat_file, &
1616  "The file with the sensible heat flux, in "//&
1617  "variable shflx.", fail_if_missing=.true.)
1618  call get_param(param_file, mdl, "SHORTWAVEUP_FILE", cs%shortwaveup_file, &
1619  "The file with the upward shortwave heat flux.", &
1620  fail_if_missing=.true.)
1621  call get_param(param_file, mdl, "SHORTWAVEDOWN_FILE", cs%shortwave_file, &
1622  "The file with the downward shortwave heat flux.", &
1623  fail_if_missing=.true.)
1624  call get_param(param_file, mdl, "SNOW_FILE", cs%snow_file, &
1625  "The file with the downward frozen precip flux, in "//&
1626  "variable snow.", fail_if_missing=.true.)
1627  call get_param(param_file, mdl, "PRECIP_FILE", cs%rain_file, &
1628  "The file with the downward total precip flux, in "//&
1629  "variable precip.", fail_if_missing=.true.)
1630  call get_param(param_file, mdl, "FRESHDISCHARGE_FILE", cs%runoff_file, &
1631  "The file with the fresh and frozen runoff/calving fluxes, "//&
1632  "invariables disch_w and disch_s.", fail_if_missing=.true.)
1633 
1634  ! These variable names are hard-coded, per the archaic OMIP conventions.
1635  cs%latentheat_file = cs%evaporation_file ; cs%latent_var = "evap"
1636  cs%LW_var = "lwdn_sfc"; cs%SW_var = "swdn_sfc"; cs%sens_var = "shflx"
1637  cs%evap_var = "evap"; cs%rain_var = "precip"; cs%snow_var = "snow"
1638  cs%lrunoff_var = "disch_w"; cs%frunoff_var = "disch_s"
1639 
1640  else
1641  call get_param(param_file, mdl, "LONGWAVE_FILE", cs%longwave_file, &
1642  "The file with the longwave heat flux, in the variable "//&
1643  "given by LONGWAVE_FORCING_VAR.", fail_if_missing=.true.)
1644  call get_param(param_file, mdl, "LONGWAVE_FORCING_VAR", cs%LW_var, &
1645  "The variable with the longwave forcing field.", default="LW")
1646 
1647  call get_param(param_file, mdl, "SHORTWAVE_FILE", cs%shortwave_file, &
1648  "The file with the shortwave heat flux, in the variable "//&
1649  "given by SHORTWAVE_FORCING_VAR.", fail_if_missing=.true.)
1650  call get_param(param_file, mdl, "SHORTWAVE_FORCING_VAR", cs%SW_var, &
1651  "The variable with the shortwave forcing field.", default="SW")
1652 
1653  call get_param(param_file, mdl, "EVAPORATION_FILE", cs%evaporation_file, &
1654  "The file with the evaporative moisture flux, in the "//&
1655  "variable given by EVAP_FORCING_VAR.", fail_if_missing=.true.)
1656  call get_param(param_file, mdl, "EVAP_FORCING_VAR", cs%evap_var, &
1657  "The variable with the evaporative moisture flux.", &
1658  default="evap")
1659 
1660  call get_param(param_file, mdl, "LATENTHEAT_FILE", cs%latentheat_file, &
1661  "The file with the latent heat flux, in the variable "//&
1662  "given by LATENT_FORCING_VAR.", fail_if_missing=.true.)
1663  call get_param(param_file, mdl, "LATENT_FORCING_VAR", cs%latent_var, &
1664  "The variable with the latent heat flux.", default="latent")
1665 
1666  call get_param(param_file, mdl, "SENSIBLEHEAT_FILE", cs%sensibleheat_file, &
1667  "The file with the sensible heat flux, in the variable "//&
1668  "given by SENSIBLE_FORCING_VAR.", fail_if_missing=.true.)
1669  call get_param(param_file, mdl, "SENSIBLE_FORCING_VAR", cs%sens_var, &
1670  "The variable with the sensible heat flux.", default="sensible")
1671 
1672  call get_param(param_file, mdl, "RAIN_FILE", cs%rain_file, &
1673  "The file with the liquid precipitation flux, in the "//&
1674  "variable given by RAIN_FORCING_VAR.", fail_if_missing=.true.)
1675  call get_param(param_file, mdl, "RAIN_FORCING_VAR", cs%rain_var, &
1676  "The variable with the liquid precipitation flux.", &
1677  default="liq_precip")
1678  call get_param(param_file, mdl, "SNOW_FILE", cs%snow_file, &
1679  "The file with the frozen precipitation flux, in the "//&
1680  "variable given by SNOW_FORCING_VAR.", fail_if_missing=.true.)
1681  call get_param(param_file, mdl, "SNOW_FORCING_VAR", cs%snow_var, &
1682  "The variable with the frozen precipitation flux.", &
1683  default="froz_precip")
1684 
1685  call get_param(param_file, mdl, "RUNOFF_FILE", cs%runoff_file, &
1686  "The file with the fresh and frozen runoff/calving "//&
1687  "fluxes, in variables given by LIQ_RUNOFF_FORCING_VAR "//&
1688  "and FROZ_RUNOFF_FORCING_VAR.", fail_if_missing=.true.)
1689  call get_param(param_file, mdl, "LIQ_RUNOFF_FORCING_VAR", cs%lrunoff_var, &
1690  "The variable with the liquid runoff flux.", &
1691  default="liq_runoff")
1692  call get_param(param_file, mdl, "FROZ_RUNOFF_FORCING_VAR", cs%frunoff_var, &
1693  "The variable with the frozen runoff flux.", &
1694  default="froz_runoff")
1695  endif
1696 
1697  call get_param(param_file, mdl, "SSTRESTORE_FILE", cs%SSTrestore_file, &
1698  "The file with the SST toward which to restore in the "//&
1699  "variable given by SST_RESTORE_VAR.", fail_if_missing=.true.)
1700  call get_param(param_file, mdl, "SALINITYRESTORE_FILE", cs%salinityrestore_file, &
1701  "The file with the surface salinity toward which to "//&
1702  "restore in the variable given by SSS_RESTORE_VAR.", &
1703  fail_if_missing=.true.)
1704 
1705  if (cs%archaic_OMIP_file) then
1706  cs%SST_restore_var = "TEMP" ; cs%SSS_restore_var = "SALT"
1707  else
1708  call get_param(param_file, mdl, "SST_RESTORE_VAR", cs%SST_restore_var, &
1709  "The variable with the SST toward which to restore.", &
1710  default="SST")
1711  call get_param(param_file, mdl, "SSS_RESTORE_VAR", cs%SSS_restore_var, &
1712  "The variable with the SSS toward which to restore.", &
1713  default="SSS")
1714  endif
1715 
1716  ! Add inputdir to the file names.
1717  cs%shortwave_file = trim(cs%inputdir)//trim(cs%shortwave_file)
1718  cs%longwave_file = trim(cs%inputdir)//trim(cs%longwave_file)
1719  cs%sensibleheat_file = trim(cs%inputdir)//trim(cs%sensibleheat_file)
1720  cs%latentheat_file = trim(cs%inputdir)//trim(cs%latentheat_file)
1721  cs%evaporation_file = trim(cs%inputdir)//trim(cs%evaporation_file)
1722  cs%snow_file = trim(cs%inputdir)//trim(cs%snow_file)
1723  cs%rain_file = trim(cs%inputdir)//trim(cs%rain_file)
1724  cs%runoff_file = trim(cs%inputdir)//trim(cs%runoff_file)
1725 
1726  cs%shortwaveup_file = trim(cs%inputdir)//trim(cs%shortwaveup_file)
1727  cs%longwaveup_file = trim(cs%inputdir)//trim(cs%longwaveup_file)
1728 
1729  cs%SSTrestore_file = trim(cs%inputdir)//trim(cs%SSTrestore_file)
1730  cs%salinityrestore_file = trim(cs%inputdir)//trim(cs%salinityrestore_file)
1731  elseif (trim(cs%buoy_config) == "const") then
1732  call get_param(param_file, mdl, "SENSIBLE_HEAT_FLUX", cs%constantHeatForcing, &
1733  "A constant heat forcing (positive into ocean) applied "//&
1734  "through the sensible heat flux field. ", &
1735  units='W/m2', scale=us%W_m2_to_QRZ_T, fail_if_missing=.true.)
1736  endif
1737  call get_param(param_file, mdl, "WIND_CONFIG", cs%wind_config, &
1738  "The character string that indicates how wind forcing "//&
1739  "is specified. Valid options include (file), (2gyre), "//&
1740  "(1gyre), (gyres), (zero), and (USER).", default="zero")
1741  if (trim(cs%wind_config) == "file") then
1742  call get_param(param_file, mdl, "WIND_FILE", cs%wind_file, &
1743  "The file in which the wind stresses are found in "//&
1744  "variables STRESS_X and STRESS_Y.", fail_if_missing=.true.)
1745  call get_param(param_file, mdl, "WINDSTRESS_X_VAR",cs%stress_x_var, &
1746  "The name of the x-wind stress variable in WIND_FILE.", &
1747  default="STRESS_X")
1748  call get_param(param_file, mdl, "WINDSTRESS_Y_VAR", cs%stress_y_var, &
1749  "The name of the y-wind stress variable in WIND_FILE.", &
1750  default="STRESS_Y")
1751  call get_param(param_file, mdl, "WIND_STAGGER",cs%wind_stagger, &
1752  "A character indicating how the wind stress components "//&
1753  "are staggered in WIND_FILE. This may be A or C for now.", &
1754  default="C")
1755  call get_param(param_file, mdl, "WINDSTRESS_SCALE", cs%wind_scale, &
1756  "A value by which the wind stresses in WIND_FILE are rescaled.", &
1757  default=1.0, units="nondim")
1758  call get_param(param_file, mdl, "USTAR_FORCING_VAR", cs%ustar_var, &
1759  "The name of the friction velocity variable in WIND_FILE "//&
1760  "or blank to get ustar from the wind stresses plus the "//&
1761  "gustiness.", default=" ", units="nondim")
1762  cs%wind_file = trim(cs%inputdir) // trim(cs%wind_file)
1763  endif
1764  if (trim(cs%wind_config) == "gyres") then
1765  call get_param(param_file, mdl, "TAUX_CONST", cs%gyres_taux_const, &
1766  "With the gyres wind_config, the constant offset in the "//&
1767  "zonal wind stress profile: "//&
1768  " A in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1769  units="Pa", default=0.0)
1770  call get_param(param_file, mdl, "TAUX_SIN_AMP",cs%gyres_taux_sin_amp, &
1771  "With the gyres wind_config, the sine amplitude in the "//&
1772  "zonal wind stress profile: "//&
1773  " B in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1774  units="Pa", default=0.0)
1775  call get_param(param_file, mdl, "TAUX_COS_AMP",cs%gyres_taux_cos_amp, &
1776  "With the gyres wind_config, the cosine amplitude in "//&
1777  "the zonal wind stress profile: "//&
1778  " C in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1779  units="Pa", default=0.0)
1780  call get_param(param_file, mdl, "TAUX_N_PIS",cs%gyres_taux_n_pis, &
1781  "With the gyres wind_config, the number of gyres in "//&
1782  "the zonal wind stress profile: "//&
1783  " n in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1784  units="nondim", default=0.0)
1785  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
1786  "This sets the default value for the various _2018_ANSWERS parameters.", &
1787  default=.false.)
1788  call get_param(param_file, mdl, "WIND_GYRES_2018_ANSWERS", cs%answers_2018, &
1789  "If true, use the order of arithmetic and expressions that recover the answers "//&
1790  "from the end of 2018. Otherwise, use expressions for the gyre friction velocities "//&
1791  "that are rotationally invariant and more likely to be the same between compilers.", &
1792  default=default_2018_answers)
1793  else
1794  cs%answers_2018 = .false.
1795  endif
1796  if (trim(cs%wind_config) == "scurves") then
1797  call get_param(param_file, mdl, "WIND_SCURVES_LATS", cs%scurves_ydata, &
1798  "A list of latitudes defining a piecewise scurve profile "//&
1799  "for zonal wind stress.", &
1800  units="degrees N", fail_if_missing=.true.)
1801  call get_param(param_file, mdl, "WIND_SCURVES_TAUX", cs%scurves_taux, &
1802  "A list of zonal wind stress values at latitudes "//&
1803  "WIND_SCURVES_LATS defining a piecewise scurve profile.", &
1804  units="Pa", fail_if_missing=.true.)
1805  endif
1806  if ((trim(cs%wind_config) == "2gyre") .or. &
1807  (trim(cs%wind_config) == "1gyre") .or. &
1808  (trim(cs%wind_config) == "gyres") .or. &
1809  (trim(cs%buoy_config) == "linear")) then
1810  cs%south_lat = g%south_lat
1811  cs%len_lat = g%len_lat
1812  endif
1813 
1814  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
1815  "The mean ocean density used with BOUSSINESQ true to "//&
1816  "calculate accelerations and the mass for conservation "//&
1817  "properties, or with BOUSSINSEQ false to convert some "//&
1818  "parameters from vertical units of m to kg m-2.", &
1819  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
1820  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
1821  "If true, the buoyancy fluxes drive the model back "//&
1822  "toward some specified surface state with a rate "//&
1823  "given by FLUXCONST.", default= .false.)
1824  call get_param(param_file, mdl, "LATENT_HEAT_FUSION", cs%latent_heat_fusion, &
1825  "The latent heat of fusion.", default=hlf, &
1826  units="J/kg", scale=us%J_kg_to_Q)
1827  call get_param(param_file, mdl, "LATENT_HEAT_VAPORIZATION", cs%latent_heat_vapor, &
1828  "The latent heat of fusion.", default=hlv, units="J/kg", scale=us%J_kg_to_Q)
1829  if (cs%restorebuoy) then
1830  ! These three variables use non-standard time units, but are rescaled as they are read.
1831  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
1832  "The constant that relates the restoring surface fluxes to the relative "//&
1833  "surface anomalies (akin to a piston velocity). Note the non-MKS units.", &
1834  default=0.0, units="m day-1", scale=us%m_to_Z*us%T_to_s/86400.0, &
1835  unscaled=flux_const_default)
1836 
1837  if (cs%use_temperature) then
1838  call get_param(param_file, mdl, "FLUXCONST_T", cs%Flux_const_T, &
1839  "The constant that relates the restoring surface temperature "//&
1840  "flux to the relative surface anomaly (akin to a piston "//&
1841  "velocity). Note the non-MKS units.", &
1842  units="m day-1", scale=us%m_to_Z*us%T_to_s/86400.0, &
1843  default=flux_const_default)
1844  call get_param(param_file, mdl, "FLUXCONST_S", cs%Flux_const_S, &
1845  "The constant that relates the restoring surface salinity "//&
1846  "flux to the relative surface anomaly (akin to a piston "//&
1847  "velocity). Note the non-MKS units.", &
1848  units="m day-1", scale=us%m_to_Z*us%T_to_s/86400.0, &
1849  default=flux_const_default)
1850  endif
1851 
1852  if (trim(cs%buoy_config) == "linear") then
1853  call get_param(param_file, mdl, "SST_NORTH", cs%T_north, &
1854  "With buoy_config linear, the sea surface temperature "//&
1855  "at the northern end of the domain toward which to "//&
1856  "to restore.", units="deg C", default=0.0)
1857  call get_param(param_file, mdl, "SST_SOUTH", cs%T_south, &
1858  "With buoy_config linear, the sea surface temperature "//&
1859  "at the southern end of the domain toward which to "//&
1860  "to restore.", units="deg C", default=0.0)
1861  call get_param(param_file, mdl, "SSS_NORTH", cs%S_north, &
1862  "With buoy_config linear, the sea surface salinity "//&
1863  "at the northern end of the domain toward which to "//&
1864  "to restore.", units="PSU", default=35.0)
1865  call get_param(param_file, mdl, "SSS_SOUTH", cs%S_south, &
1866  "With buoy_config linear, the sea surface salinity "//&
1867  "at the southern end of the domain toward which to "//&
1868  "to restore.", units="PSU", default=35.0)
1869  endif
1870  endif
1871  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
1872  "The gravitational acceleration of the Earth.", &
1873  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
1874 
1875  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
1876  "The background gustiness in the winds.", &
1877  units="Pa", default=0.0, scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1878  call get_param(param_file, mdl, "FIX_USTAR_GUSTLESS_BUG", cs%fix_ustar_gustless_bug, &
1879  "If true correct a bug in the time-averaging of the gustless wind friction velocity", &
1880  default=.true.)
1881  call get_param(param_file, mdl, "READ_GUST_2D", cs%read_gust_2d, &
1882  "If true, use a 2-dimensional gustiness supplied from "//&
1883  "an input file", default=.false.)
1884  if (cs%read_gust_2d) then
1885  call get_param(param_file, mdl, "GUST_2D_FILE", gust_file, &
1886  "The file in which the wind gustiness is found in "//&
1887  "variable gustiness.", fail_if_missing=.true.)
1888  call safe_alloc_ptr(cs%gust,g%isd,g%ied,g%jsd,g%jed)
1889  filename = trim(cs%inputdir) // trim(gust_file)
1890  call mom_read_data(filename,'gustiness',cs%gust,g%domain, timelevel=1, &
1891  scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z) ! units in file should be Pa
1892  endif
1893 
1894 ! All parameter settings are now known.
1895 
1896  if (trim(cs%wind_config) == "USER" .or. trim(cs%buoy_config) == "USER" ) then
1897  call user_surface_forcing_init(time, g, us, param_file, diag, cs%user_forcing_CSp)
1898  elseif (trim(cs%buoy_config) == "BFB" ) then
1899  call bfb_surface_forcing_init(time, g, us, param_file, diag, cs%BFB_forcing_CSp)
1900  elseif (trim(cs%buoy_config) == "dumbbell" ) then
1901  call dumbbell_surface_forcing_init(time, g, us, param_file, diag, cs%dumbbell_forcing_CSp)
1902  elseif (trim(cs%wind_config) == "MESO" .or. trim(cs%buoy_config) == "MESO" ) then
1903  call meso_surface_forcing_init(time, g, us, param_file, diag, cs%MESO_forcing_CSp)
1904  elseif (trim(cs%wind_config) == "ideal_hurr" .or.&
1905  trim(cs%wind_config) == "SCM_ideal_hurr") then
1906  call idealized_hurricane_wind_init(time, g, us, param_file, cs%idealized_hurricane_CSp)
1907  elseif (trim(cs%wind_config) == "const") then
1908  call get_param(param_file, mdl, "CONST_WIND_TAUX", cs%tau_x0, &
1909  "With wind_config const, this is the constant zonal "//&
1910  "wind-stress", units="Pa", fail_if_missing=.true.)
1911  call get_param(param_file, mdl, "CONST_WIND_TAUY", cs%tau_y0, &
1912  "With wind_config const, this is the constant meridional "//&
1913  "wind-stress", units="Pa", fail_if_missing=.true.)
1914  elseif (trim(cs%wind_config) == "SCM_CVmix_tests" .or. &
1915  trim(cs%buoy_config) == "SCM_CVmix_tests") then
1916  call scm_cvmix_tests_surface_forcing_init(time, g, param_file, cs%SCM_CVmix_tests_CSp)
1917  endif
1918 
1919  call register_forcing_type_diags(time, diag, us, cs%use_temperature, cs%handles)
1920 
1921  ! Set up any restart fields associated with the forcing.
1922  call restart_init(param_file, cs%restart_CSp, "MOM_forcing.res")
1923 !#CTRL# call register_ctrl_forcing_restarts(G, param_file, CS%ctrl_forcing_CSp, &
1924 !#CTRL# CS%restart_CSp)
1925  call restart_init_end(cs%restart_CSp)
1926 
1927  if (associated(cs%restart_CSp)) then
1928  call get_mom_input(dirs=dirs)
1929 
1930  new_sim = .false.
1931  if ((dirs%input_filename(1:1) == 'n') .and. &
1932  (len_trim(dirs%input_filename) == 1)) new_sim = .true.
1933  if (.not.new_sim) then
1934  call restore_state(dirs%input_filename, dirs%restart_input_dir, time_frc, &
1935  g, cs%restart_CSp)
1936  endif
1937  endif
1938 
1939  ! Determine how many time levels are in each forcing variable.
1940  if (trim(cs%buoy_config) == "file") then
1941  cs%SW_nlev = num_timelevels(cs%shortwave_file, cs%SW_var, min_dims=3)
1942  cs%LW_nlev = num_timelevels(cs%longwave_file, cs%LW_var, min_dims=3)
1943  cs%latent_nlev = num_timelevels(cs%latentheat_file, cs%latent_var, 3)
1944  cs%sens_nlev = num_timelevels(cs%sensibleheat_file, cs%sens_var, min_dims=3)
1945 
1946  cs%evap_nlev = num_timelevels(cs%evaporation_file, cs%evap_var, min_dims=3)
1947  cs%precip_nlev = num_timelevels(cs%rain_file, cs%rain_var, min_dims=3)
1948  cs%runoff_nlev = num_timelevels(cs%runoff_file, cs%lrunoff_var, 3)
1949 
1950  cs%SST_nlev = num_timelevels(cs%SSTrestore_file, cs%SST_restore_var, 3)
1951  cs%SSS_nlev = num_timelevels(cs%salinityrestore_file, cs%SSS_restore_var, 3)
1952  endif
1953 
1954  if (trim(cs%wind_config) == "file") &
1955  cs%wind_nlev = num_timelevels(cs%wind_file, cs%stress_x_var, min_dims=3)
1956 
1957 !#CTRL# call controlled_forcing_init(Time, G, param_file, diag, CS%ctrl_forcing_CSp)
1958 
1959  call user_revise_forcing_init(param_file, cs%urf_CS)
1960 
1961  call cpu_clock_end(id_clock_forcing)
1962 end subroutine surface_forcing_init
1963 
1964 
1965 !> Deallocate memory associated with the surface forcing module
1966 subroutine surface_forcing_end(CS, fluxes)
1967  type(surface_forcing_cs), pointer :: CS !< pointer to control struct returned by
1968  !! a previous surface_forcing_init call
1969  type(forcing), optional, intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1970 ! Arguments: CS - A pointer to the control structure returned by a previous
1971 ! call to surface_forcing_init, it will be deallocated here.
1972 ! (inout) fluxes - A structure containing pointers to any possible
1973 ! forcing fields. Unused fields have NULL ptrs.
1974 
1975  if (present(fluxes)) call deallocate_forcing_type(fluxes)
1976 
1977 !#CTRL# call controlled_forcing_end(CS%ctrl_forcing_CSp)
1978 
1979  if (associated(cs)) deallocate(cs)
1980  cs => null()
1981 
1982 end subroutine surface_forcing_end
1983 
1984 end module mom_surface_forcing
Control structure for BFB_surface_forcing.
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Container for parameters describing idealized wind structure.
Wraps the FMS time manager functions.
This module implements boundary forcing for MOM6.
Control structure for user_revise_forcing.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
Structure containing pointers to the forcing fields that may be used to drive MOM. All fluxes are positive into the ocean.
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Surface forcing for the dumbbell test case.
Reads the only Fortran name list needed to boot-strap the model.
Allocate the fields of a mechanical forcing type, based on either a set of input flags for each group...
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Control structure for the dumbbell test case forcing.
Wraps the MPP cpu clock functions.
Register fields for restarts.
Functions that calculate the surface wind stresses and fluxes of buoyancy or temperature/salinity and...
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
Do a set of halo updates that fill in the values at the duplicated edges of a staggered symmetric mem...
Definition: MOM_domains.F90:93
Do a halo update on a pair of arrays representing the two components of a vector. ...
Definition: MOM_domains.F90:59
Read a pair of data fields representing the two components of a vector from a file.
Definition: MOM_io.F90:82
Provides a few physical constants.
Orchestrates the registration and calling of tracer packages.
Structure that contains pointers to the mechanical forcing at the surface used to drive the liquid oc...
Allocate the fields of a (flux) forcing type, based on either a set of input flags for each group of ...
This control structure is used to store parameters associated with the MESO forcing.
Make a diagnostic available for averaging or output.
This control structure should be used to store any run-time variables associated with the user-specif...
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:75
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
Provides a template for users to code updating the forcing fluxes.
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Routines for error handling and I/O management.
The MOM6 facility for reading and writing restart files, and querying what has been read...
Definition: MOM_restart.F90:2
The control structure for orchestrating the calling of tracer packages.
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Surface forcing for the boundary-forced-basin (BFB) configuration.
An overloaded interface to log version information about modules.
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
Container for paths and parameter file names.
Initial conditions and forcing for the single column model (SCM) CVMix test set.
Handy functions for manipulating strings.
Template for user to code up surface forcing.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Sets forcing for the MESO configuration.
Structure that defines the id handles for the forcing type.
Forcing for the idealized hurricane and SCM_idealized_hurricane examples.
Do a halo update on an array.
Definition: MOM_domains.F90:54
Read a data field from a file.
Definition: MOM_io.F90:74
Container for surface forcing parameters.
An overloaded interface to read and log the values of various types of parameters.