MOM6
mom_checksums::uchksum Interface Reference

Detailed Description

Checksums an array (2d or 3d) staggered at C-grid u points.

Definition at line 34 of file MOM_checksums.F90.

Private functions

subroutine chksum_u_2d (array_m, mesg, HI_m, haloshift, symmetric, omit_corners, scale, logunit)
 Checksums a 2d array staggered at C-grid u points. More...
 
subroutine chksum_u_3d (array_m, mesg, HI_m, haloshift, symmetric, omit_corners, scale, logunit)
 Checksums a 3d array staggered at C-grid u points. More...
 

Detailed Description

Checksums an array (2d or 3d) staggered at C-grid u points.

Definition at line 34 of file MOM_checksums.F90.

Functions and subroutines

◆ chksum_u_2d()

subroutine mom_checksums::uchksum::chksum_u_2d ( real, dimension(hi_m%isdb:,hi_m%jsd:), intent(in), target  array_m,
character(len=*), intent(in)  mesg,
type(hor_index_type), intent(in), target  HI_m,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric,
logical, intent(in), optional  omit_corners,
real, intent(in), optional  scale,
integer, intent(in), optional  logunit 
)
private

Checksums a 2d array staggered at C-grid u points.

Parameters
[in]hi_mA horizontal index type
[in]array_mThe array to be checksummed
[in]mesgAn identifying message
[in]haloshiftThe width of halos to check (default 0)
[in]symmetricIf true, do the checksums on the full symmetric computational domain.
[in]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.
[in]logunitIO unit for checksum logging

Definition at line 845 of file MOM_checksums.F90.

