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...
 

Detailed Description

Make a diagnostic available for averaging or output.

Definition at line 70 of file MOM_diag_mediator.F90.

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 1257 of file MOM_diag_mediator.F90.

1257  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1258  !! previous call to register_diag_field.
1259  real, intent(in) :: field !< real value being offered for output or averaging
1260  type(diag_ctrl), target, intent(in) :: diag_cs !< Structure used to regulate diagnostic output
1261  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1262 
1263  ! Local variables
1264  real :: locfield
1265  logical :: used, is_stat
1266  type(diag_type), pointer :: diag => null()
1267 
1268  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1269  is_stat = .false. ; if (present(is_static)) is_stat = is_static
1270 
1271  ! Iterate over list of diag 'variants', e.g. CMOR aliases, call send_data
1272  ! for each one.
1273  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1274  'post_data_0d: Unregistered diagnostic id')
1275  diag => diag_cs%diags(diag_field_id)
1276 
1277  do while (associated(diag))
1278  locfield = field
1279  if (diag%conversion_factor /= 0.) &
1280  locfield = locfield * diag%conversion_factor
1281 
1282  if (diag_cs%diag_as_chksum) then
1283  call chksum0(locfield, diag%debug_str, logunit=diag_cs%chksum_iounit)
1284  else if (is_stat) then
1285  used = send_data(diag%fms_diag_id, locfield)
1286  elseif (diag_cs%ave_enabled) then
1287  used = send_data(diag%fms_diag_id, locfield, diag_cs%time_end)
1288  endif
1289  diag => diag%next
1290  enddo
1291 
1292  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 1297 of file MOM_diag_mediator.F90.

1297  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1298  !! previous call to register_diag_field.
1299  real, target, intent(in) :: field(:) !< 1-d array being offered for output or averaging
1300  type(diag_ctrl), target, intent(in) :: diag_cs !< Structure used to regulate diagnostic output
1301  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1302 
1303  ! Local variables
1304  logical :: used ! The return value of send_data is not used for anything.
1305  real, dimension(:), pointer :: locfield => null()
1306  logical :: is_stat
1307  integer :: k, ks, ke
1308  type(diag_type), pointer :: diag => null()
1309 
1310  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1311  is_stat = .false. ; if (present(is_static)) is_stat = is_static
1312 
1313  ! Iterate over list of diag 'variants', e.g. CMOR aliases.
1314  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1315  'post_data_1d_k: Unregistered diagnostic id')
1316  diag => diag_cs%diags(diag_field_id)
1317  do while (associated(diag))
1318 
1319  if ((diag%conversion_factor /= 0.) .and. (diag%conversion_factor /= 1.)) then
1320  ks = lbound(field,1) ; ke = ubound(field,1)
1321  allocate( locfield( ks:ke ) )
1322 
1323  do k=ks,ke
1324  if (field(k) == diag_cs%missing_value) then
1325  locfield(k) = diag_cs%missing_value
1326  else
1327  locfield(k) = field(k) * diag%conversion_factor
1328  endif
1329  enddo
1330  else
1331  locfield => field
1332  endif
1333 
1334  if (diag_cs%diag_as_chksum) then
1335  call zchksum(locfield, diag%debug_str, logunit=diag_cs%chksum_iounit)
1336  else if (is_stat) then
1337  used = send_data(diag%fms_diag_id, locfield)
1338  elseif (diag_cs%ave_enabled) then
1339  used = send_data(diag%fms_diag_id, locfield, diag_cs%time_end, weight=diag_cs%time_int)
1340  endif
1341  if ((diag%conversion_factor /= 0.) .and. (diag%conversion_factor /= 1.)) deallocate( locfield )
1342 
1343  diag => diag%next
1344  enddo
1345 
1346  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 1351 of file MOM_diag_mediator.F90.

1351  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1352  !! previous call to register_diag_field.
1353  real, intent(in) :: field(:,:) !< 2-d array being offered for output or averaging
1354  type(diag_ctrl), target, intent(in) :: diag_cs !< Structure used to regulate diagnostic output
1355  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1356  real, optional, intent(in) :: mask(:,:) !< If present, use this real array as the data mask.
1357 
1358  ! Local variables
1359  type(diag_type), pointer :: diag => null()
1360 
1361  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1362 
1363  ! Iterate over list of diag 'variants' (e.g. CMOR aliases) and post each.
1364  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1365  'post_data_2d: Unregistered diagnostic id')
1366  diag => diag_cs%diags(diag_field_id)
1367  do while (associated(diag))
1368  call post_data_2d_low(diag, field, diag_cs, is_static, mask)
1369  diag => diag%next
1370  enddo
1371 
1372  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 1522 of file MOM_diag_mediator.F90.

