MOM6
mom_diag_mediator::post_data Interface Reference

Detailed Description

Definition at line 80 of file MOM_diag_mediator.F90.

Private functions

subroutine post_data_3d (diag_field_id, field, diag_cs, is_static, mask, alt_h)
 
subroutine post_data_2d (diag_field_id, field, diag_cs, is_static, mask)
 
subroutine post_data_0d (diag_field_id, field, diag_cs, is_static)
 

Functions and subroutines

◆ post_data_0d()

subroutine mom_diag_mediator::post_data::post_data_0d ( integer, intent(in)  diag_field_id,
real, intent(in)  field,
type(diag_ctrl), intent(in), target  diag_cs,
logical, intent(in), optional  is_static 
)
private

Definition at line 516 of file MOM_diag_mediator.F90.

516  integer, intent(in) :: diag_field_id
517  real, intent(in) :: field
518  type(diag_ctrl), target, intent(in) :: diag_cs
519  logical, optional, intent(in) :: is_static
520 
521 ! Arguments:
522 ! (in) diag_field_id - the id for an output variable returned by a
523 ! previous call to register_diag_field.
524 ! (in) field - 0-d array being offered for output or averaging.
525 ! (inout) diag_cs - structure used to regulate diagnostic output.
526 ! (in,opt) is_static - If true, this is a static field that is always offered.
527 ! (in,opt) mask - If present, use this real array as the data mask.
528 
529  logical :: used, is_stat
530  type(diag_type), pointer :: diag => null()
531 
532  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
533  is_stat = .false. ; if (present(is_static)) is_stat = is_static
534 
535  ! Iterate over list of diag 'variants', e.g. CMOR aliases, call send_data
536  ! for each one.
537  call assert(diag_field_id < diag_cs%next_free_diag_id, &
538  'post_data_0d: Unregistered diagnostic id')
539  diag => diag_cs%diags(diag_field_id)
540  do while (associated(diag))
541  if (is_stat) then
542  used = send_data(diag%fms_diag_id, field)
543  elseif (diag_cs%ave_enabled) then
544  used = send_data(diag%fms_diag_id, field, diag_cs%time_end)
545  endif
546  diag => diag%next
547  enddo
548 
549  if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)

◆ post_data_2d()

subroutine mom_diag_mediator::post_data::post_data_2d ( integer, intent(in)  diag_field_id,
real, dimension(:,:), intent(in)  field,
type(diag_ctrl), intent(in), target  diag_cs,
logical, intent(in), optional  is_static,
real, dimension(:,:), intent(in), optional  mask 
)
private

Definition at line 590 of file MOM_diag_mediator.F90.

590  integer, intent(in) :: diag_field_id
591  real, intent(in) :: field(:,:)
592  type(diag_ctrl), target, intent(in) :: diag_cs
593  logical, optional, intent(in) :: is_static
594  real, optional, intent(in) :: mask(:,:)
595 
596 ! Arguments:
597 ! (in) diag_field_id - id for an output variable returned by a
598 ! previous call to register_diag_field.
599 ! (in) field - 2-d array being offered for output or averaging.
600 ! (inout) diag_cs - structure used to regulate diagnostic output.
601 ! (in,opt) is_static - If true, this is a static field that is always offered.
602 ! (in,opt) mask - If present, use this real array as the data mask.
603 
604  type(diag_type), pointer :: diag => null()
605 
606  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
607 
608  ! Iterate over list of diag 'variants' (e.g. CMOR aliases) and post each.
609  call assert(diag_field_id < diag_cs%next_free_diag_id, &
610  'post_data_2d: Unregistered diagnostic id')
611  diag => diag_cs%diags(diag_field_id)
612  do while (associated(diag))
613  call post_data_2d_low(diag, field, diag_cs, is_static, mask)
614  diag => diag%next
615  enddo
616 
617  if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)

◆ post_data_3d()

subroutine mom_diag_mediator::post_data::post_data_3d ( integer, intent(in)  diag_field_id,
real, dimension(:,:,:), intent(in)  field,
type(diag_ctrl), intent(in), target  diag_cs,
logical, intent(in), optional  is_static,
real, dimension(:,:,:), intent(in), optional  mask,
real, dimension(:,:,:), intent(in), optional, target  alt_h 
)
private

Definition at line 712 of file MOM_diag_mediator.F90.

