MOM6
mom_checksums::uchksum Interface Reference

Detailed Description

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

Definition at line 33 of file MOM_checksums.F90.

Private functions

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

Functions and subroutines

◆ chksum_u_2d()

subroutine mom_checksums::uchksum::chksum_u_2d ( real, dimension(hi%isdb:,hi%jsd:), 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 u 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 645 of file MOM_checksums.F90.

645  type(hor_index_type), intent(in) :: HI !< A horizontal index type
646  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array !< The array to be checksummed
647  character(len=*), intent(in) :: mesg !< An identifying message
648  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
649  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
650  !! symmetric computational domain.
651  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
652  real, optional, intent(in) :: scale !< A scaling factor for this array.
653  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
654 
655  real, allocatable, dimension(:,:) :: rescaled_array
656  real :: scaling
657  integer :: iounit !< Log IO unit
658  integer :: i, j, Is
659  real :: aMean, aMin, aMax
660  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
661  integer :: bcN, bcS, bcE, bcW
662  logical :: do_corners, sym, sym_stats
663 
664  if (checkfornans) then
665  if (is_nan(array(hi%IscB:hi%IecB,hi%jsc:hi%jec))) &
666  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
667 ! if (is_NaN(array)) &
668 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
669  endif
670 
671  scaling = 1.0 ; if (present(scale)) scaling = scale
672  iounit = error_unit; if(present(logunit)) iounit = logunit
673  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
674  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
675 
676  if (calculatestatistics) then
677  if (present(scale)) then
678  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
679  lbound(array,2):ubound(array,2)) )
680  rescaled_array(:,:) = 0.0
681  is = hi%isc ; if (sym_stats) is = hi%isc-1
682  do j=hi%jsc,hi%jec ; do i=is,hi%IecB
683  rescaled_array(i,j) = scale*array(i,j)
684  enddo ; enddo
685  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
686  deallocate(rescaled_array)
687  else
688  call substats(hi, array, sym_stats, amean, amin, amax)
689  endif
690 
691  if (is_root_pe()) &
692  call chk_sum_msg("u-point:", amean, amin, amax, mesg, iounit)
693  endif
694 
695  if (.not.writechksums) return
696 
697  hshift = default_shift
698  if (present(haloshift)) hshift = haloshift
699  if (hshift<0) hshift = hi%iedB-hi%iecB
700 
701  if ( hi%iscB-hshift<hi%isdB .or. hi%iecB+hshift>hi%iedB .or. &
702  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
703  write(0,*) 'chksum_u_2d: haloshift =',hshift
704  write(0,*) 'chksum_u_2d: isd,isc,iec,ied=',hi%isdB,hi%iscB,hi%iecB,hi%iedB
705  write(0,*) 'chksum_u_2d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
706  call chksum_error(fatal,'Error in chksum_u_2d '//trim(mesg))
707  endif
708 
709  bc0 = subchk(array, hi, 0, 0, scaling)
710 
711  sym = .false. ; if (present(symmetric)) sym = symmetric
712 
713  if ((hshift==0) .and. .not.sym) then
714  if (is_root_pe()) call chk_sum_msg("u-point:", bc0, mesg, iounit)
715  return
716  endif
717 
718  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
719 
720  if (hshift==0) then
721  bcw = subchk(array, hi, -hshift-1, 0, scaling)
722  if (is_root_pe()) call chk_sum_msg_w("u-point:", bc0, bcw, mesg, iounit)
723  elseif (do_corners) then
724  if (sym) then
725  bcsw = subchk(array, hi, -hshift-1, -hshift, scaling)
726  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
727  else
728  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
729  bcnw = subchk(array, hi, -hshift, hshift, scaling)
730  endif
731  bcse = subchk(array, hi, hshift, -hshift, scaling)
732  bcne = subchk(array, hi, hshift, hshift, scaling)
733 
734  if (is_root_pe()) &
735  call chk_sum_msg("u-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
736  else
737  bcs = subchk(array, hi, 0, -hshift, scaling)
738  bce = subchk(array, hi, hshift, 0, scaling)
739  if (sym) then
740  bcw = subchk(array, hi, -hshift-1, 0, scaling)
741  else
742  bcw = subchk(array, hi, -hshift, 0, scaling)
743  endif
744  bcn = subchk(array, hi, 0, hshift, scaling)
745 
746  if (is_root_pe()) &
747  call chk_sum_msg_nsew("u-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
748  endif
749 
750  contains
751 
752  integer function subchk(array, HI, di, dj, scale)
753  type(hor_index_type), intent(in) :: HI !< A horizontal index type
754  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array !< The array to be checksummed
755  integer, intent(in) :: di !< i- direction array shift for this checksum
756  integer, intent(in) :: dj !< j- direction array shift for this checksum
757  real, intent(in) :: scale !< A scaling factor for this array.
758  integer :: i, j, bc
759  subchk = 0
760  ! This line deliberately uses the h-point computational domain.
761  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
762  bc = bitcount(abs(scale*array(i,j)))
763  subchk = subchk + bc
764  enddo ; enddo
765  call sum_across_pes(subchk)
766  subchk=mod(subchk, bc_modulus)
767  end function subchk
768 
769  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
770  type(hor_index_type), intent(in) :: HI !< A horizontal index type
771  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array !< The array to be checksummed
772  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
773  !! full symmetric computational domain.
774  real, intent(out) :: aMean, aMin, aMax
775 
776  integer :: i, j, n, IsB
777 
778  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
779 
780  amin = array(hi%isc,hi%jsc) ; amax = amin
781  do j=hi%jsc,hi%jec ; do i=isb,hi%IecB
782  amin = min(amin, array(i,j))
783  amax = max(amax, array(i,j))
784  enddo ; enddo
785  ! This line deliberately uses the h-point computational domain.
786  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
787  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
788  call sum_across_pes(n)
789  call min_across_pes(amin)
790  call max_across_pes(amax)
791  amean = amean / real(n)
792  end subroutine substats
793 

◆ chksum_u_3d()

subroutine mom_checksums::uchksum::chksum_u_3d ( real, dimension(hi%isdb:,hi%jsd:,:), 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 u 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 1240 of file MOM_checksums.F90.

1240  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1241  real, dimension(HI%isdB:,HI%Jsd:,:), intent(in) :: array !< The array to be checksummed
1242  character(len=*), intent(in) :: mesg !< An identifying message
1243  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1244  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1245  !! symmetric computational domain.
1246  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1247  real, optional, intent(in) :: scale !< A scaling factor for this array.
1248  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1249 
1250  real, allocatable, dimension(:,:,:) :: rescaled_array
1251  real :: scaling
1252  integer :: iounit !< Log IO unit
1253  integer :: i, j, k, Is
1254  real :: aMean, aMin, aMax
1255  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1256  integer :: bcN, bcS, bcE, bcW
1257  logical :: do_corners, sym, sym_stats
1258 
1259  if (checkfornans) then
1260  if (is_nan(array(hi%IscB:hi%IecB,hi%jsc:hi%jec,:))) &
1261  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1262 ! if (is_NaN(array)) &
1263 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1264  endif
1265 
1266  scaling = 1.0 ; if (present(scale)) scaling = scale
1267  iounit = error_unit; if(present(logunit)) iounit = logunit
1268  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1269  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1270 
1271  if (calculatestatistics) then
1272  if (present(scale)) then
1273  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1274  lbound(array,2):ubound(array,2), &
1275  lbound(array,3):ubound(array,3)) )
1276  rescaled_array(:,:,:) = 0.0
1277  is = hi%isc ; if (sym_stats) is = hi%isc-1
1278  do k=1,size(array,3) ; do j=hi%jsc,hi%jec ; do i=is,hi%IecB
1279  rescaled_array(i,j,k) = scale*array(i,j,k)
1280  enddo ; enddo ; enddo
1281  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1282  deallocate(rescaled_array)
1283  else
1284  call substats(hi, array, sym_stats, amean, amin, amax)
1285  endif
1286  if (is_root_pe()) &
1287  call chk_sum_msg("u-point:", amean, amin, amax, mesg, iounit)
1288  endif
1289 
1290  if (.not.writechksums) return
1291 
1292  hshift = default_shift
1293  if (present(haloshift)) hshift = haloshift
1294  if (hshift<0) hshift = hi%ied-hi%iec
1295 
1296  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1297  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1298  write(0,*) 'chksum_u_3d: haloshift =',hshift
1299  write(0,*) 'chksum_u_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1300  write(0,*) 'chksum_u_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1301  call chksum_error(fatal,'Error in chksum_u_3d '//trim(mesg))
1302  endif
1303 
1304  bc0 = subchk(array, hi, 0, 0, scaling)
1305 
1306  sym = .false. ; if (present(symmetric)) sym = symmetric
1307 
1308  if ((hshift==0) .and. .not.sym) then
1309  if (is_root_pe()) call chk_sum_msg("u-point:", bc0, mesg, iounit)
1310  return
1311  endif
1312 
1313  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1314 
1315  if (hshift==0) then
1316  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1317  if (is_root_pe()) call chk_sum_msg_w("u-point:", bc0, bcw, mesg, iounit)
1318  elseif (do_corners) then
1319  if (sym) then
1320  bcsw = subchk(array, hi, -hshift-1, -hshift, scaling)
1321  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
1322  else
1323  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1324  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1325  endif
1326  bcse = subchk(array, hi, hshift, -hshift, scaling)
1327  bcne = subchk(array, hi, hshift, hshift, scaling)
1328 
1329  if (is_root_pe()) &
1330  call chk_sum_msg("u-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1331  else
1332  bcs = subchk(array, hi, 0, -hshift, scaling)
1333  bce = subchk(array, hi, hshift, 0, scaling)
1334  if (sym) then
1335  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1336  else
1337  bcw = subchk(array, hi, -hshift, 0, scaling)
1338  endif
1339  bcn = subchk(array, hi, 0, hshift, scaling)
1340 
1341  if (is_root_pe()) &
1342  call chk_sum_msg_nsew("u-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1343  endif
1344 
1345  contains
1346 
1347  integer function subchk(array, HI, di, dj, scale)
1348  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1349  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
1350  integer, intent(in) :: di !< i- direction array shift for this checksum
1351  integer, intent(in) :: dj !< j- direction array shift for this checksum
1352  real, intent(in) :: scale !< A scaling factor for this array.
1353  integer :: i, j, k, bc
1354  subchk = 0
1355  ! This line deliberately uses the h-point computational domain.
1356  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1357  bc = bitcount(abs(scale*array(i,j,k)))
1358  subchk = subchk + bc
1359  enddo ; enddo ; enddo
1360  call sum_across_pes(subchk)
1361  subchk=mod(subchk, bc_modulus)
1362  end function subchk
1363 
1364  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1365  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1366  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
1367  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1368  !! full symmetric computational domain.
1369  real, intent(out) :: aMean, aMin, aMax
1370 
1371  integer :: i, j, k, n, IsB
1372 
1373  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
1374 
1375  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1376  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc,hi%jec ; do i=isb,hi%IecB
1377  amin = min(amin, array(i,j,k))
1378  amax = max(amax, array(i,j,k))
1379  enddo ; enddo ; enddo
1380  ! This line deliberately uses the h-point computational domain.
1381  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1382  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1383  call sum_across_pes(n)
1384  call min_across_pes(amin)
1385  call max_across_pes(amax)
1386  amean = amean / real(n)
1387  end subroutine substats
1388 

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