MOM6
mom_coms::reproducing_sum Interface Reference

Detailed Description

Definition at line 62 of file MOM_coms.F90.

Private functions

real function reproducing_sum_2d (array, isr, ier, jsr, jer, EFP_sum, reproducing, overflow_check, err)
 
real function reproducing_sum_3d (array, isr, ier, jsr, jer, sums, EFP_sum, err)
 

Functions and subroutines

◆ reproducing_sum_2d()

real function mom_coms::reproducing_sum::reproducing_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,
type(efp_type), intent(out), optional  EFP_sum,
logical, intent(in), optional  reproducing,
logical, intent(in), optional  overflow_check,
integer, intent(out), optional  err 
)
private

Definition at line 80 of file MOM_coms.F90.

80  real, dimension(:,:), intent(in) :: array
81  integer, optional, intent(in) :: isr, ier, jsr, jer
82  type(efp_type), optional, intent(out) :: efp_sum
83  logical, optional, intent(in) :: reproducing
84  logical, optional, intent(in) :: overflow_check
85  integer, optional, intent(out) :: err
86  real :: sum ! Result
87 
88  ! This subroutine uses a conversion to an integer representation
89  ! of real numbers to give order-invariant sums that will reproduce
90  ! across PE count. This idea comes from R. Hallberg and A. Adcroft.
91 
92  integer(kind=8), dimension(ni) :: ints_sum
93  integer(kind=8) :: ival, prec_error
94  real :: rsum(1), rs
95  real :: max_mag_term
96  logical :: repro, over_check
97  character(len=256) :: mesg
98  integer :: i, j, n, is, ie, js, je, sgn
99 
100  if (num_pes() > max_count_prec) call mom_error(fatal, &
101  "reproducing_sum: Too many processors are being used for the value of "//&
102  "prec. Reduce prec to (2^63-1)/num_PEs.")
103 
104  prec_error = (2_8**62 + (2_8**62 - 1)) / num_pes()
105 
106  is = 1 ; ie = size(array,1) ; js = 1 ; je = size(array,2 )
107  if (present(isr)) then
108  if (isr < is) call mom_error(fatal, &
109  "Value of isr too small in reproducing_sum_2d.")
110  is = isr
111  endif
112  if (present(ier)) then
113  if (ier > ie) call mom_error(fatal, &
114  "Value of ier too large in reproducing_sum_2d.")
115  ie = ier
116  endif
117  if (present(jsr)) then
118  if (jsr < js) call mom_error(fatal, &
119  "Value of jsr too small in reproducing_sum_2d.")
120  js = jsr
121  endif
122  if (present(jer)) then
123  if (jer > je) call mom_error(fatal, &
124  "Value of jer too large in reproducing_sum_2d.")
125  je = jer
126  endif
127 
128  repro = .true. ; if (present(reproducing)) repro = reproducing
129  over_check = .true. ; if (present(overflow_check)) over_check = overflow_check
130 
131  if (repro) then
132  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
133  ints_sum(:) = 0
134  if (over_check) then
135  if ((je+1-js)*(ie+1-is) < max_count_prec) then
136  do j=js,je ; do i=is,ie
137  call increment_ints_faster(ints_sum, array(i,j), max_mag_term);
138  enddo ; enddo
139  call carry_overflow(ints_sum, prec_error)
140  elseif ((ie+1-is) < max_count_prec) then
141  do j=js,je
142  do i=is,ie
143  call increment_ints_faster(ints_sum, array(i,j), max_mag_term);
144  enddo
145  call carry_overflow(ints_sum, prec_error)
146  enddo
147  else
148  do j=js,je ; do i=is,ie
149  call increment_ints(ints_sum, real_to_ints(array(i,j), prec_error), &
150  prec_error);
151  enddo ; enddo
152  endif
153  else
154  do j=js,je ; do i=is,ie
155  sgn = 1 ; if (array(i,j)<0.0) sgn = -1
156  rs = abs(array(i,j))
157  do n=1,ni
158  ival = int(rs*i_pr(n), 8)
159  rs = rs - ival*pr(n)
160  ints_sum(n) = ints_sum(n) + sgn*ival
161  enddo
162  enddo ; enddo
163  call carry_overflow(ints_sum, prec_error)
164  endif
165 
166  if (present(err)) then
167  err = 0
168  if (overflow_error) &
169  err = err+2
170  if (nan_error) &
171  err = err+4
172  if (err > 0) then ; do n=1,ni ; ints_sum(n) = 0 ; enddo ; endif
173  else
174  if (nan_error) then
175  call mom_error(fatal, "NaN in input field of reproducing_sum(_2d).")
176  endif
177  if (abs(max_mag_term) >= prec_error*pr(1)) then
178  write(mesg, '(ES13.5)') max_mag_term
179  call mom_error(fatal,"Overflow in reproducing_sum(_2d) conversion of "//trim(mesg))
180  endif
181  if (overflow_error) then
182  call mom_error(fatal, "Overflow in reproducing_sum(_2d).")
183  endif
184  endif
185 
186  call sum_across_pes(ints_sum, ni)
187 
188  call regularize_ints(ints_sum)
189  sum = ints_to_real(ints_sum)
190  else
191  rsum(1) = 0.0
192  do j=js,je ; do i=is,ie
193  rsum(1) = rsum(1) + array(i,j);
194  enddo ; enddo
195  call sum_across_pes(rsum,1)
196  sum = rsum(1)
197 
198  if (present(err)) then ; err = 0 ; endif
199 
200  if (debug .or. present(efp_sum)) then
201  overflow_error = .false.
202  ints_sum = real_to_ints(sum, prec_error, overflow_error)
203  if (overflow_error) then
204  if (present(err)) then
205  err = err + 2
206  else
207  write(mesg, '(ES13.5)') sum
208  call mom_error(fatal,"Repro_sum_2d: Overflow in real_to_ints conversion of "//trim(mesg))
209  endif
210  endif
211  endif
212  endif
213 
214  if (present(efp_sum)) efp_sum%v(:) = ints_sum(:)
215 
216  if (debug) then
217  write(mesg,'("2d RS: ", ES24.16, 6 Z17.16)') sum, ints_sum(1:ni)
218  call mom_mesg(mesg, 3)
219  endif
220 

◆ reproducing_sum_3d()

real function mom_coms::reproducing_sum::reproducing_sum_3d ( 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,
real, dimension(:), intent(out), optional  sums,
type(efp_type), intent(out), optional  EFP_sum,
integer, intent(out), optional  err 
)
private

Definition at line 225 of file MOM_coms.F90.

225  real, dimension(:,:,:), intent(in) :: array
226  integer, optional, intent(in) :: isr, ier, jsr, jer
227  real, dimension(:), optional, intent(out) :: sums
228  type(efp_type), optional, intent(out) :: efp_sum
229  integer, optional, intent(out) :: err
230  real :: sum ! Result
231 
232  ! This subroutine uses a conversion to an integer representation
233  ! of real numbers to give order-invariant sums that will reproduce
234  ! across PE count. This idea comes from R. Hallberg and A. Adcroft.
235 
236  real :: max_mag_term
237  integer(kind=8), dimension(ni) :: ints_sum
238  integer(kind=8), dimension(ni,size(array,3)) :: ints_sums
239  integer(kind=8) :: prec_error
240  character(len=256) :: mesg
241  integer :: i, j, k, is, ie, js, je, ke, isz, jsz, n
242 
243  if (num_pes() > max_count_prec) call mom_error(fatal, &
244  "reproducing_sum: Too many processors are being used for the value of "//&
245  "prec. Reduce prec to (2^63-1)/num_PEs.")
246 
247  prec_error = (2_8**62 + (2_8**62 - 1)) / num_pes()
248  max_mag_term = 0.0
249 
250  is = 1 ; ie = size(array,1) ; js = 1 ; je = size(array,2) ; ke = size(array,3)
251  if (present(isr)) then
252  if (isr < is) call mom_error(fatal, &
253  "Value of isr too small in reproducing_sum(_3d).")
254  is = isr
255  endif
256  if (present(ier)) then
257  if (ier > ie) call mom_error(fatal, &
258  "Value of ier too large in reproducing_sum(_3d).")
259  ie = ier
260  endif
261  if (present(jsr)) then
262  if (jsr < js) call mom_error(fatal, &
263  "Value of jsr too small in reproducing_sum(_3d).")
264  js = jsr
265  endif
266  if (present(jer)) then
267  if (jer > je) call mom_error(fatal, &
268  "Value of jer too large in reproducing_sum(_3d).")
269  je = jer
270  endif
271  jsz = je+1-js; isz = ie+1-is
272 
273  if (present(sums)) then
274  if (size(sums) > ke) call mom_error(fatal, "Sums is smaller than "//&
275  "the vertical extent of array in reproducing_sum(_3d).")
276  ints_sums(:,:) = 0
277  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
278  if (jsz*isz < max_count_prec) then
279  do k=1,ke
280  do j=js,je ; do i=is,ie
281  call increment_ints_faster(ints_sums(:,k), array(i,j,k), max_mag_term);
282  enddo ; enddo
283  call carry_overflow(ints_sums(:,k), prec_error)
284  enddo
285  elseif (isz < max_count_prec) then
286  do k=1,ke ; do j=js,je
287  do i=is,ie
288  call increment_ints_faster(ints_sums(:,k), array(i,j,k), max_mag_term);
289  enddo
290  call carry_overflow(ints_sums(:,k), prec_error)
291  enddo ; enddo
292  else
293  do k=1,ke ; do j=js,je ; do i=is,ie
294  call increment_ints(ints_sums(:,k), &
295  real_to_ints(array(i,j,k), prec_error), prec_error);
296  enddo ; enddo ; enddo
297  endif
298  if (present(err)) then
299  err = 0
300  if (abs(max_mag_term) >= prec_error*pr(1)) err = err+1
301  if (overflow_error) err = err+2
302  if (nan_error) err = err+2
303  if (err > 0) then ; do k=1,ke ; do n=1,ni ; ints_sums(n,k) = 0 ; enddo ; enddo ; endif
304  else
305  if (nan_error) call mom_error(fatal, "NaN in input field of reproducing_sum(_3d).")
306  if (abs(max_mag_term) >= prec_error*pr(1)) then
307  write(mesg, '(ES13.5)') max_mag_term
308  call mom_error(fatal,"Overflow in reproducing_sum(_3d) conversion of "//trim(mesg))
309  endif
310  if (overflow_error) call mom_error(fatal, "Overflow in reproducing_sum(_3d).")
311  endif
312 
313  call sum_across_pes(ints_sums(:,1:ke), ni*ke)
314 
315  sum = 0.0
316  do k=1,ke
317  call regularize_ints(ints_sums(:,k))
318  sums(k) = ints_to_real(ints_sums(:,k))
319  sum = sum + sums(k)
320  enddo
321 
322  if (present(efp_sum)) then
323  efp_sum%v(:) = 0
324  do k=1,ke ; call increment_ints(efp_sum%v(:), ints_sums(:,k)) ; enddo
325  endif
326 
327  if (debug) then
328  do n=1,ni ; ints_sum(n) = 0 ; enddo
329  do k=1,ke ; do n=1,ni ; ints_sum(n) = ints_sum(n) + ints_sums(n,k) ; enddo ; enddo
330  write(mesg,'("3D RS: ", ES24.16, 6 Z17.16)') sum, ints_sum(1:ni)
331  call mom_mesg(mesg, 3)
332  endif
333  else
334  ints_sum(:) = 0
335  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
336  if (jsz*isz < max_count_prec) then
337  do k=1,ke
338  do j=js,je ; do i=is,ie
339  call increment_ints_faster(ints_sum, array(i,j,k), max_mag_term);
340  enddo ; enddo
341  call carry_overflow(ints_sum, prec_error)
342  enddo
343  elseif (isz < max_count_prec) then
344  do k=1,ke ; do j=js,je
345  do i=is,ie
346  call increment_ints_faster(ints_sum, array(i,j,k), max_mag_term);
347  enddo
348  call carry_overflow(ints_sum, prec_error)
349  enddo ; enddo
350  else
351  do k=1,ke ; do j=js,je ; do i=is,ie
352  call increment_ints(ints_sum, real_to_ints(array(i,j,k), prec_error), &
353  prec_error);
354  enddo ; enddo ; enddo
355  endif
356  if (present(err)) then
357  err = 0
358  if (abs(max_mag_term) >= prec_error*pr(1)) err = err+1
359  if (overflow_error) err = err+2
360  if (nan_error) err = err+2
361  if (err > 0) then ; do n=1,ni ; ints_sum(n) = 0 ; enddo ; endif
362  else
363  if (nan_error) call mom_error(fatal, "NaN in input field of reproducing_sum(_3d).")
364  if (abs(max_mag_term) >= prec_error*pr(1)) then
365  write(mesg, '(ES13.5)') max_mag_term
366  call mom_error(fatal,"Overflow in reproducing_sum(_3d) conversion of "//trim(mesg))
367  endif
368  if (overflow_error) call mom_error(fatal, "Overflow in reproducing_sum(_3d).")
369  endif
370 
371  call sum_across_pes(ints_sum, ni)
372 
373  call regularize_ints(ints_sum)
374  sum = ints_to_real(ints_sum)
375 
376  if (present(efp_sum)) efp_sum%v(:) = ints_sum(:)
377 
378  if (debug) then
379  write(mesg,'("3d RS: ", ES24.16, 6 Z17.16)') sum, ints_sum(1:ni)
380  call mom_mesg(mesg, 3)
381  endif
382  endif
383 

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