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

Contents of /trunk/dyn3d/etat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 99 - (show annotations)
Wed Jul 2 18:39:15 2014 UTC (9 years, 10 months ago) by guez
File size: 11743 byte(s)
Created procedure test_disvert (following LMDZ). Added procedures
hybrid and funcd in module disvert_m. Upgraded compute_ab from
internal procedure of disvert to module procedure. Added variables y,
ya in module disvert_m. Upgraded s from local variable of procedure
disvert to module variable.

Renamed allowed value of variable vert_sampling in procedure disvert
from "read" to "read_hybrid". Added possibility to read pressure
values, value "read_pressure". Replaced vertical distribution for
value "param" by the distribution "strato_correct" from LMDZ (but kept
the value "param"). In case "tropo", replaced 1 by dsigmin (following
LMDZ). In case "strato", replaced 0.3 by dsigmin (following LMDZ).

Changed computation of bp in procedure compute_ab.

Removed debugindex case in clmain. Removed useless argument rlon of
procedure clmain. Removed useless variables ytaux, ytauy of procedure
clmain.

Removed intermediary variables tsol, qsol, tsolsrf, tslab in procedure
etat0.

Removed variable ok_veget:. coupling with the model Orchid is not
possible. Removed variable ocean: modeling an ocean slab is not
possible.

Removed useless variables tmp_rriv and tmp_rcoa from module
interface_surf.

