Cmake build action (#865)
[FMS.git] / mpp / mpp_utilities.F90
blobf0313fbeb0fe93cb6c233cc07ac707fd8ea692d3
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.
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.
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
23 !> Currently only holds one routine for finding global min and max
25 !> @file
26 !> @brief File for @ref mpp_utilities_mod
28 !> @addtogroup mpp_utilities_mod
29 !> @{
30 module mpp_utilities_mod
32 !-----------------------------------------------------------------------
33 ! Include variable "version" to be written to log file.
34 #include<file_version.h>
35 !-----------------------------------------------------------------------
37   public :: mpp_array_global_min_max
39 contains
41 !#######################################################################
42 !> @brief Compute and return the global min and max of an array
43 !! and the corresponding lat-lon-depth locations .
45 !> This algorithm works only for an input array that has a unique global
46 !! max and min location. This is assured by introducing a factor that distinguishes
47 !! the values of extrema at each processor.
49 !! Vectorized using maxloc() and minloc() intrinsic functions by
50 !! Russell.Fiedler@csiro.au (May 2005).
52 !! Modified by Zhi.Liang@noaa.gov (July 2005)
54 !! Modified by Niki.Zadeh@noaa.gov (Feb. 2009)
56 subroutine mpp_array_global_min_max(in_array, tmask,isd,jsd,isc,iec,jsc,jec,nk, g_min, g_max, &
57                                     geo_x,geo_y,geo_z, xgmin, ygmin, zgmin, xgmax, ygmax, zgmax)
59   use mpp_mod,           only: mpp_min, mpp_max, mpp_pe, mpp_sum
61   real, dimension(isd:,jsd:,:), intent(in) :: in_array
62   real, dimension(isd:,jsd:,:), intent(in) :: tmask
63   integer,                      intent(in) :: isd,jsd,isc,iec,jsc,jec,nk
64   real,                         intent(out):: g_min, g_max
65   real, dimension(isd:,jsd:),   intent(in) :: geo_x,geo_y
66   real, dimension(:),           intent(in) :: geo_z
67   real,                         intent(out):: xgmin, ygmin, zgmin, xgmax, ygmax, zgmax
71   real    :: tmax, tmin, tmax0, tmin0
72   integer :: itmax, jtmax, ktmax, itmin, jtmin, ktmin
73   integer :: igmax, jgmax, kgmax, igmin, jgmin, kgmin
74   real    :: fudge
76   ! arrays to enable vectorization
77   integer :: iminarr(3),imaxarr(3)
79   g_min=-88888888888.0 ; g_max=-999999999.0
81   tmax=-1.e10;tmin=1.e10
82   itmax=0;jtmax=0;ktmax=0
83   itmin=0;jtmin=0;ktmin=0
85   if(ANY(tmask(isc:iec,jsc:jec,:) > 0.)) then
86      iminarr=minloc(in_array(isc:iec,jsc:jec,:),tmask(isc:iec,jsc:jec,:) > 0.)
87      imaxarr=maxloc(in_array(isc:iec,jsc:jec,:),tmask(isc:iec,jsc:jec,:) > 0.)
88      itmin=iminarr(1)+isc-1
89      jtmin=iminarr(2)+jsc-1
90      ktmin=iminarr(3)
91      itmax=imaxarr(1)+isc-1
92      jtmax=imaxarr(2)+jsc-1
93      ktmax=imaxarr(3)
94      tmin=in_array(itmin,jtmin,ktmin)
95      tmax=in_array(itmax,jtmax,ktmax)
96   end if
98   ! use "fudge" to distinguish processors when tracer extreme is independent of processor
99   fudge = 1.0 + 1.e-12*real(mpp_pe() )
100   tmax = tmax*fudge
101   tmin = tmin*fudge
102   if(tmax == 0.0) then
103     tmax = tmax + 1.e-12*real(mpp_pe() )
104   endif
105   if(tmin == 0.0) then
106     tmin = tmin + 1.e-12*real(mpp_pe() )
107   endif
110   tmax0=tmax;tmin0=tmin
112   call mpp_max(tmax)
113   call mpp_min(tmin)
115   g_max = tmax
116   g_min = tmin
118   !Now find the location of the global extrema.
119   !
120   !Note that the fudge factor above guarantees that the location of max (min) is uinque,
121   ! since tmax0 (tmin0) has slightly different values on each processor.
122   !Otherwise, the function in_array(i,j,k) could be equal to global max (min) at more
123   ! than one point in space and this would be a much more difficult problem to solve.
124   !
125   !mpp_max trick
126   !-999 on all current PE's
127   xgmax=-999.; ygmax=-999.; zgmax=-999.
128   xgmin=-999.; ygmin=-999.; zgmin=-999.
131   !except when
132   if (tmax0 == tmax) then !This happens ONLY on ONE processor because of fudge factor above.
133      xgmax=geo_x(itmax,jtmax)
134      ygmax=geo_y(itmax,jtmax)
135      zgmax=geo_z(ktmax)
136   endif
138   call mpp_max(xgmax)
139   call mpp_max(ygmax)
140   call mpp_max(zgmax)
142   if (tmin0 == tmin) then !This happens ONLY on ONE processor because of fudge factor above.
143      xgmin=geo_x(itmin,jtmin)
144      ygmin=geo_y(itmin,jtmin)
145      zgmin=geo_z(ktmin)
146   endif
148   call mpp_max(xgmin)
149   call mpp_max(ygmin)
150   call mpp_max(zgmin)
152   return
155 end subroutine mpp_array_global_min_max
156 ! </SUBROUTINE>  NAME="mpp_array_global_min_max"
160 end module mpp_utilities_mod
161 !> @}
162 ! close documentation grouping