MOM6
mom_checksums::hchksum Interface Reference

Detailed Description

Definition at line 56 of file MOM_checksums.F90.

Private functions

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_h_3d (array, mesg, HI, haloshift, omit_corners, scale)
 chksum_h_3d performs checksums on a 3d array staggered at tracer points. More...
 

Functions and subroutines

◆ chksum_h_2d()

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

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)

◆ chksum_h_3d()

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

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)

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