MOM6
mom_coms::reproducing_sum_efp Interface Reference

Detailed Description

Find an accurate and order-invariant sum of a distributed 2d field, returning the result in the form of an extended fixed point value that can be converted back with EFP_to_real.

Definition at line 59 of file MOM_coms.F90.

Private functions

type(efp_type) function reproducing_efp_sum_2d (array, isr, ier, jsr, jer, overflow_check, err, only_on_PE)
 This subroutine uses a conversion to an integer representation of real numbers to give an order-invariant sum of distributed 2-D arrays that reproduces across domain decomposition, with the result returned as an extended fixed point type that can be converted back to a real number using EFP_to_real. This technique is described in Hallberg & Adcroft, 2014, Parallel Computing, doi:10.1016/j.parco.2014.04.007. More...
 

Detailed Description

Find an accurate and order-invariant sum of a distributed 2d field, returning the result in the form of an extended fixed point value that can be converted back with EFP_to_real.

Definition at line 59 of file MOM_coms.F90.

Functions and subroutines

◆ reproducing_efp_sum_2d()

type(efp_type) function mom_coms::reproducing_sum_efp::reproducing_efp_sum_2d ( real, dimension(:,:), intent(in)  array,
integer, intent(in), optional  isr,
integer, intent(in), optional  ier,
integer, intent(in), optional  jsr,
integer, intent(in), optional  jer,
logical, intent(in), optional  overflow_check,
integer, intent(out), optional  err,
logical, intent(in), optional  only_on_PE 
)
private

This subroutine uses a conversion to an integer representation of real numbers to give an order-invariant sum of distributed 2-D arrays that reproduces across domain decomposition, with the result returned as an extended fixed point type that can be converted back to a real number using EFP_to_real. This technique is described in Hallberg & Adcroft, 2014, Parallel Computing, doi:10.1016/j.parco.2014.04.007.

Parameters
[in]arrayThe array to be summed
[in]isrThe starting i-index of the sum, noting that the array indices starts at 1
[in]ierThe ending i-index of the sum, noting that the array indices starts at 1
[in]jsrThe starting j-index of the sum, noting that the array indices starts at 1
[in]jerThe ending j-index of the sum, noting that the array indices starts at 1
[in]overflow_checkIf present and false, disable checking for overflows in incremental results. This can speed up calculations if the number of values being summed is small enough
[out]errIf present, return an error code instead of triggering any fatal errors directly from this routine.
[in]only_on_peIf present and true, do not do the sum across processors, only reporting the local sum
Returns
The result in extended fixed point format

Definition at line 93 of file MOM_coms.F90.

