MOM6
mom_checksums Module Reference

Data Types

interface  bchksum
 
interface  bchksum_pair
 
interface  chk_sum_msg
 
interface  chksum
 
interface  hchksum
 
interface  hchksum_pair
 
interface  is_nan
 
interface  qchksum
 
interface  uchksum
 
interface  uvchksum
 
interface  vchksum
 

Functions/Subroutines

subroutine chksum_pair_h_2d (mesg, arrayA, arrayB, HI, haloshift, omit_corners, scale)
 
subroutine chksum_pair_h_3d (mesg, arrayA, arrayB, HI, haloshift, omit_corners, scale)
 
subroutine chksum_h_2d (array, mesg, HI, haloshift, omit_corners, scale)
 chksum_h_2d performs checksums on a 2d array staggered at tracer points. More...
 
subroutine chksum_pair_b_2d (mesg, arrayA, arrayB, HI, haloshift, symmetric, omit_corners, scale)
 
subroutine chksum_pair_b_3d (mesg, arrayA, arrayB, HI, haloshift, symmetric, omit_corners, scale)
 
subroutine chksum_b_2d (array, mesg, HI, haloshift, symmetric, omit_corners, scale)
 chksum_B_2d performs checksums on a 2d array staggered at corner points. More...
 
subroutine chksum_uv_2d (mesg, arrayU, arrayV, HI, haloshift, symmetric, omit_corners, scale)
 
subroutine chksum_uv_3d (mesg, arrayU, arrayV, HI, haloshift, symmetric, omit_corners, scale)
 
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_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_h_3d (array, mesg, HI, haloshift, omit_corners, scale)
 chksum_h_3d performs checksums on a 3d array staggered at tracer points. More...
 
subroutine chksum_b_3d (array, mesg, HI, haloshift, symmetric, omit_corners, scale)
 chksum_B_3d performs checksums on a 3d array staggered at corner 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...
 
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...
 
subroutine chksum1d (array, mesg, start_i, end_i, compare_PEs)
 chksum1d does a checksum of a 1-dimensional array. More...
 
subroutine chksum2d (array, mesg)
 chksum2d does a checksum of all data in a 2-d array. More...
 
subroutine chksum3d (array, mesg)
 chksum3d does a checksum of all data in a 2-d array. More...
 
logical function is_nan_0d (x)
 This function returns .true. if x is a NaN, and .false. otherwise. More...
 
logical function is_nan_1d (x, skip_mpp)
 This function returns .true. if any element of x is a NaN, and .false. otherwise. More...
 
logical function is_nan_2d (x)
 This function returns .true. if any element of x is a NaN, and .false. otherwise. More...
 
logical function is_nan_3d (x)
 This function returns .true. if any element of x is a NaN, and .false. otherwise. More...
 
subroutine chk_sum_msg1 (fmsg, bc0, mesg)
 
subroutine chk_sum_msg5 (fmsg, bc0, bcSW, bcSE, bcNW, bcNE, mesg)
 
subroutine chk_sum_msg_nsew (fmsg, bc0, bcN, bcS, bcE, bcW, mesg)
 
subroutine chk_sum_msg_s (fmsg, bc0, bcS, mesg)
 
subroutine chk_sum_msg_w (fmsg, bc0, bcW, mesg)
 
subroutine chk_sum_msg2 (fmsg, bc0, bcSW, mesg)
 
subroutine chk_sum_msg3 (fmsg, aMean, aMin, aMax, mesg)
 
subroutine, public mom_checksums_init (param_file)
 MOM_checksums_init initializes the MOM_checksums module. As it happens, the only thing that it does is to log the version of this module. More...
 
subroutine chksum_error (signal, message)
 

Variables

integer, parameter default_shift =0
 
logical calculatestatistics =.true.
 
logical writechksums =.true.
 
logical checkfornans =.true.
 

Function/Subroutine Documentation

◆ chk_sum_msg1()

