/[lmdze]/trunk/Sources/dyn3d/etat0.f
ViewVC logotype

Contents of /trunk/Sources/dyn3d/etat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 156 - (show annotations)
Thu Jul 16 17:39:10 2015 UTC (8 years, 10 months ago) by guez
File size: 11906 byte(s)
In procedure cltracrn, no need for local variable zx_trs, use directly
local_trs.

In (re)startphy.nc, agglomerate variables for different surface types
into a single variable with an added dimension.

In phyredem, bring together all definitions, do not use redef.

1 module etat0_mod
2
3 use indicesol, only: nbsrf
4 use dimphy, only: klon
5
6 IMPLICIT NONE
7
8 REAL pctsrf(klon, nbsrf)
9 ! ("pctsrf(i, :)" is the composition of the surface at horizontal
10 ! position "i")
11
12 private nbsrf, klon
13
14 contains
15
16 SUBROUTINE etat0(phis)
17
18 ! From "etat0_netcdf.F", version 1.3, 2005/05/25 13:10:09
19
20 use caldyn0_m, only: caldyn0
21 use comconst, only: cpp, kappa, iniconst
22 use comgeom, only: aire_2d, apoln, apols, cu_2d, cv_2d, inigeom
23 use dimens_m, only: iim, jjm, llm, nqmx
24 use dimphy, only: zmasq
25 use dimsoil, only: nsoilmx
26 use disvert_m, only: ap, bp, preff, pa, disvert
27 use dynetat0_m, only: day_ref, annee_ref, xprimp025, xprimm025, rlatu1, &
28 rlatu2, rlatu, rlatv, yprimu1, yprimu2, rlonu, rlonv, xprimu, xprimv
29 use dynredem0_m, only: dynredem0
30 use dynredem1_m, only: dynredem1
31 use exner_hyb_m, only: exner_hyb
32 use fxhyp_m, only: fxhyp
33 use fyhyp_m, only: fyhyp
34 use geopot_m, only: geopot
35 use grid_atob, only: grille_m
36 use grid_change, only: init_dyn_phy, dyn_phy
37 use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra
38 use iniadvtrac_m, only: iniadvtrac
39 use inifilr_m, only: inifilr
40 use massdair_m, only: massdair
41 use netcdf, only: nf90_nowrite
42 use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, &
43 nf95_inq_varid, nf95_open
44 use nr_util, only: pi, assert
45 use paramet_m, only: ip1jm, ip1jmp1
46 use phyetat0_m, only: rlat, rlon
47 use phyredem_m, only: phyredem
48 use q_sat_m, only: q_sat
49 use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz
50 use regr_pr_o3_m, only: regr_pr_o3
51 use startdyn, only: start_init_dyn
52 USE start_init_orog_m, only: start_init_orog, mask
53 use start_init_phys_m, only: start_init_phys
54 use start_inter_3d_m, only: start_inter_3d
55 use temps, only: itau_phy
56 use test_disvert_m, only: test_disvert
57 use unit_nml_m, only: unit_nml
58
59 REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1)
60 ! surface geopotential, in m2 s-2
61
62 ! Local:
63
64 REAL, dimension(iim + 1, jjm + 1, llm):: ucov, t3d, teta
65 REAL vcov(iim + 1, jjm, llm)
66
67 REAL q(iim + 1, jjm + 1, llm, nqmx)
68 ! (mass fractions of trace species
69 ! "q(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)"
70 ! and pressure level "pls(i, j, l)".)
71
72 real qsat(iim + 1, jjm + 1, llm) ! mass fraction of saturating water vapor
73 REAL sn(klon)
74 REAL qsolsrf(klon, nbsrf), snsrf(klon, nbsrf)
75 REAL albe(klon, nbsrf), evap(klon, nbsrf)
76 REAL tsoil(klon, nsoilmx, nbsrf)
77 REAL radsol(klon), rain_fall(klon), snow_fall(klon)
78 REAL solsw(klon), sollw(klon), fder(klon)
79 !IM "slab" ocean
80 real seaice(klon) ! kg m-2
81 REAL frugs(klon, nbsrf), agesno(klon, nbsrf)
82 REAL rugmer(klon)
83 real, dimension(iim + 1, jjm + 1):: zmea_2d, zstd_2d, zsig_2d, zgam_2d
84 real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d
85 real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, ps
86 REAL zmea(klon), zstd(klon)
87 REAL zsig(klon), zgam(klon)
88 REAL zthe(klon)
89 REAL zpic(klon), zval(klon)
90 REAL t_ancien(klon, llm), q_ancien(klon, llm)
91 REAL run_off_lic_0(klon)
92 real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
93
94 ! D\'eclarations pour lecture glace de mer :
95 INTEGER iml_lic, jml_lic
96 INTEGER ncid, varid
97 REAL, pointer:: dlon_lic(:), dlat_lic(:)
98 REAL, ALLOCATABLE:: fraclic(:, :) ! fraction land ice
99 REAL flic_tmp(iim + 1, jjm + 1) ! fraction land ice temporary
100
101 INTEGER l, ji
102
103 REAL pk(iim + 1, jjm + 1, llm) ! fonction d'Exner aux milieux des couches
104 real pks(iim + 1, jjm + 1)
105
106 REAL masse(iim + 1, jjm + 1, llm)
107 REAL phi(iim + 1, jjm + 1, llm)
108 REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm)
109 REAL w(iim + 1, jjm + 1, llm)
110
111 real sig1(klon, llm) ! section adiabatic updraft
112 real w01(klon, llm) ! vertical velocity within adiabatic updraft
113
114 real pls(iim + 1, jjm + 1, llm)
115 ! (pressure at mid-layer of LMDZ grid, in Pa)
116 ! "pls(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
117 ! for layer "l")
118
119 REAL p3d(iim + 1, jjm + 1, llm+1) ! pressure at layer interfaces, in Pa
120 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
121 ! for interface "l")
122
123 namelist /etat0_nml/ day_ref, annee_ref
124
125 !---------------------------------
126
127 print *, "Call sequence information: etat0"
128
129 print *, "Enter namelist 'etat0_nml'."
130 read(unit=*, nml=etat0_nml)
131 write(unit_nml, nml=etat0_nml)
132
133 CALL iniconst
134
135 ! Construct a grid:
136
137 pa = 5e4
138 CALL disvert
139 call test_disvert
140
141 CALL fyhyp(rlatu, rlatv, rlatu2, yprimu2, rlatu1, yprimu1)
142 CALL fxhyp(xprimm025, rlonv, xprimv, rlonu, xprimu, xprimp025)
143
144 rlatu(1) = pi / 2.
145 rlatu(jjm + 1) = -rlatu(1)
146
147 CALL inigeom
148 CALL inifilr
149
150 rlat(1) = 90.
151 rlat(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * 180. / pi
152 ! (with conversion to degrees)
153 rlat(klon) = - 90.
154
155 rlon(1) = 0.
156 rlon(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * 180. / pi
157 ! (with conversion to degrees)
158 rlon(klon) = 0.
159
160 call start_init_orog(phis, zmea_2d, zstd_2d, zsig_2d, zgam_2d, zthe_2d, &
161 zpic_2d, zval_2d) ! also compute "mask"
162 call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points
163 zmasq = pack(mask, dyn_phy)
164 PRINT *, 'Masque construit'
165
166 call start_init_phys(tsol_2d, qsol_2d)
167 CALL start_init_dyn(tsol_2d, phis, ps)
168
169 ! Compute pressure on intermediate levels:
170 forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
171 CALL exner_hyb(ps, p3d, pks, pk)
172 call assert(MINVAL(pk) /= MAXVAL(pk), '"pk" should not be constant')
173
174 pls = preff * (pk / cpp)**(1. / kappa)
175 PRINT *, "minval(pls) = ", minval(pls)
176 print *, "maxval(pls) = ", maxval(pls)
177
178 call start_inter_3d('U', rlonv, rlatv, pls, ucov)
179 forall (l = 1: llm) ucov(:iim, :, l) = ucov(:iim, :, l) * cu_2d(:iim, :)
180 ucov(iim+1, :, :) = ucov(1, :, :)
181
182 call start_inter_3d('V', rlonu, rlatu(:jjm), pls(:, :jjm, :), vcov)
183 forall (l = 1: llm) vcov(:iim, :, l) = vcov(:iim, :, l) * cv_2d(:iim, :)
184 vcov(iim + 1, :, :) = vcov(1, :, :)
185
186 call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
187 PRINT *, 'minval(t3d) = ', minval(t3d)
188 print *, "maxval(t3d) = ", maxval(t3d)
189
190 teta(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)
191 teta(iim + 1, :, :) = teta(1, :, :)
192 DO l = 1, llm
193 teta(:, 1, l) = SUM(aire_2d(:, 1) * teta(:, 1, l)) / apoln
194 teta(:, jjm + 1, l) = SUM(aire_2d(:, jjm + 1) * teta(:, jjm + 1, l)) &
195 / apols
196 ENDDO
197
198 ! Calcul de l'humidit\'e \`a saturation :
199 qsat = q_sat(t3d, pls)
200 PRINT *, "minval(qsat) = ", minval(qsat)
201 print *, "maxval(qsat) = ", maxval(qsat)
202 IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'
203
204 ! Water vapor:
205 call start_inter_3d('R', rlonu, rlatv, pls, q(:, :, :, 1))
206 q(:, :, :, 1) = 0.01 * q(:, :, :, 1) * qsat
207 WHERE (q(:, :, :, 1) < 0.) q(:, :, :, 1) = 1E-10
208 DO l = 1, llm
209 q(:, 1, l, 1) = SUM(aire_2d(:, 1) * q(:, 1, l, 1)) / apoln
210 q(:, jjm + 1, l, 1) &
211 = SUM(aire_2d(:, jjm + 1) * q(:, jjm + 1, l, 1)) / apols
212 ENDDO
213
214 q(:, :, :, 2:4) = 0. ! liquid water, radon and lead
215
216 if (nqmx >= 5) then
217 ! Ozone:
218 call regr_lat_time_coefoz
219 call regr_pr_o3(p3d, q(:, :, :, 5))
220 ! Convert from mole fraction to mass fraction:
221 q(:, :, :, 5) = q(:, :, :, 5) * 48. / 29.
222 end if
223
224 sn = 0. ! snow
225 radsol = 0.
226 seaice = 0.
227 rugmer = 0.001
228 zmea = pack(zmea_2d, dyn_phy)
229 zstd = pack(zstd_2d, dyn_phy)
230 zsig = pack(zsig_2d, dyn_phy)
231 zgam = pack(zgam_2d, dyn_phy)
232 zthe = pack(zthe_2d, dyn_phy)
233 zpic = pack(zpic_2d, dyn_phy)
234 zval = pack(zval_2d, dyn_phy)
235
236 ! On initialise les sous-surfaces.
237 ! Lecture du fichier glace de terre pour fixer la fraction de terre
238 ! et de glace de terre :
239
240 call nf95_open("landiceref.nc", nf90_nowrite, ncid)
241
242 call nf95_inq_varid(ncid, 'longitude', varid)
243 call nf95_gw_var(ncid, varid, dlon_lic)
244 iml_lic = size(dlon_lic)
245
246 call nf95_inq_varid(ncid, 'latitude', varid)
247 call nf95_gw_var(ncid, varid, dlat_lic)
248 jml_lic = size(dlat_lic)
249
250 call nf95_inq_varid(ncid, 'landice', varid)
251 ALLOCATE(fraclic(iml_lic, jml_lic))
252 call nf95_get_var(ncid, varid, fraclic)
253
254 call nf95_close(ncid)
255
256 ! Interpolation sur la grille T du mod\`ele :
257 PRINT *, 'Dimensions de "landiceref.nc"'
258 print *, "iml_lic = ", iml_lic
259 print *, "jml_lic = ", jml_lic
260
261 ! Si les coordonn\'ees sont en degr\'es, on les transforme :
262 IF (MAXVAL(dlon_lic) > pi) THEN
263 dlon_lic = dlon_lic * pi / 180.
264 ENDIF
265 IF (maxval(dlat_lic) > pi) THEN
266 dlat_lic = dlat_lic * pi/ 180.
267 ENDIF
268
269 flic_tmp(:iim, :) = grille_m(dlon_lic, dlat_lic, fraclic, rlonv(:iim), &
270 rlatu)
271 flic_tmp(iim + 1, :) = flic_tmp(1, :)
272
273 deallocate(dlon_lic, dlat_lic) ! pointers
274
275 ! Passage sur la grille physique
276 pctsrf = 0.
277 pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)
278 ! Ad\'equation avec le maque terre/mer
279 WHERE (pctsrf(:, is_lic) < EPSFRA) pctsrf(:, is_lic) = 0.
280 WHERE (zmasq < EPSFRA) pctsrf(:, is_lic) = 0.
281 pctsrf(:, is_ter) = zmasq
282 where (zmasq > EPSFRA)
283 where (pctsrf(:, is_lic) >= zmasq)
284 pctsrf(:, is_lic) = zmasq
285 pctsrf(:, is_ter) = 0.
286 elsewhere
287 pctsrf(:, is_ter) = zmasq - pctsrf(:, is_lic)
288 where (pctsrf(:, is_ter) < EPSFRA)
289 pctsrf(:, is_ter) = 0.
290 pctsrf(:, is_lic) = zmasq
291 end where
292 end where
293 end where
294
295 ! Sous-surface oc\'ean et glace de mer (pour d\'emarrer on met glace
296 ! de mer \`a 0) :
297 pctsrf(:, is_oce) = 1. - zmasq
298 WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.
299
300 ! V\'erification que somme des sous-surfaces vaut 1 :
301 ji = count(abs(sum(pctsrf, dim = 2) - 1.) > EPSFRA)
302 IF (ji /= 0) then
303 PRINT *, 'Bad surface percentages for ', ji, 'points'
304 end IF
305
306 ! Calcul interm\'ediaire :
307 CALL massdair(p3d, masse)
308
309 forall (l = 1:llm)
310 masse(:, 1, l) = SUM(aire_2d(:iim, 1) * masse(:iim, 1, l)) / apoln
311 masse(:, jjm + 1, l) = &
312 SUM(aire_2d(:iim, jjm + 1) * masse(:iim, jjm + 1, l)) / apols
313 END forall
314
315 ! Initialisation pour traceurs:
316 call iniadvtrac
317 itau_phy = 0
318
319 CALL geopot(teta, pk , pks, phis, phi)
320 CALL caldyn0(ucov, vcov, teta, ps, masse, pk, phis, phi, w, pbaru, &
321 pbarv)
322 CALL dynredem0("start.nc", day_ref, phis)
323 CALL dynredem1("start.nc", vcov, ucov, teta, q, masse, ps, itau=0)
324
325 ! Initialisations :
326 snsrf(:, is_ter) = sn
327 snsrf(:, is_lic) = sn
328 snsrf(:, is_oce) = sn
329 snsrf(:, is_sic) = sn
330 albe(:, is_ter) = 0.08
331 albe(:, is_lic) = 0.6
332 albe(:, is_oce) = 0.5
333 albe(:, is_sic) = 0.6
334 evap = 0.
335 qsolsrf = 150.
336 tsoil = spread(spread(pack(tsol_2d, dyn_phy), 2, nsoilmx), 3, nbsrf)
337 rain_fall = 0.
338 snow_fall = 0.
339 solsw = 165.
340 sollw = -53.
341 t_ancien = 273.15
342 q_ancien = 0.
343 agesno = 0.
344 seaice = 0.
345
346 frugs(:, is_oce) = rugmer
347 frugs(:, is_ter) = MAX(1e-5, zstd * zsig / 2)
348 frugs(:, is_lic) = MAX(1e-5, zstd * zsig / 2)
349 frugs(:, is_sic) = 0.001
350 fder = 0.
351 clwcon = 0.
352 rnebcon = 0.
353 ratqs = 0.
354 run_off_lic_0 = 0.
355 sig1 = 0.
356 w01 = 0.
357
358 call phyredem("startphy.nc", pctsrf, tsoil(:, 1, :), tsoil, &
359 tsoil(:, 1, is_oce), seaice, qsolsrf, pack(qsol_2d, dyn_phy), snsrf, &
360 albe, evap, rain_fall, snow_fall, solsw, sollw, fder, radsol, frugs, &
361 agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, t_ancien, &
362 q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01)
363
364 END SUBROUTINE etat0
365
366 end module etat0_mod

  ViewVC Help
Powered by ViewVC 1.1.21