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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (hide 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 guez 3 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 guez 15 ! ("pctsrf(i, :)" is the composition of the surface at horizontal
10 guez 78 ! position "i")
11 guez 3
12     private nbsrf, klon
13    
14     contains
15    
16 guez 107 SUBROUTINE etat0(phis)
17 guez 3
18 guez 90 ! From "etat0_netcdf.F", version 1.3, 2005/05/25 13:10:09
19 guez 3
20 guez 27 use caldyn0_m, only: caldyn0
21 guez 79 use comconst, only: cpp, kappa, iniconst
22 guez 27 use comgeom, only: rlatu, rlonv, rlonu, rlatv, aire_2d, apoln, apols, &
23 guez 78 cu_2d, cv_2d, inigeom
24 guez 27 use dimens_m, only: iim, jjm, llm, nqmx
25 guez 3 use dimphy, only: zmasq
26     use dimsoil, only: nsoilmx
27 guez 79 use disvert_m, only: ap, bp, preff, pa, disvert
28 guez 129 use dynetat0_m, only: day_ref, annee_ref
29 guez 27 use dynredem0_m, only: dynredem0
30     use dynredem1_m, only: dynredem1
31     use exner_hyb_m, only: exner_hyb
32 guez 43 use geopot_m, only: geopot
33 guez 3 use grid_atob, only: grille_m
34     use grid_change, only: init_dyn_phy, dyn_phy
35 guez 27 use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra
36 guez 18 use iniadvtrac_m, only: iniadvtrac
37 guez 54 use inifilr_m, only: inifilr
38 guez 67 use massdair_m, only: massdair
39 guez 48 use netcdf, only: nf90_nowrite
40 guez 68 use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, &
41     nf95_inq_varid, nf95_open
42 guez 90 use nr_util, only: pi, assert
43 guez 27 use paramet_m, only: ip1jm, ip1jmp1
44 guez 138 use phyetat0_m, only: rlat, rlon
45 guez 27 use phyredem_m, only: phyredem
46     use q_sat_m, only: q_sat
47 guez 7 use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz
48     use regr_pr_o3_m, only: regr_pr_o3
49 guez 68 use startdyn, only: start_init_dyn
50 guez 78 USE start_init_orog_m, only: start_init_orog, mask
51 guez 43 use start_init_phys_m, only: start_init_phys
52 guez 54 use start_inter_3d_m, only: start_inter_3d
53 guez 129 use temps, only: itau_phy
54 guez 99 use test_disvert_m, only: test_disvert
55 guez 129 use unit_nml_m, only: unit_nml
56 guez 3
57 guez 107 REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1)
58     ! surface geopotential, in m2 s-2
59 guez 3
60 guez 107 ! Local:
61    
62 guez 73 REAL, dimension(iim + 1, jjm + 1, llm):: ucov, t3d, teta
63 guez 65 REAL vcov(iim + 1, jjm, llm)
64 guez 3
65 guez 68 REAL q(iim + 1, jjm + 1, llm, nqmx)
66 guez 3 ! (mass fractions of trace species
67 guez 68 ! "q(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)"
68 guez 3 ! and pressure level "pls(i, j, l)".)
69    
70     real qsat(iim + 1, jjm + 1, llm) ! mass fraction of saturating water vapor
71 guez 99 REAL sn(klon)
72     REAL qsolsrf(klon, nbsrf), snsrf(klon, nbsrf)
73 guez 3 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 guez 78 real, dimension(iim + 1, jjm + 1):: zmea_2d, zstd_2d, zsig_2d, zgam_2d
83 guez 3 real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d
84 guez 73 real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, ps
85 guez 3 REAL zmea(klon), zstd(klon)
86     REAL zsig(klon), zgam(klon)
87     REAL zthe(klon)
88     REAL zpic(klon), zval(klon)
89 guez 78 REAL t_ancien(klon, llm), q_ancien(klon, llm)
90 guez 3 REAL run_off_lic_0(klon)
91     real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
92 guez 68
93 guez 90 ! D\'eclarations pour lecture glace de mer :
94 guez 68 INTEGER iml_lic, jml_lic
95     INTEGER ncid, varid
96     REAL, pointer:: dlon_lic(:), dlat_lic(:)
97 guez 3 REAL, ALLOCATABLE:: fraclic(:, :) ! fraction land ice
98 guez 68 REAL flic_tmp(iim + 1, jjm + 1) ! fraction land ice temporary
99 guez 3
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 guez 70 REAL phi(iim + 1, jjm + 1, llm)
107 guez 3 REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm)
108 guez 91 REAL w(iim + 1, jjm + 1, llm)
109 guez 3
110 guez 72 real sig1(klon, llm) ! section adiabatic updraft
111     real w01(klon, llm) ! vertical velocity within adiabatic updraft
112    
113 guez 97 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 guez 129 namelist /etat0_nml/ day_ref, annee_ref
123    
124 guez 3 !---------------------------------
125    
126     print *, "Call sequence information: etat0"
127    
128 guez 129 print *, "Enter namelist 'etat0_nml'."
129     read(unit=*, nml=etat0_nml)
130     write(unit_nml, nml=etat0_nml)
131    
132 guez 79 CALL iniconst
133 guez 3
134 guez 36 ! Construct a grid:
135    
136 guez 3 pa = 5e4
137 guez 79 CALL disvert
138 guez 99 call test_disvert
139 guez 3 CALL inigeom
140     CALL inifilr
141    
142 guez 138 rlat(1) = 90.
143     rlat(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * 180. / pi
144 guez 3 ! (with conversion to degrees)
145 guez 138 rlat(klon) = - 90.
146 guez 3
147 guez 138 rlon(1) = 0.
148     rlon(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * 180. / pi
149 guez 3 ! (with conversion to degrees)
150 guez 138 rlon(klon) = 0.
151 guez 3
152 guez 78 call start_init_orog(phis, zmea_2d, zstd_2d, zsig_2d, zgam_2d, zthe_2d, &
153     zpic_2d, zval_2d) ! also compute "mask"
154 guez 3 call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points
155 guez 15 zmasq = pack(mask, dyn_phy)
156 guez 3 PRINT *, 'Masque construit'
157    
158 guez 43 call start_init_phys(tsol_2d, qsol_2d)
159 guez 78 CALL start_init_dyn(tsol_2d, phis, ps)
160 guez 3
161     ! Compute pressure on intermediate levels:
162 guez 73 forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
163     CALL exner_hyb(ps, p3d, pks, pk)
164 guez 90 call assert(MINVAL(pk) /= MAXVAL(pk), '"pk" should not be constant')
165 guez 3
166 guez 36 pls = preff * (pk / cpp)**(1. / kappa)
167     PRINT *, "minval(pls) = ", minval(pls)
168     print *, "maxval(pls) = ", maxval(pls)
169 guez 3
170 guez 65 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 guez 3
174 guez 65 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 guez 3
178 guez 23 call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
179 guez 78 PRINT *, 'minval(t3d) = ', minval(t3d)
180 guez 36 print *, "maxval(t3d) = ", maxval(t3d)
181 guez 3
182 guez 73 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 guez 3 / apols
188     ENDDO
189    
190 guez 90 ! Calcul de l'humidit\'e \`a saturation :
191 guez 36 qsat = q_sat(t3d, pls)
192     PRINT *, "minval(qsat) = ", minval(qsat)
193     print *, "maxval(qsat) = ", maxval(qsat)
194 guez 3 IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'
195    
196     ! Water vapor:
197 guez 68 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 guez 3 DO l = 1, llm
201 guez 68 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 guez 3 ENDDO
205    
206 guez 68 q(:, :, :, 2:4) = 0. ! liquid water, radon and lead
207 guez 3
208 guez 5 if (nqmx >= 5) then
209     ! Ozone:
210 guez 7 call regr_lat_time_coefoz
211 guez 97 call regr_pr_o3(p3d, q(:, :, :, 5))
212 guez 7 ! Convert from mole fraction to mass fraction:
213 guez 78 q(:, :, :, 5) = q(:, :, :, 5) * 48. / 29.
214 guez 5 end if
215 guez 3
216 guez 13 sn = 0. ! snow
217     radsol = 0.
218     seaice = 0.
219     rugmer = 0.001
220 guez 78 zmea = pack(zmea_2d, dyn_phy)
221 guez 13 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 guez 3
228 guez 49 ! On initialise les sous-surfaces.
229 guez 3 ! Lecture du fichier glace de terre pour fixer la fraction de terre
230 guez 49 ! et de glace de terre :
231 guez 68
232 guez 48 call nf95_open("landiceref.nc", nf90_nowrite, ncid)
233 guez 68
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 guez 48 call nf95_inq_varid(ncid, 'landice', varid)
243 guez 68 ALLOCATE(fraclic(iml_lic, jml_lic))
244 guez 48 call nf95_get_var(ncid, varid, fraclic)
245 guez 68
246 guez 48 call nf95_close(ncid)
247 guez 3
248 guez 90 ! Interpolation sur la grille T du mod\`ele :
249 guez 68 PRINT *, 'Dimensions de "landiceref.nc"'
250 guez 3 print *, "iml_lic = ", iml_lic
251     print *, "jml_lic = ", jml_lic
252    
253 guez 90 ! Si les coordonn\'ees sont en degr\'es, on les transforme :
254 guez 78 IF (MAXVAL(dlon_lic) > pi) THEN
255 guez 68 dlon_lic = dlon_lic * pi / 180.
256 guez 3 ENDIF
257 guez 73 IF (maxval(dlat_lic) > pi) THEN
258 guez 68 dlat_lic = dlat_lic * pi/ 180.
259 guez 3 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 guez 68 deallocate(dlon_lic, dlat_lic) ! pointers
266    
267 guez 3 ! Passage sur la grille physique
268 guez 15 pctsrf = 0.
269 guez 3 pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)
270 guez 90 ! Ad\'equation avec le maque terre/mer
271 guez 73 WHERE (pctsrf(:, is_lic) < EPSFRA) pctsrf(:, is_lic) = 0.
272 guez 13 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 guez 3 pctsrf(:, is_ter) = 0.
278     elsewhere
279 guez 13 pctsrf(:, is_ter) = zmasq - pctsrf(:, is_lic)
280 guez 3 where (pctsrf(:, is_ter) < EPSFRA)
281     pctsrf(:, is_ter) = 0.
282 guez 13 pctsrf(:, is_lic) = zmasq
283 guez 3 end where
284     end where
285     end where
286    
287 guez 90 ! Sous-surface oc\'ean et glace de mer (pour d\'emarrer on met glace
288     ! de mer \`a 0) :
289 guez 13 pctsrf(:, is_oce) = 1. - zmasq
290 guez 3 WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.
291    
292 guez 90 ! V\'erification que somme des sous-surfaces vaut 1 :
293 guez 15 ji = count(abs(sum(pctsrf, dim = 2) - 1.) > EPSFRA)
294     IF (ji /= 0) then
295 guez 91 PRINT *, 'Bad surface percentages for ', ji, 'points'
296 guez 15 end IF
297 guez 3
298 guez 90 ! Calcul interm\'ediaire :
299 guez 3 CALL massdair(p3d, masse)
300    
301 guez 78 forall (l = 1:llm)
302 guez 3 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 guez 5 call iniadvtrac
309 guez 3 itau_phy = 0
310    
311 guez 78 CALL geopot(teta, pk , pks, phis, phi)
312 guez 73 CALL caldyn0(ucov, vcov, teta, ps, masse, pk, phis, phi, w, pbaru, &
313 guez 23 pbarv)
314 guez 129 CALL dynredem0("start.nc", day_ref, phis)
315 guez 73 CALL dynredem1("start.nc", vcov, ucov, teta, q, masse, ps, itau=0)
316 guez 3
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 guez 15 evap = 0.
328 guez 99 qsolsrf = 150.
329     tsoil = spread(spread(pack(tsol_2d, dyn_phy), 2, nsoilmx), 3, nbsrf)
330 guez 3 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 guez 13 frugs(:, is_oce) = rugmer
340 guez 78 frugs(:, is_ter) = MAX(1e-5, zstd * zsig / 2)
341     frugs(:, is_lic) = MAX(1e-5, zstd * zsig / 2)
342 guez 3 frugs(:, is_sic) = 0.001
343     fder = 0.
344     clwcon = 0.
345     rnebcon = 0.
346     ratqs = 0.
347     run_off_lic_0 = 0.
348 guez 72 sig1 = 0.
349     w01 = 0.
350 guez 3
351 guez 138 call phyredem("startphy.nc", pctsrf, tsoil(:, 1, :), tsoil, &
352 guez 99 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 guez 3
357     END SUBROUTINE etat0
358    
359     end module etat0_mod

  ViewVC Help
Powered by ViewVC 1.1.21