subroutine mom_checksums::chk_sum_msg1 ( character(len=*), intent(in)  fmsg,
integer, intent(in)  bc0,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1403 of file MOM_checksums.F90.

1403  character(len=*), intent(in) :: fmsg, mesg
1404  integer, intent(in) :: bc0
1405  if (is_root_pe()) write(0,'(A,1(A,I10,X),A)') fmsg," c=",bc0,trim(mesg)

◆ chk_sum_msg2()

subroutine mom_checksums::chk_sum_msg2 ( character(len=*), intent(in)  fmsg,
integer, intent(in)  bc0,
integer, intent(in)  bcSW,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1447 of file MOM_checksums.F90.

1447  character(len=*), intent(in) :: fmsg, mesg
1448  integer, intent(in) :: bc0,bcsw
1449  if (is_root_pe()) write(0,'(A,2(A,I9,1X),A)') &
1450  fmsg," c=",bc0,"s/w=",bcsw,trim(mesg)

◆ chk_sum_msg3()

subroutine mom_checksums::chk_sum_msg3 ( character(len=*), intent(in)  fmsg,
real, intent(in)  aMean,
real, intent(in)  aMin,
real, intent(in)  aMax,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1456 of file MOM_checksums.F90.

1456  character(len=*), intent(in) :: fmsg, mesg
1457  real, intent(in) :: amean,amin,amax
1458  if (is_root_pe()) write(0,'(A,3(A,ES25.16,1X),A)') &
1459  fmsg," mean=",amean,"min=",amin,"max=",amax,trim(mesg)

◆ chk_sum_msg5()

subroutine mom_checksums::chk_sum_msg5 ( character(len=*), intent(in)  fmsg,
integer, intent(in)  bc0,
integer, intent(in)  bcSW,
integer, intent(in)  bcSE,
integer, intent(in)  bcNW,
integer, intent(in)  bcNE,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1411 of file MOM_checksums.F90.

1411  character(len=*), intent(in) :: fmsg, mesg
1412  integer, intent(in) :: bc0,bcsw,bcse,bcnw,bcne
1413  if (is_root_pe()) write(0,'(A,5(A,I10,1X),A)') &
1414  fmsg," c=",bc0,"sw=",bcsw,"se=",bcse,"nw=",bcnw,"ne=",bcne,trim(mesg)

◆ chk_sum_msg_nsew()

subroutine mom_checksums::chk_sum_msg_nsew ( character(len=*), intent(in)  fmsg,
integer, intent(in)  bc0,
integer, intent(in)  bcN,
integer, intent(in)  bcS,
integer, intent(in)  bcE,
integer, intent(in)  bcW,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1420 of file MOM_checksums.F90.

Referenced by chksum_b_2d(), chksum_b_3d(), chksum_h_2d(), chksum_h_3d(), chksum_u_2d(), chksum_u_3d(), chksum_v_2d(), and chksum_v_3d().

1420  character(len=*), intent(in) :: fmsg, mesg
1421  integer, intent(in) :: bc0, bcn, bcs, bce, bcw
1422  if (is_root_pe()) write(0,'(A,5(A,I10,1X),A)') &
1423  fmsg," c=",bc0,"N=",bcn,"S=",bcs,"E=",bce,"W=",bcw,trim(mesg)
Here is the caller graph for this function:

◆ chk_sum_msg_s()

subroutine mom_checksums::chk_sum_msg_s ( character(len=*), intent(in)  fmsg,
integer, intent(in)  bc0,
integer, intent(in)  bcS,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1429 of file MOM_checksums.F90.

Referenced by chksum_v_2d(), and chksum_v_3d().

1429  character(len=*), intent(in) :: fmsg, mesg
1430  integer, intent(in) :: bc0, bcs
1431  if (is_root_pe()) write(0,'(A,2(A,I10,1X),A)') &
1432  fmsg," c=",bc0,"S=",bcs,trim(mesg)
Here is the caller graph for this function:

◆ chk_sum_msg_w()

subroutine mom_checksums::chk_sum_msg_w ( character(len=*), intent(in)  fmsg,
integer, intent(in)  bc0,
integer, intent(in)  bcW,
character(len=*), intent(in)  mesg 
)
private

Definition at line 1438 of file MOM_checksums.F90.

Referenced by chksum_u_2d(), and chksum_u_3d().

1438  character(len=*), intent(in) :: fmsg, mesg
1439  integer, intent(in) :: bc0, bcw
1440  if (is_root_pe()) write(0,'(A,2(A,I10,1X),A)') &
1441  fmsg," c=",bc0,"W=",bcw,trim(mesg)
Here is the caller graph for this function:

◆ chksum1d()

subroutine mom_checksums::chksum1d ( real, dimension(:), intent(in)  array,
character(len=*), intent(in)  mesg,
integer, intent(in), optional  start_i,
integer, intent(in), optional  end_i,
logical, intent(in), optional  compare_PEs 
)
private

chksum1d does a checksum of a 1-dimensional array.

Parameters
[in]arrayThe array to be summed (index starts at 1).
[in]mesgAn identifying message.
[in]start_iThe starting index for the sum (default 1)
[in]end_iThe ending index for the sum (default all)
[in]compare_pesIf true, compare across PEs instead of summing and list the root_PE value (default true)

Definition at line 1201 of file MOM_checksums.F90.

1201  real, dimension(:), intent(in) :: array !< The array to be summed (index starts at 1).
1202  character(len=*), intent(in) :: mesg !< An identifying message.
1203  integer, optional, intent(in) :: start_i !< The starting index for the sum (default 1)
1204  integer, optional, intent(in) :: end_i !< The ending index for the sum (default all)
1205  logical, optional, intent(in) :: compare_pes !< If true, compare across PEs instead of summing
1206  !! and list the root_PE value (default true)
1207 
1208  integer :: is, ie, i, bc, sum1, sum_bc
1209  integer :: bitcount
1210  real :: sum
1211  real, allocatable :: sum_here(:)
1212  logical :: compare
1213  integer :: pe_num ! pe number of the data
1214  integer :: npes ! Total number of processsors
1215 
1216  is = lbound(array,1) ; ie = ubound(array,1)
1217  if (present(start_i)) is = start_i
1218  if (present(end_i)) ie = end_i
1219  compare = .true. ; if (present(compare_pes)) compare = compare_pes
1220 
1221  sum = 0.0 ; sum_bc = 0
1222  do i=is,ie
1223  sum = sum + array(i)
1224  bc = bitcount(abs(array(i)))
1225  sum_bc = sum_bc + bc
1226  enddo
1227 
1228  pe_num = pe_here() + 1 - root_pe() ; npes = num_pes()
1229  allocate(sum_here(npes)) ; sum_here(:) = 0.0 ; sum_here(pe_num) = sum
1230  call sum_across_pes(sum_here,npes)
1231 
1232  sum1 = sum_bc
1233  call sum_across_pes(sum1)
1234 
1235  if (.not.compare) then
1236  sum = 0.0
1237  do i=1,npes ; sum = sum + sum_here(i) ; enddo
1238  sum_bc = sum1
1239  elseif (is_root_pe()) then
1240  if (sum1 /= npes*sum_bc) &
1241  write(0, '(A40," bitcounts do not match across PEs: ",I12,1X,I12)') &
1242  mesg, sum1, npes*sum_bc
1243  do i=1,npes ; if (sum /= sum_here(i)) then
1244  write(0, '(A40," PE ",i4," sum mismatches root_PE: ",3(ES22.13,1X))') &
1245  mesg, i, sum_here(i), sum, sum_here(i)-sum
1246  endif ; enddo
1247  endif
1248  deallocate(sum_here)
1249 
1250  if (is_root_pe()) &
1251  write(0,'(A50,1X,ES25.16,1X,I12)') mesg, sum, sum_bc
1252 
int bitcount(double *x)
Definition: bitcount.c:22

◆ chksum2d()

subroutine mom_checksums::chksum2d ( real, dimension(:,:)  array,
character(len=*)  mesg 
)
private

chksum2d does a checksum of all data in a 2-d array.

Definition at line 1261 of file MOM_checksums.F90.

1261 
1262  real, dimension(:,:) :: array
1263  character(len=*) :: mesg
1264 
1265  integer :: bitcount
1266  integer :: xs,xe,ys,ye,i,j,sum1,bc
1267  real :: sum
1268 
1269  xs = lbound(array,1) ; xe = ubound(array,1)
1270  ys = lbound(array,2) ; ye = ubound(array,2)
1271 
1272  sum = 0.0 ; sum1 = 0
1273  do i=xs,xe ; do j=ys,ye
1274  bc = bitcount(abs(array(i,j)))
1275  sum1 = sum1 + bc
1276  enddo ; enddo
1277  call sum_across_pes(sum1)
1278 
1279  sum = reproducing_sum(array(:,:))
1280 
1281  if (is_root_pe()) &
1282  write(0,'(A50,1X,ES25.16,1X,I12)') mesg, sum, sum1
1283 ! write(0,'(A40,1X,Z16.16,1X,Z16.16,1X,ES25.16,1X,I12)') &
1284 ! mesg, sum, sum1, sum, sum1
1285 
int bitcount(double *x)
Definition: bitcount.c:22

◆ chksum3d()

subroutine mom_checksums::chksum3d ( real, dimension(:,:,:)  array,
character(len=*)  mesg 
)
private

chksum3d does a checksum of all data in a 2-d array.

Definition at line 1290 of file MOM_checksums.F90.

1290 
1291  real, dimension(:,:,:) :: array
1292  character(len=*) :: mesg
1293 
1294  integer :: bitcount
1295  integer :: xs,xe,ys,ye,zs,ze,i,j,k, bc,sum1
1296  real :: sum
1297 
1298  xs = lbound(array,1) ; xe = ubound(array,1)
1299  ys = lbound(array,2) ; ye = ubound(array,2)
1300  zs = lbound(array,3) ; ze = ubound(array,3)
1301 
1302  sum = 0.0 ; sum1 = 0
1303  do i=xs,xe ; do j=ys,ye ; do k=zs,ze
1304  bc = bitcount(abs(array(i,j,k)))
1305  sum1 = sum1 + bc
1306  enddo ; enddo ; enddo
1307 
1308  call sum_across_pes(sum1)
1309  sum = reproducing_sum(array(:,:,:))
1310 
1311  if (is_root_pe()) &
1312  write(0,'(A50,1X,ES25.16,1X,I12)') mesg, sum, sum1
1313 ! write(0,'(A40,1X,Z16.16,1X,Z16.16,1X,ES25.16,1X,I12)') &
1314 ! mesg, sum, sum1, sum, sum1
1315 
int bitcount(double *x)
Definition: bitcount.c:22

◆ chksum_b_2d()

subroutine mom_checksums::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 
)
private

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

Definition at line 285 of file MOM_checksums.F90.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_pair_b_2d().

285  type(hor_index_type), intent(in) :: hi !< A horizontal index type
286  real, dimension(HI%IsdB:,HI%JsdB:), &
287  intent(in) :: array !< The array to be checksummed
288  character(len=*), intent(in) :: mesg !< An identifying message
289  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
290  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the
291  !! full symmetric computational domain.
292  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
293  real, optional, intent(in) :: scale !< A scaling factor for this array.
294 
295  real :: scaling
296  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
297  integer :: bcn, bcs, bce, bcw
298  logical :: do_corners, sym, sym_stats
299 
300  if (checkfornans) then
301  if (is_nan(array(hi%IscB:hi%IecB,hi%JscB:hi%JecB))) &
302  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
303 ! if (is_NaN(array)) &
304 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
305  endif
306  scaling = 1.0 ; if (present(scale)) scaling = scale
307 
308  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
309  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
310  if (calculatestatistics) call substats(hi, array, mesg, sym_stats, scaling)
311 
312  if (.not.writechksums) return
313 
314  hshift = default_shift
315  if (present(haloshift)) hshift = haloshift
316  if (hshift<0) hshift = hi%ied-hi%iec
317 
318  if ( hi%iscB-hshift<hi%isdB .or. hi%iecB+hshift>hi%iedB .or. &
319  hi%jscB-hshift<hi%jsdB .or. hi%jecB+hshift>hi%jedB ) then
320  write(0,*) 'chksum_B_2d: haloshift =',hshift
321  write(0,*) 'chksum_B_2d: isd,isc,iec,ied=',hi%isdB,hi%iscB,hi%iecB,hi%iedB
322  write(0,*) 'chksum_B_2d: jsd,jsc,jec,jed=',hi%jsdB,hi%jscB,hi%jecB,hi%jedB
323  call chksum_error(fatal,'Error in chksum_B_2d '//trim(mesg))
324  endif
325 
326  bc0 = subchk(array, hi, 0, 0, scaling)
327 
328  sym = .false. ; if (present(symmetric)) sym = symmetric
329 
330  if ((hshift==0) .and. .not.sym) then
331  if (is_root_pe()) call chk_sum_msg("B-point:",bc0,mesg)
332  return
333  endif
334 
335  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
336 
337  if (do_corners) then
338  if (sym) then
339  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
340  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
341  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
342  else
343  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
344  bcse = subchk(array, hi, hshift, -hshift, scaling)
345  bcnw = subchk(array, hi, -hshift, hshift, scaling)
346  endif
347  bcne = subchk(array, hi, hshift, hshift, scaling)
348 
349  if (is_root_pe()) call chk_sum_msg("B-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
350  else
351  bcs = subchk(array, hi, 0, -hshift, scaling)
352  bce = subchk(array, hi, hshift, 0, scaling)
353  bcw = subchk(array, hi, -hshift, 0, scaling)
354  bcn = subchk(array, hi, 0, hshift, scaling)
355 
356  if (is_root_pe()) call chk_sum_msg_nsew("B-point:",bc0,bcn,bcs,bce,bcw,mesg)
357  endif
358 
359  contains
360 
361  integer function subchk(array, HI, di, dj, scale)
362  type(hor_index_type), intent(in) :: hi
363  real, dimension(HI%IsdB:,HI%JsdB:), intent(in) :: array
364  integer, intent(in) :: di, dj
365  real, intent(in) :: scale
366  integer :: bitcount, i, j, bc
367  subchk = 0
368  ! This line deliberately uses the h-point computational domain.
369  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
370  bc = bitcount(abs(scale*array(i,j)))
371  subchk = subchk + bc
372  enddo; enddo
373  call sum_across_pes(subchk)
374  subchk=mod(subchk,1000000000)
375  end function subchk
376 
377  subroutine substats(HI, array, mesg, sym_stats, scale)
378  type(hor_index_type), intent(in) :: hi
379  real, dimension(HI%IsdB:,HI%JsdB:), intent(in) :: array
380  character(len=*), intent(in) :: mesg
381  logical, intent(in) :: sym_stats
382  real, intent(in) :: scale
383  integer :: i, j, n, isb, jsb
384  real :: amean, amin, amax
385 
386  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
387  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
388 
389  amin = array(hi%isc,hi%jsc) ; amax = amin
390  do j=jsb,hi%JecB ; do i=isb,hi%IecB
391  amin = min(amin, array(i,j))
392  amax = max(amax, array(i,j))
393  enddo ; enddo
394  ! This line deliberately uses the h-point computational domain.
395  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
396  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
397  call sum_across_pes(n)
398  call min_across_pes(amin)
399  call max_across_pes(amax)
400  amean = amean / real(n)
401  if (is_root_pe()) call chk_sum_msg("B-point:",amean*scale,amin*scale,amax*scale,mesg)
402  end subroutine substats
403 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_b_3d()

subroutine mom_checksums::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 
)
private

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

Definition at line 814 of file MOM_checksums.F90.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_pair_b_3d().

814  type(hor_index_type), intent(in) :: hi !< A horizontal index type
815  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
816  character(len=*), intent(in) :: mesg !< An identifying message
817  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
818  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
819  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
820  real, optional, intent(in) :: scale !< A scaling factor for this array.
821 
822  real :: scaling
823  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
824  integer :: bcn, bcs, bce, bcw
825  logical :: do_corners, sym, sym_stats
826 
827  if (checkfornans) then
828  if (is_nan(array(hi%IscB:hi%IecB,hi%JscB:hi%JecB,:))) &
829  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
830 ! if (is_NaN(array)) &
831 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
832  endif
833  scaling = 1.0 ; if (present(scale)) scaling = scale
834 
835  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
836  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
837  if (calculatestatistics) call substats(hi, array, mesg, sym_stats, scaling)
838 
839  if (.not.writechksums) return
840 
841  hshift = default_shift
842  if (present(haloshift)) hshift = haloshift
843  if (hshift<0) hshift = hi%ied-hi%iec
844 
845  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
846  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
847  write(0,*) 'chksum_B_3d: haloshift =',hshift
848  write(0,*) 'chksum_B_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
849  write(0,*) 'chksum_B_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
850  call chksum_error(fatal,'Error in chksum_B_3d '//trim(mesg))
851  endif
852 
853  bc0 = subchk(array, hi, 0, 0, scaling)
854 
855  sym = .false. ; if (present(symmetric)) sym = symmetric
856 
857  if ((hshift==0) .and. .not.sym) then
858  if (is_root_pe()) call chk_sum_msg("B-point:",bc0,mesg)
859  return
860  endif
861 
862  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
863 
864  if (do_corners) then
865  if (sym) then
866  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
867  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
868  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
869  else
870  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
871  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
872  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
873  endif
874  bcne = subchk(array, hi, hshift, hshift, scaling)
875 
876  if (is_root_pe()) call chk_sum_msg("B-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
877  else
878  if (sym) then
879  bcs = subchk(array, hi, 0, -hshift-1, scaling)
880  bcw = subchk(array, hi, -hshift-1, 0, scaling)
881  else
882  bcs = subchk(array, hi, 0, -hshift, scaling)
883  bcw = subchk(array, hi, -hshift, 0, scaling)
884  endif
885  bce = subchk(array, hi, hshift, 0, scaling)
886  bcn = subchk(array, hi, 0, hshift, scaling)
887 
888  if (is_root_pe()) call chk_sum_msg_nsew("B-point:",bc0,bcn,bcs,bce,bcw,mesg)
889  endif
890 
891  contains
892 
893  integer function subchk(array, HI, di, dj, scale)
894  type(hor_index_type), intent(in) :: hi
895  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array
896  integer, intent(in) :: di, dj
897  real, intent(in) :: scale
898  integer :: bitcount, i, j, k, bc
899  subchk = 0
900  ! This line deliberately uses the h-point computational domain.
901  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
902  bc = bitcount(abs(scale*array(i,j,k)))
903  subchk = subchk + bc
904  enddo ; enddo ; enddo
905  call sum_across_pes(subchk)
906  subchk=mod(subchk,1000000000)
907  end function subchk
908 
909  subroutine substats(HI, array, mesg, sym_stats, scale)
910  type(hor_index_type), intent(in) :: hi
911  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array
912  character(len=*), intent(in) :: mesg
913  logical, intent(in) :: sym_stats
914  real, intent(in) :: scale
915  integer :: i, j, k, n, isb, jsb
916  real :: amean, amin, amax
917 
918  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
919  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
920 
921  amin = array(hi%isc,hi%jsc,1) ; amax = amin
922  do k=lbound(array,3),ubound(array,3) ; do j=jsb,hi%JecB ; do i=isb,hi%IecB
923  amin = min(amin, array(i,j,k))
924  amax = max(amax, array(i,j,k))
925  enddo ; enddo ; enddo
926  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
927  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
928  call sum_across_pes(n)
929  call min_across_pes(amin)
930  call max_across_pes(amax)
931  amean = amean / real(n)
932  if (is_root_pe()) call chk_sum_msg("B-point:",amean*scale,amin*scale,amax*scale,mesg)
933  end subroutine substats
934 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_error()

subroutine mom_checksums::chksum_error ( integer, intent(in)  signal,
character(len=*), intent(in)  message 
)
private

Definition at line 1479 of file MOM_checksums.F90.

Referenced by chksum_b_2d(), chksum_b_3d(), chksum_h_2d(), chksum_h_3d(), chksum_u_2d(), chksum_u_3d(), chksum_v_2d(), and chksum_v_3d().

1479  ! Wrapper for MOM_error to help place specific break points in
1480  ! debuggers
1481  integer, intent(in) :: signal
1482  character(len=*), intent(in) :: message
1483  call mom_error(signal, message)
Here is the caller graph for this function:

◆ chksum_h_2d()

subroutine mom_checksums::chksum_h_2d ( real, dimension(hi%isd:,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  omit_corners,
real, intent(in), optional  scale 
)
private

chksum_h_2d performs checksums on a 2d array staggered at tracer 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]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.

Definition at line 129 of file MOM_checksums.F90.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chksum_error(), default_shift, mom_error_handler::is_root_pe(), subchk(), substats(), and writechksums.

Referenced by chksum_pair_h_2d().

129  type(hor_index_type), intent(in) :: hi !< A horizontal index type
130  real, dimension(HI%isd:,HI%jsd:), intent(in) :: array !< The array to be checksummed
131  character(len=*), intent(in) :: mesg !< An identifying message
132  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
133  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
134  real, optional, intent(in) :: scale !< A scaling factor for this array.
135 
136  real :: scaling
137  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
138  integer :: bcn, bcs, bce, bcw
139  logical :: do_corners
140 
141  if (checkfornans) then
142  if (is_nan(array(hi%isc:hi%iec,hi%jsc:hi%jec))) &
143  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
144 ! if (is_NaN(array)) &
145 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
146  endif
147  scaling = 1.0 ; if (present(scale)) scaling = scale
148 
149  if (calculatestatistics) call substats(hi, array, mesg, scaling)
150 
151  if (.not.writechksums) return
152 
153  hshift = default_shift
154  if (present(haloshift)) hshift = haloshift
155  if (hshift<0) hshift = hi%ied-hi%iec
156 
157  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
158  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
159  write(0,*) 'chksum_h_2d: haloshift =',hshift
160  write(0,*) 'chksum_h_2d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
161  write(0,*) 'chksum_h_2d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
162  call chksum_error(fatal,'Error in chksum_h_2d '//trim(mesg))
163  endif
164 
165  bc0 = subchk(array, hi, 0, 0, scaling)
166 
167  if (hshift==0) then
168  if (is_root_pe()) call chk_sum_msg("h-point:",bc0,mesg)
169  return
170  endif
171 
172  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
173 
174  if (do_corners) then
175  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
176  bcse = subchk(array, hi, hshift, -hshift, scaling)
177  bcnw = subchk(array, hi, -hshift, hshift, scaling)
178  bcne = subchk(array, hi, hshift, hshift, scaling)
179 
180  if (is_root_pe()) call chk_sum_msg("h-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
181  else
182  bcs = subchk(array, hi, 0, -hshift, scaling)
183  bce = subchk(array, hi, hshift, 0, scaling)
184  bcw = subchk(array, hi, -hshift, 0, scaling)
185  bcn = subchk(array, hi, 0, hshift, scaling)
186 
187  if (is_root_pe()) call chk_sum_msg_nsew("h-point:",bc0,bcn,bcs,bce,bcw,mesg)
188  endif
189 
190  contains
191 
192  integer function subchk(array, HI, di, dj, scale)
193  type(hor_index_type), intent(in) :: hi
194  real, dimension(HI%isd:,HI%jsd:), intent(in) :: array
195  integer, intent(in) :: di, dj
196  real, intent(in) :: scale
197  integer :: bitcount, i, j, bc
198  subchk = 0
199  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
200  bc = bitcount(abs(scale*array(i,j)))
201  subchk = subchk + bc
202  enddo; enddo
203  call sum_across_pes(subchk)
204  subchk=mod(subchk,1000000000)
205  end function subchk
206 
207  subroutine substats(HI, array, mesg, scale)
208  type(hor_index_type), intent(in) :: hi
209  real, dimension(HI%isd:,HI%jsd:), intent(in) :: array
210  character(len=*), intent(in) :: mesg
211  real, intent(in) :: scale
212  integer :: i, j, n
213  real :: amean, amin, amax
214 
215  amin = array(hi%isc,hi%jsc)
216  amax = array(hi%isc,hi%jsc)
217  n = 0
218  do j=hi%jsc,hi%jec ; do i=hi%isc,hi%iec
219  amin = min(amin, array(i,j))
220  amax = max(amax, array(i,j))
221  n = n + 1
222  enddo ; enddo
223  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
224  call sum_across_pes(n)
225  call min_across_pes(amin)
226  call max_across_pes(amax)
227  amean = amean / real(n)
228  if (is_root_pe()) call chk_sum_msg("h-point:",amean*scale,amin*scale,amax*scale,mesg)
229  end subroutine substats
230 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_h_3d()

subroutine mom_checksums::chksum_h_3d ( real, dimension(hi%isd:,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  omit_corners,
real, intent(in), optional  scale 
)
private

chksum_h_3d performs checksums on a 3d array staggered at tracer 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]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.

Definition at line 706 of file MOM_checksums.F90.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_pair_h_3d().

706  type(hor_index_type), intent(in) :: hi !< A horizontal index type
707  real, dimension(HI%isd:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
708  character(len=*), intent(in) :: mesg !< An identifying message
709  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
710  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
711  real, optional, intent(in) :: scale !< A scaling factor for this array.
712 
713  real :: scaling
714  integer :: bc0, bcsw, bcse, bcnw, bcne, hshift
715  integer :: bcn, bcs, bce, bcw
716  logical :: do_corners
717 
718  if (checkfornans) then
719  if (is_nan(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))) &
720  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
721 ! if (is_NaN(array)) &
722 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
723  endif
724  scaling = 1.0 ; if (present(scale)) scaling = scale
725 
726  if (calculatestatistics) call substats(hi, array, mesg, scaling)
727 
728  if (.not.writechksums) return
729 
730  hshift = default_shift
731  if (present(haloshift)) hshift = haloshift
732  if (hshift<0) hshift = hi%ied-hi%iec
733 
734  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
735  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
736  write(0,*) 'chksum_h_3d: haloshift =',hshift
737  write(0,*) 'chksum_h_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
738  write(0,*) 'chksum_h_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
739  call chksum_error(fatal,'Error in chksum_h_3d '//trim(mesg))
740  endif
741 
742  bc0 = subchk(array, hi, 0, 0, scaling)
743 
744  if (hshift==0) then
745  if (is_root_pe()) call chk_sum_msg("h-point:",bc0,mesg)
746  return
747  endif
748 
749  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
750 
751  if (do_corners) then
752  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
753  bcse = subchk(array, hi, hshift, -hshift, scaling)
754  bcnw = subchk(array, hi, -hshift, hshift, scaling)
755  bcne = subchk(array, hi, hshift, hshift, scaling)
756 
757  if (is_root_pe()) call chk_sum_msg("h-point:",bc0,bcsw,bcse,bcnw,bcne,mesg)
758  else
759  bcs = subchk(array, hi, 0, -hshift, scaling)
760  bce = subchk(array, hi, hshift, 0, scaling)
761  bcw = subchk(array, hi, -hshift, 0, scaling)
762  bcn = subchk(array, hi, 0, hshift, scaling)
763 
764  if (is_root_pe()) call chk_sum_msg_nsew("h-point:",bc0,bcn,bcs,bce,bcw,mesg)
765  endif
766 
767  contains
768 
769  integer function subchk(array, HI, di, dj, scale)
770  type(hor_index_type), intent(in) :: hi
771  real, dimension(HI%isd:,HI%jsd:,:), intent(in) :: array
772  integer, intent(in) :: di, dj
773  real, intent(in) :: scale
774  integer :: bitcount, i, j, k, bc
775  subchk = 0
776  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
777  bc = bitcount(abs(scale*array(i,j,k)))
778  subchk = subchk + bc
779  enddo ; enddo ; enddo
780  call sum_across_pes(subchk)
781  subchk=mod(subchk,1000000000)
782  end function subchk
783 
784  subroutine substats(HI, array, mesg, scale)
785  type(hor_index_type), intent(in) :: hi
786  real, dimension(HI%isd:,HI%jsd:,:), intent(in) :: array
787  character(len=*), intent(in) :: mesg
788  real, intent(in) :: scale
789  integer :: i, j, k, n
790  real :: amean, amin, amax
791 
792  amin = array(hi%isc,hi%jsc,1)
793  amax = array(hi%isc,hi%jsc,1)
794  n = 0
795  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc,hi%jec ; do i=hi%isc,hi%iec
796  amin = min(amin, array(i,j,k))
797  amax = max(amax, array(i,j,k))
798  n = n + 1
799  enddo ; enddo ; enddo
800  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
801  call sum_across_pes(n)
802  call min_across_pes(amin)
803  call max_across_pes(amax)
804  amean = amean / real(n)
805  if (is_root_pe()) call chk_sum_msg("h-point:",amean*scale,amin*scale,amax*scale,mesg)
806  end subroutine substats
807 
subroutine substats(HI, array, mesg, scale)
int bitcount(double *x)
Definition: bitcount.c:22
integer function subchk(array, HI, di, dj, scale)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_pair_b_2d()

subroutine mom_checksums::chksum_pair_b_2d ( character(len=*), intent(in)  mesg,
real, dimension(hi%isd:,hi%jsd:), intent(in)  arrayA,
real, dimension(hi%isd:,hi%jsd:), intent(in)  arrayB,
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
Parameters
[in]mesgIdentifying messages
[in]hiA horizontal index type
[in]arraybThe arrays to be checksummed
[in]symmetricIf true, do the checksums on the full symmetric computational domain.
[in]haloshiftThe width of halos to check (default 0)
[in]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.

Definition at line 236 of file MOM_checksums.F90.

References chksum_b_2d().

236  character(len=*), intent(in) :: mesg !< Identifying messages
237  type(hor_index_type), intent(in) :: hi !< A horizontal index type
238  real, dimension(HI%isd:,HI%jsd:), intent(in) :: arraya, arrayb !< The arrays to be checksummed
239  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
240  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
241  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
242  real, optional, intent(in) :: scale !< A scaling factor for this array.
243 
244  logical :: sym
245 
246  sym = .false. ; if (present(symmetric)) sym = symmetric
247 
248  if (present(haloshift)) then
249  call chksum_b_2d(arraya, 'x '//mesg, hi, haloshift, symmetric=sym, &
250  omit_corners=omit_corners, scale=scale)
251  call chksum_b_2d(arrayb, 'y '//mesg, hi, haloshift, symmetric=sym, &
252  omit_corners=omit_corners, scale=scale)
253  else
254  call chksum_b_2d(arraya, 'x '//mesg, hi, symmetric=sym, scale=scale)
255  call chksum_b_2d(arrayb, 'y '//mesg, hi, symmetric=sym, scale=scale)
256  endif
257 
Here is the call graph for this function:

◆ chksum_pair_b_3d()

subroutine mom_checksums::chksum_pair_b_3d ( character(len=*), intent(in)  mesg,
real, dimension(hi%isdb:,hi%jsdb:, :), intent(in)  arrayA,
real, dimension(hi%isdb:,hi%jsdb:, :), intent(in)  arrayB,
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
Parameters
[in]mesgIdentifying messages
[in]hiA horizontal index type
[in]arraybThe arrays to be checksummed
[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 261 of file MOM_checksums.F90.

References chksum_b_3d().

261  character(len=*), intent(in) :: mesg !< Identifying messages
262  type(hor_index_type), intent(in) :: hi !< A horizontal index type
263  real, dimension(HI%IsdB:,HI%JsdB:, :), intent(in) :: arraya, arrayb !< The arrays to be checksummed
264  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
265  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
266  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
267  real, optional, intent(in) :: scale !< A scaling factor for this array.
268 
269  logical :: sym
270 
271  if (present(haloshift)) then
272  call chksum_b_3d(arraya, 'x '//mesg, hi, haloshift, symmetric, &
273  omit_corners, scale=scale)
274  call chksum_b_3d(arrayb, 'y '//mesg, hi, haloshift, symmetric, &
275  omit_corners, scale=scale)
276  else
277  call chksum_b_3d(arraya, 'x '//mesg, hi, symmetric=symmetric, scale=scale)
278  call chksum_b_3d(arrayb, 'y '//mesg, hi, symmetric=symmetric, scale=scale)
279  endif
280 
Here is the call graph for this function:

◆ chksum_pair_h_2d()

subroutine mom_checksums::chksum_pair_h_2d ( character(len=*), intent(in)  mesg,
real, dimension(hi%isd:,hi%jsd:), intent(in)  arrayA,
real, dimension(hi%isd:,hi%jsd:), intent(in)  arrayB,
type(hor_index_type), intent(in)  HI,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  omit_corners,
real, intent(in), optional  scale 
)
private
Parameters
[in]mesgIdentifying messages
[in]hiA horizontal index type
[in]arraybThe arrays to be checksummed
[in]haloshiftThe width of halos to check (default 0)
[in]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.

Definition at line 92 of file MOM_checksums.F90.

References chksum_h_2d().

92  character(len=*), intent(in) :: mesg !< Identifying messages
93  type(hor_index_type), intent(in) :: hi !< A horizontal index type
94  real, dimension(HI%isd:,HI%jsd:), intent(in) :: arraya, arrayb !< The arrays to be checksummed
95  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
96  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
97  real, optional, intent(in) :: scale !< A scaling factor for this array.
98 
99  if (present(haloshift)) then
100  call chksum_h_2d(arraya, 'x '//mesg, hi, haloshift, omit_corners, scale=scale)
101  call chksum_h_2d(arrayb, 'y '//mesg, hi, haloshift, omit_corners, scale=scale)
102  else
103  call chksum_h_2d(arraya, 'x '//mesg, hi, scale=scale)
104  call chksum_h_2d(arrayb, 'y '//mesg, hi, scale=scale)
105  endif
106 
Here is the call graph for this function:

◆ chksum_pair_h_3d()

subroutine mom_checksums::chksum_pair_h_3d ( character(len=*), intent(in)  mesg,
real, dimension(hi%isd:,hi%jsd:, :), intent(in)  arrayA,
real, dimension(hi%isd:,hi%jsd:, :), intent(in)  arrayB,
type(hor_index_type), intent(in)  HI,
integer, intent(in), optional  haloshift,
logical, intent(in), optional  omit_corners,
real, intent(in), optional  scale 
)
private
Parameters
[in]mesgIdentifying messages
[in]hiA horizontal index type
[in]arraybThe arrays to be checksummed
[in]haloshiftThe width of halos to check (default 0)
[in]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.

Definition at line 110 of file MOM_checksums.F90.

References chksum_h_3d().

110  character(len=*), intent(in) :: mesg !< Identifying messages
111  type(hor_index_type), intent(in) :: hi !< A horizontal index type
112  real, dimension(HI%isd:,HI%jsd:, :), intent(in) :: arraya, arrayb !< The arrays to be checksummed
113  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
114  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
115  real, optional, intent(in) :: scale !< A scaling factor for this array.
116 
117  if (present(haloshift)) then
118  call chksum_h_3d(arraya, 'x '//mesg, hi, haloshift, omit_corners, scale=scale)
119  call chksum_h_3d(arrayb, 'y '//mesg, hi, haloshift, omit_corners, scale=scale)
120  else
121  call chksum_h_3d(arraya, 'x '//mesg, hi, scale=scale)
122  call chksum_h_3d(arrayb, 'y '//mesg, hi, scale=scale)
123  endif
124 
Here is the call graph for this function:

◆ chksum_u_2d()

subroutine mom_checksums::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.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chk_sum_msg_w(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_uv_2d().

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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_u_3d()

subroutine mom_checksums::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.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chk_sum_msg_w(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_uv_3d().

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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_uv_2d()

subroutine mom_checksums::chksum_uv_2d ( character(len=*), intent(in)  mesg,
real, dimension(hi%isdb:,hi%jsd:), intent(in)  arrayU,
real, dimension(hi%isd:,hi%jsdb:), intent(in)  arrayV,
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
Parameters
[in]mesgIdentifying messages
[in]hiA horizontal index type
[in]arrayuThe u-component array to be checksummed
[in]arrayvThe v-component array to be checksummed
[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 these arrays.

Definition at line 409 of file MOM_checksums.F90.

References chksum_u_2d(), and chksum_v_2d().

409  character(len=*), intent(in) :: mesg !< Identifying messages
410  type(hor_index_type), intent(in) :: hi !< A horizontal index type
411  real, dimension(HI%IsdB:,HI%jsd:), intent(in) :: arrayu !< The u-component array to be checksummed
412  real, dimension(HI%isd:,HI%JsdB:), intent(in) :: arrayv !< The v-component array to be checksummed
413  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
414  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
415  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
416  real, optional, intent(in) :: scale !< A scaling factor for these arrays.
417 
418  if (present(haloshift)) then
419  call chksum_u_2d(arrayu, 'u '//mesg, hi, haloshift, symmetric, omit_corners, scale)
420  call chksum_v_2d(arrayv, 'v '//mesg, hi, haloshift, symmetric, omit_corners, scale)
421  else
422  call chksum_u_2d(arrayu, 'u '//mesg, hi, symmetric=symmetric)
423  call chksum_v_2d(arrayv, 'v '//mesg, hi, symmetric=symmetric)
424  endif
425 
Here is the call graph for this function:

◆ chksum_uv_3d()

subroutine mom_checksums::chksum_uv_3d ( character(len=*), intent(in)  mesg,
real, dimension(hi%isdb:,hi%jsd:,:), intent(in)  arrayU,
real, dimension(hi%isd:,hi%jsdb:,:), intent(in)  arrayV,
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
Parameters
[in]mesgIdentifying messages
[in]hiA horizontal index type
[in]arrayuThe u-component array to be checksummed
[in]arrayvThe v-component array to be checksummed
[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 these arrays.

Definition at line 429 of file MOM_checksums.F90.

References chksum_u_3d(), and chksum_v_3d().

429  character(len=*), intent(in) :: mesg !< Identifying messages
430  type(hor_index_type), intent(in) :: hi !< A horizontal index type
431  real, dimension(HI%IsdB:,HI%jsd:,:), intent(in) :: arrayu !< The u-component array to be checksummed
432  real, dimension(HI%isd:,HI%JsdB:,:), intent(in) :: arrayv !< The v-component array to be checksummed
433  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
434  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full symmetric computational domain.
435  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
436  real, optional, intent(in) :: scale !< A scaling factor for these arrays.
437 
438  if (present(haloshift)) then
439  call chksum_u_3d(arrayu, 'u '//mesg, hi, haloshift, symmetric, omit_corners, scale)
440  call chksum_v_3d(arrayv, 'v '//mesg, hi, haloshift, symmetric, omit_corners, scale)
441  else
442  call chksum_u_3d(arrayu, 'u '//mesg, hi, symmetric=symmetric)
443  call chksum_v_3d(arrayv, 'v '//mesg, hi, symmetric=symmetric)
444  endif
445 
Here is the call graph for this function:

◆ chksum_v_2d()

subroutine mom_checksums::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.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chk_sum_msg_s(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_uv_2d().

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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chksum_v_3d()

subroutine mom_checksums::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.

References calculatestatistics, checkfornans, chk_sum_msg_nsew(), chk_sum_msg_s(), chksum_error(), default_shift, subchk(), substats(), and writechksums.

Referenced by chksum_uv_3d().

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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ is_nan_0d()

logical function mom_checksums::is_nan_0d ( real, intent(in)  x)
private

This function returns .true. if x is a NaN, and .false. otherwise.

Parameters
[in]xThe value to be checked for NaNs.

Definition at line 1322 of file MOM_checksums.F90.

Referenced by is_nan_1d(), is_nan_2d(), and is_nan_3d().

1322  real, intent(in) :: x !< The value to be checked for NaNs.
1323  logical :: is_nan_0d
1324 
1325  !is_NaN_0d = (((x < 0.0) .and. (x >= 0.0)) .or. &
1326  ! (.not.(x < 0.0) .and. .not.(x >= 0.0)))
1327  if (((x < 0.0) .and. (x >= 0.0)) .or. &
1328  (.not.(x < 0.0) .and. .not.(x >= 0.0))) then
1329  is_nan_0d = .true.
1330  else
1331  is_nan_0d = .false.
1332  endif
1333 
Here is the caller graph for this function:

◆ is_nan_1d()

logical function mom_checksums::is_nan_1d ( real, dimension(:), intent(in)  x,
logical, optional  skip_mpp 
)
private

This function returns .true. if any element of x is a NaN, and .false. otherwise.

Parameters
[in]xThe array to be checked for NaNs.
skip_mppIf true, only check this array only on the local PE (default false).

Definition at line 1340 of file MOM_checksums.F90.

References is_nan_0d().

1340  real, dimension(:), intent(in) :: x !< The array to be checked for NaNs.
1341  logical :: is_nan_1d
1342  logical, optional :: skip_mpp !< If true, only check this array only on the local PE (default false).
1343 
1344  integer :: i, n
1345  logical :: call_mpp
1346 
1347  n = 0
1348  do i = lbound(x,1), ubound(x,1)
1349  if (is_nan_0d(x(i))) n = n + 1
1350  enddo
1351  call_mpp = .true.
1352  if (present(skip_mpp)) call_mpp = .not.skip_mpp
1353 
1354  if (call_mpp) call sum_across_pes(n)
1355  is_nan_1d = .false.
1356  if (n>0) is_nan_1d = .true.
1357 
Here is the call graph for this function:

◆ is_nan_2d()

logical function mom_checksums::is_nan_2d ( real, dimension(:,:), intent(in)  x)
private

This function returns .true. if any element of x is a NaN, and .false. otherwise.

Parameters
[in]xThe array to be checked for NaNs.

Definition at line 1364 of file MOM_checksums.F90.

References is_nan_0d().

1364  real, dimension(:,:), intent(in) :: x !< The array to be checked for NaNs.
1365  logical :: is_nan_2d
1366 
1367  integer :: i, j, n
1368 
1369  n = 0
1370  do j = lbound(x,2), ubound(x,2) ; do i = lbound(x,1), ubound(x,1)
1371  if (is_nan_0d(x(i,j))) n = n + 1
1372  enddo ; enddo
1373  call sum_across_pes(n)
1374  is_nan_2d = .false.
1375  if (n>0) is_nan_2d = .true.
1376 
Here is the call graph for this function:

◆ is_nan_3d()

logical function mom_checksums::is_nan_3d ( real, dimension(:,:,:), intent(in)  x)
private

This function returns .true. if any element of x is a NaN, and .false. otherwise.

Parameters
[in]xThe array to be checked for NaNs.

Definition at line 1383 of file MOM_checksums.F90.

References is_nan_0d().

1383  real, dimension(:,:,:), intent(in) :: x !< The array to be checked for NaNs.
1384  logical :: is_nan_3d
1385 
1386  integer :: i, j, k, n
1387 
1388  n = 0
1389  do k = lbound(x,3), ubound(x,3)
1390  do j = lbound(x,2), ubound(x,2) ; do i = lbound(x,1), ubound(x,1)
1391  if (is_nan_0d(x(i,j,k))) n = n + 1
1392  enddo ; enddo
1393  enddo
1394  call sum_across_pes(n)
1395  is_nan_3d = .false.
1396  if (n>0) is_nan_3d = .true.
1397 
Here is the call graph for this function:

◆ mom_checksums_init()

subroutine, public mom_checksums::mom_checksums_init ( type(param_file_type), intent(in)  param_file)

MOM_checksums_init initializes the MOM_checksums module. As it happens, the only thing that it does is to log the version of this module.

Parameters
[in]param_fileA structure to parse for run-time parameters

Definition at line 1467 of file MOM_checksums.F90.

Referenced by mom_debugging::mom_debugging_init().

1467  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1468 ! This include declares and sets the variable "version".
1469 #include "version_variable.h"
1470  character(len=40) :: mdl = "MOM_checksums" ! This module's name.
1471 
1472  call log_version(param_file, mdl, version)
1473 
Here is the caller graph for this function:

Variable Documentation

◆ calculatestatistics

logical mom_checksums::calculatestatistics =.true.
private

Definition at line 82 of file MOM_checksums.F90.

Referenced by chksum_b_2d(), chksum_b_3d(), chksum_h_2d(), chksum_h_3d(), chksum_u_2d(), chksum_u_3d(), chksum_v_2d(), and chksum_v_3d().

82 logical :: calculatestatistics=.true. ! If true, report min, max and mean.

◆ checkfornans

logical mom_checksums::checkfornans =.true.
private

Definition at line 84 of file MOM_checksums.F90.

Referenced by chksum_b_2d(), chksum_b_3d(), chksum_h_2d(), chksum_h_3d(), chksum_u_2d(), chksum_u_3d(), chksum_v_2d(), and chksum_v_3d().

84 logical :: checkfornans=.true. ! If true, checks array for NaNs and cause

◆ default_shift

integer, parameter mom_checksums::default_shift =0
private

Definition at line 81 of file MOM_checksums.F90.

Referenced by chksum_b_2d(), chksum_b_3d(), chksum_h_2d(), chksum_h_3d(), chksum_u_2d(), chksum_u_3d(), chksum_v_2d(), and chksum_v_3d().

81 integer, parameter :: default_shift=0

◆ writechksums

logical mom_checksums::writechksums =.true.
private

Definition at line 83 of file MOM_checksums.F90.

Referenced by chksum_b_2d(), chksum_b_3d(), chksum_h_2d(), chksum_h_3d(), chksum_u_2d(), chksum_u_3d(), chksum_v_2d(), and chksum_v_3d().

83 logical :: writechksums=.true. ! If true, report the bitcount checksum