712 
713  integer, intent(in) :: diag_field_id
714  real, intent(in) :: field(:,:,:)
715  type(diag_ctrl), target, intent(in) :: diag_cs
716  logical, optional, intent(in) :: is_static
717  real, optional, intent(in) :: mask(:,:,:)
718  real, target, optional, intent(in) :: alt_h(:,:,:)
719 
720 ! Arguments:
721 ! (in) diag_field_id - id for an output variable returned by a
722 ! previous call to register_diag_field.
723 ! (in) field - 3-d array being offered for output or averaging
724 ! (inout) diag - structure used to regulate diagnostic output
725 ! (in) static - If true, this is a static field that is always offered.
726 ! (in,opt) mask - If present, use this real array as the data mask.
727 
728  type(diag_type), pointer :: diag => null()
729  integer :: nz, i, j, k
730  real, dimension(:,:,:), allocatable :: remapped_field
731  logical :: staggered_in_x, staggered_in_y
732  real, dimension(:,:,:), pointer :: h_diag
733 
734  if(present(alt_h)) then
735  h_diag => alt_h
736  else
737  h_diag => diag_cs%h
738  endif
739 
740  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
741 
742  ! Iterate over list of diag 'variants', e.g. CMOR aliases, different vertical
743  ! grids, and post each.
744  call assert(diag_field_id < diag_cs%next_free_diag_id, &
745  'post_data_3d: Unregistered diagnostic id')
746  diag => diag_cs%diags(diag_field_id)
747  do while (associated(diag))
748  call assert(associated(diag%axes), 'post_data_3d: axes is not associated')
749 
750  staggered_in_x = diag%axes%is_u_point .or. diag%axes%is_q_point
751  staggered_in_y = diag%axes%is_v_point .or. diag%axes%is_q_point
752 
753  if (diag%v_extensive .and. .not.diag%axes%is_native) then
754  ! The field is vertically integrated and needs to be re-gridded
755  if (present(mask)) then
756  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
757  endif
758 
759  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
760  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz))
761  call vertically_reintegrate_diag_field( &
762  diag_cs%diag_remap_cs(diag%axes%vertical_coordinate_number), &
763  diag_cs%G, h_diag, staggered_in_x, staggered_in_y, &
764  diag%mask3d, diag_cs%missing_value, field, remapped_field)
765  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
766  if (associated(diag%mask3d)) then
767  ! Since 3d masks do not vary in the vertical, just use as much as is
768  ! needed.
769  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
770  mask=diag%mask3d(:,:,:diag%axes%nz))
771  else
772  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
773  endif
774  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
775  deallocate(remapped_field)
776  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
777  elseif (diag%axes%needs_remapping) then
778  ! Remap this field to another vertical coordinate.
779  if (present(mask)) then
780  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
781  endif
782 
783  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
784  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz))
785  call diag_remap_do_remap(diag_cs%diag_remap_cs( &
786  diag%axes%vertical_coordinate_number), &
787  diag_cs%G, h_diag, staggered_in_x, staggered_in_y, &
788  diag%mask3d, diag_cs%missing_value, field, remapped_field)
789  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
790  if (associated(diag%mask3d)) then
791  ! Since 3d masks do not vary in the vertical, just use as much as is
792  ! needed.
793  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
794  mask=diag%mask3d(:,:,:diag%axes%nz))
795  else
796  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
797  endif
798  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
799  deallocate(remapped_field)
800  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
801  elseif (diag%axes%needs_interpolating) then
802  ! Interpolate this field to another vertical coordinate.
803  if (present(mask)) then
804  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
805  endif
806 
807  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
808  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz+1))
809  call vertically_interpolate_diag_field(diag_cs%diag_remap_cs( &
810  diag%axes%vertical_coordinate_number), &
811  diag_cs%G, h_diag, staggered_in_x, staggered_in_y, &
812  diag%mask3d, diag_cs%missing_value, field, remapped_field)
813  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
814  if (associated(diag%mask3d)) then
815  ! Since 3d masks do not vary in the vertical, just use as much as is
816  ! needed.
817  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
818  mask=diag%mask3d(:,:,:diag%axes%nz+1))
819  else
820  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
821  endif
822  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
823  deallocate(remapped_field)
824  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
825  else
826  call post_data_3d_low(diag, field, diag_cs, is_static, mask)
827  endif
828  diag => diag%next
829  enddo
830  if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)
831 

The documentation for this interface was generated from the following file: