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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 225 - (show annotations)
Mon Oct 16 12:35:41 2017 UTC (6 years, 6 months ago) by guez
File size: 11547 byte(s)
LMDZE is now in Fortran 2003 (use of allocatable arguments).

gradsdef was not used.

Change names: [uv]10m to [uv]10m_srf in clmain, y[uv]1 to
[uv]1lay. Remove useless complication: zx_alf[12]. Do not modify
[uv]1lay after initial definition from [uv].

Add [uv]10m_srf to output.

Change names in physiq: [uv]10m to [uv]10m_srf, z[uv]10m to [uv]10m,
corresponding to NetCDF output names.

Remove unused complication couchelimite and useless variable inirnpb
in phytrac.

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

  ViewVC Help
Powered by ViewVC 1.1.21