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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 278 - (show annotations)
Thu Jul 12 17:53:18 2018 UTC (5 years, 10 months ago) by guez
File size: 12452 byte(s)
Created procedure phyetat0_new to avoid side effects on variables of
module phyetat0_m. Had then to change test_ozonecm: read rlat from
netCDF file instead of specifying it in test_ozonecm.

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\c{}cois 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 dimensions, only: iim, jjm
32 use mva9_m, only: mva9
33 use nr_util, only: assert, pi
34
35 ! Coordinates of input field:
36 REAL, intent(in):: xdata(:) ! (iusn)
37 REAL, intent(in):: ydata(:) ! (jusn)
38
39 REAL, intent(in):: zdata(:, :) ! (iusn, jusn) input field, in m
40 REAL, intent(in):: x(:), y(:) ! coordinates of output field
41
42 ! Correlations of US Navy orography gradients:
43
44 REAL, intent(out):: zphi(:, :) ! (iim + 1, jjm + 1)
45 ! geoptential height of orography, not smoothed, in m
46
47 real, intent(out):: zmea(:, :) ! (iim + 1, jjm + 1) smoothed orography
48 real, intent(out):: zstd(:, :) ! (iim + 1, jjm + 1) Standard deviation
49 REAL, intent(out):: zsig(:, :) ! (iim + 1, jjm + 1) Slope
50 real, intent(out):: zgam(:, :) ! (iim + 1, jjm + 1) Anisotropy
51
52 real, intent(out):: zthe(:, :) ! (iim + 1, jjm + 1)
53 ! Orientation of the small axis
54
55 REAL, intent(out):: zpic(:, :) ! (iim + 1, jjm + 1) Maximum altitude
56 real, intent(out):: zval(:, :) ! (iim + 1, jjm + 1) Minimum altitude
57
58 real, intent(out):: mask(:, :) ! (iim + 1, jjm + 1) fraction of land
59
60 ! Local:
61
62 ! In this version it is assumed that the input data come from
63 ! the US Navy dataset:
64 integer, parameter:: iusn = 2160, jusn = 1080
65 integer, parameter:: iext = 216
66 REAL xusn(iusn + 2 * iext), yusn(jusn + 2)
67 REAL zusn(iusn + 2 * iext, jusn + 2) ! in m
68
69 ! Intermediate fields (correlations of orography gradient)
70 REAL, dimension(iim + 1, jjm + 1):: ztz, zxtzx, zytzy, zxtzy, weight
71
72 ! Correlations of US Navy orography gradients:
73 REAL, dimension(iusn + 2 * iext, jusn + 2):: zxtzxusn, zytzyusn, zxtzyusn
74
75 real, dimension(iim + 1, jjm + 1):: mask_tmp, num_tot, num_lan
76 REAL a(iim + 1), b(iim + 1), c(jjm + 1), d(jjm + 1)
77 real weighx, weighy, xincr, xk, xp, xm, xw, xq, xl
78 real zbordnor, zdeltax, zbordsud, zdeltay, zbordoue, zlenx, zleny, zmeasud
79 real zweinor, zweisud, zmeanor, zbordest
80 integer ii, i, jj, j
81 real, parameter:: rad = 6371229.
82
83 !-------------------------------
84
85 print *, "Call sequence information: grid_noro"
86
87 call assert((/size(xdata), size(zdata, 1)/) == iusn, "grid_noro iusn")
88 call assert((/size(ydata), size(zdata, 2)/) == jusn, "grid_noro jusn")
89
90 call assert((/size(x), size(zphi, 1), size(zmea, 1), size(zstd, 1), &
91 size(zsig, 1), size(zgam, 1), size(zthe, 1), size(zpic, 1), &
92 size(zval, 1), size(mask, 1)/) == iim + 1, "grid_noro iim")
93
94 call assert((/size(y), size(zphi, 2), size(zmea, 2), size(zstd, 2), &
95 size(zsig, 2), size(zgam, 2), size(zthe, 2), size(zpic, 2), &
96 size(zval, 2), size(mask, 2)/) == jjm + 1, "grid_noro jjm")
97
98 print *, "Parameters of subgrid-scale orography"
99 zdeltay = 2. * pi / real(jusn) * rad
100
101 ! Extension of the US Navy database for computations at boundaries:
102
103 DO j = 1, jusn
104 yusn(j + 1) = ydata(j)
105 DO i = 1, iusn
106 zusn(i + iext, j + 1) = zdata(i, j)
107 xusn(i + iext) = xdata(i)
108 ENDDO
109 DO i = 1, iext
110 zusn(i, j + 1) = zdata(iusn - iext + i, j)
111 xusn(i) = xdata(iusn - iext + i) - 2. * pi
112 zusn(iusn + iext + i, j + 1) = zdata(i, j)
113 xusn(iusn + iext + i) = xdata(i) + 2. * pi
114 ENDDO
115 ENDDO
116
117 yusn(1) = ydata(1) + (ydata(1) - ydata(2))
118 yusn(jusn + 2) = ydata(jusn) + (ydata(jusn) - ydata(jusn - 1))
119 DO i = 1, iusn / 2 + iext
120 zusn(i, 1) = zusn(i + iusn / 2, 2)
121 zusn(i + iusn / 2 + iext, 1) = zusn(i, 2)
122 zusn(i, jusn + 2) = zusn(i + iusn / 2, jusn + 1)
123 zusn(i + iusn / 2 + iext, jusn + 2) = zusn(i, jusn + 1)
124 ENDDO
125
126 ! Compute limits of model gridpoint area (regular grid)
127
128 a(1) = x(1) - (x(2) - x(1)) / 2.0
129 b(1) = (x(1) + x(2)) / 2.0
130 DO i = 2, iim
131 a(i) = b(i - 1)
132 b(i) = (x(i) + x(i + 1)) / 2.0
133 ENDDO
134 a(iim + 1) = b(iim)
135 b(iim + 1) = x(iim + 1) + (x(iim + 1) - x(iim)) / 2.0
136
137 c(1) = y(1) - (y(2) - y(1)) / 2.0
138 d(1) = (y(1) + y(2)) / 2.0
139 DO j = 2, jjm
140 c(j) = d(j - 1)
141 d(j) = (y(j) + y(j + 1)) / 2.0
142 ENDDO
143 c(jjm + 1) = d(jjm)
144 d(jjm + 1) = y(jjm + 1) + (y(jjm + 1) - y(jjm)) / 2.0
145
146 ! Initialisations :
147 weight = 0.
148 zxtzx = 0.
149 zytzy = 0.
150 zxtzy = 0.
151 ztz = 0.
152 zmea = 0.
153 zpic = - 1E10
154 zval = 1E10
155
156 ! Compute slopes correlations on US Navy grid
157
158 zytzyusn = 0.
159 zxtzxusn = 0.
160 zxtzyusn = 0.
161
162 DO j = 2, jusn + 1
163 zdeltax = zdeltay * cos(yusn(j))
164 DO i = 2, iusn + 2 * iext - 1
165 zytzyusn(i, j) = (zusn(i, j + 1) - zusn(i, j - 1))**2 / zdeltay**2
166 zxtzxusn(i, j) = (zusn(i + 1, j) - zusn(i - 1, j))**2 / zdeltax**2
167 zxtzyusn(i, j) = (zusn(i, j + 1) - zusn(i, j - 1)) / zdeltay &
168 * (zusn(i + 1, j) - zusn(i - 1, j)) / zdeltax
169 ENDDO
170 ENDDO
171
172 ! Summation over gridpoint area
173
174 zleny = pi / real(jusn) * rad
175 xincr = pi / 2. / real(jusn)
176 DO ii = 1, iim + 1
177 DO jj = 1, jjm + 1
178 num_tot(ii, jj) = 0.
179 num_lan(ii, jj) = 0.
180 DO j = 2, jusn + 1
181 zlenx = zleny * cos(yusn(j))
182 zdeltax = zdeltay * cos(yusn(j))
183 zbordnor = (c(jj) - yusn(j) + xincr) * rad
184 zbordsud = (yusn(j) - d(jj) + xincr) * rad
185 weighy = MAX(0., min(zbordnor, zbordsud, zleny))
186 IF (weighy /= 0) THEN
187 DO i = 2, iusn + 2 * iext - 1
188 zbordest = (xusn(i) - a(ii) + xincr) * rad * cos(yusn(j))
189 zbordoue = (b(ii) + xincr - xusn(i)) * rad * cos(yusn(j))
190 weighx = MAX(0., min(zbordest, zbordoue, zlenx))
191 IF (weighx /= 0) THEN
192 num_tot(ii, jj) = num_tot(ii, jj) + 1.
193 if (zusn(i, j) >= 1.) then
194 num_lan(ii, jj) = num_lan(ii, jj) + 1.
195 end if
196 weight(ii, jj) = weight(ii, jj) + weighx * weighy
197 zxtzx(ii, jj) = zxtzx(ii, jj) &
198 + zxtzxusn(i, j) * weighx * weighy
199 zytzy(ii, jj) = zytzy(ii, jj) &
200 + zytzyusn(i, j) * weighx * weighy
201 zxtzy(ii, jj) = zxtzy(ii, jj) &
202 + zxtzyusn(i, j) * weighx * weighy
203 ztz(ii, jj) = ztz(ii, jj) &
204 + zusn(i, j) * zusn(i, j) * weighx * weighy
205 ! mean
206 zmea(ii, jj) = zmea(ii, jj) + zusn(i, j) * weighx * weighy
207 ! peacks
208 zpic(ii, jj) = max(zpic(ii, jj), zusn(i, j))
209 ! valleys
210 zval(ii, jj) = min(zval(ii, jj), zusn(i, j))
211 ENDIF
212 ENDDO
213 ENDIF
214 ENDDO
215 ENDDO
216 ENDDO
217
218 if (any(weight == 0.)) then
219 print *, "zero weight in grid_noro"
220 stop 1
221 end if
222
223 ! Compute parameters needed by the Lott & Miller (1997) and Lott
224 ! (1999) subgrid-scale orographic scheme.
225
226 DO ii = 1, iim + 1
227 DO jj = 1, jjm + 1
228 mask(ii, jj) = num_lan(ii, jj) / num_tot(ii, jj)
229 ! Mean orography:
230 zmea (ii, jj) = zmea (ii, jj) / weight(ii, jj)
231 zxtzx(ii, jj) = zxtzx(ii, jj) / weight(ii, jj)
232 zytzy(ii, jj) = zytzy(ii, jj) / weight(ii, jj)
233 zxtzy(ii, jj) = zxtzy(ii, jj) / weight(ii, jj)
234 ztz(ii, jj) = ztz(ii, jj) / weight(ii, jj)
235 ! Standard deviation:
236 zstd(ii, jj) = sqrt(MAX(0., ztz(ii, jj) - zmea(ii, jj)**2))
237 ENDDO
238 ENDDO
239
240 ! Correct values of horizontal slope near the poles:
241 DO ii = 1, iim + 1
242 zxtzx(ii, 1) = zxtzx(ii, 2)
243 zxtzx(ii, jjm + 1) = zxtzx(ii, jjm)
244 zxtzy(ii, 1) = zxtzy(ii, 2)
245 zxtzy(ii, jjm + 1) = zxtzy(ii, jjm)
246 zytzy(ii, 1) = zytzy(ii, 2)
247 zytzy(ii, jjm + 1) = zytzy(ii, jjm)
248 ENDDO
249
250 ! Masque prenant en compte maximum de terre. On met un seuil \`a 10
251 ! % de terre car en dessous les param\`etres de surface n'ont pas de
252 ! sens.
253 mask_tmp = merge(1., 0., mask >= 0.1)
254
255 zphi(:iim, :) = zmea(:iim, :) * mask_tmp(:iim, :)
256 ! (zmea is not yet smoothed)
257
258 ! Filters to smooth out fields for input into subgrid-scale
259 ! orographic scheme.
260
261 ! First filter, moving average over 9 points.
262 CALL MVA9(zmea)
263 CALL MVA9(zstd)
264 CALL MVA9(zpic)
265 CALL MVA9(zval)
266 CALL MVA9(zxtzx)
267 CALL MVA9(zxtzy)
268 CALL MVA9(zytzy)
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 <= xw) xp = 0.
280 if (xq <= xw) xq = xw
281 if (abs(xm) <= 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 ENDDO
290 ENDDO
291
292 zmea(:iim, :) = zmea(:iim, :) * mask_tmp(:iim, :)
293 zpic(:iim, :) = zpic(:iim, :) * mask_tmp(:iim, :)
294 zval(:iim, :) = zval(:iim, :) * mask_tmp(:iim, :)
295 zstd(:iim, :) = zstd(:iim, :) * mask_tmp(:iim, :)
296
297 print *, 'MEAN ORO: ', MAXVAL(zmea(:iim, :))
298 print *, 'ST. DEV.: ', MAXVAL(zstd(:iim, :))
299 print *, 'PENTE: ', MAXVAL(zsig(:iim, :))
300 print *, 'ANISOTROP: ', MAXVAL(zgam(:iim, :))
301 print *, 'ANGLE: ', minval(zthe(:iim, :)), MAXVAL(zthe(:iim, :))
302 print *, 'pic: ', MAXVAL(zpic(:iim, :))
303 print *, 'val: ', MAXVAL(zval(:iim, :))
304
305 ! gamma and theta at 1. and 0. at poles
306 zmea(iim + 1, :) = zmea(1, :)
307 zphi(iim + 1, :) = zphi(1, :)
308 zpic(iim + 1, :) = zpic(1, :)
309 zval(iim + 1, :) = zval(1, :)
310 zstd(iim + 1, :) = zstd(1, :)
311 zsig(iim + 1, :) = zsig(1, :)
312 zgam(iim + 1, :) = zgam(1, :)
313 zthe(iim + 1, :) = zthe(1, :)
314
315 zweinor = sum(weight(:iim, 1))
316 zweisud = sum(weight(:iim, jjm + 1))
317 zmeanor = sum(zmea(:iim, 1) * weight(:iim, 1))
318 zmeasud = sum(zmea(:iim, jjm + 1) * weight(:iim, jjm + 1))
319
320 zmea(:, 1) = zmeanor / zweinor
321 zmea(:, jjm + 1) = zmeasud / zweisud
322
323 zphi(:, 1) = zmeanor / zweinor
324 zphi(:, jjm + 1) = zmeasud / zweisud
325
326 zpic(:, 1) = sum(zpic(:iim, 1) * weight(:iim, 1)) / zweinor
327 zpic(:, jjm + 1) = sum(zpic(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
328 / zweisud
329
330 zval(:, 1) = sum(zval(:iim, 1) * weight(:iim, 1)) / zweinor
331 zval(:, jjm + 1) = sum(zval(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
332 / zweisud
333
334 zstd(:, 1) = sum(zstd(:iim, 1) * weight(:iim, 1)) / zweinor
335 zstd(:, jjm + 1) = sum(zstd(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
336 / zweisud
337
338 zsig(:, 1) = sum(zsig(:iim, 1) * weight(:iim, 1)) / zweinor
339 zsig(:, jjm + 1) = sum(zsig(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
340 / zweisud
341
342 zgam(:, 1) = 1.
343 zgam(:, jjm + 1) = 1.
344
345 zthe(:, 1) = 0.
346 zthe(:, jjm + 1) = 0.
347
348 END SUBROUTINE grid_noro
349
350 end module grid_noro_m

  ViewVC Help
Powered by ViewVC 1.1.21