MOM6
mom_checksums::vchksum Interface Reference

Detailed Description

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

Definition at line 39 of file MOM_checksums.F90.

Private functions

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

Detailed Description

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

Definition at line 39 of file MOM_checksums.F90.

Functions and subroutines

◆ chksum_v_2d()

subroutine mom_checksums::vchksum::chksum_v_2d ( real, dimension(hi_m%isd:,hi_m%jsdb:), 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 v 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 1026 of file MOM_checksums.F90.

1026  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
1027  real, dimension(HI_m%isd:,HI_m%JsdB:), target, intent(in) :: array_m !< The array to be checksummed
1028  character(len=*), intent(in) :: mesg !< An identifying message
1029  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1030  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1031  !! symmetric computational domain.
1032  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1033  real, optional, intent(in) :: scale !< A scaling factor for this array.
1034  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1035 
1036  real, pointer :: array(:,:) ! Field array on the input grid
1037  real, allocatable, dimension(:,:) :: rescaled_array
1038  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
1039  real :: scaling
1040  integer :: iounit !< Log IO unit
1041  integer :: i, j, Js
1042  real :: aMean, aMin, aMax
1043  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1044  integer :: bcN, bcS, bcE, bcW
1045  logical :: do_corners, sym, sym_stats
1046  integer :: turns ! Quarter turns from input to model grid
1047 
1048  ! Rotate array to the input grid
1049  turns = hi_m%turns
1050  if (modulo(turns, 4) /= 0) then
1051  allocate(hi)
1052  call rotate_hor_index(hi_m, -turns, hi)
1053  if (modulo(turns, 2) /= 0) then
1054  ! Arrays originating from u-points must be handled by uchksum
1055  allocate(array(hi%IsdB:hi%IedB, hi%jsd:hi%jed))
1056  call rotate_array(array_m, -turns, array)
1057  call uchksum(array, mesg, hi, haloshift, symmetric, omit_corners, scale, logunit)
1058  return
1059  else
1060  allocate(array(hi%isd:hi%ied, hi%JsdB:hi%JedB))
1061  call rotate_array(array_m, -turns, array)
1062  endif
1063  else
1064  hi => hi_m
1065  array => array_m
1066  endif
1067 
1068  if (checkfornans) then
1069  if (is_nan(array(hi%isc:hi%iec,hi%JscB:hi%JecB))) &
1070  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1071 ! if (is_NaN(array)) &
1072 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1073  endif
1074 
1075  scaling = 1.0 ; if (present(scale)) scaling = scale
1076  iounit = error_unit; if(present(logunit)) iounit = logunit
1077  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1078  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1079 
1080  if (calculatestatistics) then
1081  if (present(scale)) then
1082  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1083  lbound(array,2):ubound(array,2)) )
1084  rescaled_array(:,:) = 0.0
1085  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
1086  do j=js,hi%JecB ; do i=hi%isc,hi%iec
1087  rescaled_array(i,j) = scale*array(i,j)
1088  enddo ; enddo
1089  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1090  deallocate(rescaled_array)
1091  else
1092  call substats(hi, array, sym_stats, amean, amin, amax)
1093  endif
1094 
1095  if (is_root_pe()) &
1096  call chk_sum_msg("v-point:", amean, amin, amax, mesg, iounit)
1097  endif
1098 
1099  if (.not.writechksums) return
1100 
1101  hshift = default_shift
1102  if (present(haloshift)) hshift = haloshift
1103  if (hshift<0) hshift = hi%ied-hi%iec
1104 
1105  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1106  hi%jscB-hshift<hi%jsdB .or. hi%jecB+hshift>hi%jedB ) then
1107  write(0,*) 'chksum_v_2d: haloshift =',hshift
1108  write(0,*) 'chksum_v_2d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1109  write(0,*) 'chksum_v_2d: jsd,jsc,jec,jed=',hi%jsdB,hi%jscB,hi%jecB,hi%jedB
1110  call chksum_error(fatal,'Error in chksum_v_2d '//trim(mesg))
1111  endif
1112 
1113  bc0 = subchk(array, hi, 0, 0, scaling)
1114 
1115  sym = .false. ; if (present(symmetric)) sym = symmetric
1116 
1117  if ((hshift==0) .and. .not.sym) then
1118  if (is_root_pe()) call chk_sum_msg("v-point:", bc0, mesg, iounit)
1119  return
1120  endif
1121 
1122  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1123 
1124  if (hshift==0) then
1125  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1126  if (is_root_pe()) call chk_sum_msg_s("v-point:", bc0, bcs, mesg, iounit)
1127  elseif (do_corners) then
1128  if (sym) then
1129  bcsw = subchk(array, hi, -hshift, -hshift-1, scaling)
1130  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1131  else
1132  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1133  bcse = subchk(array, hi, hshift, -hshift, scaling)
1134  endif
1135  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1136  bcne = subchk(array, hi, hshift, hshift, scaling)
1137 
1138  if (is_root_pe()) &
1139  call chk_sum_msg("v-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1140  else
1141  if (sym) then
1142  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1143  else
1144  bcs = subchk(array, hi, 0, -hshift, scaling)
1145  endif
1146  bce = subchk(array, hi, hshift, 0, scaling)
1147  bcw = subchk(array, hi, -hshift, 0, scaling)
1148  bcn = subchk(array, hi, 0, hshift, scaling)
1149 
1150  if (is_root_pe()) &
1151  call chk_sum_msg_nsew("v-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1152  endif
1153 
1154  contains
1155 
1156  integer function subchk(array, HI, di, dj, scale)
1157  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1158  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
1159  integer, intent(in) :: di !< i- direction array shift for this checksum
1160  integer, intent(in) :: dj !< j- direction array shift for this checksum
1161  real, intent(in) :: scale !< A scaling factor for this array.
1162  integer :: i, j, bc
1163  subchk = 0
1164  ! This line deliberately uses the h-point computational domain.
1165  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
1166  bc = bitcount(abs(scale*array(i,j)))
1167  subchk = subchk + bc
1168  enddo ; enddo
1169  call sum_across_pes(subchk)
1170  subchk=mod(subchk, bc_modulus)
1171  end function subchk
1172 
1173  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1174  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1175  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
1176  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1177  !! full symmetric computational domain.
1178  real, intent(out) :: aMean !< Array mean
1179  real, intent(out) :: aMin !< Array minimum
1180  real, intent(out) :: aMax !< Array maximum
1181 
1182  integer :: i, j, n, JsB
1183 
1184  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
1185 
1186  amin = array(hi%isc,hi%jsc) ; amax = amin
1187  do j=jsb,hi%JecB ; do i=hi%isc,hi%iec
1188  amin = min(amin, array(i,j))
1189  amax = max(amax, array(i,j))
1190  enddo ; enddo
1191  ! This line deliberately uses the h-computational domain.
1192  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
1193  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
1194  call sum_across_pes(n)
1195  call min_across_pes(amin)
1196  call max_across_pes(amax)
1197  amean = amean / real(n)
1198  end subroutine substats
1199 

◆ chksum_v_3d()

subroutine mom_checksums::vchksum::chksum_v_3d ( real, dimension(hi_m%isd:,hi_m%jsdb:,:), 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 v 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 1705 of file MOM_checksums.F90.

1705  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
1706  real, dimension(HI_m%isd:,HI_m%JsdB:,:), target, intent(in) :: array_m !< The array to be checksummed
1707  character(len=*), intent(in) :: mesg !< An identifying message
1708  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1709  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1710  !! symmetric computational domain.
1711  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1712  real, optional, intent(in) :: scale !< A scaling factor for this array.
1713  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1714 
1715  real, pointer :: array(:,:,:) ! Field array on the input grid
1716  real, allocatable, dimension(:,:,:) :: rescaled_array
1717  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
1718  real :: scaling
1719  integer :: iounit !< Log IO unit
1720  integer :: i, j, k, Js
1721  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1722  integer :: bcN, bcS, bcE, bcW
1723  real :: aMean, aMin, aMax
1724  logical :: do_corners, sym, sym_stats
1725  integer :: turns ! Quarter turns from input to model grid
1726 
1727  ! Rotate array to the input grid
1728  turns = hi_m%turns
1729  if (modulo(turns, 4) /= 0) then
1730  allocate(hi)
1731  call rotate_hor_index(hi_m, -turns, hi)
1732  if (modulo(turns, 2) /= 0) then
1733  ! Arrays originating from u-points must be handled by uchksum
1734  allocate(array(hi%IsdB:hi%IedB, hi%jsd:hi%jed, size(array_m, 3)))
1735  call rotate_array(array_m, -turns, array)
1736  call uchksum(array, mesg, hi, haloshift, symmetric, omit_corners, scale, logunit)
1737  return
1738  else
1739  allocate(array(hi%isd:hi%ied, hi%JsdB:hi%JedB, size(array_m, 3)))
1740  call rotate_array(array_m, -turns, array)
1741  endif
1742  else
1743  hi => hi_m
1744  array => array_m
1745  endif
1746 
1747  if (checkfornans) then
1748  if (is_nan(array(hi%isc:hi%iec,hi%JscB:hi%JecB,:))) &
1749  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1750 ! if (is_NaN(array)) &
1751 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1752  endif
1753 
1754  scaling = 1.0 ; if (present(scale)) scaling = scale
1755  iounit = error_unit; if(present(logunit)) iounit = logunit
1756  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1757  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1758 
1759  if (calculatestatistics) then
1760  if (present(scale)) then
1761  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1762  lbound(array,2):ubound(array,2), &
1763  lbound(array,3):ubound(array,3)) )
1764  rescaled_array(:,:,:) = 0.0
1765  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
1766  do k=1,size(array,3) ; do j=js,hi%JecB ; do i=hi%isc,hi%iec
1767  rescaled_array(i,j,k) = scale*array(i,j,k)
1768  enddo ; enddo ; enddo
1769  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1770  deallocate(rescaled_array)
1771  else
1772  call substats(hi, array, sym_stats, amean, amin, amax)
1773  endif
1774  if (is_root_pe()) &
1775  call chk_sum_msg("v-point:", amean, amin, amax, mesg, iounit)
1776  endif
1777 
1778  if (.not.writechksums) return
1779 
1780  hshift = default_shift
1781  if (present(haloshift)) hshift = haloshift
1782  if (hshift<0) hshift = hi%ied-hi%iec
1783 
1784  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1785  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1786  write(0,*) 'chksum_v_3d: haloshift =',hshift
1787  write(0,*) 'chksum_v_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1788  write(0,*) 'chksum_v_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1789  call chksum_error(fatal,'Error in chksum_v_3d '//trim(mesg))
1790  endif
1791 
1792  bc0 = subchk(array, hi, 0, 0, scaling)
1793 
1794  sym = .false. ; if (present(symmetric)) sym = symmetric
1795 
1796  if ((hshift==0) .and. .not.sym) then
1797  if (is_root_pe()) call chk_sum_msg("v-point:", bc0, mesg, iounit)
1798  return
1799  endif
1800 
1801  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1802 
1803  if (hshift==0) then
1804  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1805  if (is_root_pe()) call chk_sum_msg_s("v-point:", bc0, bcs, mesg, iounit)
1806  elseif (do_corners) then
1807  if (sym) then
1808  bcsw = subchk(array, hi, -hshift, -hshift-1, scaling)
1809  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1810  else
1811  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1812  bcse = subchk(array, hi, hshift, -hshift, scaling)
1813  endif
1814  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1815  bcne = subchk(array, hi, hshift, hshift, scaling)
1816 
1817  if (is_root_pe()) &
1818  call chk_sum_msg("v-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1819  else
1820  if (sym) then
1821  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1822  else
1823  bcs = subchk(array, hi, 0, -hshift, scaling)
1824  endif
1825  bce = subchk(array, hi, hshift, 0, scaling)
1826  bcw = subchk(array, hi, -hshift, 0, scaling)
1827  bcn = subchk(array, hi, 0, hshift, scaling)
1828 
1829  if (is_root_pe()) &
1830  call chk_sum_msg_nsew("v-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1831  endif
1832 
1833  contains
1834 
1835  integer function subchk(array, HI, di, dj, scale)
1836  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1837  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1838  integer, intent(in) :: di !< i- direction array shift for this checksum
1839  integer, intent(in) :: dj !< j- direction array shift for this checksum
1840  real, intent(in) :: scale !< A scaling factor for this array.
1841  integer :: i, j, k, bc
1842  subchk = 0
1843  ! This line deliberately uses the h-point computational domain.
1844  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1845  bc = bitcount(abs(scale*array(i,j,k)))
1846  subchk = subchk + bc
1847  enddo ; enddo ; enddo
1848  call sum_across_pes(subchk)
1849  subchk=mod(subchk, bc_modulus)
1850  end function subchk
1851 
1852  !subroutine subStats(HI, array, mesg, sym_stats)
1853  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1854  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1855  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1856  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1857  !! full symmetric computational domain.
1858  real, intent(out) :: aMean !< Mean of array over domain
1859  real, intent(out) :: aMin !< Minimum of array over domain
1860  real, intent(out) :: aMax !< Maximum of array over domain
1861 
1862  integer :: i, j, k, n, JsB
1863 
1864  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
1865 
1866  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1867  do k=lbound(array,3),ubound(array,3) ; do j=jsb,hi%JecB ; do i=hi%isc,hi%iec
1868  amin = min(amin, array(i,j,k))
1869  amax = max(amax, array(i,j,k))
1870  enddo ; enddo ; enddo
1871  ! This line deliberately uses the h-point computational domain.
1872  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1873  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1874  call sum_across_pes(n)
1875  call min_across_pes(amin)
1876  call max_across_pes(amax)
1877  amean = amean / real(n)
1878  end subroutine substats
1879 

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