/[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 326 - (show annotations)
Mon Jun 10 00:29:10 2019 UTC (4 years, 11 months ago) by guez
File size: 12422 byte(s)
GNUmakefile does not depend on `settings.mk`.

Rename argument precip1 of `cv_driver` to rain, matching the name of the
corresponding actual argument in concvl.

Add Netcdf variable pmflxr to history file, following LMDZ.

In procedure phyetat0, no need for intermediary variable fractint.

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.) &
194 num_lan(ii, jj) = num_lan(ii, jj) + 1.
195 weight(ii, jj) = weight(ii, jj) + weighx * weighy
196 zxtzx(ii, jj) = zxtzx(ii, jj) &
197 + zxtzxusn(i, j) * weighx * weighy
198 zytzy(ii, jj) = zytzy(ii, jj) &
199 + zytzyusn(i, j) * weighx * weighy
200 zxtzy(ii, jj) = zxtzy(ii, jj) &
201 + zxtzyusn(i, j) * weighx * weighy
202 ztz(ii, jj) = ztz(ii, jj) &
203 + zusn(i, j) * zusn(i, j) * weighx * weighy
204 ! mean
205 zmea(ii, jj) = zmea(ii, jj) + zusn(i, j) * weighx * weighy
206 ! peacks
207 zpic(ii, jj) = max(zpic(ii, jj), zusn(i, j))
208 ! valleys
209 zval(ii, jj) = min(zval(ii, jj), zusn(i, j))
210 ENDIF
211 ENDDO
212 ENDIF
213 ENDDO
214 ENDDO
215 ENDDO
216
217 if (any(weight == 0.)) then
218 print *, "zero weight in grid_noro"
219 stop 1
220 end if
221
222 ! Compute parameters needed by the Lott & Miller (1997) and Lott
223 ! (1999) subgrid-scale orographic scheme.
224
225 DO ii = 1, iim + 1
226 DO jj = 1, jjm + 1
227 mask(ii, jj) = num_lan(ii, jj) / num_tot(ii, jj)
228 ! Mean orography:
229 zmea (ii, jj) = zmea (ii, jj) / weight(ii, jj)
230 zxtzx(ii, jj) = zxtzx(ii, jj) / weight(ii, jj)
231 zytzy(ii, jj) = zytzy(ii, jj) / weight(ii, jj)
232 zxtzy(ii, jj) = zxtzy(ii, jj) / weight(ii, jj)
233 ztz(ii, jj) = ztz(ii, jj) / weight(ii, jj)
234 ! Standard deviation:
235 zstd(ii, jj) = sqrt(MAX(0., ztz(ii, jj) - zmea(ii, jj)**2))
236 ENDDO
237 ENDDO
238
239 ! Correct values of horizontal slope near the poles:
240 DO ii = 1, iim + 1
241 zxtzx(ii, 1) = zxtzx(ii, 2)
242 zxtzx(ii, jjm + 1) = zxtzx(ii, jjm)
243 zxtzy(ii, 1) = zxtzy(ii, 2)
244 zxtzy(ii, jjm + 1) = zxtzy(ii, jjm)
245 zytzy(ii, 1) = zytzy(ii, 2)
246 zytzy(ii, jjm + 1) = zytzy(ii, jjm)
247 ENDDO
248
249 ! Masque prenant en compte maximum de terre. On met un seuil \`a 10
250 ! % de terre car en dessous les param\`etres de surface n'ont pas de
251 ! sens.
252 mask_tmp = merge(1., 0., mask >= 0.1)
253
254 zphi(:iim, :) = zmea(:iim, :) * mask_tmp(:iim, :)
255 ! (zmea is not yet smoothed)
256
257 ! Filters to smooth out fields for input into subgrid-scale
258 ! orographic scheme.
259
260 ! First filter, moving average over 9 points.
261 CALL MVA9(zmea)
262 CALL MVA9(zstd)
263 CALL MVA9(zpic)
264 CALL MVA9(zval)
265 CALL MVA9(zxtzx)
266 CALL MVA9(zxtzy)
267 CALL MVA9(zytzy)
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 <= xw) xp = 0.
279 if (xq <= xw) xq = xw
280 if (abs(xm) <= 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 ENDDO
289 ENDDO
290
291 zmea(:iim, :) = zmea(:iim, :) * mask_tmp(:iim, :)
292 zpic(:iim, :) = zpic(:iim, :) * mask_tmp(:iim, :)
293 zval(:iim, :) = zval(:iim, :) * mask_tmp(:iim, :)
294 zstd(:iim, :) = zstd(:iim, :) * mask_tmp(:iim, :)
295
296 print *, 'MEAN ORO: ', MAXVAL(zmea(:iim, :))
297 print *, 'ST. DEV.: ', MAXVAL(zstd(:iim, :))
298 print *, 'PENTE: ', MAXVAL(zsig(:iim, :))
299 print *, 'ANISOTROP: ', MAXVAL(zgam(:iim, :))
300 print *, 'ANGLE: ', minval(zthe(:iim, :)), MAXVAL(zthe(:iim, :))
301 print *, 'pic: ', MAXVAL(zpic(:iim, :))
302 print *, 'val: ', MAXVAL(zval(:iim, :))
303
304 ! gamma and theta at 1. and 0. at poles
305 zmea(iim + 1, :) = zmea(1, :)
306 zphi(iim + 1, :) = zphi(1, :)
307 zpic(iim + 1, :) = zpic(1, :)
308 zval(iim + 1, :) = zval(1, :)
309 zstd(iim + 1, :) = zstd(1, :)
310 zsig(iim + 1, :) = zsig(1, :)
311 zgam(iim + 1, :) = zgam(1, :)
312 zthe(iim + 1, :) = zthe(1, :)
313
314 zweinor = sum(weight(:iim, 1))
315 zweisud = sum(weight(:iim, jjm + 1))
316 zmeanor = sum(zmea(:iim, 1) * weight(:iim, 1))
317 zmeasud = sum(zmea(:iim, jjm + 1) * weight(:iim, jjm + 1))
318
319 zmea(:, 1) = zmeanor / zweinor
320 zmea(:, jjm + 1) = zmeasud / zweisud
321
322 zphi(:, 1) = zmeanor / zweinor
323 zphi(:, jjm + 1) = zmeasud / zweisud
324
325 zpic(:, 1) = sum(zpic(:iim, 1) * weight(:iim, 1)) / zweinor
326 zpic(:, jjm + 1) = sum(zpic(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
327 / zweisud
328
329 zval(:, 1) = sum(zval(:iim, 1) * weight(:iim, 1)) / zweinor
330 zval(:, jjm + 1) = sum(zval(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
331 / zweisud
332
333 zstd(:, 1) = sum(zstd(:iim, 1) * weight(:iim, 1)) / zweinor
334 zstd(:, jjm + 1) = sum(zstd(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
335 / zweisud
336
337 zsig(:, 1) = sum(zsig(:iim, 1) * weight(:iim, 1)) / zweinor
338 zsig(:, jjm + 1) = sum(zsig(:iim, jjm + 1) * weight(:iim, jjm + 1)) &
339 / zweisud
340
341 zgam(:, 1) = 1.
342 zgam(:, jjm + 1) = 1.
343
344 zthe(:, 1) = 0.
345 zthe(:, jjm + 1) = 0.
346
347 END SUBROUTINE grid_noro
348
349 end module grid_noro_m

  ViewVC Help
Powered by ViewVC 1.1.21