1 |
module grid_noro_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
SUBROUTINE grid_noro(xdata, ydata, zdata, x, y, zphi, zmea, zstd, zsig, & |
8 |
zgam, zthe, zpic, zval, mask) |
9 |
|
10 |
! From dyn3d/grid_noro.F, version 1.1.1.1 2004/05/19 12:53:06 |
11 |
|
12 |
! Authors: F. Lott, Z. X. Li, A. Harzallah and L. Fairhead |
13 |
|
14 |
! Compute the parameters of the SSO scheme as described in |
15 |
! Lott and Miller (1997) and Lott (1999). |
16 |
! Target points are on a rectangular grid: |
17 |
! jjm + 1 latitudes including North and South Poles; |
18 |
! iim + 1 longitudes, with periodicity: longitude(iim + 1) = longitude(1) |
19 |
! At the poles the field value is repeated iim + 1 times. |
20 |
|
21 |
! The parameters a, b, c, d represent the limite of the target |
22 |
! gridpoint region. The means over this region are calculated from |
23 |
! USN data, ponderated by a weight proportional to the surface |
24 |
! occupied by the data inside the model gridpoint area. In most |
25 |
! circumstances, this weight is the ratio between the surface of |
26 |
! the USN gridpoint area and the surface of the model gridpoint |
27 |
! area. See "grid_noto.txt". |
28 |
|
29 |
use dimens_m, only: iim, jjm |
30 |
use nr_util, only: assert, pi |
31 |
use mva9_m, only: mva9 |
32 |
|
33 |
REAL, intent(in):: xdata(:), ydata(:) ! coordinates of input field |
34 |
REAL, intent(in):: zdata(:, :) ! input field |
35 |
REAL, intent(in):: x(:), y(:) ! ccordinates output field |
36 |
|
37 |
! Correlations of USN orography gradients: |
38 |
|
39 |
REAL, intent(out):: zphi(:, :) |
40 |
real, intent(out):: zmea(:, :) ! Mean orography |
41 |
real, intent(out):: zstd(:, :) ! Standard deviation |
42 |
REAL zsig(:, :) ! Slope |
43 |
real zgam(:, :) ! Anisotropy |
44 |
real zthe(:, :) ! Orientation of the small axis |
45 |
REAL, intent(out):: zpic(:, :) ! Maximum altitude |
46 |
real, intent(out):: zval(:, :) ! Minimum altitude |
47 |
|
48 |
real, intent(out):: mask(:, :) ! fraction of land |
49 |
|
50 |
! Variables local to the procedure: |
51 |
|
52 |
! In this version it is assumed that the input data come from |
53 |
! the US Navy dataset: |
54 |
integer, parameter:: iusn = 2160, jusn = 1080 |
55 |
integer, parameter:: iext = 216 |
56 |
REAL xusn(iusn + 2 * iext), yusn(jusn + 2) |
57 |
REAL zusn(iusn + 2 * iext, jusn + 2) |
58 |
|
59 |
! Intermediate fields (correlations of orography gradient) |
60 |
|
61 |
REAL ztz(iim + 1, jjm + 1), zxtzx(iim + 1, jjm + 1) |
62 |
REAL zytzy(iim + 1, jjm + 1), zxtzy(iim + 1, jjm + 1) |
63 |
REAL weight(iim + 1, jjm + 1) |
64 |
|
65 |
! Correlations of USN orography gradients: |
66 |
REAL, dimension(iusn + 2 * iext, jusn + 2):: zxtzxusn, zytzyusn, zxtzyusn |
67 |
|
68 |
real mask_tmp(size(x), size(y)) |
69 |
real num_tot(iim + 1, jjm + 1), num_lan(iim + 1, jjm + 1) |
70 |
|
71 |
REAL a(iim + 1), b(iim + 1), c(jjm + 1), d(jjm + 1) |
72 |
real rad, weighx, weighy, xincr, xk, xp, xm, xw, xq, xl |
73 |
real zbordnor, zdeltax, zbordsud, zdeltay, zbordoue, zlenx, zleny, zmeasud |
74 |
real zllmpic, zllmmea, zllmgam, zllmthe, zllmstd, zllmsig, zllmval |
75 |
real zpicnor, zminthe, zsigsud, zstdnor, zstdsud, zvalsud, zvalnor |
76 |
real zweinor, zweisud, zsignor, zpicsud, zmeanor, zbordest |
77 |
integer ii, i, jj, j |
78 |
|
79 |
!------------------------------- |
80 |
|
81 |
print *, "Call sequence information: grid_noro" |
82 |
|
83 |
call assert((/size(xdata), size(zdata, 1)/) == iusn, "grid_noro iusn") |
84 |
call assert((/size(ydata), size(zdata, 2)/) == jusn, "grid_noro jusn") |
85 |
|
86 |
call assert((/size(x), size(zphi, 1), size(zmea, 1), size(zstd, 1), & |
87 |
size(zsig, 1), size(zgam, 1), size(zthe, 1), size(zpic, 1), & |
88 |
size(zval, 1), size(mask, 1)/) == iim + 1, "grid_noro iim") |
89 |
|
90 |
call assert((/size(y), size(zphi, 2), size(zmea, 2), size(zstd, 2), & |
91 |
size(zsig, 2), size(zgam, 2), size(zthe, 2), size(zpic, 2), & |
92 |
size(zval, 2), size(mask, 2)/) == jjm + 1, "grid_noro jjm") |
93 |
|
94 |
print *, "Paramètres de l'orographie à l'échelle sous-maille" |
95 |
rad = 6371229. |
96 |
zdeltay = 2. * pi / real(jusn) * rad |
97 |
|
98 |
! Extension of the USN database to POCEED computations at boundaries: |
99 |
|
100 |
DO j = 1, jusn |
101 |
yusn(j + 1) = ydata(j) |
102 |
DO i = 1, iusn |
103 |
zusn(i + iext, j + 1) = zdata(i, j) |
104 |
xusn(i + iext) = xdata(i) |
105 |
ENDDO |
106 |
DO i = 1, iext |
107 |
zusn(i, j + 1) = zdata(iusn - iext + i, j) |
108 |
xusn(i) = xdata(iusn - iext + i) - 2. * pi |
109 |
zusn(iusn + iext + i, j + 1) = zdata(i, j) |
110 |
xusn(iusn + iext + i) = xdata(i) + 2. * pi |
111 |
ENDDO |
112 |
ENDDO |
113 |
|
114 |
yusn(1) = ydata(1) + (ydata(1) - ydata(2)) |
115 |
yusn(jusn + 2) = ydata(jusn) + (ydata(jusn) - ydata(jusn - 1)) |
116 |
DO i = 1, iusn / 2 + iext |
117 |
zusn(i, 1) = zusn(i + iusn / 2, 2) |
118 |
zusn(i + iusn / 2 + iext, 1) = zusn(i, 2) |
119 |
zusn(i, jusn + 2) = zusn(i + iusn / 2, jusn + 1) |
120 |
zusn(i + iusn / 2 + iext, jusn + 2) = zusn(i, jusn + 1) |
121 |
ENDDO |
122 |
|
123 |
! COMPUTE LIMITS OF MODEL GRIDPOINT AREA (REGULAR GRID) |
124 |
|
125 |
a(1) = x(1) - (x(2) - x(1)) / 2.0 |
126 |
b(1) = (x(1) + x(2)) / 2.0 |
127 |
DO i = 2, iim |
128 |
a(i) = b(i - 1) |
129 |
b(i) = (x(i) + x(i + 1)) / 2.0 |
130 |
ENDDO |
131 |
a(iim + 1) = b(iim) |
132 |
b(iim + 1) = x(iim + 1) + (x(iim + 1) - x(iim)) / 2.0 |
133 |
|
134 |
c(1) = y(1) - (y(2) - y(1)) / 2.0 |
135 |
d(1) = (y(1) + y(2)) / 2.0 |
136 |
DO j = 2, jjm |
137 |
c(j) = d(j - 1) |
138 |
d(j) = (y(j) + y(j + 1)) / 2.0 |
139 |
ENDDO |
140 |
c(jjm + 1) = d(jjm) |
141 |
d(jjm + 1) = y(jjm + 1) + (y(jjm + 1) - y(jjm)) / 2.0 |
142 |
|
143 |
! Initialisations : |
144 |
weight = 0. |
145 |
zxtzx = 0. |
146 |
zytzy = 0. |
147 |
zxtzy = 0. |
148 |
ztz = 0. |
149 |
zmea = 0. |
150 |
zpic = - 1E10 |
151 |
zval = 1E10 |
152 |
|
153 |
! COMPUTE SLOPES CORRELATIONS ON USN GRID |
154 |
|
155 |
zytzyusn = 0. |
156 |
zxtzxusn = 0. |
157 |
zxtzyusn = 0. |
158 |
|
159 |
DO j = 2, jusn + 1 |
160 |
zdeltax = zdeltay * cos(yusn(j)) |
161 |
DO i = 2, iusn + 2 * iext - 1 |
162 |
zytzyusn(i, j) = (zusn(i, j + 1) - zusn(i, j - 1))**2 / zdeltay**2 |
163 |
zxtzxusn(i, j) = (zusn(i + 1, j) - zusn(i - 1, j))**2 / zdeltax**2 |
164 |
zxtzyusn(i, j) = (zusn(i, j + 1) - zusn(i, j - 1)) / zdeltay & |
165 |
* (zusn(i + 1, j) - zusn(i - 1, j)) / zdeltax |
166 |
ENDDO |
167 |
ENDDO |
168 |
|
169 |
! SUMMATION OVER GRIDPOINT AREA |
170 |
|
171 |
zleny = pi / real(jusn) * rad |
172 |
xincr = pi / 2. / real(jusn) |
173 |
DO ii = 1, iim + 1 |
174 |
DO jj = 1, jjm + 1 |
175 |
num_tot(ii, jj) = 0. |
176 |
num_lan(ii, jj) = 0. |
177 |
DO j = 2, jusn + 1 |
178 |
zlenx = zleny * cos(yusn(j)) |
179 |
zdeltax = zdeltay * cos(yusn(j)) |
180 |
zbordnor = (c(jj) - yusn(j) + xincr) * rad |
181 |
zbordsud = (yusn(j) - d(jj) + xincr) * rad |
182 |
weighy = AMAX1(0., amin1(zbordnor, zbordsud, zleny)) |
183 |
IF (weighy /= 0) THEN |
184 |
DO i = 2, iusn + 2 * iext - 1 |
185 |
zbordest = (xusn(i) - a(ii) + xincr) * rad * cos(yusn(j)) |
186 |
zbordoue = (b(ii) + xincr - xusn(i)) * rad * cos(yusn(j)) |
187 |
weighx = AMAX1(0., amin1(zbordest, zbordoue, zlenx)) |
188 |
IF (weighx /= 0) THEN |
189 |
num_tot(ii, jj) = num_tot(ii, jj) + 1. |
190 |
if (zusn(i, j) >= 1.) then |
191 |
num_lan(ii, jj) = num_lan(ii, jj) + 1. |
192 |
end if |
193 |
weight(ii, jj) = weight(ii, jj) + weighx * weighy |
194 |
zxtzx(ii, jj) = zxtzx(ii, jj) & |
195 |
+ zxtzxusn(i, j) * weighx * weighy |
196 |
zytzy(ii, jj) = zytzy(ii, jj) & |
197 |
+ zytzyusn(i, j) * weighx * weighy |
198 |
zxtzy(ii, jj) = zxtzy(ii, jj) & |
199 |
+ zxtzyusn(i, j) * weighx * weighy |
200 |
ztz(ii, jj) = ztz(ii, jj) & |
201 |
+ zusn(i, j) * zusn(i, j) * weighx * weighy |
202 |
! mean |
203 |
zmea(ii, jj) = zmea(ii, jj) + zusn(i, j) * weighx * weighy |
204 |
! peacks |
205 |
zpic(ii, jj) = amax1(zpic(ii, jj), zusn(i, j)) |
206 |
! valleys |
207 |
zval(ii, jj) = amin1(zval(ii, jj), zusn(i, j)) |
208 |
ENDIF |
209 |
ENDDO |
210 |
ENDIF |
211 |
ENDDO |
212 |
ENDDO |
213 |
ENDDO |
214 |
|
215 |
if (any(weight == 0.)) stop "zero weight in grid_noro" |
216 |
|
217 |
! COMPUTE PARAMETERS NEEDED BY THE LOTT & MILLER (1997) AND |
218 |
! LOTT (1999) SSO SCHEME. |
219 |
|
220 |
zllmmea = 0. |
221 |
zllmstd = 0. |
222 |
zllmsig = 0. |
223 |
zllmgam = 0. |
224 |
zllmpic = 0. |
225 |
zllmval = 0. |
226 |
zllmthe = 0. |
227 |
zminthe = 0. |
228 |
DO ii = 1, iim + 1 |
229 |
DO jj = 1, jjm + 1 |
230 |
mask(ii, jj) = num_lan(ii, jj) / num_tot(ii, jj) |
231 |
! Mean Orography: |
232 |
zmea (ii, jj) = zmea (ii, jj) / weight(ii, jj) |
233 |
zxtzx(ii, jj) = zxtzx(ii, jj) / weight(ii, jj) |
234 |
zytzy(ii, jj) = zytzy(ii, jj) / weight(ii, jj) |
235 |
zxtzy(ii, jj) = zxtzy(ii, jj) / weight(ii, jj) |
236 |
ztz(ii, jj) = ztz(ii, jj) / weight(ii, jj) |
237 |
! Standard deviation: |
238 |
zstd(ii, jj) = sqrt(MAX(0., ztz(ii, jj) - zmea(ii, jj)**2)) |
239 |
ENDDO |
240 |
ENDDO |
241 |
|
242 |
! CORRECT VALUES OF HORIZONTAL SLOPE NEAR THE POLES: |
243 |
DO ii = 1, iim + 1 |
244 |
zxtzx(ii, 1) = zxtzx(ii, 2) |
245 |
zxtzx(ii, jjm + 1) = zxtzx(ii, jjm) |
246 |
zxtzy(ii, 1) = zxtzy(ii, 2) |
247 |
zxtzy(ii, jjm + 1) = zxtzy(ii, jjm) |
248 |
zytzy(ii, 1) = zytzy(ii, 2) |
249 |
zytzy(ii, jjm + 1) = zytzy(ii, jjm) |
250 |
ENDDO |
251 |
|
252 |
! FILTERS TO SMOOTH OUT FIELDS FOR INPUT INTO SSO SCHEME. |
253 |
|
254 |
! FIRST FILTER, MOVING AVERAGE OVER 9 POINTS. |
255 |
CALL MVA9(zmea) |
256 |
CALL MVA9(zstd) |
257 |
CALL MVA9(zpic) |
258 |
CALL MVA9(zval) |
259 |
CALL MVA9(zxtzx) |
260 |
CALL MVA9(zxtzy) |
261 |
CALL MVA9(zytzy) |
262 |
|
263 |
! Masque prenant en compte maximum de terre. On seuille à 10 % de |
264 |
! terre car en dessous les paramètres de surface n'ont pas de |
265 |
! sens. |
266 |
mask_tmp = 0. |
267 |
WHERE (mask >= 0.1) mask_tmp = 1. |
268 |
|
269 |
DO ii = 1, iim |
270 |
DO jj = 1, jjm + 1 |
271 |
! Coefficients K, L et M: |
272 |
xk = (zxtzx(ii, jj) + zytzy(ii, jj)) / 2. |
273 |
xl = (zxtzx(ii, jj) - zytzy(ii, jj)) / 2. |
274 |
xm = zxtzy(ii, jj) |
275 |
xp = xk - sqrt(xl**2 + xm**2) |
276 |
xq = xk + sqrt(xl**2 + xm**2) |
277 |
xw = 1e-8 |
278 |
if(xp.le.xw) xp = 0. |
279 |
if(xq.le.xw) xq = xw |
280 |
if(abs(xm).le.xw) xm = xw * sign(1., xm) |
281 |
! modification pour masque de terre fractionnaire |
282 |
! slope: |
283 |
zsig(ii, jj) = sqrt(xq) * mask_tmp(ii, jj) |
284 |
! isotropy: |
285 |
zgam(ii, jj) = xp / xq * mask_tmp(ii, jj) |
286 |
! angle theta: |
287 |
zthe(ii, jj) = 57.29577951 * atan2(xm, xl) / 2. * mask_tmp(ii, jj) |
288 |
zphi(ii, jj) = zmea(ii, jj) * mask_tmp(ii, jj) |
289 |
zmea(ii, jj) = zmea(ii, jj) * mask_tmp(ii, jj) |
290 |
zpic(ii, jj) = zpic(ii, jj) * mask_tmp(ii, jj) |
291 |
zval(ii, jj) = zval(ii, jj) * mask_tmp(ii, jj) |
292 |
zstd(ii, jj) = zstd(ii, jj) * mask_tmp(ii, jj) |
293 |
zllmmea = AMAX1(zmea(ii, jj), zllmmea) |
294 |
zllmstd = AMAX1(zstd(ii, jj), zllmstd) |
295 |
zllmsig = AMAX1(zsig(ii, jj), zllmsig) |
296 |
zllmgam = AMAX1(zgam(ii, jj), zllmgam) |
297 |
zllmthe = AMAX1(zthe(ii, jj), zllmthe) |
298 |
zminthe = amin1(zthe(ii, jj), zminthe) |
299 |
zllmpic = AMAX1(zpic(ii, jj), zllmpic) |
300 |
zllmval = AMAX1(zval(ii, jj), zllmval) |
301 |
ENDDO |
302 |
ENDDO |
303 |
|
304 |
print *, 'MEAN ORO: ', zllmmea |
305 |
print *, 'ST. DEV.: ', zllmstd |
306 |
print *, 'PENTE: ', zllmsig |
307 |
print *, 'ANISOTROP: ', zllmgam |
308 |
print *, 'ANGLE: ', zminthe, zllmthe |
309 |
print *, 'pic: ', zllmpic |
310 |
print *, 'val: ', zllmval |
311 |
|
312 |
! gamma and theta at 1. and 0. at poles |
313 |
zmea(iim + 1, :) = zmea(1, :) |
314 |
zphi(iim + 1, :) = zphi(1, :) |
315 |
zpic(iim + 1, :) = zpic(1, :) |
316 |
zval(iim + 1, :) = zval(1, :) |
317 |
zstd(iim + 1, :) = zstd(1, :) |
318 |
zsig(iim + 1, :) = zsig(1, :) |
319 |
zgam(iim + 1, :) = zgam(1, :) |
320 |
zthe(iim + 1, :) = zthe(1, :) |
321 |
|
322 |
zmeanor = 0. |
323 |
zmeasud = 0. |
324 |
zstdnor = 0. |
325 |
zstdsud = 0. |
326 |
zsignor = 0. |
327 |
zsigsud = 0. |
328 |
zweinor = 0. |
329 |
zweisud = 0. |
330 |
zpicnor = 0. |
331 |
zpicsud = 0. |
332 |
zvalnor = 0. |
333 |
zvalsud = 0. |
334 |
|
335 |
DO ii = 1, iim |
336 |
zweinor = zweinor + weight(ii, 1) |
337 |
zweisud = zweisud + weight(ii, jjm + 1) |
338 |
zmeanor = zmeanor + zmea(ii, 1) * weight(ii, 1) |
339 |
zmeasud = zmeasud + zmea(ii, jjm + 1) * weight(ii, jjm + 1) |
340 |
zstdnor = zstdnor + zstd(ii, 1) * weight(ii, 1) |
341 |
zstdsud = zstdsud + zstd(ii, jjm + 1) * weight(ii, jjm + 1) |
342 |
zsignor = zsignor + zsig(ii, 1) * weight(ii, 1) |
343 |
zsigsud = zsigsud + zsig(ii, jjm + 1) * weight(ii, jjm + 1) |
344 |
zpicnor = zpicnor + zpic(ii, 1) * weight(ii, 1) |
345 |
zpicsud = zpicsud + zpic(ii, jjm + 1) * weight(ii, jjm + 1) |
346 |
zvalnor = zvalnor + zval(ii, 1) * weight(ii, 1) |
347 |
zvalsud = zvalsud + zval(ii, jjm + 1) * weight(ii, jjm + 1) |
348 |
ENDDO |
349 |
|
350 |
zmea(:, 1) = zmeanor / zweinor |
351 |
zmea(:, jjm + 1) = zmeasud / zweisud |
352 |
|
353 |
zphi(:, 1) = zmeanor / zweinor |
354 |
zphi(:, jjm + 1) = zmeasud / zweisud |
355 |
|
356 |
zpic(:, 1) = zpicnor / zweinor |
357 |
zpic(:, jjm + 1) = zpicsud / zweisud |
358 |
|
359 |
zval(:, 1) = zvalnor / zweinor |
360 |
zval(:, jjm + 1) = zvalsud / zweisud |
361 |
|
362 |
zstd(:, 1) = zstdnor / zweinor |
363 |
zstd(:, jjm + 1) = zstdsud / zweisud |
364 |
|
365 |
zsig(:, 1) = zsignor / zweinor |
366 |
zsig(:, jjm + 1) = zsigsud / zweisud |
367 |
|
368 |
zgam(:, 1) = 1. |
369 |
zgam(:, jjm + 1) = 1. |
370 |
|
371 |
zthe(:, 1) = 0. |
372 |
zthe(:, jjm + 1) = 0. |
373 |
|
374 |
END SUBROUTINE grid_noro |
375 |
|
376 |
end module grid_noro_m |