MOM6
mom_checksums::hchksum Interface Reference

Detailed Description

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

Definition at line 49 of file MOM_checksums.F90.

Private functions

subroutine chksum_h_2d (array_m, mesg, HI_m, haloshift, omit_corners, scale, logunit)
 Checksums a 2d array staggered at tracer points. More...
 
subroutine chksum_h_3d (array_m, mesg, HI_m, haloshift, omit_corners, scale, logunit)
 Checksums a 3d array staggered at tracer points. More...
 

Detailed Description

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

Definition at line 49 of file MOM_checksums.F90.

Functions and subroutines

◆ chksum_h_2d()

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

Checksums a 2d array staggered at tracer points.

Parameters
[in]hi_mHorizontal index bounds of the model grid
[in]array_mField array on the model grid
[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.
[in]logunitIO unit for checksum logging

Definition at line 305 of file MOM_checksums.F90.

306  type(hor_index_type), target, intent(in) :: HI_m !< Horizontal index bounds of the model grid
307  real, dimension(HI_m%isd:,HI_m%jsd:), target, intent(in) :: array_m !< Field array on the model grid
308  character(len=*), intent(in) :: mesg !< An identifying message
309  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
310  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
311  real, optional, intent(in) :: scale !< A scaling factor for this array.
312  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
313 
314  real, pointer :: array(:,:) ! Field array on the input grid
315  real, allocatable, dimension(:,:) :: rescaled_array
316  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
317  real :: scaling
318  integer :: iounit !< Log IO unit
319  integer :: i, j
320  real :: aMean, aMin, aMax
321  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
322  integer :: bcN, bcS, bcE, bcW
323  logical :: do_corners
324  integer :: turns ! Quarter turns from input to model grid
325 
326  ! Rotate array to the input grid
327  turns = hi_m%turns
328  if (modulo(turns, 4) /= 0) then
329  allocate(hi)
330  call rotate_hor_index(hi_m, -turns, hi)
331  allocate(array(hi%isd:hi%ied, hi%jsd:hi%jed))
332  call rotate_array(array_m, -turns, array)
333  else
334  hi => hi_m
335  array => array_m
336  endif
337 
338  if (checkfornans) then
339  if (is_nan(array(hi%isc:hi%iec,hi%jsc:hi%jec))) &
340  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
341 ! if (is_NaN(array)) &
342 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
343  endif
344 
345  scaling = 1.0 ; if (present(scale)) scaling = scale
346  iounit = error_unit; if(present(logunit)) iounit = logunit
347 
348  if (calculatestatistics) then
349  if (present(scale)) then
350  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
351  lbound(array,2):ubound(array,2)) )
352  rescaled_array(:,:) = 0.0
353  do j=hi%jsc,hi%jec ; do i=hi%isc,hi%iec
354  rescaled_array(i,j) = scale*array(i,j)
355  enddo ; enddo
356  call substats(hi, rescaled_array, amean, amin, amax)
357  deallocate(rescaled_array)
358  else
359  call substats(hi, array, amean, amin, amax)
360  endif
361 
362  if (is_root_pe()) &
363  call chk_sum_msg("h-point:", amean, amin, amax, mesg, iounit)
364  endif
365 
366  if (.not.writechksums) return
367 
368  hshift = default_shift
369  if (present(haloshift)) hshift = haloshift
370  if (hshift<0) hshift = hi%ied-hi%iec
371 
372  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
373  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
374  write(0,*) 'chksum_h_2d: haloshift =',hshift
375  write(0,*) 'chksum_h_2d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
376  write(0,*) 'chksum_h_2d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
377  call chksum_error(fatal,'Error in chksum_h_2d '//trim(mesg))
378  endif
379 
380  bc0 = subchk(array, hi, 0, 0, scaling)
381 
382  if (hshift==0) then
383  if (is_root_pe()) call chk_sum_msg("h-point:", bc0, mesg, iounit)
384  return
385  endif
386 
387  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
388 
389  if (do_corners) then
390  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
391  bcse = subchk(array, hi, hshift, -hshift, scaling)
392  bcnw = subchk(array, hi, -hshift, hshift, scaling)
393  bcne = subchk(array, hi, hshift, hshift, scaling)
394 
395  if (is_root_pe()) &
396  call chk_sum_msg("h-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
397  else
398  bcs = subchk(array, hi, 0, -hshift, scaling)
399  bce = subchk(array, hi, hshift, 0, scaling)
400  bcw = subchk(array, hi, -hshift, 0, scaling)
401  bcn = subchk(array, hi, 0, hshift, scaling)
402 
403  if (is_root_pe()) &
404  call chk_sum_msg_nsew("h-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
405  endif
406 
407  contains
408  integer function subchk(array, HI, di, dj, scale)
409  type(hor_index_type), intent(in) :: HI !< A horizontal index type
410  real, dimension(HI%isd:,HI%jsd:), intent(in) :: array !< The array to be checksummed
411  integer, intent(in) :: di !< i- direction array shift for this checksum
412  integer, intent(in) :: dj !< j- direction array shift for this checksum
413  real, intent(in) :: scale !< A scaling factor for this array.
414  integer :: i, j, bc
415  subchk = 0
416  do j=hi%jsc+dj,hi%jec+dj; do i=hi%isc+di,hi%iec+di
417  bc = bitcount(abs(scale*array(i,j)))
418  subchk = subchk + bc
419  enddo ; enddo
420  call sum_across_pes(subchk)
421  subchk=mod(subchk, bc_modulus)
422  end function subchk
423 
424  subroutine substats(HI, array, aMean, aMin, aMax)
425  type(hor_index_type), intent(in) :: HI !< A horizontal index type
426  real, dimension(HI%isd:,HI%jsd:), intent(in) :: array !< The array to be checksummed
427  real, intent(out) :: aMean !< Array mean
428  real, intent(out) :: aMin !< Array minimum
429  real, intent(out) :: aMax !< Array maximum
430 
431  integer :: i, j, n
432 
433  amin = array(hi%isc,hi%jsc)
434  amax = array(hi%isc,hi%jsc)
435  n = 0
436  do j=hi%jsc,hi%jec ; do i=hi%isc,hi%iec
437  amin = min(amin, array(i,j))
438  amax = max(amax, array(i,j))
439  n = n + 1
440  enddo ; enddo
441  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec))
442  call sum_across_pes(n)
443  call min_across_pes(amin)
444  call max_across_pes(amax)
445  amean = amean / real(n)
446  end subroutine substats
447 

◆ chksum_h_3d()

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

Checksums a 3d array staggered at tracer 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]omit_cornersIf true, avoid checking diagonal shifts
[in]scaleA scaling factor for this array.
[in]logunitIO unit for checksum logging

