MOM6
MOM_diag_remap.F90
Go to the documentation of this file.
1 !> This module is used for runtime remapping of diagnostics to z star, sigma and
2 !! rho vertical coordinates. It defines the diag_remap_ctrl type which
3 !! represents a remapping of diagnostics to a particular vertical coordinate.
4 !! The module is used by the diag mediator module in the following way:
5 !! 1) _init() is called to initialise a diag_remap_ctrl instance.
6 !! 2) _configure_axes() is called to read the configuration file and set up the
7 !! vertical coordinate / axes definitions.
8 !! 3) _get_axes_info() returns information needed for the diag mediator to
9 !! define new axes for the remapped diagnostics.
10 !! 4) _update() is called periodically (whenever h, T or S change) to either
11 !! create or update the target remapping grids.
12 !! 5) _do_remap() is called from within a diag post() to do the remapping before
13 !! the diagnostic is written out.
14 
16 
17 ! This file is part of MOM6. See LICENSE.md for the license.
18 
19 use mom_coms, only : sum_across_pes
20 use mom_error_handler, only : mom_error, fatal, assert
23 use mom_io, only : slasher, mom_read_data
24 use mom_io, only : file_exists, field_size
26 use mom_grid, only : ocean_grid_type
28 use mom_eos, only : eos_type
35 use regrid_consts, only : coordinatemode
38 use coord_rho, only : build_rho_column
39 
40 use diag_axis_mod, only : get_diag_axis_name
41 use diag_manager_mod, only : diag_axis_init
42 
43 implicit none ; private
44 
45 public diag_remap_ctrl
53 
54 !> This type represents remapping of diagnostics to a particular vertical
55 !! coordinate.
56 !! There is one of these types for each vertical coordinate. The vertical axes
57 !! of a diagnostic will reference an instance of this type indicating how (or
58 !! if) the diagnostic should be vertically remapped when being posted.
60  logical :: configured = .false. !< Whether vertical coordinate has been configured
61  logical :: initialized = .false. !< Whether remappping initialized
62  logical :: used = .false. !< Whether this coordinate actually gets used.
63  integer :: vertical_coord = 0 !< The vertical coordinate that we remap to
64  character(len=10) :: vertical_coord_name ='' !< The coordinate name as understood by ALE
65  character(len=16) :: diag_coord_name = '' !< A name for the purpose of run-time parameters
66  character(len=8) :: diag_module_suffix = '' !< The suffix for the module to appear in diag_table
67  type(remapping_cs) :: remap_cs !< Remapping control structure use for this axes
68  type(regridding_cs) :: regrid_cs !< Regridding control structure that defines the coordiantes for this axes
69  integer :: nz = 0 !< Number of vertical levels used for remapping
70  real, dimension(:,:,:), allocatable :: h !< Remap grid thicknesses
71  real, dimension(:), allocatable :: dz !< Nominal layer thicknesses
72  integer :: interface_axes_id = 0 !< Vertical axes id for remapping at interfaces
73  integer :: layer_axes_id = 0 !< Vertical axes id for remapping on layers
74 end type diag_remap_ctrl
75 
76 contains
77 
78 !> Initialize a diagnostic remapping type with the given vertical coordinate.
79 subroutine diag_remap_init(remap_cs, coord_tuple)
80  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
81  character(len=*), intent(in) :: coord_tuple !< A string in form of
82  !! MODULE_SUFFIX PARAMETER_SUFFIX COORDINATE_NAME
83 
84  remap_cs%diag_module_suffix = trim(extractword(coord_tuple, 1))
85  remap_cs%diag_coord_name = trim(extractword(coord_tuple, 2))
86  remap_cs%vertical_coord_name = trim(extractword(coord_tuple, 3))
87  remap_cs%vertical_coord = coordinatemode(remap_cs%vertical_coord_name)
88  remap_cs%configured = .false.
89  remap_cs%initialized = .false.
90  remap_cs%used = .false.
91  remap_cs%nz = 0
92 
93 end subroutine diag_remap_init
94 
95 !> De-init a diagnostic remapping type.
96 !! Free allocated memory.
97 subroutine diag_remap_end(remap_cs)
98  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
99 
100  if (allocated(remap_cs%h)) deallocate(remap_cs%h)
101  if (allocated(remap_cs%dz)) deallocate(remap_cs%dz)
102  remap_cs%configured = .false.
103  remap_cs%initialized = .false.
104  remap_cs%used = .false.
105  remap_cs%nz = 0
106 
107 end subroutine diag_remap_end
108 
109 !> Inform that all diagnostics have been registered.
110 !! If _set_active() has not been called on the remapping control structure
111 !! will be disabled. This saves time in the case that a vertical coordinate was
112 !! configured but no diagnostics which use the coordinate appeared in the
113 !! diag_table.
114 subroutine diag_remap_diag_registration_closed(remap_cs)
115  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
116 
117  if (.not. remap_cs%used) then
118  call diag_remap_end(remap_cs)
119  endif
120 
122 
123 !> Indicate that this remapping type is actually used by the diag manager.
124 !! If this is never called then the type will be disabled to save time.
125 !! See further explanation with diag_remap_registration_closed.
126 subroutine diag_remap_set_active(remap_cs)
127  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
128 
129  remap_cs%used = .true.
130 
131 end subroutine diag_remap_set_active
132 
133 !> Configure the vertical axes for a diagnostic remapping control structure.
134 !! Reads a configuration parameters to determine coordinate generation.
135 subroutine diag_remap_configure_axes(remap_cs, GV, param_file)
136  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remap control structure
137  type(verticalgrid_type), intent(in) :: GV !< ocean vertical grid structure
138  type(param_file_type), intent(in) :: param_file !< Parameter file structure
139  ! Local variables
140  integer :: nzi(4), nzl(4), k
141  character(len=200) :: inputdir, string, filename, int_varname, layer_varname
142  character(len=40) :: mod = "MOM_diag_remap" ! This module's name.
143  character(len=8) :: units, expected_units
144  character(len=34) :: longname, string2
145 
146  character(len=256) :: err_msg
147  logical :: ierr
148 
149  real, allocatable, dimension(:) :: interfaces, layers
150 
151  call initialize_regridding(remap_cs%regrid_cs, gv, gv%max_depth, param_file, mod, &
152  trim(remap_cs%vertical_coord_name), "DIAG_COORD", trim(remap_cs%diag_coord_name))
153  call set_regrid_params(remap_cs%regrid_cs, min_thickness=0., integrate_downward_for_e=.false.)
154 
155  remap_cs%nz = get_regrid_size(remap_cs%regrid_cs)
156 
157  if (remap_cs%vertical_coord == coordinatemode('SIGMA')) then
158  units = 'nondim'
159  longname = 'Fraction'
160  elseif (remap_cs%vertical_coord == coordinatemode('RHO')) then
161  units = 'kg m-3'
162  longname = 'Target Potential Density'
163  else
164  units = 'meters'
165  longname = 'Depth'
166  endif
167 
168  ! Make axes objects
169  allocate(interfaces(remap_cs%nz+1))
170  allocate(layers(remap_cs%nz))
171 
172  interfaces(:) = getcoordinateinterfaces(remap_cs%regrid_cs)
173  layers(:) = 0.5 * ( interfaces(1:remap_cs%nz) + interfaces(2:remap_cs%nz+1) )
174 
175  remap_cs%interface_axes_id = diag_axis_init(lowercase(trim(remap_cs%diag_coord_name))//'_i', &
176  interfaces, trim(units), 'z', &
177  trim(longname)//' at interface', direction=-1)
178  remap_cs%layer_axes_id = diag_axis_init(lowercase(trim(remap_cs%diag_coord_name))//'_l', &
179  layers, trim(units), 'z', &
180  trim(longname)//' at cell center', direction=-1, &
181  edges=remap_cs%interface_axes_id)
182 
183  ! Axes have now been configured.
184  remap_cs%configured = .true.
185 
186  deallocate(interfaces)
187  deallocate(layers)
188 
189 end subroutine diag_remap_configure_axes
190 
191 !> Get layer and interface axes ids for this coordinate
192 !! Needed when defining axes groups.
193 subroutine diag_remap_get_axes_info(remap_cs, nz, id_layer, id_interface)
194  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coordinate control structure
195  integer, intent(out) :: nz !< Number of vertical levels for the coordinate
196  integer, intent(out) :: id_layer !< 1D-axes id for layer points
197  integer, intent(out) :: id_interface !< 1D-axes id for interface points
198 
199  nz = remap_cs%nz
200  id_layer = remap_cs%layer_axes_id
201  id_interface = remap_cs%interface_axes_id
202 
203 end subroutine diag_remap_get_axes_info
204 
205 
206 !> Whether or not the axes for this vertical coordinated has been configured.
207 !! Configuration is complete when diag_remap_configure_axes() has been
208 !! successfully called.
209 function diag_remap_axes_configured(remap_cs)
210  type(diag_remap_ctrl), intent(in) :: remap_cs
211  logical :: diag_remap_axes_configured
212 
213  diag_remap_axes_configured = remap_cs%configured
214 
215 end function
216 
217 !> Build/update target vertical grids for diagnostic remapping.
218 !! \note The target grids need to be updated whenever sea surface
219 !! height or layer thicknesses changes. In the case of density-based
220 !! coordinates then technically we should also regenerate the
221 !! target grid whenever T/S change.
222 subroutine diag_remap_update(remap_cs, G, h, T, S, eqn_of_state)
223  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diagnostic coordinate control structure
224  type(ocean_grid_type), pointer :: G !< The ocean's grid type
225  real, dimension(:, :, :), intent(in) :: h, T, S !< New thickness, T and S
226  type(eos_type), pointer, intent(in) :: eqn_of_state !< A pointer to the equation of state
227 
228  ! Local variables
229  integer :: i, j, k, nz
230  real, dimension(remap_cs%nz + 1) :: zInterfaces
231  real, dimension(remap_cs%nz) :: resolution
232 
233  ! Note that coordinateMode('LAYER') is never 'configured' so will
234  ! always return here.
235  if (.not. remap_cs%configured) then
236  return
237  endif
238 
239  nz = remap_cs%nz
240 
241  if (.not. remap_cs%initialized) then
242  ! Initialize remapping and regridding on the first call
243  call initialize_remapping(remap_cs%remap_cs, 'PPM_IH4', boundary_extrapolation=.false.)
244  allocate(remap_cs%h(g%isd:g%ied,g%jsd:g%jed, nz))
245  remap_cs%initialized = .true.
246  endif
247 
248  ! Calculate remapping thicknesses for different target grids based on
249  ! nominal/target interface locations. This happens for every call on the
250  ! assumption that h, T, S has changed.
251  do j=g%jsd, g%jed
252  do i=g%isd, g%ied
253  if (g%mask2dT(i,j)==0.) then
254  remap_cs%h(i,j,:) = 0.
255  cycle
256  endif
257 
258  if (remap_cs%vertical_coord == coordinatemode('ZSTAR')) then
259  call build_zstar_column(get_zlike_cs(remap_cs%regrid_cs), nz, &
260  g%bathyT(i,j), sum(h(i,j,:)), zinterfaces)
261  elseif (remap_cs%vertical_coord == coordinatemode('SIGMA')) then
262  call build_sigma_column(get_sigma_cs(remap_cs%regrid_cs), nz, &
263  g%bathyT(i,j), sum(h(i,j,:)), zinterfaces)
264  elseif (remap_cs%vertical_coord == coordinatemode('RHO')) then
265  call build_rho_column(get_rho_cs(remap_cs%regrid_cs), remap_cs%remap_cs, nz, &
266  g%bathyT(i,j), h(i,j,:), t(i, j, :), s(i, j, :), &
267  eqn_of_state, zinterfaces)
268  elseif (remap_cs%vertical_coord == coordinatemode('SLIGHT')) then
269 ! call build_slight_column(remap_cs%regrid_cs,remap_cs%remap_cs, nz, &
270 ! G%bathyT(i,j), sum(h(i,j,:)), zInterfaces)
271  call mom_error(fatal,"diag_remap_update: SLIGHT coordinate not coded for diagnostics yet!")
272  elseif (remap_cs%vertical_coord == coordinatemode('HYCOM1')) then
273 ! call build_hycom1_column(remap_cs%regrid_cs, nz, &
274 ! G%bathyT(i,j), sum(h(i,j,:)), zInterfaces)
275  call mom_error(fatal,"diag_remap_update: HYCOM1 coordinate not coded for diagnostics yet!")
276  endif
277  remap_cs%h(i,j,:) = zinterfaces(1:nz) - zinterfaces(2:nz+1)
278  enddo
279  enddo
280 
281 end subroutine diag_remap_update
282 
283 !> Remap diagnostic field to alternative vertical grid.
284 subroutine diag_remap_do_remap(remap_cs, G, h, staggered_in_x, staggered_in_y, &
285  mask, missing_value, field, remapped_field)
286  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
287  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
288  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
289  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
290  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
291  real, dimension(:,:,:), pointer :: mask !< A mask for the field
292  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
293  real, dimension(:,:,:), intent(in) :: field(:,:,:) !< The diagnostic field to be remapped
294  real, dimension(:,:,:), intent(inout) :: remapped_field !< Field remapped to new coordinate
295  ! Local variables
296  real, dimension(remap_cs%nz) :: h_dest
297  real, dimension(size(h,3)) :: h_src
298  logical :: mask_vanished_layers
299  integer :: nz_src, nz_dest
300  integer :: i, j, k
301 
302  call assert(remap_cs%initialized, 'diag_remap_do_remap: remap_cs not initialized.')
303  call assert(size(field, 3) == size(h, 3), &
304  'diag_remap_do_remap: Remap field and thickness z-axes do not match.')
305 
306  nz_src = size(field,3)
307  nz_dest = remap_cs%nz
308  mask_vanished_layers = (remap_cs%vertical_coord == coordinatemode('ZSTAR'))
309  remapped_field(:,:,:) = missing_value
310 
311  if (staggered_in_x .and. .not. staggered_in_y) then
312  ! U-points
313  do j=g%jsc, g%jec
314  do i=g%iscB, g%iecB
315  if (associated(mask)) then
316  if (mask(i,j,1) == 0.) cycle
317  endif
318  h_src(:) = 0.5 * (h(i,j,:) + h(i+1,j,:))
319  h_dest(:) = 0.5 * (remap_cs%h(i,j,:) + remap_cs%h(i+1,j,:))
320  call remapping_core_h(remap_cs%remap_cs, nz_src, h_src(:), field(i,j,:), &
321  nz_dest, h_dest(:), remapped_field(i,j,:))
322  if (mask_vanished_layers) then ! This only works for z-like output
323  do k=1, nz_dest
324  if (h_dest(k) == 0.) remapped_field(i, j, k:nz_dest) = missing_value
325  enddo
326  endif
327  enddo
328  enddo
329  elseif (staggered_in_y .and. .not. staggered_in_x) then
330  ! V-points
331  do j=g%jscB, g%jecB
332  do i=g%isc, g%iec
333  if (associated(mask)) then
334  if (mask(i,j,1) == 0.) cycle
335  endif
336  h_src(:) = 0.5 * (h(i,j,:) + h(i,j+1,:))
337  h_dest(:) = 0.5 * (remap_cs%h(i,j,:) + remap_cs%h(i,j+1,:) )
338  call remapping_core_h(remap_cs%remap_cs, nz_src, h_src(:), field(i,j,:), &
339  nz_dest, h_dest(:), remapped_field(i,j,:))
340  if (mask_vanished_layers) then ! This only works for z-like output
341  do k=1, nz_dest
342  if (h_dest(k) == 0.) remapped_field(i,j,k) = missing_value
343  enddo
344  endif
345  enddo
346  enddo
347  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
348  ! H-points
349  do j=g%jsc, g%jec
350  do i=g%isc, g%iec
351  if (associated(mask)) then
352  if (mask(i,j, 1) == 0.) cycle
353  endif
354  h_dest(:) = remap_cs%h(i,j,:)
355  call remapping_core_h(remap_cs%remap_cs, nz_src, h(i,j,:), field(i,j,:), &
356  nz_dest, h_dest(:), remapped_field(i,j,:))
357  if (mask_vanished_layers) then ! This only works for z-like output
358  do k=1, nz_dest
359  if (h_dest(k)==0.) remapped_field(i,j,k) = missing_value
360  enddo
361  endif
362  enddo
363  enddo
364  else
365  call assert(.false., 'diag_remap_do_remap: Unsupported axis combination')
366  endif
367 
368 end subroutine diag_remap_do_remap
369 
370 !> Vertically re-grid an already vertically-integrated diagnostic field to alternative vertical grid.
371 subroutine vertically_reintegrate_diag_field(remap_cs, G, h, staggered_in_x, staggered_in_y, &
372  mask, missing_value, field, reintegrated_field)
373  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
374  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
375  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
376  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
377  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
378  real, dimension(:,:,:), pointer :: mask !< A mask for the field
379  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
380  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped
381  real, dimension(:,:,:), intent(inout) :: reintegrated_field !< Field argument remapped to alternative coordinate
382  ! Local variables
383  real, dimension(remap_cs%nz) :: h_dest
384  real, dimension(size(h,3)) :: h_src
385  integer :: nz_src, nz_dest
386  integer :: i, j, k
387 
388  call assert(remap_cs%initialized, 'vertically_reintegrate_diag_field: remap_cs not initialized.')
389  call assert(size(field, 3) == size(h, 3), &
390  'vertically_reintegrate_diag_field: Remap field and thickness z-axes do not match.')
391 
392  nz_src = size(field,3)
393  nz_dest = remap_cs%nz
394  reintegrated_field(:,:,:) = missing_value
395 
396  if (staggered_in_x .and. .not. staggered_in_y) then
397  ! U-points
398  do j=g%jsc, g%jec
399  do i=g%iscB, g%iecB
400  if (associated(mask)) then
401  if (mask(i,j,1) == 0.) cycle
402  endif
403  h_src(:) = 0.5 * (h(i,j,:) + h(i+1,j,:))
404  h_dest(:) = 0.5 * ( remap_cs%h(i,j,:) + remap_cs%h(i+1,j,:) )
405  call reintegrate_column(nz_src, h_src, field(i,j,:), &
406  nz_dest, h_dest, missing_value, reintegrated_field(i,j,:))
407  enddo
408  enddo
409  elseif (staggered_in_y .and. .not. staggered_in_x) then
410  ! V-points
411  do j=g%jscB, g%jecB
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,:) + h(i,j+1,:))
417  h_dest(:) = 0.5 * ( remap_cs%h(i,j,:) + remap_cs%h(i,j+1,:) )
418  call reintegrate_column(nz_src, h_src, field(i,j,:), &
419  nz_dest, h_dest, missing_value, reintegrated_field(i,j,:))
420  enddo
421  enddo
422  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
423  ! H-points
424  do j=g%jsc, g%jec
425  do i=g%isc, g%iec
426  if (associated(mask)) then
427  if (mask(i,j, 1) == 0.) cycle
428  endif
429  h_src(:) = h(i,j,:)
430  h_dest(:) = remap_cs%h(i,j,:)
431  call reintegrate_column(nz_src, h_src, field(i,j,:), &
432  nz_dest, h_dest, missing_value, reintegrated_field(i,j,:))
433  enddo
434  enddo
435  else
436  call assert(.false., 'vertically_reintegrate_diag_field: Q point remapping is not coded yet.')
437  endif
438 
440 
441 !> Vertically interpolate diagnostic field to alternative vertical grid.
442 subroutine vertically_interpolate_diag_field(remap_cs, G, h, staggered_in_x, staggered_in_y, &
443  mask, missing_value, field, interpolated_field)
444  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
445  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
446  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
447  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
448  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
449  real, dimension(:,:,:), pointer :: mask !< A mask for the field
450  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
451  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped
452  real, dimension(:,:,:), intent(inout) :: interpolated_field !< Field argument remapped to alternative coordinate
453  ! Local variables
454  real, dimension(remap_cs%nz) :: h_dest
455  real, dimension(size(h,3)) :: h_src
456  integer :: nz_src, nz_dest
457  integer :: i, j, k
458 
459  call assert(remap_cs%initialized, 'vertically_interpolate_diag_field: remap_cs not initialized.')
460  call assert(size(field, 3) == size(h, 3)+1, &
461  'vertically_interpolate_diag_field: Remap field and thickness z-axes do not match.')
462 
463  interpolated_field(:,:,:) = missing_value
464 
465  nz_src = size(h,3)
466  nz_dest = remap_cs%nz
467 
468  if (staggered_in_x .and. .not. staggered_in_y) then
469  ! U-points
470  do j=g%jsc, g%jec
471  do i=g%iscB, g%iecB
472  if (associated(mask)) then
473  if (mask(i,j,1) == 0.) cycle
474  endif
475  h_src(:) = 0.5 * (h(i,j,:) + h(i+1,j,:))
476  h_dest(:) = 0.5 * ( remap_cs%h(i,j,:) + remap_cs%h(i+1,j,:) )
477  call interpolate_column(nz_src, h_src, field(i,j,:), &
478  nz_dest, h_dest, missing_value, interpolated_field(i,j,:))
479  enddo
480  enddo
481  elseif (staggered_in_y .and. .not. staggered_in_x) then
482  ! V-points
483  do j=g%jscB, g%jecB
484  do i=g%isc, g%iec
485  if (associated(mask)) then
486  if (mask(i,j,1) == 0.) cycle
487  endif
488  h_src(:) = 0.5 * (h(i,j,:) + h(i,j+1,:))
489  h_dest(:) = 0.5 * ( remap_cs%h(i,j,:) + remap_cs%h(i,j+1,:) )
490  call interpolate_column(nz_src, h_src, field(i,j,:), &
491  nz_dest, h_dest, missing_value, interpolated_field(i,j,:))
492  enddo
493  enddo
494  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
495  ! H-points
496  do j=g%jsc, g%jec
497  do i=g%isc, g%iec
498  if (associated(mask)) then
499  if (mask(i,j, 1) == 0.) cycle
500  endif
501  h_src(:) = h(i,j,:)
502  h_dest(:) = remap_cs%h(i,j,:)
503  call interpolate_column(nz_src, h_src, field(i,j,:), &
504  nz_dest, h_dest, missing_value, interpolated_field(i,j,:))
505  enddo
506  enddo
507  else
508  call assert(.false., 'vertically_interpolate_diag_field: Q point remapping is not coded yet.')
509  endif
510 
512 
513 !> Horizontally average field
514 subroutine horizontally_average_diag_field(G, h, staggered_in_x, staggered_in_y, &
515  is_layer, is_extensive, &
516  missing_value, field, averaged_field)
517  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
518  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
519  logical, intent(in) :: staggered_in_x !< True if the x-axis location is at u or q points
520  logical, intent(in) :: staggered_in_y !< True if the y-axis location is at v or q points
521  logical, intent(in) :: is_layer !< True if the z-axis location is at h points
522  logical, intent(in) :: is_extensive !< True if the z-direction is spatially integrated (over layers)
523  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
524  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped
525  real, dimension(:), intent(inout) :: averaged_field !< Field argument horizontally averaged
526  ! Local variables
527  real, dimension(size(field, 3)) :: vol_sum, stuff_sum ! nz+1 is needed for interface averages
528  real :: v1, v2, total_volume, total_stuff, val
529  integer :: i, j, k, nz
530 
531  nz = size(field, 3)
532 
533  if (staggered_in_x .and. .not. staggered_in_y) then
534  if (is_layer) then
535  ! U-points
536  do k=1,nz
537  vol_sum(k) = 0.
538  stuff_sum(k) = 0.
539  if (is_extensive) then
540  do j=g%jsc, g%jec ; do i=g%isc, g%iec
541  v1 = g%areaCu(i,j)
542  v2 = g%areaCu(i-1,j)
543  vol_sum(k) = vol_sum(k) + 0.5 * ( v1 + v2 ) * g%mask2dT(i,j)
544  stuff_sum(k) = stuff_sum(k) + 0.5 * ( v1 * field(i,j,k) + v2 * field(i-1,j,k) ) * g%mask2dT(i,j)
545  enddo ; enddo
546  else ! Intensive
547  do j=g%jsc, g%jec ; do i=g%isc, g%iec
548  v1 = g%areaCu(i,j) * 0.5 * ( h(i,j,k) + h(i+1,j,k) )
549  v2 = g%areaCu(i-1,j) * 0.5 * ( h(i,j,k) + h(i-1,j,k) )
550  vol_sum(k) = vol_sum(k) + 0.5 * ( v1 + v2 ) * g%mask2dT(i,j)
551  stuff_sum(k) = stuff_sum(k) + 0.5 * ( v1 * field(i,j,k) + v2 * field(i-1,j,k) ) * g%mask2dT(i,j)
552  enddo ; enddo
553  endif
554  enddo
555  else ! Interface
556  do k=1,nz
557  vol_sum(k) = 0.
558  stuff_sum(k) = 0.
559  do j=g%jsc, g%jec ; do i=g%isc, g%iec
560  v1 = g%areaCu(i,j)
561  v2 = g%areaCu(i-1,j)
562  vol_sum(k) = vol_sum(k) + 0.5 * ( v1 + v2 ) * g%mask2dT(i,j)
563  stuff_sum(k) = stuff_sum(k) + 0.5 * ( v1 * field(i,j,k) + v2 * field(i-1,j,k) ) * g%mask2dT(i,j)
564  enddo ; enddo
565  enddo
566  endif
567  elseif (staggered_in_y .and. .not. staggered_in_x) then
568  if (is_layer) then
569  ! V-points
570  do k=1,nz
571  vol_sum(k) = 0.
572  stuff_sum(k) = 0.
573  if (is_extensive) then
574  do j=g%jsc, g%jec ; do i=g%isc, g%iec
575  v1 = g%areaCv(i,j)
576  v2 = g%areaCv(i,j-1)
577  vol_sum(k) = vol_sum(k) + 0.5 * ( v1 + v2 ) * g%mask2dT(i,j)
578  stuff_sum(k) = stuff_sum(k) + 0.5 * ( v1 * field(i,j,k) + v2 * field(i,j-1,k) ) * g%mask2dT(i,j)
579  enddo ; enddo
580  else ! Intensive
581  do j=g%jsc, g%jec ; do i=g%isc, g%iec
582  v1 = g%areaCv(i,j) * 0.5 * ( h(i,j,k) + h(i,j+1,k) )
583  v2 = g%areaCv(i,j-1) * 0.5 * ( h(i,j,k) + h(i,j-1,k) )
584  vol_sum(k) = vol_sum(k) + 0.5 * ( v1 + v2 ) * g%mask2dT(i,j)
585  stuff_sum(k) = stuff_sum(k) + 0.5 * ( v1 * field(i,j,k) + v2 * field(i,j-1,k) ) * g%mask2dT(i,j)
586  enddo ; enddo
587  endif
588  enddo
589  else ! Interface
590  do k=1,nz
591  vol_sum(k) = 0.
592  stuff_sum(k) = 0.
593  do j=g%jsc, g%jec ; do i=g%isc, g%iec
594  v1 = g%areaCv(i,j)
595  v2 = g%areaCv(i,j-1)
596  vol_sum(k) = vol_sum(k) + 0.5 * ( v1 + v2 ) * g%mask2dT(i,j)
597  stuff_sum(k) = stuff_sum(k) + 0.5 * ( v1 * field(i,j,k) + v2 * field(i,j-1,k) ) * g%mask2dT(i,j)
598  enddo ; enddo
599  enddo
600  endif
601  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
602  if (is_layer) then
603  ! H-points
604  do k=1,nz
605  vol_sum(k) = 0.
606  stuff_sum(k) = 0.
607  if (is_extensive) then
608  do j=g%jsc, g%jec ; do i=g%isc, g%iec
609  if (g%mask2dT(i,j)>0. .and. h(i,j,k)>0.) then
610  v1 = g%areaT(i,j)
611  vol_sum(k) = vol_sum(k) + v1
612  stuff_sum(k) = stuff_sum(k) + v1 * field(i,j,k)
613  endif
614  enddo ; enddo
615  else ! Intensive
616  do j=g%jsc, g%jec ; do i=g%isc, g%iec
617  if (g%mask2dT(i,j)>0. .and. h(i,j,k)>0.) then
618  v1 = g%areaT(i,j) * h(i,j,k)
619  vol_sum(k) = vol_sum(k) + v1
620  stuff_sum(k) = stuff_sum(k) + v1 * field(i,j,k)
621  endif
622  enddo ; enddo
623  endif
624  enddo
625  else ! Interface
626  do k=1,nz
627  vol_sum(k) = 0.
628  stuff_sum(k) = 0.
629  do j=g%jsc, g%jec ; do i=g%isc, g%iec
630  val = field(i,j,k)
631  if (g%mask2dT(i,j)>0. .and. val/=missing_value) then
632  v1 = g%areaT(i,j)
633  vol_sum(k) = vol_sum(k) + v1
634  stuff_sum(k) = stuff_sum(k) + v1 * field(i,j,k)
635  endif
636  enddo ; enddo
637  enddo
638  endif
639  else
640  call assert(.false., 'horizontally_average_diag_field: Q point averaging is not coded yet.')
641  endif
642 
643  call sum_across_pes(vol_sum, nz)
644  call sum_across_pes(stuff_sum, nz)
645 
646  do k=1,nz
647  if (vol_sum(k)>0.) then
648  averaged_field(k) = stuff_sum(k) / vol_sum(k)
649  else
650  averaged_field(k) = missing_value
651  endif
652  enddo
653 
654 end subroutine horizontally_average_diag_field
655 
656 end module mom_diag_remap
integer function coordinatemode(string)
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer...
Provides kernels for single-column interpolation, re-integration (re-mapping of integrated quantities...
real function, dimension(nk), public uniformresolution(nk, coordMode, maxDepth, rhoLight, rhoHeavy)
type(rho_cs) function, public get_rho_cs(CS)
subroutine, public diag_remap_set_active(remap_cs)
Indicate that this remapping type is actually used by the diag manager. If this is never called then ...
subroutine, public diag_remap_init(remap_cs, coord_tuple)
Initialize a diagnostic remapping type with the given vertical coordinate.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
logical function, public diag_remap_axes_configured(remap_cs)
Whether or not the axes for this vertical coordinated has been configured. Configuration is complete ...
type(zlike_cs) function, public get_zlike_cs(CS)
subroutine, public diag_remap_diag_registration_closed(remap_cs)
Inform that all diagnostics have been registered. If _set_active() has not been called on the remappi...
Provides the ocean grid type.
Definition: MOM_grid.F90:2
subroutine, public diag_remap_end(remap_cs)
De-init a diagnostic remapping type. Free allocated memory.
subroutine, public vertically_interpolate_diag_field(remap_cs, G, h, staggered_in_x, staggered_in_y, mask, missing_value, field, interpolated_field)
Vertically interpolate diagnostic field to alternative vertical grid.
Generates vertical grids as part of the ALE algorithm.
Provides column-wise vertical remapping functions.
This module contains I/O framework code.
Definition: MOM_io.F90:2
integer function, public get_regrid_size(CS)
Returns the number of levels/layers in the regridding control structure.
subroutine, public reintegrate_column(nsrc, h_src, uh_src, ndest, h_dest, missing_value, uh_dest)
Conservatively calculate integrated data, uh_dest, on grid h_dest, from layer-integrated data...
Container for remapping parameters.
subroutine, public diag_remap_get_axes_info(remap_cs, nz, id_layer, id_interface)
Get layer and interface axes ids for this coordinate Needed when defining axes groups.
character(len=len(input_string)) function, public lowercase(input_string)
subroutine, public build_rho_column(CS, remapCS, nz, depth, h, T, S, eqn_of_state, zInterface)
Definition: coord_rho.F90:86
Regridding control structure.
This type represents remapping of diagnostics to a particular vertical coordinate. There is one of these types for each vertical coordinate. The vertical axes of a diagnostic will reference an instance of this type indicating how (or if) the diagnostic should be vertically remapped when being posted.
subroutine, public set_regrid_params(CS, boundary_extrapolation, min_thickness, old_grid_weight, interp_scheme, depth_of_time_filter_shallow, depth_of_time_filter_deep, compress_fraction, dz_min_surface, nz_fixed_surface, Rho_ML_avg_depth, nlay_ML_to_interior, fix_haloclines, halocline_filt_len, halocline_strat_tol, integrate_downward_for_e, adaptTimeRatio, adaptZoom, adaptZoomCoeff, adaptBuoyCoeff, adaptAlpha, adaptDoMin)
Can be used to set any of the parameters for MOM_regridding.
subroutine, public assert(logical_arg, msg)
Issues a FATAL error if the assertion fails, i.e. the first argument is false.
character(len=120) function, public extractword(string, n)
Returns the string corresponding to the nth word in the argument or "" if the string is not long enou...
subroutine, public diag_remap_update(remap_cs, G, h, T, S, eqn_of_state)
Build/update target vertical grids for diagnostic remapping.
Regrid columns for the sigma coordinate.
Definition: coord_sigma.F90:2
real function, dimension(cs%nk+1), public getcoordinateinterfaces(CS)
Query the target coordinate interface positions.
subroutine, public build_zstar_column(CS, nz, depth, total_thickness, zInterface, z_rigid_top, eta_orig)
Builds a z* coordinate with a minimum thickness.
Definition: coord_zlike.F90:60
Regrid columns for a z-like coordinate (z-star, z-level)
Definition: coord_zlike.F90:2
subroutine, public horizontally_average_diag_field(G, h, staggered_in_x, staggered_in_y, is_layer, is_extensive, missing_value, field, averaged_field)
Horizontally average field.
This module is used for runtime remapping of diagnostics to z star, sigma and rho vertical coordinate...
type(sigma_cs) function, public get_sigma_cs(CS)
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
subroutine, public remapping_core_h(CS, n0, h0, u0, n1, h1, u1)
Remaps column of values u0 on grid h0 to grid h1 assuming the top edge is aligned.
subroutine, public diag_remap_do_remap(remap_cs, G, h, staggered_in_x, staggered_in_y, mask, missing_value, field, remapped_field)
Remap diagnostic field to alternative vertical grid.
Regrid columns for the continuous isopycnal (rho) coordinate.
Definition: coord_rho.F90:2
subroutine, public vertically_reintegrate_diag_field(remap_cs, G, h, staggered_in_x, staggered_in_y, mask, missing_value, field, reintegrated_field)
Vertically re-grid an already vertically-integrated diagnostic field to alternative vertical grid...
Contains constants for interpreting input parameters that control regridding.
subroutine, public initialize_remapping(CS, remapping_scheme, boundary_extrapolation, check_reconstruction, check_remapping, force_bounds_in_subcell)
Constructor for remapping control structure.
subroutine, public mom_error(level, message, all_print)
subroutine, public build_sigma_column(CS, nz, depth, totalThickness, zInterface)
Build a sigma coordinate column.
Definition: coord_sigma.F90:61
subroutine, public initialize_regridding(CS, GV, max_depth, param_file, mod, coord_mode, param_prefix, param_suffix)
Initialization and configures a regridding control structure based on customizable run-time parameter...
subroutine, public diag_remap_configure_axes(remap_cs, GV, param_file)
Configure the vertical axes for a diagnostic remapping control structure. Reads a configuration param...
A control structure for the equation of state.
Definition: MOM_EOS.F90:55
subroutine, public interpolate_column(nsrc, h_src, u_src, ndest, h_dest, missing_value, u_dest)
Linearly interpolate interface data, u_src, from grid h_src to a grid h_dest.