Moved initialization of variables da, mp, phi in procedure physiq to
to inside the test iflag_con >= 3.

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 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 serre, only: alphax
50 use startdyn, only: start_init_dyn
51 USE start_init_orog_m, only: start_init_orog, mask
52 use start_init_phys_m, only: start_init_phys
53 use start_inter_3d_m, only: start_inter_3d
54 use temps, only: itau_phy, annee_ref, day_ref
55 use test_disvert_m, only: test_disvert
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 sn(klon)
73 REAL 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 seaice(klon) ! kg m-2
81 REAL frugs(klon, nbsrf), agesno(klon, nbsrf)
82 REAL rugmer(klon)
83 REAL phis(iim + 1, jjm + 1) ! surface geopotential, in m2 s-2
84 real, dimension(iim + 1, jjm + 1):: zmea_2d, zstd_2d, zsig_2d, zgam_2d
85 real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d
86 real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, ps
87 REAL zmea(klon), zstd(klon)
88 REAL zsig(klon), zgam(klon)
89 REAL zthe(klon)
90 REAL zpic(klon), zval(klon)
91 REAL t_ancien(klon, llm), q_ancien(klon, llm)
92 REAL run_off_lic_0(klon)
93 real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
94
95 ! D\'eclarations pour lecture glace de mer :
96 INTEGER iml_lic, jml_lic
97 INTEGER ncid, varid
98 REAL, pointer:: dlon_lic(:), dlat_lic(:)
99 REAL, ALLOCATABLE:: fraclic(:, :) ! fraction land ice
100 REAL flic_tmp(iim + 1, jjm + 1) ! fraction land ice temporary
101
102 INTEGER l, ji
103
104 REAL pk(iim + 1, jjm + 1, llm) ! fonction d'Exner aux milieux des couches
105 real pks(iim + 1, jjm + 1)
106
107 REAL masse(iim + 1, jjm + 1, llm)
108 REAL phi(iim + 1, jjm + 1, llm)
109 REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm)
110 REAL w(iim + 1, jjm + 1, llm)
111
112 real sig1(klon, llm) ! section adiabatic updraft
113 real w01(klon, llm) ! vertical velocity within adiabatic updraft
114
115 real pls(iim + 1, jjm + 1, llm)
116 ! (pressure at mid-layer of LMDZ grid, in Pa)
117 ! "pls(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
118 ! for layer "l")
119
120 REAL p3d(iim + 1, jjm + 1, llm+1) ! pressure at layer interfaces, in Pa
121 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
122 ! for interface "l")
123
124 !---------------------------------
125
126 print *, "Call sequence information: etat0"
127
128 CALL iniconst
129
130 ! Construct a grid:
131
132 pa = 5e4
133 CALL disvert
134 call test_disvert
135 CALL inigeom
136 CALL inifilr
137
138 latfi(1) = 90.
139 latfi(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * 180. / pi
140 ! (with conversion to degrees)
141 latfi(klon) = - 90.
142
143 lonfi(1) = 0.
144 lonfi(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * 180. / pi
145 ! (with conversion to degrees)
146 lonfi(klon) = 0.
147
148 call start_init_orog(phis, zmea_2d, zstd_2d, zsig_2d, zgam_2d, zthe_2d, &
149 zpic_2d, zval_2d) ! also compute "mask"
150 call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points
151 zmasq = pack(mask, dyn_phy)
152 PRINT *, 'Masque construit'
153
154 call start_init_phys(tsol_2d, qsol_2d)
155 CALL start_init_dyn(tsol_2d, phis, ps)
156
157 ! Compute pressure on intermediate levels:
158 forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
159 CALL exner_hyb(ps, p3d, pks, pk)
160 call assert(MINVAL(pk) /= MAXVAL(pk), '"pk" should not be constant')
161
162 pls = preff * (pk / cpp)**(1. / kappa)
163 PRINT *, "minval(pls) = ", minval(pls)
164 print *, "maxval(pls) = ", maxval(pls)
165
166 call start_inter_3d('U', rlonv, rlatv, pls, ucov)
167 forall (l = 1: llm) ucov(:iim, :, l) = ucov(:iim, :, l) * cu_2d(:iim, :)
168 ucov(iim+1, :, :) = ucov(1, :, :)
169
170 call start_inter_3d('V', rlonu, rlatu(:jjm), pls(:, :jjm, :), vcov)
171 forall (l = 1: llm) vcov(:iim, :, l) = vcov(:iim, :, l) * cv_2d(:iim, :)
172 vcov(iim + 1, :, :) = vcov(1, :, :)
173
174 call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
175 PRINT *, 'minval(t3d) = ', minval(t3d)
176 print *, "maxval(t3d) = ", maxval(t3d)
177
178 teta(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)
179 teta(iim + 1, :, :) = teta(1, :, :)
180 DO l = 1, llm
181 teta(:, 1, l) = SUM(aire_2d(:, 1) * teta(:, 1, l)) / apoln
182 teta(:, jjm + 1, l) = SUM(aire_2d(:, jjm + 1) * teta(:, jjm + 1, l)) &
183 / apols
184 ENDDO
185
186 ! Calcul de l'humidit\'e \`a saturation :
187 qsat = q_sat(t3d, pls)
188 PRINT *, "minval(qsat) = ", minval(qsat)
189 print *, "maxval(qsat) = ", maxval(qsat)
190 IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'
191
192 ! Water vapor:
193 call start_inter_3d('R', rlonu, rlatv, pls, q(:, :, :, 1))
194 q(:, :, :, 1) = 0.01 * q(:, :, :, 1) * qsat
195 WHERE (q(:, :, :, 1) < 0.) q(:, :, :, 1) = 1E-10
196 DO l = 1, llm
197 q(:, 1, l, 1) = SUM(aire_2d(:, 1) * q(:, 1, l, 1)) / apoln
198 q(:, jjm + 1, l, 1) &
199 = SUM(aire_2d(:, jjm + 1) * q(:, jjm + 1, l, 1)) / apols
200 ENDDO
201
202 q(:, :, :, 2:4) = 0. ! liquid water, radon and lead
203
204 if (nqmx >= 5) then
205 ! Ozone:
206 call regr_lat_time_coefoz
207 call regr_pr_o3(p3d, q(:, :, :, 5))
208 ! Convert from mole fraction to mass fraction:
209 q(:, :, :, 5) = q(:, :, :, 5) * 48. / 29.
210 end if
211
212 sn = 0. ! snow
213 radsol = 0.
214 seaice = 0.
215 rugmer = 0.001
216 zmea = pack(zmea_2d, dyn_phy)
217 zstd = pack(zstd_2d, dyn_phy)
218 zsig = pack(zsig_2d, dyn_phy)
219 zgam = pack(zgam_2d, dyn_phy)
220 zthe = pack(zthe_2d, dyn_phy)
221 zpic = pack(zpic_2d, dyn_phy)
222 zval = pack(zval_2d, dyn_phy)
223
224 ! On initialise les sous-surfaces.
225 ! Lecture du fichier glace de terre pour fixer la fraction de terre
226 ! et de glace de terre :
227
228 call nf95_open("landiceref.nc", nf90_nowrite, ncid)
229
230 call nf95_inq_varid(ncid, 'longitude', varid)
231 call nf95_gw_var(ncid, varid, dlon_lic)
232 iml_lic = size(dlon_lic)
233
234 call nf95_inq_varid(ncid, 'latitude', varid)
235 call nf95_gw_var(ncid, varid, dlat_lic)
236 jml_lic = size(dlat_lic)
237
238 call nf95_inq_varid(ncid, 'landice', varid)
239 ALLOCATE(fraclic(iml_lic, jml_lic))
240 call nf95_get_var(ncid, varid, fraclic)
241
242 call nf95_close(ncid)
243
244 ! Interpolation sur la grille T du mod\`ele :
245 PRINT *, 'Dimensions de "landiceref.nc"'
246 print *, "iml_lic = ", iml_lic
247 print *, "jml_lic = ", jml_lic
248
249 ! Si les coordonn\'ees sont en degr\'es, on les transforme :
250 IF (MAXVAL(dlon_lic) > pi) THEN
251 dlon_lic = dlon_lic * pi / 180.
252 ENDIF
253 IF (maxval(dlat_lic) > pi) THEN
254 dlat_lic = dlat_lic * pi/ 180.
255 ENDIF
256
257 flic_tmp(:iim, :) = grille_m(dlon_lic, dlat_lic, fraclic, rlonv(:iim), &
258 rlatu)
259 flic_tmp(iim + 1, :) = flic_tmp(1, :)
260
261 deallocate(dlon_lic, dlat_lic) ! pointers
262
263 ! Passage sur la grille physique
264 pctsrf = 0.
265 pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)
266 ! Ad\'equation avec le maque terre/mer
267 WHERE (pctsrf(:, is_lic) < EPSFRA) pctsrf(:, is_lic) = 0.
268 WHERE (zmasq < EPSFRA) pctsrf(:, is_lic) = 0.
269 pctsrf(:, is_ter) = zmasq
270 where (zmasq > EPSFRA)
271 where (pctsrf(:, is_lic) >= zmasq)
272 pctsrf(:, is_lic) = zmasq
273 pctsrf(:, is_ter) = 0.
274 elsewhere
275 pctsrf(:, is_ter) = zmasq - pctsrf(:, is_lic)
276 where (pctsrf(:, is_ter) < EPSFRA)
277 pctsrf(:, is_ter) = 0.
278 pctsrf(:, is_lic) = zmasq
279 end where
280 end where
281 end where
282
283 ! Sous-surface oc\'ean et glace de mer (pour d\'emarrer on met glace
284 ! de mer \`a 0) :
285 pctsrf(:, is_oce) = 1. - zmasq
286 WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.
287
288 ! V\'erification que somme des sous-surfaces vaut 1 :
289 ji = count(abs(sum(pctsrf, dim = 2) - 1.) > EPSFRA)
290 IF (ji /= 0) then
291 PRINT *, 'Bad surface percentages for ', ji, 'points'
292 end IF
293
294 ! Calcul interm\'ediaire :
295 CALL massdair(p3d, masse)
296
297 print *, 'ALPHAX = ', alphax
298
299 forall (l = 1:llm)
300 masse(:, 1, l) = SUM(aire_2d(:iim, 1) * masse(:iim, 1, l)) / apoln
301 masse(:, jjm + 1, l) = &
302 SUM(aire_2d(:iim, jjm + 1) * masse(:iim, jjm + 1, l)) / apols
303 END forall
304
305 ! Initialisation pour traceurs:
306 call iniadvtrac
307 itau_phy = 0
308 day_ref = dayref
309 annee_ref = anneeref
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", dayref, 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", latfi, lonfi, 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 CALL histclo
357
358 END SUBROUTINE etat0
359
360 end module etat0_mod

  ViewVC Help
Powered by ViewVC 1.1.21