MOM6
mom_diag_remap Module Reference

Detailed Description

provides runtime remapping of diagnostics to z star, sigma and rho vertical coordinates.

The diag_remap_ctrl type represents a remapping of diagnostics to a particular vertical coordinate. The module is used by the diag mediator module in the following way:

  1. diag_remap_init() is called to initialize a diag_remap_ctrl instance.
  2. diag_remap_configure_axes() is called to read the configuration file and set up the vertical coordinate / axes definitions.
  3. diag_remap_get_axes_info() returns information needed for the diag mediator to define new axes for the remapped diagnostics.
  4. diag_remap_update() is called periodically (whenever h, T or S change) to either create or update the target remapping grids.
  5. diag_remap_do_remap() is called from within a diag post() to do the remapping before the diagnostic is written out.

Data Types

type  diag_remap_ctrl
 Represents remapping of diagnostics to a particular vertical coordinate. More...
 

Functions/Subroutines

subroutine, public diag_remap_init (remap_cs, coord_tuple)
 Initialize a diagnostic remapping type with the given vertical coordinate. More...
 
subroutine, public diag_remap_end (remap_cs)
 De-init a diagnostic remapping type. Free allocated memory. More...
 
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 remapping control structure will be disabled. This saves time in the case that a vertical coordinate was configured but no diagnostics which use the coordinate appeared in the diag_table. More...
 
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 the type will be disabled to save time. See further explanation with diag_remap_registration_closed. More...
 
subroutine, public diag_remap_configure_axes (remap_cs, GV, US, param_file)
 Configure the vertical axes for a diagnostic remapping control structure. Reads a configuration parameters to determine coordinate generation. More...
 
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. More...
 
logical function, public diag_remap_axes_configured (remap_cs)
 Whether or not the axes for this vertical coordinated has been configured. Configuration is complete when diag_remap_configure_axes() has been successfully called. More...
 
subroutine, public diag_remap_update (remap_cs, G, GV, US, h, T, S, eqn_of_state)
 Build/update target vertical grids for diagnostic remapping. More...
 
subroutine, public diag_remap_do_remap (remap_cs, G, GV, h, staggered_in_x, staggered_in_y, mask, missing_value, field, remapped_field)
 Remap diagnostic field to alternative vertical grid. More...
 
subroutine, public diag_remap_calc_hmask (remap_cs, G, mask)
 Calculate masks for target grid. More...
 
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. More...
 
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. More...
 
subroutine, public horizontally_average_diag_field (G, GV, h, staggered_in_x, staggered_in_y, is_layer, is_extensive, missing_value, field, averaged_field, averaged_mask)
 Horizontally average field. More...
 

Function/Subroutine Documentation

◆ diag_remap_axes_configured()

logical function, public mom_diag_remap::diag_remap_axes_configured ( type(diag_remap_ctrl), intent(in)  remap_cs)

Whether or not the axes for this vertical coordinated has been configured. Configuration is complete when diag_remap_configure_axes() has been successfully called.

Parameters
[in]remap_csDiagnostic coordinate control structure

Definition at line 255 of file MOM_diag_remap.F90.

255  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coordinate control structure
256  logical :: diag_remap_axes_configured
257 
258  diag_remap_axes_configured = remap_cs%configured
259 

◆ diag_remap_calc_hmask()

subroutine, public mom_diag_remap::diag_remap_calc_hmask ( type(diag_remap_ctrl), intent(in)  remap_cs,
type(ocean_grid_type), intent(in)  G,
real, dimension(:,:,:), intent(out)  mask 
)

Calculate masks for target grid.

Parameters
[in]remap_csDiagnostic coodinate control structure
[in]gOcean grid structure
[out]maskh-point mask for target grid

Definition at line 435 of file MOM_diag_remap.F90.

435  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
436  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
437  real, dimension(:,:,:), intent(out) :: mask !< h-point mask for target grid
438  ! Local variables
439  real, dimension(remap_cs%nz) :: h_dest
440  integer :: i, j, k
441  logical :: mask_vanished_layers
442  real :: h_tot, h_err
443 
444  call assert(remap_cs%initialized, 'diag_remap_calc_hmask: remap_cs not initialized.')
445 
446  ! Only z*-like diagnostic coordinates should have a 3d mask
447  mask_vanished_layers = (remap_cs%vertical_coord == coordinatemode('ZSTAR'))
448  mask(:,:,:) = 0.
449 
450  do j=g%jsc-1, g%jec+1 ; do i=g%isc-1, g%iec+1
451  if (g%mask2dT(i,j)>0.) then
452  if (mask_vanished_layers) then
453  h_dest(:) = remap_cs%h(i,j,:)
454  h_tot = 0.
455  h_err = 0.
456  do k=1, remap_cs%nz
457  h_tot = h_tot + h_dest(k)
458  ! This is an overestimate of how thick a vanished layer might be, that
459  ! appears due to round-off.
460  h_err = h_err + epsilon(h_tot) * h_tot
461  ! Mask out vanished layers
462  if (h_dest(k)<=8.*h_err) then
463  mask(i,j,k) = 0.
464  else
465  mask(i,j,k) = 1.
466  endif
467  enddo
468  else ! all layers might contain data
469  mask(i,j,:) = 1.
470  endif
471  endif
472  enddo ; enddo
473 