93  real, dimension(:,:), intent(in) :: array !< The array to be summed
94  integer, optional, intent(in) :: isr !< The starting i-index of the sum, noting
95  !! that the array indices starts at 1
96  integer, optional, intent(in) :: ier !< The ending i-index of the sum, noting
97  !! that the array indices starts at 1
98  integer, optional, intent(in) :: jsr !< The starting j-index of the sum, noting
99  !! that the array indices starts at 1
100  integer, optional, intent(in) :: jer !< The ending j-index of the sum, noting
101  !! that the array indices starts at 1
102  logical, optional, intent(in) :: overflow_check !< If present and false, disable
103  !! checking for overflows in incremental results.
104  !! This can speed up calculations if the number
105  !! of values being summed is small enough
106  integer, optional, intent(out) :: err !< If present, return an error code instead of
107  !! triggering any fatal errors directly from
108  !! this routine.
109  logical, optional, intent(in) :: only_on_PE !< If present and true, do not do the sum
110  !! across processors, only reporting the local sum
111  type(EFP_type) :: EFP_sum !< The result in extended fixed point format
112 
113  ! This subroutine uses a conversion to an integer representation
114  ! of real numbers to give order-invariant sums that will reproduce
115  ! across PE count. This idea comes from R. Hallberg and A. Adcroft.
116 
117  integer(kind=8), dimension(ni) :: ints_sum
118  integer(kind=8) :: ival, prec_error
119  real :: rs
120  real :: max_mag_term
121  logical :: over_check, do_sum_across_PEs
122  character(len=256) :: mesg
123  integer :: i, j, n, is, ie, js, je, sgn
124 
125  if (num_pes() > max_count_prec) call mom_error(fatal, &
126  "reproducing_sum: Too many processors are being used for the value of "//&
127  "prec. Reduce prec to (2^63-1)/num_PEs.")
128 
129  prec_error = (2_8**62 + (2_8**62 - 1)) / num_pes()
130 
131  is = 1 ; ie = size(array,1) ; js = 1 ; je = size(array,2 )
132  if (present(isr)) then
133  if (isr < is) call mom_error(fatal, "Value of isr too small in reproducing_EFP_sum_2d.")
134  is = isr
135  endif
136  if (present(ier)) then
137  if (ier > ie) call mom_error(fatal, "Value of ier too large in reproducing_EFP_sum_2d.")
138  ie = ier
139  endif
140  if (present(jsr)) then
141  if (jsr < js) call mom_error(fatal, "Value of jsr too small in reproducing_EFP_sum_2d.")
142  js = jsr
143  endif
144  if (present(jer)) then
145  if (jer > je) call mom_error(fatal, "Value of jer too large in reproducing_EFP_sum_2d.")
146  je = jer
147  endif
148 
149  over_check = .true. ; if (present(overflow_check)) over_check = overflow_check
150  do_sum_across_pes = .true. ; if (present(only_on_pe)) do_sum_across_pes = .not.only_on_pe
151 
152  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
153  ints_sum(:) = 0
154  if (over_check) then
155  if ((je+1-js)*(ie+1-is) < max_count_prec) then
156  do j=js,je ; do i=is,ie
157  call increment_ints_faster(ints_sum, array(i,j), max_mag_term)
158  enddo ; enddo
159  call carry_overflow(ints_sum, prec_error)
160  elseif ((ie+1-is) < max_count_prec) then
161  do j=js,je
162  do i=is,ie
163  call increment_ints_faster(ints_sum, array(i,j), max_mag_term)
164  enddo
165  call carry_overflow(ints_sum, prec_error)
166  enddo
167  else
168  do j=js,je ; do i=is,ie
169  call increment_ints(ints_sum, real_to_ints(array(i,j), prec_error), &
170  prec_error)
171  enddo ; enddo
172  endif
173  else
174  do j=js,je ; do i=is,ie
175  sgn = 1 ; if (array(i,j)<0.0) sgn = -1
176  rs = abs(array(i,j))
177  do n=1,ni
178  ival = int(rs*i_pr(n), 8)
179  rs = rs - ival*pr(n)
180  ints_sum(n) = ints_sum(n) + sgn*ival
181  enddo
182  enddo ; enddo
183  call carry_overflow(ints_sum, prec_error)
184  endif
185 
186  if (present(err)) then
187  err = 0
188  if (overflow_error) &
189  err = err+2
190  if (nan_error) &
191  err = err+4
192  if (err > 0) then ; do n=1,ni ; ints_sum(n) = 0 ; enddo ; endif
193  else
194  if (nan_error) then
195  call mom_error(fatal, "NaN in input field of reproducing_EFP_sum(_2d).")
196  endif
197  if (abs(max_mag_term) >= prec_error*pr(1)) then
198  write(mesg, '(ES13.5)') max_mag_term
199  call mom_error(fatal,"Overflow in reproducing_EFP_sum(_2d) conversion of "//trim(mesg))
200  endif
201  if (overflow_error) then
202  call mom_error(fatal, "Overflow in reproducing_EFP_sum(_2d).")
203  endif
204  endif
205 
206  if (do_sum_across_pes) call sum_across_pes(ints_sum, ni)
207 
208  call regularize_ints(ints_sum)
209 
210  efp_sum%v(:) = ints_sum(:)
211 

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