FMS  2024.03
Flexible Modeling System
mpp_utilities.F90
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 !***********************************************************************
19 !> @defgroup mpp_utilities_mod mpp_utilities_mod
20 !> @ingroup mpp
21 !> @brief Module for utiltity routines to be used in @ref mpp modules
22 !!
23 !> Currently only holds one routine for finding global min and max
24 
25 !> @addtogroup mpp_utilities_mod
26 !> @{
27 module mpp_utilities_mod
28 
29 implicit none
30 !-----------------------------------------------------------------------
31 ! Include variable "version" to be written to log file.
32 #include<file_version.h>
33 !-----------------------------------------------------------------------
34 
35  public :: mpp_array_global_min_max
36 
37 contains
38 
39 !#######################################################################
40 !> @brief Compute and return the global min and max of an array
41 !! and the corresponding lat-lon-depth locations .
42 !!
43 !> This algorithm works only for an input array that has a unique global
44 !! max and min location. This is assured by introducing a factor that distinguishes
45 !! the values of extrema at each processor.
46 !!
47 !! Vectorized using maxloc() and minloc() intrinsic functions by
48 !! Russell.Fiedler@csiro.au (May 2005).
49 !!
50 !! Modified by Zhi.Liang@noaa.gov (July 2005)
51 !!
52 !! Modified by Niki.Zadeh@noaa.gov (Feb. 2009)
53 !!
54 subroutine mpp_array_global_min_max(in_array, tmask,isd,jsd,isc,iec,jsc,jec,nk, g_min, g_max, &
55  geo_x,geo_y,geo_z, xgmin, ygmin, zgmin, xgmax, ygmax, zgmax)
56 
57  use mpp_mod, only: mpp_min, mpp_max, mpp_pe, mpp_sum
58 
59  integer, intent(in) :: isd,jsd,isc,iec,jsc,jec,nk
60  real, dimension(isd:,jsd:,:), intent(in) :: in_array
61  real, dimension(isd:,jsd:,:), intent(in) :: tmask
62  real, intent(out):: g_min, g_max
63  real, dimension(isd:,jsd:), intent(in) :: geo_x,geo_y
64  real, dimension(:), intent(in) :: geo_z
65  real, intent(out):: xgmin, ygmin, zgmin, xgmax, ygmax, zgmax
66 
67  real :: tmax, tmin, tmax0, tmin0
68  integer :: itmax, jtmax, ktmax, itmin, jtmin, ktmin
69  real :: fudge
70 
71  ! arrays to enable vectorization
72  integer :: iminarr(3),imaxarr(3)
73 
74  g_min=-88888888888.0 ; g_max=-999999999.0
75 
76  tmax=-1.e10;tmin=1.e10
77  itmax=0;jtmax=0;ktmax=0
78  itmin=0;jtmin=0;ktmin=0
79 
80  if(any(tmask(isc:iec,jsc:jec,:) > 0.)) then
81  iminarr=minloc(in_array(isc:iec,jsc:jec,:),tmask(isc:iec,jsc:jec,:) > 0.)
82  imaxarr=maxloc(in_array(isc:iec,jsc:jec,:),tmask(isc:iec,jsc:jec,:) > 0.)
83  itmin=iminarr(1)+isc-1
84  jtmin=iminarr(2)+jsc-1
85  ktmin=iminarr(3)
86  itmax=imaxarr(1)+isc-1
87  jtmax=imaxarr(2)+jsc-1
88  ktmax=imaxarr(3)
89  tmin=in_array(itmin,jtmin,ktmin)
90  tmax=in_array(itmax,jtmax,ktmax)
91  end if
92 
93  ! use "fudge" to distinguish processors when tracer extreme is independent of processor
94  fudge = 1.0 + 1.e-12*real(mpp_pe() )
95  tmax = tmax*fudge
96  tmin = tmin*fudge
97  if(tmax == 0.0) then
98  tmax = tmax + 1.e-12*real(mpp_pe() )
99  endif
100  if(tmin == 0.0) then
101  tmin = tmin + 1.e-12*real(mpp_pe() )
102  endif
103 
104 
105  tmax0=tmax;tmin0=tmin
106 
107  call mpp_max(tmax)
108  call mpp_min(tmin)
109 
110  g_max = tmax
111  g_min = tmin
112 
113  !Now find the location of the global extrema.
114  !
115  !Note that the fudge factor above guarantees that the location of max (min) is uinque,
116  ! since tmax0 (tmin0) has slightly different values on each processor.
117  !Otherwise, the function in_array(i,j,k) could be equal to global max (min) at more
118  ! than one point in space and this would be a much more difficult problem to solve.
119  !
120  !mpp_max trick
121  !-999 on all current PE's
122  xgmax=-999.; ygmax=-999.; zgmax=-999.
123  xgmin=-999.; ygmin=-999.; zgmin=-999.
124 
125 
126  !except when
127  if (tmax0 == tmax) then !This happens ONLY on ONE processor because of fudge factor above.
128  xgmax=geo_x(itmax,jtmax)
129  ygmax=geo_y(itmax,jtmax)
130  zgmax=geo_z(ktmax)
131  endif
132 
133  call mpp_max(xgmax)
134  call mpp_max(ygmax)
135  call mpp_max(zgmax)
136 
137  if (tmin0 == tmin) then !This happens ONLY on ONE processor because of fudge factor above.
138  xgmin=geo_x(itmin,jtmin)
139  ygmin=geo_y(itmin,jtmin)
140  zgmin=geo_z(ktmin)
141  endif
142 
143  call mpp_max(xgmin)
144  call mpp_max(ygmin)
145  call mpp_max(zgmin)
146 
147  return
148 
149 
150 end subroutine mpp_array_global_min_max
151 ! </SUBROUTINE> NAME="mpp_array_global_min_max"
152 
153 
154 
155 end module mpp_utilities_mod
156 !> @}
157 ! close documentation grouping
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Reduction operations. Find the max of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:538
Reduction operations. Find the min of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:560
Reduction operation.
Definition: mpp.F90:597
subroutine, public mpp_array_global_min_max(in_array, tmask, isd, jsd, isc, iec, jsc, jec, nk, g_min, g_max, geo_x, geo_y, geo_z, xgmin, ygmin, zgmin, xgmax, ygmax, zgmax)
Compute and return the global min and max of an array and the corresponding lat-lon-depth locations .