References mom_error_handler::assert(), and regrid_consts::coordinatemode().

Here is the call graph for this function:

◆ diag_remap_configure_axes()

subroutine, public mom_diag_remap::diag_remap_configure_axes ( type(diag_remap_ctrl), intent(inout)  remap_cs,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file 
)

Configure the vertical axes for a diagnostic remapping control structure. Reads a configuration parameters to determine coordinate generation.

Parameters
[in,out]remap_csDiag remap control structure
[in]gvocean vertical grid structure
[in]usA dimensional unit scaling type
[in]param_fileParameter file structure

Definition at line 180 of file MOM_diag_remap.F90.

180  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remap control structure
181  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
182  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
183  type(param_file_type), intent(in) :: param_file !< Parameter file structure
184  ! Local variables
185  integer :: nzi(4), nzl(4), k
186  character(len=200) :: inputdir, string, filename, int_varname, layer_varname
187  character(len=40) :: mod = "MOM_diag_remap" ! This module's name.
188  character(len=8) :: units, expected_units
189  character(len=34) :: longname, string2
190 
191  character(len=256) :: err_msg
192  logical :: ierr
193 
194  real, allocatable, dimension(:) :: interfaces, layers
195 
196  call initialize_regridding(remap_cs%regrid_cs, gv, us, gv%max_depth, param_file, mod, &
197  trim(remap_cs%vertical_coord_name), "DIAG_COORD", trim(remap_cs%diag_coord_name))
198  call set_regrid_params(remap_cs%regrid_cs, min_thickness=0., integrate_downward_for_e=.false.)
199 
200  remap_cs%nz = get_regrid_size(remap_cs%regrid_cs)
201 
202  if (remap_cs%vertical_coord == coordinatemode('SIGMA')) then
203  units = 'nondim'
204  longname = 'Fraction'
205  elseif (remap_cs%vertical_coord == coordinatemode('RHO')) then
206  units = 'kg m-3'
207  longname = 'Target Potential Density'
208  else
209  units = 'meters'
210  longname = 'Depth'
211  endif
212 
213  ! Make axes objects
214  allocate(interfaces(remap_cs%nz+1))
215  allocate(layers(remap_cs%nz))
216 
217  interfaces(:) = getcoordinateinterfaces(remap_cs%regrid_cs)
218  layers(:) = 0.5 * ( interfaces(1:remap_cs%nz) + interfaces(2:remap_cs%nz+1) )
219 
220  remap_cs%interface_axes_id = diag_axis_init(lowercase(trim(remap_cs%diag_coord_name))//'_i', &
221  interfaces, trim(units), 'z', &
222  trim(longname)//' at interface', direction=-1)
223  remap_cs%layer_axes_id = diag_axis_init(lowercase(trim(remap_cs%diag_coord_name))//'_l', &
224  layers, trim(units), 'z', &
225  trim(longname)//' at cell center', direction=-1, &
226  edges=remap_cs%interface_axes_id)
227 
228  ! Axes have now been configured.
229  remap_cs%configured = .true.
230 
231  deallocate(interfaces)
232  deallocate(layers)
233 

References regrid_consts::coordinatemode(), and mom_string_functions::lowercase().

Here is the call graph for this function:

◆ diag_remap_diag_registration_closed()

subroutine, public mom_diag_remap::diag_remap_diag_registration_closed ( type(diag_remap_ctrl), intent(inout)  remap_cs)

Inform that all diagnostics have been registered. If _set_active() has not been called on the remapping control structure will be disabled. This saves time in the case that a vertical coordinate was configured but no diagnostics which use the coordinate appeared in the diag_table.

Parameters
[in,out]remap_csDiag remapping control structure

Definition at line 159 of file MOM_diag_remap.F90.

159  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
160 
161  if (.not. remap_cs%used) then
162  call diag_remap_end(remap_cs)
163  endif
164 

References diag_remap_end().

Here is the call graph for this function:

◆ diag_remap_do_remap()

subroutine, public mom_diag_remap::diag_remap_do_remap ( type(diag_remap_ctrl), intent(in)  remap_cs,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(:,:,:), intent(in)  h,
logical, intent(in)  staggered_in_x,
logical, intent(in)  staggered_in_y,
real, dimension(:,:,:), pointer  mask,
real, intent(in)  missing_value,
real, dimension(:,:,:), intent(in)  field,
real, dimension(:,:,:), intent(inout)  remapped_field 
)

Remap diagnostic field to alternative vertical grid.

Parameters
[in]remap_csDiagnostic coodinate control structure
[in]gOcean grid structure
[in]gvocean vertical grid structure
[in]hThe current thicknesses
[in]staggered_in_xTrue is the x-axis location is at u or q points
[in]staggered_in_yTrue is the y-axis location is at v or q points
maskA mask for the field
[in]missing_valueA missing_value to assign land/vanished points
[in]fieldThe diagnostic field to be remapped
[in,out]remapped_fieldField remapped to new coordinate

Definition at line 340 of file MOM_diag_remap.F90.

340  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
341  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
342  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
343  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
344  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
345  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
346  real, dimension(:,:,:), pointer :: mask !< A mask for the field
347  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
348  real, dimension(:,:,:), intent(in) :: field(:,:,:) !< The diagnostic field to be remapped
349  real, dimension(:,:,:), intent(inout) :: remapped_field !< Field remapped to new coordinate
350  ! Local variables
351  real, dimension(remap_cs%nz) :: h_dest
352  real, dimension(size(h,3)) :: h_src
353  real :: h_neglect, h_neglect_edge
354  integer :: nz_src, nz_dest
355  integer :: i, j, k !< Grid index
356  integer :: i1, j1 !< 1-based index
357  integer :: i_lo, i_hi, j_lo, j_hi !< (uv->h) interpolation indices
358  integer :: shift !< Symmetric offset for 1-based indexing
359 
360  call assert(remap_cs%initialized, 'diag_remap_do_remap: remap_cs not initialized.')
361  call assert(size(field, 3) == size(h, 3), &
362  'diag_remap_do_remap: Remap field and thickness z-axes do not match.')
363 
364  !### Try replacing both of these with GV%H_subroundoff
365  if (gv%Boussinesq) then
366  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
367  else
368  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
369  endif
370 
371  nz_src = size(field,3)
372  nz_dest = remap_cs%nz
373  remapped_field(:,:,:) = 0.
374 
375  ! Symmetric grid offset under 1-based indexing; see header for details.
376  shift = 0; if (g%symmetric) shift = 1
377 
378  if (staggered_in_x .and. .not. staggered_in_y) then
379  ! U-points
380  do j=g%jsc, g%jec
381  do i=g%iscB, g%iecB
382  i1 = i - g%isdB + 1
383  i_lo = i1 - shift; i_hi = i_lo + 1
384  if (associated(mask)) then
385  if (mask(i,j,1) == 0.) cycle
386  endif
387  h_src(:) = 0.5 * (h(i_lo,j,:) + h(i_hi,j,:))
388  h_dest(:) = 0.5 * (remap_cs%h(i_lo,j,:) + remap_cs%h(i_hi,j,:))
389  call remapping_core_h(remap_cs%remap_cs, &
390  nz_src, h_src(:), field(i1,j,:), &
391  nz_dest, h_dest(:), remapped_field(i1,j,:), &
392  h_neglect, h_neglect_edge)
393  enddo
394  enddo
395  elseif (staggered_in_y .and. .not. staggered_in_x) then
396  ! V-points
397  do j=g%jscB, g%jecB
398  j1 = j - g%jsdB + 1
399  j_lo = j1 - shift; j_hi = j_lo + 1
400  do i=g%isc, g%iec
401  if (associated(mask)) then
402  if (mask(i,j,1) == 0.) cycle
403  endif
404  h_src(:) = 0.5 * (h(i,j_lo,:) + h(i,j_hi,:))
405  h_dest(:) = 0.5 * (remap_cs%h(i,j_lo,:) + remap_cs%h(i,j_hi,:))
406  call remapping_core_h(remap_cs%remap_cs, &
407  nz_src, h_src(:), field(i,j1,:), &
408  nz_dest, h_dest(:), remapped_field(i,j1,:), &
409  h_neglect, h_neglect_edge)
410  enddo
411  enddo
412  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
413  ! H-points
414  do j=g%jsc, g%jec
415  do i=g%isc, g%iec
416  if (associated(mask)) then
417  if (mask(i,j,1) == 0.) cycle
418  endif
419  h_src(:) = h(i,j,:)
420  h_dest(:) = remap_cs%h(i,j,:)
421  call remapping_core_h(remap_cs%remap_cs, &
422  nz_src, h_src(:), field(i,j,:), &
423  nz_dest, h_dest(:), remapped_field(i,j,:), &
424  h_neglect, h_neglect_edge)
425  enddo
426  enddo
427  else
428  call assert(.false., 'diag_remap_do_remap: Unsupported axis combination')
429  endif
430 

References mom_error_handler::assert(), and mom_remapping::remapping_core_h().

Here is the call graph for this function:

◆ diag_remap_end()

subroutine, public mom_diag_remap::diag_remap_end ( type(diag_remap_ctrl), intent(inout)  remap_cs)

De-init a diagnostic remapping type. Free allocated memory.

Parameters
[in,out]remap_csDiag remapping control structure

Definition at line 142 of file MOM_diag_remap.F90.

142  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
143 
144  if (allocated(remap_cs%h)) deallocate(remap_cs%h)
145  if (allocated(remap_cs%dz)) deallocate(remap_cs%dz)
146  remap_cs%configured = .false.
147  remap_cs%initialized = .false.
148  remap_cs%used = .false.
149  remap_cs%nz = 0
150 

Referenced by diag_remap_diag_registration_closed().

Here is the caller graph for this function:

◆ diag_remap_get_axes_info()

subroutine, public mom_diag_remap::diag_remap_get_axes_info ( type(diag_remap_ctrl), intent(in)  remap_cs,
integer, intent(out)  nz,
integer, intent(out)  id_layer,
integer, intent(out)  id_interface 
)

Get layer and interface axes ids for this coordinate Needed when defining axes groups.

Parameters
[in]remap_csDiagnostic coordinate control structure
[out]nzNumber of vertical levels for the coordinate
[out]id_layer1D-axes id for layer points
[out]id_interface1D-axes id for interface points

Definition at line 239 of file MOM_diag_remap.F90.

239  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coordinate control structure
240  integer, intent(out) :: nz !< Number of vertical levels for the coordinate
241  integer, intent(out) :: id_layer !< 1D-axes id for layer points
242  integer, intent(out) :: id_interface !< 1D-axes id for interface points
243 
244  nz = remap_cs%nz
245  id_layer = remap_cs%layer_axes_id
246  id_interface = remap_cs%interface_axes_id
247 

◆ diag_remap_init()

subroutine, public mom_diag_remap::diag_remap_init ( type(diag_remap_ctrl), intent(inout)  remap_cs,
character(len=*), intent(in)  coord_tuple 
)

Initialize a diagnostic remapping type with the given vertical coordinate.

Parameters
[in,out]remap_csDiag remapping control structure
[in]coord_tupleA string in form of MODULE_SUFFIX PARAMETER_SUFFIX COORDINATE_NAME

Definition at line 124 of file MOM_diag_remap.F90.

124  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
125  character(len=*), intent(in) :: coord_tuple !< A string in form of
126  !! MODULE_SUFFIX PARAMETER_SUFFIX COORDINATE_NAME
127 
128  remap_cs%diag_module_suffix = trim(extractword(coord_tuple, 1))
129  remap_cs%diag_coord_name = trim(extractword(coord_tuple, 2))
130  remap_cs%vertical_coord_name = trim(extractword(coord_tuple, 3))
131  remap_cs%vertical_coord = coordinatemode(remap_cs%vertical_coord_name)
132  remap_cs%configured = .false.
133  remap_cs%initialized = .false.
134  remap_cs%used = .false.
135  remap_cs%nz = 0
136 

References regrid_consts::coordinatemode(), and mom_string_functions::extractword().

Here is the call graph for this function:

◆ diag_remap_set_active()

subroutine, public mom_diag_remap::diag_remap_set_active ( type(diag_remap_ctrl), intent(inout)  remap_cs)

Indicate that this remapping type is actually used by the diag manager. If this is never called then the type will be disabled to save time. See further explanation with diag_remap_registration_closed.

Parameters
[in,out]remap_csDiag remapping control structure

Definition at line 171 of file MOM_diag_remap.F90.

171  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diag remapping control structure
172 
173  remap_cs%used = .true.
174 

◆ diag_remap_update()

subroutine, public mom_diag_remap::diag_remap_update ( type(diag_remap_ctrl), intent(inout)  remap_cs,
type(ocean_grid_type), pointer  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension(:, :, :), intent(in)  h,
real, dimension(:, :, :), intent(in)  T,
real, dimension(:, :, :), intent(in)  S,
type(eos_type), pointer  eqn_of_state 
)

Build/update target vertical grids for diagnostic remapping.

Note
The target grids need to be updated whenever sea surface height or layer thicknesses changes. In the case of density-based coordinates then technically we should also regenerate the target grid whenever T/S change.
Parameters
[in,out]remap_csDiagnostic coordinate control structure
gThe ocean's grid type
[in]gvocean vertical grid structure
[in]usA dimensional unit scaling type
[in]hNew thickness
[in]tNew T
[in]sNew S
eqn_of_stateA pointer to the equation of state

Definition at line 268 of file MOM_diag_remap.F90.

268  type(diag_remap_ctrl), intent(inout) :: remap_cs !< Diagnostic coordinate control structure
269  type(ocean_grid_type), pointer :: G !< The ocean's grid type
270  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
271  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
272  real, dimension(:, :, :), intent(in) :: h !< New thickness
273  real, dimension(:, :, :), intent(in) :: T !< New T
274  real, dimension(:, :, :), intent(in) :: S !< New S
275  type(EOS_type), pointer :: eqn_of_state !< A pointer to the equation of state
276 
277  ! Local variables
278  real, dimension(remap_cs%nz + 1) :: zInterfaces
279  real :: h_neglect, h_neglect_edge
280  integer :: i, j, k, nz
281 
282  ! Note that coordinateMode('LAYER') is never 'configured' so will
283  ! always return here.
284  if (.not. remap_cs%configured) then
285  return
286  endif
287 
288  !### Try replacing both of these with GV%H_subroundoff
289  if (gv%Boussinesq) then
290  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
291  else
292  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
293  endif
294  nz = remap_cs%nz
295 
296  if (.not. remap_cs%initialized) then
297  ! Initialize remapping and regridding on the first call
298  call initialize_remapping(remap_cs%remap_cs, 'PPM_IH4', boundary_extrapolation=.false.)
299  allocate(remap_cs%h(g%isd:g%ied,g%jsd:g%jed, nz))
300  remap_cs%initialized = .true.
301  endif
302 
303  ! Calculate remapping thicknesses for different target grids based on
304  ! nominal/target interface locations. This happens for every call on the
305  ! assumption that h, T, S has changed.
306  do j=g%jsc-1, g%jec+1 ; do i=g%isc-1, g%iec+1
307  if (g%mask2dT(i,j)==0.) then
308  remap_cs%h(i,j,:) = 0.
309  cycle
310  endif
311 
312  if (remap_cs%vertical_coord == coordinatemode('ZSTAR')) then
313  call build_zstar_column(get_zlike_cs(remap_cs%regrid_cs), &
314  gv%Z_to_H*g%bathyT(i,j), sum(h(i,j,:)), &
315  zinterfaces, zscale=gv%Z_to_H)
316  elseif (remap_cs%vertical_coord == coordinatemode('SIGMA')) then
317  call build_sigma_column(get_sigma_cs(remap_cs%regrid_cs), &
318  gv%Z_to_H*g%bathyT(i,j), sum(h(i,j,:)), zinterfaces)
319  elseif (remap_cs%vertical_coord == coordinatemode('RHO')) then
320  call build_rho_column(get_rho_cs(remap_cs%regrid_cs), g%ke, &
321  us%Z_to_m*g%bathyT(i,j), h(i,j,:), t(i,j,:), s(i,j,:), &
322  eqn_of_state, zinterfaces, h_neglect, h_neglect_edge)
323  elseif (remap_cs%vertical_coord == coordinatemode('SLIGHT')) then
324 ! call build_slight_column(remap_cs%regrid_cs,remap_cs%remap_cs, nz, &
325 ! US%Z_to_m*G%bathyT(i,j), sum(h(i,j,:)), zInterfaces)
326  call mom_error(fatal,"diag_remap_update: SLIGHT coordinate not coded for diagnostics yet!")
327  elseif (remap_cs%vertical_coord == coordinatemode('HYCOM1')) then
328 ! call build_hycom1_column(remap_cs%regrid_cs, nz, &
329 ! US%Z_to_m*G%bathyT(i,j), sum(h(i,j,:)), zInterfaces)
330  call mom_error(fatal,"diag_remap_update: HYCOM1 coordinate not coded for diagnostics yet!")
331  endif
332  remap_cs%h(i,j,:) = zinterfaces(1:nz) - zinterfaces(2:nz+1)
333  enddo ; enddo
334 

References coord_rho::build_rho_column(), coord_sigma::build_sigma_column(), coord_zlike::build_zstar_column(), regrid_consts::coordinatemode(), mom_regridding::get_rho_cs(), mom_regridding::get_sigma_cs(), mom_regridding::get_zlike_cs(), and mom_error_handler::mom_error().

Here is the call graph for this function:

◆ horizontally_average_diag_field()

subroutine, public mom_diag_remap::horizontally_average_diag_field ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(:,:,:), intent(in)  h,
logical, intent(in)  staggered_in_x,
logical, intent(in)  staggered_in_y,
logical, intent(in)  is_layer,
logical, intent(in)  is_extensive,
real, intent(in)  missing_value,
real, dimension(:,:,:), intent(in)  field,
real, dimension(:), intent(inout)  averaged_field,
logical, dimension(:), intent(inout)  averaged_mask 
)

