MOM6
mom_checksums::uchksum Interface Reference

Detailed Description

Definition at line 44 of file MOM_checksums.F90.

Private functions

subroutine chksum_u_2d (array, mesg, HI, haloshift, symmetric, omit_corners, scale)
 chksum_u_2d performs checksums on a 2d array staggered at C-grid u points. More...
 
subroutine chksum_u_3d (array, mesg, HI, haloshift, symmetric, omit_corners, scale)
 chksum_u_3d performs checksums on 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 
)
private

chksum_u_2d performs checksums on 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.

Definition at line 450 of file MOM_checksums.F90.

450  type(hor_index_type), intent(in) :: hi !< A horizontal index type
451  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array !< The array to be checksummed
452  character(len=*), intent(in) :: mesg !< An identifying message
453  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
454  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
455  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
456  real, optional, intent(in) :: scale !< A scaling factor for this array.
457 
458  real :: scaling
459  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
460  integer :: bcn, bcs, bce, bcw
461  logical :: do_corners, sym, sym_stats
462 
463  if (checkfornans) then
464  if (is_nan(array(hi%IscB:hi%IecB,hi%jsc:hi%jec))) &
465  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
466 ! if (is_NaN(array)) &
467 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
468  endif
469  scaling = 1.0 ; if (present(scale)) scaling = scale
470 
471  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
472  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
473  if (calculatestatistics) call substats(hi, array, mesg, sym_stats, scaling)
474 
475  if (.not.writechksums) return
476 
477  hshift = default_shift
478  if (present(haloshift)) hshift = haloshift
479  if (hshift<0) hshift = hi%iedB-hi%iecB
480 
481  if ( hi%iscB-hshift<hi%isdB .or. hi%iecB+hshift>hi%iedB .or. &
482  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
483  write(0,*) 'chksum_u_2d: haloshift =',hshift
484  write(0,*) 'chksum_u_2d: isd,isc,iec,ied=',hi%isdB,hi%iscB,hi%iecB,hi%iedB
485  write(0,*) 'chksum_u_2d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
486  call chksum_error(fatal,'Error in chksum_u_2d '//trim(mesg))
487  endif
488 
489  bc0 = subchk(array, hi, 0, 0, scaling)
490 
491  sym = .false. ; if (present(symmetric)) sym = symmetric
492 
493  if ((hshift==0) .and. .not.sym) then
494  if (is_root_pe()) call chk_sum_msg("u-point:",bc0,mesg)
495  return
496  endif
497 
498  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
499 
500  if (hshift==0) then
501  bcw = subchk(array, hi, -hshift-1, 0, scaling)
502  if (is_root_pe()) call chk_sum_msg_w("u-point:",bc0,bcw,mesg)
503  elseif (do_corners) then
504  if (sym) then
505  bcsw = subchk(array, hi, -hshift-1, -hshift, scaling)
506  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
507  else
508  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
509  bcnw = subchk(array, hi, -hshift, hshift, scaling)
510  endif
511  bcse = subchk(array, hi, hshift, -hshift, scaling)
512  bcne = subchk(array, hi, hshift, hshift, scaling)
513 
514  if (is_root_pe()) call chk_sum_msg("u-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
515  else
516  bcs = subchk(array, hi, 0, -hshift, scaling)
517  bce = subchk(array, hi, hshift, 0, scaling)
518  if (sym) then
519  bcw = subchk(array, hi, -hshift-1, 0, scaling)
520  else
521  bcw = subchk(array, hi, -hshift, 0, scaling)
522  endif
523  bcn = subchk(array, hi, 0, hshift, scaling)
524 
525  if (is_root_pe()) call chk_sum_msg_nsew("u-point:",bc0,bcn,bcs,bce,bcw,mesg)
526  endif
527 
528  contains
529 
530  integer function subchk(array, HI, di, dj, scale)
531  type(hor_index_type), intent(in) :: hi
532  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array
533  integer, intent(in) :: di, dj
534  real, intent(in) :: scale
535  integer :: bitcount, i, j, bc
536  subchk = 0
537  ! This line deliberately uses the h-point computational domain.
538  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
539  bc = bitcount(abs(scale*array(i,j)))
540  subchk = subchk + bc
541  enddo; enddo
542  call sum_across_pes(subchk)
543  subchk=mod(subchk,1000000000)
544  end function subchk
545 
546  subroutine substats(HI, array, mesg, sym_stats, scale)
547  type(hor_index_type), intent(in) :: hi
548  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: array
549  character(len=*), intent(in) :: mesg
550  logical, intent(in) :: sym_stats
551  real, intent(in) :: scale
552  integer :: i, j, n, isb
553  real :: amean, amin, amax
554 
555  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
556 
557  amin = array(hi%isc,hi%jsc) ; amax = amin
558  do j=hi%jsc,hi%jec ; do i=isb,hi%IecB
559  amin = min(amin, array(i,j))
560  amax = max(amax, array(i,j))
561  enddo ; enddo
562  ! This line deliberately uses the h-point computational domain.
563  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
564  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
565  call sum_across_pes(n)
566  call min_across_pes(amin)
567  call max_across_pes(amax)
568  amean = amean / real(n)
569  if (is_root_pe()) call chk_sum_msg("u-point:",amean*scale,amin*scale,amax*scale,mesg)
570  end subroutine substats
571 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)

◆ 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 
)
private

chksum_u_3d performs checksums on 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.

Definition at line 941 of file MOM_checksums.F90.

941  type(hor_index_type), intent(in) :: hi !< A horizontal index type
942  real, dimension(HI%isdB:,HI%Jsd:,:), intent(in) :: array !< The array to be checksummed
943  character(len=*), intent(in) :: mesg !< An identifying message
944  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
945  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
946  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
947  real, optional, intent(in) :: scale !< A scaling factor for this array.
948 
949  real :: scaling
950  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
951  integer :: bcn, bcs, bce, bcw
952  logical :: do_corners, sym, sym_stats
953 
954  if (checkfornans) then
955  if (is_nan(array(hi%IscB:hi%IecB,hi%jsc:hi%jec,:))) &
956  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
957 ! if (is_NaN(array)) &
958 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
959  endif
960  scaling = 1.0 ; if (present(scale)) scaling = scale
961 
962  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
963  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
964  if (calculatestatistics) call substats(hi, array, mesg, sym_stats, scaling)
965 
966  if (.not.writechksums) return
967 
968  hshift = default_shift
969  if (present(haloshift)) hshift = haloshift
970  if (hshift<0) hshift = hi%ied-hi%iec
971 
972  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
973  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
974  write(0,*) 'chksum_u_3d: haloshift =',hshift
975  write(0,*) 'chksum_u_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
976  write(0,*) 'chksum_u_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
977  call chksum_error(fatal,'Error in chksum_u_3d '//trim(mesg))
978  endif
979 
980  bc0 = subchk(array, hi, 0, 0, scaling)
981 
982  sym = .false. ; if (present(symmetric)) sym = symmetric
983 
984  if ((hshift==0) .and. .not.sym) then
985  if (is_root_pe()) call chk_sum_msg("u-point:",bc0,mesg)
986  return
987  endif
988 
989  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
990 
991  if (hshift==0) then
992  bcw = subchk(array, hi, -hshift-1, 0, scaling)
993  if (is_root_pe()) call chk_sum_msg_w("u-point:",bc0,bcw,mesg)
994  elseif (do_corners) then
995  if (sym) then
996  bcsw = subchk(array, hi, -hshift-1, -hshift, scaling)
997  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
998  else
999  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1000  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1001  endif
1002  bcse = subchk(array, hi, hshift, -hshift, scaling)
1003  bcne = subchk(array, hi, hshift, hshift, scaling)
1004 
1005  if (is_root_pe()) call chk_sum_msg("u-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
1006  else
1007  bcs = subchk(array, hi, 0, -hshift, scaling)
1008  bce = subchk(array, hi, hshift, 0, scaling)
1009  if (sym) then
1010  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1011  else
1012  bcw = subchk(array, hi, -hshift, 0, scaling)
1013  endif
1014  bcn = subchk(array, hi, 0, hshift, scaling)
1015 
1016  if (is_root_pe()) call chk_sum_msg_nsew("u-point:",bc0,bcn,bcs,bce,bcw,mesg)
1017  endif
1018 
1019  contains
1020 
1021  integer function subchk(array, HI, di, dj, scale)
1022  type(hor_index_type), intent(in) :: hi
1023  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: array
1024  integer, intent(in) :: di, dj
1025  real, intent(in) :: scale
1026  integer :: bitcount, i, j, k, bc
1027  subchk = 0
1028  ! This line deliberately uses the h-point computational domain.
1029  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1030  bc = bitcount(abs(scale*array(i,j,k)))
1031  subchk = subchk + bc
1032  enddo ; enddo ; enddo
1033  call sum_across_pes(subchk)
1034  subchk=mod(subchk,1000000000)
1035  end function subchk
1036 
1037  subroutine substats(HI, array, mesg, sym_stats, scale)
1038  type(hor_index_type), intent(in) :: hi
1039  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: array
1040  character(len=*), intent(in) :: mesg
1041  logical, intent(in) :: sym_stats
1042  real, intent(in) :: scale
1043  integer :: i, j, k, n, isb
1044  real :: amean, amin, amax
1045 
1046  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
1047 
1048  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1049  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc,hi%jec ; do i=isb,hi%IecB
1050  amin = min(amin, array(i,j,k))
1051  amax = max(amax, array(i,j,k))
1052  enddo ; enddo ; enddo
1053  ! This line deliberately uses the h-point computational domain.
1054  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1055  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1056  call sum_across_pes(n)
1057  call min_across_pes(amin)
1058  call max_across_pes(amax)
1059  amean = amean / real(n)
1060  if (is_root_pe()) call chk_sum_msg("u-point:",amean*scale,amin*scale,amax*scale,mesg)
1061  end subroutine substats
1062 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)

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