MOM6
mom_surface_forcing_nuopc.F90
Go to the documentation of this file.
1 !> Converts the input ESMF data (import data) to a MOM-specific data type (surface_forcing_CS).
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_coms, only : reproducing_sum
7 use mom_constants, only : hlv, hlf
8 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
9 use mom_cpu_clock, only : clock_subcomponent
10 use mom_diag_mediator, only : diag_ctrl
11 use mom_diag_mediator, only : safe_alloc_ptr, time_type
13 use mom_domains, only : global_field_sum, bitwise_exact_sum
14 use mom_domains, only : agrid, bgrid_ne, cgrid_ne, to_all
15 use mom_domains, only : to_north, to_east, omit_corners
16 use mom_error_handler, only : mom_error, warning, fatal, is_root_pe, mom_mesg
23 use mom_grid, only : ocean_grid_type
24 use mom_io, only : slasher, write_version_number, mom_read_data
30 use mom_variables, only : surface
33 
34 use coupler_types_mod, only : coupler_2d_bc_type, coupler_type_write_chksums
35 use coupler_types_mod, only : coupler_type_initialized, coupler_type_spawn
36 use coupler_types_mod, only : coupler_type_copy_data
37 use data_override_mod, only : data_override_init, data_override
38 use fms_mod, only : stdout
39 use mpp_mod, only : mpp_chksum
40 use time_interp_external_mod, only : init_external_field, time_interp_external
41 use time_interp_external_mod, only : time_interp_external_init
42 
43 implicit none ; private
44 
45 #include <MOM_memory.h>
46 
52 
55 private surface_forcing_end
56 
57 !> Contains pointers to the forcing fields which may be used to drive MOM.
58 !! All fluxes are positive downward.
59 type, public :: surface_forcing_cs ; private
60  integer :: wind_stagger !< AGRID, BGRID_NE, or CGRID_NE (integer values
61  !! from MOM_domains) to indicate the staggering of
62  !! the winds that are being provided in calls to
63  !! update_ocean_model.
64  logical :: use_temperature !! If true, temp and saln used as state variables
65  real :: wind_stress_multiplier !< A multiplier applied to incoming wind stress (nondim).
66 
67  real :: rho0 !< Boussinesq reference density [R ~> kg m-3]
68  real :: area_surf = -1.0 !< total ocean surface area [m^2]
69  real :: latent_heat_fusion !< latent heat of fusion [J/kg]
70  real :: latent_heat_vapor !< latent heat of vaporization [J/kg]
71 
72  real :: max_p_surf !< maximum surface pressure that can be
73  !! exerted by the atmosphere and floating sea-ice,
74  !! in Pa. This is needed because the FMS coupling
75  !! structure does not limit the water that can be
76  !! frozen out of the ocean and the ice-ocean heat
77  !! fluxes are treated explicitly.
78  logical :: use_limited_p_ssh !< If true, return the sea surface height with
79  !! the correction for the atmospheric (and sea-ice)
80  !! pressure limited by max_p_surf instead of the
81  !! full atmospheric pressure. The default is true.
82 
83  real :: gust_const !< constant unresolved background gustiness for ustar [R L Z T-1 ~> Pa]
84  logical :: read_gust_2d !< If true, use a 2-dimensional gustiness supplied
85  !! from an input file.
86  real, pointer, dimension(:,:) :: &
87  tke_tidal => null(), & !< turbulent kinetic energy introduced to the
88  !! bottom boundary layer by drag on the tidal flows [R Z3 T-3 ~> W m-2]
89  gust => null(), & !< spatially varying unresolved background
90  !! gustiness that contributes to ustar [R L Z T-1 ~> Pa].
91  !! gust is used when read_gust_2d is true.
92  ustar_tidal => null() !< tidal contribution to the bottom friction velocity [Z T-1 ~> m s-1]
93  real :: cd_tides !< drag coefficient that applies to the tides (nondimensional)
94  real :: utide !< constant tidal velocity to use if read_tideamp
95  !! is false [Z T-1 ~> m s-1]
96  logical :: read_tideamp !< If true, spatially varying tidal amplitude read from a file.
97 
98  logical :: rigid_sea_ice !< If true, sea-ice exerts a rigidity that acts
99  !! to damp surface deflections (especially surface
100  !! gravity waves). The default is false.
101  real :: g_earth !< Gravitational acceleration [m s-2]
102  real :: kv_sea_ice !! viscosity in sea-ice that resists sheared vertical motions [m2 s-1]
103  real :: density_sea_ice !< typical density of sea-ice [kg m-3]. The value is
104  !! only used to convert the ice pressure into
105  !! appropriate units for use with Kv_sea_ice.
106  real :: rigid_sea_ice_mass !< A mass per unit area of sea-ice beyond which
107  !! sea-ice viscosity becomes effective, in kg m-2,
108  !! typically of order 1000 [kg m-2].
109  logical :: allow_flux_adjustments !< If true, use data_override to obtain flux adjustments
110  logical :: liquid_runoff_from_data !< If true, use data_override to obtain liquid runoff
111 
112  real :: flux_const !< piston velocity for surface restoring [m/s]
113  logical :: salt_restore_as_sflux !< If true, SSS restore as salt flux instead of water flux
114  logical :: adjust_net_srestore_to_zero !< adjust srestore to zero (for both salt_flux or vprec)
115  logical :: adjust_net_srestore_by_scaling !< adjust srestore w/o moving zero contour
116  logical :: adjust_net_fresh_water_to_zero !< adjust net surface fresh-water (w/ restoring) to zero
117  logical :: use_net_fw_adjustment_sign_bug !< use the wrong sign when adjusting net FW
118  logical :: adjust_net_fresh_water_by_scaling !< adjust net surface fresh-water w/o moving zero contour
119  logical :: mask_srestore_under_ice !< If true, use an ice mask defined by frazil
120  !< criteria for salinity restoring.
121  real :: ice_salt_concentration !< salt concentration for sea ice [kg/kg]
122  logical :: mask_srestore_marginal_seas !< if true, then mask SSS restoring in marginal seas
123  real :: max_delta_srestore !< maximum delta salinity used for restoring
124  real :: max_delta_trestore !< maximum delta sst used for restoring
125  real, pointer, dimension(:,:) :: basin_mask => null() !< mask for SSS restoring by basin
126 
127  type(diag_ctrl), pointer :: diag !< structure to regulate diagnostic output timing
128  character(len=200) :: inputdir !< directory where NetCDF input files are
129  character(len=200) :: salt_restore_file !< filename for salt restoring data
130  character(len=30) :: salt_restore_var_name !< name of surface salinity in salt_restore_file
131  logical :: mask_srestore !< if true, apply a 2-dimensional mask to the surface
132  !< salinity restoring fluxes. The masking file should be
133  !< in inputdir/salt_restore_mask.nc and the field should
134  !! be named 'mask'
135  real, pointer, dimension(:,:) :: srestore_mask => null() !< mask for SSS restoring
136  character(len=200) :: temp_restore_file !< filename for sst restoring data
137  character(len=30) :: temp_restore_var_name !< name of surface temperature in temp_restore_file
138  logical :: mask_trestore !< if true, apply a 2-dimensional mask to the surface
139  !! temperature restoring fluxes. The masking file should be
140  !! in inputdir/temp_restore_mask.nc and the field should
141  !! be named 'mask'
142  real, pointer, dimension(:,:) :: trestore_mask => null() !< mask for SST restoring
143  integer :: id_srestore = -1 !< id number for time_interp_external.
144  integer :: id_trestore = -1 !< id number for time_interp_external.
145 
146  ! Diagnostics handles
147  type(forcing_diags), public :: handles
148 
149  type(mom_restart_cs), pointer :: restart_csp => null()
150  type(user_revise_forcing_cs), pointer :: urf_cs => null()
151 end type surface_forcing_cs
152 
153 !> Structure corresponding to forcing, but with the elements, units, and conventions
154 !! that exactly conform to the use for MOM-based coupled models.
155 type, public :: ice_ocean_boundary_type
156  real, pointer, dimension(:,:) :: lrunoff =>null() !< liquid runoff [kg/m2/s]
157  real, pointer, dimension(:,:) :: frunoff =>null() !< ice runoff [kg/m2/s]
158  real, pointer, dimension(:,:) :: u_flux =>null() !< i-direction wind stress [Pa]
159  real, pointer, dimension(:,:) :: v_flux =>null() !< j-direction wind stress [Pa]
160  real, pointer, dimension(:,:) :: t_flux =>null() !< sensible heat flux [W/m2]
161  real, pointer, dimension(:,:) :: q_flux =>null() !< specific humidity flux [kg/m2/s]
162  real, pointer, dimension(:,:) :: salt_flux =>null() !< salt flux [kg/m2/s]
163  real, pointer, dimension(:,:) :: seaice_melt_heat =>null() !< sea ice and snow melt heat flux [W/m2]
164  real, pointer, dimension(:,:) :: seaice_melt =>null() !< water flux due to sea ice and snow melting [kg/m2/s]
165  real, pointer, dimension(:,:) :: lw_flux =>null() !< long wave radiation [W/m2]
166  real, pointer, dimension(:,:) :: sw_flux_vis_dir =>null() !< direct visible sw radiation [W/m2]
167  real, pointer, dimension(:,:) :: sw_flux_vis_dif =>null() !< diffuse visible sw radiation [W/m2]
168  real, pointer, dimension(:,:) :: sw_flux_nir_dir =>null() !< direct Near InfraRed sw radiation [W/m2]
169  real, pointer, dimension(:,:) :: sw_flux_nir_dif =>null() !< diffuse Near InfraRed sw radiation [W/m2]
170  real, pointer, dimension(:,:) :: lprec =>null() !< mass flux of liquid precip [kg/m2/s]
171  real, pointer, dimension(:,:) :: fprec =>null() !< mass flux of frozen precip [kg/m2/s]
172  real, pointer, dimension(:,:) :: ustar_berg =>null() !< frictional velocity beneath icebergs [m/s]
173  real, pointer, dimension(:,:) :: area_berg =>null() !< area covered by icebergs[m2/m2]
174  real, pointer, dimension(:,:) :: mass_berg =>null() !< mass of icebergs(kg/m2)
175  real, pointer, dimension(:,:) :: lrunoff_hflx =>null() !< heat content of liquid runoff [W/m2]
176  real, pointer, dimension(:,:) :: frunoff_hflx =>null() !< heat content of frozen runoff [W/m2]
177  real, pointer, dimension(:,:) :: p =>null() !< pressure of overlying ice and atmosphere
178  !< on ocean surface [Pa]
179  real, pointer, dimension(:,:) :: mi =>null() !< mass of ice [kg/m2]
180  real, pointer, dimension(:,:) :: ice_rigidity =>null() !< rigidity of the sea ice, sea-ice and
181  !! ice-shelves, expressed as a coefficient
182  !! for divergence damping, as determined
183  !! outside of the ocean model in [m3/s]
184  integer :: xtype !< The type of the exchange - REGRID, REDIST or DIRECT
185  type(coupler_2d_bc_type) :: fluxes !< A structure that may contain an array of
186  !! named fields used for passive tracer fluxes.
187  integer :: wind_stagger = -999 !< A flag indicating the spatial discretization of
188  !! wind stresses. This flag may be set by the
189  !! flux-exchange code, based on what the sea-ice
190  !! model is providing. Otherwise, the value from
191  !! the surface_forcing_CS is used.
193 
195 
196 contains
197 
198 !> This subroutine translates the Ice_ocean_boundary_type into a MOM
199 !! thermodynamic forcing type, including changes of units, sign conventions,
200 !! and putting the fields into arrays with MOM-standard halos.
201 subroutine convert_iob_to_fluxes(IOB, fluxes, index_bounds, Time, valid_time, G, US, CS, &
202  sfc_state, restore_salt, restore_temp)
204  target, intent(in) :: iob !< An ice-ocean boundary type with fluxes to drive
205  !! the ocean in a coupled model
206  type(forcing), intent(inout) :: fluxes !< A structure containing pointers to all
207  !! possible mass, heat or salt flux forcing fields.
208  !! Unused fields have NULL ptrs.
209  integer, dimension(4), intent(in) :: index_bounds !< The i- and j- size of the arrays in IOB.
210  type(time_type), intent(in) :: time !< The time of the fluxes, used for interpolating the
211  !! salinity to the right time, when it is being restored.
212  real, intent(in) :: valid_time !< The amount of time over which these fluxes
213  !! should be applied [s].
214  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
215  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
216  type(surface_forcing_cs),pointer :: cs !< A pointer to the control structure returned by a
217  !! previous call to surface_forcing_init.
218  type(surface), intent(in) :: sfc_state !< A structure containing fields that describe the
219  !! surface state of the ocean.
220  logical, optional, intent(in) :: restore_salt !< If true, salinity is restored to a target value.
221  logical, optional, intent(in) :: restore_temp !< If true, temperature is restored to a target value.
222 
223  ! local variables
224  real, dimension(SZI_(G),SZJ_(G)) :: &
225  data_restore, & !< The surface value toward which to restore [g/kg or degC]
226  sst_anom, & !< Instantaneous sea surface temperature anomalies from a target value [deg C]
227  sss_anom, & !< Instantaneous sea surface salinity anomalies from a target value [g/kg]
228  sss_mean, & !< A (mean?) salinity about which to normalize local salinity
229  !! anomalies when calculating restorative precipitation anomalies [g/kg]
230  pme_adj, & !< The adjustment to PminusE that will cause the salinity
231  !! to be restored toward its target value [kg/(m^2 * s)]
232  net_fw, & !< The area integrated net freshwater flux into the ocean [kg/s]
233  net_fw2, & !< The area integrated net freshwater flux into the ocean [kg/s]
234  work_sum, & !< A 2-d array that is used as the work space for a global
235  !! sum, used with units of m2 or [kg/s]
236  open_ocn_mask !< a binary field indicating where ice is present based on frazil criteria
237 
238  integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq, i0, j0
239  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, isr, ier, jsr, jer
240  integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
241 
242  logical :: restore_salinity !< local copy of the argument restore_salt, if it
243  !! is present, or false (no restoring) otherwise.
244  logical :: restore_sst !< local copy of the argument restore_temp, if it
245  !! is present, or false (no restoring) otherwise.
246  real :: delta_sss !< temporary storage for sss diff from restoring value
247  real :: delta_sst !< temporary storage for sst diff from restoring value
248  real :: kg_m2_s_conversion !< A combination of unit conversion factors for rescaling
249  !! mass fluxes [R Z s m2 kg-1 T-1 ~> 1].
250 
251  real :: c_p !< heat capacity of seawater ( J/(K kg) )
252  real :: sign_for_net_fw_bug !< Should be +1. but an old bug can be recovered by using -1.
253 
254  call cpu_clock_begin(id_clock_forcing)
255 
256  isc_bnd = index_bounds(1) ; iec_bnd = index_bounds(2)
257  jsc_bnd = index_bounds(3) ; jec_bnd = index_bounds(4)
258  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
259  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
260  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
261  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
262  isr = is-isd+1 ; ier = ie-isd+1 ; jsr = js-jsd+1 ; jer = je-jsd+1
263 
264  kg_m2_s_conversion = us%kg_m3_to_R*us%m_to_Z*us%T_to_s
265  c_p = fluxes%C_p
266  open_ocn_mask(:,:) = 1.0
267  pme_adj(:,:) = 0.0
268  fluxes%vPrecGlobalAdj = 0.0
269  fluxes%vPrecGlobalScl = 0.0
270  fluxes%saltFluxGlobalAdj = 0.0
271  fluxes%saltFluxGlobalScl = 0.0
272  fluxes%netFWGlobalAdj = 0.0
273  fluxes%netFWGlobalScl = 0.0
274 
275  restore_salinity = .false.
276  if (present(restore_salt)) restore_salinity = restore_salt
277  restore_sst = .false.
278  if (present(restore_temp)) restore_sst = restore_temp
279 
280  ! allocation and initialization if this is the first time that this
281  ! flux type has been used.
282  if (fluxes%dt_buoy_accum < 0) then
283  call allocate_forcing_type(g, fluxes, water=.true., heat=.true., &
284  ustar=.true., press=.true.)
285 
286  call safe_alloc_ptr(fluxes%sw_vis_dir,isd,ied,jsd,jed)
287  call safe_alloc_ptr(fluxes%sw_vis_dif,isd,ied,jsd,jed)
288  call safe_alloc_ptr(fluxes%sw_nir_dir,isd,ied,jsd,jed)
289  call safe_alloc_ptr(fluxes%sw_nir_dif,isd,ied,jsd,jed)
290 
291  call safe_alloc_ptr(fluxes%p_surf ,isd,ied,jsd,jed)
292  call safe_alloc_ptr(fluxes%p_surf_full,isd,ied,jsd,jed)
293  if (cs%use_limited_P_SSH) then
294  fluxes%p_surf_SSH => fluxes%p_surf
295  else
296  fluxes%p_surf_SSH => fluxes%p_surf_full
297  endif
298 
299  call safe_alloc_ptr(fluxes%salt_flux,isd,ied,jsd,jed)
300  call safe_alloc_ptr(fluxes%salt_flux_in,isd,ied,jsd,jed)
301  call safe_alloc_ptr(fluxes%salt_flux_added,isd,ied,jsd,jed)
302 
303  call safe_alloc_ptr(fluxes%TKE_tidal,isd,ied,jsd,jed)
304  call safe_alloc_ptr(fluxes%ustar_tidal,isd,ied,jsd,jed)
305 
306  if (cs%allow_flux_adjustments) then
307  call safe_alloc_ptr(fluxes%heat_added,isd,ied,jsd,jed)
308  call safe_alloc_ptr(fluxes%salt_flux_added,isd,ied,jsd,jed)
309  endif
310 
311  do j=js-2,je+2 ; do i=is-2,ie+2
312  fluxes%TKE_tidal(i,j) = cs%TKE_tidal(i,j)
313  fluxes%ustar_tidal(i,j) = cs%ustar_tidal(i,j)
314  enddo ; enddo
315 
316  if (restore_temp) call safe_alloc_ptr(fluxes%heat_added,isd,ied,jsd,jed)
317 
318  endif ! endif for allocation and initialization
319 
320  if (((associated(iob%ustar_berg) .and. (.not.associated(fluxes%ustar_berg))) &
321  .or. (associated(iob%area_berg) .and. (.not.associated(fluxes%area_berg)))) &
322  .or. (associated(iob%mass_berg) .and. (.not.associated(fluxes%mass_berg)))) &
323  call allocate_forcing_type(g, fluxes, iceberg=.true.)
324 
325  if ((.not.coupler_type_initialized(fluxes%tr_fluxes)) .and. &
326  coupler_type_initialized(iob%fluxes)) &
327  call coupler_type_spawn(iob%fluxes, fluxes%tr_fluxes, &
328  (/is,is,ie,ie/), (/js,js,je,je/))
329  ! It might prove valuable to use the same array extents as the rest of the
330  ! ocean model, rather than using haloless arrays, in which case the last line
331  ! would be: ( (/isd,is,ie,ied/), (/jsd,js,je,jed/))
332 
333  ! allocation and initialization on first call to this routine
334  if (cs%area_surf < 0.0) then
335  do j=js,je ; do i=is,ie
336  work_sum(i,j) = us%L_to_m**2*g%areaT(i,j) * g%mask2dT(i,j)
337  enddo ; enddo
338  cs%area_surf = reproducing_sum(work_sum, isr, ier, jsr, jer)
339  endif ! endif for allocation and initialization
340 
341 
342  ! Indicate that there are new unused fluxes.
343  fluxes%fluxes_used = .false.
344  fluxes%dt_buoy_accum = us%s_to_T*valid_time
345 
346  if (cs%allow_flux_adjustments) then
347  fluxes%heat_added(:,:)=0.0
348  fluxes%salt_flux_added(:,:)=0.0
349  endif
350 
351  do j=js,je ; do i=is,ie
352  fluxes%salt_flux(i,j) = 0.0
353  fluxes%vprec(i,j) = 0.0
354  enddo ; enddo
355 
356  ! Salinity restoring logic
357  if (restore_salinity) then
358  call time_interp_external(cs%id_srestore,time,data_restore)
359  ! open_ocn_mask indicates where to restore salinity (1 means restore, 0 does not)
360  open_ocn_mask(:,:) = 1.0
361  if (cs%mask_srestore_under_ice) then ! Do not restore under sea-ice
362  do j=js,je ; do i=is,ie
363  if (sfc_state%SST(i,j) <= -0.0539*sfc_state%SSS(i,j)) open_ocn_mask(i,j)=0.0
364  enddo ; enddo
365  endif
366  if (cs%salt_restore_as_sflux) then
367  do j=js,je ; do i=is,ie
368  delta_sss = data_restore(i,j)- sfc_state%SSS(i,j)
369  delta_sss = sign(1.0,delta_sss)*min(abs(delta_sss),cs%max_delta_srestore)
370  fluxes%salt_flux(i,j) = 1.e-3*g%mask2dT(i,j) * (cs%Rho0*us%m_to_Z*us%T_to_s*cs%Flux_const)* &
371  (cs%basin_mask(i,j)*open_ocn_mask(i,j)*cs%srestore_mask(i,j)) *delta_sss ! kg Salt m-2 s-1
372  enddo ; enddo
373  if (cs%adjust_net_srestore_to_zero) then
374  if (cs%adjust_net_srestore_by_scaling) then
375  call adjust_area_mean_to_zero(fluxes%salt_flux, g, fluxes%saltFluxGlobalScl, &
376  unit_scale=us%R_to_kg_m3*us%Z_to_m*us%s_to_T)
377  fluxes%saltFluxGlobalAdj = 0.
378  else
379  work_sum(is:ie,js:je) = us%L_to_m**2*us%R_to_kg_m3*us%Z_to_m*us%s_to_T * &
380  g%areaT(is:ie,js:je)*fluxes%salt_flux(is:ie,js:je)
381  fluxes%saltFluxGlobalAdj = reproducing_sum(work_sum(:,:), isr,ier, jsr,jer)/cs%area_surf
382  fluxes%salt_flux(is:ie,js:je) = fluxes%salt_flux(is:ie,js:je) - kg_m2_s_conversion * fluxes%saltFluxGlobalAdj
383  endif
384  endif
385  fluxes%salt_flux_added(is:ie,js:je) = fluxes%salt_flux(is:ie,js:je) ! Diagnostic
386  else
387  do j=js,je ; do i=is,ie
388  if (g%mask2dT(i,j) > 0.5) then
389  delta_sss = sfc_state%SSS(i,j) - data_restore(i,j)
390  delta_sss = sign(1.0,delta_sss)*min(abs(delta_sss),cs%max_delta_srestore)
391  fluxes%vprec(i,j) = (cs%basin_mask(i,j)*open_ocn_mask(i,j)*cs%srestore_mask(i,j))* &
392  (us%m_to_Z*us%T_to_s * cs%Rho0*cs%Flux_const) * &
393  delta_sss / (0.5*(sfc_state%SSS(i,j) + data_restore(i,j)))
394  endif
395  enddo ; enddo
396  if (cs%adjust_net_srestore_to_zero) then
397  if (cs%adjust_net_srestore_by_scaling) then
398  call adjust_area_mean_to_zero(fluxes%vprec, g, fluxes%vPrecGlobalScl, &
399  unit_scale=us%R_to_kg_m3*us%Z_to_m*us%s_to_T)
400  fluxes%vPrecGlobalAdj = 0.
401  else
402  work_sum(is:ie,js:je) = us%L_to_m**2*g%areaT(is:ie,js:je) * &
403  us%R_to_kg_m3*us%Z_to_m*us%s_to_T*fluxes%vprec(is:ie,js:je)
404  fluxes%vPrecGlobalAdj = reproducing_sum(work_sum(:,:), isr, ier, jsr, jer) / cs%area_surf
405  do j=js,je ; do i=is,ie
406  fluxes%vprec(i,j) = ( fluxes%vprec(i,j) - kg_m2_s_conversion * fluxes%vPrecGlobalAdj ) * g%mask2dT(i,j)
407  enddo ; enddo
408  endif
409  endif
410  endif
411  endif
412 
413  ! SST restoring logic
414  if (restore_sst) then
415  call time_interp_external(cs%id_trestore,time,data_restore)
416  do j=js,je ; do i=is,ie
417  delta_sst = data_restore(i,j)- sfc_state%SST(i,j)
418  delta_sst = sign(1.0,delta_sst)*min(abs(delta_sst),cs%max_delta_trestore)
419  fluxes%heat_added(i,j) = g%mask2dT(i,j) * cs%trestore_mask(i,j) * &
420  (us%R_to_kg_m3*cs%Rho0*fluxes%C_p) * delta_sst * cs%Flux_const ! W m-2
421  enddo ; enddo
422  endif
423 
424  ! Check that liquid runoff has a place to go
425  if (cs%liquid_runoff_from_data .and. .not. associated(iob%lrunoff)) then
426  call mom_error(fatal, "liquid runoff is being added via data_override but "// &
427  "there is no associated runoff in the IOB")
428  return
429  end if
430 
431  ! obtain fluxes from IOB; note the staggering of indices
432  i0 = is - isc_bnd ; j0 = js - jsc_bnd
433  do j=js,je ; do i=is,ie
434 
435  if (associated(iob%lprec)) &
436  fluxes%lprec(i,j) = kg_m2_s_conversion * iob%lprec(i-i0,j-j0) * g%mask2dT(i,j)
437 
438  if (associated(iob%fprec)) &
439  fluxes%fprec(i,j) = kg_m2_s_conversion * iob%fprec(i-i0,j-j0) * g%mask2dT(i,j)
440 
441  if (associated(iob%q_flux)) &
442  fluxes%evap(i,j) = kg_m2_s_conversion * iob%q_flux(i-i0,j-j0) * g%mask2dT(i,j)
443 
444  ! liquid runoff flux
445  if (associated(iob%lrunoff)) then
446  if(cs%liquid_runoff_from_data)call data_override('OCN', 'runoff', iob%lrunoff, time)
447  fluxes%lrunoff(i,j) = kg_m2_s_conversion * iob%lrunoff(i-i0,j-j0) * g%mask2dT(i,j)
448  endif
449 
450  ! ice runoff flux
451  if (associated(iob%frunoff)) then
452  fluxes%frunoff(i,j) = kg_m2_s_conversion * iob%frunoff(i-i0,j-j0) * g%mask2dT(i,j)
453  endif
454 
455  if (associated(iob%ustar_berg)) &
456  fluxes%ustar_berg(i,j) = us%m_to_Z*us%T_to_s * iob%ustar_berg(i-i0,j-j0) * g%mask2dT(i,j)
457 
458  if (associated(iob%area_berg)) &
459  fluxes%area_berg(i,j) = iob%area_berg(i-i0,j-j0) * g%mask2dT(i,j)
460 
461  if (associated(iob%mass_berg)) &
462  fluxes%mass_berg(i,j) = iob%mass_berg(i-i0,j-j0) * g%mask2dT(i,j)
463 
464  if (associated(iob%lrunoff_hflx)) &
465  fluxes%heat_content_lrunoff(i,j) = iob%lrunoff_hflx(i-i0,j-j0) * g%mask2dT(i,j)
466 
467  if (associated(iob%frunoff_hflx)) &
468  fluxes%heat_content_frunoff(i,j) = kg_m2_s_conversion * iob%frunoff_hflx(i-i0,j-j0) * g%mask2dT(i,j)
469 
470  if (associated(iob%lw_flux)) &
471  fluxes%LW(i,j) = iob%lw_flux(i-i0,j-j0) * g%mask2dT(i,j)
472 
473  if (associated(iob%t_flux)) &
474  fluxes%sens(i,j) = iob%t_flux(i-i0,j-j0) * g%mask2dT(i,j)
475 
476  ! sea ice and snow melt heat flux [W/m2]
477  if (associated(iob%seaice_melt_heat)) &
478  fluxes%seaice_melt_heat(i,j) = g%mask2dT(i,j) * iob%seaice_melt_heat(i-i0,j-j0)
479 
480  ! water flux due to sea ice and snow melt [kg/m2/s]
481  if (associated(iob%seaice_melt)) &
482  fluxes%seaice_melt(i,j) = kg_m2_s_conversion * g%mask2dT(i,j) * iob%seaice_melt(i-i0,j-j0)
483 
484  fluxes%latent(i,j) = 0.0
485  if (associated(iob%fprec)) then
486  fluxes%latent(i,j) = fluxes%latent(i,j) + iob%fprec(i-i0,j-j0)*cs%latent_heat_fusion
487  fluxes%latent_fprec_diag(i,j) = g%mask2dT(i,j) * iob%fprec(i-i0,j-j0)*cs%latent_heat_fusion
488  endif
489  if (associated(iob%frunoff)) then
490  fluxes%latent(i,j) = fluxes%latent(i,j) + iob%frunoff(i-i0,j-j0)*cs%latent_heat_fusion
491  fluxes%latent_frunoff_diag(i,j) = g%mask2dT(i,j) * iob%frunoff(i-i0,j-j0)*cs%latent_heat_fusion
492  endif
493  if (associated(iob%q_flux)) then
494  fluxes%latent(i,j) = fluxes%latent(i,j) + iob%q_flux(i-i0,j-j0)*cs%latent_heat_vapor
495  fluxes%latent_evap_diag(i,j) = g%mask2dT(i,j) * iob%q_flux(i-i0,j-j0)*cs%latent_heat_vapor
496  endif
497  fluxes%latent(i,j) = g%mask2dT(i,j) * fluxes%latent(i,j)
498 
499  if (associated(iob%sw_flux_vis_dir)) &
500  fluxes%sw_vis_dir(i,j) = g%mask2dT(i,j) * iob%sw_flux_vis_dir(i-i0,j-j0)
501 
502  if (associated(iob%sw_flux_vis_dif)) &
503  fluxes%sw_vis_dif(i,j) = g%mask2dT(i,j) * iob%sw_flux_vis_dif(i-i0,j-j0)
504 
505  if (associated(iob%sw_flux_nir_dir)) &
506  fluxes%sw_nir_dir(i,j) = g%mask2dT(i,j) * iob%sw_flux_nir_dir(i-i0,j-j0)
507 
508  if (associated(iob%sw_flux_nir_dif)) &
509  fluxes%sw_nir_dif(i,j) = g%mask2dT(i,j) * iob%sw_flux_nir_dif(i-i0,j-j0)
510 
511  fluxes%sw(i,j) = fluxes%sw_vis_dir(i,j) + fluxes%sw_vis_dif(i,j) + &
512  fluxes%sw_nir_dir(i,j) + fluxes%sw_nir_dif(i,j)
513 
514  enddo ; enddo
515 
516  ! applied surface pressure from atmosphere and cryosphere
517  if (associated(iob%p)) then
518  if (cs%max_p_surf >= 0.0) then
519  do j=js,je ; do i=is,ie
520  fluxes%p_surf_full(i,j) = g%mask2dT(i,j) * iob%p(i-i0,j-j0)
521  fluxes%p_surf(i,j) = min(fluxes%p_surf_full(i,j),cs%max_p_surf)
522  enddo; enddo
523  else
524  do j=js,je ; do i=is,ie
525  fluxes%p_surf_full(i,j) = g%mask2dT(i,j) * iob%p(i-i0,j-j0)
526  fluxes%p_surf(i,j) = fluxes%p_surf_full(i,j)
527  enddo; enddo
528  endif
529  fluxes%accumulate_p_surf = .true. ! Multiple components may contribute to surface pressure.
530  endif
531 
532  if (associated(iob%salt_flux)) then
533  do j=js,je ; do i=is,ie
534  fluxes%salt_flux(i,j) = g%mask2dT(i,j)*(fluxes%salt_flux(i,j) + kg_m2_s_conversion*iob%salt_flux(i-i0,j-j0))
535  fluxes%salt_flux_in(i,j) = g%mask2dT(i,j)*( kg_m2_s_conversion*iob%salt_flux(i-i0,j-j0) )
536  enddo ; enddo
537  endif
538 
539  ! adjust the NET fresh-water flux to zero, if flagged
540  if (cs%adjust_net_fresh_water_to_zero) then
541  sign_for_net_fw_bug = 1.
542  if (cs%use_net_FW_adjustment_sign_bug) sign_for_net_fw_bug = -1.
543  do j=js,je ; do i=is,ie
544  net_fw(i,j) = us%R_to_kg_m3*us%Z_to_m*us%s_to_T * &
545  (((fluxes%lprec(i,j) + fluxes%fprec(i,j) + fluxes%seaice_melt(i,j)) + &
546  (fluxes%lrunoff(i,j) + fluxes%frunoff(i,j))) + &
547  (fluxes%evap(i,j) + fluxes%vprec(i,j)) ) * us%L_to_m**2*g%areaT(i,j)
548  net_fw2(i,j) = net_fw(i,j) / (us%L_to_m**2*g%areaT(i,j))
549  enddo ; enddo
550 
551  if (cs%adjust_net_fresh_water_by_scaling) then
552  call adjust_area_mean_to_zero(net_fw2, g, fluxes%netFWGlobalScl)
553  do j=js,je ; do i=is,ie
554  fluxes%vprec(i,j) = fluxes%vprec(i,j) + us%kg_m3_to_R*us%m_to_Z*us%T_to_s * &
555  (net_fw2(i,j) - net_fw(i,j)/(us%L_to_m**2*g%areaT(i,j))) * g%mask2dT(i,j)
556  enddo ; enddo
557  else
558  fluxes%netFWGlobalAdj = reproducing_sum(net_fw(:,:), isr, ier, jsr, jer) / cs%area_surf
559  do j=js,je ; do i=is,ie
560  fluxes%vprec(i,j) = ( fluxes%vprec(i,j) - kg_m2_s_conversion * fluxes%netFWGlobalAdj ) * g%mask2dT(i,j)
561  enddo ; enddo
562  endif
563 
564  endif
565 
566  if (coupler_type_initialized(fluxes%tr_fluxes) .and. &
567  coupler_type_initialized(iob%fluxes)) &
568  call coupler_type_copy_data(iob%fluxes, fluxes%tr_fluxes)
569 
570  if (cs%allow_flux_adjustments) then
571  ! Apply adjustments to fluxes
572  call apply_flux_adjustments(g, us, cs, time, fluxes)
573  endif
574 
575  ! Allow for user-written code to alter fluxes after all the above
576  call user_alter_forcing(sfc_state, fluxes, time, g, cs%urf_CS)
577 
578  call cpu_clock_end(id_clock_forcing)
579 
580 end subroutine convert_iob_to_fluxes
581 
582 !> This subroutine translates the Ice_ocean_boundary_type into a MOM
583 !! mechanical forcing type, including changes of units, sign conventions,
584 !! and putting the fields into arrays with MOM-standard halos.
585 subroutine convert_iob_to_forces(IOB, forces, index_bounds, Time, G, US, CS)
587  target, intent(in) :: iob !< An ice-ocean boundary type with fluxes to drive
588  !! the ocean in a coupled model
589  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
590  integer, dimension(4), intent(in) :: index_bounds !< The i- and j- size of the arrays in IOB.
591  type(time_type), intent(in) :: time !< The time of the fluxes, used for interpolating the
592  !! salinity to the right time, when it is being restored.
593  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
594  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
595  type(surface_forcing_cs),pointer :: cs !< A pointer to the control structure returned by a
596  !! previous call to surface_forcing_init.
597 
598  ! local variables
599  real, dimension(SZIB_(G),SZJB_(G)) :: &
600  taux_at_q, & !< Zonal wind stresses at q points [Pa]
601  tauy_at_q !< Meridional wind stresses at q points [Pa]
602 
603  real, dimension(SZI_(G),SZJ_(G)) :: &
604  rigidity_at_h, & !< Ice rigidity at tracer points (m3 s-1)
605  taux_at_h, & !< Zonal wind stresses at h points [Pa]
606  tauy_at_h !< Meridional wind stresses at h points [Pa]
607 
608  real :: gustiness !< unresolved gustiness that contributes to ustar [R Z L T-2 ~> Pa]
609  real :: irho0 !< inverse of the mean density in [Z L-1 R-1 ~> m3 kg-1]
610  real :: taux2, tauy2 !< squared wind stresses [R2 Z2 L2 T-4 ~> Pa2]
611  real :: tau_mag !< magnitude of the wind stress [R Z L T-2 ~> Pa]
612  real :: pa_conversion ! A unit conversion factor from Pa to the internal wind stress units [R Z L T-2 Pa-1 ~> 1]
613  real :: stress_conversion ! A unit conversion factor from Pa times any stress multiplier [R Z L T-2 Pa-1 ~> 1]
614  real :: i_gearth !< 1.0 / G_Earth [s2 m-1]
615  real :: kv_rho_ice !< (CS%kv_sea_ice / CS%density_sea_ice) ( m^5/(s*kg) )
616  real :: mass_ice !< mass of sea ice at a face (kg/m^2)
617  real :: mass_eff !< effective mass of sea ice for rigidity (kg/m^2)
618 
619  integer :: wind_stagger !< AGRID, BGRID_NE, or CGRID_NE (integers from MOM_domains)
620  integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq, i0, j0
621  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, isr, ier, jsr, jer
622  integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
623 
624  call cpu_clock_begin(id_clock_forcing)
625 
626  isc_bnd = index_bounds(1) ; iec_bnd = index_bounds(2)
627  jsc_bnd = index_bounds(3) ; jec_bnd = index_bounds(4)
628  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
629  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
630  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
631  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
632  isr = is-isd+1 ; ier = ie-isd+1 ; jsr = js-jsd+1 ; jer = je-jsd+1
633  i0 = is - isc_bnd ; j0 = js - jsc_bnd
634 
635  irho0 = us%L_to_Z / cs%Rho0
636  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
637  stress_conversion = pa_conversion * cs%wind_stress_multiplier
638 
639  ! allocation and initialization if this is the first time that this
640  ! mechanical forcing type has been used.
641  if (.not.forces%initialized) then
642 
643  call allocate_mech_forcing(g, forces, stress=.true., ustar=.true., press=.true.)
644 
645  call safe_alloc_ptr(forces%p_surf,isd,ied,jsd,jed)
646  call safe_alloc_ptr(forces%p_surf_full,isd,ied,jsd,jed)
647 
648  if (cs%rigid_sea_ice) then
649  call safe_alloc_ptr(forces%rigidity_ice_u,isdb,iedb,jsd,jed)
650  call safe_alloc_ptr(forces%rigidity_ice_v,isd,ied,jsdb,jedb)
651  endif
652 
653  forces%initialized = .true.
654  endif
655 
656  if ( (associated(iob%area_berg) .and. (.not. associated(forces%area_berg))) .or. &
657  (associated(iob%mass_berg) .and. (.not. associated(forces%mass_berg))) ) &
658  call allocate_mech_forcing(g, forces, iceberg=.true.)
659  if (associated(iob%ice_rigidity)) then
660  rigidity_at_h(:,:) = 0.0
661  call safe_alloc_ptr(forces%rigidity_ice_u,isdb,iedb,jsd,jed)
662  call safe_alloc_ptr(forces%rigidity_ice_v,isd,ied,jsdb,jedb)
663  endif
664 
665  forces%accumulate_rigidity = .true. ! Multiple components may contribute to rigidity.
666  if (associated(forces%rigidity_ice_u)) forces%rigidity_ice_u(:,:) = 0.0
667  if (associated(forces%rigidity_ice_v)) forces%rigidity_ice_v(:,:) = 0.0
668 
669  ! applied surface pressure from atmosphere and cryosphere
670  if (cs%use_limited_P_SSH) then
671  forces%p_surf_SSH => forces%p_surf
672  else
673  forces%p_surf_SSH => forces%p_surf_full
674  endif
675  if (associated(iob%p)) then
676  if (cs%max_p_surf >= 0.0) then
677  do j=js,je ; do i=is,ie
678  forces%p_surf_full(i,j) = g%mask2dT(i,j) * iob%p(i-i0,j-j0)
679  forces%p_surf(i,j) = min(forces%p_surf_full(i,j),cs%max_p_surf)
680  enddo ; enddo
681  else
682  do j=js,je ; do i=is,ie
683  forces%p_surf_full(i,j) = g%mask2dT(i,j) * iob%p(i-i0,j-j0)
684  forces%p_surf(i,j) = forces%p_surf_full(i,j)
685  enddo ; enddo
686  endif
687  else
688  do j=js,je ; do i=is,ie
689  forces%p_surf_full(i,j) = 0.0
690  forces%p_surf(i,j) = 0.0
691  enddo ; enddo
692  endif
693  forces%accumulate_p_surf = .true. ! Multiple components may contribute to surface pressure.
694 
695  ! TODO: this does not seem correct for NEMS
696 #ifdef CESMCOUPLED
697  wind_stagger = agrid
698 #else
699  wind_stagger = cs%wind_stagger
700 #endif
701 
702  if (wind_stagger == bgrid_ne) then
703  ! This is necessary to fill in the halo points.
704  taux_at_q(:,:) = 0.0 ; tauy_at_q(:,:) = 0.0
705  endif
706  if (wind_stagger == agrid) then
707  ! This is necessary to fill in the halo points.
708  taux_at_h(:,:) = 0.0 ; tauy_at_h(:,:) = 0.0
709  endif
710 
711  ! obtain fluxes from IOB; note the staggering of indices
712  do j=js,je ; do i=is,ie
713  if (associated(iob%area_berg)) &
714  forces%area_berg(i,j) = iob%area_berg(i-i0,j-j0) * g%mask2dT(i,j)
715 
716  if (associated(iob%mass_berg)) &
717  forces%mass_berg(i,j) = iob%mass_berg(i-i0,j-j0) * g%mask2dT(i,j)
718 
719  if (associated(iob%ice_rigidity)) &
720  rigidity_at_h(i,j) = iob%ice_rigidity(i-i0,j-j0) * g%mask2dT(i,j)
721 
722  if (wind_stagger == bgrid_ne) then
723  if (associated(iob%u_flux)) taux_at_q(i,j) = iob%u_flux(i-i0,j-j0) * stress_conversion
724  if (associated(iob%v_flux)) tauy_at_q(i,j) = iob%v_flux(i-i0,j-j0) * stress_conversion
725  elseif (wind_stagger == agrid) then
726  if (associated(iob%u_flux)) taux_at_h(i,j) = iob%u_flux(i-i0,j-j0) * stress_conversion
727  if (associated(iob%v_flux)) tauy_at_h(i,j) = iob%v_flux(i-i0,j-j0) * stress_conversion
728  else ! C-grid wind stresses.
729  if (associated(iob%u_flux)) forces%taux(i,j) = iob%u_flux(i-i0,j-j0) * stress_conversion
730  if (associated(iob%v_flux)) forces%tauy(i,j) = iob%v_flux(i-i0,j-j0) * stress_conversion
731  endif
732 
733  enddo ; enddo
734 
735  ! surface momentum stress related fields as function of staggering
736  if (wind_stagger == bgrid_ne) then
737  if (g%symmetric) &
738  call fill_symmetric_edges(taux_at_q, tauy_at_q, g%Domain, stagger=bgrid_ne)
739  call pass_vector(taux_at_q, tauy_at_q, g%Domain, stagger=bgrid_ne, halo=1)
740 
741  do j=js,je ; do i=isq,ieq
742  forces%taux(i,j) = 0.0
743  if ((g%mask2dBu(i,j) + g%mask2dBu(i,j-1)) > 0) &
744  forces%taux(i,j) = (g%mask2dBu(i,j)*taux_at_q(i,j) + &
745  g%mask2dBu(i,j-1)*taux_at_q(i,j-1)) / &
746  (g%mask2dBu(i,j) + g%mask2dBu(i,j-1))
747  enddo ; enddo
748 
749  do j=jsq,jeq ; do i=is,ie
750  forces%tauy(i,j) = 0.0
751  if ((g%mask2dBu(i,j) + g%mask2dBu(i-1,j)) > 0) &
752  forces%tauy(i,j) = (g%mask2dBu(i,j)*tauy_at_q(i,j) + &
753  g%mask2dBu(i-1,j)*tauy_at_q(i-1,j)) / &
754  (g%mask2dBu(i,j) + g%mask2dBu(i-1,j))
755  enddo ; enddo
756 
757  ! ustar is required for the bulk mixed layer formulation. The background value
758  ! of 0.02 Pa is a relatively small value intended to give reasonable behavior
759  ! in regions of very weak winds.
760 
761  do j=js,je ; do i=is,ie
762  tau_mag = 0.0 ; gustiness = cs%gust_const
763  if (((g%mask2dBu(i,j) + g%mask2dBu(i-1,j-1)) + &
764  (g%mask2dBu(i,j-1) + g%mask2dBu(i-1,j))) > 0) then
765  tau_mag = sqrt(((g%mask2dBu(i,j)*(taux_at_q(i,j)**2 + tauy_at_q(i,j)**2) + &
766  g%mask2dBu(i-1,j-1)*(taux_at_q(i-1,j-1)**2 + tauy_at_q(i-1,j-1)**2)) + &
767  (g%mask2dBu(i,j-1)*(taux_at_q(i,j-1)**2 + tauy_at_q(i,j-1)**2) + &
768  g%mask2dBu(i-1,j)*(taux_at_q(i-1,j)**2 + tauy_at_q(i-1,j)**2)) ) / &
769  ((g%mask2dBu(i,j) + g%mask2dBu(i-1,j-1)) + (g%mask2dBu(i,j-1) + g%mask2dBu(i-1,j))) )
770  if (cs%read_gust_2d) gustiness = cs%gust(i,j)
771  endif
772  forces%ustar(i,j) = sqrt(gustiness*irho0 + irho0*tau_mag)
773  enddo ; enddo
774 
775  elseif (wind_stagger == agrid) then
776  call pass_vector(taux_at_h, tauy_at_h, g%Domain, to_all+omit_corners, stagger=agrid, halo=1)
777 
778  do j=js,je ; do i=isq,ieq
779  forces%taux(i,j) = 0.0
780  if ((g%mask2dT(i,j) + g%mask2dT(i+1,j)) > 0) &
781  forces%taux(i,j) = (g%mask2dT(i,j)*taux_at_h(i,j) + &
782  g%mask2dT(i+1,j)*taux_at_h(i+1,j)) / &
783  (g%mask2dT(i,j) + g%mask2dT(i+1,j))
784  enddo ; enddo
785 
786  do j=jsq,jeq ; do i=is,ie
787  forces%tauy(i,j) = 0.0
788  if ((g%mask2dT(i,j) + g%mask2dT(i,j+1)) > 0) &
789  forces%tauy(i,j) = (g%mask2dT(i,j)*tauy_at_h(i,j) + &
790  g%mask2dT(i,j+1)*tauy_at_h(i,j+1)) / &
791  (g%mask2dT(i,j) + g%mask2dT(i,j+1))
792  enddo ; enddo
793 
794  do j=js,je ; do i=is,ie
795  gustiness = cs%gust_const
796  if (cs%read_gust_2d .and. (g%mask2dT(i,j) > 0)) gustiness = cs%gust(i,j)
797  forces%ustar(i,j) = sqrt(gustiness*irho0 + irho0 * g%mask2dT(i,j) * &
798  sqrt(taux_at_h(i,j)**2 + tauy_at_h(i,j)**2))
799  enddo ; enddo
800 
801  else ! C-grid wind stresses.
802  if (g%symmetric) &
803  call fill_symmetric_edges(forces%taux, forces%tauy, g%Domain)
804  call pass_vector(forces%taux, forces%tauy, g%Domain, halo=1)
805 
806  do j=js,je ; do i=is,ie
807  taux2 = 0.0
808  if ((g%mask2dCu(i-1,j) + g%mask2dCu(i,j)) > 0) &
809  taux2 = (g%mask2dCu(i-1,j)*forces%taux(i-1,j)**2 + &
810  g%mask2dCu(i,j)*forces%taux(i,j)**2) / (g%mask2dCu(i-1,j) + g%mask2dCu(i,j))
811 
812  tauy2 = 0.0
813  if ((g%mask2dCv(i,j-1) + g%mask2dCv(i,j)) > 0) &
814  tauy2 = (g%mask2dCv(i,j-1)*forces%tauy(i,j-1)**2 + &
815  g%mask2dCv(i,j)*forces%tauy(i,j)**2) / (g%mask2dCv(i,j-1) + g%mask2dCv(i,j))
816 
817  if (cs%read_gust_2d) then
818  forces%ustar(i,j) = sqrt(cs%gust(i,j)*irho0 + irho0*sqrt(taux2 + tauy2))
819  else
820  forces%ustar(i,j) = sqrt(cs%gust_const*irho0 + irho0*sqrt(taux2 + tauy2))
821  endif
822  enddo ; enddo
823 
824  endif ! endif for wind related fields
825 
826  ! sea ice related dynamic fields
827  if (associated(iob%ice_rigidity)) then
828  call pass_var(rigidity_at_h, g%Domain, halo=1)
829  do i=is-1,ie ; do j=js,je
830  forces%rigidity_ice_u(i,j) = forces%rigidity_ice_u(i,j) + &
831  min(rigidity_at_h(i,j), rigidity_at_h(i+1,j))
832  enddo ; enddo
833  do i=is,ie ; do j=js-1,je
834  forces%rigidity_ice_v(i,j) = forces%rigidity_ice_v(i,j) + &
835  min(rigidity_at_h(i,j), rigidity_at_h(i,j+1))
836  enddo ; enddo
837  endif
838 
839  if (cs%rigid_sea_ice) then
840  call pass_var(forces%p_surf_full, g%Domain, halo=1)
841  i_gearth = 1.0 / cs%g_Earth
842  kv_rho_ice = (cs%kv_sea_ice / cs%density_sea_ice)
843  do i=is-1,ie ; do j=js,je
844  mass_ice = min(forces%p_surf_full(i,j), forces%p_surf_full(i+1,j)) * i_gearth
845  mass_eff = 0.0
846  if (mass_ice > cs%rigid_sea_ice_mass) then
847  mass_eff = (mass_ice - cs%rigid_sea_ice_mass) **2 / &
848  (mass_ice + cs%rigid_sea_ice_mass)
849  endif
850  forces%rigidity_ice_u(i,j) = forces%rigidity_ice_u(i,j) + kv_rho_ice * mass_eff
851  enddo ; enddo
852  do i=is,ie ; do j=js-1,je
853  mass_ice = min(forces%p_surf_full(i,j), forces%p_surf_full(i,j+1)) * i_gearth
854  mass_eff = 0.0
855  if (mass_ice > cs%rigid_sea_ice_mass) then
856  mass_eff = (mass_ice - cs%rigid_sea_ice_mass) **2 / &
857  (mass_ice + cs%rigid_sea_ice_mass)
858  endif
859  forces%rigidity_ice_v(i,j) = forces%rigidity_ice_v(i,j) + kv_rho_ice * mass_eff
860  enddo ; enddo
861  endif
862 
863  if (cs%allow_flux_adjustments) then
864  ! Apply adjustments to forces
865  call apply_force_adjustments(g, us, cs, time, forces)
866  endif
867 
868 !### ! Allow for user-written code to alter fluxes after all the above
869 !### call user_alter_mech_forcing(forces, Time, G, CS%urf_CS)
870 
871  call cpu_clock_end(id_clock_forcing)
872 end subroutine convert_iob_to_forces
873 
874 !> Adds thermodynamic flux adjustments obtained via data_override
875 !! Component name is 'OCN'
876 !! Available adjustments are:
877 !! - hflx_adj (Heat flux into the ocean, in W m-2)
878 !! - sflx_adj (Salt flux into the ocean, in kg salt m-2 s-1)
879 !! - prcme_adj (Fresh water flux into the ocean, in kg m-2 s-1)
880 subroutine apply_flux_adjustments(G, US, CS, Time, fluxes)
881  type(ocean_grid_type), intent(inout) :: g !< Ocean grid structure
882  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
883  type(surface_forcing_cs), pointer :: cs !< Surface forcing control structure
884  type(time_type), intent(in) :: time !< Model time structure
885  type(forcing), intent(inout) :: fluxes !< Surface fluxes structure
886 
887  ! Local variables
888  real, dimension(SZI_(G),SZJ_(G)) :: temp_at_h !< Fluxes at h points [W m-2 or kg m-2 s-1]
889 
890  integer :: isc, iec, jsc, jec, i, j
891  logical :: overrode_h
892 
893  isc = g%isc; iec = g%iec ; jsc = g%jsc; jec = g%jec
894 
895  overrode_h = .false.
896  call data_override('OCN', 'hflx_adj', temp_at_h(isc:iec,jsc:jec), time, override=overrode_h)
897 
898  if (overrode_h) then ; do j=jsc,jec ; do i=isc,iec
899  fluxes%heat_added(i,j) = fluxes%heat_added(i,j) + temp_at_h(i,j)* g%mask2dT(i,j)
900  enddo ; enddo ; endif
901  ! Not needed? ! if (overrode_h) call pass_var(fluxes%heat_added, G%Domain)
902 
903  overrode_h = .false.
904  call data_override('OCN', 'sflx_adj', temp_at_h(isc:iec,jsc:jec), time, override=overrode_h)
905 
906  if (overrode_h) then ; do j=jsc,jec ; do i=isc,iec
907  fluxes%salt_flux_added(i,j) = fluxes%salt_flux_added(i,j) + &
908  us%kg_m3_to_R*us%m_to_Z*us%T_to_s * temp_at_h(i,j)* g%mask2dT(i,j)
909  enddo ; enddo ; endif
910  ! Not needed? ! if (overrode_h) call pass_var(fluxes%salt_flux_added, G%Domain)
911 
912  overrode_h = .false.
913  call data_override('OCN', 'prcme_adj', temp_at_h(isc:iec,jsc:jec), time, override=overrode_h)
914 
915  if (overrode_h) then ; do j=jsc,jec ; do i=isc,iec
916  fluxes%vprec(i,j) = fluxes%vprec(i,j) + us%kg_m3_to_R*us%m_to_Z*us%T_to_s * temp_at_h(i,j)* g%mask2dT(i,j)
917  enddo ; enddo ; endif
918  ! Not needed? ! if (overrode_h) call pass_var(fluxes%vprec, G%Domain)
919 end subroutine apply_flux_adjustments
920 
921 !> Adds mechanical forcing adjustments obtained via data_override
922 !! Component name is 'OCN'
923 !! Available adjustments are:
924 !! - taux_adj (Zonal wind stress delta, positive to the east, in Pa)
925 !! - tauy_adj (Meridional wind stress delta, positive to the north, in Pa)
926 subroutine apply_force_adjustments(G, US, CS, Time, forces)
927  type(ocean_grid_type), intent(inout) :: g !< Ocean grid structure
928  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
929  type(surface_forcing_cs), pointer :: cs !< Surface forcing control structure
930  type(time_type), intent(in) :: time !< Model time structure
931  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
932 
933  ! Local variables
934  real, dimension(SZI_(G),SZJ_(G)) :: tempx_at_h !< Delta to zonal wind stress at h points [R Z L T-2 ~> Pa]
935  real, dimension(SZI_(G),SZJ_(G)) :: tempy_at_h !< Delta to meridional wind stress at h points [R Z L T-2 ~> Pa]
936 
937  integer :: isc, iec, jsc, jec, i, j
938  real :: dlondx, dlondy, rdlon, cosa, sina, zonal_tau, merid_tau
939  real :: pa_conversion ! A unit conversion factor from Pa to the internal units [R Z L T-2 Pa-1 ~> 1]
940  logical :: overrode_x, overrode_y
941 
942  isc = g%isc; iec = g%iec ; jsc = g%jsc; jec = g%jec
943  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
944 
945  tempx_at_h(:,:) = 0.0 ; tempy_at_h(:,:) = 0.0
946  ! Either reads data or leaves contents unchanged
947  overrode_x = .false. ; overrode_y = .false.
948  call data_override('OCN', 'taux_adj', tempx_at_h(isc:iec,jsc:jec), time, override=overrode_x)
949  call data_override('OCN', 'tauy_adj', tempy_at_h(isc:iec,jsc:jec), time, override=overrode_y)
950 
951  if (overrode_x .or. overrode_y) then
952  if (.not. (overrode_x .and. overrode_y)) call mom_error(fatal,"apply_flux_adjustments: "//&
953  "Both taux_adj and tauy_adj must be specified, or neither, in data_table")
954 
955  ! Rotate winds
956  call pass_vector(tempx_at_h, tempy_at_h, g%Domain, to_all, agrid, halo=1)
957  do j=jsc-1,jec+1 ; do i=isc-1,iec+1
958  dlondx = g%geoLonCu(i,j) - g%geoLonCu(i-1,j)
959  dlondy = g%geoLonCv(i,j) - g%geoLonCv(i,j-1)
960  rdlon = sqrt( dlondx * dlondx + dlondy * dlondy )
961  if (rdlon > 0.) rdlon = 1. / rdlon
962  cosa = dlondx * rdlon
963  sina = dlondy * rdlon
964  zonal_tau = pa_conversion * tempx_at_h(i,j)
965  merid_tau = pa_conversion * tempy_at_h(i,j)
966  tempx_at_h(i,j) = cosa * zonal_tau - sina * merid_tau
967  tempy_at_h(i,j) = sina * zonal_tau + cosa * merid_tau
968  enddo ; enddo
969 
970  ! Average to C-grid locations
971  do j=jsc,jec ; do i=isc-1,iec
972  forces%taux(i,j) = forces%taux(i,j) + 0.5 * ( tempx_at_h(i,j) + tempx_at_h(i+1,j) )
973  enddo ; enddo
974 
975  do j=jsc-1,jec ; do i=isc,iec
976  forces%tauy(i,j) = forces%tauy(i,j) + 0.5 * ( tempy_at_h(i,j) + tempy_at_h(i,j+1) )
977  enddo ; enddo
978  endif ! overrode_x .or. overrode_y
979 
980 end subroutine apply_force_adjustments
981 
982 !> Save any restart files associated with the surface forcing.
983 subroutine forcing_save_restart(CS, G, Time, directory, time_stamped, &
984  filename_suffix)
985  type(surface_forcing_cs), pointer :: cs !< A pointer to the control structure returned
986  !! by a previous call to surface_forcing_init
987  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
988  type(time_type), intent(in) :: time !< The current model time
989  character(len=*), intent(in) :: directory !< The directory into which to write the
990  !! restart files
991  logical, optional, intent(in) :: time_stamped !< If true, the restart file names include
992  !! a unique time stamp. The default is false.
993  character(len=*), optional, intent(in) :: filename_suffix !< An optional suffix (e.g., a time-
994  !! stamp) to append to the restart file names.
995 
996  if (.not.associated(cs)) return
997  if (.not.associated(cs%restart_CSp)) return
998  call save_restart(directory, time, g, cs%restart_CSp, time_stamped)
999 
1000 end subroutine forcing_save_restart
1001 
1002 !> Initialize the surface forcing, including setting parameters and allocating permanent memory.
1003 subroutine surface_forcing_init(Time, G, US, param_file, diag, CS, restore_salt, restore_temp)
1004  type(time_type), intent(in) :: time !< The current model time
1005  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
1006  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1007  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1008  type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate
1009  !! diagnostic output
1010  type(surface_forcing_cs), pointer :: cs !< A pointer that is set to point to the control
1011  !! structure for this module
1012  logical, optional, intent(in) :: restore_salt !< If present and true surface salinity
1013  !! restoring will be applied in this model.
1014  logical, optional, intent(in) :: restore_temp !< If present and true surface temperature
1015  !! restoring will be applied in this model.
1016 
1017  ! Local variables
1018  real :: utide ! The RMS tidal velocity [Z T-1 ~> m s-1].
1019  type(directories) :: dirs
1020  logical :: new_sim, iceberg_flux_diags
1021  type(time_type) :: time_frc
1022  character(len=200) :: tideamp_file, gust_file, salt_file, temp_file ! Input file names.
1023 ! This include declares and sets the variable "version".
1024 #include "version_variable.h"
1025  character(len=40) :: mdl = "MOM_surface_forcing_nuopc" ! This module's name.
1026  character(len=48) :: stagger
1027  character(len=48) :: flnam
1028  character(len=240) :: basin_file
1029  integer :: i, j, isd, ied, jsd, jed
1030 
1031  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1032 
1033  if (associated(cs)) then
1034  call mom_error(warning, "surface_forcing_init called with an associated "// &
1035  "control structure.")
1036  return
1037  endif
1038  allocate(cs)
1039 
1040  id_clock_forcing=cpu_clock_id('Ocean surface forcing', grain=clock_subcomponent)
1041  call cpu_clock_begin(id_clock_forcing)
1042 
1043  cs%diag => diag
1044 
1045  call write_version_number(version)
1046  ! Read all relevant parameters and write them to the model log.
1047  call log_version(param_file, mdl, version, "")
1048 
1049  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, &
1050  "The directory in which all input files are found.", &
1051  default=".")
1052  cs%inputdir = slasher(cs%inputdir)
1053  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
1054  "If true, Temperature and salinity are used as state "//&
1055  "variables.", default=.true.)
1056  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
1057  "The mean ocean density used with BOUSSINESQ true to "//&
1058  "calculate accelerations and the mass for conservation "//&
1059  "properties, or with BOUSSINSEQ false to convert some "//&
1060  "parameters from vertical units of m to kg m-2.", &
1061  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
1062  call get_param(param_file, mdl, "LATENT_HEAT_FUSION", cs%latent_heat_fusion, &
1063  "The latent heat of fusion.", units="J/kg", default=hlf)
1064  call get_param(param_file, mdl, "LATENT_HEAT_VAPORIZATION", cs%latent_heat_vapor, &
1065  "The latent heat of fusion.", units="J/kg", default=hlv)
1066  call get_param(param_file, mdl, "MAX_P_SURF", cs%max_p_surf, &
1067  "The maximum surface pressure that can be exerted by the "//&
1068  "atmosphere and floating sea-ice or ice shelves. This is "//&
1069  "needed because the FMS coupling structure does not "//&
1070  "limit the water that can be frozen out of the ocean and "//&
1071  "the ice-ocean heat fluxes are treated explicitly. No "//&
1072  "limit is applied if a negative value is used.", units="Pa", &
1073  default=-1.0)
1074  call get_param(param_file, mdl, "ADJUST_NET_SRESTORE_TO_ZERO", &
1075  cs%adjust_net_srestore_to_zero, &
1076  "If true, adjusts the salinity restoring seen to zero "//&
1077  "whether restoring is via a salt flux or virtual precip.",&
1078  default=restore_salt)
1079  call get_param(param_file, mdl, "ADJUST_NET_SRESTORE_BY_SCALING", &
1080  cs%adjust_net_srestore_by_scaling, &
1081  "If true, adjustments to salt restoring to achieve zero net are "//&
1082  "made by scaling values without moving the zero contour.",&
1083  default=.false.)
1084  call get_param(param_file, mdl, "ADJUST_NET_FRESH_WATER_TO_ZERO", &
1085  cs%adjust_net_fresh_water_to_zero, &
1086  "If true, adjusts the net fresh-water forcing seen "//&
1087  "by the ocean (including restoring) to zero.", default=.false.)
1088  if (cs%adjust_net_fresh_water_to_zero) &
1089  call get_param(param_file, mdl, "USE_NET_FW_ADJUSTMENT_SIGN_BUG", &
1090  cs%use_net_FW_adjustment_sign_bug, &
1091  "If true, use the wrong sign for the adjustment to "//&
1092  "the net fresh-water.", default=.false.)
1093  call get_param(param_file, mdl, "ADJUST_NET_FRESH_WATER_BY_SCALING", &
1094  cs%adjust_net_fresh_water_by_scaling, &
1095  "If true, adjustments to net fresh water to achieve zero net are "//&
1096  "made by scaling values without moving the zero contour.",&
1097  default=.false.)
1098  call get_param(param_file, mdl, "ICE_SALT_CONCENTRATION", &
1099  cs%ice_salt_concentration, &
1100  "The assumed sea-ice salinity needed to reverse engineer the "//&
1101  "melt flux (or ice-ocean fresh-water flux).", &
1102  units="kg/kg", default=0.005)
1103  call get_param(param_file, mdl, "USE_LIMITED_PATM_SSH", cs%use_limited_P_SSH, &
1104  "If true, return the sea surface height with the "//&
1105  "correction for the atmospheric (and sea-ice) pressure "//&
1106  "limited by max_p_surf instead of the full atmospheric "//&
1107  "pressure.", default=.true.)
1108 
1109  call get_param(param_file, mdl, "WIND_STAGGER", stagger, &
1110  "A case-insensitive character string to indicate the "//&
1111  "staggering of the input wind stress field. Valid "//&
1112  "values are 'A', 'B', or 'C'.", default="C")
1113  if (uppercase(stagger(1:1)) == 'A') then ; cs%wind_stagger = agrid
1114  elseif (uppercase(stagger(1:1)) == 'B') then ; cs%wind_stagger = bgrid_ne
1115  elseif (uppercase(stagger(1:1)) == 'C') then ; cs%wind_stagger = cgrid_ne
1116  else ; call mom_error(fatal,"surface_forcing_init: WIND_STAGGER = "// &
1117  trim(stagger)//" is invalid.") ; endif
1118  call get_param(param_file, mdl, "WIND_STRESS_MULTIPLIER", cs%wind_stress_multiplier, &
1119  "A factor multiplying the wind-stress given to the ocean by the "//&
1120  "coupler. This is used for testing and should be =1.0 for any "//&
1121  "production runs.", default=1.0)
1122 
1123  if (restore_salt) then
1124  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
1125  "The constant that relates the restoring surface fluxes "//&
1126  "to the relative surface anomalies (akin to a piston "//&
1127  "velocity). Note the non-MKS units.", units="m day-1", &
1128  fail_if_missing=.true.)
1129  call get_param(param_file, mdl, "SALT_RESTORE_FILE", cs%salt_restore_file, &
1130  "A file in which to find the surface salinity to use for restoring.", &
1131  default="salt_restore.nc")
1132  call get_param(param_file, mdl, "SALT_RESTORE_VARIABLE", cs%salt_restore_var_name, &
1133  "The name of the surface salinity variable to read from "//&
1134  "SALT_RESTORE_FILE for restoring salinity.", &
1135  default="salt")
1136 ! Convert CS%Flux_const from m day-1 to m s-1.
1137  cs%Flux_const = cs%Flux_const / 86400.0
1138 
1139  call get_param(param_file, mdl, "SRESTORE_AS_SFLUX", cs%salt_restore_as_sflux, &
1140  "If true, the restoring of salinity is applied as a salt "//&
1141  "flux instead of as a freshwater flux.", default=.false.)
1142  call get_param(param_file, mdl, "MAX_DELTA_SRESTORE", cs%max_delta_srestore, &
1143  "The maximum salinity difference used in restoring terms.", &
1144  units="PSU or g kg-1", default=999.0)
1145  call get_param(param_file, mdl, "MASK_SRESTORE_UNDER_ICE", &
1146  cs%mask_srestore_under_ice, &
1147  "If true, disables SSS restoring under sea-ice based on a frazil "//&
1148  "criteria (SST<=Tf). Only used when RESTORE_SALINITY is True.", &
1149  default=.false.)
1150  call get_param(param_file, mdl, "MASK_SRESTORE_MARGINAL_SEAS", &
1151  cs%mask_srestore_marginal_seas, &
1152  "If true, disable SSS restoring in marginal seas. Only used when "//&
1153  "RESTORE_SALINITY is True.", default=.false.)
1154  call get_param(param_file, mdl, "BASIN_FILE", basin_file, &
1155  "A file in which to find the basin masks, in variable 'basin'.", &
1156  default="basin.nc")
1157  basin_file = trim(cs%inputdir) // trim(basin_file)
1158  call safe_alloc_ptr(cs%basin_mask,isd,ied,jsd,jed) ; cs%basin_mask(:,:) = 1.0
1159  if (cs%mask_srestore_marginal_seas) then
1160  call mom_read_data(basin_file,'basin',cs%basin_mask,g%domain, timelevel=1)
1161  do j=jsd,jed ; do i=isd,ied
1162  if (cs%basin_mask(i,j) >= 6.0) then ; cs%basin_mask(i,j) = 0.0
1163  else ; cs%basin_mask(i,j) = 1.0 ; endif
1164  enddo ; enddo
1165  endif
1166  call get_param(param_file, mdl, "MASK_SRESTORE", cs%mask_srestore, &
1167  "If true, read a file (salt_restore_mask) containing "//&
1168  "a mask for SSS restoring.", default=.false.)
1169  endif
1170 
1171  if (restore_temp) then
1172  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
1173  "The constant that relates the restoring surface fluxes "//&
1174  "to the relative surface anomalies (akin to a piston "//&
1175  "velocity). Note the non-MKS units.", units="m day-1", &
1176  fail_if_missing=.true.)
1177  call get_param(param_file, mdl, "SST_RESTORE_FILE", cs%temp_restore_file, &
1178  "A file in which to find the surface temperature to use for restoring.", &
1179  default="temp_restore.nc")
1180  call get_param(param_file, mdl, "SST_RESTORE_VARIABLE", cs%temp_restore_var_name, &
1181  "The name of the surface temperature variable to read from "//&
1182  "SST_RESTORE_FILE for restoring sst.", &
1183  default="temp")
1184  ! Convert CS%Flux_const from m day-1 to m s-1.
1185  cs%Flux_const = cs%Flux_const / 86400.0
1186 
1187  call get_param(param_file, mdl, "MAX_DELTA_TRESTORE", cs%max_delta_trestore, &
1188  "The maximum sst difference used in restoring terms.", &
1189  units="degC ", default=999.0)
1190  call get_param(param_file, mdl, "MASK_TRESTORE", cs%mask_trestore, &
1191  "If true, read a file (temp_restore_mask) containing "//&
1192  "a mask for SST restoring.", default=.false.)
1193 
1194  endif
1195 
1196 ! Optionally read tidal amplitude from input file (m s-1) on model grid.
1197 ! Otherwise use default tidal amplitude for bottom frictionally-generated
1198 ! dissipation. Default cd_tides is chosen to yield approx 1 TWatt of
1199 ! work done against tides globally using OSU tidal amplitude.
1200  call get_param(param_file, mdl, "CD_TIDES", cs%cd_tides, &
1201  "The drag coefficient that applies to the tides.", &
1202  units="nondim", default=1.0e-4)
1203  call get_param(param_file, mdl, "READ_TIDEAMP", cs%read_TIDEAMP, &
1204  "If true, read a file (given by TIDEAMP_FILE) containing "//&
1205  "the tidal amplitude with INT_TIDE_DISSIPATION.", default=.false.)
1206  if (cs%read_TIDEAMP) then
1207  call get_param(param_file, mdl, "TIDEAMP_FILE", tideamp_file, &
1208  "The path to the file containing the spatially varying "//&
1209  "tidal amplitudes with INT_TIDE_DISSIPATION.", &
1210  default="tideamp.nc")
1211  cs%utide=0.0
1212  else
1213  call get_param(param_file, mdl, "UTIDE", cs%utide, &
1214  "The constant tidal amplitude used with INT_TIDE_DISSIPATION.", &
1215  units="m s-1", default=0.0, scale=us%m_to_Z*us%T_to_s)
1216  endif
1217 
1218  call safe_alloc_ptr(cs%TKE_tidal,isd,ied,jsd,jed)
1219  call safe_alloc_ptr(cs%ustar_tidal,isd,ied,jsd,jed)
1220 
1221  if (cs%read_TIDEAMP) then
1222  tideamp_file = trim(cs%inputdir) // trim(tideamp_file)
1223  call mom_read_data(tideamp_file,'tideamp',cs%TKE_tidal,g%domain,timelevel=1, scale=us%m_to_Z*us%T_to_s)
1224  do j=jsd, jed; do i=isd, ied
1225  utide = cs%TKE_tidal(i,j)
1226  cs%TKE_tidal(i,j) = g%mask2dT(i,j)*cs%Rho0*cs%cd_tides*(utide*utide*utide)
1227  cs%ustar_tidal(i,j) = sqrt(cs%cd_tides)*utide
1228  enddo ; enddo
1229  else
1230  do j=jsd,jed; do i=isd,ied
1231  utide = cs%utide
1232  cs%TKE_tidal(i,j) = cs%Rho0*cs%cd_tides*(utide*utide*utide)
1233  cs%ustar_tidal(i,j) = sqrt(cs%cd_tides)*utide
1234  enddo ; enddo
1235  endif
1236 
1237  call time_interp_external_init
1238 
1239 ! Optionally read a x-y gustiness field in place of a global
1240 ! constant.
1241 
1242  call get_param(param_file, mdl, "READ_GUST_2D", cs%read_gust_2d, &
1243  "If true, use a 2-dimensional gustiness supplied from "//&
1244  "an input file", default=.false.)
1245  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
1246  "The background gustiness in the winds.", &
1247  units="Pa", default=0.02, scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1248  if (cs%read_gust_2d) then
1249  call get_param(param_file, mdl, "GUST_2D_FILE", gust_file, &
1250  "The file in which the wind gustiness is found in "//&
1251  "variable gustiness.")
1252 
1253  call safe_alloc_ptr(cs%gust,isd,ied,jsd,jed)
1254  gust_file = trim(cs%inputdir) // trim(gust_file)
1255  call mom_read_data(gust_file,'gustiness',cs%gust,g%domain, timelevel=1, &
1256  scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z) ! units in file should be Pa
1257  endif
1258 
1259 ! See whether sufficiently thick sea ice should be treated as rigid.
1260  call get_param(param_file, mdl, "USE_RIGID_SEA_ICE", cs%rigid_sea_ice, &
1261  "If true, sea-ice is rigid enough to exert a "//&
1262  "nonhydrostatic pressure that resist vertical motion.", &
1263  default=.false.)
1264  if (cs%rigid_sea_ice) then
1265  call get_param(param_file, mdl, "G_EARTH", cs%g_Earth, &
1266  "The gravitational acceleration of the Earth.", &
1267  units="m s-2", default = 9.80)
1268  call get_param(param_file, mdl, "SEA_ICE_MEAN_DENSITY", cs%density_sea_ice, &
1269  "A typical density of sea ice, used with the kinematic "//&
1270  "viscosity, when USE_RIGID_SEA_ICE is true.", units="kg m-3", &
1271  default=900.0)
1272  call get_param(param_file, mdl, "SEA_ICE_VISCOSITY", cs%Kv_sea_ice, &
1273  "The kinematic viscosity of sufficiently thick sea ice "//&
1274  "for use in calculating the rigidity of sea ice.", &
1275  units="m2 s-1", default=1.0e9)
1276  call get_param(param_file, mdl, "SEA_ICE_RIGID_MASS", cs%rigid_sea_ice_mass, &
1277  "The mass of sea-ice per unit area at which the sea-ice "//&
1278  "starts to exhibit rigidity", units="kg m-2", default=1000.0)
1279  endif
1280 
1281  call get_param(param_file, mdl, "ALLOW_ICEBERG_FLUX_DIAGNOSTICS", iceberg_flux_diags, &
1282  "If true, makes available diagnostics of fluxes from icebergs "//&
1283  "as seen by MOM6.", default=.false.)
1284  call register_forcing_type_diags(time, diag, us, cs%use_temperature, cs%handles, &
1285  use_berg_fluxes=iceberg_flux_diags)
1286 
1287  call get_param(param_file, mdl, "ALLOW_FLUX_ADJUSTMENTS", cs%allow_flux_adjustments, &
1288  "If true, allows flux adjustments to specified via the "//&
1289  "data_table using the component name 'OCN'.", default=.false.)
1290 
1291  call get_param(param_file, mdl, "LIQUID_RUNOFF_FROM_DATA", cs%liquid_runoff_from_data, &
1292  "If true, allows liquid river runoff to be specified via the "//&
1293  "data_table using the component name 'OCN'.", default=.false.)
1294 
1295  if (cs%allow_flux_adjustments .or. cs%liquid_runoff_from_data) then
1296  call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
1297  endif
1298 
1299  if (present(restore_salt)) then ; if (restore_salt) then
1300  salt_file = trim(cs%inputdir) // trim(cs%salt_restore_file)
1301  cs%id_srestore = init_external_field(salt_file, cs%salt_restore_var_name, domain=g%Domain%mpp_domain)
1302  call safe_alloc_ptr(cs%srestore_mask,isd,ied,jsd,jed); cs%srestore_mask(:,:) = 1.0
1303  if (cs%mask_srestore) then ! read a 2-d file containing a mask for restoring fluxes
1304  flnam = trim(cs%inputdir) // 'salt_restore_mask.nc'
1305  call mom_read_data(flnam,'mask', cs%srestore_mask, g%domain, timelevel=1)
1306  endif
1307  endif ; endif
1308 
1309  if (present(restore_temp)) then ; if (restore_temp) then
1310  temp_file = trim(cs%inputdir) // trim(cs%temp_restore_file)
1311  cs%id_trestore = init_external_field(temp_file, cs%temp_restore_var_name, domain=g%Domain%mpp_domain)
1312  call safe_alloc_ptr(cs%trestore_mask,isd,ied,jsd,jed); cs%trestore_mask(:,:) = 1.0
1313  if (cs%mask_trestore) then ! read a 2-d file containing a mask for restoring fluxes
1314  flnam = trim(cs%inputdir) // 'temp_restore_mask.nc'
1315  call mom_read_data(flnam, 'mask', cs%trestore_mask, g%domain, timelevel=1)
1316  endif
1317  endif ; endif
1318 
1319  ! Set up any restart fields associated with the forcing.
1320  call restart_init(param_file, cs%restart_CSp, "MOM_forcing.res")
1321  call restart_init_end(cs%restart_CSp)
1322 
1323  if (associated(cs%restart_CSp)) then
1324  call get_mom_input(dirs=dirs)
1325 
1326  new_sim = .false.
1327  if ((dirs%input_filename(1:1) == 'n') .and. &
1328  (len_trim(dirs%input_filename) == 1)) new_sim = .true.
1329  if (.not.new_sim) then
1330  call restore_state(dirs%input_filename, dirs%restart_input_dir, time_frc, &
1331  g, cs%restart_CSp)
1332  endif
1333  endif
1334 
1335  call user_revise_forcing_init(param_file, cs%urf_CS)
1336 
1337  call cpu_clock_end(id_clock_forcing)
1338 end subroutine surface_forcing_init
1339 
1340 !> Clean up and deallocate any memory associated with this module and its children.
1341 subroutine surface_forcing_end(CS, fluxes)
1342  type(surface_forcing_cs), pointer :: cs !< A pointer to the control structure returned by
1343  !! a previous call to surface_forcing_init, it will
1344  !! be deallocated here.
1345  type(forcing), optional, intent(inout) :: fluxes !< A structure containing pointers to all
1346  !! possible mass, heat or salt flux forcing fields.
1347  !! If present, it will be deallocated here.
1348 
1349  if (present(fluxes)) call deallocate_forcing_type(fluxes)
1350 
1351  if (associated(cs)) deallocate(cs)
1352  cs => null()
1353 
1354 end subroutine surface_forcing_end
1355 
1356 !> Write out a set of messages with checksums of the fields in an ice_ocen_boundary type
1357 subroutine ice_ocn_bnd_type_chksum(id, timestep, iobt)
1359  character(len=*), intent(in) :: id !< An identifying string for this call
1360  integer, intent(in) :: timestep !< The number of elapsed timesteps
1361  type(ice_ocean_boundary_type), &
1362  intent(in) :: iobt !< An ice-ocean boundary type with fluxes to drive the
1363  !! ocean in a coupled model whose checksums are reported
1364 
1365  ! local variables
1366  integer :: n,m, outunit
1367 
1368  outunit = stdout()
1369 
1370  write(outunit,*) "BEGIN CHECKSUM(ice_ocean_boundary_type):: ", id, timestep
1371  write(outunit,100) 'iobt%u_flux ' , mpp_chksum( iobt%u_flux )
1372  write(outunit,100) 'iobt%v_flux ' , mpp_chksum( iobt%v_flux )
1373  write(outunit,100) 'iobt%t_flux ' , mpp_chksum( iobt%t_flux )
1374  write(outunit,100) 'iobt%q_flux ' , mpp_chksum( iobt%q_flux )
1375  write(outunit,100) 'iobt%salt_flux ' , mpp_chksum( iobt%salt_flux )
1376  write(outunit,100) 'iobt%seaice_melt_heat' , mpp_chksum( iobt%seaice_melt_heat)
1377  write(outunit,100) 'iobt%seaice_melt ' , mpp_chksum( iobt%seaice_melt )
1378  write(outunit,100) 'iobt%lw_flux ' , mpp_chksum( iobt%lw_flux )
1379  write(outunit,100) 'iobt%sw_flux_vis_dir' , mpp_chksum( iobt%sw_flux_vis_dir)
1380  write(outunit,100) 'iobt%sw_flux_vis_dif' , mpp_chksum( iobt%sw_flux_vis_dif)
1381  write(outunit,100) 'iobt%sw_flux_nir_dir' , mpp_chksum( iobt%sw_flux_nir_dir)
1382  write(outunit,100) 'iobt%sw_flux_nir_dif' , mpp_chksum( iobt%sw_flux_nir_dif)
1383  write(outunit,100) 'iobt%lprec ' , mpp_chksum( iobt%lprec )
1384  write(outunit,100) 'iobt%fprec ' , mpp_chksum( iobt%fprec )
1385  write(outunit,100) 'iobt%lrunoff ' , mpp_chksum( iobt%lrunoff )
1386  write(outunit,100) 'iobt%frunoff ' , mpp_chksum( iobt%frunoff )
1387  write(outunit,100) 'iobt%p ' , mpp_chksum( iobt%p )
1388  if (associated(iobt%ustar_berg)) &
1389  write(outunit,100) 'iobt%ustar_berg ' , mpp_chksum( iobt%ustar_berg )
1390  if (associated(iobt%area_berg)) &
1391  write(outunit,100) 'iobt%area_berg ' , mpp_chksum( iobt%area_berg )
1392  if (associated(iobt%mass_berg)) &
1393  write(outunit,100) 'iobt%mass_berg ' , mpp_chksum( iobt%mass_berg )
1394 100 FORMAT(" CHECKSUM::",a20," = ",z20)
1395 
1396  call coupler_type_write_chksums(iobt%fluxes, outunit, 'iobt%')
1397 
1398 end subroutine ice_ocn_bnd_type_chksum
1399 
1400 end module mom_surface_forcing_nuopc
user_revise_forcing::user_revise_forcing_cs
Control structure for user_revise_forcing.
Definition: user_revise_forcing.F90:23
mom_restart::restore_state
subroutine, public restore_state(filename, directory, day, G, CS)
restore_state reads the model state from previously generated files. All restart variables are read f...
Definition: MOM_restart.F90:984
mom_forcing_type::allocate_forcing_type
subroutine, public allocate_forcing_type(G, fluxes, water, heat, ustar, press, shelf, iceberg, salt)
Conditionally allocate fields within the forcing type.
Definition: MOM_forcing_type.F90:2811
mom_forcing_type::mech_forcing
Structure that contains pointers to the mechanical forcing at the surface used to drive the liquid oc...
Definition: MOM_forcing_type.F90:187
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_spatial_means
Functions and routines to take area, volume, mass-weighted, layerwise, zonal or meridional means.
Definition: MOM_spatial_means.F90:2
mom_spatial_means::adjust_area_mean_to_zero
subroutine, public adjust_area_mean_to_zero(array, G, scaling, unit_scale)
Adjust 2d array such that area mean is zero without moving the zero contour.
Definition: MOM_spatial_means.F90:357
mom_surface_forcing_nuopc::ice_ocean_boundary_type
Structure corresponding to forcing, but with the elements, units, and conventions that exactly confor...
Definition: mom_surface_forcing_nuopc.F90:155
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_constants
Provides a few physical constants.
Definition: MOM_constants.F90:2
user_revise_forcing
Provides a template for users to code updating the forcing fluxes.
Definition: user_revise_forcing.F90:2
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_forcing_type::deallocate_forcing_type
subroutine, public deallocate_forcing_type(fluxes)
Deallocate the forcing type.
Definition: MOM_forcing_type.F90:2931
mom_get_input::directories
Container for paths and parameter file names.
Definition: MOM_get_input.F90:20
mom_forcing_type::allocate_mech_forcing
subroutine, public allocate_mech_forcing(G, forces, stress, ustar, shelf, press, iceberg)
Conditionally allocate fields within the mechanical forcing type.
Definition: MOM_forcing_type.F90:2879
mom_surface_forcing_nuopc::apply_flux_adjustments
subroutine, private apply_flux_adjustments(G, US, CS, Time, fluxes)
Adds thermodynamic flux adjustments obtained via data_override Component name is 'OCN' Available adju...
Definition: mom_surface_forcing_nuopc.F90:881
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_surface_forcing_nuopc
Converts the input ESMF data (import data) to a MOM-specific data type (surface_forcing_CS).
Definition: mom_surface_forcing_nuopc.F90:2
mom_error_handler::mom_mesg
subroutine, public mom_mesg(message, verb, all_print)
This provides a convenient interface for writing an informative comment.
Definition: MOM_error_handler.F90:53
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_string_functions::uppercase
character(len=len(input_string)) function, public uppercase(input_string)
Return a string in which all uppercase letters have been replaced by their lowercase counterparts.
Definition: MOM_string_functions.F90:42
mom_restart::mom_restart_cs
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:72
mom_surface_forcing_nuopc::surface_forcing_end
subroutine, private surface_forcing_end(CS, fluxes)
Clean up and deallocate any memory associated with this module and its children.
Definition: mom_surface_forcing_nuopc.F90:1342
mom_forcing_type::mech_forcing_diags
subroutine, public mech_forcing_diags(forces, dt, G, time_end, diag, handles)
Offer mechanical forcing fields for diagnostics for those fields registered as part of register_forci...
Definition: MOM_forcing_type.F90:2201
mom_get_input
Reads the only Fortran name list needed to boot-strap the model.
Definition: MOM_get_input.F90:6
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_surface_forcing_nuopc::ice_ocn_bnd_type_chksum
subroutine, public ice_ocn_bnd_type_chksum(id, timestep, iobt)
Write out a set of messages with checksums of the fields in an ice_ocen_boundary type.
Definition: mom_surface_forcing_nuopc.F90:1358
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_domains::to_all
integer, parameter, public to_all
A flag for passing in all directions.
Definition: MOM_domains.F90:132
mom_restart
The MOM6 facility for reading and writing restart files, and querying what has been read.
Definition: MOM_restart.F90:2
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_surface_forcing_nuopc::convert_iob_to_fluxes
subroutine, public convert_iob_to_fluxes(IOB, fluxes, index_bounds, Time, valid_time, G, US, CS, sfc_state, restore_salt, restore_temp)
This subroutine translates the Ice_ocean_boundary_type into a MOM thermodynamic forcing type,...
Definition: mom_surface_forcing_nuopc.F90:203
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
user_revise_forcing::user_alter_forcing
subroutine, public user_alter_forcing(state, fluxes, day, G, CS)
This subroutine sets the surface wind stresses.
Definition: user_revise_forcing.F90:34
mom_cpu_clock
Wraps the MPP cpu clock functions.
Definition: MOM_cpu_clock.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_get_input::get_mom_input
subroutine, public get_mom_input(param_file, dirs, check_params, default_input_filename, ensemble_num)
Get the names of the I/O directories and initialization file. Also calls the subroutine that opens ru...
Definition: MOM_get_input.F90:34
mom_surface_forcing_nuopc::surface_forcing_cs
Contains pointers to the forcing fields which may be used to drive MOM. All fluxes are positive downw...
Definition: mom_surface_forcing_nuopc.F90:59
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
user_revise_forcing::user_revise_forcing_init
subroutine, public user_revise_forcing_init(param_file, CS)
Initialize the user_revise_forcing control structure.
Definition: user_revise_forcing.F90:49
mom_restart::register_restart_field
Register fields for restarts.
Definition: MOM_restart.F90:107
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_error_handler::is_root_pe
logical function, public is_root_pe()
This returns .true. if the current PE is the root PE.
Definition: MOM_error_handler.F90:44
mom_restart::restart_init
subroutine, public restart_init(param_file, CS, restart_root)
Initialize this module and set up a restart control structure.
Definition: MOM_restart.F90:1421
mom_surface_forcing_nuopc::forcing_save_restart
subroutine, public forcing_save_restart(CS, G, Time, directory, time_stamped, filename_suffix)
Save any restart files associated with the surface forcing.
Definition: mom_surface_forcing_nuopc.F90:985
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:50
mom_surface_forcing_nuopc::id_clock_forcing
integer id_clock_forcing
Definition: mom_surface_forcing_nuopc.F90:194
mom_restart::restart_init_end
subroutine, public restart_init_end(CS)
Indicate that all variables have now been registered.
Definition: MOM_restart.F90:1479
mom_surface_forcing_nuopc::apply_force_adjustments
subroutine, private apply_force_adjustments(G, US, CS, Time, forces)
Adds mechanical forcing adjustments obtained via data_override Component name is 'OCN' Available adju...
Definition: mom_surface_forcing_nuopc.F90:927
mom_error_handler::mom_error
subroutine, public mom_error(level, message, all_print)
This provides a convenient interface for writing an mpp_error message with run-time filter based on a...
Definition: MOM_error_handler.F90:72
mom_restart::save_restart
subroutine, public save_restart(directory, time, G, CS, time_stamped, filename, GV)
save_restart saves all registered variables to restart files.
Definition: MOM_restart.F90:781
mom_forcing_type::forcing_diags
Structure that defines the id handles for the forcing type.
Definition: MOM_forcing_type.F90:237
mom_coms::reproducing_sum
Find an accurate and order-invariant sum of distributed 2d or 3d fields.
Definition: MOM_coms.F90:54
mom_surface_forcing_nuopc::convert_iob_to_forces
subroutine, public convert_iob_to_forces(IOB, forces, index_bounds, Time, G, US, CS)
This subroutine translates the Ice_ocean_boundary_type into a MOM mechanical forcing type,...
Definition: mom_surface_forcing_nuopc.F90:586
mom_forcing_type::deallocate_mech_forcing
subroutine, public deallocate_mech_forcing(forces)
Deallocate the mechanical forcing type.
Definition: MOM_forcing_type.F90:2983
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_surface_forcing_nuopc::surface_forcing_init
subroutine, public surface_forcing_init(Time, G, US, param_file, diag, CS, restore_salt, restore_temp)
Initialize the surface forcing, including setting parameters and allocating permanent memory.
Definition: mom_surface_forcing_nuopc.F90:1004
mom_forcing_type::register_forcing_type_diags
subroutine, public register_forcing_type_diags(Time, diag, US, use_temperature, handles, use_berg_fluxes)
Register members of the forcing type for diagnostics.
Definition: MOM_forcing_type.F90:1221
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
mom_domains::fill_symmetric_edges
Do a set of halo updates that fill in the values at the duplicated edges of a staggered symmetric mem...
Definition: MOM_domains.F90:88
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239