Horizontally average field.

Parameters
[in]gOcean grid structure
[in]gvThe ocean vertical grid structure
[in]hThe current thicknesses
[in]staggered_in_xTrue if the x-axis location is at u or q points
[in]staggered_in_yTrue if the y-axis location is at v or q points
[in]is_layerTrue if the z-axis location is at h points
[in]is_extensiveTrue if the z-direction is spatially integrated (over layers)
[in]missing_valueA missing_value to assign land/vanished points
[in]fieldThe diagnostic field to be remapped
[in,out]averaged_fieldField argument horizontally averaged
[in,out]averaged_maskMask for horizontally averaged field

Definition at line 644 of file MOM_diag_remap.F90.

644  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
645  type(verticalGrid_type), intent(in) :: GV !< The ocean vertical grid structure
646  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
647  logical, intent(in) :: staggered_in_x !< True if the x-axis location is at u or q points
648  logical, intent(in) :: staggered_in_y !< True if the y-axis location is at v or q points
649  logical, intent(in) :: is_layer !< True if the z-axis location is at h points
650  logical, intent(in) :: is_extensive !< True if the z-direction is spatially integrated (over layers)
651  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
652  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped
653  real, dimension(:), intent(inout) :: averaged_field !< Field argument horizontally averaged
654  logical, dimension(:), intent(inout) :: averaged_mask !< Mask for horizontally averaged field
655 
656  ! Local variables
657  real, dimension(G%isc:G%iec, G%jsc:G%jec, size(field,3)) :: volume, stuff
658  real, dimension(size(field, 3)) :: vol_sum, stuff_sum ! nz+1 is needed for interface averages
659  real :: height
660  integer :: i, j, k, nz
661  integer :: i1, j1 !< 1-based index
662 
663  nz = size(field, 3)
664 
665  ! TODO: These averages could potentially be modified to use the function in
666  ! the MOM_spatial_means module.
667  ! NOTE: Reproducible sums must be computed in the original MKS units
668 
669  if (staggered_in_x .and. .not. staggered_in_y) then
670  if (is_layer) then
671  ! U-points
672  do k=1,nz
673  vol_sum(k) = 0.
674  stuff_sum(k) = 0.
675  if (is_extensive) then
676  do j=g%jsc, g%jec ; do i=g%isc, g%iec
677  i1 = i - g%isdB + 1
678  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCu(i,j)) * g%mask2dCu(i,j)
679  stuff(i,j,k) = volume(i,j,k) * field(i1,j,k)
680  enddo ; enddo
681  else ! Intensive
682  do j=g%jsc, g%jec ; do i=g%isc, g%iec
683  i1 = i - g%isdB + 1
684  height = 0.5 * (h(i,j,k) + h(i+1,j,k))
685  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCu(i,j)) &
686  * (gv%H_to_m * height) * g%mask2dCu(i,j)
687  stuff(i,j,k) = volume(i,j,k) * field(i1,j,k)
688  enddo ; enddo
689  endif
690  enddo
691  else ! Interface
692  do k=1,nz
693  do j=g%jsc, g%jec ; do i=g%isc, g%iec
694  i1 = i - g%isdB + 1
695  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCu(i,j)) * g%mask2dCu(i,j)
696  stuff(i,j,k) = volume(i,j,k) * field(i1,j,k)
697  enddo ; enddo
698  enddo
699  endif
700  elseif (staggered_in_y .and. .not. staggered_in_x) then
701  if (is_layer) then
702  ! V-points
703  do k=1,nz
704  if (is_extensive) then
705  do j=g%jsc, g%jec ; do i=g%isc, g%iec
706  j1 = j - g%jsdB + 1
707  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCv(i,j)) * g%mask2dCv(i,j)
708  stuff(i,j,k) = volume(i,j,k) * field(i,j1,k)
709  enddo ; enddo
710  else ! Intensive
711  do j=g%jsc, g%jec ; do i=g%isc, g%iec
712  j1 = j - g%jsdB + 1
713  height = 0.5 * (h(i,j,k) + h(i,j+1,k))
714  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCv(i,j)) &
715  * (gv%H_to_m * height) * g%mask2dCv(i,j)
716  stuff(i,j,k) = volume(i,j,k) * field(i,j1,k)
717  enddo ; enddo
718  endif
719  enddo
720  else ! Interface
721  do k=1,nz
722  do j=g%jsc, g%jec ; do i=g%isc, g%iec
723  j1 = j - g%jsdB + 1
724  volume(i,j,k) = (g%US%L_to_m**2 * g%areaCv(i,j)) * g%mask2dCv(i,j)
725  stuff(i,j,k) = volume(i,j,k) * field(i,j1,k)
726  enddo ; enddo
727  enddo
728  endif
729  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
730  if (is_layer) then
731  ! H-points
732  do k=1,nz
733  if (is_extensive) then
734  do j=g%jsc, g%jec ; do i=g%isc, g%iec
735  if (h(i,j,k) > 0.) then
736  volume(i,j,k) = (g%US%L_to_m**2 * g%areaT(i,j)) * g%mask2dT(i,j)
737  stuff(i,j,k) = volume(i,j,k) * field(i,j,k)
738  else
739  volume(i,j,k) = 0.
740  stuff(i,j,k) = 0.
741  endif
742  enddo ; enddo
743  else ! Intensive
744  do j=g%jsc, g%jec ; do i=g%isc, g%iec
745  volume(i,j,k) = (g%US%L_to_m**2 * g%areaT(i,j)) &
746  * (gv%H_to_m * h(i,j,k)) * g%mask2dT(i,j)
747  stuff(i,j,k) = volume(i,j,k) * field(i,j,k)
748  enddo ; enddo
749  endif
750  enddo
751  else ! Interface
752  do k=1,nz
753  do j=g%jsc, g%jec ; do i=g%isc, g%iec
754  volume(i,j,k) = (g%US%L_to_m**2 * g%areaT(i,j)) * g%mask2dT(i,j)
755  stuff(i,j,k) = volume(i,j,k) * field(i,j,k)
756  enddo ; enddo
757  enddo
758  endif
759  else
760  call assert(.false., 'horizontally_average_diag_field: Q point averaging is not coded yet.')
761  endif
762 
763  do k = 1,nz
764  vol_sum(k) = reproducing_sum(volume(:,:,k))
765  stuff_sum(k) = reproducing_sum(stuff(:,:,k))
766  enddo
767 
768  averaged_mask(:) = .true.
769  do k=1,nz
770  if (vol_sum(k) > 0.) then
771  averaged_field(k) = stuff_sum(k) / vol_sum(k)
772  else
773  averaged_field(k) = 0.
774  averaged_mask(k) = .false.
775  endif
776  enddo
777 

References mom_error_handler::assert().

Referenced by mom_diag_mediator::post_xy_average().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ vertically_interpolate_diag_field()

subroutine, public mom_diag_remap::vertically_interpolate_diag_field ( type(diag_remap_ctrl), intent(in)  remap_cs,
type(ocean_grid_type), intent(in)  G,
real, dimension(:,:,:), intent(in)  h,
logical, intent(in)  staggered_in_x,
logical, intent(in)  staggered_in_y,
real, dimension(:,:,:), pointer  mask,
real, intent(in)  missing_value,
real, dimension(:,:,:), intent(in)  field,
real, dimension(:,:,:), intent(inout)  interpolated_field 
)

Vertically interpolate diagnostic field to alternative vertical grid.

Parameters
[in]remap_csDiagnostic coodinate control structure
[in]gOcean grid structure
[in]hThe current thicknesses
[in]staggered_in_xTrue is the x-axis location is at u or q points
[in]staggered_in_yTrue is the y-axis location is at v or q points
maskA mask for the field
[in]missing_valueA missing_value to assign land/vanished points
[in]fieldThe diagnostic field to be remapped
[in,out]interpolated_fieldField argument remapped to alternative coordinate

Definition at line 560 of file MOM_diag_remap.F90.

560  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
561  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
562  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
563  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
564  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
565  real, dimension(:,:,:), pointer :: mask !< A mask for the field
566  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
567  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped
568  real, dimension(:,:,:), intent(inout) :: interpolated_field !< Field argument remapped to alternative coordinate
569  ! Local variables
570  real, dimension(remap_cs%nz) :: h_dest
571  real, dimension(size(h,3)) :: h_src
572  integer :: nz_src, nz_dest
573  integer :: i, j, k !< Grid index
574  integer :: i1, j1 !< 1-based index
575  integer :: i_lo, i_hi, j_lo, j_hi !< (uv->h) interpolation indices
576  integer :: shift !< Symmetric offset for 1-based indexing
577 
578  call assert(remap_cs%initialized, 'vertically_interpolate_diag_field: remap_cs not initialized.')
579  call assert(size(field, 3) == size(h, 3)+1, &
580  'vertically_interpolate_diag_field: Remap field and thickness z-axes do not match.')
581 
582  interpolated_field(:,:,:) = 0.
583 
584  nz_src = size(h,3)
585  nz_dest = remap_cs%nz
586 
587  ! Symmetric grid offset under 1-based indexing; see header for details.
588  shift = 0; if (g%symmetric) shift = 1
589 
590  if (staggered_in_x .and. .not. staggered_in_y) then
591  ! U-points
592  do j=g%jsc, g%jec
593  do i=g%iscB, g%iecB
594  i1 = i - g%isdB + 1
595  i_lo = i1 - shift; i_hi = i_lo + 1
596  if (associated(mask)) then
597  if (mask(i,j,1) == 0.) cycle
598  endif
599  h_src(:) = 0.5 * (h(i_lo,j,:) + h(i_hi,j,:))
600  h_dest(:) = 0.5 * (remap_cs%h(i_lo,j,:) + remap_cs%h(i_hi,j,:))
601  call interpolate_column(nz_src, h_src, field(i1,j,:), &
602  nz_dest, h_dest, 0., interpolated_field(i1,j,:))
603  enddo
604  enddo
605  elseif (staggered_in_y .and. .not. staggered_in_x) then
606  ! V-points
607  do j=g%jscB, g%jecB
608  j1 = j - g%jsdB + 1
609  j_lo = j1 - shift; j_hi = j_lo + 1
610  do i=g%isc, g%iec
611  if (associated(mask)) then
612  if (mask(i,j,1) == 0.) cycle
613  endif
614  h_src(:) = 0.5 * (h(i,j_lo,:) + h(i,j_hi,:))
615  h_dest(:) = 0.5 * (remap_cs%h(i,j_lo,:) + remap_cs%h(i,j_hi,:))
616  call interpolate_column(nz_src, h_src, field(i,j1,:), &
617  nz_dest, h_dest, 0., interpolated_field(i,j1,:))
618  enddo
619  enddo
620  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
621  ! H-points
622  do j=g%jsc, g%jec
623  do i=g%isc, g%iec
624  if (associated(mask)) then
625  if (mask(i,j,1) == 0.) cycle
626  endif
627  h_src(:) = h(i,j,:)
628  h_dest(:) = remap_cs%h(i,j,:)
629  call interpolate_column(nz_src, h_src, field(i,j,:), &
630  nz_dest, h_dest, 0., interpolated_field(i,j,:))
631  enddo
632  enddo
633  else
634  call assert(.false., 'vertically_interpolate_diag_field: Q point remapping is not coded yet.')
635  endif
636 

