/[lmdze]/trunk/phylmd/Orography/grid_noro_m.f
ViewVC logotype

Contents of /trunk/phylmd/Orography/grid_noro_m.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 95 - (show annotations)
Wed Apr 2 12:59:54 2014 UTC (10 years, 2 months ago) by guez
File size: 13188 byte(s)
Removed argument ps of calfis (was not done in revision 91, error in
log message).

Removed optional actual argument pkf of the call to exner_hyb before
calfis, in leapfrog. pkf was not used before the next call to
exner_hyb.

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

  ViewVC Help
Powered by ViewVC 1.1.21