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

Contents of /trunk/dyn3d/etat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (show annotations)
Wed Mar 12 21:16:36 2014 UTC (10 years, 2 months ago) by guez
File size: 11748 byte(s)
Removed procedures ini_histday, ini_histhf, write_histday and
write_histhf.

Divided file regr_pr_coefoz.f into regr_pr_av.f and
regr_pr_int.f. (Following LMDZ.) Divided module regr_pr_coefoz into
modules regr_pr_av_m and regr_pr_int_m. Renamed regr_pr_av_coefoz to
regr_pr_av and regr_pr_int_coefoz to regr_pr_int. The idea is that
those procedures are more general than Mobidic.

Removed argument dudyn of calfis and physiq. dudyn is not used either
in LMDZ. Removed computation in calfis of unused variable zpsrf (not
used either in LMDZ). Removed useless computation of dqfi in calfis
(part 62): the results were overwritten. (Same in LMDZ.)

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

  ViewVC Help
Powered by ViewVC 1.1.21