MOM6
mom_checksums::bchksum Interface Reference

Detailed Description

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

Definition at line 54 of file MOM_checksums.F90.

Private functions

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

Detailed Description

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

Definition at line 54 of file MOM_checksums.F90.

Functions and subroutines

◆ chksum_b_2d()

subroutine mom_checksums::bchksum::chksum_b_2d ( real, dimension(hi_m%isdb:,hi_m%jsdb:), intent(in), target  array_m,
character(len=*), intent(in)  mesg,
type(hor_index_type), intent(in), target  HI_m,
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]hi_mA horizontal index type
[in]array_mThe 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 570 of file MOM_checksums.F90.

570  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
571  real, dimension(HI_m%IsdB:,HI_m%JsdB:), &
572  target, intent(in) :: array_m !< The array to be checksummed
573  character(len=*), intent(in) :: mesg !< An identifying message
574  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
575  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the
576  !! full symmetric computational domain.
577  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
578  real, optional, intent(in) :: scale !< A scaling factor for this array.
579  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
580 
581  real, pointer :: array(:,:) ! Field array on the input grid
582  real, allocatable, dimension(:,:) :: rescaled_array
583  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
584  real :: scaling
585  integer :: iounit !< Log IO unit
586  integer :: i, j, Is, Js
587  real :: aMean, aMin, aMax
588  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
589  integer :: bcN, bcS, bcE, bcW
590  logical :: do_corners, sym, sym_stats
591  integer :: turns ! Quarter turns from input to model grid
592 
593  ! Rotate array to the input grid
594  turns = hi_m%turns
595  if (modulo(turns, 4) /= 0) then
596  allocate(hi)
597  call rotate_hor_index(hi_m, -turns, hi)
598  allocate(array(hi%IsdB:hi%IedB, hi%JsdB:hi%JedB))
599  call rotate_array(array_m, -turns, array)
600  else
601  hi => hi_m
602  array => array_m
603  endif
604 
605  if (checkfornans) then
606  if (is_nan(array(hi%IscB:hi%IecB,hi%JscB:hi%JecB))) &
607  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
608 ! if (is_NaN(array)) &
609 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
610  endif
611 
612  scaling = 1.0 ; if (present(scale)) scaling = scale
613  iounit = error_unit; if(present(logunit)) iounit = logunit
614  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
615  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
616 
617  if (calculatestatistics) then
618  if (present(scale)) then
619  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
620  lbound(array,2):ubound(array,2)) )
621  rescaled_array(:,:) = 0.0
622  is = hi%isc ; if (sym_stats) is = hi%isc-1
623  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
624  do j=js,hi%JecB ; do i=is,hi%IecB
625  rescaled_array(i,j) = scale*array(i,j)
626  enddo ; enddo
627  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
628  deallocate(rescaled_array)
629  else
630  call substats(hi, array, sym_stats, amean, amin, amax)
631  endif
632  if (is_root_pe()) &
633  call chk_sum_msg("B-point:", amean, amin, amax, mesg, iounit)
634  endif
635 
636  if (.not.writechksums) return
637 
638  hshift = default_shift
639  if (present(haloshift)) hshift = haloshift
640  if (hshift<0) hshift = hi%ied-hi%iec
641 
642  if ( hi%iscB-hshift<hi%isdB .or. hi%iecB+hshift>hi%iedB .or. &
643  hi%jscB-hshift<hi%jsdB .or. hi%jecB+hshift>hi%jedB ) then
644  write(0,*) 'chksum_B_2d: haloshift =',hshift
645  write(0,*) 'chksum_B_2d: isd,isc,iec,ied=',hi%isdB,hi%iscB,hi%iecB,hi%iedB
646  write(0,*) 'chksum_B_2d: jsd,jsc,jec,jed=',hi%jsdB,hi%jscB,hi%jecB,hi%jedB
647  call chksum_error(fatal,'Error in chksum_B_2d '//trim(mesg))
648  endif
649 
650  bc0 = subchk(array, hi, 0, 0, scaling)
651 
652  sym = .false. ; if (present(symmetric)) sym = symmetric
653 
654  if ((hshift==0) .and. .not.sym) then
655  if (is_root_pe()) call chk_sum_msg("B-point:", bc0, mesg, iounit)
656  return
657  endif
658 
659  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
660 
661  if (do_corners) then
662  if (sym) then
663  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
664  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
665  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
666  else
667  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
668  bcse = subchk(array, hi, hshift, -hshift, scaling)
669  bcnw = subchk(array, hi, -hshift, hshift, scaling)
670  endif
671  bcne = subchk(array, hi, hshift, hshift, scaling)
672 
673  if (is_root_pe()) &
674  call chk_sum_msg("B-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
675  else
676  bcs = subchk(array, hi, 0, -hshift, scaling)
677  bce = subchk(array, hi, hshift, 0, scaling)
678  bcw = subchk(array, hi, -hshift, 0, scaling)
679  bcn = subchk(array, hi, 0, hshift, scaling)
680 
681  if (is_root_pe()) &
682  call chk_sum_msg_nsew("B-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
683  endif
684 
685  contains
686 
687  integer function subchk(array, HI, di, dj, scale)
688  type(hor_index_type), intent(in) :: HI !< A horizontal index type
689  real, dimension(HI%IsdB:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
690  integer, intent(in) :: di !< i- direction array shift for this checksum
691  integer, intent(in) :: dj !< j- direction array shift for this checksum
692  real, intent(in) :: scale !< A scaling factor for this array.
693  integer :: i, j, bc
694  subchk = 0
695  ! This line deliberately uses the h-point computational domain.
696  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
697  bc = bitcount(abs(scale*array(i,j)))
698  subchk = subchk + bc
699  enddo ; enddo
700  call sum_across_pes(subchk)
701  subchk=mod(subchk, bc_modulus)
702  end function subchk
703 
704  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
705  type(hor_index_type), intent(in) :: HI !< A horizontal index type
706  real, dimension(HI%IsdB:,HI%JsdB:), intent(in) :: array !< The array to be checksummed
707  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
708  !! full symmetric computational domain.
709  real, intent(out) :: aMean !< Array mean
710  real, intent(out) :: aMin !< Array minimum
711  real, intent(out) :: aMax !< Array maximum
712 
713  integer :: i, j, n, IsB, JsB
714 
715  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
716  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
717 
718  amin = array(hi%isc,hi%jsc) ; amax = amin
719  do j=jsb,hi%JecB ; do i=isb,hi%IecB
720  amin = min(amin, array(i,j))
721  amax = max(amax, array(i,j))
722  enddo ; enddo
723  ! This line deliberately uses the h-point computational domain.
724  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
725  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc)
726  call sum_across_pes(n)
727  call min_across_pes(amin)
728  call max_across_pes(amax)
729  amean = amean / real(n)
730  end subroutine substats
731 

◆ chksum_b_3d()

subroutine mom_checksums::bchksum::chksum_b_3d ( real, dimension(hi_m%isdb:,hi_m%jsdb:,:), intent(in), target  array_m,
character(len=*), intent(in)  mesg,
type(hor_index_type), intent(in), target  HI_m,
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]hi_mA horizontal index type
[in]array_mThe 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 1354 of file MOM_checksums.F90.

1354  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
1355  real, dimension(HI_m%IsdB:,HI_m%JsdB:,:), target, intent(in) :: array_m !< The array to be checksummed
1356  character(len=*), intent(in) :: mesg !< An identifying message
1357  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1358  logical, optional, intent(in) :: symmetric !< If true, do the checksums on the full
1359  !! symmetric computational domain.
1360  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1361  real, optional, intent(in) :: scale !< A scaling factor for this array.
1362  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1363 
1364  real, pointer :: array(:,:,:) ! Field array on the input grid
1365  real, allocatable, dimension(:,:,:) :: rescaled_array
1366  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
1367  real :: scaling
1368  integer :: iounit !< Log IO unit
1369  integer :: i, j, k, Is, Js
1370  real :: aMean, aMin, aMax
1371  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1372  integer :: bcN, bcS, bcE, bcW
1373  logical :: do_corners, sym, sym_stats
1374  integer :: turns ! Quarter turns from input to model grid
1375 
1376  ! Rotate array to the input grid
1377  turns = hi_m%turns
1378  if (modulo(turns, 4) /= 0) then
1379  allocate(hi)
1380  call rotate_hor_index(hi_m, -turns, hi)
1381  allocate(array(hi%IsdB:hi%IedB, hi%JsdB:hi%JedB, size(array_m, 3)))
1382  call rotate_array(array_m, -turns, array)
1383  else
1384  hi => hi_m
1385  array => array_m
1386  endif
1387 
1388  if (checkfornans) then
1389  if (is_nan(array(hi%IscB:hi%IecB,hi%JscB:hi%JecB,:))) &
1390  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1391 ! if (is_NaN(array)) &
1392 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1393  endif
1394 
1395  scaling = 1.0 ; if (present(scale)) scaling = scale
1396  iounit = error_unit; if(present(logunit)) iounit = logunit
1397  sym_stats = .false. ; if (present(symmetric)) sym_stats = symmetric
1398  if (present(haloshift)) then ; if (haloshift > 0) sym_stats = .true. ; endif
1399 
1400  if (calculatestatistics) then
1401  if (present(scale)) then
1402  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1403  lbound(array,2):ubound(array,2), &
1404  lbound(array,3):ubound(array,3)) )
1405  rescaled_array(:,:,:) = 0.0
1406  is = hi%isc ; if (sym_stats) is = hi%isc-1
1407  js = hi%jsc ; if (sym_stats) js = hi%jsc-1
1408  do k=1,size(array,3) ; do j=js,hi%JecB ; do i=is,hi%IecB
1409  rescaled_array(i,j,k) = scale*array(i,j,k)
1410  enddo ; enddo ; enddo
1411  call substats(hi, rescaled_array, sym_stats, amean, amin, amax)
1412  deallocate(rescaled_array)
1413  else
1414  call substats(hi, array, sym_stats, amean, amin, amax)
1415  endif
1416 
1417  if (is_root_pe()) &
1418  call chk_sum_msg("B-point:", amean, amin, amax, mesg, iounit)
1419  endif
1420 
1421  if (.not.writechksums) return
1422 
1423  hshift = default_shift
1424  if (present(haloshift)) hshift = haloshift
1425  if (hshift<0) hshift = hi%ied-hi%iec
1426 
1427  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1428  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1429  write(0,*) 'chksum_B_3d: haloshift =',hshift
1430  write(0,*) 'chksum_B_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1431  write(0,*) 'chksum_B_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1432  call chksum_error(fatal,'Error in chksum_B_3d '//trim(mesg))
1433  endif
1434 
1435  bc0 = subchk(array, hi, 0, 0, scaling)
1436 
1437  sym = .false. ; if (present(symmetric)) sym = symmetric
1438 
1439  if ((hshift==0) .and. .not.sym) then
1440  if (is_root_pe()) call chk_sum_msg("B-point:", bc0, mesg, iounit)
1441  return
1442  endif
1443 
1444  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1445 
1446  if (do_corners) then
1447  if (sym) then
1448  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
1449  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1450  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
1451  else
1452  bcsw = subchk(array, hi, -hshift-1, -hshift-1, scaling)
1453  bcse = subchk(array, hi, hshift, -hshift-1, scaling)
1454  bcnw = subchk(array, hi, -hshift-1, hshift, scaling)
1455  endif
1456  bcne = subchk(array, hi, hshift, hshift, scaling)
1457 
1458  if (is_root_pe()) &
1459  call chk_sum_msg("B-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1460  else
1461  if (sym) then
1462  bcs = subchk(array, hi, 0, -hshift-1, scaling)
1463  bcw = subchk(array, hi, -hshift-1, 0, scaling)
1464  else
1465  bcs = subchk(array, hi, 0, -hshift, scaling)
1466  bcw = subchk(array, hi, -hshift, 0, scaling)
1467  endif
1468  bce = subchk(array, hi, hshift, 0, scaling)
1469  bcn = subchk(array, hi, 0, hshift, scaling)
1470 
1471  if (is_root_pe()) &
1472  call chk_sum_msg_nsew("B-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1473  endif
1474 
1475  contains
1476 
1477  integer function subchk(array, HI, di, dj, scale)
1478  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1479  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1480  integer, intent(in) :: di !< i- direction array shift for this checksum
1481  integer, intent(in) :: dj !< j- direction array shift for this checksum
1482  real, intent(in) :: scale !< A scaling factor for this array.
1483  integer :: i, j, k, bc
1484  subchk = 0
1485  ! This line deliberately uses the h-point computational domain.
1486  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1487  bc = bitcount(abs(scale*array(i,j,k)))
1488  subchk = subchk + bc
1489  enddo ; enddo ; enddo
1490  call sum_across_pes(subchk)
1491  subchk=mod(subchk, bc_modulus)
1492  end function subchk
1493 
1494  subroutine substats(HI, array, sym_stats, aMean, aMin, aMax)
1495  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1496  real, dimension(HI%IsdB:,HI%JsdB:,:), intent(in) :: array !< The array to be checksummed
1497  logical, intent(in) :: sym_stats !< If true, evaluate the statistics on the
1498  !! full symmetric computational domain.
1499  real, intent(out) :: aMean !< Array mean
1500  real, intent(out) :: aMin !< Array minimum
1501  real, intent(out) :: aMax !< Array maximum
1502 
1503  integer :: i, j, k, n, IsB, JsB
1504 
1505  isb = hi%isc ; if (sym_stats) isb = hi%isc-1
1506  jsb = hi%jsc ; if (sym_stats) jsb = hi%jsc-1
1507 
1508  amin = array(hi%isc,hi%jsc,1) ; amax = amin
1509  do k=lbound(array,3),ubound(array,3) ; do j=jsb,hi%JecB ; do i=isb,hi%IecB
1510  amin = min(amin, array(i,j,k))
1511  amax = max(amax, array(i,j,k))
1512  enddo ; enddo ; enddo
1513  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1514  n = (1 + hi%jec - hi%jsc) * (1 + hi%iec - hi%isc) * size(array,3)
1515  call sum_across_pes(n)
1516  call min_across_pes(amin)
1517  call max_across_pes(amax)
1518  amean = amean / real(n)
1519  end subroutine substats
1520 

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