1522 
1523  integer, intent(in) :: diag_field_id !< The id for an output variable returned by a
1524  !! previous call to register_diag_field.
1525  real, intent(in) :: field(:,:,:) !< 3-d array being offered for output or averaging
1526  type(diag_ctrl), target, intent(in) :: diag_cs !< Structure used to regulate diagnostic output
1527  logical, optional, intent(in) :: is_static !< If true, this is a static field that is always offered.
1528  real, optional, intent(in) :: mask(:,:,:) !< If present, use this real array as the data mask.
1529  real, dimension(:,:,:), &
1530  target, optional, intent(in) :: alt_h !< An alternate thickness to use for vertically
1531  !! remapping this diagnostic [H ~> m or kg m-2].
1532 
1533  ! Local variables
1534  type(diag_type), pointer :: diag => null()
1535  integer :: nz, i, j, k
1536  real, dimension(:,:,:), allocatable :: remapped_field
1537  logical :: staggered_in_x, staggered_in_y
1538  real, dimension(:,:,:), pointer :: h_diag => null()
1539 
1540  if (id_clock_diag_mediator>0) call cpu_clock_begin(id_clock_diag_mediator)
1541 
1542  ! For intensive variables only, we can choose to use a different diagnostic grid
1543  ! to map to
1544  if (present(alt_h)) then
1545  h_diag => alt_h
1546  else
1547  h_diag => diag_cs%h
1548  endif
1549 
1550  ! Iterate over list of diag 'variants', e.g. CMOR aliases, different vertical
1551  ! grids, and post each.
1552  call assert(diag_field_id < diag_cs%next_free_diag_id, &
1553  'post_data_3d: Unregistered diagnostic id')
1554  diag => diag_cs%diags(diag_field_id)
1555  do while (associated(diag))
1556  call assert(associated(diag%axes), 'post_data_3d: axes is not associated')
1557 
1558  staggered_in_x = diag%axes%is_u_point .or. diag%axes%is_q_point
1559  staggered_in_y = diag%axes%is_v_point .or. diag%axes%is_q_point
1560 
1561  if (diag%v_extensive .and. .not.diag%axes%is_native) then
1562  ! The field is vertically integrated and needs to be re-gridded
1563  if (present(mask)) then
1564  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
1565  endif
1566 
1567  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1568  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz))
1569  call vertically_reintegrate_diag_field( &
1570  diag_cs%diag_remap_cs(diag%axes%vertical_coordinate_number), diag_cs%G, &
1571  diag_cs%h_begin, &
1572  diag_cs%diag_remap_cs(diag%axes%vertical_coordinate_number)%h_extensive, &
1573  staggered_in_x, staggered_in_y, diag%axes%mask3d, field, remapped_field)
1574  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1575  if (associated(diag%axes%mask3d)) then
1576  ! Since 3d masks do not vary in the vertical, just use as much as is
1577  ! needed.
1578  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
1579  mask=diag%axes%mask3d)
1580  else
1581  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
1582  endif
1583  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1584  deallocate(remapped_field)
1585  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1586  elseif (diag%axes%needs_remapping) then
1587  ! Remap this field to another vertical coordinate.
1588  if (present(mask)) then
1589  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
1590  endif
1591 
1592  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1593  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz))
1594  call diag_remap_do_remap(diag_cs%diag_remap_cs(diag%axes%vertical_coordinate_number), &
1595  diag_cs%G, diag_cs%GV, h_diag, staggered_in_x, staggered_in_y, &
1596  diag%axes%mask3d, field, remapped_field)
1597  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1598  if (associated(diag%axes%mask3d)) then
1599  ! Since 3d masks do not vary in the vertical, just use as much as is
1600  ! needed.
1601  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
1602  mask=diag%axes%mask3d)
1603  else
1604  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
1605  endif
1606  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1607  deallocate(remapped_field)
1608  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1609  elseif (diag%axes%needs_interpolating) then
1610  ! Interpolate this field to another vertical coordinate.
1611  if (present(mask)) then
1612  call mom_error(fatal,"post_data_3d: no mask for regridded field.")
1613  endif
1614 
1615  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1616  allocate(remapped_field(size(field,1), size(field,2), diag%axes%nz+1))
1617  call vertically_interpolate_diag_field(diag_cs%diag_remap_cs( &
1618  diag%axes%vertical_coordinate_number), &
1619  diag_cs%G, h_diag, staggered_in_x, staggered_in_y, &
1620  diag%axes%mask3d, field, remapped_field)
1621  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1622  if (associated(diag%axes%mask3d)) then
1623  ! Since 3d masks do not vary in the vertical, just use as much as is
1624  ! needed.
1625  call post_data_3d_low(diag, remapped_field, diag_cs, is_static, &
1626  mask=diag%axes%mask3d)
1627  else
1628  call post_data_3d_low(diag, remapped_field, diag_cs, is_static)
1629  endif
1630  if (id_clock_diag_remap>0) call cpu_clock_begin(id_clock_diag_remap)
1631  deallocate(remapped_field)
1632  if (id_clock_diag_remap>0) call cpu_clock_end(id_clock_diag_remap)
1633  else
1634  call post_data_3d_low(diag, field, diag_cs, is_static, mask)
1635  endif
1636  diag => diag%next
1637  enddo
1638  if (id_clock_diag_mediator>0) call cpu_clock_end(id_clock_diag_mediator)
1639 

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