MOM6
mom_checksums::bchksum Interface Reference

Detailed Description

Checksums an array (2d or 3d) staggered at corner points.

Definition at line 53 of file MOM_checksums.F90.

Private functions

subroutine chksum_b_2d (array, mesg, HI, haloshift, symmetric, omit_corners, scale, logunit)
 Checksums a 2d array staggered at corner points. More...
 
subroutine chksum_b_3d (array, mesg, HI, haloshift, symmetric, omit_corners, scale, logunit)
 Checksums a 3d array staggered at corner points. More...
 

Functions and subroutines

◆ chksum_b_2d()

subroutine mom_checksums::bchksum::chksum_b_2d ( real, dimension(hi%isdb:,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 corner 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 439 of file MOM_checksums.F90.

439  type(hor_index_type), intent(in) :: HI !< A horizontal index type
440  real, dimension(HI%IsdB:,HI%JsdB:), &
441  intent(in) :: array !< The array to be checksummed
442  character(len=*), intent(in) :: mesg !< An identifying message
443  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
444  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the
445  !! full symmetric computational domain.
446  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
447  real, optional, intent(in) :: scale !< A scaling factor for this array.
448  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
449 
450  real, allocatable, dimension(:,:) :: rescaled_array
451  real :: scaling
452  integer :: iounit !< Log IO unit
453  integer :: i, j, Is, Js
454  real :: aMean, aMin, aMax
455  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
456  integer :: bcN, bcS, bcE, bcW
457  logical :: do_corners, sym, sym_stats
458 
459  if (checkfornans) then
460  if (is_nan(array(hi%IscB:hi%IecB,hi%JscB:hi%JecB))) &
461  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
462 ! if (is_NaN(array)) &
463 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
464  endif
465 
466  scaling = 1.0 ; if (present(scale)) scaling = scale
467  iounit = error_unit; if(present(logunit)) iounit = logunit
468  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
469  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
470 
471  if (calculatestatistics) then
472  if (present(scale)) then
473  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
474  lbound(array,2):ubound(array,2)) )
475  rescaled_array(:,:) = 0.0
476  is = hi%isc ; if (sym_stats) is = hi%isc-1
477  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
478  do j=js,hi%JecB ; do i=is,hi%IecB
479  rescaled_array(i,j) = scale*array(i,j)
480  enddo ; enddo
481  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
482  deallocate(rescaled_array)
483  else
484  call substats(hi, array, sym_stats, amean, amin, amax)
485  endif
486  if (is_root_pe()) &
487  call chk_sum_msg("B-point:", amean, amin, amax, mesg, iounit)
488  endif
489 
490  if (.not.writechksums) return
491 
492  hshift = default_shift
493  if (present(haloshift)) hshift = haloshift
494  if (hshift<0) hshift = hi%ied-hi%iec
495 
496  if ( hi%iscB-hshift<hi%isdB .or. hi%iecB+hshift>hi%iedB .or. &
497  hi%jscB-hshift<hi%jsdB .or. hi%jecB+hshift>hi%jedB ) then
498  write(0,*) 'chksum_B_2d: haloshift =',hshift
499  write(0,*) 'chksum_B_2d: isd,isc,iec,ied=',hi%isdB,hi%iscB,hi%iecB,hi%iedB
500  write(0,*) 'chksum_B_2d: jsd,jsc,jec,jed=',hi%jsdB,hi%jscB,hi%jecB,hi%jedB
501  call chksum_error(fatal,'Error in chksum_B_2d '//trim(mesg))
502  endif
503 
504  bc0 = subchk(array, hi, 0, 0, scaling)
505 
506  sym = .false. ; if (present(symmetric)) sym = symmetric
507 
508  if ((hshift==0) .and. .not.sym) then
509  if (is_root_pe()) call chk_sum_msg("B-point:", bc0, mesg, iounit)
510  return
511  endif
512 
513  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
514 
515  if (do_corners) then
516  if (sym) then
517  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
518  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
519  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
520  else
521  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
522  bcse = subchk(array, hi, hshift, -hshift, scaling)
523  bcnw = subchk(array, hi, -hshift, hshift, scaling)
524  endif
525  bcne = subchk(array, hi, hshift, hshift, scaling)
526 
527  if (is_root_pe()) &
528  call chk_sum_msg("B-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
529  else
530  bcs = subchk(array, hi, 0, -hshift, scaling)
531  bce = subchk(array, hi, hshift, 0, scaling)
532  bcw = subchk(array, hi, -hshift, 0, scaling)
533  bcn = subchk(array, hi, 0, hshift, scaling)
534 
535  if (is_root_pe()) &
536  call chk_sum_msg_nsew("B-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
537  endif
538 
539  contains
540 
541  integer function subchk(array, HI, di, dj, scale)
542  type(hor_index_type), intent(in) :: HI !< A horizontal index type
543  real, dimension(HI%IsdB:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
544  integer, intent(in) :: di !< i- direction array shift for this checksum
545  integer, intent(in) :: dj !< j- direction array shift for this checksum
546  real, intent(in) :: scale !< A scaling factor for this array.
547  integer :: i, j, bc
548  subchk = 0
549  ! This line deliberately uses the h-point computational domain.
550  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
551  bc = bitcount(abs(scale*array(i,j)))
552  subchk = subchk + bc
553  enddo ; enddo
554  call sum_across_pes(subchk)
555  subchk=mod(subchk, bc_modulus)
556  end function subchk
557 
558  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
559  type(hor_index_type), intent(in) :: HI !< A horizontal index type
560  real, dimension(HI%IsdB:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
561  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
562  !! full symmetric computational domain.
563  real, intent(out) :: aMean, aMin, aMax
564 
565  integer :: i, j, n, IsB, JsB
566 
567  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
568  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
569 
570  amin = array(hi%isc,hi%jsc) ; amax = amin
571  do j=jsb,hi%JecB ; do i=isb,hi%IecB
572  amin = min(amin, array(i,j))
573  amax = max(amax, array(i,j))
574  enddo ; enddo
575  ! This line deliberately uses the h-point computational domain.
576  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
577  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
578  call sum_across_pes(n)
579  call min_across_pes(amin)
580  call max_across_pes(amax)
581  amean = amean / real(n)
582  end subroutine substats
583 

◆ chksum_b_3d()

subroutine mom_checksums::bchksum::chksum_b_3d ( real, dimension(hi%isdb:,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 corner 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 1085 of file MOM_checksums.F90.

1085  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1086  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1087  character(len=*), intent(in) :: mesg !< An identifying message
1088  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1089  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1090  !! symmetric computational domain.
1091  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1092  real, optional, intent(in) :: scale !< A scaling factor for this array.
1093  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1094 
1095  real, allocatable, dimension(:,:,:) :: rescaled_array
1096  real :: scaling
1097  integer :: iounit !< Log IO unit
1098  integer :: i, j, k, Is, Js
1099  real :: aMean, aMin, aMax
1100  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1101  integer :: bcN, bcS, bcE, bcW
1102  logical :: do_corners, sym, sym_stats
1103 
1104  if (checkfornans) then
1105  if (is_nan(array(hi%IscB:hi%IecB,hi%JscB:hi%JecB,:))) &
1106  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1107 ! if (is_NaN(array)) &
1108 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1109  endif
1110 
1111  scaling = 1.0 ; if (present(scale)) scaling = scale
1112  iounit = error_unit; if(present(logunit)) iounit = logunit
1113  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1114  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1115 
1116  if (calculatestatistics) then
1117  if (present(scale)) then
1118  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1119  lbound(array,2):ubound(array,2), &
1120  lbound(array,3):ubound(array,3)) )
1121  rescaled_array(:,:,:) = 0.0
1122  is = hi%isc ; if (sym_stats) is = hi%isc-1
1123  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
1124  do k=1,size(array,3) ; do j=js,hi%JecB ; do i=is,hi%IecB
1125  rescaled_array(i,j,k) = scale*array(i,j,k)
1126  enddo ; enddo ; enddo
1127  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1128  deallocate(rescaled_array)
1129  else
1130  call substats(hi, array, sym_stats, amean, amin, amax)
1131  endif
1132 
1133  if (is_root_pe()) &
1134  call chk_sum_msg("B-point:", amean, amin, amax, mesg, iounit)
1135  endif
1136 
1137  if (.not.writechksums) return
1138 
1139  hshift = default_shift
1140  if (present(haloshift)) hshift = haloshift
1141  if (hshift<0) hshift = hi%ied-hi%iec
1142 
1143  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1144  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1145  write(0,*) 'chksum_B_3d: haloshift =',hshift
1146  write(0,*) 'chksum_B_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1147  write(0,*) 'chksum_B_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1148  call chksum_error(fatal,'Error in chksum_B_3d '//trim(mesg))
1149  endif
1150 
1151  bc0 = subchk(array, hi, 0, 0, scaling)
1152 
1153  sym = .false. ; if (present(symmetric)) sym = symmetric
1154 
1155  if ((hshift==0) .and. .not.sym) then
1156  if (is_root_pe()) call chk_sum_msg("B-point:", bc0, mesg, iounit)
1157  return
1158  endif
1159 
1160  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1161 
1162  if (do_corners) then
1163  if (sym) then
1164  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
1165  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1166  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
1167  else
1168  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
1169  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1170  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
1171  endif
1172  bcne = subchk(array, hi, hshift, hshift, scaling)
1173 
1174  if (is_root_pe()) &
1175  call chk_sum_msg("B-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1176  else
1177  if (sym) then
1178  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1179  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1180  else
1181  bcs = subchk(array, hi, 0, -hshift, scaling)
1182  bcw = subchk(array, hi, -hshift, 0, scaling)
1183  endif
1184  bce = subchk(array, hi, hshift, 0, scaling)
1185  bcn = subchk(array, hi, 0, hshift, scaling)
1186 
1187  if (is_root_pe()) &
1188  call chk_sum_msg_nsew("B-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1189  endif
1190 
1191  contains
1192 
1193  integer function subchk(array, HI, di, dj, scale)
1194  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1195  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1196  integer, intent(in) :: di !< i- direction array shift for this checksum
1197  integer, intent(in) :: dj !< j- direction array shift for this checksum
1198  real, intent(in) :: scale !< A scaling factor for this array.
1199  integer :: i, j, k, bc
1200  subchk = 0
1201  ! This line deliberately uses the h-point computational domain.
1202  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1203  bc = bitcount(abs(scale*array(i,j,k)))
1204  subchk = subchk + bc
1205  enddo ; enddo ; enddo
1206  call sum_across_pes(subchk)
1207  subchk=mod(subchk, bc_modulus)
1208  end function subchk
1209 
1210  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1211  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1212  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1213  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1214  !! full symmetric computational domain.
1215  real, intent(out) :: aMean, aMin, aMax
1216 
1217  integer :: i, j, k, n, IsB, JsB
1218 
1219  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
1220  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
1221 
1222  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1223  do k=lbound(array,3),ubound(array,3) ; do j=jsb,hi%JecB ; do i=isb,hi%IecB
1224  amin = min(amin, array(i,j,k))
1225  amax = max(amax, array(i,j,k))
1226  enddo ; enddo ; enddo
1227  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1228  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1229  call sum_across_pes(n)
1230  call min_across_pes(amin)
1231  call max_across_pes(amax)
1232  amean = amean / real(n)
1233  end subroutine substats
1234 

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