83 use diag_axis_mod,
only : get_diag_axis_name
84 use diag_manager_mod,
only : diag_axis_init
87 implicit none ;
private 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
105 logical :: configured = .false.
106 logical :: initialized = .false.
107 logical :: used = .false.
108 integer :: vertical_coord = 0
109 character(len=10) :: vertical_coord_name =
'' 110 character(len=16) :: diag_coord_name =
'' 111 character(len=8) :: diag_module_suffix =
'' 115 real,
dimension(:,:,:),
allocatable :: h
116 real,
dimension(:,:,:),
allocatable :: h_extensive
118 integer :: interface_axes_id = 0
119 integer :: layer_axes_id = 0
120 logical :: answers_2018
128 subroutine diag_remap_init(remap_cs, coord_tuple, answers_2018)
130 character(len=*),
intent(in) :: coord_tuple
132 logical,
intent(in) :: answers_2018
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
146 end subroutine diag_remap_init
150 subroutine diag_remap_end(remap_cs)
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.
159 end subroutine diag_remap_end
166 subroutine diag_remap_diag_registration_closed(remap_cs)
169 if (.not. remap_cs%used)
then 170 call diag_remap_end(remap_cs)
173 end subroutine diag_remap_diag_registration_closed
178 subroutine diag_remap_set_active(remap_cs)
181 remap_cs%used = .true.
183 end subroutine diag_remap_set_active
187 subroutine diag_remap_configure_axes(remap_cs, GV, US, param_file)
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" 196 character(len=8) :: units, expected_units
197 character(len=34) :: longname, string2
199 character(len=256) :: err_msg
202 real,
allocatable,
dimension(:) :: interfaces, layers
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.)
208 remap_cs%nz = get_regrid_size(remap_cs%regrid_cs)
210 if (remap_cs%vertical_coord == coordinatemode(
'SIGMA'))
then 212 longname =
'Fraction' 213 elseif (remap_cs%vertical_coord == coordinatemode(
'RHO'))
then 215 longname =
'Target Potential Density' 222 allocate(interfaces(remap_cs%nz+1))
223 allocate(layers(remap_cs%nz))
225 interfaces(:) = getcoordinateinterfaces(remap_cs%regrid_cs)
226 layers(:) = 0.5 * ( interfaces(1:remap_cs%nz) + interfaces(2:remap_cs%nz+1) )
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)
237 remap_cs%configured = .true.
239 deallocate(interfaces)
242 end subroutine diag_remap_configure_axes
246 subroutine diag_remap_get_axes_info(remap_cs, nz, id_layer, id_interface)
248 integer,
intent(out) :: nz
249 integer,
intent(out) :: id_layer
250 integer,
intent(out) :: id_interface
253 id_layer = remap_cs%layer_axes_id
254 id_interface = remap_cs%interface_axes_id
256 end subroutine diag_remap_get_axes_info
262 function diag_remap_axes_configured(remap_cs)
264 logical :: diag_remap_axes_configured
266 diag_remap_axes_configured = remap_cs%configured
275 subroutine diag_remap_update(remap_cs, G, GV, US, h, T, S, eqn_of_state, h_target)
280 real,
dimension(:,:,:),
intent(in) :: h
281 real,
dimension(:,:,:),
intent(in) :: t
282 real,
dimension(:,:,:),
intent(in) :: s
283 type(
eos_type),
pointer :: eqn_of_state
284 real,
dimension(:,:,:),
intent(inout) :: h_target
287 real,
dimension(remap_cs%nz + 1) :: zinterfaces
288 real :: h_neglect, h_neglect_edge
289 integer :: i, j, k, nz
293 if (.not. remap_cs%configured)
then 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
302 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
306 if (.not. remap_cs%initialized)
then 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.
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 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 336 call mom_error(fatal,
"diag_remap_update: SLIGHT coordinate not coded for diagnostics yet!")
337 elseif (remap_cs%vertical_coord == coordinatemode(
'HYCOM1'))
then 340 call mom_error(fatal,
"diag_remap_update: HYCOM1 coordinate not coded for diagnostics yet!")
343 h_target(i,j,k) = zinterfaces(k) - zinterfaces(k+1)
347 end subroutine diag_remap_update
350 subroutine diag_remap_do_remap(remap_cs, G, GV, h, staggered_in_x, staggered_in_y, &
351 mask, field, remapped_field)
355 real,
dimension(:,:,:),
intent(in) :: h
356 logical,
intent(in) :: staggered_in_x
357 logical,
intent(in) :: staggered_in_y
358 real,
dimension(:,:,:),
pointer :: mask
359 real,
dimension(:,:,:),
intent(in) :: field(:,:,:)
360 real,
dimension(:,:,:),
intent(inout) :: remapped_field
362 real,
dimension(remap_cs%nz) :: h_dest
363 real,
dimension(size(h,3)) :: h_src
364 real :: h_neglect, h_neglect_edge
365 integer :: nz_src, nz_dest
368 integer :: i_lo, i_hi, j_lo, j_hi
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.')
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
380 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
383 nz_src =
size(field,3)
384 nz_dest = remap_cs%nz
385 remapped_field(:,:,:) = 0.
388 shift = 0;
if (g%symmetric) shift = 1
390 if (staggered_in_x .and. .not. staggered_in_y)
then 395 i_lo = i1 - shift; i_hi = i_lo + 1
396 if (
associated(mask))
then 397 if (mask(i,j,1) == 0.) cycle
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)
407 elseif (staggered_in_y .and. .not. staggered_in_x)
then 411 j_lo = j1 - shift; j_hi = j_lo + 1
413 if (
associated(mask))
then 414 if (mask(i,j,1) == 0.) cycle
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)
424 elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y))
then 428 if (
associated(mask))
then 429 if (mask(i,j,1) == 0.) cycle
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)
440 call assert(.false.,
'diag_remap_do_remap: Unsupported axis combination')
443 end subroutine diag_remap_do_remap
446 subroutine diag_remap_calc_hmask(remap_cs, G, mask)
449 real,
dimension(:,:,:),
intent(out) :: mask
451 real,
dimension(remap_cs%nz) :: h_dest
453 logical :: mask_vanished_layers
457 call assert(remap_cs%initialized,
'diag_remap_calc_hmask: remap_cs not initialized.')
460 mask_vanished_layers = (remap_cs%vertical_coord == coordinatemode(
'ZSTAR'))
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,:)
470 h_tot = h_tot + h_dest(k)
473 h_err = h_err + epsilon(h_tot) * h_tot
475 if (h_dest(k)<=8.*h_err)
then 487 end subroutine diag_remap_calc_hmask
490 subroutine vertically_reintegrate_diag_field(remap_cs, G, h, h_target, staggered_in_x, staggered_in_y, &
491 mask, field, reintegrated_field)
494 real,
dimension(:,:,:),
intent(in) :: h
495 real,
dimension(:,:,:),
intent(in) :: h_target
496 logical,
intent(in) :: staggered_in_x
497 logical,
intent(in) :: staggered_in_y
498 real,
dimension(:,:,:),
pointer :: mask
499 real,
dimension(:,:,:),
intent(in) :: field
500 real,
dimension(:,:,:),
intent(inout) :: reintegrated_field
502 real,
dimension(remap_cs%nz) :: h_dest
503 real,
dimension(size(h,3)) :: h_src
504 integer :: nz_src, nz_dest
507 integer :: i_lo, i_hi, j_lo, j_hi
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.')
514 nz_src =
size(field,3)
515 nz_dest = remap_cs%nz
516 reintegrated_field(:,:,:) = 0.
519 shift = 0;
if (g%symmetric) shift = 1
521 if (staggered_in_x .and. .not. staggered_in_y)
then 526 i_lo = i1 - shift; i_hi = i_lo + 1
527 if (
associated(mask))
then 528 if (mask(i,j,1) == 0.) cycle
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,:))
536 elseif (staggered_in_y .and. .not. staggered_in_x)
then 540 j_lo = j1 - shift; j_hi = j_lo + 1
542 if (
associated(mask))
then 543 if (mask(i,j,1) == 0.) cycle
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,:))
551 elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y))
then 555 if (
associated(mask))
then 556 if (mask(i,j,1) == 0.) cycle
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,:))
565 call assert(.false.,
'vertically_reintegrate_diag_field: Q point remapping is not coded yet.')
568 end subroutine vertically_reintegrate_diag_field
571 subroutine vertically_interpolate_diag_field(remap_cs, G, h, staggered_in_x, staggered_in_y, &
572 mask, field, interpolated_field)
575 real,
dimension(:,:,:),
intent(in) :: h
576 logical,
intent(in) :: staggered_in_x
577 logical,
intent(in) :: staggered_in_y
578 real,
dimension(:,:,:),
pointer :: mask
579 real,
dimension(:,:,:),
intent(in) :: field
580 real,
dimension(:,:,:),
intent(inout) :: interpolated_field
582 real,
dimension(remap_cs%nz) :: h_dest
583 real,
dimension(size(h,3)) :: h_src
584 integer :: nz_src, nz_dest
587 integer :: i_lo, i_hi, j_lo, j_hi
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.')
594 interpolated_field(:,:,:) = 0.
597 nz_dest = remap_cs%nz
600 shift = 0;
if (g%symmetric) shift = 1
602 if (staggered_in_x .and. .not. staggered_in_y)
then 607 i_lo = i1 - shift; i_hi = i_lo + 1
608 if (
associated(mask))
then 609 if (mask(i,j,1) == 0.) cycle
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,:))
617 elseif (staggered_in_y .and. .not. staggered_in_x)
then 621 j_lo = j1 - shift; j_hi = j_lo + 1
623 if (
associated(mask))
then 624 if (mask(i,j,1) == 0.) cycle
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,:))
632 elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y))
then 636 if (
associated(mask))
then 637 if (mask(i,j,1) == 0.) cycle
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,:))
646 call assert(.false.,
'vertically_interpolate_diag_field: Q point remapping is not coded yet.')
649 end subroutine vertically_interpolate_diag_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, &
658 real,
dimension(:,:,:),
intent(in) :: h
659 logical,
intent(in) :: staggered_in_x
660 logical,
intent(in) :: staggered_in_y
661 logical,
intent(in) :: is_layer
662 logical,
intent(in) :: is_extensive
663 real,
dimension(:,:,:),
intent(in) :: field
664 real,
dimension(:),
intent(inout) :: averaged_field
665 logical,
dimension(:),
intent(inout) :: averaged_mask
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
670 type(
efp_type),
dimension(2*size(field,3)) :: sums_efp
672 integer :: i, j, k, nz
681 if (staggered_in_x .and. .not. staggered_in_y)
then 687 if (is_extensive)
then 688 do j=g%jsc, g%jec ;
do i=g%isc, g%iec
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)
694 do j=g%jsc, g%jec ;
do i=g%isc, g%iec
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)
705 do j=g%jsc, g%jec ;
do i=g%isc, g%iec
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)
712 elseif (staggered_in_y .and. .not. staggered_in_x)
then 716 if (is_extensive)
then 717 do j=g%jsc, g%jec ;
do i=g%isc, g%iec
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)
723 do j=g%jsc, g%jec ;
do i=g%isc, g%iec
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)
734 do j=g%jsc, g%jec ;
do i=g%isc, g%iec
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)
741 elseif ((.not. staggered_in_x) .and. (.not. staggered_in_y))
then 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)
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)
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)
772 call assert(.false.,
'horizontally_average_diag_field: Q point averaging is not coded yet.')
783 vol_sum(k) = efp_to_real(sums_efp(2*k-1))
784 stuff_sum(k) = efp_to_real(sums_efp(2*k))
787 averaged_mask(:) = .true.
789 if (vol_sum(k) > 0.)
then 790 averaged_field(k) = stuff_sum(k) / vol_sum(k)
792 averaged_field(k) = 0.
793 averaged_mask(k) = .false.
797 end subroutine horizontally_average_diag_field
Provides kernels for single-column interpolation, re-integration (re-mapping of integrated quantities...
Ocean grid type. See mom_grid for details.
Sum a value or 1-d array of values across processors, returning the sums in place.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Generates vertical grids as part of the ALE algorithm.
Provides column-wise vertical remapping functions.
This module contains I/O framework code.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Container for remapping parameters.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Describes various unit conversion factors.
Regridding control structure.
Represents remapping of diagnostics to a particular vertical coordinate.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
Regrid columns for the sigma coordinate.
Regrid columns for a z-like coordinate (z-star, z-level)
provides runtime remapping of diagnostics to z star, sigma and rho vertical coordinates.
Provides subroutines for quantities specific to the equation of state.
Describes the vertical ocean grid, including unit conversion factors.
Indicate whether a file exists, perhaps with domain decomposition.
Regrid columns for the continuous isopycnal (rho) coordinate.
Contains constants for interpreting input parameters that control regridding.
Find an accurate and order-invariant sum of a distributed 2d field, returning the result in the form ...
Handy functions for manipulating strings.
Provides a transparent vertical ocean grid type and supporting routines.
Read a data field from a file.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.
A control structure for the equation of state.
The Extended Fixed Point (EFP) type provides a public interface for doing sums and taking differences...