MOM6
MOM_diag_remap.F90
1 !> provides runtime remapping of diagnostics to z star, sigma and
2 !! rho vertical coordinates.
3 !!
4 !! The diag_remap_ctrl type represents a remapping of diagnostics to a particular
5 !! vertical coordinate. The module is used by the diag mediator module in the
6 !! following way:
7 !! 1. diag_remap_init() is called to initialize a diag_remap_ctrl instance.
8 !! 2. diag_remap_configure_axes() is called to read the configuration file and set up the
9 !! vertical coordinate / axes definitions.
10 !! 3. diag_remap_get_axes_info() returns information needed for the diag mediator to
11 !! define new axes for the remapped diagnostics.
12 !! 4. diag_remap_update() is called periodically (whenever h, T or S change) to either
13 !! create or update the target remapping grids.
14 !! 5. diag_remap_do_remap() is called from within a diag post() to do the remapping before
15 !! the diagnostic is written out.
16 
17 
18 ! NOTE: In the following functions, the fields are passed using 1-based
19 ! indexing, which requires special handling within the grid index loops.
20 !
21 ! * diag_remap_do_remap
22 ! * vertically_reintegrate_diag_field
23 ! * vertically_interpolate_diag_field
24 ! * horizontally_average_diag_field
25 !
26 ! Symmetric grids add an additional row of western and southern points to u-
27 ! and v-grids. Non-symmetric grids are 1-based and symmetric grids are
28 ! zero-based, allowing the same expressions to be used when accessing the
29 ! fields. But if u- or v-points become 1-indexed, as in these functions, then
30 ! the stencils must be re-assessed.
31 !
32 ! For interpolation between h and u grids, we use the following relations:
33 !
34 ! h->u: f_u(ig) = 0.5 * (f_h( ig ) + f_h(ig+1))
35 ! f_u(i1) = 0.5 * (f_h(i1-1) + f_h( i1 ))
36 !
37 ! u->h: f_h(ig) = 0.5 * (f_u(ig-1) + f_u( ig ))
38 ! f_h(i1) = 0.5 * (f_u( i1 ) + f_u(i1+1))
39 !
40 ! where ig is the grid index and i1 is the 1-based index. That is, a 1-based
41 ! u-point is ahead of its matching h-point in non-symmetric mode, but behind
42 ! its matching h-point in non-symmetric mode.
43 !
44 ! We can combine these expressions by applying to ig a -1 shift on u-grids and
45 ! a +1 shift on h-grids in symmetric mode.
46 !
47 ! We do not adjust the h-point indices, since they are assumed to be 1-based.
48 ! This is only correct when global indexing is disabled. If global indexing is
49 ! enabled, then all indices will need to be defined relative to the data
50 ! domain.
51 !
52 ! Finally, note that the mask input fields are pointers to arrays which are
53 ! zero-indexed, and do not need any corrections over grid index loops.
54 
55 
57 
58 ! This file is part of MOM6. See LICENSE.md for the license.
59 
60 use mom_coms, only : reproducing_sum_efp, efp_to_real
61 use mom_coms, only : efp_type, assignment(=), efp_sum_across_pes
62 use mom_error_handler, only : mom_error, fatal, assert, warning
63 use mom_diag_vkernels, only : interpolate_column, reintegrate_column
65 use mom_io, only : slasher, mom_read_data
66 use mom_io, only : file_exists, field_size
67 use mom_string_functions, only : lowercase, extractword
68 use mom_grid, only : ocean_grid_type
71 use mom_eos, only : eos_type
72 use mom_remapping, only : remapping_cs, initialize_remapping
73 use mom_remapping, only : remapping_core_h
74 use mom_regridding, only : regridding_cs, initialize_regridding
75 use mom_regridding, only : set_regrid_params, get_regrid_size
76 use mom_regridding, only : getcoordinateinterfaces
77 use mom_regridding, only : get_zlike_cs, get_sigma_cs, get_rho_cs
78 use regrid_consts, only : coordinatemode
79 use coord_zlike, only : build_zstar_column
80 use coord_sigma, only : build_sigma_column
81 use coord_rho, only : build_rho_column
82 
83 use diag_axis_mod, only : get_diag_axis_name
84 use diag_manager_mod, only : diag_axis_init
85 
86 use mom_debugging, only : check_column_integrals
87 implicit none ; private
88 
89 public diag_remap_ctrl
90 public diag_remap_init, diag_remap_end, diag_remap_update, diag_remap_do_remap
91 public diag_remap_configure_axes, diag_remap_axes_configured
92 public diag_remap_calc_hmask
93 public diag_remap_get_axes_info, diag_remap_set_active
94 public diag_remap_diag_registration_closed
95 public vertically_reintegrate_diag_field
96 public vertically_interpolate_diag_field
97 public horizontally_average_diag_field
98 
99 !> Represents remapping of diagnostics to a particular vertical coordinate.
100 !!
101 !! There is one of these types for each vertical coordinate. The vertical axes
102 !! of a diagnostic will reference an instance of this type indicating how (or
103 !! if) the diagnostic should be vertically remapped when being posted.
105  logical :: configured = .false. !< Whether vertical coordinate has been configured
106  logical :: initialized = .false. !< Whether remappping initialized
107  logical :: used = .false. !< Whether this coordinate actually gets used.
108  integer :: vertical_coord = 0 !< The vertical coordinate that we remap to
109  character(len=10) :: vertical_coord_name ='' !< The coordinate name as understood by ALE
110  character(len=16) :: diag_coord_name = '' !< A name for the purpose of run-time parameters
111  character(len=8) :: diag_module_suffix = '' !< The suffix for the module to appear in diag_table
112  type(remapping_cs) :: remap_cs !< Remapping control structure use for this axes
113  type(regridding_cs) :: regrid_cs !< Regridding control structure that defines the coordinates for this axes
114  integer :: nz = 0 !< Number of vertical levels used for remapping
115  real, dimension(:,:,:), allocatable :: h !< Remap grid thicknesses [H ~> m or kg m-2]
116  real, dimension(:,:,:), allocatable :: h_extensive !< Remap grid thicknesses for extensive
117  !! variables [H ~> m or kg m-2]
118  integer :: interface_axes_id = 0 !< Vertical axes id for remapping at interfaces
119  integer :: layer_axes_id = 0 !< Vertical axes id for remapping on layers
120  logical :: answers_2018 !< If true, use the order of arithmetic and expressions for remapping
121  !! that recover the answers from the end of 2018. Otherwise, use
122  !! updated more robust forms of the same expressions.
123 end type diag_remap_ctrl
124 
125 contains
126 
127 !> Initialize a diagnostic remapping type with the given vertical coordinate.
128 subroutine diag_remap_init(remap_cs, coord_tuple, answers_2018)
129  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
130  character(len=*), intent(in) :: coord_tuple !< A string in form of
131  !! MODULE_SUFFIX PARAMETER_SUFFIX COORDINATE_NAME
132  logical, intent(in) :: answers_2018 !< If true, use the order of arithmetic and expressions
133  !! for remapping that recover the answers from the end of 2018.
134  !! Otherwise, use more robust forms of the same expressions.
135 
136  remap_cs%diag_module_suffix = trim(extractword(coord_tuple, 1))
137  remap_cs%diag_coord_name = trim(extractword(coord_tuple, 2))
138  remap_cs%vertical_coord_name = trim(extractword(coord_tuple, 3))
139  remap_cs%vertical_coord = coordinatemode(remap_cs%vertical_coord_name)
140  remap_cs%configured = .false.
141  remap_cs%initialized = .false.
142  remap_cs%used = .false.
143  remap_cs%answers_2018 = answers_2018
144  remap_cs%nz = 0
145 
146 end subroutine diag_remap_init
147 
148 !> De-init a diagnostic remapping type.
149 !! Free allocated memory.
150 subroutine diag_remap_end(remap_cs)
151  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
152 
153  if (allocated(remap_cs%h)) deallocate(remap_cs%h)
154  remap_cs%configured = .false.
155  remap_cs%initialized = .false.
156  remap_cs%used = .false.
157  remap_cs%nz = 0
158 
159 end subroutine diag_remap_end
160 
161 !> Inform that all diagnostics have been registered.
162 !! If _set_active() has not been called on the remapping control structure
163 !! will be disabled. This saves time in the case that a vertical coordinate was
164 !! configured but no diagnostics which use the coordinate appeared in the
165 !! diag_table.
166 subroutine diag_remap_diag_registration_closed(remap_cs)
167  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
168 
169  if (.not. remap_cs%used) then
170  call diag_remap_end(remap_cs)
171  endif
172 
173 end subroutine diag_remap_diag_registration_closed
174 
175 !> Indicate that this remapping type is actually used by the diag manager.
176 !! If this is never called then the type will be disabled to save time.
177 !! See further explanation with diag_remap_registration_closed.
178 subroutine diag_remap_set_active(remap_cs)
179  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
180 
181  remap_cs%used = .true.
182 
183 end subroutine diag_remap_set_active
184 
185 !> Configure the vertical axes for a diagnostic remapping control structure.
186 !! Reads a configuration parameters to determine coordinate generation.
187 subroutine diag_remap_configure_axes(remap_cs, GV, US, param_file)
188  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remap control structure
189  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
190  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
191  type(param_file_type), intent(in) :: param_file !< Parameter file structure
192  ! Local variables
193  integer :: nzi(4), nzl(4), k
194  character(len=200) :: inputdir, string, filename, int_varname, layer_varname
195  character(len=40) :: mod = "MOM_diag_remap" ! This module's name.
196  character(len=8) :: units, expected_units
197  character(len=34) :: longname, string2
198 
199  character(len=256) :: err_msg
200  logical :: ierr
201 
202  real, allocatable, dimension(:) :: interfaces, layers
203 
204  call initialize_regridding(remap_cs%regrid_cs, gv, us, gv%max_depth, param_file, mod, &
205  trim(remap_cs%vertical_coord_name), "DIAG_COORD", trim(remap_cs%diag_coord_name))
206  call set_regrid_params(remap_cs%regrid_cs, min_thickness=0., integrate_downward_for_e=.false.)
207 
208  remap_cs%nz = get_regrid_size(remap_cs%regrid_cs)
209 
210  if (remap_cs%vertical_coord == coordinatemode('SIGMA')) then
211  units = 'nondim'
212  longname = 'Fraction'
213  elseif (remap_cs%vertical_coord == coordinatemode('RHO')) then
214  units = 'kg m-3'
215  longname = 'Target Potential Density'
216  else
217  units = 'meters'
218  longname = 'Depth'
219  endif
220 
221  ! Make axes objects
222  allocate(interfaces(remap_cs%nz+1))
223  allocate(layers(remap_cs%nz))
224 
225  interfaces(:) = getcoordinateinterfaces(remap_cs%regrid_cs)
226  layers(:) = 0.5 * ( interfaces(1:remap_cs%nz) + interfaces(2:remap_cs%nz+1) )
227 
228  remap_cs%interface_axes_id = diag_axis_init(lowercase(trim(remap_cs%diag_coord_name))//'_i', &
229  interfaces, trim(units), 'z', &
230  trim(longname)//' at interface', direction=-1)
231  remap_cs%layer_axes_id = diag_axis_init(lowercase(trim(remap_cs%diag_coord_name))//'_l', &
232  layers, trim(units), 'z', &
233  trim(longname)//' at cell center', direction=-1, &
234  edges=remap_cs%interface_axes_id)
235 
236  ! Axes have now been configured.
237  remap_cs%configured = .true.
238 
239  deallocate(interfaces)
240  deallocate(layers)
241 
242 end subroutine diag_remap_configure_axes
243 
244 !> Get layer and interface axes ids for this coordinate
245 !! Needed when defining axes groups.
246 subroutine diag_remap_get_axes_info(remap_cs, nz, id_layer, id_interface)
247  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coordinate control structure
248  integer, intent(out) :: nz !< Number of vertical levels for the coordinate
249  integer, intent(out) :: id_layer !< 1D-axes id for layer points
250  integer, intent(out) :: id_interface !< 1D-axes id for interface points
251 
252  nz = remap_cs%nz
253  id_layer = remap_cs%layer_axes_id
254  id_interface = remap_cs%interface_axes_id
255 
256 end subroutine diag_remap_get_axes_info
257 
258 
259 !> Whether or not the axes for this vertical coordinated has been configured.
260 !! Configuration is complete when diag_remap_configure_axes() has been
261 !! successfully called.
262 function diag_remap_axes_configured(remap_cs)
263  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coordinate control structure
264  logical :: diag_remap_axes_configured
265 
266  diag_remap_axes_configured = remap_cs%configured
267 
268 end function
269 
270 !> Build/update target vertical grids for diagnostic remapping.
271 !! \note The target grids need to be updated whenever sea surface
272 !! height or layer thicknesses changes. In the case of density-based
273 !! coordinates then technically we should also regenerate the
274 !! target grid whenever T/S change.
275 subroutine diag_remap_update(remap_cs, G, GV, US, h, T, S, eqn_of_state, h_target)
276  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diagnostic coordinate control structure
277  type(ocean_grid_type), pointer :: g !< The ocean's grid type
278  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
279  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
280  real, dimension(:,:,:), intent(in) :: h !< New thickness [H ~> m or kg m-2]
281  real, dimension(:,:,:), intent(in) :: t !< New temperatures [degC]
282  real, dimension(:,:,:), intent(in) :: s !< New salinities [ppt]
283  type(eos_type), pointer :: eqn_of_state !< A pointer to the equation of state
284  real, dimension(:,:,:), intent(inout) :: h_target !< The new diagnostic thicknesses [H ~> m or kg m-2]
285 
286  ! Local variables
287  real, dimension(remap_cs%nz + 1) :: zinterfaces ! Interface positions [H ~> m or kg m-2]
288  real :: h_neglect, h_neglect_edge ! Negligible thicknesses [H ~> m or kg m-2]
289  integer :: i, j, k, nz
290 
291  ! Note that coordinateMode('LAYER') is never 'configured' so will
292  ! always return here.
293  if (.not. remap_cs%configured) then
294  return
295  endif
296 
297  if (.not.remap_cs%answers_2018) then
298  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
299  elseif (gv%Boussinesq) then
300  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
301  else
302  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
303  endif
304  nz = remap_cs%nz
305 
306  if (.not. remap_cs%initialized) then
307  ! Initialize remapping and regridding on the first call
308  call initialize_remapping(remap_cs%remap_cs, 'PPM_IH4', boundary_extrapolation=.false., &
309  answers_2018=remap_cs%answers_2018)
310  remap_cs%initialized = .true.
311  endif
312 
313  ! Calculate remapping thicknesses for different target grids based on
314  ! nominal/target interface locations. This happens for every call on the
315  ! assumption that h, T, S has changed.
316  do j=g%jsc-1, g%jec+1 ; do i=g%isc-1, g%iec+1
317  if (g%mask2dT(i,j)==0.) then
318  h_target(i,j,:) = 0.
319  cycle
320  endif
321 
322  if (remap_cs%vertical_coord == coordinatemode('ZSTAR')) then
323  call build_zstar_column(get_zlike_cs(remap_cs%regrid_cs), &
324  gv%Z_to_H*g%bathyT(i,j), sum(h(i,j,:)), &
325  zinterfaces, zscale=gv%Z_to_H)
326  elseif (remap_cs%vertical_coord == coordinatemode('SIGMA')) then
327  call build_sigma_column(get_sigma_cs(remap_cs%regrid_cs), &
328  gv%Z_to_H*g%bathyT(i,j), sum(h(i,j,:)), zinterfaces)
329  elseif (remap_cs%vertical_coord == coordinatemode('RHO')) then
330  call build_rho_column(get_rho_cs(remap_cs%regrid_cs), g%ke, &
331  gv%Z_to_H*g%bathyT(i,j), h(i,j,:), t(i,j,:), s(i,j,:), &
332  eqn_of_state, zinterfaces, h_neglect, h_neglect_edge)
333  elseif (remap_cs%vertical_coord == coordinatemode('SLIGHT')) then
334 ! call build_slight_column(remap_cs%regrid_cs,remap_cs%remap_cs, nz, &
335 ! GV%Z_to_H*G%bathyT(i,j), sum(h(i,j,:)), zInterfaces)
336  call mom_error(fatal,"diag_remap_update: SLIGHT coordinate not coded for diagnostics yet!")
337  elseif (remap_cs%vertical_coord == coordinatemode('HYCOM1')) then
338 ! call build_hycom1_column(remap_cs%regrid_cs, nz, &
339 ! GV%Z_to_H*G%bathyT(i,j), sum(h(i,j,:)), zInterfaces)
340  call mom_error(fatal,"diag_remap_update: HYCOM1 coordinate not coded for diagnostics yet!")
341  endif
342  do k = 1,nz
343  h_target(i,j,k) = zinterfaces(k) - zinterfaces(k+1)
344  enddo
345  enddo ; enddo
346 
347 end subroutine diag_remap_update
348 
349 !> Remap diagnostic field to alternative vertical grid.
350 subroutine diag_remap_do_remap(remap_cs, G, GV, h, staggered_in_x, staggered_in_y, &
351  mask, field, remapped_field)
352  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
353  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
354  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
355  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses [H ~> m or kg m-2]
356  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
357  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
358  real, dimension(:,:,:), pointer :: mask !< A mask for the field [nondim]
359  real, dimension(:,:,:), intent(in) :: field(:,:,:) !< The diagnostic field to be remapped [A]
360  real, dimension(:,:,:), intent(inout) :: remapped_field !< Field remapped to new coordinate [A]
361  ! Local variables
362  real, dimension(remap_cs%nz) :: h_dest ! Destination thicknesses [H ~> m or kg m-2]
363  real, dimension(size(h,3)) :: h_src ! A column of source thicknesses [H ~> m or kg m-2]
364  real :: h_neglect, h_neglect_edge ! Negligible thicknesses [H ~> m or kg m-2]
365  integer :: nz_src, nz_dest
366  integer :: i, j, k !< Grid index
367  integer :: i1, j1 !< 1-based index
368  integer :: i_lo, i_hi, j_lo, j_hi !< (uv->h) interpolation indices
369  integer :: shift !< Symmetric offset for 1-based indexing
370 
371  call assert(remap_cs%initialized, 'diag_remap_do_remap: remap_cs not initialized.')
372  call assert(size(field, 3) == size(h, 3), &
373  'diag_remap_do_remap: Remap field and thickness z-axes do not match.')
374 
375  if (.not.remap_cs%answers_2018) then
376  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
377  elseif (gv%Boussinesq) then
378  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
379  else
380  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
381  endif
382 
383  nz_src = size(field,3)
384  nz_dest = remap_cs%nz
385  remapped_field(:,:,:) = 0.
386 
387  ! Symmetric grid offset under 1-based indexing; see header for details.
388  shift = 0; if (g%symmetric) shift = 1
389 
390  if (staggered_in_x .and. .not. staggered_in_y) then
391  ! U-points
392  do j=g%jsc, g%jec
393  do i=g%iscB, g%iecB
394  i1 = i - g%isdB + 1
395  i_lo = i1 - shift; i_hi = i_lo + 1
396  if (associated(mask)) then
397  if (mask(i,j,1) == 0.) cycle
398  endif
399  h_src(:) = 0.5 * (h(i_lo,j,:) + h(i_hi,j,:))
400  h_dest(:) = 0.5 * (remap_cs%h(i_lo,j,:) + remap_cs%h(i_hi,j,:))
401  call remapping_core_h(remap_cs%remap_cs, &
402  nz_src, h_src(:), field(i1,j,:), &
403  nz_dest, h_dest(:), remapped_field(i1,j,:), &
404  h_neglect, h_neglect_edge)
405  enddo
406  enddo
407  elseif (staggered_in_y .and. .not. staggered_in_x) then
408  ! V-points
409  do j=g%jscB, g%jecB
410  j1 = j - g%jsdB + 1
411  j_lo = j1 - shift; j_hi = j_lo + 1
412  do i=g%isc, g%iec
413  if (associated(mask)) then
414  if (mask(i,j,1) == 0.) cycle
415  endif
416  h_src(:) = 0.5 * (h(i,j_lo,:) + h(i,j_hi,:))
417  h_dest(:) = 0.5 * (remap_cs%h(i,j_lo,:) + remap_cs%h(i,j_hi,:))
418  call remapping_core_h(remap_cs%remap_cs, &
419  nz_src, h_src(:), field(i,j1,:), &
420  nz_dest, h_dest(:), remapped_field(i,j1,:), &
421  h_neglect, h_neglect_edge)
422  enddo
423  enddo
424  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
425  ! H-points
426  do j=g%jsc, g%jec
427  do i=g%isc, g%iec
428  if (associated(mask)) then
429  if (mask(i,j,1) == 0.) cycle
430  endif
431  h_src(:) = h(i,j,:)
432  h_dest(:) = remap_cs%h(i,j,:)
433  call remapping_core_h(remap_cs%remap_cs, &
434  nz_src, h_src(:), field(i,j,:), &
435  nz_dest, h_dest(:), remapped_field(i,j,:), &
436  h_neglect, h_neglect_edge)
437  enddo
438  enddo
439  else
440  call assert(.false., 'diag_remap_do_remap: Unsupported axis combination')
441  endif
442 
443 end subroutine diag_remap_do_remap
444 
445 !> Calculate masks for target grid
446 subroutine diag_remap_calc_hmask(remap_cs, G, mask)
447  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
448  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
449  real, dimension(:,:,:), intent(out) :: mask !< h-point mask for target grid [nondim]
450  ! Local variables
451  real, dimension(remap_cs%nz) :: h_dest ! Destination thicknesses [H ~> m or kg m-2]
452  integer :: i, j, k
453  logical :: mask_vanished_layers
454  real :: h_tot ! Sum of all thicknesses [H ~> m or kg m-2]
455  real :: h_err ! An estimate of a negligible thickness [H ~> m or kg m-2]
456 
457  call assert(remap_cs%initialized, 'diag_remap_calc_hmask: remap_cs not initialized.')
458 
459  ! Only z*-like diagnostic coordinates should have a 3d mask
460  mask_vanished_layers = (remap_cs%vertical_coord == coordinatemode('ZSTAR'))
461  mask(:,:,:) = 0.
462 
463  do j=g%jsc-1, g%jec+1 ; do i=g%isc-1, g%iec+1
464  if (g%mask2dT(i,j)>0.) then
465  if (mask_vanished_layers) then
466  h_dest(:) = remap_cs%h(i,j,:)
467  h_tot = 0.
468  h_err = 0.
469  do k=1, remap_cs%nz
470  h_tot = h_tot + h_dest(k)
471  ! This is an overestimate of how thick a vanished layer might be, that
472  ! appears due to round-off.
473  h_err = h_err + epsilon(h_tot) * h_tot
474  ! Mask out vanished layers
475  if (h_dest(k)<=8.*h_err) then
476  mask(i,j,k) = 0.
477  else
478  mask(i,j,k) = 1.
479  endif
480  enddo
481  else ! all layers might contain data
482  mask(i,j,:) = 1.
483  endif
484  endif
485  enddo ; enddo
486 
487 end subroutine diag_remap_calc_hmask
488 
489 !> Vertically re-grid an already vertically-integrated diagnostic field to alternative vertical grid.
490 subroutine vertically_reintegrate_diag_field(remap_cs, G, h, h_target, staggered_in_x, staggered_in_y, &
491  mask, field, reintegrated_field)
492  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
493  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
494  real, dimension(:,:,:), intent(in) :: h !< The thicknesses of the source grid [H ~> m or kg m-2]
495  real, dimension(:,:,:), intent(in) :: h_target !< The thicknesses of the target grid [H ~> m or kg m-2]
496  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
497  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
498  real, dimension(:,:,:), pointer :: mask !< A mask for the field [nondim]
499  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped [A]
500  real, dimension(:,:,:), intent(inout) :: reintegrated_field !< Field argument remapped to alternative coordinate [A]
501  ! Local variables
502  real, dimension(remap_cs%nz) :: h_dest ! Destination thicknesses [H ~> m or kg m-2]
503  real, dimension(size(h,3)) :: h_src ! A column of source thicknesses [H ~> m or kg m-2]
504  integer :: nz_src, nz_dest
505  integer :: i, j, k !< Grid index
506  integer :: i1, j1 !< 1-based index
507  integer :: i_lo, i_hi, j_lo, j_hi !< (uv->h) interpolation indices
508  integer :: shift !< Symmetric offset for 1-based indexing
509 
510  call assert(remap_cs%initialized, 'vertically_reintegrate_diag_field: remap_cs not initialized.')
511  call assert(size(field, 3) == size(h, 3), &
512  'vertically_reintegrate_diag_field: Remap field and thickness z-axes do not match.')
513 
514  nz_src = size(field,3)
515  nz_dest = remap_cs%nz
516  reintegrated_field(:,:,:) = 0.
517 
518  ! Symmetric grid offset under 1-based indexing; see header for details.
519  shift = 0; if (g%symmetric) shift = 1
520 
521  if (staggered_in_x .and. .not. staggered_in_y) then
522  ! U-points
523  do j=g%jsc, g%jec
524  do i=g%iscB, g%iecB
525  i1 = i - g%isdB + 1
526  i_lo = i1 - shift; i_hi = i_lo + 1
527  if (associated(mask)) then
528  if (mask(i,j,1) == 0.) cycle
529  endif
530  h_src(:) = 0.5 * (h(i_lo,j,:) + h(i_hi,j,:))
531  h_dest(:) = 0.5 * (h_target(i_lo,j,:) + h_target(i_hi,j,:))
532  call reintegrate_column(nz_src, h_src, field(i1,j,:), &
533  nz_dest, h_dest, 0., reintegrated_field(i1,j,:))
534  enddo
535  enddo
536  elseif (staggered_in_y .and. .not. staggered_in_x) then
537  ! V-points
538  do j=g%jscB, g%jecB
539  j1 = j - g%jsdB + 1
540  j_lo = j1 - shift; j_hi = j_lo + 1
541  do i=g%isc, g%iec
542  if (associated(mask)) then
543  if (mask(i,j,1) == 0.) cycle
544  endif
545  h_src(:) = 0.5 * (h(i,j_lo,:) + h(i,j_hi,:))
546  h_dest(:) = 0.5 * (h_target(i,j_lo,:) + h_target(i,j_hi,:))
547  call reintegrate_column(nz_src, h_src, field(i,j1,:), &
548  nz_dest, h_dest, 0., reintegrated_field(i,j1,:))
549  enddo
550  enddo
551  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
552  ! H-points
553  do j=g%jsc, g%jec
554  do i=g%isc, g%iec
555  if (associated(mask)) then
556  if (mask(i,j,1) == 0.) cycle
557  endif
558  h_src(:) = h(i,j,:)
559  h_dest(:) = h_target(i,j,:)
560  call reintegrate_column(nz_src, h_src, field(i,j,:), &
561  nz_dest, h_dest, 0., reintegrated_field(i,j,:))
562  enddo
563  enddo
564  else
565  call assert(.false., 'vertically_reintegrate_diag_field: Q point remapping is not coded yet.')
566  endif
567 
568 end subroutine vertically_reintegrate_diag_field
569 
570 !> Vertically interpolate diagnostic field to alternative vertical grid.
571 subroutine vertically_interpolate_diag_field(remap_cs, G, h, staggered_in_x, staggered_in_y, &
572  mask, field, interpolated_field)
573  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
574  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
575  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses [H ~> m or kg m-2]
576  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
577  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
578  real, dimension(:,:,:), pointer :: mask !< A mask for the field [nondim]
579  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped [A]
580  real, dimension(:,:,:), intent(inout) :: interpolated_field !< Field argument remapped to alternative coordinate [A]
581  ! Local variables
582  real, dimension(remap_cs%nz) :: h_dest ! Destination thicknesses [H ~> m or kg m-2]
583  real, dimension(size(h,3)) :: h_src ! A column of source thicknesses [H ~> m or kg m-2]
584  integer :: nz_src, nz_dest
585  integer :: i, j, k !< Grid index
586  integer :: i1, j1 !< 1-based index
587  integer :: i_lo, i_hi, j_lo, j_hi !< (uv->h) interpolation indices
588  integer :: shift !< Symmetric offset for 1-based indexing
589 
590  call assert(remap_cs%initialized, 'vertically_interpolate_diag_field: remap_cs not initialized.')
591  call assert(size(field, 3) == size(h, 3)+1, &
592  'vertically_interpolate_diag_field: Remap field and thickness z-axes do not match.')
593 
594  interpolated_field(:,:,:) = 0.
595 
596  nz_src = size(h,3)
597  nz_dest = remap_cs%nz
598 
599  ! Symmetric grid offset under 1-based indexing; see header for details.
600  shift = 0; if (g%symmetric) shift = 1
601 
602  if (staggered_in_x .and. .not. staggered_in_y) then
603  ! U-points
604  do j=g%jsc, g%jec
605  do i=g%iscB, g%iecB
606  i1 = i - g%isdB + 1
607  i_lo = i1 - shift; i_hi = i_lo + 1
608  if (associated(mask)) then
609  if (mask(i,j,1) == 0.) cycle
610  endif
611  h_src(:) = 0.5 * (h(i_lo,j,:) + h(i_hi,j,:))
612  h_dest(:) = 0.5 * (remap_cs%h(i_lo,j,:) + remap_cs%h(i_hi,j,:))
613  call interpolate_column(nz_src, h_src, field(i1,j,:), &
614  nz_dest, h_dest, 0., interpolated_field(i1,j,:))
615  enddo
616  enddo
617  elseif (staggered_in_y .and. .not. staggered_in_x) then
618  ! V-points
619  do j=g%jscB, g%jecB
620  j1 = j - g%jsdB + 1
621  j_lo = j1 - shift; j_hi = j_lo + 1
622  do i=g%isc, g%iec
623  if (associated(mask)) then
624  if (mask(i,j,1) == 0.) cycle
625  endif
626  h_src(:) = 0.5 * (h(i,j_lo,:) + h(i,j_hi,:))
627  h_dest(:) = 0.5 * (remap_cs%h(i,j_lo,:) + remap_cs%h(i,j_hi,:))
628  call interpolate_column(nz_src, h_src, field(i,j1,:), &
629  nz_dest, h_dest, 0., interpolated_field(i,j1,:))
630  enddo
631  enddo
632  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
633  ! H-points
634  do j=g%jsc, g%jec
635  do i=g%isc, g%iec
636  if (associated(mask)) then
637  if (mask(i,j,1) == 0.) cycle
638  endif
639  h_src(:) = h(i,j,:)
640  h_dest(:) = remap_cs%h(i,j,:)
641  call interpolate_column(nz_src, h_src, field(i,j,:), &
642  nz_dest, h_dest, 0., interpolated_field(i,j,:))
643  enddo
644  enddo
645  else
646  call assert(.false., 'vertically_interpolate_diag_field: Q point remapping is not coded yet.')
647  endif
648 
649 end subroutine vertically_interpolate_diag_field
650 
651 !> Horizontally average field
652 subroutine horizontally_average_diag_field(G, GV, h, staggered_in_x, staggered_in_y, &
653  is_layer, is_extensive, &
654  field, averaged_field, &
655  averaged_mask)
656  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
657  type(verticalgrid_type), intent(in) :: gv !< The ocean vertical grid structure
658  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses [H ~> m or kg m-2]
659  logical, intent(in) :: staggered_in_x !< True if the x-axis location is at u or q points
660  logical, intent(in) :: staggered_in_y !< True if the y-axis location is at v or q points
661  logical, intent(in) :: is_layer !< True if the z-axis location is at h points
662  logical, intent(in) :: is_extensive !< True if the z-direction is spatially integrated (over layers)
663  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped [A]
664  real, dimension(:), intent(inout) :: averaged_field !< Field argument horizontally averaged [A]
665  logical, dimension(:), intent(inout) :: averaged_mask !< Mask for horizontally averaged field [nondim]
666 
667  ! Local variables
668  real, dimension(G%isc:G%iec, G%jsc:G%jec, size(field,3)) :: volume, stuff
669  real, dimension(size(field, 3)) :: vol_sum, stuff_sum ! nz+1 is needed for interface averages
670  type(efp_type), dimension(2*size(field,3)) :: sums_efp ! Sums of volume or stuff by layer
671  real :: height ! An average thickness attributed to an velocity point [H ~> m or kg m-2]
672  integer :: i, j, k, nz
673  integer :: i1, j1 !< 1-based index
674 
675  nz = size(field, 3)
676 
677  ! TODO: These averages could potentially be modified to use the function in
678  ! the MOM_spatial_means module.
679  ! NOTE: Reproducible sums must be computed in the original MKS units
680 
681  if (staggered_in_x .and. .not. staggered_in_y) then
682  if (is_layer) then
683  ! U-points
684  do k=1,nz
685  vol_sum(k) = 0.
686  stuff_sum(k) = 0.
687  if (is_extensive) then
688  do j=g%jsc, g%jec ; do i=g%isc, g%iec
689  i1 = i - g%isdB + 1
690  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCu(i,j)) * g%mask2dCu(i,j)
691  stuff(i,j,k) = volume(i,j,k) * field(i1,j,k)
692  enddo ; enddo
693  else ! Intensive
694  do j=g%jsc, g%jec ; do i=g%isc, g%iec
695  i1 = i - g%isdB + 1
696  height = 0.5 * (h(i,j,k) + h(i+1,j,k))
697  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCu(i,j)) &
698  * (gv%H_to_m * height) * g%mask2dCu(i,j)
699  stuff(i,j,k) = volume(i,j,k) * field(i1,j,k)
700  enddo ; enddo
701  endif
702  enddo
703  else ! Interface
704  do k=1,nz
705  do j=g%jsc, g%jec ; do i=g%isc, g%iec
706  i1 = i - g%isdB + 1
707  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCu(i,j)) * g%mask2dCu(i,j)
708  stuff(i,j,k) = volume(i,j,k) * field(i1,j,k)
709  enddo ; enddo
710  enddo
711  endif
712  elseif (staggered_in_y .and. .not. staggered_in_x) then
713  if (is_layer) then
714  ! V-points
715  do k=1,nz
716  if (is_extensive) then
717  do j=g%jsc, g%jec ; do i=g%isc, g%iec
718  j1 = j - g%jsdB + 1
719  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCv(i,j)) * g%mask2dCv(i,j)
720  stuff(i,j,k) = volume(i,j,k) * field(i,j1,k)
721  enddo ; enddo
722  else ! Intensive
723  do j=g%jsc, g%jec ; do i=g%isc, g%iec
724  j1 = j - g%jsdB + 1
725  height = 0.5 * (h(i,j,k) + h(i,j+1,k))
726  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCv(i,j)) &
727  * (gv%H_to_m * height) * g%mask2dCv(i,j)
728  stuff(i,j,k) = volume(i,j,k) * field(i,j1,k)
729  enddo ; enddo
730  endif
731  enddo
732  else ! Interface
733  do k=1,nz
734  do j=g%jsc, g%jec ; do i=g%isc, g%iec
735  j1 = j - g%jsdB + 1
736  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCv(i,j)) * g%mask2dCv(i,j)
737  stuff(i,j,k) = volume(i,j,k) * field(i,j1,k)
738  enddo ; enddo
739  enddo
740  endif
741  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
742  if (is_layer) then
743  ! H-points
744  do k=1,nz
745  if (is_extensive) then
746  do j=g%jsc, g%jec ; do i=g%isc, g%iec
747  if (h(i,j,k) > 0.) then
748  volume(i,j,k) = (g%US%L_to_m**2 * g%areaT(i,j)) * g%mask2dT(i,j)
749  stuff(i,j,k) = volume(i,j,k) * field(i,j,k)
750  else
751  volume(i,j,k) = 0.
752  stuff(i,j,k) = 0.
753  endif
754  enddo ; enddo
755  else ! Intensive
756  do j=g%jsc, g%jec ; do i=g%isc, g%iec
757  volume(i,j,k) = (g%US%L_to_m**2 * g%areaT(i,j)) &
758  * (gv%H_to_m * h(i,j,k)) * g%mask2dT(i,j)
759  stuff(i,j,k) = volume(i,j,k) * field(i,j,k)
760  enddo ; enddo
761  endif
762  enddo
763  else ! Interface
764  do k=1,nz
765  do j=g%jsc, g%jec ; do i=g%isc, g%iec
766  volume(i,j,k) = (g%US%L_to_m**2 * g%areaT(i,j)) * g%mask2dT(i,j)
767  stuff(i,j,k) = volume(i,j,k) * field(i,j,k)
768  enddo ; enddo
769  enddo
770  endif
771  else
772  call assert(.false., 'horizontally_average_diag_field: Q point averaging is not coded yet.')
773  endif
774 
775  ! Packing the sums into a single array with a single call to sum across PEs saves reduces
776  ! the costs of communication.
777  do k=1,nz
778  sums_efp(2*k-1) = reproducing_sum_efp(volume(:,:,k), only_on_pe=.true.)
779  sums_efp(2*k) = reproducing_sum_efp(stuff(:,:,k), only_on_pe=.true.)
780  enddo
781  call efp_sum_across_pes(sums_efp, 2*nz)
782  do k=1,nz
783  vol_sum(k) = efp_to_real(sums_efp(2*k-1))
784  stuff_sum(k) = efp_to_real(sums_efp(2*k))
785  enddo
786 
787  averaged_mask(:) = .true.
788  do k=1,nz
789  if (vol_sum(k) > 0.) then
790  averaged_field(k) = stuff_sum(k) / vol_sum(k)
791  else
792  averaged_field(k) = 0.
793  averaged_mask(k) = .false.
794  endif
795  enddo
796 
797 end subroutine horizontally_average_diag_field
798 
799 end module mom_diag_remap
coord_rho
Regrid columns for the continuous isopycnal (rho) coordinate.
Definition: coord_rho.F90:2
mom_diag_vkernels
Provides kernels for single-column interpolation, re-integration (re-mapping of integrated quantities...
Definition: MOM_diag_vkernels.F90:3
mom_regridding::regridding_cs
Regridding control structure.
Definition: MOM_regridding.F90:45
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_coms::efp_type
The Extended Fixed Point (EFP) type provides a public interface for doing sums and taking differences...
Definition: MOM_coms.F90:74
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_remapping::remapping_cs
Container for remapping parameters.
Definition: MOM_remapping.F90:24
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_coms::reproducing_sum_efp
Find an accurate and order-invariant sum of a distributed 2d field, returning the result in the form ...
Definition: MOM_coms.F90:59
mom_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
regrid_consts
Contains constants for interpreting input parameters that control regridding.
Definition: regrid_consts.F90:2
mom_remapping
Provides column-wise vertical remapping functions.
Definition: MOM_remapping.F90:2
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:108
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
coord_zlike
Regrid columns for a z-like coordinate (z-star, z-level)
Definition: coord_zlike.F90:2
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_diag_remap
provides runtime remapping of diagnostics to z star, sigma and rho vertical coordinates.
Definition: MOM_diag_remap.F90:56
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_regridding
Generates vertical grids as part of the ALE algorithm.
Definition: MOM_regridding.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_coms::efp_sum_across_pes
Sum a value or 1-d array of values across processors, returning the sums in place.
Definition: MOM_coms.F90:64
mom_diag_remap::diag_remap_ctrl
Represents remapping of diagnostics to a particular vertical coordinate.
Definition: MOM_diag_remap.F90:104
coord_sigma
Regrid columns for the sigma coordinate.
Definition: coord_sigma.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26