MOM6
MOM_diagnostics.F90
1 !> Calculates any requested diagnostic quantities
2 !! that are not calculated in the various subroutines.
3 !! Diagnostic quantities are requested by allocating them memory.
5 
6 ! This file is part of MOM6. See LICENSE.md for the license.
7 
8 use mom_coms, only : reproducing_sum
9 use mom_density_integrals, only : int_density_dz
10 use mom_diag_mediator, only : post_data, get_diag_time_end
11 use mom_diag_mediator, only : register_diag_field, register_scalar_field
12 use mom_diag_mediator, only : register_static_field, diag_register_area_ids
13 use mom_diag_mediator, only : diag_ctrl, time_type, safe_alloc_ptr
14 use mom_diag_mediator, only : diag_get_volume_cell_measure_dm_id
16 use mom_diag_mediator, only : diag_save_grids, diag_restore_grids, diag_copy_storage_to_diag
17 use mom_domains, only : create_group_pass, do_group_pass, group_pass_type
18 use mom_domains, only : to_north, to_east
19 use mom_eos, only : calculate_density, calculate_density_derivs, eos_domain
20 use mom_eos, only : gsw_sp_from_sr, gsw_pt_from_ct
21 use mom_error_handler, only : mom_error, fatal, warning
23 use mom_grid, only : ocean_grid_type
25 use mom_spatial_means, only : global_area_mean, global_layer_mean
26 use mom_spatial_means, only : global_volume_mean, global_area_integral
27 use mom_tracer_registry, only : tracer_registry_type, post_tracer_transport_diagnostics
31 use mom_verticalgrid, only : verticalgrid_type, get_thickness_units
32 use mom_wave_speed, only : wave_speed, wave_speed_cs, wave_speed_init
33 use coupler_types_mod, only : coupler_type_send_data
34 
35 implicit none ; private
36 
37 #include <MOM_memory.h>
38 
39 public calculate_diagnostic_fields, register_time_deriv, write_static_fields
40 public find_eta
41 public register_surface_diags, post_surface_dyn_diags, post_surface_thermo_diags
42 public register_transport_diags, post_transport_diagnostics
43 public mom_diagnostics_init, mom_diagnostics_end
44 
45 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
46 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
47 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
48 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
49 
50 !> The control structure for the MOM_diagnostics module
51 type, public :: diagnostics_cs ; private
52  real :: mono_n2_column_fraction = 0. !< The lower fraction of water column over which N2 is limited as
53  !! monotonic for the purposes of calculating the equivalent
54  !! barotropic wave speed.
55  real :: mono_n2_depth = -1. !< The depth below which N2 is limited as monotonic for the purposes of
56  !! calculating the equivalent barotropic wave speed [Z ~> m].
57 
58  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
59  !! regulate the timing of diagnostic output.
60 
61  ! following arrays store diagnostics calculated here and unavailable outside.
62 
63  ! following fields have nz+1 levels.
64  real, pointer, dimension(:,:,:) :: &
65  e => null(), & !< interface height [Z ~> m]
66  e_d => null() !< interface height above bottom [Z ~> m]
67 
68  ! following fields have nz layers.
69  real, pointer, dimension(:,:,:) :: &
70  du_dt => null(), & !< net i-acceleration [L T-2 ~> m s-2]
71  dv_dt => null(), & !< net j-acceleration [L T-2 ~> m s-2]
72  dh_dt => null(), & !< thickness rate of change [H T-1 ~> m s-1 or kg m-2 s-1]
73  p_ebt => null() !< Equivalent barotropic modal structure [nondim]
74  ! hf_du_dt => NULL(), hf_dv_dt => NULL() !< du_dt, dv_dt x fract. thickness [L T-2 ~> m s-2].
75  ! 3D diagnostics hf_du(dv)_dt are commented because there is no clarity on proper remapping grid option.
76  ! The code is retained for degugging purposes in the future.
77 
78  real, pointer, dimension(:,:,:) :: h_rlay => null() !< Layer thicknesses in potential density
79  !! coordinates [H ~> m or kg m-2]
80  real, pointer, dimension(:,:,:) :: uh_rlay => null() !< Zonal transports in potential density
81  !! coordinates [H m2 s-1 ~> m3 s-1 or kg s-1]
82  real, pointer, dimension(:,:,:) :: vh_rlay => null() !< Meridional transports in potential density
83  !! coordinates [H m2 s-1 ~> m3 s-1 or kg s-1]
84  real, pointer, dimension(:,:,:) :: uhgm_rlay => null() !< Zonal Gent-McWilliams transports in potential density
85  !! coordinates [H m2 s-1 ~> m3 s-1 or kg s-1]
86  real, pointer, dimension(:,:,:) :: vhgm_rlay => null() !< Meridional Gent-McWilliams transports in potential density
87  !! coordinates [H m2 s-1 ~> m3 s-1 or kg s-1]
88 
89  ! following fields are 2-D.
90  real, pointer, dimension(:,:) :: &
91  cg1 => null(), & !< First baroclinic gravity wave speed [L T-1 ~> m s-1]
92  rd1 => null(), & !< First baroclinic deformation radius [L ~> m]
93  cfl_cg1 => null(), & !< CFL for first baroclinic gravity wave speed [nondim]
94  cfl_cg1_x => null(), & !< i-component of CFL for first baroclinic gravity wave speed [nondim]
95  cfl_cg1_y => null() !< j-component of CFL for first baroclinic gravity wave speed [nondim]
96 
97  ! The following arrays hold diagnostics in the layer-integrated energy budget.
98  real, pointer, dimension(:,:,:) :: &
99  ke => null(), & !< KE per unit mass [L2 T-2 ~> m2 s-2]
100  dke_dt => null(), & !< time derivative of the layer KE [H L2 T-3 ~> m3 s-3]
101  pe_to_ke => null(), & !< potential energy to KE term [m3 s-3]
102  ke_coradv => null(), & !< KE source from the combined Coriolis and
103  !! advection terms [H L2 T-3 ~> m3 s-3].
104  !! The Coriolis source should be zero, but is not due to truncation
105  !! errors. There should be near-cancellation of the global integral
106  !! of this spurious Coriolis source.
107  ke_adv => null(), & !< KE source from along-layer advection [H L2 T-3 ~> m3 s-3]
108  ke_visc => null(), & !< KE source from vertical viscosity [H L2 T-3 ~> m3 s-3]
109  ke_horvisc => null(), & !< KE source from horizontal viscosity [H L2 T-3 ~> m3 s-3]
110  ke_dia => null() !< KE source from diapycnal diffusion [H L2 T-3 ~> m3 s-3]
111 
112  !>@{ Diagnostic IDs
113  integer :: id_u = -1, id_v = -1, id_h = -1
114  integer :: id_e = -1, id_e_d = -1
115  integer :: id_du_dt = -1, id_dv_dt = -1
116  ! integer :: id_hf_du_dt = -1, id_hf_dv_dt = -1
117  integer :: id_hf_du_dt_2d = -1, id_hf_dv_dt_2d = -1
118  integer :: id_col_ht = -1, id_dh_dt = -1
119  integer :: id_ke = -1, id_dkedt = -1
120  integer :: id_pe_to_ke = -1, id_ke_coradv = -1
121  integer :: id_ke_adv = -1, id_ke_visc = -1
122  integer :: id_ke_horvisc = -1, id_ke_dia = -1
123  integer :: id_uh_rlay = -1, id_vh_rlay = -1
124  integer :: id_uhgm_rlay = -1, id_vhgm_rlay = -1
125  integer :: id_h_rlay = -1, id_rd1 = -1
126  integer :: id_rml = -1, id_rcv = -1
127  integer :: id_cg1 = -1, id_cfl_cg1 = -1
128  integer :: id_cfl_cg1_x = -1, id_cfl_cg1_y = -1
129  integer :: id_cg_ebt = -1, id_rd_ebt = -1
130  integer :: id_p_ebt = -1
131  integer :: id_temp_int = -1, id_salt_int = -1
132  integer :: id_mass_wt = -1, id_col_mass = -1
133  integer :: id_masscello = -1, id_masso = -1
134  integer :: id_volcello = -1
135  integer :: id_tpot = -1, id_sprac = -1
136  integer :: id_tob = -1, id_sob = -1
137  integer :: id_thetaoga = -1, id_soga = -1
138  integer :: id_sosga = -1, id_tosga = -1
139  integer :: id_temp_layer_ave = -1, id_salt_layer_ave = -1
140  integer :: id_pbo = -1
141  integer :: id_thkcello = -1, id_rhoinsitu = -1
142  integer :: id_rhopot0 = -1, id_rhopot2 = -1
143  integer :: id_drho_dt = -1, id_drho_ds = -1
144  integer :: id_h_pre_sync = -1
145  !>@}
146  !> The control structure for calculating wave speed.
147  type(wave_speed_cs), pointer :: wave_speed_csp => null()
148 
149  type(p3d) :: var_ptr(max_fields_) !< pointers to variables used in the calculation
150  !! of time derivatives
151  type(p3d) :: deriv(max_fields_) !< Time derivatives of various fields
152  type(p3d) :: prev_val(max_fields_) !< Previous values of variables used in the calculation
153  !! of time derivatives
154  !< previous values of variables used in calculation of time derivatives
155  integer :: nlay(max_fields_) !< The number of layers in each diagnostics
156  integer :: num_time_deriv = 0 !< The number of time derivative diagnostics
157 
158  type(group_pass_type) :: pass_ke_uv !< A handle used for group halo passes
159 
160 end type diagnostics_cs
161 
162 
163 !> A structure with diagnostic IDs of the surface and integrated variables
164 type, public :: surface_diag_ids ; private
165  !>@{ Diagnostic IDs for 2-d surface and bottom flux and state fields
166  !Diagnostic IDs for 2-d surface and bottom fields
167  integer :: id_zos = -1, id_zossq = -1
168  integer :: id_volo = -1, id_speed = -1
169  integer :: id_ssh = -1, id_ssh_ga = -1
170  integer :: id_sst = -1, id_sst_sq = -1, id_sstcon = -1
171  integer :: id_sss = -1, id_sss_sq = -1, id_sssabs = -1
172  integer :: id_ssu = -1, id_ssv = -1
173 
174  ! Diagnostic IDs for heat and salt flux fields
175  integer :: id_fraz = -1
176  integer :: id_salt_deficit = -1
177  integer :: id_heat_pme = -1
178  integer :: id_intern_heat = -1
179  !>@}
180 end type surface_diag_ids
181 
182 
183 !> A structure with diagnostic IDs of mass transport related diagnostics
184 type, public :: transport_diag_ids ; private
185  !>@{ Diagnostics for tracer horizontal transport
186  integer :: id_uhtr = -1, id_umo = -1, id_umo_2d = -1
187  integer :: id_vhtr = -1, id_vmo = -1, id_vmo_2d = -1
188  integer :: id_dynamics_h = -1, id_dynamics_h_tendency = -1
189  !>@}
190 end type transport_diag_ids
191 
192 
193 contains
194 !> Diagnostics not more naturally calculated elsewhere are computed here.
195 subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
196  dt, diag_pre_sync, G, GV, US, CS, eta_bt)
197  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
198  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
199  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
200  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
201  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
202  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
203  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
204  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
205  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
206  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
207  intent(in) :: uh !< Transport through zonal faces = u*h*dy,
208  !! [H L2 T-1 ~> m3 s-1 or kg s-1].
209  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
210  intent(in) :: vh !< Transport through meridional faces = v*h*dx,
211  !! [H L2 T-1 ~> m3 s-1 or kg s-1].
212  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
213  !! thermodynamic variables.
214  type(accel_diag_ptrs), intent(in) :: adp !< structure with pointers to
215  !! accelerations in momentum equation.
216  type(cont_diag_ptrs), intent(in) :: cdp !< structure with pointers to
217  !! terms in continuity equation.
218  real, dimension(:,:), pointer :: p_surf !< A pointer to the surface pressure [R L2 T-2 ~> Pa].
219  !! If p_surf is not associated, it is the same
220  !! as setting the surface pressure to 0.
221  real, intent(in) :: dt !< The time difference since the last
222  !! call to this subroutine [T ~> s].
223  type(diag_grid_storage), intent(in) :: diag_pre_sync !< Target grids from previous timestep
224  type(diagnostics_cs), intent(inout) :: cs !< Control structure returned by a
225  !! previous call to diagnostics_init.
226  real, dimension(SZI_(G),SZJ_(G)), &
227  optional, intent(in) :: eta_bt !< An optional barotropic
228  !! variable that gives the "correct" free surface height (Boussinesq) or total water column
229  !! mass per unit area (non-Boussinesq). This is used to dilate the layer thicknesses when
230  !! calculating interface heights [H ~> m or kg m-2].
231 
232  ! Local variables
233  integer, dimension(2) :: eosdom ! The i-computational domain for the equation of state
234  integer i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz, nkmb
235 
236  real :: rcv(szi_(g),szj_(g),szk_(g)) ! Coordinate variable potential density [R ~> kg m-3].
237  real :: work_3d(szi_(g),szj_(g),szk_(g)) ! A 3-d temporary work array.
238  real :: work_2d(szi_(g),szj_(g)) ! A 2-d temporary work array.
239  real :: rho_in_situ(szi_(g)) ! In situ density [R ~> kg m-3]
240 
241  real, allocatable, dimension(:,:) :: &
242  hf_du_dt_2d, hf_dv_dt_2d ! z integeral of hf_du_dt, hf_dv_dt [L T-2 ~> m s-2].
243 
244  ! tmp array for surface properties
245  real :: surface_field(szi_(g),szj_(g))
246  real :: pressure_1d(szi_(g)) ! Temporary array for pressure when calling EOS [R L2 T-2 ~> Pa]
247  real :: wt, wt_p
248 
249  real :: f2_h ! Squared Coriolis parameter at to h-points [T-2 ~> s-2]
250  real :: mag_beta ! Magnitude of the gradient of f [T-1 L-1 ~> s-1 m-1]
251  real :: absurdly_small_freq2 ! Frequency squared used to avoid division by 0 [T-2 ~> s-2]
252 
253  integer :: k_list
254 
255  real, dimension(SZK_(G)) :: temp_layer_ave, salt_layer_ave
256  real :: thetaoga, soga, masso, tosga, sosga
257 
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  nz = g%ke ; nkmb = gv%nk_rho_varies
261 
262  ! This value is roughly (pi / (the age of the universe) )^2.
263  absurdly_small_freq2 = 1e-34*us%T_to_s**2
264 
265  if (loc(cs)==0) call mom_error(fatal, &
266  "calculate_diagnostic_fields: Module must be initialized before used.")
267 
268  call calculate_derivs(dt, g, cs)
269 
270  if (dt > 0.0) then
271  call diag_save_grids(cs%diag)
272  call diag_copy_storage_to_diag(cs%diag, diag_pre_sync)
273 
274  if (cs%id_h_pre_sync > 0) &
275  call post_data(cs%id_h_pre_sync, diag_pre_sync%h_state, cs%diag, alt_h = diag_pre_sync%h_state)
276 
277  if (cs%id_du_dt>0) call post_data(cs%id_du_dt, cs%du_dt, cs%diag, alt_h = diag_pre_sync%h_state)
278 
279  if (cs%id_dv_dt>0) call post_data(cs%id_dv_dt, cs%dv_dt, cs%diag, alt_h = diag_pre_sync%h_state)
280 
281  if (cs%id_dh_dt>0) call post_data(cs%id_dh_dt, cs%dh_dt, cs%diag, alt_h = diag_pre_sync%h_state)
282 
283  !! Diagnostics for terms multiplied by fractional thicknesses
284 
285  ! 3D diagnostics hf_du(dv)_dt are commented because there is no clarity on proper remapping grid option.
286  ! The code is retained for degugging purposes in the future.
287  !if (CS%id_hf_du_dt > 0) then
288  ! do k=1,nz ; do j=js,je ; do I=Isq,Ieq
289  ! CS%hf_du_dt(I,j,k) = CS%du_dt(I,j,k) * ADp%diag_hfrac_u(I,j,k)
290  ! enddo ; enddo ; enddo
291  ! call post_data(CS%id_hf_du_dt, CS%hf_du_dt, CS%diag, alt_h = diag_pre_sync%h_state)
292  !endif
293 
294  !if (CS%id_hf_dv_dt > 0) then
295  ! do k=1,nz ; do J=Jsq,Jeq ; do i=is,ie
296  ! CS%hf_dv_dt(i,J,k) = CS%dv_dt(i,J,k) * ADp%diag_hfrac_v(i,J,k)
297  ! enddo ; enddo ; enddo
298  ! call post_data(CS%id_hf_dv_dt, CS%hf_dv_dt, CS%diag, alt_h = diag_pre_sync%h_state)
299  !endif
300 
301  if (cs%id_hf_du_dt_2d > 0) then
302  allocate(hf_du_dt_2d(g%IsdB:g%IedB,g%jsd:g%jed))
303  hf_du_dt_2d(:,:) = 0.0
304  do k=1,nz ; do j=js,je ; do i=isq,ieq
305  hf_du_dt_2d(i,j) = hf_du_dt_2d(i,j) + cs%du_dt(i,j,k) * adp%diag_hfrac_u(i,j,k)
306  enddo ; enddo ; enddo
307  call post_data(cs%id_hf_du_dt_2d, hf_du_dt_2d, cs%diag)
308  deallocate(hf_du_dt_2d)
309  endif
310 
311  if (cs%id_hf_dv_dt_2d > 0) then
312  allocate(hf_dv_dt_2d(g%isd:g%ied,g%JsdB:g%JedB))
313  hf_dv_dt_2d(:,:) = 0.0
314  do k=1,nz ; do j=jsq,jeq ; do i=is,ie
315  hf_dv_dt_2d(i,j) = hf_dv_dt_2d(i,j) + cs%dv_dt(i,j,k) * adp%diag_hfrac_v(i,j,k)
316  enddo ; enddo ; enddo
317  call post_data(cs%id_hf_dv_dt_2d, hf_dv_dt_2d, cs%diag)
318  deallocate(hf_dv_dt_2d)
319  endif
320 
321  call diag_restore_grids(cs%diag)
322 
323  call calculate_energy_diagnostics(u, v, h, uh, vh, adp, cdp, g, gv, us, cs)
324  endif
325 
326  ! smg: is the following robust to ALE? It seems a bit opaque.
327  ! If the model is NOT in isopycnal mode then nkmb=0. But we need all the
328  ! following diagnostics to treat all layers as variable density, so we set
329  ! nkmb = nz, on the expectation that loops nkmb+1,nz will not iterate.
330  ! This behavior is ANSI F77 but some compiler options can force at least
331  ! one iteration that would break the following one-line workaround!
332  if (nkmb==0 .and. nz > 1) nkmb = nz
333 
334  if (cs%id_u > 0) call post_data(cs%id_u, u, cs%diag)
335 
336  if (cs%id_v > 0) call post_data(cs%id_v, v, cs%diag)
337 
338  if (cs%id_h > 0) call post_data(cs%id_h, h, cs%diag)
339 
340  if (associated(cs%e)) then
341  call find_eta(h, tv, g, gv, us, cs%e, eta_bt)
342  if (cs%id_e > 0) call post_data(cs%id_e, cs%e, cs%diag)
343  endif
344 
345  if (associated(cs%e_D)) then
346  if (associated(cs%e)) then
347  do k=1,nz+1 ; do j=js,je ; do i=is,ie
348  cs%e_D(i,j,k) = cs%e(i,j,k) + g%bathyT(i,j)
349  enddo ; enddo ; enddo
350  else
351  call find_eta(h, tv, g, gv, us, cs%e_D, eta_bt)
352  do k=1,nz+1 ; do j=js,je ; do i=is,ie
353  cs%e_D(i,j,k) = cs%e_D(i,j,k) + g%bathyT(i,j)
354  enddo ; enddo ; enddo
355  endif
356 
357  if (cs%id_e_D > 0) call post_data(cs%id_e_D, cs%e_D, cs%diag)
358  endif
359 
360  ! mass per area of grid cell (for Bouss, use Rho0)
361  if (cs%id_masscello > 0) then
362  do k=1,nz; do j=js,je ; do i=is,ie
363  work_3d(i,j,k) = gv%H_to_kg_m2*h(i,j,k)
364  enddo ; enddo ; enddo
365  call post_data(cs%id_masscello, work_3d, cs%diag)
366  endif
367 
368  ! mass of liquid ocean (for Bouss, use Rho0). The reproducing sum requires the use of MKS units.
369  if (cs%id_masso > 0) then
370  work_2d(:,:) = 0.0
371  do k=1,nz ; do j=js,je ; do i=is,ie
372  work_2d(i,j) = work_2d(i,j) + (gv%H_to_kg_m2*h(i,j,k)) * us%L_to_m**2*g%areaT(i,j)
373  enddo ; enddo ; enddo
374  masso = reproducing_sum(work_2d)
375  call post_data(cs%id_masso, masso, cs%diag)
376  endif
377 
378  ! diagnose thickness/volumes of grid cells [m]
379  if (cs%id_thkcello>0 .or. cs%id_volcello>0) then
380  if (gv%Boussinesq) then ! thkcello = h for Boussinesq
381  if (cs%id_thkcello > 0) then ; if (gv%H_to_Z == 1.0) then
382  call post_data(cs%id_thkcello, h, cs%diag)
383  else
384  do k=1,nz; do j=js,je ; do i=is,ie
385  work_3d(i,j,k) = gv%H_to_Z*h(i,j,k)
386  enddo ; enddo ; enddo
387  call post_data(cs%id_thkcello, work_3d, cs%diag)
388  endif ; endif
389  if (cs%id_volcello > 0) then ! volcello = h*area for Boussinesq
390  do k=1,nz; do j=js,je ; do i=is,ie
391  work_3d(i,j,k) = ( gv%H_to_Z*h(i,j,k) ) * us%Z_to_m*us%L_to_m**2*g%areaT(i,j)
392  enddo ; enddo ; enddo
393  call post_data(cs%id_volcello, work_3d, cs%diag)
394  endif
395  else ! thkcello = dp/(rho*g) for non-Boussinesq
396  eosdom(:) = eos_domain(g%HI)
397  do j=js,je
398  if (associated(p_surf)) then ! Pressure loading at top of surface layer [R L2 T-2 ~> Pa]
399  do i=is,ie
400  pressure_1d(i) = p_surf(i,j)
401  enddo
402  else
403  do i=is,ie
404  pressure_1d(i) = 0.0
405  enddo
406  endif
407  do k=1,nz ! Integrate vertically downward for pressure
408  do i=is,ie ! Pressure for EOS at the layer center [R L2 T-2 ~> Pa]
409  pressure_1d(i) = pressure_1d(i) + 0.5*(gv%H_to_RZ*gv%g_Earth)*h(i,j,k)
410  enddo
411  ! Store in-situ density [R ~> kg m-3] in work_3d
412  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pressure_1d, rho_in_situ, &
413  tv%eqn_of_state, eosdom)
414  do i=is,ie ! Cell thickness = dz = dp/(g*rho) (meter); store in work_3d
415  work_3d(i,j,k) = (gv%H_to_RZ*h(i,j,k)) / rho_in_situ(i)
416  enddo
417  do i=is,ie ! Pressure for EOS at the bottom interface [R L2 T-2 ~> Pa]
418  pressure_1d(i) = pressure_1d(i) + 0.5*(gv%H_to_RZ*gv%g_Earth)*h(i,j,k)
419  enddo
420  enddo ! k
421  enddo ! j
422  if (cs%id_thkcello > 0) call post_data(cs%id_thkcello, work_3d, cs%diag)
423  if (cs%id_volcello > 0) then
424  do k=1,nz; do j=js,je ; do i=is,ie ! volcello = dp/(rho*g)*area for non-Boussinesq
425  work_3d(i,j,k) = us%Z_to_m*us%L_to_m**2*g%areaT(i,j) * work_3d(i,j,k)
426  enddo ; enddo ; enddo
427  call post_data(cs%id_volcello, work_3d, cs%diag)
428  endif
429  endif
430  endif
431 
432  ! Calculate additional, potentially derived temperature diagnostics
433  if (tv%T_is_conT) then
434  ! Internal T&S variables are conservative temperature & absolute salinity,
435  ! so they need to converted to potential temperature and practical salinity
436  ! for some diagnostics using TEOS-10 function calls.
437  if ((cs%id_Tpot > 0) .or. (cs%id_tob > 0)) then
438  do k=1,nz ; do j=js,je ; do i=is,ie
439  work_3d(i,j,k) = gsw_pt_from_ct(tv%S(i,j,k),tv%T(i,j,k))
440  enddo ; enddo ; enddo
441  if (cs%id_Tpot > 0) call post_data(cs%id_Tpot, work_3d, cs%diag)
442  if (cs%id_tob > 0) call post_data(cs%id_tob, work_3d(:,:,nz), cs%diag, mask=g%mask2dT)
443  endif
444  else
445  ! Internal T&S variables are potential temperature & practical salinity
446  if (cs%id_tob > 0) call post_data(cs%id_tob, tv%T(:,:,nz), cs%diag, mask=g%mask2dT)
447  endif
448 
449  ! Calculate additional, potentially derived salinity diagnostics
450  if (tv%S_is_absS) then
451  ! Internal T&S variables are conservative temperature & absolute salinity,
452  ! so they need to converted to potential temperature and practical salinity
453  ! for some diagnostics using TEOS-10 function calls.
454  if ((cs%id_Sprac > 0) .or. (cs%id_sob > 0)) then
455  do k=1,nz ; do j=js,je ; do i=is,ie
456  work_3d(i,j,k) = gsw_sp_from_sr(tv%S(i,j,k))
457  enddo ; enddo ; enddo
458  if (cs%id_Sprac > 0) call post_data(cs%id_Sprac, work_3d, cs%diag)
459  if (cs%id_sob > 0) call post_data(cs%id_sob, work_3d(:,:,nz), cs%diag, mask=g%mask2dT)
460  endif
461  else
462  ! Internal T&S variables are potential temperature & practical salinity
463  if (cs%id_sob > 0) call post_data(cs%id_sob, tv%S(:,:,nz), cs%diag, mask=g%mask2dT)
464  endif
465 
466  ! volume mean potential temperature
467  if (cs%id_thetaoga>0) then
468  thetaoga = global_volume_mean(tv%T, h, g, gv)
469  call post_data(cs%id_thetaoga, thetaoga, cs%diag)
470  endif
471 
472  ! area mean SST
473  if (cs%id_tosga > 0) then
474  do j=js,je ; do i=is,ie
475  surface_field(i,j) = tv%T(i,j,1)
476  enddo ; enddo
477  tosga = global_area_mean(surface_field, g)
478  call post_data(cs%id_tosga, tosga, cs%diag)
479  endif
480 
481  ! volume mean salinity
482  if (cs%id_soga>0) then
483  soga = global_volume_mean(tv%S, h, g, gv)
484  call post_data(cs%id_soga, soga, cs%diag)
485  endif
486 
487  ! area mean SSS
488  if (cs%id_sosga > 0) then
489  do j=js,je ; do i=is,ie
490  surface_field(i,j) = tv%S(i,j,1)
491  enddo ; enddo
492  sosga = global_area_mean(surface_field, g)
493  call post_data(cs%id_sosga, sosga, cs%diag)
494  endif
495 
496  ! layer mean potential temperature
497  if (cs%id_temp_layer_ave>0) then
498  temp_layer_ave = global_layer_mean(tv%T, h, g, gv)
499  call post_data(cs%id_temp_layer_ave, temp_layer_ave, cs%diag)
500  endif
501 
502  ! layer mean salinity
503  if (cs%id_salt_layer_ave>0) then
504  salt_layer_ave = global_layer_mean(tv%S, h, g, gv)
505  call post_data(cs%id_salt_layer_ave, salt_layer_ave, cs%diag)
506  endif
507 
508  call calculate_vertical_integrals(h, tv, p_surf, g, gv, us, cs)
509 
510  if ((cs%id_Rml > 0) .or. (cs%id_Rcv > 0) .or. associated(cs%h_Rlay) .or. &
511  associated(cs%uh_Rlay) .or. associated(cs%vh_Rlay) .or. &
512  associated(cs%uhGM_Rlay) .or. associated(cs%vhGM_Rlay)) then
513 
514  if (associated(tv%eqn_of_state)) then
515  eosdom(:) = eos_domain(g%HI, halo=1)
516  pressure_1d(:) = tv%P_Ref
517  !$OMP parallel do default(shared)
518  do k=1,nz ; do j=js-1,je+1
519  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pressure_1d, rcv(:,j,k), tv%eqn_of_state, &
520  eosdom)
521  enddo ; enddo
522  else ! Rcv should not be used much in this case, so fill in sensible values.
523  do k=1,nz ; do j=js-1,je+1 ; do i=is-1,ie+1
524  rcv(i,j,k) = gv%Rlay(k)
525  enddo ; enddo ; enddo
526  endif
527  if (cs%id_Rml > 0) call post_data(cs%id_Rml, rcv, cs%diag)
528  if (cs%id_Rcv > 0) call post_data(cs%id_Rcv, rcv, cs%diag)
529 
530  if (associated(cs%h_Rlay)) then
531  k_list = nz/2
532 !$OMP parallel do default(none) shared(is,ie,js,je,nz,nkmb,CS,Rcv,h,GV) &
533 !$OMP private(wt,wt_p) firstprivate(k_list)
534  do j=js,je
535  do k=1,nkmb ; do i=is,ie
536  cs%h_Rlay(i,j,k) = 0.0
537  enddo ; enddo
538  do k=nkmb+1,nz ; do i=is,ie
539  cs%h_Rlay(i,j,k) = h(i,j,k)
540  enddo ; enddo
541  do k=1,nkmb ; do i=is,ie
542  call find_weights(gv%Rlay, rcv(i,j,k), k_list, nz, wt, wt_p)
543  cs%h_Rlay(i,j,k_list) = cs%h_Rlay(i,j,k_list) + h(i,j,k)*wt
544  cs%h_Rlay(i,j,k_list+1) = cs%h_Rlay(i,j,k_list+1) + h(i,j,k)*wt_p
545  enddo ; enddo
546  enddo
547 
548  if (cs%id_h_Rlay > 0) call post_data(cs%id_h_Rlay, cs%h_Rlay, cs%diag)
549  endif
550 
551  if (associated(cs%uh_Rlay)) then
552  k_list = nz/2
553 !$OMP parallel do default(none) shared(Isq,Ieq,js,je,nz,nkmb,Rcv,CS,GV,uh) &
554 !$OMP private(wt,wt_p) firstprivate(k_list)
555  do j=js,je
556  do k=1,nkmb ; do i=isq,ieq
557  cs%uh_Rlay(i,j,k) = 0.0
558  enddo ; enddo
559  do k=nkmb+1,nz ; do i=isq,ieq
560  cs%uh_Rlay(i,j,k) = uh(i,j,k)
561  enddo ; enddo
562  k_list = nz/2
563  do k=1,nkmb ; do i=isq,ieq
564  call find_weights(gv%Rlay, 0.5*(rcv(i,j,k)+rcv(i+1,j,k)), k_list, nz, wt, wt_p)
565  cs%uh_Rlay(i,j,k_list) = cs%uh_Rlay(i,j,k_list) + uh(i,j,k)*wt
566  cs%uh_Rlay(i,j,k_list+1) = cs%uh_Rlay(i,j,k_list+1) + uh(i,j,k)*wt_p
567  enddo ; enddo
568  enddo
569 
570  if (cs%id_uh_Rlay > 0) call post_data(cs%id_uh_Rlay, cs%uh_Rlay, cs%diag)
571  endif
572 
573  if (associated(cs%vh_Rlay)) then
574  k_list = nz/2
575 !$OMP parallel do default(none) shared(Jsq,Jeq,is,ie,nz,nkmb,Rcv,CS,GV,vh) &
576 !$OMP private(wt,wt_p) firstprivate(k_list)
577  do j=jsq,jeq
578  do k=1,nkmb ; do i=is,ie
579  cs%vh_Rlay(i,j,k) = 0.0
580  enddo ; enddo
581  do k=nkmb+1,nz ; do i=is,ie
582  cs%vh_Rlay(i,j,k) = vh(i,j,k)
583  enddo ; enddo
584  do k=1,nkmb ; do i=is,ie
585  call find_weights(gv%Rlay, 0.5*(rcv(i,j,k)+rcv(i,j+1,k)), k_list, nz, wt, wt_p)
586  cs%vh_Rlay(i,j,k_list) = cs%vh_Rlay(i,j,k_list) + vh(i,j,k)*wt
587  cs%vh_Rlay(i,j,k_list+1) = cs%vh_Rlay(i,j,k_list+1) + vh(i,j,k)*wt_p
588  enddo ; enddo
589  enddo
590 
591  if (cs%id_vh_Rlay > 0) call post_data(cs%id_vh_Rlay, cs%vh_Rlay, cs%diag)
592  endif
593 
594  if (associated(cs%uhGM_Rlay) .and. associated(cdp%uhGM)) then
595  k_list = nz/2
596 !$OMP parallel do default(none) shared(Isq,Ieq,js,je,nz,nkmb,Rcv,CDP,CS,GV) &
597 !$OMP private(wt,wt_p) firstprivate(k_list)
598  do j=js,je
599  do k=1,nkmb ; do i=isq,ieq
600  cs%uhGM_Rlay(i,j,k) = 0.0
601  enddo ; enddo
602  do k=nkmb+1,nz ; do i=isq,ieq
603  cs%uhGM_Rlay(i,j,k) = cdp%uhGM(i,j,k)
604  enddo ; enddo
605  do k=1,nkmb ; do i=isq,ieq
606  call find_weights(gv%Rlay, 0.5*(rcv(i,j,k)+rcv(i+1,j,k)), k_list, nz, wt, wt_p)
607  cs%uhGM_Rlay(i,j,k_list) = cs%uhGM_Rlay(i,j,k_list) + cdp%uhGM(i,j,k)*wt
608  cs%uhGM_Rlay(i,j,k_list+1) = cs%uhGM_Rlay(i,j,k_list+1) + cdp%uhGM(i,j,k)*wt_p
609  enddo ; enddo
610  enddo
611 
612  if (cs%id_uhGM_Rlay > 0) call post_data(cs%id_uhGM_Rlay, cs%uhGM_Rlay, cs%diag)
613  endif
614 
615  if (associated(cs%vhGM_Rlay) .and. associated(cdp%vhGM)) then
616  k_list = nz/2
617 !$OMP parallel do default(none) shared(is,ie,Jsq,Jeq,nz,nkmb,CS,CDp,Rcv,GV) &
618 !$OMP private(wt,wt_p) firstprivate(k_list)
619  do j=jsq,jeq
620  do k=1,nkmb ; do i=is,ie
621  cs%vhGM_Rlay(i,j,k) = 0.0
622  enddo ; enddo
623  do k=nkmb+1,nz ; do i=is,ie
624  cs%vhGM_Rlay(i,j,k) = cdp%vhGM(i,j,k)
625  enddo ; enddo
626  do k=1,nkmb ; do i=is,ie
627  call find_weights(gv%Rlay, 0.5*(rcv(i,j,k)+rcv(i,j+1,k)), k_list, nz, wt, wt_p)
628  cs%vhGM_Rlay(i,j,k_list) = cs%vhGM_Rlay(i,j,k_list) + cdp%vhGM(i,j,k)*wt
629  cs%vhGM_Rlay(i,j,k_list+1) = cs%vhGM_Rlay(i,j,k_list+1) + cdp%vhGM(i,j,k)*wt_p
630  enddo ; enddo
631  enddo
632 
633  if (cs%id_vhGM_Rlay > 0) call post_data(cs%id_vhGM_Rlay, cs%vhGM_Rlay, cs%diag)
634  endif
635  endif
636 
637  if (associated(tv%eqn_of_state)) then
638  eosdom(:) = eos_domain(g%HI)
639  if (cs%id_rhopot0 > 0) then
640  pressure_1d(:) = 0.
641  !$OMP parallel do default(shared)
642  do k=1,nz ; do j=js,je
643  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pressure_1d, rcv(:,j,k), &
644  tv%eqn_of_state, eosdom)
645  enddo ; enddo
646  if (cs%id_rhopot0 > 0) call post_data(cs%id_rhopot0, rcv, cs%diag)
647  endif
648  if (cs%id_rhopot2 > 0) then
649  pressure_1d(:) = 2.0e7*us%kg_m3_to_R*us%m_s_to_L_T**2 ! 2000 dbars
650  !$OMP parallel do default(shared)
651  do k=1,nz ; do j=js,je
652  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pressure_1d, rcv(:,j,k), &
653  tv%eqn_of_state, eosdom)
654  enddo ; enddo
655  if (cs%id_rhopot2 > 0) call post_data(cs%id_rhopot2, rcv, cs%diag)
656  endif
657  if (cs%id_rhoinsitu > 0) then
658  !$OMP parallel do default(shared) private(pressure_1d)
659  do j=js,je
660  pressure_1d(:) = 0. ! Start at p=0 Pa at surface
661  do k=1,nz
662  pressure_1d(:) = pressure_1d(:) + 0.5 * h(:,j,k) * (gv%H_to_RZ*gv%g_Earth) ! Pressure in middle of layer k
663  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pressure_1d, rcv(:,j,k), &
664  tv%eqn_of_state, eosdom)
665  pressure_1d(:) = pressure_1d(:) + 0.5 * h(:,j,k) * (gv%H_to_RZ*gv%g_Earth) ! Pressure at bottom of layer k
666  enddo
667  enddo
668  if (cs%id_rhoinsitu > 0) call post_data(cs%id_rhoinsitu, rcv, cs%diag)
669  endif
670 
671  if (cs%id_drho_dT > 0 .or. cs%id_drho_dS > 0) then
672  !$OMP parallel do default(shared) private(pressure_1d)
673  do j=js,je
674  pressure_1d(:) = 0. ! Start at p=0 Pa at surface
675  do k=1,nz
676  pressure_1d(:) = pressure_1d(:) + 0.5 * h(:,j,k) * gv%H_to_Pa ! Pressure in middle of layer k
677  ! To avoid storing more arrays, put drho_dT into Rcv, and drho_dS into work3d
678  call calculate_density_derivs(tv%T(:,j,k),tv%S(:,j,k),pressure_1d, &
679  rcv(:,j,k),work_3d(:,j,k),is,ie-is+1, tv%eqn_of_state)
680  pressure_1d(:) = pressure_1d(:) + 0.5 * h(:,j,k) * gv%H_to_Pa ! Pressure at bottom of layer k
681  enddo
682  enddo
683  if (cs%id_drho_dT > 0) call post_data(cs%id_drho_dT, rcv, cs%diag)
684  if (cs%id_drho_dS > 0) call post_data(cs%id_drho_dS, work_3d, cs%diag)
685  endif
686  endif
687 
688  if ((cs%id_cg1>0) .or. (cs%id_Rd1>0) .or. (cs%id_cfl_cg1>0) .or. &
689  (cs%id_cfl_cg1_x>0) .or. (cs%id_cfl_cg1_y>0)) then
690  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp)
691  if (cs%id_cg1>0) call post_data(cs%id_cg1, cs%cg1, cs%diag)
692  if (cs%id_Rd1>0) then
693  !$OMP parallel do default(shared) private(f2_h,mag_beta)
694  do j=js,je ; do i=is,ie
695  ! Blend the equatorial deformation radius with the standard one.
696  f2_h = absurdly_small_freq2 + 0.25 * &
697  ((g%CoriolisBu(i,j)**2 + g%CoriolisBu(i-1,j-1)**2) + &
698  (g%CoriolisBu(i-1,j)**2 + g%CoriolisBu(i,j-1)**2))
699  mag_beta = sqrt(0.5 * ( &
700  (((g%CoriolisBu(i,j)-g%CoriolisBu(i-1,j)) * g%IdxCv(i,j))**2 + &
701  ((g%CoriolisBu(i,j-1)-g%CoriolisBu(i-1,j-1)) * g%IdxCv(i,j-1))**2) + &
702  (((g%CoriolisBu(i,j)-g%CoriolisBu(i,j-1)) * g%IdyCu(i,j))**2 + &
703  ((g%CoriolisBu(i-1,j)-g%CoriolisBu(i-1,j-1)) * g%IdyCu(i-1,j))**2) ))
704  cs%Rd1(i,j) = cs%cg1(i,j) / sqrt(f2_h + cs%cg1(i,j) * mag_beta)
705 
706  enddo ; enddo
707  call post_data(cs%id_Rd1, cs%Rd1, cs%diag)
708  endif
709  if (cs%id_cfl_cg1>0) then
710  do j=js,je ; do i=is,ie
711  cs%cfl_cg1(i,j) = (dt*cs%cg1(i,j)) * (g%IdxT(i,j) + g%IdyT(i,j))
712  enddo ; enddo
713  call post_data(cs%id_cfl_cg1, cs%cfl_cg1, cs%diag)
714  endif
715  if (cs%id_cfl_cg1_x>0) then
716  do j=js,je ; do i=is,ie
717  cs%cfl_cg1_x(i,j) = (dt*cs%cg1(i,j)) * g%IdxT(i,j)
718  enddo ; enddo
719  call post_data(cs%id_cfl_cg1_x, cs%cfl_cg1_x, cs%diag)
720  endif
721  if (cs%id_cfl_cg1_y>0) then
722  do j=js,je ; do i=is,ie
723  cs%cfl_cg1_y(i,j) = (dt*cs%cg1(i,j)) * g%IdyT(i,j)
724  enddo ; enddo
725  call post_data(cs%id_cfl_cg1_y, cs%cfl_cg1_y, cs%diag)
726  endif
727  endif
728  if ((cs%id_cg_ebt>0) .or. (cs%id_Rd_ebt>0) .or. (cs%id_p_ebt>0)) then
729  if (cs%id_p_ebt>0) then
730  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp, use_ebt_mode=.true., &
731  mono_n2_column_fraction=cs%mono_N2_column_fraction, &
732  mono_n2_depth=cs%mono_N2_depth, modal_structure=cs%p_ebt)
733  call post_data(cs%id_p_ebt, cs%p_ebt, cs%diag)
734  else
735  call wave_speed(h, tv, g, gv, us, cs%cg1, cs%wave_speed_CSp, use_ebt_mode=.true., &
736  mono_n2_column_fraction=cs%mono_N2_column_fraction, &
737  mono_n2_depth=cs%mono_N2_depth)
738  endif
739  if (cs%id_cg_ebt>0) call post_data(cs%id_cg_ebt, cs%cg1, cs%diag)
740  if (cs%id_Rd_ebt>0) then
741  !$OMP parallel do default(shared) private(f2_h,mag_beta)
742  do j=js,je ; do i=is,ie
743  ! Blend the equatorial deformation radius with the standard one.
744  f2_h = absurdly_small_freq2 + 0.25 * &
745  ((g%CoriolisBu(i,j)**2 + g%CoriolisBu(i-1,j-1)**2) + &
746  (g%CoriolisBu(i-1,j)**2 + g%CoriolisBu(i,j-1)**2))
747  mag_beta = sqrt(0.5 * ( &
748  (((g%CoriolisBu(i,j)-g%CoriolisBu(i-1,j)) * g%IdxCv(i,j))**2 + &
749  ((g%CoriolisBu(i,j-1)-g%CoriolisBu(i-1,j-1)) * g%IdxCv(i,j-1))**2) + &
750  (((g%CoriolisBu(i,j)-g%CoriolisBu(i,j-1)) * g%IdyCu(i,j))**2 + &
751  ((g%CoriolisBu(i-1,j)-g%CoriolisBu(i-1,j-1)) * g%IdyCu(i-1,j))**2) ))
752  cs%Rd1(i,j) = cs%cg1(i,j) / sqrt(f2_h + cs%cg1(i,j) * mag_beta)
753 
754  enddo ; enddo
755  call post_data(cs%id_Rd_ebt, cs%Rd1, cs%diag)
756  endif
757  endif
758 
759 end subroutine calculate_diagnostic_fields
760 
761 !> This subroutine finds the location of R_in in an increasing ordered
762 !! list, Rlist, returning as k the element such that
763 !! Rlist(k) <= R_in < Rlist(k+1), and where wt and wt_p are the linear
764 !! weights that should be assigned to elements k and k+1.
765 subroutine find_weights(Rlist, R_in, k, nz, wt, wt_p)
766  real, dimension(:), &
767  intent(in) :: Rlist !< The list of target densities [R ~> kg m-3]
768  real, intent(in) :: R_in !< The density being inserted into Rlist [R ~> kg m-3]
769  integer, intent(inout) :: k !< The value of k such that Rlist(k) <= R_in < Rlist(k+1)
770  !! The input value is a first guess
771  integer, intent(in) :: nz !< The number of layers in Rlist
772  real, intent(out) :: wt !< The weight of layer k for interpolation, nondim
773  real, intent(out) :: wt_p !< The weight of layer k+1 for interpolation, nondim
774 
775  ! This subroutine finds location of R_in in an increasing ordered
776  ! list, Rlist, returning as k the element such that
777  ! Rlist(k) <= R_in < Rlist(k+1), and where wt and wt_p are the linear
778  ! weights that should be assigned to elements k and k+1.
779 
780  integer :: k_upper, k_lower, k_new, inc
781 
782  ! First, bracket the desired point.
783  if ((k < 1) .or. (k > nz)) k = nz/2
784 
785  k_upper = k ; k_lower = k ; inc = 1
786  if (r_in < rlist(k)) then
787  do
788  k_lower = max(k_lower-inc, 1)
789  if ((k_lower == 1) .or. (r_in >= rlist(k_lower))) exit
790  k_upper = k_lower
791  inc = inc*2
792  enddo
793  else
794  do
795  k_upper = min(k_upper+inc, nz)
796  if ((k_upper == nz) .or. (r_in < rlist(k_upper))) exit
797  k_lower = k_upper
798  inc = inc*2
799  enddo
800  endif
801 
802  if ((k_lower == 1) .and. (r_in <= rlist(k_lower))) then
803  k = 1 ; wt = 1.0 ; wt_p = 0.0
804  elseif ((k_upper == nz) .and. (r_in >= rlist(k_upper))) then
805  k = nz-1 ; wt = 0.0 ; wt_p = 1.0
806  else
807  do
808  if (k_upper <= k_lower+1) exit
809  k_new = (k_upper + k_lower) / 2
810  if (r_in < rlist(k_new)) then
811  k_upper = k_new
812  else
813  k_lower = k_new
814  endif
815  enddo
816 
817 ! Uncomment this as a code check
818 ! if ((R_in < Rlist(k_lower)) .or. (R_in >= Rlist(k_upper)) .or. (k_upper-k_lower /= 1)) &
819 ! write (*,*) "Error: ",R_in," is not between R(",k_lower,") = ", &
820 ! Rlist(k_lower)," and R(",k_upper,") = ",Rlist(k_upper),"."
821  k = k_lower
822  wt = (rlist(k_upper) - r_in) / (rlist(k_upper) - rlist(k_lower))
823  wt_p = 1.0 - wt
824 
825  endif
826 
827 end subroutine find_weights
828 
829 !> This subroutine calculates vertical integrals of several tracers, along
830 !! with the mass-weight of these tracers, the total column mass, and the
831 !! carefully calculated column height.
832 subroutine calculate_vertical_integrals(h, tv, p_surf, G, GV, US, CS)
833  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
834  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
835  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
836  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
837  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
838  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
839  !! thermodynamic variables.
840  real, dimension(:,:), pointer :: p_surf !< A pointer to the surface pressure [R L2 T-2 ~> Pa].
841  !! If p_surf is not associated, it is the same
842  !! as setting the surface pressure to 0.
843  type(diagnostics_CS), intent(inout) :: CS !< Control structure returned by a
844  !! previous call to diagnostics_init.
845 
846  real, dimension(SZI_(G), SZJ_(G)) :: &
847  z_top, & ! Height of the top of a layer or the ocean [Z ~> m].
848  z_bot, & ! Height of the bottom of a layer (for id_mass) or the
849  ! (positive) depth of the ocean (for id_col_ht) [Z ~> m].
850  mass, & ! integrated mass of the water column [R Z ~> kg m-2]. For
851  ! non-Boussinesq models this is rho*dz. For Boussinesq
852  ! models, this is either the integral of in-situ density
853  ! (rho*dz for col_mass) or reference density (Rho_0*dz for mass_wt).
854  btm_pres,&! The pressure at the ocean bottom, or CMIP variable 'pbo'.
855  ! This is the column mass multiplied by gravity plus the pressure
856  ! at the ocean surface [R L2 T-2 ~> Pa].
857  dpress, & ! Change in hydrostatic pressure across a layer [R L2 T-2 ~> Pa].
858  tr_int ! vertical integral of a tracer times density,
859  ! (Rho_0 in a Boussinesq model) [TR kg m-2].
860  real :: IG_Earth ! Inverse of gravitational acceleration [T2 Z L-2 ~> s2 m-1].
861 
862  integer :: i, j, k, is, ie, js, je, nz
863  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
864 
865  if (cs%id_mass_wt > 0) then
866  do j=js,je ; do i=is,ie ; mass(i,j) = 0.0 ; enddo ; enddo
867  do k=1,nz ; do j=js,je ; do i=is,ie
868  mass(i,j) = mass(i,j) + gv%H_to_RZ*h(i,j,k)
869  enddo ; enddo ; enddo
870  call post_data(cs%id_mass_wt, mass, cs%diag)
871  endif
872 
873  if (cs%id_temp_int > 0) then
874  do j=js,je ; do i=is,ie ; tr_int(i,j) = 0.0 ; enddo ; enddo
875  do k=1,nz ; do j=js,je ; do i=is,ie
876  tr_int(i,j) = tr_int(i,j) + (gv%H_to_RZ*h(i,j,k))*tv%T(i,j,k)
877  enddo ; enddo ; enddo
878  call post_data(cs%id_temp_int, tr_int, cs%diag)
879  endif
880 
881  if (cs%id_salt_int > 0) then
882  do j=js,je ; do i=is,ie ; tr_int(i,j) = 0.0 ; enddo ; enddo
883  do k=1,nz ; do j=js,je ; do i=is,ie
884  tr_int(i,j) = tr_int(i,j) + (gv%H_to_RZ*h(i,j,k))*tv%S(i,j,k)
885  enddo ; enddo ; enddo
886  call post_data(cs%id_salt_int, tr_int, cs%diag)
887  endif
888 
889  if (cs%id_col_ht > 0) then
890  call find_eta(h, tv, g, gv, us, z_top)
891  do j=js,je ; do i=is,ie
892  z_bot(i,j) = z_top(i,j) + g%bathyT(i,j)
893  enddo ; enddo
894  call post_data(cs%id_col_ht, z_bot, cs%diag)
895  endif
896 
897  ! NOTE: int_density_z expects z_top and z_btm values from [ij]sq to [ij]eq+1
898  if (cs%id_col_mass > 0 .or. cs%id_pbo > 0) then
899  do j=js,je ; do i=is,ie ; mass(i,j) = 0.0 ; enddo ; enddo
900  if (gv%Boussinesq) then
901  if (associated(tv%eqn_of_state)) then
902  ig_earth = 1.0 / gv%g_Earth
903 ! do j=js,je ; do i=is,ie ; z_bot(i,j) = -P_SURF(i,j)/GV%H_to_Pa ; enddo ; enddo
904  do j=g%jscB,g%jecB+1 ; do i=g%iscB,g%iecB+1
905  z_bot(i,j) = 0.0
906  enddo ; enddo
907  do k=1,nz
908  do j=g%jscB,g%jecB+1 ; do i=g%iscB,g%iecB+1
909  z_top(i,j) = z_bot(i,j)
910  z_bot(i,j) = z_top(i,j) - gv%H_to_Z*h(i,j,k)
911  enddo ; enddo
912  call int_density_dz(tv%T(:,:,k), tv%S(:,:,k), z_top, z_bot, 0.0, gv%Rho0, gv%g_Earth, &
913  g%HI, tv%eqn_of_state, us, dpress)
914  do j=js,je ; do i=is,ie
915  mass(i,j) = mass(i,j) + dpress(i,j) * ig_earth
916  enddo ; enddo
917  enddo
918  else
919  do k=1,nz ; do j=js,je ; do i=is,ie
920  mass(i,j) = mass(i,j) + (gv%H_to_Z*gv%Rlay(k))*h(i,j,k)
921  enddo ; enddo ; enddo
922  endif
923  else
924  do k=1,nz ; do j=js,je ; do i=is,ie
925  mass(i,j) = mass(i,j) + gv%H_to_RZ*h(i,j,k)
926  enddo ; enddo ; enddo
927  endif
928  if (cs%id_col_mass > 0) then
929  call post_data(cs%id_col_mass, mass, cs%diag)
930  endif
931  if (cs%id_pbo > 0) then
932  do j=js,je ; do i=is,ie ; btm_pres(i,j) = 0.0 ; enddo ; enddo
933  ! 'pbo' is defined as the sea water pressure at the sea floor
934  ! pbo = (mass * g) + p_surf
935  ! where p_surf is the sea water pressure at sea water surface.
936  do j=js,je ; do i=is,ie
937  btm_pres(i,j) = gv%g_Earth * mass(i,j)
938  if (associated(p_surf)) then
939  btm_pres(i,j) = btm_pres(i,j) + p_surf(i,j)
940  endif
941  enddo ; enddo
942  call post_data(cs%id_pbo, btm_pres, cs%diag)
943  endif
944  endif
945 
946 end subroutine calculate_vertical_integrals
947 
948 !> This subroutine calculates terms in the mechanical energy budget.
949 subroutine calculate_energy_diagnostics(u, v, h, uh, vh, ADp, CDp, G, GV, US, CS)
950  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
951  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
952  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
953  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
954  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
955  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
956  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
957  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
958  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
959  intent(in) :: uh !< Transport through zonal faces=u*h*dy,
960  !! [H L2 T-1 ~> m3 s-1 or kg s-1].
961  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
962  intent(in) :: vh !< Transport through merid faces=v*h*dx,
963  !! [H L2 T-1 ~> m3 s-1 or kg s-1].
964  type(accel_diag_ptrs), intent(in) :: ADp !< Structure pointing to accelerations in momentum equation.
965  type(cont_diag_ptrs), intent(in) :: CDp !< Structure pointing to terms in continuity equations.
966  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
967  type(diagnostics_CS), intent(inout) :: CS !< Control structure returned by a previous call to
968  !! diagnostics_init.
969 
970 ! This subroutine calculates terms in the mechanical energy budget.
971 
972  ! Local variables
973  real :: KE_u(SZIB_(G),SZJ_(G))
974  real :: KE_v(SZI_(G),SZJB_(G))
975  real :: KE_h(SZI_(G),SZJ_(G))
976 
977  integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
978  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
979  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
980 
981  do j=js-1,je ; do i=is-1,ie
982  ke_u(i,j) = 0.0 ; ke_v(i,j) = 0.0
983  enddo ; enddo
984 
985  if (associated(cs%KE)) then
986  do k=1,nz ; do j=js,je ; do i=is,ie
987  cs%KE(i,j,k) = ((u(i,j,k) * u(i,j,k) + u(i-1,j,k) * u(i-1,j,k)) &
988  + (v(i,j,k) * v(i,j,k) + v(i,j-1,k) * v(i,j-1,k))) * 0.25
989  ! DELETE THE FOLLOWING... Make this 0 to test the momentum balance,
990  ! or a huge number to test the continuity balance.
991  ! CS%KE(i,j,k) *= 1e20
992  enddo ; enddo ; enddo
993  if (cs%id_KE > 0) call post_data(cs%id_KE, cs%KE, cs%diag)
994  endif
995 
996  if (.not.g%symmetric) then
997  if (associated(cs%dKE_dt) .OR. associated(cs%PE_to_KE) .OR. associated(cs%KE_CorAdv) .OR. &
998  associated(cs%KE_adv) .OR. associated(cs%KE_visc) .OR. associated(cs%KE_horvisc).OR. &
999  associated(cs%KE_dia) ) then
1000  call create_group_pass(cs%pass_KE_uv, ke_u, ke_v, g%Domain, to_north+to_east)
1001  endif
1002  endif
1003 
1004  if (associated(cs%dKE_dt)) then
1005  do k=1,nz
1006  do j=js,je ; do i=isq,ieq
1007  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * cs%du_dt(i,j,k)
1008  enddo ; enddo
1009  do j=jsq,jeq ; do i=is,ie
1010  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * cs%dv_dt(i,j,k)
1011  enddo ; enddo
1012  do j=js,je ; do i=is,ie
1013  ke_h(i,j) = cs%KE(i,j,k) * cs%dh_dt(i,j,k)
1014  enddo ; enddo
1015  if (.not.g%symmetric) &
1016  call do_group_pass(cs%pass_KE_uv, g%domain)
1017  do j=js,je ; do i=is,ie
1018  cs%dKE_dt(i,j,k) = ke_h(i,j) + 0.5 * g%IareaT(i,j) &
1019  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1020  enddo ; enddo
1021  enddo
1022  if (cs%id_dKEdt > 0) call post_data(cs%id_dKEdt, cs%dKE_dt, cs%diag)
1023  endif
1024 
1025  if (associated(cs%PE_to_KE)) then
1026  do k=1,nz
1027  do j=js,je ; do i=isq,ieq
1028  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * adp%PFu(i,j,k)
1029  enddo ; enddo
1030  do j=jsq,jeq ; do i=is,ie
1031  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * adp%PFv(i,j,k)
1032  enddo ; enddo
1033  if (.not.g%symmetric) &
1034  call do_group_pass(cs%pass_KE_uv, g%domain)
1035  do j=js,je ; do i=is,ie
1036  cs%PE_to_KE(i,j,k) = 0.5 * g%IareaT(i,j) &
1037  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1038  enddo ; enddo
1039  enddo
1040  if (cs%id_PE_to_KE > 0) call post_data(cs%id_PE_to_KE, cs%PE_to_KE, cs%diag)
1041  endif
1042 
1043  if (associated(cs%KE_CorAdv)) then
1044  do k=1,nz
1045  do j=js,je ; do i=isq,ieq
1046  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * adp%CAu(i,j,k)
1047  enddo ; enddo
1048  do j=jsq,jeq ; do i=is,ie
1049  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * adp%CAv(i,j,k)
1050  enddo ; enddo
1051  do j=js,je ; do i=is,ie
1052  ke_h(i,j) = -cs%KE(i,j,k) * g%IareaT(i,j) &
1053  * (uh(i,j,k) - uh(i-1,j,k) + vh(i,j,k) - vh(i,j-1,k))
1054  enddo ; enddo
1055  if (.not.g%symmetric) &
1056  call do_group_pass(cs%pass_KE_uv, g%domain)
1057  do j=js,je ; do i=is,ie
1058  cs%KE_CorAdv(i,j,k) = ke_h(i,j) + 0.5 * g%IareaT(i,j) &
1059  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1060  enddo ; enddo
1061  enddo
1062  if (cs%id_KE_Coradv > 0) call post_data(cs%id_KE_Coradv, cs%KE_Coradv, cs%diag)
1063  endif
1064 
1065  if (associated(cs%KE_adv)) then
1066  ! NOTE: All terms in KE_adv are multipled by -1, which can easily produce
1067  ! negative zeros and may signal a reproducibility issue over land.
1068  ! We resolve this by re-initializing and only evaluating over water points.
1069  ke_u(:,:) = 0. ; ke_v(:,:) = 0.
1070  do k=1,nz
1071  do j=js,je ; do i=isq,ieq
1072  if (g%mask2dCu(i,j) /= 0.) &
1073  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * adp%gradKEu(i,j,k)
1074  enddo ; enddo
1075  do j=jsq,jeq ; do i=is,ie
1076  if (g%mask2dCv(i,j) /= 0.) &
1077  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * adp%gradKEv(i,j,k)
1078  enddo ; enddo
1079  do j=js,je ; do i=is,ie
1080  ke_h(i,j) = -cs%KE(i,j,k) * g%IareaT(i,j) &
1081  * (uh(i,j,k) - uh(i-1,j,k) + vh(i,j,k) - vh(i,j-1,k))
1082  enddo ; enddo
1083  if (.not.g%symmetric) &
1084  call do_group_pass(cs%pass_KE_uv, g%domain)
1085  do j=js,je ; do i=is,ie
1086  cs%KE_adv(i,j,k) = ke_h(i,j) + 0.5 * g%IareaT(i,j) &
1087  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1088  enddo ; enddo
1089  enddo
1090  if (cs%id_KE_adv > 0) call post_data(cs%id_KE_adv, cs%KE_adv, cs%diag)
1091  endif
1092 
1093  if (associated(cs%KE_visc)) then
1094  do k=1,nz
1095  do j=js,je ; do i=isq,ieq
1096  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * adp%du_dt_visc(i,j,k)
1097  enddo ; enddo
1098  do j=jsq,jeq ; do i=is,ie
1099  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * adp%dv_dt_visc(i,j,k)
1100  enddo ; enddo
1101  if (.not.g%symmetric) &
1102  call do_group_pass(cs%pass_KE_uv, g%domain)
1103  do j=js,je ; do i=is,ie
1104  cs%KE_visc(i,j,k) = 0.5 * g%IareaT(i,j) &
1105  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1106  enddo ; enddo
1107  enddo
1108  if (cs%id_KE_visc > 0) call post_data(cs%id_KE_visc, cs%KE_visc, cs%diag)
1109  endif
1110 
1111  if (associated(cs%KE_horvisc)) then
1112  do k=1,nz
1113  do j=js,je ; do i=isq,ieq
1114  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * adp%diffu(i,j,k)
1115  enddo ; enddo
1116  do j=jsq,jeq ; do i=is,ie
1117  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * adp%diffv(i,j,k)
1118  enddo ; enddo
1119  if (.not.g%symmetric) &
1120  call do_group_pass(cs%pass_KE_uv, g%domain)
1121  do j=js,je ; do i=is,ie
1122  cs%KE_horvisc(i,j,k) = 0.5 * g%IareaT(i,j) &
1123  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1124  enddo ; enddo
1125  enddo
1126  if (cs%id_KE_horvisc > 0) call post_data(cs%id_KE_horvisc, cs%KE_horvisc, cs%diag)
1127  endif
1128 
1129  if (associated(cs%KE_dia)) then
1130  do k=1,nz
1131  do j=js,je ; do i=isq,ieq
1132  ke_u(i,j) = uh(i,j,k) * g%dxCu(i,j) * adp%du_dt_dia(i,j,k)
1133  enddo ; enddo
1134  do j=jsq,jeq ; do i=is,ie
1135  ke_v(i,j) = vh(i,j,k) * g%dyCv(i,j) * adp%dv_dt_dia(i,j,k)
1136  enddo ; enddo
1137  do j=js,je ; do i=is,ie
1138  ke_h(i,j) = cs%KE(i,j,k) &
1139  * (us%T_to_s * (cdp%diapyc_vel(i,j,k) - cdp%diapyc_vel(i,j,k+1)))
1140  enddo ; enddo
1141  if (.not.g%symmetric) &
1142  call do_group_pass(cs%pass_KE_uv, g%domain)
1143  do j=js,je ; do i=is,ie
1144  cs%KE_dia(i,j,k) = ke_h(i,j) + 0.5 * g%IareaT(i,j) &
1145  * (ke_u(i,j) + ke_u(i-1,j) + ke_v(i,j) + ke_v(i,j-1))
1146  enddo ; enddo
1147  enddo
1148  if (cs%id_KE_dia > 0) call post_data(cs%id_KE_dia, cs%KE_dia, cs%diag)
1149  endif
1150 
1151 end subroutine calculate_energy_diagnostics
1152 
1153 !> This subroutine registers fields to calculate a diagnostic time derivative.
1154 subroutine register_time_deriv(lb, f_ptr, deriv_ptr, CS)
1155  integer, intent(in), dimension(3) :: lb !< Lower index bound of f_ptr
1156  real, dimension(lb(1):,lb(2):,:), target :: f_ptr
1157  !< Time derivative operand
1158  real, dimension(lb(1):,lb(2):,:), target :: deriv_ptr
1159  !< Time derivative of f_ptr
1160  type(diagnostics_cs), pointer :: cs !< Control structure returned by previous call to
1161  !! diagnostics_init.
1162 
1163  ! This subroutine registers fields to calculate a diagnostic time derivative.
1164  ! NOTE: Lower bound is required for grid indexing in calculate_derivs().
1165  ! We assume that the vertical axis is 1-indexed.
1166 
1167  integer :: m !< New index of deriv_ptr in CS%deriv
1168  integer :: ub(3) !< Upper index bound of f_ptr, based on shape.
1169 
1170  if (.not.associated(cs)) call mom_error(fatal, &
1171  "register_time_deriv: Module must be initialized before it is used.")
1172 
1173  if (cs%num_time_deriv >= max_fields_) then
1174  call mom_error(warning,"MOM_diagnostics: Attempted to register more than " // &
1175  "MAX_FIELDS_ diagnostic time derivatives via register_time_deriv.")
1176  return
1177  endif
1178 
1179  m = cs%num_time_deriv+1 ; cs%num_time_deriv = m
1180 
1181  ub(:) = lb(:) + shape(f_ptr) - 1
1182 
1183  cs%nlay(m) = size(f_ptr, 3)
1184  cs%deriv(m)%p => deriv_ptr
1185  allocate(cs%prev_val(m)%p(lb(1):ub(1), lb(2):ub(2), cs%nlay(m)))
1186 
1187  cs%var_ptr(m)%p => f_ptr
1188  cs%prev_val(m)%p(:,:,:) = f_ptr(:,:,:)
1189 
1190 end subroutine register_time_deriv
1191 
1192 !> This subroutine calculates all registered time derivatives.
1193 subroutine calculate_derivs(dt, G, CS)
1194  real, intent(in) :: dt !< The time interval over which differences occur [T ~> s].
1195  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
1196  type(diagnostics_CS), intent(inout) :: CS !< Control structure returned by previous call to
1197  !! diagnostics_init.
1198 
1199 ! This subroutine calculates all registered time derivatives.
1200  real :: Idt ! The inverse timestep [T-1 ~> s-1]
1201  integer :: i, j, k, m
1202 
1203  if (dt > 0.0) then ; idt = 1.0/dt
1204  else ; return ; endif
1205 
1206  ! Because the field is unknown, its grid index bounds are also unknown.
1207  ! Additionally, two of the fields (dudt, dvdt) require calculation of spatial
1208  ! derivatives when computing d(KE)/dt. This raises issues in non-symmetric
1209  ! mode, where the symmetric boundaries (west, south) may not be updated.
1210 
1211  ! For this reason, we explicitly loop from isc-1:iec and jsc-1:jec, in order
1212  ! to force boundary value updates, even though it may not be strictly valid
1213  ! for all fields. Note this assumes a halo, and that it has been updated.
1214 
1215  do m=1,cs%num_time_deriv
1216  do k=1,cs%nlay(m) ; do j=g%jsc-1,g%jec ; do i=g%isc-1,g%iec
1217  cs%deriv(m)%p(i,j,k) = (cs%var_ptr(m)%p(i,j,k) - cs%prev_val(m)%p(i,j,k)) * idt
1218  cs%prev_val(m)%p(i,j,k) = cs%var_ptr(m)%p(i,j,k)
1219  enddo ; enddo ; enddo
1220  enddo
1221 
1222 end subroutine calculate_derivs
1223 
1224 !> This routine posts diagnostics of various dynamic ocean surface quantities,
1225 !! including velocities, speed and sea surface height, at the time the ocean
1226 !! state is reported back to the caller
1227 subroutine post_surface_dyn_diags(IDs, G, diag, sfc_state, ssh)
1228  type(surface_diag_ids), intent(in) :: ids !< A structure with the diagnostic IDs.
1229  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
1230  type(diag_ctrl), intent(in) :: diag !< regulates diagnostic output
1231  type(surface), intent(in) :: sfc_state !< structure describing the ocean surface state
1232  real, dimension(SZI_(G),SZJ_(G)), &
1233  intent(in) :: ssh !< Time mean surface height without corrections for ice displacement [m]
1234 
1235  ! Local variables
1236  real, dimension(SZI_(G),SZJ_(G)) :: speed ! The surface speed [L T-1 ~> m s-1]
1237  integer :: i, j, is, ie, js, je
1238 
1239  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1240 
1241  if (ids%id_ssh > 0) &
1242  call post_data(ids%id_ssh, ssh, diag, mask=g%mask2dT)
1243 
1244  if (ids%id_ssu > 0) &
1245  call post_data(ids%id_ssu, sfc_state%u, diag, mask=g%mask2dCu)
1246 
1247  if (ids%id_ssv > 0) &
1248  call post_data(ids%id_ssv, sfc_state%v, diag, mask=g%mask2dCv)
1249 
1250  if (ids%id_speed > 0) then
1251  do j=js,je ; do i=is,ie
1252  speed(i,j) = sqrt(0.5*(sfc_state%u(i-1,j)**2 + sfc_state%u(i,j)**2) + &
1253  0.5*(sfc_state%v(i,j-1)**2 + sfc_state%v(i,j)**2))
1254  enddo ; enddo
1255  call post_data(ids%id_speed, speed, diag, mask=g%mask2dT)
1256  endif
1257 
1258 end subroutine post_surface_dyn_diags
1259 
1260 
1261 !> This routine posts diagnostics of various ocean surface and integrated
1262 !! quantities at the time the ocean state is reported back to the caller
1263 subroutine post_surface_thermo_diags(IDs, G, GV, US, diag, dt_int, sfc_state, tv, &
1264  ssh, ssh_ibc)
1265  type(surface_diag_ids), intent(in) :: ids !< A structure with the diagnostic IDs.
1266  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
1267  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
1268  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1269  type(diag_ctrl), intent(in) :: diag !< regulates diagnostic output
1270  real, intent(in) :: dt_int !< total time step associated with these diagnostics [T ~> s].
1271  type(surface), intent(in) :: sfc_state !< structure describing the ocean surface state
1272  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
1273  real, dimension(SZI_(G),SZJ_(G)), &
1274  intent(in) :: ssh !< Time mean surface height without corrections for ice displacement [m]
1275  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: ssh_ibc !< Time mean surface height with corrections
1276  !! for ice displacement and the inverse barometer [m]
1277 
1278  real, dimension(SZI_(G),SZJ_(G)) :: work_2d ! A 2-d work array
1279  real, dimension(SZI_(G),SZJ_(G)) :: &
1280  zos ! dynamic sea lev (zero area mean) from inverse-barometer adjusted ssh [m]
1281  real :: i_time_int ! The inverse of the time interval [T-1 ~> s-1].
1282  real :: zos_area_mean, volo, ssh_ga
1283  integer :: i, j, is, ie, js, je
1284 
1285  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1286 
1287  ! area mean SSH
1288  if (ids%id_ssh_ga > 0) then
1289  ssh_ga = global_area_mean(ssh, g)
1290  call post_data(ids%id_ssh_ga, ssh_ga, diag)
1291  endif
1292 
1293  ! post the dynamic sea level, zos, and zossq.
1294  ! zos is ave_ssh with sea ice inverse barometer removed,
1295  ! and with zero global area mean.
1296  if (ids%id_zos > 0 .or. ids%id_zossq > 0) then
1297  zos(:,:) = 0.0
1298  do j=js,je ; do i=is,ie
1299  zos(i,j) = ssh_ibc(i,j)
1300  enddo ; enddo
1301  zos_area_mean = global_area_mean(zos, g)
1302  do j=js,je ; do i=is,ie
1303  zos(i,j) = zos(i,j) - g%mask2dT(i,j)*zos_area_mean
1304  enddo ; enddo
1305  if (ids%id_zos > 0) call post_data(ids%id_zos, zos, diag, mask=g%mask2dT)
1306  if (ids%id_zossq > 0) then
1307  do j=js,je ; do i=is,ie
1308  work_2d(i,j) = zos(i,j)*zos(i,j)
1309  enddo ; enddo
1310  call post_data(ids%id_zossq, work_2d, diag, mask=g%mask2dT)
1311  endif
1312  endif
1313 
1314  ! post total volume of the liquid ocean
1315  if (ids%id_volo > 0) then
1316  do j=js,je ; do i=is,ie
1317  work_2d(i,j) = g%mask2dT(i,j)*(ssh(i,j) + us%Z_to_m*g%bathyT(i,j))
1318  enddo ; enddo
1319  volo = global_area_integral(work_2d, g)
1320  call post_data(ids%id_volo, volo, diag)
1321  endif
1322 
1323  ! Use Adcroft's rule of reciprocals; it does the right thing here.
1324  i_time_int = 0.0 ; if (dt_int > 0.0) i_time_int = 1.0 / dt_int
1325 
1326  ! post time-averaged rate of frazil formation
1327  if (associated(tv%frazil) .and. (ids%id_fraz > 0)) then
1328  do j=js,je ; do i=is,ie
1329  work_2d(i,j) = tv%frazil(i,j) * i_time_int
1330  enddo ; enddo
1331  call post_data(ids%id_fraz, work_2d, diag, mask=g%mask2dT)
1332  endif
1333 
1334  ! post time-averaged salt deficit
1335  if (associated(tv%salt_deficit) .and. (ids%id_salt_deficit > 0)) then
1336  do j=js,je ; do i=is,ie
1337  work_2d(i,j) = tv%salt_deficit(i,j) * i_time_int
1338  enddo ; enddo
1339  call post_data(ids%id_salt_deficit, work_2d, diag, mask=g%mask2dT)
1340  endif
1341 
1342  ! post temperature of P-E+R
1343  if (associated(tv%TempxPmE) .and. (ids%id_Heat_PmE > 0)) then
1344  do j=js,je ; do i=is,ie
1345  work_2d(i,j) = tv%TempxPmE(i,j) * (tv%C_p * i_time_int)
1346  enddo ; enddo
1347  call post_data(ids%id_Heat_PmE, work_2d, diag, mask=g%mask2dT)
1348  endif
1349 
1350  ! post geothermal heating or internal heat source/sinks
1351  if (associated(tv%internal_heat) .and. (ids%id_intern_heat > 0)) then
1352  do j=js,je ; do i=is,ie
1353  work_2d(i,j) = tv%internal_heat(i,j) * (tv%C_p * i_time_int)
1354  enddo ; enddo
1355  call post_data(ids%id_intern_heat, work_2d, diag, mask=g%mask2dT)
1356  endif
1357 
1358  if (tv%T_is_conT) then
1359  ! Internal T&S variables are conservative temperature & absolute salinity
1360  if (ids%id_sstcon > 0) call post_data(ids%id_sstcon, sfc_state%SST, diag, mask=g%mask2dT)
1361  ! Use TEOS-10 function calls convert T&S diagnostics from conservative temp
1362  ! to potential temperature.
1363  do j=js,je ; do i=is,ie
1364  work_2d(i,j) = gsw_pt_from_ct(sfc_state%SSS(i,j), sfc_state%SST(i,j))
1365  enddo ; enddo
1366  if (ids%id_sst > 0) call post_data(ids%id_sst, work_2d, diag, mask=g%mask2dT)
1367  else
1368  ! Internal T&S variables are potential temperature & practical salinity
1369  if (ids%id_sst > 0) call post_data(ids%id_sst, sfc_state%SST, diag, mask=g%mask2dT)
1370  endif
1371 
1372  if (tv%S_is_absS) then
1373  ! Internal T&S variables are conservative temperature & absolute salinity
1374  if (ids%id_sssabs > 0) call post_data(ids%id_sssabs, sfc_state%SSS, diag, mask=g%mask2dT)
1375  ! Use TEOS-10 function calls convert T&S diagnostics from absolute salinity
1376  ! to practical salinity.
1377  do j=js,je ; do i=is,ie
1378  work_2d(i,j) = gsw_sp_from_sr(sfc_state%SSS(i,j))
1379  enddo ; enddo
1380  if (ids%id_sss > 0) call post_data(ids%id_sss, work_2d, diag, mask=g%mask2dT)
1381  else
1382  ! Internal T&S variables are potential temperature & practical salinity
1383  if (ids%id_sss > 0) call post_data(ids%id_sss, sfc_state%SSS, diag, mask=g%mask2dT)
1384  endif
1385 
1386  if (ids%id_sst_sq > 0) then
1387  do j=js,je ; do i=is,ie
1388  work_2d(i,j) = sfc_state%SST(i,j)*sfc_state%SST(i,j)
1389  enddo ; enddo
1390  call post_data(ids%id_sst_sq, work_2d, diag, mask=g%mask2dT)
1391  endif
1392  if (ids%id_sss_sq > 0) then
1393  do j=js,je ; do i=is,ie
1394  work_2d(i,j) = sfc_state%SSS(i,j)*sfc_state%SSS(i,j)
1395  enddo ; enddo
1396  call post_data(ids%id_sss_sq, work_2d, diag, mask=g%mask2dT)
1397  endif
1398 
1399  call coupler_type_send_data(sfc_state%tr_fields, get_diag_time_end(diag))
1400 
1401 end subroutine post_surface_thermo_diags
1402 
1403 
1404 !> This routine posts diagnostics of the transports, including the subgridscale
1405 !! contributions.
1406 subroutine post_transport_diagnostics(G, GV, US, uhtr, vhtr, h, IDs, diag_pre_dyn, &
1407  diag, dt_trans, Reg)
1408  type(ocean_grid_type), intent(inout) :: g !< ocean grid structure
1409  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
1410  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1411  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(in) :: uhtr !< Accumulated zonal thickness fluxes
1412  !! used to advect tracers [H L2 ~> m3 or kg]
1413  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(in) :: vhtr !< Accumulated meridional thickness fluxes
1414  !! used to advect tracers [H L2 ~> m3 or kg]
1415  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1416  intent(in) :: h !< The updated layer thicknesses [H ~> m or kg m-2]
1417  type(transport_diag_ids), intent(in) :: ids !< A structure with the diagnostic IDs.
1418  type(diag_grid_storage), intent(inout) :: diag_pre_dyn !< Stored grids from before dynamics
1419  type(diag_ctrl), intent(inout) :: diag !< regulates diagnostic output
1420  real, intent(in) :: dt_trans !< total time step associated with the transports [T ~> s].
1421  type(tracer_registry_type), pointer :: reg !< Pointer to the tracer registry
1422 
1423  ! Local variables
1424  real, dimension(SZIB_(G), SZJ_(G)) :: umo2d ! Diagnostics of integrated mass transport [R Z L2 T-1 ~> kg s-1]
1425  real, dimension(SZI_(G), SZJB_(G)) :: vmo2d ! Diagnostics of integrated mass transport [R Z L2 T-1 ~> kg s-1]
1426  real, dimension(SZIB_(G), SZJ_(G), SZK_(G)) :: umo ! Diagnostics of layer mass transport [R Z L2 T-1 ~> kg s-1]
1427  real, dimension(SZI_(G), SZJB_(G), SZK_(G)) :: vmo ! Diagnostics of layer mass transport [R Z L2 T-1 ~> kg s-1]
1428  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_tend ! Change in layer thickness due to dynamics
1429  ! [H s-1 ~> m s-1 or kg m-2 s-1].
1430  real :: idt ! The inverse of the time interval [T-1 ~> s-1]
1431  real :: h_to_rz_dt ! A conversion factor from accumulated transports to fluxes
1432  ! [R Z H-1 T-1 ~> kg m-3 s-1 or s-1].
1433  integer :: i, j, k, is, ie, js, je, nz
1434  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1435 
1436  idt = 1. / dt_trans
1437  h_to_rz_dt = gv%H_to_RZ * idt
1438 
1439  call diag_save_grids(diag)
1440  call diag_copy_storage_to_diag(diag, diag_pre_dyn)
1441 
1442  if (ids%id_umo_2d > 0) then
1443  umo2d(:,:) = 0.0
1444  do k=1,nz ; do j=js,je ; do i=is-1,ie
1445  umo2d(i,j) = umo2d(i,j) + uhtr(i,j,k) * h_to_rz_dt
1446  enddo ; enddo ; enddo
1447  call post_data(ids%id_umo_2d, umo2d, diag)
1448  endif
1449  if (ids%id_umo > 0) then
1450  ! Convert to kg/s.
1451  do k=1,nz ; do j=js,je ; do i=is-1,ie
1452  umo(i,j,k) = uhtr(i,j,k) * h_to_rz_dt
1453  enddo ; enddo ; enddo
1454  call post_data(ids%id_umo, umo, diag, alt_h = diag_pre_dyn%h_state)
1455  endif
1456  if (ids%id_vmo_2d > 0) then
1457  vmo2d(:,:) = 0.0
1458  do k=1,nz ; do j=js-1,je ; do i=is,ie
1459  vmo2d(i,j) = vmo2d(i,j) + vhtr(i,j,k) * h_to_rz_dt
1460  enddo ; enddo ; enddo
1461  call post_data(ids%id_vmo_2d, vmo2d, diag)
1462  endif
1463  if (ids%id_vmo > 0) then
1464  ! Convert to kg/s.
1465  do k=1,nz ; do j=js-1,je ; do i=is,ie
1466  vmo(i,j,k) = vhtr(i,j,k) * h_to_rz_dt
1467  enddo ; enddo ; enddo
1468  call post_data(ids%id_vmo, vmo, diag, alt_h = diag_pre_dyn%h_state)
1469  endif
1470 
1471  if (ids%id_uhtr > 0) call post_data(ids%id_uhtr, uhtr, diag, alt_h = diag_pre_dyn%h_state)
1472  if (ids%id_vhtr > 0) call post_data(ids%id_vhtr, vhtr, diag, alt_h = diag_pre_dyn%h_state)
1473  if (ids%id_dynamics_h > 0) call post_data(ids%id_dynamics_h, diag_pre_dyn%h_state, diag, &
1474  alt_h = diag_pre_dyn%h_state)
1475  ! Post the change in thicknesses
1476  if (ids%id_dynamics_h_tendency > 0) then
1477  h_tend(:,:,:) = 0.
1478  do k=1,nz ; do j=js,je ; do i=is,ie
1479  h_tend(i,j,k) = (h(i,j,k) - diag_pre_dyn%h_state(i,j,k))*idt
1480  enddo ; enddo ; enddo
1481  call post_data(ids%id_dynamics_h_tendency, h_tend, diag, alt_h = diag_pre_dyn%h_state)
1482  endif
1483 
1484  call post_tracer_transport_diagnostics(g, gv, reg, diag_pre_dyn%h_state, diag)
1485 
1486  call diag_restore_grids(diag)
1487 
1488 end subroutine post_transport_diagnostics
1489 
1490 !> This subroutine registers various diagnostics and allocates space for fields
1491 !! that other diagnostis depend upon.
1492 subroutine mom_diagnostics_init(MIS, ADp, CDp, Time, G, GV, US, param_file, diag, CS, tv)
1493  type(ocean_internal_state), intent(in) :: mis !< For "MOM Internal State" a set of pointers to
1494  !! the fields and accelerations that make up the
1495  !! ocean's internal physical state.
1496  type(accel_diag_ptrs), intent(inout) :: adp !< Structure with pointers to momentum equation
1497  !! terms.
1498  type(cont_diag_ptrs), intent(inout) :: cdp !< Structure with pointers to continuity
1499  !! equation terms.
1500  type(time_type), intent(in) :: time !< Current model time.
1501  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
1502  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
1503  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1504  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
1505  !! parameters.
1506  type(diag_ctrl), target, intent(inout) :: diag !< Structure to regulate diagnostic output.
1507  type(diagnostics_cs), pointer :: cs !< Pointer set to point to control structure
1508  !! for this module.
1509  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
1510  !! thermodynamic variables.
1511 
1512  ! Local variables
1513  real :: omega, f2_min, convert_h
1514  ! This include declares and sets the variable "version".
1515 # include "version_variable.h"
1516  character(len=40) :: mdl = "MOM_diagnostics" ! This module's name.
1517  character(len=48) :: thickness_units, flux_units
1518  real :: wave_speed_min ! A floor in the first mode speed below which 0 is returned [L T-1 ~> m s-1]
1519  real :: wave_speed_tol ! The fractional tolerance for finding the wave speeds [nondim]
1520  logical :: better_speed_est ! If true, use a more robust estimate of the first
1521  ! mode wave speed as the starting point for iterations.
1522  logical :: use_temperature, adiabatic
1523  logical :: default_2018_answers, remap_answers_2018
1524  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, nz, nkml, nkbl
1525  integer :: is, ie, js, je, isq, ieq, jsq, jeq, i, j
1526 
1527  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1528  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1529  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
1530  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1531 
1532  if (associated(cs)) then
1533  call mom_error(warning, "MOM_diagnostics_init called with an associated "// &
1534  "control structure.")
1535  return
1536  endif
1537  allocate(cs)
1538 
1539  cs%diag => diag
1540  use_temperature = associated(tv%T)
1541  call get_param(param_file, mdl, "ADIABATIC", adiabatic, default=.false., &
1542  do_not_log=.true.)
1543 
1544  ! Read all relevant parameters and write them to the model log.
1545  call log_version(param_file, mdl, version, "")
1546  call get_param(param_file, mdl, "DIAG_EBT_MONO_N2_COLUMN_FRACTION", cs%mono_N2_column_fraction, &
1547  "The lower fraction of water column over which N2 is limited as monotonic "// &
1548  "for the purposes of calculating the equivalent barotropic wave speed.", &
1549  units='nondim', default=0.)
1550  call get_param(param_file, mdl, "DIAG_EBT_MONO_N2_DEPTH", cs%mono_N2_depth, &
1551  "The depth below which N2 is limited as monotonic for the "// &
1552  "purposes of calculating the equivalent barotropic wave speed.", &
1553  units='m', scale=us%m_to_Z, default=-1.)
1554  call get_param(param_file, mdl, "INTERNAL_WAVE_SPEED_TOL", wave_speed_tol, &
1555  "The fractional tolerance for finding the wave speeds.", &
1556  units="nondim", default=0.001)
1557  !### Set defaults so that wave_speed_min*wave_speed_tol >= 1e-9 m s-1
1558  call get_param(param_file, mdl, "INTERNAL_WAVE_SPEED_MIN", wave_speed_min, &
1559  "A floor in the first mode speed below which 0 used instead.", &
1560  units="m s-1", default=0.0, scale=us%m_s_to_L_T)
1561  call get_param(param_file, mdl, "INTERNAL_WAVE_SPEED_BETTER_EST", better_speed_est, &
1562  "If true, use a more robust estimate of the first mode wave speed as the "//&
1563  "starting point for iterations.", default=.false.) !### Change the default.
1564  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
1565  "This sets the default value for the various _2018_ANSWERS parameters.", &
1566  default=.false.)
1567  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", remap_answers_2018, &
1568  "If true, use the order of arithmetic and expressions that recover the "//&
1569  "answers from the end of 2018. Otherwise, use updated and more robust "//&
1570  "forms of the same expressions.", default=default_2018_answers)
1571 
1572  if (gv%Boussinesq) then
1573  thickness_units = "m" ; flux_units = "m3 s-1" ; convert_h = gv%H_to_m
1574  else
1575  thickness_units = "kg m-2" ; flux_units = "kg s-1" ; convert_h = gv%H_to_kg_m2
1576  endif
1577 
1578  cs%id_masscello = register_diag_field('ocean_model', 'masscello', diag%axesTL,&
1579  time, 'Mass per unit area of liquid ocean grid cell', 'kg m-2', &
1580  standard_name='sea_water_mass_per_unit_area', v_extensive=.true.)
1581 
1582  cs%id_masso = register_scalar_field('ocean_model', 'masso', time, &
1583  diag, 'Mass of liquid ocean', 'kg', standard_name='sea_water_mass')
1584 
1585  cs%id_thkcello = register_diag_field('ocean_model', 'thkcello', diag%axesTL, time, &
1586  long_name = 'Cell Thickness', standard_name='cell_thickness', &
1587  units='m', conversion=us%Z_to_m, v_extensive=.true.)
1588  cs%id_h_pre_sync = register_diag_field('ocean_model', 'h_pre_sync', diag%axesTL, time, &
1589  long_name = 'Cell thickness from the previous timestep', &
1590  units='m', conversion=gv%H_to_m, v_extensive=.true.)
1591 
1592  ! Note that CS%id_volcello would normally be registered here but because it is a "cell measure" and
1593  ! must be registered first. We earlier stored the handle of volcello but need it here for posting
1594  ! by this module.
1595  cs%id_volcello = diag_get_volume_cell_measure_dm_id(diag)
1596 
1597  if (use_temperature) then
1598  if (tv%T_is_conT) then
1599  cs%id_Tpot = register_diag_field('ocean_model', 'temp', diag%axesTL, &
1600  time, 'Potential Temperature', 'degC')
1601  endif
1602  if (tv%S_is_absS) then
1603  cs%id_Sprac = register_diag_field('ocean_model', 'salt', diag%axesTL, &
1604  time, 'Salinity', 'psu')
1605  endif
1606 
1607  cs%id_tob = register_diag_field('ocean_model','tob', diag%axesT1, time, &
1608  long_name='Sea Water Potential Temperature at Sea Floor', &
1609  standard_name='sea_water_potential_temperature_at_sea_floor', units='degC')
1610  cs%id_sob = register_diag_field('ocean_model','sob',diag%axesT1, time, &
1611  long_name='Sea Water Salinity at Sea Floor', &
1612  standard_name='sea_water_salinity_at_sea_floor', units='psu')
1613 
1614  cs%id_temp_layer_ave = register_diag_field('ocean_model', 'temp_layer_ave', &
1615  diag%axesZL, time, 'Layer Average Ocean Temperature', 'degC')
1616  cs%id_salt_layer_ave = register_diag_field('ocean_model', 'salt_layer_ave', &
1617  diag%axesZL, time, 'Layer Average Ocean Salinity', 'psu')
1618 
1619  cs%id_thetaoga = register_scalar_field('ocean_model', 'thetaoga', &
1620  time, diag, 'Global Mean Ocean Potential Temperature', 'degC',&
1621  standard_name='sea_water_potential_temperature')
1622  cs%id_soga = register_scalar_field('ocean_model', 'soga', &
1623  time, diag, 'Global Mean Ocean Salinity', 'psu', &
1624  standard_name='sea_water_salinity')
1625 
1626  cs%id_tosga = register_scalar_field('ocean_model', 'sst_global', time, diag,&
1627  long_name='Global Area Average Sea Surface Temperature', &
1628  units='degC', standard_name='sea_surface_temperature', &
1629  cmor_field_name='tosga', cmor_standard_name='sea_surface_temperature', &
1630  cmor_long_name='Sea Surface Temperature')
1631  cs%id_sosga = register_scalar_field('ocean_model', 'sss_global', time, diag,&
1632  long_name='Global Area Average Sea Surface Salinity', &
1633  units='psu', standard_name='sea_surface_salinity', &
1634  cmor_field_name='sosga', cmor_standard_name='sea_surface_salinity', &
1635  cmor_long_name='Sea Surface Salinity')
1636  endif
1637 
1638  cs%id_u = register_diag_field('ocean_model', 'u', diag%axesCuL, time, &
1639  'Zonal velocity', 'm s-1', conversion=us%L_T_to_m_s, cmor_field_name='uo', &
1640  cmor_standard_name='sea_water_x_velocity', cmor_long_name='Sea Water X Velocity')
1641  cs%id_v = register_diag_field('ocean_model', 'v', diag%axesCvL, time, &
1642  'Meridional velocity', 'm s-1', conversion=us%L_T_to_m_s, cmor_field_name='vo', &
1643  cmor_standard_name='sea_water_y_velocity', cmor_long_name='Sea Water Y Velocity')
1644  cs%id_h = register_diag_field('ocean_model', 'h', diag%axesTL, time, &
1645  'Layer Thickness', thickness_units, v_extensive=.true., conversion=convert_h)
1646 
1647  cs%id_e = register_diag_field('ocean_model', 'e', diag%axesTi, time, &
1648  'Interface Height Relative to Mean Sea Level', 'm', conversion=us%Z_to_m)
1649  if (cs%id_e>0) call safe_alloc_ptr(cs%e,isd,ied,jsd,jed,nz+1)
1650 
1651  cs%id_e_D = register_diag_field('ocean_model', 'e_D', diag%axesTi, time, &
1652  'Interface Height above the Seafloor', 'm', conversion=us%Z_to_m)
1653  if (cs%id_e_D>0) call safe_alloc_ptr(cs%e_D,isd,ied,jsd,jed,nz+1)
1654 
1655  cs%id_Rml = register_diag_field('ocean_model', 'Rml', diag%axesTL, time, &
1656  'Mixed Layer Coordinate Potential Density', 'kg m-3', conversion=us%R_to_kg_m3)
1657 
1658  cs%id_Rcv = register_diag_field('ocean_model', 'Rho_cv', diag%axesTL, time, &
1659  'Coordinate Potential Density', 'kg m-3', conversion=us%R_to_kg_m3)
1660 
1661  cs%id_rhopot0 = register_diag_field('ocean_model', 'rhopot0', diag%axesTL, time, &
1662  'Potential density referenced to surface', 'kg m-3', conversion=us%R_to_kg_m3)
1663  cs%id_rhopot2 = register_diag_field('ocean_model', 'rhopot2', diag%axesTL, time, &
1664  'Potential density referenced to 2000 dbar', 'kg m-3', conversion=us%R_to_kg_m3)
1665  cs%id_rhoinsitu = register_diag_field('ocean_model', 'rhoinsitu', diag%axesTL, time, &
1666  'In situ density', 'kg m-3', conversion=us%R_to_kg_m3)
1667  cs%id_drho_dT = register_diag_field('ocean_model', 'drho_dT', diag%axesTL, time, &
1668  'Partial derivative of rhoinsitu with respect to temperature (alpha)', 'kg m-3 degC-1')
1669  cs%id_drho_dS = register_diag_field('ocean_model', 'drho_dS', diag%axesTL, time, &
1670  'Partial derivative of rhoinsitu with respect to salinity (beta)', 'kg^2 g-1 m-3')
1671 
1672  cs%id_du_dt = register_diag_field('ocean_model', 'dudt', diag%axesCuL, time, &
1673  'Zonal Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
1674  if ((cs%id_du_dt>0) .and. .not.associated(cs%du_dt)) then
1675  call safe_alloc_ptr(cs%du_dt,isdb,iedb,jsd,jed,nz)
1676  call register_time_deriv(lbound(mis%u), mis%u, cs%du_dt, cs)
1677  endif
1678 
1679  cs%id_dv_dt = register_diag_field('ocean_model', 'dvdt', diag%axesCvL, time, &
1680  'Meridional Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
1681  if ((cs%id_dv_dt>0) .and. .not.associated(cs%dv_dt)) then
1682  call safe_alloc_ptr(cs%dv_dt,isd,ied,jsdb,jedb,nz)
1683  call register_time_deriv(lbound(mis%v), mis%v, cs%dv_dt, cs)
1684  endif
1685 
1686  cs%id_dh_dt = register_diag_field('ocean_model', 'dhdt', diag%axesTL, time, &
1687  'Thickness tendency', trim(thickness_units)//" s-1", conversion=convert_h*us%s_to_T, v_extensive=.true.)
1688  if ((cs%id_dh_dt>0) .and. .not.associated(cs%dh_dt)) then
1689  call safe_alloc_ptr(cs%dh_dt,isd,ied,jsd,jed,nz)
1690  call register_time_deriv(lbound(mis%h), mis%h, cs%dh_dt, cs)
1691  endif
1692 
1693  !CS%id_hf_du_dt = register_diag_field('ocean_model', 'hf_dudt', diag%axesCuL, Time, &
1694  ! 'Fractional Thickness-weighted Zonal Acceleration', 'm s-2', v_extensive=.true., &
1695  ! conversion=US%L_T2_to_m_s2)
1696  !if (CS%id_hf_du_dt > 0) then
1697  ! call safe_alloc_ptr(CS%hf_du_dt,IsdB,IedB,jsd,jed,nz)
1698  ! if (.not.associated(CS%du_dt)) then
1699  ! call safe_alloc_ptr(CS%du_dt,IsdB,IedB,jsd,jed,nz)
1700  ! call register_time_deriv(lbound(MIS%u), MIS%u, CS%du_dt, CS)
1701  ! endif
1702  ! call safe_alloc_ptr(ADp%diag_hfrac_u,IsdB,IedB,jsd,jed,nz)
1703  !endif
1704 
1705  !CS%id_hf_dv_dt = register_diag_field('ocean_model', 'hf_dvdt', diag%axesCvL, Time, &
1706  ! 'Fractional Thickness-weighted Meridional Acceleration', 'm s-2', v_extensive=.true., &
1707  ! conversion=US%L_T2_to_m_s2)
1708  !if (CS%id_hf_dv_dt > 0) then
1709  ! call safe_alloc_ptr(CS%hf_dv_dt,isd,ied,JsdB,JedB,nz)
1710  ! if (.not.associated(CS%dv_dt)) then
1711  ! call safe_alloc_ptr(CS%dv_dt,isd,ied,JsdB,JedB,nz)
1712  ! call register_time_deriv(lbound(MIS%v), MIS%v, CS%dv_dt, CS)
1713  ! endif
1714  ! call safe_alloc_ptr(ADp%diag_hfrac_v,isd,ied,Jsd,JedB,nz)
1715  !endif
1716 
1717  cs%id_hf_du_dt_2d = register_diag_field('ocean_model', 'hf_dudt_2d', diag%axesCu1, time, &
1718  'Depth-sum Fractional Thickness-weighted Zonal Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
1719  if (cs%id_hf_du_dt_2d > 0) then
1720  if (.not.associated(cs%du_dt)) then
1721  call safe_alloc_ptr(cs%du_dt,isdb,iedb,jsd,jed,nz)
1722  call register_time_deriv(lbound(mis%u), mis%u, cs%du_dt, cs)
1723  endif
1724  call safe_alloc_ptr(adp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
1725  endif
1726 
1727  cs%id_hf_dv_dt_2d = register_diag_field('ocean_model', 'hf_dvdt_2d', diag%axesCv1, time, &
1728  'Depth-sum Fractional Thickness-weighted Meridional Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
1729  if (cs%id_hf_dv_dt_2d > 0) then
1730  if (.not.associated(cs%dv_dt)) then
1731  call safe_alloc_ptr(cs%dv_dt,isd,ied,jsdb,jedb,nz)
1732  call register_time_deriv(lbound(mis%v), mis%v, cs%dv_dt, cs)
1733  endif
1734  call safe_alloc_ptr(adp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
1735  endif
1736 
1737  ! layer thickness variables
1738  !if (GV%nk_rho_varies > 0) then
1739  cs%id_h_Rlay = register_diag_field('ocean_model', 'h_rho', diag%axesTL, time, &
1740  'Layer thicknesses in pure potential density coordinates', &
1741  thickness_units, conversion=convert_h)
1742  if (cs%id_h_Rlay>0) call safe_alloc_ptr(cs%h_Rlay,isd,ied,jsd,jed,nz)
1743 
1744  cs%id_uh_Rlay = register_diag_field('ocean_model', 'uh_rho', diag%axesCuL, time, &
1745  'Zonal volume transport in pure potential density coordinates', &
1746  flux_units, conversion=us%L_to_m**2*us%s_to_T*convert_h)
1747  if (cs%id_uh_Rlay>0) call safe_alloc_ptr(cs%uh_Rlay,isdb,iedb,jsd,jed,nz)
1748 
1749  cs%id_vh_Rlay = register_diag_field('ocean_model', 'vh_rho', diag%axesCvL, time, &
1750  'Meridional volume transport in pure potential density coordinates', &
1751  flux_units, conversion=us%L_to_m**2*us%s_to_T*convert_h)
1752  if (cs%id_vh_Rlay>0) call safe_alloc_ptr(cs%vh_Rlay,isd,ied,jsdb,jedb,nz)
1753 
1754  cs%id_uhGM_Rlay = register_diag_field('ocean_model', 'uhGM_rho', diag%axesCuL, time, &
1755  'Zonal volume transport due to interface height diffusion in pure potential '//&
1756  'density coordinates', flux_units, conversion=us%L_to_m**2*us%s_to_T*convert_h)
1757  if (cs%id_uhGM_Rlay>0) call safe_alloc_ptr(cs%uhGM_Rlay,isdb,iedb,jsd,jed,nz)
1758 
1759  cs%id_vhGM_Rlay = register_diag_field('ocean_model', 'vhGM_rho', diag%axesCvL, time, &
1760  'Meridional volume transport due to interface height diffusion in pure potential '//&
1761  'density coordinates', flux_units, conversion=us%L_to_m**2*us%s_to_T*convert_h)
1762  if (cs%id_vhGM_Rlay>0) call safe_alloc_ptr(cs%vhGM_Rlay,isd,ied,jsdb,jedb,nz)
1763  !endif
1764 
1765 
1766  ! terms in the kinetic energy budget
1767  cs%id_KE = register_diag_field('ocean_model', 'KE', diag%axesTL, time, &
1768  'Layer kinetic energy per unit mass', &
1769  'm2 s-2', conversion=us%L_T_to_m_s**2)
1770  if (cs%id_KE>0) call safe_alloc_ptr(cs%KE,isd,ied,jsd,jed,nz)
1771 
1772  cs%id_dKEdt = register_diag_field('ocean_model', 'dKE_dt', diag%axesTL, time, &
1773  'Kinetic Energy Tendency of Layer', &
1774  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1775  if (cs%id_dKEdt>0) call safe_alloc_ptr(cs%dKE_dt,isd,ied,jsd,jed,nz)
1776 
1777  cs%id_PE_to_KE = register_diag_field('ocean_model', 'PE_to_KE', diag%axesTL, time, &
1778  'Potential to Kinetic Energy Conversion of Layer', &
1779  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1780  if (cs%id_PE_to_KE>0) call safe_alloc_ptr(cs%PE_to_KE,isd,ied,jsd,jed,nz)
1781 
1782  cs%id_KE_Coradv = register_diag_field('ocean_model', 'KE_Coradv', diag%axesTL, time, &
1783  'Kinetic Energy Source from Coriolis and Advection', &
1784  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1785  if (cs%id_KE_Coradv>0) call safe_alloc_ptr(cs%KE_Coradv,isd,ied,jsd,jed,nz)
1786 
1787  cs%id_KE_adv = register_diag_field('ocean_model', 'KE_adv', diag%axesTL, time, &
1788  'Kinetic Energy Source from Advection', &
1789  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1790  if (cs%id_KE_adv>0) call safe_alloc_ptr(cs%KE_adv,isd,ied,jsd,jed,nz)
1791 
1792  cs%id_KE_visc = register_diag_field('ocean_model', 'KE_visc', diag%axesTL, time, &
1793  'Kinetic Energy Source from Vertical Viscosity and Stresses', &
1794  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1795  if (cs%id_KE_visc>0) call safe_alloc_ptr(cs%KE_visc,isd,ied,jsd,jed,nz)
1796 
1797  cs%id_KE_horvisc = register_diag_field('ocean_model', 'KE_horvisc', diag%axesTL, time, &
1798  'Kinetic Energy Source from Horizontal Viscosity', &
1799  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1800  if (cs%id_KE_horvisc>0) call safe_alloc_ptr(cs%KE_horvisc,isd,ied,jsd,jed,nz)
1801 
1802  if (.not. adiabatic) then
1803  cs%id_KE_dia = register_diag_field('ocean_model', 'KE_dia', diag%axesTL, time, &
1804  'Kinetic Energy Source from Diapycnal Diffusion', &
1805  'm3 s-3', conversion=gv%H_to_m*(us%L_T_to_m_s**2)*us%s_to_T)
1806  if (cs%id_KE_dia>0) call safe_alloc_ptr(cs%KE_dia,isd,ied,jsd,jed,nz)
1807  endif
1808 
1809  ! gravity wave CFLs
1810  cs%id_cg1 = register_diag_field('ocean_model', 'cg1', diag%axesT1, time, &
1811  'First baroclinic gravity wave speed', 'm s-1', conversion=us%L_T_to_m_s)
1812  cs%id_Rd1 = register_diag_field('ocean_model', 'Rd1', diag%axesT1, time, &
1813  'First baroclinic deformation radius', 'm', conversion=us%L_to_m)
1814  cs%id_cfl_cg1 = register_diag_field('ocean_model', 'CFL_cg1', diag%axesT1, time, &
1815  'CFL of first baroclinic gravity wave = dt*cg1*(1/dx+1/dy)', 'nondim')
1816  cs%id_cfl_cg1_x = register_diag_field('ocean_model', 'CFL_cg1_x', diag%axesT1, time, &
1817  'i-component of CFL of first baroclinic gravity wave = dt*cg1*/dx', 'nondim')
1818  cs%id_cfl_cg1_y = register_diag_field('ocean_model', 'CFL_cg1_y', diag%axesT1, time, &
1819  'j-component of CFL of first baroclinic gravity wave = dt*cg1*/dy', 'nondim')
1820  cs%id_cg_ebt = register_diag_field('ocean_model', 'cg_ebt', diag%axesT1, time, &
1821  'Equivalent barotropic gravity wave speed', 'm s-1', conversion=us%L_T_to_m_s)
1822  cs%id_Rd_ebt = register_diag_field('ocean_model', 'Rd_ebt', diag%axesT1, time, &
1823  'Equivalent barotropic deformation radius', 'm', conversion=us%L_to_m)
1824  cs%id_p_ebt = register_diag_field('ocean_model', 'p_ebt', diag%axesTL, time, &
1825  'Equivalent barotropic modal strcuture', 'nondim')
1826 
1827  if ((cs%id_cg1>0) .or. (cs%id_Rd1>0) .or. (cs%id_cfl_cg1>0) .or. &
1828  (cs%id_cfl_cg1_x>0) .or. (cs%id_cfl_cg1_y>0) .or. &
1829  (cs%id_cg_ebt>0) .or. (cs%id_Rd_ebt>0) .or. (cs%id_p_ebt>0)) then
1830  call wave_speed_init(cs%wave_speed_CSp, remap_answers_2018=remap_answers_2018, &
1831  better_speed_est=better_speed_est, min_speed=wave_speed_min, &
1832  wave_speed_tol=wave_speed_tol)
1833  call wave_speed_init(cs%wave_speed_CSp, remap_answers_2018=remap_answers_2018)
1834  call safe_alloc_ptr(cs%cg1,isd,ied,jsd,jed)
1835  if (cs%id_Rd1>0) call safe_alloc_ptr(cs%Rd1,isd,ied,jsd,jed)
1836  if (cs%id_Rd_ebt>0) call safe_alloc_ptr(cs%Rd1,isd,ied,jsd,jed)
1837  if (cs%id_cfl_cg1>0) call safe_alloc_ptr(cs%cfl_cg1,isd,ied,jsd,jed)
1838  if (cs%id_cfl_cg1_x>0) call safe_alloc_ptr(cs%cfl_cg1_x,isd,ied,jsd,jed)
1839  if (cs%id_cfl_cg1_y>0) call safe_alloc_ptr(cs%cfl_cg1_y,isd,ied,jsd,jed)
1840  if (cs%id_p_ebt>0) call safe_alloc_ptr(cs%p_ebt,isd,ied,jsd,jed,nz)
1841  endif
1842 
1843  cs%id_mass_wt = register_diag_field('ocean_model', 'mass_wt', diag%axesT1, time, &
1844  'The column mass for calculating mass-weighted average properties', 'kg m-2', conversion=us%RZ_to_kg_m2)
1845 
1846  if (use_temperature) then
1847  cs%id_temp_int = register_diag_field('ocean_model', 'temp_int', diag%axesT1, time, &
1848  'Density weighted column integrated potential temperature', 'degC kg m-2', conversion=us%RZ_to_kg_m2, &
1849  cmor_field_name='opottempmint', &
1850  cmor_long_name='integral_wrt_depth_of_product_of_sea_water_density_and_potential_temperature',&
1851  cmor_standard_name='Depth integrated density times potential temperature')
1852 
1853  cs%id_salt_int = register_diag_field('ocean_model', 'salt_int', diag%axesT1, time, &
1854  'Density weighted column integrated salinity', 'psu kg m-2', conversion=us%RZ_to_kg_m2, &
1855  cmor_field_name='somint', &
1856  cmor_long_name='integral_wrt_depth_of_product_of_sea_water_density_and_salinity',&
1857  cmor_standard_name='Depth integrated density times salinity')
1858  endif
1859 
1860  cs%id_col_mass = register_diag_field('ocean_model', 'col_mass', diag%axesT1, time, &
1861  'The column integrated in situ density', 'kg m-2', conversion=us%RZ_to_kg_m2)
1862 
1863  cs%id_col_ht = register_diag_field('ocean_model', 'col_height', diag%axesT1, time, &
1864  'The height of the water column', 'm', conversion=us%Z_to_m)
1865  cs%id_pbo = register_diag_field('ocean_model', 'pbo', diag%axesT1, time, &
1866  long_name='Sea Water Pressure at Sea Floor', standard_name='sea_water_pressure_at_sea_floor', &
1867  units='Pa', conversion=us%RL2_T2_to_Pa)
1868 
1869  call set_dependent_diagnostics(mis, adp, cdp, g, cs)
1870 
1871 end subroutine mom_diagnostics_init
1872 
1873 
1874 !> Register diagnostics of the surface state and integrated quantities
1875 subroutine register_surface_diags(Time, G, US, IDs, diag, tv)
1876  type(time_type), intent(in) :: time !< current model time
1877  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
1878  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1879  type(surface_diag_ids), intent(inout) :: ids !< A structure with the diagnostic IDs.
1880  type(diag_ctrl), intent(inout) :: diag !< regulates diagnostic output
1881  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
1882 
1883  ! Vertically integrated, budget, and surface state diagnostics
1884  ids%id_volo = register_scalar_field('ocean_model', 'volo', time, diag,&
1885  long_name='Total volume of liquid ocean', units='m3', &
1886  standard_name='sea_water_volume')
1887  ids%id_zos = register_diag_field('ocean_model', 'zos', diag%axesT1, time,&
1888  standard_name = 'sea_surface_height_above_geoid', &
1889  long_name= 'Sea surface height above geoid', units='m')
1890  ids%id_zossq = register_diag_field('ocean_model', 'zossq', diag%axesT1, time,&
1891  standard_name='square_of_sea_surface_height_above_geoid', &
1892  long_name='Square of sea surface height above geoid', units='m2')
1893  ids%id_ssh = register_diag_field('ocean_model', 'SSH', diag%axesT1, time, &
1894  'Sea Surface Height', 'm')
1895  ids%id_ssh_ga = register_scalar_field('ocean_model', 'ssh_ga', time, diag,&
1896  long_name='Area averaged sea surface height', units='m', &
1897  standard_name='area_averaged_sea_surface_height')
1898  ids%id_ssu = register_diag_field('ocean_model', 'SSU', diag%axesCu1, time, &
1899  'Sea Surface Zonal Velocity', 'm s-1', conversion=us%L_T_to_m_s)
1900  ids%id_ssv = register_diag_field('ocean_model', 'SSV', diag%axesCv1, time, &
1901  'Sea Surface Meridional Velocity', 'm s-1', conversion=us%L_T_to_m_s)
1902  ids%id_speed = register_diag_field('ocean_model', 'speed', diag%axesT1, time, &
1903  'Sea Surface Speed', 'm s-1', conversion=us%L_T_to_m_s)
1904 
1905  if (associated(tv%T)) then
1906  ids%id_sst = register_diag_field('ocean_model', 'SST', diag%axesT1, time, &
1907  'Sea Surface Temperature', 'degC', cmor_field_name='tos', &
1908  cmor_long_name='Sea Surface Temperature', &
1909  cmor_standard_name='sea_surface_temperature')
1910  ids%id_sst_sq = register_diag_field('ocean_model', 'SST_sq', diag%axesT1, time, &
1911  'Sea Surface Temperature Squared', 'degC2', cmor_field_name='tossq', &
1912  cmor_long_name='Square of Sea Surface Temperature ', &
1913  cmor_standard_name='square_of_sea_surface_temperature')
1914  ids%id_sss = register_diag_field('ocean_model', 'SSS', diag%axesT1, time, &
1915  'Sea Surface Salinity', 'psu', cmor_field_name='sos', &
1916  cmor_long_name='Sea Surface Salinity', &
1917  cmor_standard_name='sea_surface_salinity')
1918  ids%id_sss_sq = register_diag_field('ocean_model', 'SSS_sq', diag%axesT1, time, &
1919  'Sea Surface Salinity Squared', 'psu', cmor_field_name='sossq', &
1920  cmor_long_name='Square of Sea Surface Salinity ', &
1921  cmor_standard_name='square_of_sea_surface_salinity')
1922  if (tv%T_is_conT) then
1923  ids%id_sstcon = register_diag_field('ocean_model', 'conSST', diag%axesT1, time, &
1924  'Sea Surface Conservative Temperature', 'Celsius')
1925  endif
1926  if (tv%S_is_absS) then
1927  ids%id_sssabs = register_diag_field('ocean_model', 'absSSS', diag%axesT1, time, &
1928  'Sea Surface Absolute Salinity', 'g kg-1')
1929  endif
1930  if (associated(tv%frazil)) then
1931  ids%id_fraz = register_diag_field('ocean_model', 'frazil', diag%axesT1, time, &
1932  'Heat from frazil formation', 'W m-2', conversion=us%QRZ_T_to_W_m2, &
1933  cmor_field_name='hfsifrazil', &
1934  cmor_standard_name='heat_flux_into_sea_water_due_to_frazil_ice_formation', &
1935  cmor_long_name='Heat Flux into Sea Water due to Frazil Ice Formation')
1936  endif
1937  endif
1938 
1939  ids%id_salt_deficit = register_diag_field('ocean_model', 'salt_deficit', diag%axesT1, time, &
1940  'Salt sink in ocean due to ice flux', &
1941  'psu m-2 s-1', conversion=us%RZ_T_to_kg_m2s)
1942  ids%id_Heat_PmE = register_diag_field('ocean_model', 'Heat_PmE', diag%axesT1, time, &
1943  'Heat flux into ocean from mass flux into ocean', &
1944  'W m-2', conversion=us%QRZ_T_to_W_m2)
1945  ids%id_intern_heat = register_diag_field('ocean_model', 'internal_heat', diag%axesT1, time,&
1946  'Heat flux into ocean from geothermal or other internal sources', &
1947  'W m-2', conversion=us%QRZ_T_to_W_m2)
1948 
1949 end subroutine register_surface_diags
1950 
1951 !> Register certain diagnostics related to transports
1952 subroutine register_transport_diags(Time, G, GV, US, IDs, diag)
1953  type(time_type), intent(in) :: time !< current model time
1954  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
1955  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
1956  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1957  type(transport_diag_ids), intent(inout) :: ids !< A structure with the diagnostic IDs.
1958  type(diag_ctrl), intent(inout) :: diag !< regulates diagnostic output
1959 
1960  real :: h_convert
1961  character(len=48) :: thickness_units
1962 
1963  thickness_units = get_thickness_units(gv)
1964  if (gv%Boussinesq) then
1965  h_convert = gv%H_to_m
1966  else
1967  h_convert = gv%H_to_kg_m2
1968  endif
1969 
1970  ! Diagnostics related to tracer and mass transport
1971  ids%id_uhtr = register_diag_field('ocean_model', 'uhtr', diag%axesCuL, time, &
1972  'Accumulated zonal thickness fluxes to advect tracers', 'kg', &
1973  y_cell_method='sum', v_extensive=.true., conversion=h_convert*us%L_to_m**2)
1974  ids%id_vhtr = register_diag_field('ocean_model', 'vhtr', diag%axesCvL, time, &
1975  'Accumulated meridional thickness fluxes to advect tracers', 'kg', &
1976  x_cell_method='sum', v_extensive=.true., conversion=h_convert*us%L_to_m**2)
1977  ids%id_umo = register_diag_field('ocean_model', 'umo', &
1978  diag%axesCuL, time, 'Ocean Mass X Transport', &
1979  'kg s-1', conversion=us%RZ_T_to_kg_m2s*us%L_to_m**2, &
1980  standard_name='ocean_mass_x_transport', y_cell_method='sum', v_extensive=.true.)
1981  ids%id_vmo = register_diag_field('ocean_model', 'vmo', &
1982  diag%axesCvL, time, 'Ocean Mass Y Transport', &
1983  'kg s-1', conversion=us%RZ_T_to_kg_m2s*us%L_to_m**2, &
1984  standard_name='ocean_mass_y_transport', x_cell_method='sum', v_extensive=.true.)
1985  ids%id_umo_2d = register_diag_field('ocean_model', 'umo_2d', &
1986  diag%axesCu1, time, 'Ocean Mass X Transport Vertical Sum', &
1987  'kg s-1', conversion=us%RZ_T_to_kg_m2s*us%L_to_m**2, &
1988  standard_name='ocean_mass_x_transport_vertical_sum', y_cell_method='sum')
1989  ids%id_vmo_2d = register_diag_field('ocean_model', 'vmo_2d', &
1990  diag%axesCv1, time, 'Ocean Mass Y Transport Vertical Sum', &
1991  'kg s-1', conversion=us%RZ_T_to_kg_m2s*us%L_to_m**2, &
1992  standard_name='ocean_mass_y_transport_vertical_sum', x_cell_method='sum')
1993  ids%id_dynamics_h = register_diag_field('ocean_model','dynamics_h', &
1994  diag%axesTl, time, 'Change in layer thicknesses due to horizontal dynamics', &
1995  'm s-1', v_extensive=.true., conversion=gv%H_to_m)
1996  ids%id_dynamics_h_tendency = register_diag_field('ocean_model','dynamics_h_tendency', &
1997  diag%axesTl, time, 'Change in layer thicknesses due to horizontal dynamics', &
1998  'm s-1', v_extensive=.true., conversion=gv%H_to_m*us%s_to_T)
1999 
2000 end subroutine register_transport_diags
2001 
2002 !> Offers the static fields in the ocean grid type for output via the diag_manager.
2003 subroutine write_static_fields(G, GV, US, tv, diag)
2004  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
2005  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
2006  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
2007  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
2008  type(diag_ctrl), target, intent(inout) :: diag !< regulates diagnostic output
2009 
2010  ! Local variables
2011  integer :: id
2012  logical :: use_temperature
2013 
2014  id = register_static_field('ocean_model', 'geolat', diag%axesT1, &
2015  'Latitude of tracer (T) points', 'degrees_north')
2016  if (id > 0) call post_data(id, g%geoLatT, diag, .true.)
2017 
2018  id = register_static_field('ocean_model', 'geolon', diag%axesT1, &
2019  'Longitude of tracer (T) points', 'degrees_east')
2020  if (id > 0) call post_data(id, g%geoLonT, diag, .true.)
2021 
2022  id = register_static_field('ocean_model', 'geolat_c', diag%axesB1, &
2023  'Latitude of corner (Bu) points', 'degrees_north', interp_method='none')
2024  if (id > 0) call post_data(id, g%geoLatBu, diag, .true.)
2025 
2026  id = register_static_field('ocean_model', 'geolon_c', diag%axesB1, &
2027  'Longitude of corner (Bu) points', 'degrees_east', interp_method='none')
2028  if (id > 0) call post_data(id, g%geoLonBu, diag, .true.)
2029 
2030  id = register_static_field('ocean_model', 'geolat_v', diag%axesCv1, &
2031  'Latitude of meridional velocity (Cv) points', 'degrees_north', interp_method='none')
2032  if (id > 0) call post_data(id, g%geoLatCv, diag, .true.)
2033 
2034  id = register_static_field('ocean_model', 'geolon_v', diag%axesCv1, &
2035  'Longitude of meridional velocity (Cv) points', 'degrees_east', interp_method='none')
2036  if (id > 0) call post_data(id, g%geoLonCv, diag, .true.)
2037 
2038  id = register_static_field('ocean_model', 'geolat_u', diag%axesCu1, &
2039  'Latitude of zonal velocity (Cu) points', 'degrees_north', interp_method='none')
2040  if (id > 0) call post_data(id, g%geoLatCu, diag, .true.)
2041 
2042  id = register_static_field('ocean_model', 'geolon_u', diag%axesCu1, &
2043  'Longitude of zonal velocity (Cu) points', 'degrees_east', interp_method='none')
2044  if (id > 0) call post_data(id, g%geoLonCu, diag, .true.)
2045 
2046  id = register_static_field('ocean_model', 'area_t', diag%axesT1, &
2047  'Surface area of tracer (T) cells', 'm2', conversion=us%L_to_m**2, &
2048  cmor_field_name='areacello', cmor_standard_name='cell_area', &
2049  cmor_long_name='Ocean Grid-Cell Area', &
2050  x_cell_method='sum', y_cell_method='sum', area_cell_method='sum')
2051  if (id > 0) then
2052  call post_data(id, g%areaT, diag, .true.)
2053  call diag_register_area_ids(diag, id_area_t=id)
2054  endif
2055 
2056  id = register_static_field('ocean_model', 'area_u', diag%axesCu1, &
2057  'Surface area of x-direction flow (U) cells', 'm2', conversion=us%L_to_m**2, &
2058  cmor_field_name='areacello_cu', cmor_standard_name='cell_area', &
2059  cmor_long_name='Ocean Grid-Cell Area', &
2060  x_cell_method='sum', y_cell_method='sum', area_cell_method='sum')
2061  if (id > 0) call post_data(id, g%areaCu, diag, .true.)
2062 
2063  id = register_static_field('ocean_model', 'area_v', diag%axesCv1, &
2064  'Surface area of y-direction flow (V) cells', 'm2', conversion=us%L_to_m**2, &
2065  cmor_field_name='areacello_cv', cmor_standard_name='cell_area', &
2066  cmor_long_name='Ocean Grid-Cell Area', &
2067  x_cell_method='sum', y_cell_method='sum', area_cell_method='sum')
2068  if (id > 0) call post_data(id, g%areaCv, diag, .true.)
2069 
2070  id = register_static_field('ocean_model', 'area_q', diag%axesB1, &
2071  'Surface area of B-grid flow (Q) cells', 'm2', conversion=us%L_to_m**2, &
2072  cmor_field_name='areacello_bu', cmor_standard_name='cell_area', &
2073  cmor_long_name='Ocean Grid-Cell Area', &
2074  x_cell_method='sum', y_cell_method='sum', area_cell_method='sum')
2075  if (id > 0) call post_data(id, g%areaBu, diag, .true.)
2076 
2077  id = register_static_field('ocean_model', 'depth_ocean', diag%axesT1, &
2078  'Depth of the ocean at tracer points', 'm', conversion=us%Z_to_m, &
2079  standard_name='sea_floor_depth_below_geoid', &
2080  cmor_field_name='deptho', cmor_long_name='Sea Floor Depth', &
2081  cmor_standard_name='sea_floor_depth_below_geoid', area=diag%axesT1%id_area, &
2082  x_cell_method='mean', y_cell_method='mean', area_cell_method='mean')
2083  if (id > 0) call post_data(id, g%bathyT, diag, .true., mask=g%mask2dT)
2084 
2085  id = register_static_field('ocean_model', 'wet', diag%axesT1, &
2086  '0 if land, 1 if ocean at tracer points', 'none', area=diag%axesT1%id_area)
2087  if (id > 0) call post_data(id, g%mask2dT, diag, .true.)
2088 
2089  id = register_static_field('ocean_model', 'wet_c', diag%axesB1, &
2090  '0 if land, 1 if ocean at corner (Bu) points', 'none', interp_method='none')
2091  if (id > 0) call post_data(id, g%mask2dBu, diag, .true.)
2092 
2093  id = register_static_field('ocean_model', 'wet_u', diag%axesCu1, &
2094  '0 if land, 1 if ocean at zonal velocity (Cu) points', 'none', interp_method='none')
2095  if (id > 0) call post_data(id, g%mask2dCu, diag, .true.)
2096 
2097  id = register_static_field('ocean_model', 'wet_v', diag%axesCv1, &
2098  '0 if land, 1 if ocean at meridional velocity (Cv) points', 'none', interp_method='none')
2099  if (id > 0) call post_data(id, g%mask2dCv, diag, .true.)
2100 
2101  id = register_static_field('ocean_model', 'Coriolis', diag%axesB1, &
2102  'Coriolis parameter at corner (Bu) points', 's-1', interp_method='none', conversion=us%s_to_T)
2103  if (id > 0) call post_data(id, g%CoriolisBu, diag, .true.)
2104 
2105  id = register_static_field('ocean_model', 'dxt', diag%axesT1, &
2106  'Delta(x) at thickness/tracer points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2107  if (id > 0) call post_data(id, g%dxT, diag, .true.)
2108 
2109  id = register_static_field('ocean_model', 'dyt', diag%axesT1, &
2110  'Delta(y) at thickness/tracer points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2111  if (id > 0) call post_data(id, g%dyT, diag, .true.)
2112 
2113  id = register_static_field('ocean_model', 'dxCu', diag%axesCu1, &
2114  'Delta(x) at u points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2115  if (id > 0) call post_data(id, g%dxCu, diag, .true.)
2116 
2117  id = register_static_field('ocean_model', 'dyCu', diag%axesCu1, &
2118  'Delta(y) at u points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2119  if (id > 0) call post_data(id, g%dyCu, diag, .true.)
2120 
2121  id = register_static_field('ocean_model', 'dxCv', diag%axesCv1, &
2122  'Delta(x) at v points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2123  if (id > 0) call post_data(id, g%dxCv, diag, .true.)
2124 
2125  id = register_static_field('ocean_model', 'dyCv', diag%axesCv1, &
2126  'Delta(y) at v points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2127  if (id > 0) call post_data(id, g%dyCv, diag, .true.)
2128 
2129  id = register_static_field('ocean_model', 'dyCuo', diag%axesCu1, &
2130  'Open meridional grid spacing at u points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2131  if (id > 0) call post_data(id, g%dy_Cu, diag, .true.)
2132 
2133  id = register_static_field('ocean_model', 'dxCvo', diag%axesCv1, &
2134  'Open zonal grid spacing at v points (meter)', 'm', interp_method='none', conversion=us%L_to_m)
2135  if (id > 0) call post_data(id, g%dx_Cv, diag, .true.)
2136 
2137  id = register_static_field('ocean_model', 'sin_rot', diag%axesT1, &
2138  'sine of the clockwise angle of the ocean grid north to true north', 'none')
2139  if (id > 0) call post_data(id, g%sin_rot, diag, .true.)
2140 
2141  id = register_static_field('ocean_model', 'cos_rot', diag%axesT1, &
2142  'cosine of the clockwise angle of the ocean grid north to true north', 'none')
2143  if (id > 0) call post_data(id, g%cos_rot, diag, .true.)
2144 
2145 
2146  ! This static diagnostic is from CF 1.8, and is the fraction of a cell
2147  ! covered by ocean, given as a percentage (poorly named).
2148  id = register_static_field('ocean_model', 'area_t_percent', diag%axesT1, &
2149  'Percentage of cell area covered by ocean', '%', conversion=100.0, &
2150  cmor_field_name='sftof', cmor_standard_name='SeaAreaFraction', &
2151  cmor_long_name='Sea Area Fraction', &
2152  x_cell_method='mean', y_cell_method='mean', area_cell_method='mean')
2153  if (id > 0) call post_data(id, g%mask2dT, diag, .true.)
2154 
2155 
2156  id = register_static_field('ocean_model','Rho_0', diag%axesNull, &
2157  'mean ocean density used with the Boussinesq approximation', &
2158  'kg m-3', conversion=us%R_to_kg_m3, cmor_field_name='rhozero', &
2159  cmor_standard_name='reference_sea_water_density_for_boussinesq_approximation', &
2160  cmor_long_name='reference sea water density for boussinesq approximation')
2161  if (id > 0) call post_data(id, gv%Rho0, diag, .true.)
2162 
2163  use_temperature = associated(tv%T)
2164  if (use_temperature) then
2165  id = register_static_field('ocean_model','C_p', diag%axesNull, &
2166  'heat capacity of sea water', 'J kg-1 K-1', conversion=us%Q_to_J_kg, &
2167  cmor_field_name='cpocean', &
2168  cmor_standard_name='specific_heat_capacity_of_sea_water', &
2169  cmor_long_name='specific_heat_capacity_of_sea_water')
2170  if (id > 0) call post_data(id, tv%C_p, diag, .true.)
2171  endif
2172 
2173 end subroutine write_static_fields
2174 
2175 
2176 !> This subroutine sets up diagnostics upon which other diagnostics depend.
2177 subroutine set_dependent_diagnostics(MIS, ADp, CDp, G, CS)
2178  type(ocean_internal_state), intent(in) :: MIS !< For "MOM Internal State" a set of pointers to
2179  !! the fields and accelerations making up ocean
2180  !! internal physical state.
2181  type(accel_diag_ptrs), intent(inout) :: ADp !< Structure pointing to accelerations in
2182  !! momentum equation.
2183  type(cont_diag_ptrs), intent(inout) :: CDp !< Structure pointing to terms in continuity
2184  !! equation.
2185  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
2186  type(diagnostics_CS), pointer :: CS !< Pointer to the control structure for this
2187  !! module.
2188 
2189 ! This subroutine sets up diagnostics upon which other diagnostics depend.
2190  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, nz
2191  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
2192  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
2193 
2194  if (associated(cs%dKE_dt) .or. associated(cs%PE_to_KE) .or. &
2195  associated(cs%KE_CorAdv) .or. associated(cs%KE_adv) .or. &
2196  associated(cs%KE_visc) .or. associated(cs%KE_horvisc) .or. &
2197  associated(cs%KE_dia)) &
2198  call safe_alloc_ptr(cs%KE,isd,ied,jsd,jed,nz)
2199 
2200  if (associated(cs%dKE_dt)) then
2201  if (.not.associated(cs%du_dt)) then
2202  call safe_alloc_ptr(cs%du_dt,isdb,iedb,jsd,jed,nz)
2203  call register_time_deriv(lbound(mis%u), mis%u, cs%du_dt, cs)
2204  endif
2205  if (.not.associated(cs%dv_dt)) then
2206  call safe_alloc_ptr(cs%dv_dt,isd,ied,jsdb,jedb,nz)
2207  call register_time_deriv(lbound(mis%v), mis%v, cs%dv_dt, cs)
2208  endif
2209  if (.not.associated(cs%dh_dt)) then
2210  call safe_alloc_ptr(cs%dh_dt,isd,ied,jsd,jed,nz)
2211  call register_time_deriv(lbound(mis%h), mis%h, cs%dh_dt, cs)
2212  endif
2213  endif
2214 
2215  if (associated(cs%KE_adv)) then
2216  call safe_alloc_ptr(adp%gradKEu,isdb,iedb,jsd,jed,nz)
2217  call safe_alloc_ptr(adp%gradKEv,isd,ied,jsdb,jedb,nz)
2218  endif
2219 
2220  if (associated(cs%KE_visc)) then
2221  call safe_alloc_ptr(adp%du_dt_visc,isdb,iedb,jsd,jed,nz)
2222  call safe_alloc_ptr(adp%dv_dt_visc,isd,ied,jsdb,jedb,nz)
2223  endif
2224 
2225  if (associated(cs%KE_dia)) then
2226  call safe_alloc_ptr(adp%du_dt_dia,isdb,iedb,jsd,jed,nz)
2227  call safe_alloc_ptr(adp%dv_dt_dia,isd,ied,jsdb,jedb,nz)
2228  endif
2229 
2230  if (associated(cs%uhGM_Rlay)) call safe_alloc_ptr(cdp%uhGM,isdb,iedb,jsd,jed,nz)
2231  if (associated(cs%vhGM_Rlay)) call safe_alloc_ptr(cdp%vhGM,isd,ied,jsdb,jedb,nz)
2232 
2233 end subroutine set_dependent_diagnostics
2234 
2235 !> Deallocate memory associated with the diagnostics module
2236 subroutine mom_diagnostics_end(CS, ADp)
2237  type(diagnostics_cs), pointer :: cs !< Control structure returned by a
2238  !! previous call to diagnostics_init.
2239  type(accel_diag_ptrs), intent(inout) :: adp !< structure with pointers to
2240  !! accelerations in momentum equation.
2241  integer :: m
2242 
2243  if (associated(cs%e)) deallocate(cs%e)
2244  if (associated(cs%e_D)) deallocate(cs%e_D)
2245  if (associated(cs%KE)) deallocate(cs%KE)
2246  if (associated(cs%dKE_dt)) deallocate(cs%dKE_dt)
2247  if (associated(cs%PE_to_KE)) deallocate(cs%PE_to_KE)
2248  if (associated(cs%KE_Coradv)) deallocate(cs%KE_Coradv)
2249  if (associated(cs%KE_adv)) deallocate(cs%KE_adv)
2250  if (associated(cs%KE_visc)) deallocate(cs%KE_visc)
2251  if (associated(cs%KE_horvisc)) deallocate(cs%KE_horvisc)
2252  if (associated(cs%KE_dia)) deallocate(cs%KE_dia)
2253  if (associated(cs%dv_dt)) deallocate(cs%dv_dt)
2254  if (associated(cs%dh_dt)) deallocate(cs%dh_dt)
2255  if (associated(cs%du_dt)) deallocate(cs%du_dt)
2256  if (associated(cs%h_Rlay)) deallocate(cs%h_Rlay)
2257  if (associated(cs%uh_Rlay)) deallocate(cs%uh_Rlay)
2258  if (associated(cs%vh_Rlay)) deallocate(cs%vh_Rlay)
2259  if (associated(cs%uhGM_Rlay)) deallocate(cs%uhGM_Rlay)
2260  if (associated(cs%vhGM_Rlay)) deallocate(cs%vhGM_Rlay)
2261 
2262  if (associated(adp%gradKEu)) deallocate(adp%gradKEu)
2263  if (associated(adp%gradKEu)) deallocate(adp%gradKEu)
2264  if (associated(adp%du_dt_visc)) deallocate(adp%du_dt_visc)
2265  if (associated(adp%dv_dt_visc)) deallocate(adp%dv_dt_visc)
2266  if (associated(adp%du_dt_dia)) deallocate(adp%du_dt_dia)
2267  if (associated(adp%dv_dt_dia)) deallocate(adp%dv_dt_dia)
2268  if (associated(adp%du_other)) deallocate(adp%du_other)
2269  if (associated(adp%dv_other)) deallocate(adp%dv_other)
2270 
2271  if (associated(adp%diag_hfrac_u)) deallocate(adp%diag_hfrac_u)
2272  if (associated(adp%diag_hfrac_v)) deallocate(adp%diag_hfrac_v)
2273 
2274  do m=1,cs%num_time_deriv ; deallocate(cs%prev_val(m)%p) ; enddo
2275 
2276  deallocate(cs)
2277 
2278 end subroutine mom_diagnostics_end
2279 
2280 end module mom_diagnostics
Calculates the heights of sruface or all interfaces from layer thicknesses.
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Control structure for MOM_wave_speed.
Functions and routines to take area, volume, mass-weighted, layerwise, zonal or meridional means...
Functions for calculating interface heights, including free surface height.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Routines for calculating baroclinic wave speeds.
The MOM6 facility to parse input files for runtime parameters.
Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
Make a diagnostic available for averaging or output.
A structure with diagnostic IDs of the surface and integrated variables.
Describes various unit conversion factors.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
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
Type to carry basic tracer information.
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances.
Routines for error handling and I/O management.
Stores all the remapping grids and the model&#39;s native space thicknesses.
The control structure for the MOM_diagnostics module.
Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90.
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Provides integrals of density.
Set up a group of halo updates.
Definition: MOM_domains.F90:84
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Find an accurate and order-invariant sum of a distributed 2d or 3d field.
Definition: MOM_coms.F90:53
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
A structure with diagnostic IDs of mass transport related diagnostics.
An overloaded interface to read and log the values of various types of parameters.
A structure for creating arrays of pointers to 3D arrays.
Calculates any requested diagnostic quantities that are not calculated in the various subroutines...