References mom_error_handler::assert(), and mom_diag_vkernels::interpolate_column().

Here is the call graph for this function:

◆ vertically_reintegrate_diag_field()

subroutine, public mom_diag_remap::vertically_reintegrate_diag_field ( type(diag_remap_ctrl), intent(in)  remap_cs,
type(ocean_grid_type), intent(in)  G,
real, dimension(:,:,:), intent(in)  h,
logical, intent(in)  staggered_in_x,
logical, intent(in)  staggered_in_y,
real, dimension(:,:,:), pointer  mask,
real, intent(in)  missing_value,
real, dimension(:,:,:), intent(in)  field,
real, dimension(:,:,:), intent(inout)  reintegrated_field 
)

Vertically re-grid an already vertically-integrated diagnostic field to alternative vertical grid.

Parameters
[in]remap_csDiagnostic coodinate control structure
[in]gOcean grid structure
[in]hThe current thicknesses
[in]staggered_in_xTrue is the x-axis location is at u or q points
[in]staggered_in_yTrue is the y-axis location is at v or q points
maskA mask for the field
[in]missing_valueA missing_value to assign land/vanished points
[in]fieldThe diagnostic field to be remapped
[in,out]reintegrated_fieldField argument remapped to alternative coordinate

Definition at line 479 of file MOM_diag_remap.F90.

