MOM6
mom_diag_mediator::post_data Interface Reference

Detailed Description

Make a diagnostic available for averaging or output.

Definition at line 70 of file MOM_diag_mediator.F90.

Private functions

subroutine post_data_3d (diag_field_id, field, diag_cs, is_static, mask, alt_h)
 Make a real 3-d array diagnostic available for averaging or output. More...
 
subroutine post_data_2d (diag_field_id, field, diag_cs, is_static, mask)
 Make a real 2-d array diagnostic available for averaging or output. More...
 
subroutine post_data_1d_k (diag_field_id, field, diag_cs, is_static)
 Make a real 1-d array diagnostic available for averaging or output. More...
 
subroutine post_data_0d (diag_field_id, field, diag_cs, is_static)
 Make a real scalar diagnostic available for averaging or output. More...
 

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

Make a real scalar diagnostic available for averaging or output.

Parameters
[in]diag_field_idThe id for an output variable returned by a previous call to register_diag_field.
[in]fieldreal value being offered for output or averaging
[in]diag_csStructure used to regulate diagnostic output
[in]is_staticIf true, this is a static field that is always offered.

Definition at line 1199 of file MOM_diag_mediator.F90.

1199  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1200  !! previous call to register_diag_field.
1201  real, intent(in) :: field !< real value being offered for output or averaging
1202  type(diag_ctrl), target, intent(in) :: diag_CS !< Structure used to regulate diagnostic output
1203  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1204 
1205  ! Local variables
1206  logical :: used, is_stat
1207  type(diag_type), pointer :: diag => null()
1208 
1209  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1210  is_stat = .false. ; if (present(is_static)) is_stat = is_static
1211 
1212  ! Iterate over list of diag 'variants', e.g. CMOR aliases, call send_data
1213  ! for each one.
1214  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1215  'post_data_0d: Unregistered diagnostic id')
1216  diag => diag_cs%diags(diag_field_id)
1217  do while (associated(diag))
1218  if (diag_cs%diag_as_chksum) then
1219  call chksum0(field, diag%debug_str, logunit=diag_cs%chksum_iounit)
1220  else if (is_stat) then
1221  used = send_data(diag%fms_diag_id, field)
1222  elseif (diag_cs%ave_enabled) then
1223  used = send_data(diag%fms_diag_id, field, diag_cs%time_end)
1224  endif
1225  diag => diag%next
1226  enddo
1227 
1228  if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)

◆ post_data_1d_k()

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

Make a real 1-d array diagnostic available for averaging or output.

Parameters
[in]diag_field_idThe id for an output variable returned by a previous call to register_diag_field.
[in]field1-d array being offered for output or averaging
[in]diag_csStructure used to regulate diagnostic output
[in]is_staticIf true, this is a static field that is always offered.

Definition at line 1233 of file MOM_diag_mediator.F90.

1233  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1234  !! previous call to register_diag_field.
1235  real, target, intent(in) :: field(:) !< 1-d array being offered for output or averaging
1236  type(diag_ctrl), target, intent(in) :: diag_CS !< Structure used to regulate diagnostic output
1237  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1238 
1239  ! Local variables
1240  logical :: used ! The return value of send_data is not used for anything.
1241  real, dimension(:), pointer :: locfield => null()
1242  logical :: is_stat
1243  integer :: k, ks, ke
1244  type(diag_type), pointer :: diag => null()
1245 
1246  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1247  is_stat = .false. ; if (present(is_static)) is_stat = is_static
1248 
1249  ! Iterate over list of diag 'variants', e.g. CMOR aliases.
1250  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1251  'post_data_1d_k: Unregistered diagnostic id')
1252  diag => diag_cs%diags(diag_field_id)
1253  do while (associated(diag))
1254 
1255  if ((diag%conversion_factor /= 0.) .and. (diag%conversion_factor /= 1.)) then
1256  ks = lbound(field,1) ; ke = ubound(field,1)
1257  allocate( locfield( ks:ke ) )
1258 
1259  do k=ks,ke
1260  if (field(k) == diag_cs%missing_value) then
1261  locfield(k) = diag_cs%missing_value
1262  else
1263  locfield(k) = field(k) * diag%conversion_factor
1264  endif
1265  enddo
1266  else
1267  locfield => field
1268  endif
1269 
1270  if (diag_cs%diag_as_chksum) then
1271  call zchksum(locfield, diag%debug_str, logunit=diag_cs%chksum_iounit)
1272  else if (is_stat) then
1273  used = send_data(diag%fms_diag_id, locfield)
1274  elseif (diag_cs%ave_enabled) then
1275  used = send_data(diag%fms_diag_id, locfield, diag_cs%time_end, weight=diag_cs%time_int)
1276  endif
1277  if ((diag%conversion_factor /= 0.) .and. (diag%conversion_factor /= 1.)) deallocate( locfield )
1278 
1279  diag => diag%next
1280  enddo
1281 
1282  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

Make a real 2-d array diagnostic available for averaging or output.

Parameters
[in]diag_field_idThe id for an output variable returned by a previous call to register_diag_field.
[in]field2-d array being offered for output or averaging
[in]diag_csStructure used to regulate diagnostic output
[in]is_staticIf true, this is a static field that is always offered.
[in]maskIf present, use this real array as the data mask.

Definition at line 1287 of file MOM_diag_mediator.F90.

