source: LMDZ5/branches/testing/libf/phylmd/rrtm/eq_regions_mod.F90 @ 1999

Last change on this file since 1999 was 1999, checked in by Laurent Fairhead, 11 years ago

Merged trunk changes r1920:1997 into testing branch

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 11.9 KB
Line 
1module eq_regions_mod
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
63USE PARKIND1  ,ONLY : JPIM     ,JPRB
64
65IMPLICIT NONE
66
67SAVE
68
69PRIVATE
70
71PUBLIC eq_regions,l_regions_debug,n_regions_ns,n_regions_ew,n_regions,my_region_ns,my_region_ew
72
73real(kind=jprb) pi
74logical :: l_regions_debug=.false.
75integer(kind=jpim) :: n_regions_ns
76integer(kind=jpim) :: n_regions_ew
77integer(kind=jpim) :: my_region_ns
78integer(kind=jpim) :: my_region_ew
79integer(kind=jpim),allocatable :: n_regions(:)
80
81CONTAINS
82
83subroutine eq_regions(N)
84!
85! eq_regions uses the zonal equal area sphere partitioning algorithm to partition
86! the surface of a sphere into N regions of equal area and small diameter.
87!
88integer(kind=jpim),intent(in) :: N
89integer(kind=jpim) :: n_collars,j
90real(kind=jprb),allocatable :: r_regions(:)
91real(kind=jprb) :: c_polar
92
93pi=2.0_jprb*asin(1.0_jprb)
94
95n_regions(:)=0
96
97if( N == 1 )then
98
99  !
100  ! We have only one region, which must be the whole sphere.
101  !
102  n_regions(1)=1
103  n_regions_ns=1
104
105else
106
107  !
108  ! Given N, determine c_polar
109  ! the colatitude of the North polar spherical cap.
110  !
111  c_polar = polar_colat(N)
112  !
113  ! Given N, determine the ideal angle for spherical collars.
114  ! Based on N, this ideal angle, and c_polar,
115  ! determine n_collars, the number of collars between the polar caps.
116  !
117  n_collars = num_collars(N,c_polar,ideal_collar_angle(N))
118  n_regions_ns=n_collars+2
119  !
120  ! Given N, c_polar and n_collars, determine r_regions,
121  ! a list of the ideal real number of regions in each collar,
122  ! plus the polar caps.
123  ! The number of elements is n_collars+2.
124  ! r_regions[1] is 1.
125  ! r_regions[n_collars+2] is 1.
126  ! The sum of r_regions is N.
127  allocate(r_regions(n_collars+2))
128  call ideal_region_list(N,c_polar,n_collars,r_regions)
129  !
130  ! Given N and r_regions, determine n_regions, a list of the natural number
131  ! of regions in each collar and the polar caps.
132  ! This list is as close as possible to r_regions.
133  ! The number of elements is n_collars+2.
134  ! n_regions[1] is 1.
135  ! n_regions[n_collars+2] is 1.
136  ! The sum of n_regions is N.
137  !
138  call round_to_naturals(N,n_collars,r_regions)
139  deallocate(r_regions)
140  if( N /= sum(n_regions(:)) )then
141    write(*,'("eq_regions: N=",I10," sum(n_regions(:))=",I10)')N,sum(n_regions(:))
142    call abor1('eq_regions: N /= sum(n_regions)')
143  endif
144
145endif
146
147if( l_regions_debug )then
148  write(*,'("eq_regions: N=",I6," n_regions_ns=",I4)') N,n_regions_ns
149  do j=1,n_regions_ns
150    write(*,'("eq_regions: n_regions(",I4,")=",I4)') j,n_regions(j)
151  enddo
152endif
153n_regions_ew=maxval(n_regions(:))
154
155return
156end subroutine eq_regions
157
158function num_collars(N,c_polar,a_ideal) result(num_c)
159!
160!NUM_COLLARS The number of collars between the polar caps
161!
162! Given N, an ideal angle, and c_polar,
163! determine n_collars, the number of collars between the polar caps.
164!
165integer(kind=jpim),intent(in) :: N
166real(kind=jprb),intent(in) :: a_ideal,c_polar
167integer(kind=jpim) :: num_c
168logical enough
169enough = (N > 2) .and. (a_ideal > 0)
170if( enough )then
171  num_c = max(1,nint((pi-2.*c_polar)/a_ideal))
172else
173  num_c = 0
174endif
175return
176end function num_collars
177
178subroutine ideal_region_list(N,c_polar,n_collars,r_regions)
179!
180!IDEAL_REGION_LIST The ideal real number of regions in each zone
181!
182! List the ideal real number of regions in each collar, plus the polar caps.
183!
184! Given N, c_polar and n_collars, determine r_regions, a list of the ideal real
185! number of regions in each collar, plus the polar caps.
186! The number of elements is n_collars+2.
187! r_regions[1] is 1.
188! r_regions[n_collars+2] is 1.
189! The sum of r_regions is N.
190!
191integer(kind=jpim),intent(in) :: N,n_collars
192real(kind=jprb),intent(in) :: c_polar
193real(kind=jprb),intent(out) :: r_regions(n_collars+2)
194integer(kind=jpim) :: collar_n
195real(kind=jprb) :: ideal_region_area,ideal_collar_area
196real(kind=jprb) :: a_fitting
197r_regions(:)=0.0_jprb
198r_regions(1) = 1.0_jprb
199if( n_collars > 0 )then
200  !
201  ! Based on n_collars and c_polar, determine a_fitting,
202  ! the collar angle such that n_collars collars fit between the polar caps.
203  !
204  a_fitting = (pi-2.0_jprb*c_polar)/float(n_collars)
205  ideal_region_area = area_of_ideal_region(N)
206  do collar_n=1,n_collars
207    ideal_collar_area = area_of_collar(c_polar+(collar_n-1)*a_fitting, &
208     & c_polar+collar_n*a_fitting)
209    r_regions(1+collar_n) = ideal_collar_area / ideal_region_area
210  enddo
211endif
212r_regions(2+n_collars) = 1.
213return
214end subroutine ideal_region_list
215
216function ideal_collar_angle(N) result(ideal)
217!
218! IDEAL_COLLAR_ANGLE The ideal angle for spherical collars of an EQ partition
219!
220! IDEAL_COLLAR_ANGLE(N) sets ANGLE to the ideal angle for the
221! spherical collars of an EQ partition of the unit sphere S^2 into N regions.
222!
223integer(kind=jpim),intent(in) :: N
224real(kind=jprb) :: ideal
225ideal = area_of_ideal_region(N)**(0.5_jprb)
226return
227end function ideal_collar_angle
228
229subroutine round_to_naturals(N,n_collars,r_regions)
230!
231! ROUND_TO_NATURALS Round off a given list of numbers of regions
232!
233! Given N and r_regions, determine n_regions, a list of the natural number
234! of regions in each collar and the polar caps.
235! This list is as close as possible to r_regions, using rounding.
236! The number of elements is n_collars+2.
237! n_regions[1] is 1.
238! n_regions[n_collars+2] is 1.
239! The sum of n_regions is N.
240!
241integer(kind=jpim),intent(in) :: N,n_collars
242real(kind=jprb),intent(in) :: r_regions(n_collars+2)
243integer(kind=jpim) :: zone_n
244real(kind=jprb) :: discrepancy
245n_regions(1:n_collars+2) = r_regions(:)
246discrepancy = 0.0_jprb
247do zone_n = 1,n_collars+2
248    n_regions(zone_n) = nint(r_regions(zone_n)+discrepancy);
249    discrepancy = discrepancy+r_regions(zone_n)-float(n_regions(zone_n));
250enddo
251return
252end subroutine round_to_naturals
253
254function polar_colat(N) result(polar_c)
255!
256! Given N, determine the colatitude of the North polar spherical cap.
257!
258integer(kind=jpim),intent(in) :: N
259real(kind=jprb) :: area
260real(kind=jprb) :: polar_c
261if( N == 1 ) polar_c=pi
262if( N == 2 ) polar_c=pi/2.0_jprb
263if( N > 2 )then
264  area=area_of_ideal_region(N)
265  polar_c=sradius_of_cap(area)
266endif
267return
268end function polar_colat
269
270function area_of_ideal_region(N) result(area)
271!
272! AREA_OF_IDEAL_REGION(N) sets AREA to be the area of one of N equal
273! area regions on S^2, that is 1/N times AREA_OF_SPHERE.
274!
275integer(kind=jpim),intent(in) :: N
276real(kind=jprb) :: area_of_sphere
277real(kind=jprb) :: area
278area_of_sphere = (2.0_jprb*pi**1.5_jprb/gamma(1.5_jprb))
279area = area_of_sphere/float(N)
280return
281end function area_of_ideal_region
282
283function sradius_of_cap(area) result(sradius)
284!
285! SRADIUS_OF_CAP(AREA) returns the spherical radius of
286! an S^2 spherical cap of area AREA.
287!
288real(kind=jprb),intent(in) :: area
289real(kind=jprb) :: sradius
290sradius = 2.0_jprb*asin(sqrt(area/pi)/2.0_jprb)
291return
292end function sradius_of_cap
293
294function area_of_collar(a_top, a_bot) result(area)
295!
296! AREA_OF_COLLAR Area of spherical collar
297!
298! AREA_OF_COLLAR(A_TOP, A_BOT) sets AREA to be the area of an S^2 spherical
299! collar specified by A_TOP, A_BOT, where A_TOP is top (smaller) spherical radius,
300! A_BOT is bottom (larger) spherical radius.
301!
302real(kind=jprb),intent(in) :: a_top,a_bot
303real(kind=jprb) area
304area = area_of_cap(a_bot) - area_of_cap(a_top)
305return
306end function area_of_collar
307
308function area_of_cap(s_cap) result(area)
309!
310! AREA_OF_CAP Area of spherical cap
311!
312! AREA_OF_CAP(S_CAP) sets AREA to be the area of an S^2 spherical
313! cap of spherical radius S_CAP.
314!
315real(kind=jprb),intent(in) :: s_cap
316real(kind=jprb) area
317area = 4.0_jprb*pi * sin(s_cap/2.0_jprb)**2
318return
319end function area_of_cap
320
321function gamma(x) result(gamma_res)
322real(kind=jprb),intent(in) :: x
323real(kind=jprb) :: gamma_res
324real(kind=jprb) :: p0,p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13
325real(kind=jprb) :: w,y
326integer(kind=jpim) :: k,n
327parameter (&
328& p0 =   0.999999999999999990e+00_jprb,&
329& p1 =  -0.422784335098466784e+00_jprb,&
330& p2 =  -0.233093736421782878e+00_jprb,&
331& p3 =   0.191091101387638410e+00_jprb,&
332& p4 =  -0.024552490005641278e+00_jprb,&
333& p5 =  -0.017645244547851414e+00_jprb,&
334& p6 =   0.008023273027855346e+00_jprb)
335parameter (&
336& p7 =  -0.000804329819255744e+00_jprb,&
337& p8 =  -0.000360837876648255e+00_jprb,&
338& p9 =   0.000145596568617526e+00_jprb,&
339& p10 = -0.000017545539395205e+00_jprb,&
340& p11 = -0.000002591225267689e+00_jprb,&
341& p12 =  0.000001337767384067e+00_jprb,&
342& p13 = -0.000000199542863674e+00_jprb)
343n = nint(x - 2)
344w = x - (n + 2)
345y = ((((((((((((p13 * w + p12) * w + p11) * w + p10) *&
346&    w + p9) * w + p8) * w + p7) * w + p6) * w + p5) *&
347&    w + p4) * w + p3) * w + p2) * w + p1) * w + p0
348if (n .gt. 0) then
349  w = x - 1
350  do k = 2, n
351    w = w * (x - k)
352  end do
353else
354  w = 1
355  do k = 0, -n - 1
356    y = y * (x + k)
357  end do
358end if
359gamma_res = w / y
360return
361end function gamma
362
363end module eq_regions_mod
Note: See TracBrowser for help on using the repository browser.