MOM6
mom_checksums::vchksum Interface Reference

Detailed Description

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

Definition at line 38 of file MOM_checksums.F90.

Private functions

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

Functions and subroutines

◆ chksum_v_2d()

subroutine mom_checksums::vchksum::chksum_v_2d ( real, dimension(hi%isd:,hi%jsdb:), intent(in)  array,
character(len=*), intent(in)  mesg,
type(hor_index_type), intent(in)  HI,
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]hiA horizontal index type
[in]arrayThe 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 799 of file MOM_checksums.F90.

799  type(hor_index_type), intent(in) :: HI !< A horizontal index type
800  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
801  character(len=*), intent(in) :: mesg !< An identifying message
802  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
803  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
804  !! symmetric computational domain.
805  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
806  real, optional, intent(in) :: scale !< A scaling factor for this array.
807  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
808 
809  real, allocatable, dimension(:,:) :: rescaled_array
810  real :: scaling
811  integer :: iounit !< Log IO unit
812  integer :: i, j, Js
813  real :: aMean, aMin, aMax
814  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
815  integer :: bcN, bcS, bcE, bcW
816  logical :: do_corners, sym, sym_stats
817 
818  if (checkfornans) then
819  if (is_nan(array(hi%isc:hi%iec,hi%JscB:hi%JecB))) &
820  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
821 ! if (is_NaN(array)) &
822 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
823  endif
824 
825  scaling = 1.0 ; if (present(scale)) scaling = scale
826  iounit = error_unit; if(present(logunit)) iounit = logunit
827  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
828  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
829 
830  if (calculatestatistics) then
831  if (present(scale)) then
832  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
833  lbound(array,2):ubound(array,2)) )
834  rescaled_array(:,:) = 0.0
835  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
836  do j=js,hi%JecB ; do i=hi%isc,hi%iec
837  rescaled_array(i,j) = scale*array(i,j)
838  enddo ; enddo
839  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
840  deallocate(rescaled_array)
841  else
842  call substats(hi, array, sym_stats, amean, amin, amax)
843  endif
844 
845  if (is_root_pe()) &
846  call chk_sum_msg("v-point:", amean, amin, amax, mesg, iounit)
847  endif
848 
849  if (.not.writechksums) return
850 
851  hshift = default_shift
852  if (present(haloshift)) hshift = haloshift
853  if (hshift<0) hshift = hi%ied-hi%iec
854 
855  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
856  hi%jscB-hshift<hi%jsdB .or. hi%jecB+hshift>hi%jedB ) then
857  write(0,*) 'chksum_v_2d: haloshift =',hshift
858  write(0,*) 'chksum_v_2d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
859  write(0,*) 'chksum_v_2d: jsd,jsc,jec,jed=',hi%jsdB,hi%jscB,hi%jecB,hi%jedB
860  call chksum_error(fatal,'Error in chksum_v_2d '//trim(mesg))
861  endif
862 
863  bc0 = subchk(array, hi, 0, 0, scaling)
864 
865  sym = .false. ; if (present(symmetric)) sym = symmetric
866 
867  if ((hshift==0) .and. .not.sym) then
868  if (is_root_pe()) call chk_sum_msg("v-point:", bc0, mesg, iounit)
869  return
870  endif
871 
872  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
873 
874  if (hshift==0) then
875  bcs = subchk(array, hi, 0, -hshift-1, scaling)
876  if (is_root_pe()) call chk_sum_msg_s("v-point:", bc0, bcs, mesg, iounit)
877  elseif (do_corners) then
878  if (sym) then
879  bcsw = subchk(array, hi, -hshift, -hshift-1, scaling)
880  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
881  else
882  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
883  bcse = subchk(array, hi, hshift, -hshift, scaling)
884  endif
885  bcnw = subchk(array, hi, -hshift, hshift, scaling)
886  bcne = subchk(array, hi, hshift, hshift, scaling)
887 
888  if (is_root_pe()) &
889  call chk_sum_msg("v-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
890  else
891  if (sym) then
892  bcs = subchk(array, hi, 0, -hshift-1, scaling)
893  else
894  bcs = subchk(array, hi, 0, -hshift, scaling)
895  endif
896  bce = subchk(array, hi, hshift, 0, scaling)
897  bcw = subchk(array, hi, -hshift, 0, scaling)
898  bcn = subchk(array, hi, 0, hshift, scaling)
899 
900  if (is_root_pe()) &
901  call chk_sum_msg_nsew("v-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
902  endif
903 
904  contains
905 
906  integer function subchk(array, HI, di, dj, scale)
907  type(hor_index_type), intent(in) :: HI !< A horizontal index type
908  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
909  integer, intent(in) :: di !< i- direction array shift for this checksum
910  integer, intent(in) :: dj !< j- direction array shift for this checksum
911  real, intent(in) :: scale !< A scaling factor for this array.
912  integer :: i, j, bc
913  subchk = 0
914  ! This line deliberately uses the h-point computational domain.
915  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
916  bc = bitcount(abs(scale*array(i,j)))
917  subchk = subchk + bc
918  enddo ; enddo
919  call sum_across_pes(subchk)
920  subchk=mod(subchk, bc_modulus)
921  end function subchk
922 
923  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
924  type(hor_index_type), intent(in) :: HI !< A horizontal index type
925  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
926  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
927  !! full symmetric computational domain.
928  real, intent(out) :: aMean, aMin, aMax
929 
930  integer :: i, j, n, JsB
931 
932  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
933 
934  amin = array(hi%isc,hi%jsc) ; amax = amin
935  do j=jsb,hi%JecB ; do i=hi%isc,hi%iec
936  amin = min(amin, array(i,j))
937  amax = max(amax, array(i,j))
938  enddo ; enddo
939  ! This line deliberately uses the h-computational domain.
940  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
941  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
942  call sum_across_pes(n)
943  call min_across_pes(amin)
944  call max_across_pes(amax)
945  amean = amean / real(n)
946  end subroutine substats
947 

◆ chksum_v_3d()

subroutine mom_checksums::vchksum::chksum_v_3d ( real, dimension(hi%isd:,hi%jsdb:,:), intent(in)  array,
character(len=*), intent(in)  mesg,
type(hor_index_type), intent(in)  HI,
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]hiA horizontal index type
[in]arrayThe 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 1394 of file MOM_checksums.F90.

1394  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1395  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1396  character(len=*), intent(in) :: mesg !< An identifying message
1397  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1398  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1399  !! symmetric computational domain.
1400  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1401  real, optional, intent(in) :: scale !< A scaling factor for this array.
1402  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1403 
1404  real, allocatable, dimension(:,:,:) :: rescaled_array
1405  real :: scaling
1406  integer :: iounit !< Log IO unit
1407  integer :: i, j, k, Js
1408  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1409  integer :: bcN, bcS, bcE, bcW
1410  real :: aMean, aMin, aMax
1411  logical :: do_corners, sym, sym_stats
1412 
1413  if (checkfornans) then
1414  if (is_nan(array(hi%isc:hi%iec,hi%JscB:hi%JecB,:))) &
1415  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1416 ! if (is_NaN(array)) &
1417 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1418  endif
1419 
1420  scaling = 1.0 ; if (present(scale)) scaling = scale
1421  iounit = error_unit; if(present(logunit)) iounit = logunit
1422  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1423  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1424 
1425  if (calculatestatistics) then
1426  if (present(scale)) then
1427  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1428  lbound(array,2):ubound(array,2), &
1429  lbound(array,3):ubound(array,3)) )
1430  rescaled_array(:,:,:) = 0.0
1431  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
1432  do k=1,size(array,3) ; do j=js,hi%JecB ; do i=hi%isc,hi%iec
1433  rescaled_array(i,j,k) = scale*array(i,j,k)
1434  enddo ; enddo ; enddo
1435  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1436  deallocate(rescaled_array)
1437  else
1438  call substats(hi, array, sym_stats, amean, amin, amax)
1439  endif
1440  if (is_root_pe()) &
1441  call chk_sum_msg("v-point:", amean, amin, amax, mesg, iounit)
1442  endif
1443 
1444  if (.not.writechksums) return
1445 
1446  hshift = default_shift
1447  if (present(haloshift)) hshift = haloshift
1448  if (hshift<0) hshift = hi%ied-hi%iec
1449 
1450  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1451  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1452  write(0,*) 'chksum_v_3d: haloshift =',hshift
1453  write(0,*) 'chksum_v_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1454  write(0,*) 'chksum_v_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1455  call chksum_error(fatal,'Error in chksum_v_3d '//trim(mesg))
1456  endif
1457 
1458  bc0 = subchk(array, hi, 0, 0, scaling)
1459 
1460  sym = .false. ; if (present(symmetric)) sym = symmetric
1461 
1462  if ((hshift==0) .and. .not.sym) then
1463  if (is_root_pe()) call chk_sum_msg("v-point:", bc0, mesg, iounit)
1464  return
1465  endif
1466 
1467  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1468 
1469  if (hshift==0) then
1470  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1471  if (is_root_pe()) call chk_sum_msg_s("v-point:", bc0, bcs, mesg, iounit)
1472  elseif (do_corners) then
1473  if (sym) then
1474  bcsw = subchk(array, hi, -hshift, -hshift-1, scaling)
1475  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1476  else
1477  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1478  bcse = subchk(array, hi, hshift, -hshift, scaling)
1479  endif
1480  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1481  bcne = subchk(array, hi, hshift, hshift, scaling)
1482 
1483  if (is_root_pe()) &
1484  call chk_sum_msg("v-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1485  else
1486  if (sym) then
1487  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1488  else
1489  bcs = subchk(array, hi, 0, -hshift, scaling)
1490  endif
1491  bce = subchk(array, hi, hshift, 0, scaling)
1492  bcw = subchk(array, hi, -hshift, 0, scaling)
1493  bcn = subchk(array, hi, 0, hshift, scaling)
1494 
1495  if (is_root_pe()) &
1496  call chk_sum_msg_nsew("v-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1497  endif
1498 
1499  contains
1500 
1501  integer function subchk(array, HI, di, dj, scale)
1502  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1503  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1504  integer, intent(in) :: di !< i- direction array shift for this checksum
1505  integer, intent(in) :: dj !< j- direction array shift for this checksum
1506  real, intent(in) :: scale !< A scaling factor for this array.
1507  integer :: i, j, k, bc
1508  subchk = 0
1509  ! This line deliberately uses the h-point computational domain.
1510  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1511  bc = bitcount(abs(scale*array(i,j,k)))
1512  subchk = subchk + bc
1513  enddo ; enddo ; enddo
1514  call sum_across_pes(subchk)
1515  subchk=mod(subchk, bc_modulus)
1516  end function subchk
1517 
1518  !subroutine subStats(HI, array, mesg, sym_stats)
1519  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1520  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1521  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1522  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1523  !! full symmetric computational domain.
1524  real, intent(out) :: aMean, aMin, aMax !< Mean/min/max of array over domain
1525 
1526  integer :: i, j, k, n, JsB
1527 
1528  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
1529 
1530  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1531  do k=lbound(array,3),ubound(array,3) ; do j=jsb,hi%JecB ; do i=hi%isc,hi%iec
1532  amin = min(amin, array(i,j,k))
1533  amax = max(amax, array(i,j,k))
1534  enddo ; enddo ; enddo
1535  ! This line deliberately uses the h-point computational domain.
1536  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1537  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1538  call sum_across_pes(n)
1539  call min_across_pes(amin)
1540  call max_across_pes(amax)
1541  amean = amean / real(n)
1542  end subroutine substats
1543 

The documentation for this interface was generated from the following file:
substats
subroutine substats(array, aMean, aMin, aMax)
Definition: MOM_checksums.F90:179
subchk
integer function subchk(array, scale)
Definition: MOM_checksums.F90:167