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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (show annotations)
Fri May 22 23:13:19 2015 UTC (9 years ago) by guez
File size: 11641 byte(s)
Moved variable nb_files from module histcom_var to module
histbeg_totreg_m.

Removed unused argument q of writehist.

No history file is created in program ce0l so there is no need to call
histclo in etat0.

In phyredem, access variables rlat and rlon directly from module
phyetat0_m instead of having them as arguments. This is clearer for
the program gcm. There are bad side effects for the program ce0l: we
have to modify the module variables rlat and rlon in procedure etat0,
and we need the additional file phyetat0.f to compile ce0l.

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

  ViewVC Help
Powered by ViewVC 1.1.21