FV3 Bundle
mpp_efp.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
20 #include <fms_platform.h>
21 
22 use mpp_mod, only : mpp_error, fatal, warning, note
23 use mpp_mod, only : mpp_pe, mpp_root_pe, mpp_npes
24 use mpp_mod, only : mpp_sum
25 
26 implicit none ; private
27 
30 public :: operator(+), operator(-), assignment(=)
32 
33 ! This module provides interfaces to the non-domain-oriented communication
34 ! subroutines.
35 integer, parameter :: numbit = 46 ! number of bits used in the 64-bit signed integer representation.
36 integer, parameter :: numint = 6 ! The number of long integers to use to represent
37  ! a real number.
38 
39 integer(LONG_KIND), parameter :: prec=2_8**numbit ! The precision of each integer.
40 real(DOUBLE_KIND), parameter :: r_prec=2.0_8**numbit ! A real version of prec.
41 real(DOUBLE_KIND), parameter :: i_prec=1.0_8/(2.0_8**numbit) ! The inverse of prec.
42 integer, parameter :: max_count_prec=2**(63-numbit)-1
43  ! The number of values that can be added together
44  ! with the current value of prec before there will
45  ! be roundoff problems.
46 
47 real(DOUBLE_KIND), parameter, dimension(NUMINT) :: &
48  pr = (/ r_prec**2, r_prec, 1.0_8, 1.0_8/r_prec, 1.0_8/r_prec**2, 1.0_8/r_prec**3 /)
49 real(DOUBLE_KIND), parameter, dimension(NUMINT) :: &
50  i_pr = (/ 1.0_8/r_prec**2, 1.0_8/r_prec, 1.0_8, r_prec, r_prec**2, r_prec**3 /)
51 
52 logical :: overflow_error = .false., nan_error = .false.
53 logical :: debug = .false. ! Making this true enables debugging output.
54 
56  module procedure mpp_reproducing_sum_r8_2d
57  module procedure mpp_reproducing_sum_r8_3d
58  module procedure mpp_reproducing_sum_r4_2d
59 end interface mpp_reproducing_sum
60 
61 ! The Extended Fixed Point (mpp_efp) type provides a public interface for doing
62 ! sums and taking differences with this type.
63 type, public :: mpp_efp_type ; private
64  integer(kind=8), dimension(NUMINT) :: v
65 end type mpp_efp_type
66 
67 interface operator (+); module procedure mpp_efp_plus ; end interface
68 interface operator (-); module procedure mpp_efp_minus ; end interface
69 interface assignment(=); module procedure mpp_efp_assign ; end interface
70 
71 contains
72 
73 function mpp_reproducing_sum_r8_2d(array, isr, ier, jsr, jer, EFP_sum, reproducing, &
74  overflow_check, err) result(sum)
75  real(DOUBLE_KIND), dimension(:,:), intent(in) :: array
76  integer, optional, intent(in) :: isr, ier, jsr, jer
77  type(mpp_efp_type), optional, intent(out) :: efp_sum
78  logical, optional, intent(in) :: reproducing
79  logical, optional, intent(in) :: overflow_check
80  integer, optional, intent(out) :: err
81  real(DOUBLE_KIND) :: sum ! Result
82 
83  ! This subroutine uses a conversion to an integer representation
84  ! of real numbers to give order-invariant sums that will reproduce
85  ! across PE count. This idea comes from R. Hallberg and A. Adcroft.
86 
87  integer(LONG_KIND), dimension(NUMINT) :: ints_sum
88  integer(LONG_KIND) :: ival, prec_error
89  real(DOUBLE_KIND) :: rsum(1), rs
90  real(DOUBLE_KIND) :: max_mag_term
91  logical :: repro, over_check
92  character(len=256) :: mesg
93  integer :: i, j, n, is, ie, js, je, sgn
94 
95  if (mpp_npes() > max_count_prec) call mpp_error(fatal, &
96  "mpp_reproducing_sum: Too many processors are being used for the value of "//&
97  "prec. Reduce prec to (2^63-1)/mpp_npes.")
98 
99  prec_error = (2_8**62 + (2_8**62 - 1)) / mpp_npes()
100 
101  is = 1 ; ie = size(array,1) ; js = 1 ; je = size(array,2 )
102  if (present(isr)) then
103  if (isr < is) call mpp_error(fatal, &
104  "Value of isr too small in mpp_reproducing_sum_2d.")
105  is = isr
106  endif
107  if (present(ier)) then
108  if (ier > ie) call mpp_error(fatal, &
109  "Value of ier too large in mpp_reproducing_sum_2d.")
110  ie = ier
111  endif
112  if (present(jsr)) then
113  if (jsr < js) call mpp_error(fatal, &
114  "Value of jsr too small in mpp_reproducing_sum_2d.")
115  js = jsr
116  endif
117  if (present(jer)) then
118  if (jer > je) call mpp_error(fatal, &
119  "Value of jer too large in mpp_reproducing_sum_2d.")
120  je = jer
121  endif
122 
123  repro = .true. ; if (present(reproducing)) repro = reproducing
124  over_check = .true. ; if (present(overflow_check)) over_check = overflow_check
125 
126  if (repro) then
127  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
128  ints_sum(:) = 0
129  if (over_check) then
130  if ((je+1-js)*(ie+1-is) < max_count_prec) then
131  do j=js,je ; do i=is,ie
132  call increment_ints_faster(ints_sum, array(i,j), max_mag_term);
133  enddo ; enddo
134  call carry_overflow(ints_sum, prec_error)
135  elseif ((ie+1-is) < max_count_prec) then
136  do j=js,je
137  do i=is,ie
138  call increment_ints_faster(ints_sum, array(i,j), max_mag_term);
139  enddo
140  call carry_overflow(ints_sum, prec_error)
141  enddo
142  else
143  do j=js,je ; do i=is,ie
144  call increment_ints(ints_sum, real_to_ints(array(i,j), prec_error), &
145  prec_error);
146  enddo ; enddo
147  endif
148  else
149  do j=js,je ; do i=is,ie
150  sgn = 1 ; if (array(i,j)<0.0) sgn = -1
151  rs = abs(array(i,j))
152  do n=1,numint
153  ival = int(rs*i_pr(n), 8)
154  rs = rs - ival*pr(n)
155  ints_sum(n) = ints_sum(n) + sgn*ival
156  enddo
157  enddo ; enddo
158  call carry_overflow(ints_sum, prec_error)
159  endif
160 
161  if (present(err)) then
162  err = 0
163  if (overflow_error) &
164  err = err+2
165  if (nan_error) &
166  err = err+4
167  if (err > 0) then ; do n=1,numint ; ints_sum(n) = 0 ; enddo ; endif
168  else
169  if (nan_error) then
170  call mpp_error(fatal, "NaN in input field of mpp_reproducing_sum(_2d), this indicates numerical instability")
171  endif
172  if (abs(max_mag_term) >= prec_error*pr(1)) then
173  write(mesg, '(ES13.5)') max_mag_term
174  call mpp_error(fatal,"Overflow in mpp_reproducing_sum(_2d) conversion of "//trim(mesg))
175  endif
176  if (overflow_error) then
177  call mpp_error(fatal, "Overflow in mpp_reproducing_sum(_2d).")
178  endif
179  endif
180 
181  call mpp_sum(ints_sum, numint)
182 
183  call regularize_ints(ints_sum)
184  sum = ints_to_real(ints_sum)
185  else
186  rsum(1) = 0.0
187  do j=js,je ; do i=is,ie
188  rsum(1) = rsum(1) + array(i,j);
189  enddo ; enddo
190  call mpp_sum(rsum,1)
191  sum = rsum(1)
192 
193  if (present(err)) then ; err = 0 ; endif
194 
195  if (debug .or. present(efp_sum)) then
196  overflow_error = .false.
197  ints_sum = real_to_ints(sum, prec_error, overflow_error)
198  if (overflow_error) then
199  if (present(err)) then
200  err = err + 2
201  else
202  write(mesg, '(ES13.5)') sum
203  call mpp_error(fatal,"Repro_sum_2d: Overflow in real_to_ints conversion of "//trim(mesg))
204  endif
205  endif
206  endif
207  endif
208 
209  if (present(efp_sum)) efp_sum%v(:) = ints_sum(:)
210 
211  if (debug) then
212  write(mesg,'("2d RS: ", ES24.16, 6 Z17.16)') sum, ints_sum(1:numint)
213  if(mpp_pe() == mpp_root_pe()) call mpp_error(note, mesg)
214  endif
215 
216 end function mpp_reproducing_sum_r8_2d
217 
218 function mpp_reproducing_sum_r4_2d(array, isr, ier, jsr, jer, EFP_sum, reproducing, &
219  overflow_check, err) result(sum)
220  real(FLOAT_KIND), dimension(:,:), intent(in) :: array
221  integer, optional, intent(in) :: isr, ier, jsr, jer
222  type(mpp_efp_type), optional, intent(out) :: efp_sum
223  logical, optional, intent(in) :: reproducing
224  logical, optional, intent(in) :: overflow_check
225  integer, optional, intent(out) :: err
226  real(FLOAT_KIND) :: sum ! Result
227 
228  real(DOUBLE_KIND) :: array_r8(size(array,1), size(array,2))
229 
230  array_r8 = array
231 
232  sum = mpp_reproducing_sum_r8_2d(array_r8, isr, ier, jsr, jer, efp_sum, reproducing, &
233  overflow_check, err)
234 
235  return
236 
237 end function mpp_reproducing_sum_r4_2d
238 
239 
240 function mpp_reproducing_sum_r8_3d(array, isr, ier, jsr, jer, sums, EFP_sum, err) &
241  result(sum)
242  real(DOUBLE_KIND), dimension(:,:,:), intent(in) :: array
243  integer, optional, intent(in) :: isr, ier, jsr, jer
244  real(DOUBLE_KIND), dimension(:), optional, intent(out) :: sums
245  type(mpp_efp_type), optional, intent(out) :: efp_sum
246  integer, optional, intent(out) :: err
247  real(DOUBLE_KIND) :: sum ! Result
248 
249  ! This subroutine uses a conversion to an integer representation
250  ! of real numbers to give order-invariant sums that will reproduce
251  ! across PE count. This idea comes from R. Hallberg and A. Adcroft.
252 
253  real(DOUBLE_KIND) :: max_mag_term
254  integer(LONG_KIND), dimension(NUMINT) :: ints_sum
255  integer(LONG_KIND), dimension(NUMINT,size(array,3)) :: ints_sums
256  integer(LONG_KIND) :: prec_error
257  character(len=256) :: mesg
258  integer :: i, j, k, is, ie, js, je, ke, isz, jsz, n
259 
260  if (mpp_npes() > max_count_prec) call mpp_error(fatal, &
261  "mpp_reproducing_sum: Too many processors are being used for the value of "//&
262  "prec. Reduce prec to (2^63-1)/mpp_npes.")
263 
264  prec_error = (2_8**62 + (2_8**62 - 1)) / mpp_npes()
265  max_mag_term = 0.0
266 
267  is = 1 ; ie = size(array,1) ; js = 1 ; je = size(array,2) ; ke = size(array,3)
268  if (present(isr)) then
269  if (isr < is) call mpp_error(fatal, &
270  "Value of isr too small in mpp_reproducing_sum(_3d).")
271  is = isr
272  endif
273  if (present(ier)) then
274  if (ier > ie) call mpp_error(fatal, &
275  "Value of ier too large in mpp_reproducing_sum(_3d).")
276  ie = ier
277  endif
278  if (present(jsr)) then
279  if (jsr < js) call mpp_error(fatal, &
280  "Value of jsr too small in mpp_reproducing_sum(_3d).")
281  js = jsr
282  endif
283  if (present(jer)) then
284  if (jer > je) call mpp_error(fatal, &
285  "Value of jer too large in mpp_reproducing_sum(_3d).")
286  je = jer
287  endif
288  jsz = je+1-js; isz = ie+1-is
289 
290  if (present(sums)) then
291  if (size(sums) > ke) call mpp_error(fatal, "Sums is smaller than "//&
292  "the vertical extent of array in mpp_reproducing_sum(_3d).")
293  ints_sums(:,:) = 0
294  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
295  if (jsz*isz < max_count_prec) then
296  do k=1,ke
297  do j=js,je ; do i=is,ie
298  call increment_ints_faster(ints_sums(:,k), array(i,j,k), max_mag_term);
299  enddo ; enddo
300  call carry_overflow(ints_sums(:,k), prec_error)
301  enddo
302  elseif (isz < max_count_prec) then
303  do k=1,ke ; do j=js,je
304  do i=is,ie
305  call increment_ints_faster(ints_sums(:,k), array(i,j,k), max_mag_term);
306  enddo
307  call carry_overflow(ints_sums(:,k), prec_error)
308  enddo ; enddo
309  else
310  do k=1,ke ; do j=js,je ; do i=is,ie
311  call increment_ints(ints_sums(:,k), &
312  real_to_ints(array(i,j,k), prec_error), prec_error);
313  enddo ; enddo ; enddo
314  endif
315  if (present(err)) then
316  err = 0
317  if (abs(max_mag_term) >= prec_error*pr(1)) err = err+1
318  if (overflow_error) err = err+2
319  if (nan_error) err = err+2
320  if (err > 0) then ; do k=1,ke ; do n=1,numint ; ints_sums(n,k) = 0 ; enddo ; enddo ; endif
321  else
322  if (nan_error) call mpp_error(fatal, &
323  "NaN in input field of mpp_reproducing_sum(_3d), this indicates numerical instability")
324  if (abs(max_mag_term) >= prec_error*pr(1)) then
325  write(mesg, '(ES13.5)') max_mag_term
326  call mpp_error(fatal,"Overflow in mpp_reproducing_sum(_3d) conversion of "//trim(mesg))
327  endif
328  if (overflow_error) call mpp_error(fatal, "Overflow in mpp_reproducing_sum(_3d).")
329  endif
330 
331  call mpp_sum(ints_sums(:,1:ke), numint*ke)
332 
333  sum = 0.0
334  do k=1,ke
335  call regularize_ints(ints_sums(:,k))
336  sums(k) = ints_to_real(ints_sums(:,k))
337  sum = sum + sums(k)
338  enddo
339 
340  if (present(efp_sum)) then
341  efp_sum%v(:) = 0
342  do k=1,ke ; call increment_ints(efp_sum%v(:), ints_sums(:,k)) ; enddo
343  endif
344 
345  if (debug) then
346  do n=1,numint ; ints_sum(n) = 0 ; enddo
347  do k=1,ke ; do n=1,numint ; ints_sum(n) = ints_sum(n) + ints_sums(n,k) ; enddo ; enddo
348  write(mesg,'("3D RS: ", ES24.16, 6 Z17.16)') sum, ints_sum(1:numint)
349  if(mpp_pe()==mpp_root_pe()) call mpp_error(note, mesg)
350  endif
351  else
352  ints_sum(:) = 0
353  overflow_error = .false. ; nan_error = .false. ; max_mag_term = 0.0
354  if (jsz*isz < max_count_prec) then
355  do k=1,ke
356  do j=js,je ; do i=is,ie
357  call increment_ints_faster(ints_sum, array(i,j,k), max_mag_term);
358  enddo ; enddo
359  call carry_overflow(ints_sum, prec_error)
360  enddo
361  elseif (isz < max_count_prec) then
362  do k=1,ke ; do j=js,je
363  do i=is,ie
364  call increment_ints_faster(ints_sum, array(i,j,k), max_mag_term);
365  enddo
366  call carry_overflow(ints_sum, prec_error)
367  enddo ; enddo
368  else
369  do k=1,ke ; do j=js,je ; do i=is,ie
370  call increment_ints(ints_sum, real_to_ints(array(i,j,k), prec_error), &
371  prec_error);
372  enddo ; enddo ; enddo
373  endif
374  if (present(err)) then
375  err = 0
376  if (abs(max_mag_term) >= prec_error*pr(1)) err = err+1
377  if (overflow_error) err = err+2
378  if (nan_error) err = err+2
379  if (err > 0) then ; do n=1,numint ; ints_sum(n) = 0 ; enddo ; endif
380  else
381  if (nan_error) call mpp_error(fatal, &
382  "NaN in input field of mpp_reproducing_sum(_3d), this indicates numerical instability")
383  if (abs(max_mag_term) >= prec_error*pr(1)) then
384  write(mesg, '(ES13.5)') max_mag_term
385  call mpp_error(fatal,"Overflow in mpp_reproducing_sum(_3d) conversion of "//trim(mesg))
386  endif
387  if (overflow_error) call mpp_error(fatal, "Overflow in mpp_reproducing_sum(_3d).")
388  endif
389 
390  call mpp_sum(ints_sum, numint)
391 
392  call regularize_ints(ints_sum)
393  sum = ints_to_real(ints_sum)
394 
395  if (present(efp_sum)) efp_sum%v(:) = ints_sum(:)
396 
397  if (debug) then
398  write(mesg,'("3d RS: ", ES24.16, 6 Z17.16)') sum, ints_sum(1:numint)
399  if(mpp_pe()==mpp_root_pe()) call mpp_error(note, mesg)
400  endif
401  endif
402 
403 end function mpp_reproducing_sum_r8_3d
404 
405 function real_to_ints(r, prec_error, overflow) result(ints)
406  real(DOUBLE_KIND), intent(in) :: r
407  integer(LONG_KIND), optional, intent(in) :: prec_error
408  logical, optional, intent(inout) :: overflow
409  integer(LONG_KIND), dimension(NUMINT) :: ints
410  ! This subroutine converts a real number to an equivalent representation
411  ! using several long integers.
412 
413  real(DOUBLE_KIND) :: rs
414  character(len=80) :: mesg
415  integer(LONG_KIND) :: ival, prec_err
416  integer :: sgn, i
417 
418  prec_err = prec ; if (present(prec_error)) prec_err = prec_error
419  ints(:) = 0_8
420  if ((r >= 1e30) .eqv. (r < 1e30)) then ; nan_error = .true. ; return ; endif
421 
422  sgn = 1 ; if (r<0.0) sgn = -1
423  rs = abs(r)
424 
425  if (present(overflow)) then
426  if (.not.(rs < prec_err*pr(1))) overflow = .true.
427  if ((r >= 1e30) .eqv. (r < 1e30)) overflow = .true.
428  elseif (.not.(rs < prec_err*pr(1))) then
429  write(mesg, '(ES13.5)') r
430  call mpp_error(fatal,"Overflow in real_to_ints conversion of "//trim(mesg))
431  endif
432 
433  do i=1,numint
434  ival = int(rs*i_pr(i), 8)
435  rs = rs - ival*pr(i)
436  ints(i) = sgn*ival
437  enddo
438 
439 end function real_to_ints
440 
441 function ints_to_real(ints) result(r)
442  integer(LONG_KIND), dimension(NUMINT), intent(in) :: ints
443  real(DOUBLE_KIND) :: r
444  ! This subroutine reverses the conversion in real_to_ints.
445 
446  integer :: i
447 
448  r = 0.0
449  do i=1,numint ; r = r + pr(i)*ints(i) ; enddo
450 end function ints_to_real
451 
452 subroutine increment_ints(int_sum, int2, prec_error)
453  integer(LONG_KIND), dimension(NUMINT), intent(inout) :: int_sum
454  integer(LONG_KIND), dimension(NUMINT), intent(in) :: int2
455  integer(LONG_KIND), optional, intent(in) :: prec_error
456 
457  ! This subroutine increments a number with another, both using the integer
458  ! representation in real_to_ints.
459  integer :: i
460 
461  do i=numint,2,-1
462  int_sum(i) = int_sum(i) + int2(i)
463  ! Carry the local overflow.
464  if (int_sum(i) > prec) then
465  int_sum(i) = int_sum(i) - prec
466  int_sum(i-1) = int_sum(i-1) + 1
467  elseif (int_sum(i) < -prec) then
468  int_sum(i) = int_sum(i) + prec
469  int_sum(i-1) = int_sum(i-1) - 1
470  endif
471  enddo
472  int_sum(1) = int_sum(1) + int2(1)
473  if (present(prec_error)) then
474  if (abs(int_sum(1)) > prec_error) overflow_error = .true.
475  else
476  if (abs(int_sum(1)) > prec) overflow_error = .true.
477  endif
478 
479 end subroutine increment_ints
480 
481 subroutine increment_ints_faster(int_sum, r, max_mag_term)
482  integer(LONG_KIND), dimension(NUMINT), intent(inout) :: int_sum
483  real(DOUBLE_KIND), intent(in) :: r
484  real(DOUBLE_KIND), intent(inout) :: max_mag_term
485 
486  ! This subroutine increments a number with another, both using the integer
487  ! representation in real_to_ints, but without doing any carrying of overflow.
488  ! The entire operation is embedded in a single call for greater speed.
489  real(DOUBLE_KIND) :: rs
490  integer(LONG_KIND) :: ival
491  integer :: sgn, i
492 
493  if ((r >= 1e30) .eqv. (r < 1e30)) then ; nan_error = .true. ; return ; endif
494  sgn = 1 ; if (r<0.0) sgn = -1
495  rs = abs(r)
496  if (rs > abs(max_mag_term)) max_mag_term = r
497 
498  do i=1,numint
499  ival = int(rs*i_pr(i), 8)
500  rs = rs - ival*pr(i)
501  int_sum(i) = int_sum(i) + sgn*ival
502  enddo
503 
504 end subroutine increment_ints_faster
505 
506 subroutine carry_overflow(int_sum, prec_error)
507  integer(LONG_KIND), dimension(NUMINT), intent(inout) :: int_sum
508  integer(LONG_KIND), intent(in) :: prec_error
509 
510  ! This subroutine handles carrying of the overflow.
511  integer :: i, num_carry
512 
513  do i=numint,2,-1 ; if (abs(int_sum(i)) > prec) then
514  num_carry = int(int_sum(i) * i_prec)
515  int_sum(i) = int_sum(i) - num_carry*prec
516  int_sum(i-1) = int_sum(i-1) + num_carry
517  endif ; enddo
518  if (abs(int_sum(1)) > prec_error) then
519  overflow_error = .true.
520  endif
521 
522 end subroutine carry_overflow
523 
524 subroutine regularize_ints(int_sum)
525  integer(LONG_KIND), dimension(NUMINT), intent(inout) :: int_sum
526 
527  ! This subroutine carries the overflow, and then makes sure that
528  ! all integers are of the same sign as the overall value.
529  logical :: positive
530  integer :: i, num_carry
531 
532  do i=numint,2,-1 ; if (abs(int_sum(i)) > prec) then
533  num_carry = int(int_sum(i) * i_prec)
534  int_sum(i) = int_sum(i) - num_carry*prec
535  int_sum(i-1) = int_sum(i-1) + num_carry
536  endif ; enddo
537 
538  ! Determine the sign of the final number.
539  positive = .true.
540  do i=1,numint
541  if (abs(int_sum(i)) > 0) then
542  if (int_sum(i) < 0) positive = .false.
543  exit
544  endif
545  enddo
546 
547  if (positive) then
548  do i=numint,2,-1 ; if (int_sum(i) < 0) then
549  int_sum(i) = int_sum(i) + prec
550  int_sum(i-1) = int_sum(i-1) - 1
551  endif ; enddo
552  else
553  do i=numint,2,-1 ; if (int_sum(i) > 0) then
554  int_sum(i) = int_sum(i) - prec
555  int_sum(i-1) = int_sum(i-1) + 1
556  endif ; enddo
557  endif
558 
559 end subroutine regularize_ints
560 
564 end function mpp_query_efp_overflow_error
565 
566 subroutine mpp_reset_efp_overlow_error()
567  overflow_error = .false.
568 end subroutine mpp_reset_efp_overlow_error
569 
570 function mpp_efp_plus(EFP1, EFP2)
571  type(mpp_efp_type) :: mpp_efp_plus
572  type(mpp_efp_type), intent(in) :: efp1, efp2
573 
574  mpp_efp_plus = efp1
575 
576  call increment_ints(mpp_efp_plus%v(:), efp2%v(:))
577 end function mpp_efp_plus
578 
579 function mpp_efp_minus(EFP1, EFP2)
580  type(mpp_efp_type) :: mpp_efp_minus
581  type(mpp_efp_type), intent(in) :: efp1, efp2
582  integer :: i
583 
584  do i=1,numint ; mpp_efp_minus%v(i) = -1*efp2%v(i) ; enddo
585 
586  call increment_ints(mpp_efp_minus%v(:), efp1%v(:))
587 end function mpp_efp_minus
588 
589 subroutine mpp_efp_assign(EFP1, EFP2)
590  type(mpp_efp_type), intent(out) :: EFP1
591  type(mpp_efp_type), intent(in) :: EFP2
592  integer i
593  ! This subroutine assigns all components of the extended fixed point type
594  ! variable on the RHS (EFP2) to the components of the variable on the LHS
595  ! (EFP1).
596 
597  do i=1,numint ; efp1%v(i) = efp2%v(i) ; enddo
598 end subroutine mpp_efp_assign
599 
600 function mpp_efp_to_real(EFP1)
601  type(mpp_efp_type), intent(inout) :: efp1
602  real(DOUBLE_KIND) :: mpp_efp_to_real
603 
604  call regularize_ints(efp1%v)
605  mpp_efp_to_real = ints_to_real(efp1%v)
606 end function mpp_efp_to_real
607 
608 function mpp_efp_real_diff(EFP1, EFP2)
609  type(mpp_efp_type), intent(in) :: efp1, efp2
610  real(DOUBLE_KIND) :: mpp_efp_real_diff
611 
612  type(mpp_efp_type) :: efp_diff
613 
614  efp_diff = efp1 - efp2
616 
617 end function mpp_efp_real_diff
618 
619 function mpp_real_to_efp(val, overflow)
620  real(DOUBLE_KIND), intent(in) :: val
621  logical, optional, intent(inout) :: overflow
623 
624  logical :: over
625  character(len=80) :: mesg
626 
627  if (present(overflow)) then
628  mpp_real_to_efp%v(:) = real_to_ints(val, overflow=overflow)
629  else
630  over = .false.
631  mpp_real_to_efp%v(:) = real_to_ints(val, overflow=over)
632  if (over) then
633  write(mesg, '(ES13.5)') val
634  call mpp_error(fatal,"Overflow in mpp_real_to_efp conversion of "//trim(mesg))
635  endif
636  endif
637 
638 end function mpp_real_to_efp
639 
640 subroutine mpp_efp_list_sum_across_pes(EFPs, nval, errors)
641  type(mpp_efp_type), dimension(:), intent(inout) :: efps
642  integer, intent(in) :: nval
643  logical, dimension(:), optional, intent(out) :: errors
644 
645  ! This subroutine does a sum across PEs of a list of EFP variables,
646  ! returning the sums in place, with all overflows carried.
647 
648  integer(LONG_KIND), dimension(NUMINT,nval) :: ints
649  integer(LONG_KIND) :: prec_error
650  logical :: error_found
651  character(len=256) :: mesg
652  integer :: i, n
653 
654  if (mpp_npes() > max_count_prec) call mpp_error(fatal, &
655  "mpp_efp_list_sum_across_PEs: Too many processors are being used for the value of "//&
656  "prec. Reduce prec to (2^63-1)/mpp_npes.")
657 
658  prec_error = (2_8**62 + (2_8**62 - 1)) / mpp_npes()
659  ! overflow_error is an overflow error flag for the whole module.
660  overflow_error = .false. ; error_found = .false.
661 
662  do i=1,nval ; do n=1,numint ; ints(n,i) = efps(i)%v(n) ; enddo ; enddo
663 
664  call mpp_sum(ints(:,:), numint*nval)
665 
666  if (present(errors)) errors(:) = .false.
667  do i=1,nval
668  overflow_error = .false.
669  call carry_overflow(ints(:,i), prec_error)
670  do n=1,numint ; efps(i)%v(n) = ints(n,i) ; enddo
671  if (present(errors)) errors(i) = overflow_error
672  if (overflow_error) then
673  write (mesg,'("mpp_efp_list_sum_across_PEs error at ",i6," val was ",ES12.6, ", prec_error = ",ES12.6)') &
674  i, mpp_efp_to_real(efps(i)), real(prec_error)
675  if(mpp_pe()==mpp_root_pe()) call mpp_error(warning, mesg)
676  endif
677  error_found = error_found .or. overflow_error
678  enddo
679  if (error_found .and. .not.(present(errors))) then
680  call mpp_error(fatal, "Overflow in mpp_efp_list_sum_across_PEs.")
681  endif
682 
683 end subroutine mpp_efp_list_sum_across_pes
684 
685 end module mpp_efp_mod
logical function, public mpp_query_efp_overflow_error()
Definition: mpp_efp.F90:562
type(mpp_efp_type) function, public mpp_efp_plus(EFP1, EFP2)
Definition: mpp_efp.F90:571
type(mpp_efp_type) function, public mpp_efp_minus(EFP1, EFP2)
Definition: mpp_efp.F90:580
real(float_kind) function mpp_reproducing_sum_r4_2d(array, isr, ier, jsr, jer, EFP_sum, reproducing, overflow_check, err)
Definition: mpp_efp.F90:220
subroutine, public mpp_reset_efp_overlow_error()
Definition: mpp_efp.F90:567
real(double_kind), dimension(numint), parameter pr
Definition: mpp_efp.F90:47
integer, parameter numbit
Definition: mpp_efp.F90:35
type(mpp_efp_type) function, public mpp_real_to_efp(val, overflow)
Definition: mpp_efp.F90:620
Definition: mpp.F90:39
logical debug
Definition: mpp_efp.F90:53
subroutine regularize_ints(int_sum)
Definition: mpp_efp.F90:525
real(double_kind), parameter r_prec
Definition: mpp_efp.F90:40
real(double_kind) function mpp_reproducing_sum_r8_3d(array, isr, ier, jsr, jer, sums, EFP_sum, err)
Definition: mpp_efp.F90:242
integer(long_kind) function, dimension(numint) real_to_ints(r, prec_error, overflow)
Definition: mpp_efp.F90:406
logical nan_error
Definition: mpp_efp.F90:52
real(double_kind) function, public mpp_efp_real_diff(EFP1, EFP2)
Definition: mpp_efp.F90:609
subroutine increment_ints(int_sum, int2, prec_error)
Definition: mpp_efp.F90:453
real(double_kind), parameter i_prec
Definition: mpp_efp.F90:41
integer(long_kind), parameter prec
Definition: mpp_efp.F90:39
subroutine carry_overflow(int_sum, prec_error)
Definition: mpp_efp.F90:507
integer, parameter numint
Definition: mpp_efp.F90:36
subroutine increment_ints_faster(int_sum, r, max_mag_term)
Definition: mpp_efp.F90:482
integer, parameter max_count_prec
Definition: mpp_efp.F90:42
subroutine, public mpp_efp_list_sum_across_pes(EFPs, nval, errors)
Definition: mpp_efp.F90:641
real(double_kind) function ints_to_real(ints)
Definition: mpp_efp.F90:442
logical overflow_error
Definition: mpp_efp.F90:52
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
real(double_kind), dimension(numint), parameter i_pr
Definition: mpp_efp.F90:49
real(double_kind) function, public mpp_efp_to_real(EFP1)
Definition: mpp_efp.F90:601
real(double_kind) function mpp_reproducing_sum_r8_2d(array, isr, ier, jsr, jer, EFP_sum, reproducing, overflow_check, err)
Definition: mpp_efp.F90:75