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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21