MOM6
mom_checksums::vchksum Interface Reference

Detailed Description

Definition at line 48 of file MOM_checksums.F90.

Private functions

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

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

Definition at line 578 of file MOM_checksums.F90.

578  type(hor_index_type), intent(in) :: hi !< A horizontal index type
579  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
580  character(len=*), intent(in) :: mesg !< An identifying message
581  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
582  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
583  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
584  real, optional, intent(in) :: scale !< A scaling factor for this array.
585 
586  real :: scaling
587  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
588  integer :: bcn, bcs, bce, bcw
589  logical :: do_corners, sym, sym_stats
590 
591  if (checkfornans) then
592  if (is_nan(array(hi%isc:hi%iec,hi%JscB:hi%JecB))) &
593  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
594 ! if (is_NaN(array)) &
595 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
596  endif
597  scaling = 1.0 ; if (present(scale)) scaling = scale
598 
599  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
600  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
601  if (calculatestatistics) call substats(hi, array, mesg, sym_stats, scaling)
602 
603  if (.not.writechksums) return
604 
605  hshift = default_shift
606  if (present(haloshift)) hshift = haloshift
607  if (hshift<0) hshift = hi%ied-hi%iec
608 
609  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
610  hi%jscB-hshift<hi%jsdB .or. hi%jecB+hshift>hi%jedB ) then
611  write(0,*) 'chksum_v_2d: haloshift =',hshift
612  write(0,*) 'chksum_v_2d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
613  write(0,*) 'chksum_v_2d: jsd,jsc,jec,jed=',hi%jsdB,hi%jscB,hi%jecB,hi%jedB
614  call chksum_error(fatal,'Error in chksum_v_2d '//trim(mesg))
615  endif
616 
617  bc0 = subchk(array, hi, 0, 0, scaling)
618 
619  sym = .false. ; if (present(symmetric)) sym = symmetric
620 
621  if ((hshift==0) .and. .not.sym) then
622  if (is_root_pe()) call chk_sum_msg("v-point:",bc0,mesg)
623  return
624  endif
625 
626  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
627 
628  if (hshift==0) then
629  bcs = subchk(array, hi, 0, -hshift-1, scaling)
630  if (is_root_pe()) call chk_sum_msg_s("v-point:",bc0,bcs,mesg)
631  elseif (do_corners) then
632  if (sym) then
633  bcsw = subchk(array, hi, -hshift, -hshift-1, scaling)
634  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
635  else
636  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
637  bcse = subchk(array, hi, hshift, -hshift, scaling)
638  endif
639  bcnw = subchk(array, hi, -hshift, hshift, scaling)
640  bcne = subchk(array, hi, hshift, hshift, scaling)
641 
642  if (is_root_pe()) call chk_sum_msg("v-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
643  else
644  if (sym) then
645  bcs = subchk(array, hi, 0, -hshift-1, scaling)
646  else
647  bcs = subchk(array, hi, 0, -hshift, scaling)
648  endif
649  bce = subchk(array, hi, hshift, 0, scaling)
650  bcw = subchk(array, hi, -hshift, 0, scaling)
651  bcn = subchk(array, hi, 0, hshift, scaling)
652 
653  if (is_root_pe()) call chk_sum_msg_nsew("v-point:",bc0,bcn,bcs,bce,bcw,mesg)
654  endif
655 
656  contains
657 
658  integer function subchk(array, HI, di, dj, scale)
659  type(hor_index_type), intent(in) :: hi
660  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array
661  integer, intent(in) :: di, dj
662  real, intent(in) :: scale
663  integer :: bitcount, i, j, bc
664  subchk = 0
665  ! This line deliberately uses the h-point computational domain.
666  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
667  bc = bitcount(abs(scale*array(i,j)))
668  subchk = subchk + bc
669  enddo; enddo
670  call sum_across_pes(subchk)
671  subchk=mod(subchk,1000000000)
672  end function subchk
673 
674  subroutine substats(HI, array, mesg, sym_stats, scale)
675  type(hor_index_type), intent(in) :: hi
676  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: array
677  character(len=*), intent(in) :: mesg
678  logical, intent(in) :: sym_stats
679  real, intent(in) :: scale
680  integer :: i, j, n, jsb
681  real :: amean, amin, amax
682 
683  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
684 
685  amin = array(hi%isc,hi%jsc) ; amax = amin
686  do j=jsb,hi%JecB ; do i=hi%isc,hi%iec
687  amin = min(amin, array(i,j))
688  amax = max(amax, array(i,j))
689  enddo ; enddo
690  ! This line deliberately uses the h-computational domain.
691  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
692  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
693  call sum_across_pes(n)
694  call min_across_pes(amin)
695  call max_across_pes(amax)
696  amean = amean / real(n)
697  if (is_root_pe()) call chk_sum_msg("v-point:",amean*scale,amin*scale,amax*scale,mesg)
698  end subroutine substats
699 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)

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

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

Definition at line 1069 of file MOM_checksums.F90.

1069  type(hor_index_type), intent(in) :: hi !< A horizontal index type
1070  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1071  character(len=*), intent(in) :: mesg !< An identifying message
1072  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1073  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
1074  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1075  real, optional, intent(in) :: scale !< A scaling factor for this array.
1076 
1077  real :: scaling
1078  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
1079  integer :: bcn, bcs, bce, bcw
1080  logical :: do_corners, sym, sym_stats
1081 
1082  if (checkfornans) then
1083  if (is_nan(array(hi%isc:hi%iec,hi%JscB:hi%JecB,:))) &
1084  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1085 ! if (is_NaN(array)) &
1086 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1087  endif
1088  scaling = 1.0 ; if (present(scale)) scaling = scale
1089 
1090  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1091  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1092  if (calculatestatistics) call substats(hi, array, mesg, sym_stats, scaling)
1093 
1094  if (.not.writechksums) return
1095 
1096  hshift = default_shift
1097  if (present(haloshift)) hshift = haloshift
1098  if (hshift<0) hshift = hi%ied-hi%iec
1099 
1100  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1101  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1102  write(0,*) 'chksum_v_3d: haloshift =',hshift
1103  write(0,*) 'chksum_v_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1104  write(0,*) 'chksum_v_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1105  call chksum_error(fatal,'Error in chksum_v_3d '//trim(mesg))
1106  endif
1107 
1108  bc0 = subchk(array, hi, 0, 0, scaling)
1109 
1110  sym = .false. ; if (present(symmetric)) sym = symmetric
1111 
1112  if ((hshift==0) .and. .not.sym) then
1113  if (is_root_pe()) call chk_sum_msg("v-point:",bc0,mesg)
1114  return
1115  endif
1116 
1117  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1118 
1119  if (hshift==0) then
1120  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1121  if (is_root_pe()) call chk_sum_msg_s("v-point:",bc0,bcs,mesg)
1122  elseif (do_corners) then
1123  if (sym) then
1124  bcsw = subchk(array, hi, -hshift, -hshift-1, scaling)
1125  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1126  else
1127  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1128  bcse = subchk(array, hi, hshift, -hshift, scaling)
1129  endif
1130  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1131  bcne = subchk(array, hi, hshift, hshift, scaling)
1132 
1133  if (is_root_pe()) call chk_sum_msg("v-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
1134  else
1135  if (sym) then
1136  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1137  else
1138  bcs = subchk(array, hi, 0, -hshift, scaling)
1139  endif
1140  bce = subchk(array, hi, hshift, 0, scaling)
1141  bcw = subchk(array, hi, -hshift, 0, scaling)
1142  bcn = subchk(array, hi, 0, hshift, scaling)
1143 
1144  if (is_root_pe()) call chk_sum_msg_nsew("v-point:",bc0,bcn,bcs,bce,bcw,mesg)
1145  endif
1146 
1147  contains
1148 
1149  integer function subchk(array, HI, di, dj, scale)
1150  type(hor_index_type), intent(in) :: hi
1151  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array
1152  integer, intent(in) :: di, dj
1153  real, intent(in) :: scale
1154  integer :: bitcount, i, j, k, bc
1155  subchk = 0
1156  ! This line deliberately uses the h-point computational domain.
1157  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1158  bc = bitcount(abs(scale*array(i,j,k)))
1159  subchk = subchk + bc
1160  enddo ; enddo ; enddo
1161  call sum_across_pes(subchk)
1162  subchk=mod(subchk,1000000000)
1163  end function subchk
1164 
1165  subroutine substats(HI, array, mesg, sym_stats, scale)
1166  type(hor_index_type), intent(in) :: hi
1167  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: array
1168  character(len=*), intent(in) :: mesg
1169  logical, intent(in) :: sym_stats
1170  real, intent(in) :: scale
1171  integer :: i, j, k, n, jsb
1172  real :: amean, amin, amax
1173 
1174  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
1175 
1176  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1177  do k=lbound(array,3),ubound(array,3) ; do j=jsb,hi%JecB ; do i=hi%isc,hi%iec
1178  amin = min(amin, array(i,j,k))
1179  amax = max(amax, array(i,j,k))
1180  enddo ; enddo ; enddo
1181  ! This line deliberately uses the h-point computational domain.
1182  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1183  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1184  call sum_across_pes(n)
1185  call min_across_pes(amin)
1186  call max_across_pes(amax)
1187  amean = amean / real(n)
1188  if (is_root_pe()) call chk_sum_msg("v-point:",amean*scale,amin*scale,amax*scale,mesg)
1189  end subroutine substats
1190 
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: