LMDZ
eq_regions_mod.F90
Go to the documentation of this file.
2 !
3 ! Purpose.
4 ! --------
5 ! eq_regions_mod provides the code to perform a high level
6 ! partitioning of the surface of a sphere into regions of
7 ! equal area and small diameter.
8 ! the type.
9 !
10 ! Background.
11 ! -----------
12 ! This Fortran version of eq_regions is a much cut down version of the
13 ! "Recursive Zonal Equal Area (EQ) Sphere Partitioning Toolbox" of the
14 ! same name developed by Paul Leopardi at the University of New South Wales.
15 ! This version has been coded specifically for the case of partitioning the
16 ! surface of a sphere or S^dim (where dim=2) as denoted in the original code.
17 ! Only a subset of the original eq_regions package has been coded to determine
18 ! the high level distribution of regions on a sphere, as the detailed
19 ! distribution of grid points to each region is left to IFS software.
20 ! This is required to take into account the spatial distribution of grid
21 ! points in an IFS gaussian grid and provide an optimal (i.e. exact)
22 ! distribution of grid points over regions.
23 !
24 ! The following copyright notice for the eq_regions package is included from
25 ! the original MatLab release.
26 !
27 ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28 ! + Release 1.10 2005-06-26 +
29 ! + +
30 ! + Copyright (c) 2004, 2005, University of New South Wales +
31 ! + +
32 ! + Permission is hereby granted, free of charge, to any person obtaining +
33 ! + a copy of this software and associated documentation files (the +
34 ! + "Software"), to deal in the Software without restriction, including +
35 ! + without limitation the rights to use, copy, modify, merge, publish, +
36 ! + distribute, sublicense, and/or sell copies of the Software, and to +
37 ! + permit persons to whom the Software is furnished to do so, subject to +
38 ! + the following conditions: +
39 ! + +
40 ! + The above copyright notice and this permission notice shall be included +
41 ! + in all copies or substantial portions of the Software. +
42 ! + +
43 ! + THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, +
44 ! + EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF +
45 ! + MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. +
46 ! + IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY +
47 ! + CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, +
48 ! + TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE +
49 ! + SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE. +
50 ! + +
51 ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
52 !
53 ! Author.
54 ! -------
55 ! George Mozdzynski *ECMWF*
56 !
57 ! Modifications.
58 ! --------------
59 ! Original : 2006-04-15
60 !
61 !--------------------------------------------------------------------------------
62 
63 USE parkind1 ,ONLY : jpim ,jprb
64 
65 IMPLICIT NONE
66 
67 SAVE
68 
69 PRIVATE
70 
72 
73 real(kind=jprb) pi
74 logical :: l_regions_debug=.false.
75 integer(kind=jpim) :: n_regions_ns
76 integer(kind=jpim) :: n_regions_ew
77 integer(kind=jpim) :: my_region_ns
78 integer(kind=jpim) :: my_region_ew
79 integer(kind=jpim),allocatable :: n_regions(:)
80 
81 
82 !$OMP THREADPRIVATE(l_regions_debug,my_region_ew,my_region_ns,n_regions_ew,n_regions_ns,pi,n_regions)
83 
84 CONTAINS
85 
86 subroutine eq_regions(N)
87 !
88 ! eq_regions uses the zonal equal area sphere partitioning algorithm to partition
89 ! the surface of a sphere into N regions of equal area and small diameter.
90 !
91 integer(kind=jpim),intent(in) :: N
92 integer(kind=jpim) :: n_collars,j
93 real(kind=jprb),allocatable :: r_regions(:)
94 real(kind=jprb) :: c_polar
95 
96 pi=2.0_jprb*asin(1.0_jprb)
97 
98 n_regions(:)=0
99 
100 if( n == 1 )then
101 
102  !
103  ! We have only one region, which must be the whole sphere.
104  !
105  n_regions(1)=1
106  n_regions_ns=1
107 
108 else
109 
110  !
111  ! Given N, determine c_polar
112  ! the colatitude of the North polar spherical cap.
113  !
114  c_polar = polar_colat(n)
115  !
116  ! Given N, determine the ideal angle for spherical collars.
117  ! Based on N, this ideal angle, and c_polar,
118  ! determine n_collars, the number of collars between the polar caps.
119  !
120  n_collars = num_collars(n,c_polar,ideal_collar_angle(n))
121  n_regions_ns=n_collars+2
122  !
123  ! Given N, c_polar and n_collars, determine r_regions,
124  ! a list of the ideal real number of regions in each collar,
125  ! plus the polar caps.
126  ! The number of elements is n_collars+2.
127  ! r_regions[1] is 1.
128  ! r_regions[n_collars+2] is 1.
129  ! The sum of r_regions is N.
130  allocate(r_regions(n_collars+2))
131  call ideal_region_list(n,c_polar,n_collars,r_regions)
132  !
133  ! Given N and r_regions, determine n_regions, a list of the natural number
134  ! of regions in each collar and the polar caps.
135  ! This list is as close as possible to r_regions.
136  ! The number of elements is n_collars+2.
137  ! n_regions[1] is 1.
138  ! n_regions[n_collars+2] is 1.
139  ! The sum of n_regions is N.
140  !
141  call round_to_naturals(n,n_collars,r_regions)
142  deallocate(r_regions)
143  if( n /= sum(n_regions(:)) )then
144  write(*,'("eq_regions: N=",I10," sum(n_regions(:))=",I10)')n,sum(n_regions(:))
145  call abor1('eq_regions: N /= sum(n_regions)')
146  endif
147 
148 endif
149 
150 if( l_regions_debug )then
151  write(*,'("eq_regions: N=",I6," n_regions_ns=",I4)') n,n_regions_ns
152  do j=1,n_regions_ns
153  write(*,'("eq_regions: n_regions(",I4,")=",I4)') j,n_regions(j)
154  enddo
155 endif
156 n_regions_ew=maxval(n_regions(:))
157 
158 return
159 end subroutine eq_regions
160 
161 function num_collars(N,c_polar,a_ideal) result(num_c)
162 !
163 !NUM_COLLARS The number of collars between the polar caps
164 !
165 ! Given N, an ideal angle, and c_polar,
166 ! determine n_collars, the number of collars between the polar caps.
167 !
168 integer(kind=jpim),intent(in) :: N
169 real(kind=jprb),intent(in) :: a_ideal,c_polar
170 integer(kind=jpim) :: num_c
171 logical enough
172 enough = (n > 2) .and. (a_ideal > 0)
173 if( enough )then
174  num_c = max(1,nint((pi-2.*c_polar)/a_ideal))
175 else
176  num_c = 0
177 endif
178 return
179 end function num_collars
180 
181 subroutine ideal_region_list(N,c_polar,n_collars,r_regions)
182 !
183 !IDEAL_REGION_LIST The ideal real number of regions in each zone
184 !
185 ! List the ideal real number of regions in each collar, plus the polar caps.
186 !
187 ! Given N, c_polar and n_collars, determine r_regions, a list of the ideal real
188 ! number of regions in each collar, plus the polar caps.
189 ! The number of elements is n_collars+2.
190 ! r_regions[1] is 1.
191 ! r_regions[n_collars+2] is 1.
192 ! The sum of r_regions is N.
193 !
194 integer(kind=jpim),intent(in) :: N,n_collars
195 real(kind=jprb),intent(in) :: c_polar
196 real(kind=jprb),intent(out) :: r_regions(n_collars+2)
197 integer(kind=jpim) :: collar_n
198 real(kind=jprb) :: ideal_region_area,ideal_collar_area
199 real(kind=jprb) :: a_fitting
200 r_regions(:)=0.0_jprb
201 r_regions(1) = 1.0_jprb
202 if( n_collars > 0 )then
203  !
204  ! Based on n_collars and c_polar, determine a_fitting,
205  ! the collar angle such that n_collars collars fit between the polar caps.
206  !
207  a_fitting = (pi-2.0_jprb*c_polar)/float(n_collars)
208  ideal_region_area = area_of_ideal_region(n)
209  do collar_n=1,n_collars
210  ideal_collar_area = area_of_collar(c_polar+(collar_n-1)*a_fitting, &
211  & c_polar+collar_n*a_fitting)
212  r_regions(1+collar_n) = ideal_collar_area / ideal_region_area
213  enddo
214 endif
215 r_regions(2+n_collars) = 1.
216 return
217 end subroutine ideal_region_list
218 
219 function ideal_collar_angle(N) result(ideal)
220 !
221 ! IDEAL_COLLAR_ANGLE The ideal angle for spherical collars of an EQ partition
222 !
223 ! IDEAL_COLLAR_ANGLE(N) sets ANGLE to the ideal angle for the
224 ! spherical collars of an EQ partition of the unit sphere S^2 into N regions.
225 !
226 integer(kind=jpim),intent(in) :: N
227 real(kind=jprb) :: ideal
228 ideal = area_of_ideal_region(n)**(0.5_jprb)
229 return
230 end function ideal_collar_angle
231 
232 subroutine round_to_naturals(N,n_collars,r_regions)
233 !
234 ! ROUND_TO_NATURALS Round off a given list of numbers of regions
235 !
236 ! Given N and r_regions, determine n_regions, a list of the natural number
237 ! of regions in each collar and the polar caps.
238 ! This list is as close as possible to r_regions, using rounding.
239 ! The number of elements is n_collars+2.
240 ! n_regions[1] is 1.
241 ! n_regions[n_collars+2] is 1.
242 ! The sum of n_regions is N.
243 !
244 integer(kind=jpim),intent(in) :: N,n_collars
245 real(kind=jprb),intent(in) :: r_regions(n_collars+2)
246 integer(kind=jpim) :: zone_n
247 real(kind=jprb) :: discrepancy
248 n_regions(1:n_collars+2) = r_regions(:)
249 discrepancy = 0.0_jprb
250 do zone_n = 1,n_collars+2
251  n_regions(zone_n) = nint(r_regions(zone_n)+discrepancy);
252  discrepancy = discrepancy+r_regions(zone_n)-float(n_regions(zone_n));
253 enddo
254 return
255 end subroutine round_to_naturals
256 
257 function polar_colat(N) result(polar_c)
258 !
259 ! Given N, determine the colatitude of the North polar spherical cap.
260 !
261 integer(kind=jpim),intent(in) :: N
262 real(kind=jprb) :: area
263 real(kind=jprb) :: polar_c
264 if( n == 1 ) polar_c=pi
265 if( n == 2 ) polar_c=pi/2.0_jprb
266 if( n > 2 )then
267  area=area_of_ideal_region(n)
268  polar_c=sradius_of_cap(area)
269 endif
270 return
271 end function polar_colat
272 
273 function area_of_ideal_region(N) result(area)
274 !
275 ! AREA_OF_IDEAL_REGION(N) sets AREA to be the area of one of N equal
276 ! area regions on S^2, that is 1/N times AREA_OF_SPHERE.
277 !
278 integer(kind=jpim),intent(in) :: N
279 real(kind=jprb) :: area_of_sphere
280 real(kind=jprb) :: area
281 area_of_sphere = (2.0_jprb*pi**1.5_jprb/gamma(1.5_jprb))
282 area = area_of_sphere/float(n)
283 return
284 end function area_of_ideal_region
285 
286 function sradius_of_cap(area) result(sradius)
287 !
288 ! SRADIUS_OF_CAP(AREA) returns the spherical radius of
289 ! an S^2 spherical cap of area AREA.
290 !
291 real(kind=jprb),intent(in) :: area
292 real(kind=jprb) :: sradius
293 sradius = 2.0_jprb*asin(sqrt(area/pi)/2.0_jprb)
294 return
295 end function sradius_of_cap
296 
297 function area_of_collar(a_top, a_bot) result(area)
298 !
299 ! AREA_OF_COLLAR Area of spherical collar
300 !
301 ! AREA_OF_COLLAR(A_TOP, A_BOT) sets AREA to be the area of an S^2 spherical
302 ! collar specified by A_TOP, A_BOT, where A_TOP is top (smaller) spherical radius,
303 ! A_BOT is bottom (larger) spherical radius.
304 !
305 real(kind=jprb),intent(in) :: a_top,a_bot
306 real(kind=jprb) area
307 area = area_of_cap(a_bot) - area_of_cap(a_top)
308 return
309 end function area_of_collar
310 
311 function area_of_cap(s_cap) result(area)
312 !
313 ! AREA_OF_CAP Area of spherical cap
314 !
315 ! AREA_OF_CAP(S_CAP) sets AREA to be the area of an S^2 spherical
316 ! cap of spherical radius S_CAP.
317 !
318 real(kind=jprb),intent(in) :: s_cap
319 real(kind=jprb) area
320 area = 4.0_jprb*pi * sin(s_cap/2.0_jprb)**2
321 return
322 end function area_of_cap
323 
324 function gamma(x) result(gamma_res)
325 real(kind=jprb),intent(in) :: x
326 real(kind=jprb) :: gamma_res
327 real(kind=jprb) :: p0,p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13
328 real(kind=jprb) :: w,y
329 integer(kind=jpim) :: k,n
330 parameter(&
331 & p0 = 0.999999999999999990e+00_jprb,&
332 & p1 = -0.422784335098466784e+00_jprb,&
333 & p2 = -0.233093736421782878e+00_jprb,&
334 & p3 = 0.191091101387638410e+00_jprb,&
335 & p4 = -0.024552490005641278e+00_jprb,&
336 & p5 = -0.017645244547851414e+00_jprb,&
337 & p6 = 0.008023273027855346e+00_jprb)
338 parameter(&
339 & p7 = -0.000804329819255744e+00_jprb,&
340 & p8 = -0.000360837876648255e+00_jprb,&
341 & p9 = 0.000145596568617526e+00_jprb,&
342 & p10 = -0.000017545539395205e+00_jprb,&
343 & p11 = -0.000002591225267689e+00_jprb,&
344 & p12 = 0.000001337767384067e+00_jprb,&
345 & p13 = -0.000000199542863674e+00_jprb)
346 n = nint(x - 2)
347 w = x - (n + 2)
348 y = ((((((((((((p13 * w + p12) * w + p11) * w + p10) *&
349 & w + p9) * w + p8) * w + p7) * w + p6) * w + p5) *&
350 & w + p4) * w + p3) * w + p2) * w + p1) * w + p0
351 if (n .gt. 0) then
352  w = x - 1
353  do k = 2, n
354  w = w * (x - k)
355  end do
356 else
357  w = 1
358  do k = 0, -n - 1
359  y = y * (x + k)
360  end do
361 end if
362 gamma_res = w / y
363 return
364 end function gamma
365 
366 end module eq_regions_mod
real(kind=jprb) function area_of_ideal_region(N)
subroutine, public eq_regions(N)
integer(kind=jpim), public n_regions_ns
!$Id mode_top_bound COMMON comconstr && pi
Definition: comconst.h:7
subroutine abor1(CDTEXT)
Definition: abor1.F90:2
real(kind=jprb) function ideal_collar_angle(N)
real(kind=jprb) function polar_colat(N)
integer(kind=jpim), public my_region_ns
!$Id itapm1 ENDIF!IM on interpole les champs sur les niveaux STD de pression!IM a chaque pas de temps de la physique c!positionnement de l argument logique a false c!pour ne pas recalculer deux fois la meme chose!c!a cet effet un appel a plevel_new a ete deplace c!a la fin de la serie d appels c!la boucle DO nlevSTD a ete internalisee c!dans d ou la creation de cette routine c c!CALL false
Definition: calcul_STDlev.h:26
integer(kind=jpim), public my_region_ew
integer, parameter jprb
Definition: parkind1.F90:31
!$Header!integer nvarmx parameter(nfmx=10, imx=200, jmx=150, lmx=200, nvarmx=1000) real xd(imx
real(kind=jprb) function sradius_of_cap(area)
real(kind=jprb) function area_of_cap(s_cap)
integer(kind=jpim), public n_regions_ew
logical, public l_regions_debug
real(kind=jprb) function area_of_collar(a_top, a_bot)
integer, parameter jpim
Definition: parkind1.F90:13
subroutine round_to_naturals(N, n_collars, r_regions)
integer(kind=jpim), dimension(:), allocatable, public n_regions
integer(kind=jpim) function num_collars(N, c_polar, a_ideal)
subroutine ideal_region_list(N, c_polar, n_collars, r_regions)
real(kind=jprb) function gamma(x)