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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 343 - (show annotations)
Mon Oct 28 08:14:26 2019 UTC (4 years, 7 months ago) by guez
File size: 12494 byte(s)
Add output variables rld and rldcs

Add output variables rld and rldcs (following LMDZ).

In procedure cdrag, rename variables zdu2, ztsolv, ztvd, zri to du2,
tsolv, tvd, ri. Replace `exp(log(psol))` by psol.

In procedure `pbl_surface`, rename u, v to `u_seri`, `v_seri`.

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). Compute
17 ! all the parameters needed for the gravity wave drag code.
18
19 ! Target points are on a rectangular grid:
20 ! jjm + 1 latitudes including North and South Poles;
21 ! iim + 1 longitudes, with periodicity: longitude(iim + 1) = longitude(1)
22 ! At the poles the field value is repeated iim + 1 times.
23
24 ! The parameters a, b, c, d represent the limite of the target
25 ! gridpoint region. The means over this region are calculated from
26 ! US Navy data, ponderated by a weight proportional to the surface
27 ! occupied by the data inside the model gridpoint area. In most
28 ! circumstances, this weight is the ratio between the surface of
29 ! the US Navy gridpoint area and the surface of the model gridpoint
30 ! area. See "grid_noto.txt".
31
32 use dimensions, only: iim, jjm
33 use mva9_m, only: mva9
34 use nr_util, only: assert, pi
35
36 ! Coordinates of input field:
37 REAL, intent(in):: xdata(:) ! (iusn)
38 REAL, intent(in):: ydata(:) ! (jusn)
39
40 REAL, intent(in):: zdata(:, :) ! (iusn, jusn) input field, in m
41 REAL, intent(in):: x(:), y(:) ! coordinates of output field
42
43 ! Correlations of US Navy orography gradients:
44
45 REAL, intent(out):: zphi(:, :) ! (iim + 1, jjm + 1)
46 ! geoptential height of orography, not smoothed, in m
47
48 real, intent(out):: zmea(:, :) ! (iim + 1, jjm + 1) smoothed orography
49 real, intent(out):: zstd(:, :) ! (iim + 1, jjm + 1) Standard deviation
50 REAL, intent(out):: zsig(:, :) ! (iim + 1, jjm + 1) Slope
51 real, intent(out):: zgam(:, :) ! (iim + 1, jjm + 1) Anisotropy
52
53 real, intent(out):: zthe(:, :) ! (iim + 1, jjm + 1)
54 ! Orientation of the small axis
55
56 REAL, intent(out):: zpic(:, :) ! (iim + 1, jjm + 1) Maximum altitude
57 real, intent(out):: zval(:, :) ! (iim + 1, jjm + 1) Minimum altitude
58
59 real, intent(out):: mask(:, :) ! (iim + 1, jjm + 1) fraction of land
60
61 ! Local:
62
63 ! In this version it is assumed that the input data come from
64 ! the US Navy dataset:
65 integer, parameter:: iusn = 2160, jusn = 1080
66 integer, parameter:: iext = 216
67 REAL xusn(iusn + 2 * iext), yusn(jusn + 2)
68 REAL zusn(iusn + 2 * iext, jusn + 2) ! in m
69
70 ! Intermediate fields (correlations of orography gradient)
71 REAL, dimension(iim + 1, jjm + 1):: ztz, zxtzx, zytzy, zxtzy, weight
72
73 ! Correlations of US Navy orography gradients:
74 REAL, dimension(iusn + 2 * iext, jusn + 2):: zxtzxusn, zytzyusn, zxtzyusn
75
76 real, dimension(iim + 1, jjm + 1):: mask_tmp, num_tot, num_lan
77 REAL a(iim + 1), b(iim + 1), c(jjm + 1), d(jjm + 1)
78 real weighx, weighy, xincr, xk, xp, xm, xw, xq, xl
79 real zbordnor, zdeltax, zbordsud, zdeltay, zbordoue, zlenx, zleny, zmeasud
80 real zweinor, zweisud, zmeanor, zbordest
81 integer ii, i, jj, j
82 real, parameter:: rad = 6371229.
83
84 !-------------------------------
85
86 print *, "Call sequence information: grid_noro"
87
88 call assert((/size(xdata), size(zdata, 1)/) == iusn, "grid_noro iusn")
89 call assert((/size(ydata), size(zdata, 2)/) == jusn, "grid_noro jusn")
90
91 call assert((/size(x), size(zphi, 1), size(zmea, 1), size(zstd, 1), &
92 size(zsig, 1), size(zgam, 1), size(zthe, 1), size(zpic, 1), &
93 size(zval, 1), size(mask, 1)/) == iim + 1, "grid_noro iim")
94
95 call assert((/size(y), size(zphi, 2), size(zmea, 2), size(zstd, 2), &
96 size(zsig, 2), size(zgam, 2), size(zthe, 2), size(zpic, 2), &
97 size(zval, 2), size(mask, 2)/) == jjm + 1, "grid_noro jjm")
98
99 print *, "Parameters of subgrid-scale orography"
100 zdeltay = 2. * pi / real(jusn) * rad
101
102 ! Extension of the US Navy database for computations at boundaries:
103
104 DO j = 1, jusn
105 yusn(j + 1) = ydata(j)
106 DO i = 1, iusn
107 zusn(i + iext, j + 1) = zdata(i, j)
108 xusn(i + iext) = xdata(i)
109 ENDDO
110 DO i = 1, iext
111 zusn(i, j + 1) = zdata(iusn - iext + i, j)
112 xusn(i) = xdata(iusn - iext + i) - 2. * pi
113 zusn(iusn + iext + i, j + 1) = zdata(i, j)
114 xusn(iusn + iext + i) = xdata(i) + 2. * pi
115 ENDDO
116 ENDDO
117
118 yusn(1) = ydata(1) + (ydata(1) - ydata(2))
119 yusn(jusn + 2) = ydata(jusn) + (ydata(jusn) - ydata(jusn - 1))
120 DO i = 1, iusn / 2 + iext
121 zusn(i, 1) = zusn(i + iusn / 2, 2)
122 zusn(i + iusn / 2 + iext, 1) = zusn(i, 2)
123 zusn(i, jusn + 2) = zusn(i + iusn / 2, jusn + 1)
124 zusn(i + iusn / 2 + iext, jusn + 2) = zusn(i, jusn + 1)
125 ENDDO
126
127 ! Compute limits of model gridpoint area (regular grid)
128
129 a(1) = x(1) - (x(2) - x(1)) / 2.0
130 b(1) = (x(1) + x(2)) / 2.0
131 DO i = 2, iim
132 a(i) = b(i - 1)
133 b(i) = (x(i) + x(i + 1)) / 2.0
134 ENDDO
135 a(iim + 1) = b(iim)
136 b(iim + 1) = x(iim + 1) + (x(iim + 1) - x(iim)) / 2.0
137
138 c(1) = y(1) - (y(2) - y(1)) / 2.0
139 d(1) = (y(1) + y(2)) / 2.0
140 DO j = 2, jjm
141 c(j) = d(j - 1)
142 d(j) = (y(j) + y(j + 1)) / 2.0
143 ENDDO
144 c(jjm + 1) = d(jjm)
145 d(jjm + 1) = y(jjm + 1) + (y(jjm + 1) - y(jjm)) / 2.0
146
147 ! Initialisations :
148 weight = 0.
149 zxtzx = 0.
150 zytzy = 0.
151 zxtzy = 0.
152 ztz = 0.
153 zmea = 0.
154 zpic = - 1E10
155 zval = 1E10
156
157 ! Compute slopes correlations on US Navy grid
158
159 zytzyusn = 0.
160 zxtzxusn = 0.
161 zxtzyusn = 0.
162
163 DO j = 2, jusn + 1
164 zdeltax = zdeltay * cos(yusn(j))
165 DO i = 2, iusn + 2 * iext - 1
166 zytzyusn(i, j) = (zusn(i, j + 1) - zusn(i, j - 1))**2 / zdeltay**2
167 zxtzxusn(i, j) = (zusn(i + 1, j) - zusn(i - 1, j))**2 / zdeltax**2
168 zxtzyusn(i, j) = (zusn(i, j + 1) - zusn(i, j - 1)) / zdeltay &
169 * (zusn(i + 1, j) - zusn(i - 1, j)) / zdeltax
170 ENDDO
171 ENDDO
172
173 ! Summation over gridpoint area
174
175 zleny = pi / real(jusn) * rad
176 xincr = pi / 2. / real(jusn)
177 DO ii = 1, iim + 1
178 DO jj = 1, jjm + 1
179 num_tot(ii, jj) = 0.
180 num_lan(ii, jj) = 0.
181 DO j = 2, jusn + 1
182 zlenx = zleny * cos(yusn(j))
183 zdeltax = zdeltay * cos(yusn(j))
184 zbordnor = (c(jj) - yusn(j) + xincr) * rad
185 zbordsud = (yusn(j) - d(jj) + xincr) * rad
186 weighy = MAX(0., min(zbordnor, zbordsud, zleny))
187 IF (weighy /= 0) THEN
188 DO i = 2, iusn + 2 * iext - 1
189 zbordest = (xusn(i) - a(ii) + xincr) * rad * cos(yusn(j))
190 zbordoue = (b(ii) + xincr - xusn(i)) * rad * cos(yusn(j))
191 weighx = MAX(0., min(zbordest, zbordoue, zlenx))
192 IF (weighx /= 0) THEN
193 num_tot(ii, jj) = num_tot(ii, jj) + 1.
194 if (zusn(i, j) >= 1.) &
195 num_lan(ii, jj) = num_lan(ii, jj) + 1.
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