Definition at line 1203 of file MOM_checksums.F90.

1204  type(hor_index_type), target, intent(in) :: HI_m !< A horizontal index type
1205  real, dimension(HI_m%isd:,HI_m%jsd:,:), target, intent(in) :: array_m !< The array to be checksummed
1206  character(len=*), intent(in) :: mesg !< An identifying message
1207  integer, optional, intent(in) :: haloshift !< The width of halos to check (default 0)
1208  logical, optional, intent(in) :: omit_corners !< If true, avoid checking diagonal shifts
1209  real, optional, intent(in) :: scale !< A scaling factor for this array.
1210  integer, optional, intent(in) :: logunit !< IO unit for checksum logging
1211 
1212  real, pointer :: array(:,:,:) ! Field array on the input grid
1213  real, allocatable, dimension(:,:,:) :: rescaled_array
1214  type(hor_index_type), pointer :: HI ! Horizontal index bounds of the input grid
1215  real :: scaling
1216  integer :: iounit !< Log IO unit
1217  integer :: i, j, k
1218  real :: aMean, aMin, aMax
1219  integer :: bc0, bcSW, bcSE, bcNW, bcNE, hshift
1220  integer :: bcN, bcS, bcE, bcW
1221  logical :: do_corners
1222  integer :: turns ! Quarter turns from input to model grid
1223 
1224  ! Rotate array to the input grid
1225  turns = hi_m%turns
1226  if (modulo(turns, 4) /= 0) then
1227  allocate(hi)
1228  call rotate_hor_index(hi_m, -turns, hi)
1229  allocate(array(hi%isd:hi%ied, hi%jsd:hi%jed, size(array_m, 3)))
1230  call rotate_array(array_m, -turns, array)
1231  else
1232  hi => hi_m
1233  array => array_m
1234  endif
1235 
1236  if (checkfornans) then
1237  if (is_nan(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))) &
1238  call chksum_error(fatal, 'NaN detected: '//trim(mesg))
1239 ! if (is_NaN(array)) &
1240 ! call chksum_error(FATAL, 'NaN detected in halo: '//trim(mesg))
1241  endif
1242 
1243  scaling = 1.0 ; if (present(scale)) scaling = scale
1244  iounit = error_unit; if(present(logunit)) iounit = logunit
1245 
1246  if (calculatestatistics) then
1247  if (present(scale)) then
1248  allocate( rescaled_array(lbound(array,1):ubound(array,1), &
1249  lbound(array,2):ubound(array,2), &
1250  lbound(array,3):ubound(array,3)) )
1251  rescaled_array(:,:,:) = 0.0
1252  do k=1,size(array,3) ; do j=hi%jsc,hi%jec ; do i=hi%isc,hi%iec
1253  rescaled_array(i,j,k) = scale*array(i,j,k)
1254  enddo ; enddo ; enddo
1255 
1256  call substats(hi, rescaled_array, amean, amin, amax)
1257  deallocate(rescaled_array)
1258  else
1259  call substats(hi, array, amean, amin, amax)
1260  endif
1261 
1262  if (is_root_pe()) &
1263  call chk_sum_msg("h-point:", amean, amin, amax, mesg, iounit)
1264  endif
1265 
1266  if (.not.writechksums) return
1267 
1268  hshift = default_shift
1269  if (present(haloshift)) hshift = haloshift
1270  if (hshift<0) hshift = hi%ied-hi%iec
1271 
1272  if ( hi%isc-hshift<hi%isd .or. hi%iec+hshift>hi%ied .or. &
1273  hi%jsc-hshift<hi%jsd .or. hi%jec+hshift>hi%jed ) then
1274  write(0,*) 'chksum_h_3d: haloshift =',hshift
1275  write(0,*) 'chksum_h_3d: isd,isc,iec,ied=',hi%isd,hi%isc,hi%iec,hi%ied
1276  write(0,*) 'chksum_h_3d: jsd,jsc,jec,jed=',hi%jsd,hi%jsc,hi%jec,hi%jed
1277  call chksum_error(fatal,'Error in chksum_h_3d '//trim(mesg))
1278  endif
1279 
1280  bc0 = subchk(array, hi, 0, 0, scaling)
1281 
1282  if (hshift==0) then
1283  if (is_root_pe()) call chk_sum_msg("h-point:", bc0, mesg, iounit)
1284  return
1285  endif
1286 
1287  do_corners = .true. ; if (present(omit_corners)) do_corners = .not.omit_corners
1288 
1289  if (do_corners) then
1290  bcsw = subchk(array, hi, -hshift, -hshift, scaling)
1291  bcse = subchk(array, hi, hshift, -hshift, scaling)
1292  bcnw = subchk(array, hi, -hshift, hshift, scaling)
1293  bcne = subchk(array, hi, hshift, hshift, scaling)
1294 
1295  if (is_root_pe()) &
1296  call chk_sum_msg("h-point:", bc0, bcsw, bcse, bcnw, bcne, mesg, iounit)
1297  else
1298  bcs = subchk(array, hi, 0, -hshift, scaling)
1299  bce = subchk(array, hi, hshift, 0, scaling)
1300  bcw = subchk(array, hi, -hshift, 0, scaling)
1301  bcn = subchk(array, hi, 0, hshift, scaling)
1302 
1303  if (is_root_pe()) &
1304  call chk_sum_msg_nsew("h-point:", bc0, bcn, bcs, bce, bcw, mesg, iounit)
1305  endif
1306 
1307  contains
1308 
1309  integer function subchk(array, HI, di, dj, scale)
1310  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1311  real, dimension(HI%isd:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
1312  integer, intent(in) :: di !< i- direction array shift for this checksum
1313  integer, intent(in) :: dj !< j- direction array shift for this checksum
1314  real, intent(in) :: scale !< A scaling factor for this array.
1315  integer :: i, j, k, bc
1316  subchk = 0
1317  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc+dj,hi%jec+dj ; do i=hi%isc+di,hi%iec+di
1318  bc = bitcount(abs(scale*array(i,j,k)))
1319  subchk = subchk + bc
1320  enddo ; enddo ; enddo
1321  call sum_across_pes(subchk)
1322  subchk=mod(subchk, bc_modulus)
1323  end function subchk
1324 
1325  subroutine substats(HI, array, aMean, aMin, aMax)
1326  type(hor_index_type), intent(in) :: HI !< A horizontal index type
1327  real, dimension(HI%isd:,HI%jsd:,:), intent(in) :: array !< The array to be checksummed
1328  real, intent(out) :: aMean !< Array mean
1329  real, intent(out) :: aMin !< Array minimum
1330  real, intent(out) :: aMax !< Array maximum
1331 
1332  integer :: i, j, k, n
1333 
1334  amin = array(hi%isc,hi%jsc,1)
1335  amax = array(hi%isc,hi%jsc,1)
1336  n = 0
1337  do k=lbound(array,3),ubound(array,3) ; do j=hi%jsc,hi%jec ; do i=hi%isc,hi%iec
1338  amin = min(amin, array(i,j,k))
1339  amax = max(amax, array(i,j,k))
1340  n = n + 1
1341  enddo ; enddo ; enddo
1342  amean = reproducing_sum(array(hi%isc:hi%iec,hi%jsc:hi%jec,:))
1343  call sum_across_pes(n)
1344  call min_across_pes(amin)
1345  call max_across_pes(amax)
1346  amean = amean / real(n)
1347  end subroutine substats
1348 

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