1287  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1288  !! previous call to register_diag_field.
1289  real, intent(in) :: field(:,:) !< 2-d array being offered for output or averaging
1290  type(diag_ctrl), target, intent(in) :: diag_CS !< Structure used to regulate diagnostic output
1291  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1292  real, optional, intent(in) :: mask(:,:) !< If present, use this real array as the data mask.
1293 
1294  ! Local variables
1295  type(diag_type), pointer :: diag => null()
1296 
1297  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1298 
1299  ! Iterate over list of diag 'variants' (e.g. CMOR aliases) and post each.
1300  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1301  'post_data_2d: Unregistered diagnostic id')
1302  diag => diag_cs%diags(diag_field_id)
1303  do while (associated(diag))
1304  call post_data_2d_low(diag, field, diag_cs, is_static, mask)
1305  diag => diag%next
1306  enddo
1307 
1308  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

Make a real 3-d array diagnostic available for averaging or output.

Parameters
[in]diag_field_idThe id for an output variable returned by a previous call to register_diag_field.
[in]field3-d array being offered for output or averaging
[in]diag_csStructure used to regulate diagnostic output
[in]is_staticIf true, this is a static field that is always offered.
[in]maskIf present, use this real array as the data mask.
[in]alt_hAn alternate thickness to use for vertically

Definition at line 1458 of file MOM_diag_mediator.F90.

1458 
1459  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1460  !! previous call to register_diag_field.
1461  real, intent(in) :: field(:,:,:) !< 3-d array being offered for output or averaging
1462  type(diag_ctrl), target, intent(in) :: diag_CS !< Structure used to regulate diagnostic output
1463  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1464  real, optional, intent(in) :: mask(:,:,:) !< If present, use this real array as the data mask.
1465  real, dimension(:,:,:), &
1466  target, optional, intent(in) :: alt_h !< An alternate thickness to use for vertically
1467  !! remapping this diagnostic [H ~> m or kg m-2].
1468 
1469  ! Local variables
1470  type(diag_type), pointer :: diag => null()
1471  integer :: nz, i, j, k
1472  real, dimension(:,:,:), allocatable :: remapped_field
1473  logical :: staggered_in_x, staggered_in_y
1474  real, dimension(:,:,:), pointer :: h_diag => null()
1475 
1476  if (present(alt_h)) then
1477  h_diag => alt_h
1478  else
1479  h_diag => diag_cs%h
1480  endif
1481 
1482  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1483 
1484  ! Iterate over list of diag 'variants', e.g. CMOR aliases, different vertical
1485  ! grids, and post each.
1486  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1487  'post_data_3d: Unregistered diagnostic id')
1488  diag => diag_cs%diags(diag_field_id)
1489  do while (associated(diag))
1490  call assert(associated(diag%axes), 'post_data_3d: axes is not associated')
1491 
1492  staggered_in_x = diag%axes%is_u_point .or. diag%axes%is_q_point
1493  staggered_in_y = diag%axes%is_v_point .or. diag%axes%is_q_point
1494 
1495  if (diag%v_extensive .and. .not.diag%axes%is_native) then
1496  ! The field is vertically integrated and needs to be re-gridded
1497  if (present(mask)) then
1498  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
1499  endif
1500 
1501  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1502  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz))
1503  call vertically_reintegrate_diag_field( &
1504  diag_cs%diag_remap_cs(diag%axes%vertical_coordinate_number), &
1505  diag_cs%G, h_diag, staggered_in_x, staggered_in_y, &
1506  diag%axes%mask3d, diag_cs%missing_value, field, remapped_field)
1507  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1508  if (associated(diag%axes%mask3d)) then
1509  ! Since 3d masks do not vary in the vertical, just use as much as is
1510  ! needed.
1511  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
1512  mask=diag%axes%mask3d)
1513  else
1514  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
1515  endif
1516  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1517  deallocate(remapped_field)
1518  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1519  elseif (diag%axes%needs_remapping) then
1520  ! Remap this field to another vertical coordinate.
1521  if (present(mask)) then
1522  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
1523  endif
1524 
1525  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1526  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz))
1527  call diag_remap_do_remap(diag_cs%diag_remap_cs( &
1528  diag%axes%vertical_coordinate_number), &
1529  diag_cs%G, diag_cs%GV, h_diag, staggered_in_x, staggered_in_y, &
1530  diag%axes%mask3d, diag_cs%missing_value, field, remapped_field)
1531  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1532  if (associated(diag%axes%mask3d)) then
1533  ! Since 3d masks do not vary in the vertical, just use as much as is
1534  ! needed.
1535  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
1536  mask=diag%axes%mask3d)
1537  else
1538  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
1539  endif
1540  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1541  deallocate(remapped_field)
1542  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1543  elseif (diag%axes%needs_interpolating) then
1544  ! Interpolate this field to another vertical coordinate.
1545  if (present(mask)) then
1546  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
1547  endif
1548 
1549  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1550  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz+1))
1551  call vertically_interpolate_diag_field(diag_cs%diag_remap_cs( &
1552  diag%axes%vertical_coordinate_number), &
1553  diag_cs%G, h_diag, staggered_in_x, staggered_in_y, &
1554  diag%axes%mask3d, diag_cs%missing_value, field, remapped_field)
1555  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1556  if (associated(diag%axes%mask3d)) then
1557  ! Since 3d masks do not vary in the vertical, just use as much as is
1558  ! needed.
1559  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
1560  mask=diag%axes%mask3d)
1561  else
1562  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
1563  endif
1564  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1565  deallocate(remapped_field)
1566  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1567  else
1568  call post_data_3d_low(diag, field, diag_cs, is_static, mask)
1569  endif
1570  diag => diag%next
1571  enddo
1572  if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)
1573 

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