847  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
848  real, dimension(HI_m%IsdB:,HI_m%jsd:), target, intent(in) :: array_m !< The array to be checksummed
849  character(len=*), intent(in) :: mesg !< An identifying message
850  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
851  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
852  !! symmetric computational domain.
853  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
854  real, optional, intent(in) :: scale !< A scaling factor for this array.
855  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
856 
857  real, pointer :: array(:,:) ! Field array on the input grid
858  real, allocatable, dimension(:,:) :: rescaled_array
859  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
860  real :: scaling
861  integer :: iounit !< Log IO unit
862  integer :: i, j, Is
863  real :: aMean, aMin, aMax
864  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
865  integer :: bcN, bcS, bcE, bcW
866  logical :: do_corners, sym, sym_stats
867  integer :: turns ! Quarter turns from input to model grid
868 
869  ! Rotate array to the input grid
870  turns = hi_m%turns
871  if (modulo(turns, 4) /= 0) then
872  allocate(hi)
873  call rotate_hor_index(hi_m, -turns, hi)
874  if (modulo(turns, 2) /= 0) then
875  ! Arrays originating from v-points must be handled by vchksum
876  allocate(array(hi%isd:hi%ied, hi%JsdB:hi%JedB))
877  call rotate_array(array_m, -turns, array)
878  call vchksum(array, mesg, hi, haloshift, symmetric, omit_corners, scale, logunit)
879  return
880  else
881  allocate(array(hi%IsdB:hi%IedB, hi%jsd:hi%jed))
882  call rotate_array(array_m, -turns, array)
883  endif
884  else
885  hi => hi_m
886  array => array_m
887  endif
888 
889  if (checkfornans) then
890  if (is_nan(array(hi%IscB:hi%IecB,hi%jsc:hi%jec))) &
891  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
892 ! if (is_NaN(array)) &
893 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
894  endif
895 
896  scaling = 1.0 ; if (present(scale)) scaling = scale
897  iounit = error_unit; if(present(logunit)) iounit = logunit
898  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
899  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
900 
901  if (calculatestatistics) then
902  if (present(scale)) then
903  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
904  lbound(array,2):ubound(array,2)) )
905  rescaled_array(:,:) = 0.0
906  is = hi%isc ; if (sym_stats) is = hi%isc-1
907  do j=hi%jsc,hi%jec ; do i=is,hi%IecB
908  rescaled_array(i,j) = scale*array(i,j)
909  enddo ; enddo
910  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
911  deallocate(rescaled_array)
912  else
913  call substats(hi, array, sym_stats, amean, amin, amax)
914  endif
915 
916  if (is_root_pe()) &
917  call chk_sum_msg("u-point:", amean, amin, amax, mesg, iounit)
918  endif
919 
920  if (.not.writechksums) return
921 
922  hshift = default_shift
923  if (present(haloshift)) hshift = haloshift
924  if (hshift<0) hshift = hi%iedB-hi%iecB
925 
926  if ( hi%iscB-hshift<hi%isdB .or. hi%iecB+hshift>hi%iedB .or. &
927  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
928  write(0,*) 'chksum_u_2d: haloshift =',hshift
929  write(0,*) 'chksum_u_2d: isd,isc,iec,ied=',hi%isdB,hi%iscB,hi%iecB,hi%iedB
930  write(0,*) 'chksum_u_2d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
931  call chksum_error(fatal,'Error in chksum_u_2d '//trim(mesg))
932  endif
933 
934  bc0 = subchk(array, hi, 0, 0, scaling)
935 
936  sym = .false. ; if (present(symmetric)) sym = symmetric
937 
938  if ((hshift==0) .and. .not.sym) then
939  if (is_root_pe()) call chk_sum_msg("u-point:", bc0, mesg, iounit)
940  return
941  endif
942 
943  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
944 
945  if (hshift==0) then
946  bcw = subchk(array, hi, -hshift-1, 0, scaling)
947  if (is_root_pe()) call chk_sum_msg_w("u-point:", bc0, bcw, mesg, iounit)
948  elseif (do_corners) then
949  if (sym) then
950  bcsw = subchk(array, hi, -hshift-1, -hshift, scaling)
951  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
952  else
953  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
954  bcnw = subchk(array, hi, -hshift, hshift, scaling)
955  endif
956  bcse = subchk(array, hi, hshift, -hshift, scaling)
957  bcne = subchk(array, hi, hshift, hshift, scaling)
958 
959  if (is_root_pe()) &
960  call chk_sum_msg("u-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
961  else
962  bcs = subchk(array, hi, 0, -hshift, scaling)
963  bce = subchk(array, hi, hshift, 0, scaling)
964  if (sym) then
965  bcw = subchk(array, hi, -hshift-1, 0, scaling)
966  else
967  bcw = subchk(array, hi, -hshift, 0, scaling)
968  endif
969  bcn = subchk(array, hi, 0, hshift, scaling)
970 
971  if (is_root_pe()) &
972  call chk_sum_msg_nsew("u-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
973  endif
974 
975  contains
976 
977  integer function subchk(array, HI, di, dj, scale)
978  type(hor_index_type), intent(in) :: HI !< A horizontal index type
979  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array !< The array to be checksummed
980  integer, intent(in) :: di !< i- direction array shift for this checksum
981  integer, intent(in) :: dj !< j- direction array shift for this checksum
982  real, intent(in) :: scale !< A scaling factor for this array.
983  integer :: i, j, bc
984  subchk = 0
985  ! This line deliberately uses the h-point computational domain.
986  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
987  bc = bitcount(abs(scale*array(i,j)))
988  subchk = subchk + bc
989  enddo ; enddo
990  call sum_across_pes(subchk)
991  subchk=mod(subchk, bc_modulus)
992  end function subchk
993 
994  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
995  type(hor_index_type), intent(in) :: HI !< A horizontal index type
996  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array !< The array to be checksummed
997  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
998  !! full symmetric computational domain.
999  real, intent(out) :: aMean !< Array mean
1000  real, intent(out) :: aMin !< Array minimum
1001  real, intent(out) :: aMax !< Array maximum
1002 
1003  integer :: i, j, n, IsB
1004 
1005  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
1006 
1007  amin = array(hi%isc,hi%jsc) ; amax = amin
1008  do j=hi%jsc,hi%jec ; do i=isb,hi%IecB
1009  amin = min(amin, array(i,j))
1010  amax = max(amax, array(i,j))
1011  enddo ; enddo
1012  ! This line deliberately uses the h-point computational domain.
1013  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
1014  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
1015  call sum_across_pes(n)
1016  call min_across_pes(amin)
1017  call max_across_pes(amax)
1018  amean = amean / real(n)
1019  end subroutine substats
1020 

◆ chksum_u_3d()

subroutine mom_checksums::uchksum::chksum_u_3d ( real, dimension(hi_m%isdb:,hi_m%jsd:,:), intent(in), target  array_m,
character(len=*), intent(in)  mesg,
type(hor_index_type), intent(in), target  HI_m,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  symmetric,
logical, intent(in), optional  omit_corners,
real, intent(in), optional  scale,
integer, intent(in), optional  logunit 
)
private

Checksums a 3d array staggered at C-grid u points.

Parameters
[in]hi_mA horizontal index type
[in]array_mThe array to be checksummed
[in]mesgAn identifying message
[in]haloshiftThe width of halos to check (default 0)
[in]symmetricIf true, do the checksums on the full symmetric computational domain.
[in]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.
[in]logunitIO unit for checksum logging

Definition at line 1524 of file MOM_checksums.F90.

1526  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
1527  real, dimension(HI_m%isdB:,HI_m%Jsd:,:), target, intent(in) :: array_m !< The array to be checksummed
1528  character(len=*), intent(in) :: mesg !< An identifying message
1529  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1530  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1531  !! symmetric computational domain.
1532  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1533  real, optional, intent(in) :: scale !< A scaling factor for this array.
1534  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1535 
1536  real, pointer :: array(:,:,:) ! Field array on the input grid
1537  real, allocatable, dimension(:,:,:) :: rescaled_array
1538  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
1539  real :: scaling
1540  integer :: iounit !< Log IO unit
1541  integer :: i, j, k, Is
1542  real :: aMean, aMin, aMax
1543  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1544  integer :: bcN, bcS, bcE, bcW
1545  logical :: do_corners, sym, sym_stats
1546  integer :: turns ! Quarter turns from input to model grid
1547 
1548  ! Rotate array to the input grid
1549  turns = hi_m%turns
1550  if (modulo(turns, 4) /= 0) then
1551  allocate(hi)
1552  call rotate_hor_index(hi_m, -turns, hi)
1553  if (modulo(turns, 2) /= 0) then
1554  ! Arrays originating from v-points must be handled by vchksum
1555  allocate(array(hi%isd:hi%ied, hi%JsdB:hi%JedB, size(array_m, 3)))
1556  call rotate_array(array_m, -turns, array)
1557  call vchksum(array, mesg, hi, haloshift, symmetric, omit_corners, scale, logunit)
1558  return
1559  else
1560  allocate(array(hi%IsdB:hi%IedB, hi%jsd:hi%jed, size(array_m, 3)))
1561  call rotate_array(array_m, -turns, array)
1562  endif
1563  else
1564  hi => hi_m
1565  array => array_m
1566  endif
1567 
1568  if (checkfornans) then
1569  if (is_nan(array(hi%IscB:hi%IecB,hi%jsc:hi%jec,:))) &
1570  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1571 ! if (is_NaN(array)) &
1572 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1573  endif
1574 
1575  scaling = 1.0 ; if (present(scale)) scaling = scale
1576  iounit = error_unit; if(present(logunit)) iounit = logunit
1577  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1578  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1579 
1580  if (calculatestatistics) then
1581  if (present(scale)) then
1582  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1583  lbound(array,2):ubound(array,2), &
1584  lbound(array,3):ubound(array,3)) )
1585  rescaled_array(:,:,:) = 0.0
1586  is = hi%isc ; if (sym_stats) is = hi%isc-1
1587  do k=1,size(array,3) ; do j=hi%jsc,hi%jec ; do i=is,hi%IecB
1588  rescaled_array(i,j,k) = scale*array(i,j,k)
1589  enddo ; enddo ; enddo
1590  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1591  deallocate(rescaled_array)
1592  else
1593  call substats(hi, array, sym_stats, amean, amin, amax)
1594  endif
1595  if (is_root_pe()) &
1596  call chk_sum_msg("u-point:", amean, amin, amax, mesg, iounit)
1597  endif
1598 
1599  if (.not.writechksums) return
1600 
1601  hshift = default_shift
1602  if (present(haloshift)) hshift = haloshift
1603  if (hshift<0) hshift = hi%ied-hi%iec
1604 
1605  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1606  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1607  write(0,*) 'chksum_u_3d: haloshift =',hshift
1608  write(0,*) 'chksum_u_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1609  write(0,*) 'chksum_u_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1610  call chksum_error(fatal,'Error in chksum_u_3d '//trim(mesg))
1611  endif
1612 
1613  bc0 = subchk(array, hi, 0, 0, scaling)
1614 
1615  sym = .false. ; if (present(symmetric)) sym = symmetric
1616 
1617  if ((hshift==0) .and. .not.sym) then
1618  if (is_root_pe()) call chk_sum_msg("u-point:", bc0, mesg, iounit)
1619  return
1620  endif
1621 
1622  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1623 
1624  if (hshift==0) then
1625  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1626  if (is_root_pe()) call chk_sum_msg_w("u-point:", bc0, bcw, mesg, iounit)
1627  elseif (do_corners) then
1628  if (sym) then
1629  bcsw = subchk(array, hi, -hshift-1, -hshift, scaling)
1630  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
1631  else
1632  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1633  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1634  endif
1635  bcse = subchk(array, hi, hshift, -hshift, scaling)
1636  bcne = subchk(array, hi, hshift, hshift, scaling)
1637 
1638  if (is_root_pe()) &
1639  call chk_sum_msg("u-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1640  else
1641  bcs = subchk(array, hi, 0, -hshift, scaling)
1642  bce = subchk(array, hi, hshift, 0, scaling)
1643  if (sym) then
1644  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1645  else
1646  bcw = subchk(array, hi, -hshift, 0, scaling)
1647  endif
1648  bcn = subchk(array, hi, 0, hshift, scaling)
1649 
1650  if (is_root_pe()) &
1651  call chk_sum_msg_nsew("u-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1652  endif
1653 
1654  contains
1655 
1656  integer function subchk(array, HI, di, dj, scale)
1657  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1658  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
1659  integer, intent(in) :: di !< i- direction array shift for this checksum
1660  integer, intent(in) :: dj !< j- direction array shift for this checksum
1661  real, intent(in) :: scale !< A scaling factor for this array.
1662  integer :: i, j, k, bc
1663  subchk = 0
1664  ! This line deliberately uses the h-point computational domain.
1665  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1666  bc = bitcount(abs(scale*array(i,j,k)))
1667  subchk = subchk + bc
1668  enddo ; enddo ; enddo
1669  call sum_across_pes(subchk)
1670  subchk=mod(subchk, bc_modulus)
1671  end function subchk
1672 
1673  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1674  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1675  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
1676  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1677  !! full symmetric computational domain.
1678  real, intent(out) :: aMean !< Array mean
1679  real, intent(out) :: aMin !< Array minimum
1680  real, intent(out) :: aMax !< Array maximum
1681 
1682  integer :: i, j, k, n, IsB
1683 
1684  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
1685 
1686  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1687  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc,hi%jec ; do i=isb,hi%IecB
1688  amin = min(amin, array(i,j,k))
1689  amax = max(amax, array(i,j,k))
1690  enddo ; enddo ; enddo
1691  ! This line deliberately uses the h-point computational domain.
1692  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1693  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1694  call sum_across_pes(n)
1695  call min_across_pes(amin)
1696  call max_across_pes(amax)
1697  amean = amean / real(n)
1698  end subroutine substats
1699 

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