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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 163 - (show annotations)
Fri Jul 24 18:14:04 2015 UTC (8 years, 10 months ago) by guez
File size: 11693 byte(s)
In etat0, do not use masse computed by caldyn0, use masse averaged at
the poles which is computed higher in etat0 (following LMDZ). Changes
restart.nc.

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

  ViewVC Help
Powered by ViewVC 1.1.21