479  type(diag_remap_ctrl), intent(in) :: remap_cs !< Diagnostic coodinate control structure
480  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
481  real, dimension(:,:,:), intent(in) :: h !< The current thicknesses
482  logical, intent(in) :: staggered_in_x !< True is the x-axis location is at u or q points
483  logical, intent(in) :: staggered_in_y !< True is the y-axis location is at v or q points
484  real, dimension(:,:,:), pointer :: mask !< A mask for the field
485  real, intent(in) :: missing_value !< A missing_value to assign land/vanished points
486  real, dimension(:,:,:), intent(in) :: field !< The diagnostic field to be remapped
487  real, dimension(:,:,:), intent(inout) :: reintegrated_field !< Field argument remapped to alternative coordinate
488  ! Local variables
489  real, dimension(remap_cs%nz) :: h_dest
490  real, dimension(size(h,3)) :: h_src
491  integer :: nz_src, nz_dest
492  integer :: i, j, k !< Grid index
493  integer :: i1, j1 !< 1-based index
494  integer :: i_lo, i_hi, j_lo, j_hi !< (uv->h) interpolation indices
495  integer :: shift !< Symmetric offset for 1-based indexing
496 
497  call assert(remap_cs%initialized, 'vertically_reintegrate_diag_field: remap_cs not initialized.')
498  call assert(size(field, 3) == size(h, 3), &
499  'vertically_reintegrate_diag_field: Remap field and thickness z-axes do not match.')
500 
501  nz_src = size(field,3)
502  nz_dest = remap_cs%nz
503  reintegrated_field(:,:,:) = 0.
504 
505  ! Symmetric grid offset under 1-based indexing; see header for details.
506  shift = 0; if (g%symmetric) shift = 1
507 
508  if (staggered_in_x .and. .not. staggered_in_y) then
509  ! U-points
510  do j=g%jsc, g%jec
511  do i=g%iscB, g%iecB
512  i1 = i - g%isdB + 1
513  i_lo = i1 - shift; i_hi = i_lo + 1
514  if (associated(mask)) then
515  if (mask(i,j,1) == 0.) cycle
516  endif
517  h_src(:) = 0.5 * (h(i_lo,j,:) + h(i_hi,j,:))
518  h_dest(:) = 0.5 * (remap_cs%h(i_lo,j,:) + remap_cs%h(i_hi,j,:))
519  call reintegrate_column(nz_src, h_src, field(i1,j,:), &
520  nz_dest, h_dest, 0., reintegrated_field(i1,j,:))
521  enddo
522  enddo
523  elseif (staggered_in_y .and. .not. staggered_in_x) then
524  ! V-points
525  do j=g%jscB, g%jecB
526  j1 = j - g%jsdB + 1
527  j_lo = j1 - shift; j_hi = j_lo + 1
528  do i=g%isc, g%iec
529  if (associated(mask)) then
530  if (mask(i,j,1) == 0.) cycle
531  endif
532  h_src(:) = 0.5 * (h(i,j_lo,:) + h(i,j_hi,:))
533  h_dest(:) = 0.5 * (remap_cs%h(i,j_lo,:) + remap_cs%h(i,j_hi,:))
534  call reintegrate_column(nz_src, h_src, field(i,j1,:), &
535  nz_dest, h_dest, 0., reintegrated_field(i,j1,:))
536  enddo
537  enddo
538  elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y)) then
539  ! H-points
540  do j=g%jsc, g%jec
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(:) = h(i,j,:)
546  h_dest(:) = remap_cs%h(i,j,:)
547  call reintegrate_column(nz_src, h_src, field(i,j,:), &
548  nz_dest, h_dest, 0., reintegrated_field(i,j,:))
549  enddo
550  enddo
551  else
552  call assert(.false., 'vertically_reintegrate_diag_field: Q point remapping is not coded yet.')
553  endif
554 

References mom_error_handler::assert(), and mom_diag_vkernels::reintegrate_column().

Here is the call graph for this function: