20 #include <fms_platform.h> 23 use mpp_mod,
only : mpp_pe, mpp_root_pe, mpp_npes
26 implicit none ;
private 30 public ::
operator(+),
operator(-),
assignment(=)
47 real(DOUBLE_KIND),
parameter,
dimension(NUMINT) :: &
49 real(DOUBLE_KIND),
parameter,
dimension(NUMINT) :: &
64 integer(kind=8),
dimension(NUMINT) :: v
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
74 overflow_check, err)
result(sum)
75 real(DOUBLE_KIND),
dimension(:,:),
intent(in) :: array
76 integer,
optional,
intent(in) :: isr, ier, jsr, jer
78 logical,
optional,
intent(in) :: reproducing
79 logical,
optional,
intent(in) :: overflow_check
80 integer,
optional,
intent(out) :: err
81 real(DOUBLE_KIND) :: sum
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
96 "mpp_reproducing_sum: Too many processors are being used for the value of "//&
97 "prec. Reduce prec to (2^63-1)/mpp_npes.")
99 prec_error = (2_8**62 + (2_8**62 - 1)) / mpp_npes()
101 is = 1 ; ie =
size(array,1) ; js = 1 ; je =
size(array,2 )
102 if (
present(isr))
then 104 "Value of isr too small in mpp_reproducing_sum_2d.")
107 if (
present(ier))
then 109 "Value of ier too large in mpp_reproducing_sum_2d.")
112 if (
present(jsr))
then 114 "Value of jsr too small in mpp_reproducing_sum_2d.")
117 if (
present(jer))
then 119 "Value of jer too large in mpp_reproducing_sum_2d.")
123 repro = .true. ;
if (
present(reproducing)) repro = reproducing
124 over_check = .true. ;
if (
present(overflow_check)) over_check = overflow_check
131 do j=js,je ;
do i=is,ie
143 do j=js,je ;
do i=is,ie
149 do j=js,je ;
do i=is,ie
150 sgn = 1 ;
if (array(i,j)<0.0) sgn = -1
153 ival = int(rs*
i_pr(n), 8)
155 ints_sum(n) = ints_sum(n) + sgn*ival
161 if (
present(err))
then 167 if (err > 0)
then ;
do n=1,
numint ; ints_sum(n) = 0 ;
enddo ;
endif 170 call mpp_error(fatal,
"NaN in input field of mpp_reproducing_sum(_2d), this indicates numerical instability")
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))
177 call mpp_error(fatal,
"Overflow in mpp_reproducing_sum(_2d).")
187 do j=js,je ;
do i=is,ie
188 rsum(1) = rsum(1) + array(i,j);
193 if (
present(err))
then ; err = 0 ;
endif 195 if (
debug .or.
present(efp_sum))
then 199 if (
present(err))
then 202 write(mesg,
'(ES13.5)') sum
203 call mpp_error(fatal,
"Repro_sum_2d: Overflow in real_to_ints conversion of "//trim(mesg))
209 if (
present(efp_sum)) efp_sum%v(:) = ints_sum(:)
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)
219 overflow_check, err)
result(sum)
220 real(FLOAT_KIND),
dimension(:,:),
intent(in) :: array
221 integer,
optional,
intent(in) :: isr, ier, jsr, jer
223 logical,
optional,
intent(in) :: reproducing
224 logical,
optional,
intent(in) :: overflow_check
225 integer,
optional,
intent(out) :: err
226 real(FLOAT_KIND) :: sum
228 real(DOUBLE_KIND) :: array_r8(
size(array,1),
size(array,2))
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
246 integer,
optional,
intent(out) :: err
247 real(DOUBLE_KIND) :: sum
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
261 "mpp_reproducing_sum: Too many processors are being used for the value of "//&
262 "prec. Reduce prec to (2^63-1)/mpp_npes.")
264 prec_error = (2_8**62 + (2_8**62 - 1)) / mpp_npes()
267 is = 1 ; ie =
size(array,1) ; js = 1 ; je =
size(array,2) ; ke =
size(array,3)
268 if (
present(isr))
then 270 "Value of isr too small in mpp_reproducing_sum(_3d).")
273 if (
present(ier))
then 275 "Value of ier too large in mpp_reproducing_sum(_3d).")
278 if (
present(jsr))
then 280 "Value of jsr too small in mpp_reproducing_sum(_3d).")
283 if (
present(jer))
then 285 "Value of jer too large in mpp_reproducing_sum(_3d).")
288 jsz = je+1-js; isz = ie+1-is
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).")
297 do j=js,je ;
do i=is,ie
303 do k=1,ke ;
do j=js,je
310 do k=1,ke ;
do j=js,je ;
do i=is,ie
313 enddo ;
enddo ;
enddo 315 if (
present(err))
then 317 if (abs(max_mag_term) >= prec_error*
pr(1)) err = err+1
320 if (err > 0)
then ;
do k=1,ke ;
do n=1,
numint ; ints_sums(n,k) = 0 ;
enddo ;
enddo ;
endif 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))
340 if (
present(efp_sum))
then 342 do k=1,ke ;
call increment_ints(efp_sum%v(:), ints_sums(:,k)) ;
enddo 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)
356 do j=js,je ;
do i=is,ie
362 do k=1,ke ;
do j=js,je
369 do k=1,ke ;
do j=js,je ;
do i=is,ie
372 enddo ;
enddo ;
enddo 374 if (
present(err))
then 376 if (abs(max_mag_term) >= prec_error*
pr(1)) err = err+1
379 if (err > 0)
then ;
do n=1,
numint ; ints_sum(n) = 0 ;
enddo ;
endif 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))
395 if (
present(efp_sum)) efp_sum%v(:) = ints_sum(:)
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)
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
413 real(DOUBLE_KIND) :: rs
414 character(len=80) :: mesg
415 integer(LONG_KIND) :: ival, prec_err
418 prec_err =
prec ;
if (
present(prec_error)) prec_err = prec_error
420 if ((r >= 1e30) .eqv. (r < 1e30))
then ;
nan_error = .true. ;
return ;
endif 422 sgn = 1 ;
if (r<0.0) sgn = -1
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))
434 ival = int(rs*
i_pr(i), 8)
442 integer(LONG_KIND),
dimension(NUMINT),
intent(in) :: ints
443 real(DOUBLE_KIND) :: r
449 do i=1,
numint ; r = r +
pr(i)*ints(i) ;
enddo 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
462 int_sum(i) = int_sum(i) + int2(i)
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
472 int_sum(1) = int_sum(1) + int2(1)
473 if (
present(prec_error))
then 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
489 real(DOUBLE_KIND) :: rs
490 integer(LONG_KIND) :: ival
493 if ((r >= 1e30) .eqv. (r < 1e30))
then ;
nan_error = .true. ;
return ;
endif 494 sgn = 1 ;
if (r<0.0) sgn = -1
496 if (rs > abs(max_mag_term)) max_mag_term = r
499 ival = int(rs*
i_pr(i), 8)
501 int_sum(i) = int_sum(i) + sgn*ival
507 integer(LONG_KIND),
dimension(NUMINT),
intent(inout) :: int_sum
508 integer(LONG_KIND),
intent(in) :: prec_error
511 integer :: i, num_carry
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
518 if (abs(int_sum(1)) > prec_error)
then 525 integer(LONG_KIND),
dimension(NUMINT),
intent(inout) :: int_sum
530 integer :: i, num_carry
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
541 if (abs(int_sum(i)) > 0)
then 542 if (int_sum(i) < 0) positive = .false.
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
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
570 function mpp_efp_plus(EFP1, EFP2)
577 end function mpp_efp_plus
579 function mpp_efp_minus(EFP1, EFP2)
584 do i=1,
numint ; mpp_efp_minus%v(i) = -1*efp2%v(i) ;
enddo 587 end function mpp_efp_minus
589 subroutine mpp_efp_assign(EFP1, EFP2)
590 type(mpp_efp_type),
intent(out) :: EFP1
591 type(mpp_efp_type),
intent(in) :: EFP2
597 do i=1,
numint ; efp1%v(i) = efp2%v(i) ;
enddo 598 end subroutine mpp_efp_assign
614 efp_diff = efp1 - efp2
620 real(DOUBLE_KIND),
intent(in) :: val
621 logical,
optional,
intent(inout) :: overflow
625 character(len=80) :: mesg
627 if (
present(overflow))
then 633 write(mesg,
'(ES13.5)') val
634 call mpp_error(fatal,
"Overflow in mpp_real_to_efp conversion of "//trim(mesg))
642 integer,
intent(in) :: nval
643 logical,
dimension(:),
optional,
intent(out) :: errors
648 integer(LONG_KIND),
dimension(NUMINT,nval) :: ints
649 integer(LONG_KIND) :: prec_error
650 logical :: error_found
651 character(len=256) :: mesg
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.")
658 prec_error = (2_8**62 + (2_8**62 - 1)) / mpp_npes()
662 do i=1,nval ;
do n=1,
numint ; ints(n,i) = efps(i)%v(n) ;
enddo ;
enddo 666 if (
present(errors)) errors(:) = .false.
670 do n=1,
numint ; efps(i)%v(n) = ints(n,i) ;
enddo 673 write (mesg,
'("mpp_efp_list_sum_across_PEs error at ",i6," val was ",ES12.6, ", prec_error = ",ES12.6)') &
675 if(mpp_pe()==mpp_root_pe())
call mpp_error(warning, mesg)
679 if (error_found .and. .not.(
present(errors)))
then 680 call mpp_error(fatal,
"Overflow in mpp_efp_list_sum_across_PEs.")
logical function, public mpp_query_efp_overflow_error()
type(mpp_efp_type) function, public mpp_efp_plus(EFP1, EFP2)
type(mpp_efp_type) function, public mpp_efp_minus(EFP1, EFP2)
real(float_kind) function mpp_reproducing_sum_r4_2d(array, isr, ier, jsr, jer, EFP_sum, reproducing, overflow_check, err)
subroutine, public mpp_reset_efp_overlow_error()
real(double_kind), dimension(numint), parameter pr
integer, parameter numbit
type(mpp_efp_type) function, public mpp_real_to_efp(val, overflow)
subroutine regularize_ints(int_sum)
real(double_kind), parameter r_prec
real(double_kind) function mpp_reproducing_sum_r8_3d(array, isr, ier, jsr, jer, sums, EFP_sum, err)
integer(long_kind) function, dimension(numint) real_to_ints(r, prec_error, overflow)
real(double_kind) function, public mpp_efp_real_diff(EFP1, EFP2)
subroutine increment_ints(int_sum, int2, prec_error)
real(double_kind), parameter i_prec
integer(long_kind), parameter prec
subroutine carry_overflow(int_sum, prec_error)
integer, parameter numint
subroutine increment_ints_faster(int_sum, r, max_mag_term)
integer, parameter max_count_prec
subroutine, public mpp_efp_list_sum_across_pes(EFPs, nval, errors)
real(double_kind) function ints_to_real(ints)
************************************************************************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
real(double_kind) function, public mpp_efp_to_real(EFP1)
real(double_kind) function mpp_reproducing_sum_r8_2d(array, isr, ier, jsr, jer, EFP_sum, reproducing, overflow_check, err)