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