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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 175 - (hide annotations)
Fri Feb 5 16:02:34 2016 UTC (8 years, 3 months ago) by guez
File size: 11627 byte(s)
Added argument itau_phy to ini_histins, phyetat0, phytrac and
phyredem0. Removed variable itau_phy of module temps. Avoiding side
effect in etat0 and phyetat0. The procedures ini_histins, phyetat0,
phytrac and phyredem0 are all called by physiq so there is no
cascading variable penalty.

In procedure inifilr, made the condition on colat0 weaker to allow for
rounding error.

Removed arguments flux_o, flux_g and t_slab of clmain, flux_o and
flux_g of clqh and interfsurf_hq, tslab and seaice of phyetat0 and
phyredem. NetCDF variables TSLAB and SEAICE no longer in
restartphy.nc. All these variables were related to the not-implemented
slab ocean. seaice and tslab were just set to 0 in phyetat0 and never
used nor changed. flux_o and flux_g were computed in clmain but never
used in physiq.

Removed argument swnet of clqh. Was used only to compute a local
variable, swdown, which was not used.

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 139 use comgeom, only: aire_2d, apoln, apols, cu_2d, cv_2d, inigeom
23 guez 157 use conf_gcm_m, only: nday, day_step, iphysiq
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 139 use dynetat0_m, only: day_ref, annee_ref, xprimp025, xprimm025, rlatu1, &
29     rlatu2, rlatu, rlatv, yprimu1, yprimu2, rlonu, rlonv, xprimu, xprimv
30 guez 27 use dynredem0_m, only: dynredem0
31     use dynredem1_m, only: dynredem1
32     use exner_hyb_m, only: exner_hyb
33 guez 139 use fxhyp_m, only: fxhyp
34     use fyhyp_m, only: fyhyp
35 guez 43 use geopot_m, only: geopot
36 guez 3 use grid_atob, only: grille_m
37     use grid_change, only: init_dyn_phy, dyn_phy
38 guez 27 use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra
39 guez 18 use iniadvtrac_m, only: iniadvtrac
40 guez 54 use inifilr_m, only: inifilr
41 guez 67 use massdair_m, only: massdair
42 guez 48 use netcdf, only: nf90_nowrite
43 guez 157 use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, nf95_put_var, &
44 guez 68 nf95_inq_varid, nf95_open
45 guez 90 use nr_util, only: pi, assert
46 guez 27 use paramet_m, only: ip1jm, ip1jmp1
47 guez 138 use phyetat0_m, only: rlat, rlon
48 guez 157 use phyredem0_m, only: phyredem0, ncid_restartphy
49 guez 27 use phyredem_m, only: phyredem
50     use q_sat_m, only: q_sat
51 guez 7 use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz
52     use regr_pr_o3_m, only: regr_pr_o3
53 guez 68 use startdyn, only: start_init_dyn
54 guez 78 USE start_init_orog_m, only: start_init_orog, mask
55 guez 43 use start_init_phys_m, only: start_init_phys
56 guez 54 use start_inter_3d_m, only: start_inter_3d
57 guez 99 use test_disvert_m, only: test_disvert
58 guez 129 use unit_nml_m, only: unit_nml
59 guez 3
60 guez 107 REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1)
61     ! surface geopotential, in m2 s-2
62 guez 3
63 guez 107 ! Local:
64    
65 guez 73 REAL, dimension(iim + 1, jjm + 1, llm):: ucov, t3d, teta
66 guez 65 REAL vcov(iim + 1, jjm, llm)
67 guez 3
68 guez 68 REAL q(iim + 1, jjm + 1, llm, nqmx)
69 guez 3 ! (mass fractions of trace species
70 guez 68 ! "q(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)"
71 guez 3 ! and pressure level "pls(i, j, l)".)
72    
73     real qsat(iim + 1, jjm + 1, llm) ! mass fraction of saturating water vapor
74 guez 99 REAL qsolsrf(klon, nbsrf), snsrf(klon, nbsrf)
75 guez 3 REAL albe(klon, nbsrf), evap(klon, nbsrf)
76     REAL tsoil(klon, nsoilmx, nbsrf)
77 guez 157 REAL null_array(klon)
78     REAL solsw(klon), sollw(klon)
79 guez 3 !IM "slab" ocean
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 clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
91 guez 68
92 guez 90 ! D\'eclarations pour lecture glace de mer :
93 guez 68 INTEGER iml_lic, jml_lic
94     INTEGER ncid, varid
95     REAL, pointer:: dlon_lic(:), dlat_lic(:)
96 guez 3 REAL, ALLOCATABLE:: fraclic(:, :) ! fraction land ice
97 guez 68 REAL flic_tmp(iim + 1, jjm + 1) ! fraction land ice temporary
98 guez 3
99     INTEGER l, ji
100    
101     REAL pk(iim + 1, jjm + 1, llm) ! fonction d'Exner aux milieux des couches
102     real pks(iim + 1, jjm + 1)
103     REAL masse(iim + 1, jjm + 1, llm)
104 guez 70 REAL phi(iim + 1, jjm + 1, llm)
105 guez 72 real sig1(klon, llm) ! section adiabatic updraft
106     real w01(klon, llm) ! vertical velocity within adiabatic updraft
107    
108 guez 97 real pls(iim + 1, jjm + 1, llm)
109     ! (pressure at mid-layer of LMDZ grid, in Pa)
110     ! "pls(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
111     ! for layer "l")
112    
113     REAL p3d(iim + 1, jjm + 1, llm+1) ! pressure at layer interfaces, in Pa
114     ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
115     ! for interface "l")
116    
117 guez 129 namelist /etat0_nml/ day_ref, annee_ref
118    
119 guez 3 !---------------------------------
120    
121     print *, "Call sequence information: etat0"
122    
123 guez 129 print *, "Enter namelist 'etat0_nml'."
124     read(unit=*, nml=etat0_nml)
125     write(unit_nml, nml=etat0_nml)
126    
127 guez 79 CALL iniconst
128 guez 3
129 guez 36 ! Construct a grid:
130    
131 guez 3 pa = 5e4
132 guez 79 CALL disvert
133 guez 99 call test_disvert
134 guez 139
135     CALL fyhyp(rlatu, rlatv, rlatu2, yprimu2, rlatu1, yprimu1)
136     CALL fxhyp(xprimm025, rlonv, xprimv, rlonu, xprimu, xprimp025)
137    
138     rlatu(1) = pi / 2.
139     rlatu(jjm + 1) = -rlatu(1)
140    
141 guez 3 CALL inigeom
142     CALL inifilr
143    
144 guez 138 rlat(1) = 90.
145     rlat(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * 180. / pi
146 guez 3 ! (with conversion to degrees)
147 guez 138 rlat(klon) = - 90.
148 guez 3
149 guez 138 rlon(1) = 0.
150     rlon(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * 180. / pi
151 guez 3 ! (with conversion to degrees)
152 guez 138 rlon(klon) = 0.
153 guez 3
154 guez 78 call start_init_orog(phis, zmea_2d, zstd_2d, zsig_2d, zgam_2d, zthe_2d, &
155     zpic_2d, zval_2d) ! also compute "mask"
156 guez 3 call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points
157 guez 15 zmasq = pack(mask, dyn_phy)
158 guez 3 PRINT *, 'Masque construit'
159    
160 guez 43 call start_init_phys(tsol_2d, qsol_2d)
161 guez 78 CALL start_init_dyn(tsol_2d, phis, ps)
162 guez 3
163     ! Compute pressure on intermediate levels:
164 guez 73 forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
165     CALL exner_hyb(ps, p3d, pks, pk)
166 guez 90 call assert(MINVAL(pk) /= MAXVAL(pk), '"pk" should not be constant')
167 guez 3
168 guez 36 pls = preff * (pk / cpp)**(1. / kappa)
169     PRINT *, "minval(pls) = ", minval(pls)
170     print *, "maxval(pls) = ", maxval(pls)
171 guez 3
172 guez 65 call start_inter_3d('U', rlonv, rlatv, pls, ucov)
173     forall (l = 1: llm) ucov(:iim, :, l) = ucov(:iim, :, l) * cu_2d(:iim, :)
174     ucov(iim+1, :, :) = ucov(1, :, :)
175 guez 3
176 guez 65 call start_inter_3d('V', rlonu, rlatu(:jjm), pls(:, :jjm, :), vcov)
177     forall (l = 1: llm) vcov(:iim, :, l) = vcov(:iim, :, l) * cv_2d(:iim, :)
178     vcov(iim + 1, :, :) = vcov(1, :, :)
179 guez 3
180 guez 23 call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
181 guez 78 PRINT *, 'minval(t3d) = ', minval(t3d)
182 guez 36 print *, "maxval(t3d) = ", maxval(t3d)
183 guez 3
184 guez 73 teta(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)
185     teta(iim + 1, :, :) = teta(1, :, :)
186     DO l = 1, llm
187     teta(:, 1, l) = SUM(aire_2d(:, 1) * teta(:, 1, l)) / apoln
188     teta(:, jjm + 1, l) = SUM(aire_2d(:, jjm + 1) * teta(:, jjm + 1, l)) &
189 guez 3 / apols
190     ENDDO
191    
192 guez 90 ! Calcul de l'humidit\'e \`a saturation :
193 guez 36 qsat = q_sat(t3d, pls)
194     PRINT *, "minval(qsat) = ", minval(qsat)
195     print *, "maxval(qsat) = ", maxval(qsat)
196 guez 3 IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'
197    
198     ! Water vapor:
199 guez 68 call start_inter_3d('R', rlonu, rlatv, pls, q(:, :, :, 1))
200     q(:, :, :, 1) = 0.01 * q(:, :, :, 1) * qsat
201     WHERE (q(:, :, :, 1) < 0.) q(:, :, :, 1) = 1E-10
202 guez 3 DO l = 1, llm
203 guez 68 q(:, 1, l, 1) = SUM(aire_2d(:, 1) * q(:, 1, l, 1)) / apoln
204     q(:, jjm + 1, l, 1) &
205     = SUM(aire_2d(:, jjm + 1) * q(:, jjm + 1, l, 1)) / apols
206 guez 3 ENDDO
207    
208 guez 68 q(:, :, :, 2:4) = 0. ! liquid water, radon and lead
209 guez 3
210 guez 5 if (nqmx >= 5) then
211     ! Ozone:
212 guez 7 call regr_lat_time_coefoz
213 guez 97 call regr_pr_o3(p3d, q(:, :, :, 5))
214 guez 7 ! Convert from mole fraction to mass fraction:
215 guez 78 q(:, :, :, 5) = q(:, :, :, 5) * 48. / 29.
216 guez 5 end if
217 guez 3
218 guez 157 null_array = 0.
219 guez 13 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 163 ! V\'erification que la 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 guez 5 call iniadvtrac
308 guez 78 CALL geopot(teta, pk , pks, phis, phi)
309 guez 163 CALL caldyn0(ucov, vcov, teta, ps, pk, phis, phi)
310 guez 157 CALL dynredem0(day_ref, phis)
311     CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = 0)
312 guez 3
313     ! Initialisations :
314 guez 157 snsrf = 0.
315 guez 3 albe(:, is_ter) = 0.08
316     albe(:, is_lic) = 0.6
317     albe(:, is_oce) = 0.5
318     albe(:, is_sic) = 0.6
319 guez 15 evap = 0.
320 guez 99 qsolsrf = 150.
321     tsoil = spread(spread(pack(tsol_2d, dyn_phy), 2, nsoilmx), 3, nbsrf)
322 guez 3 solsw = 165.
323     sollw = -53.
324     t_ancien = 273.15
325     q_ancien = 0.
326     agesno = 0.
327    
328 guez 13 frugs(:, is_oce) = rugmer
329 guez 78 frugs(:, is_ter) = MAX(1e-5, zstd * zsig / 2)
330     frugs(:, is_lic) = MAX(1e-5, zstd * zsig / 2)
331 guez 3 frugs(:, is_sic) = 0.001
332     clwcon = 0.
333     rnebcon = 0.
334     ratqs = 0.
335 guez 72 sig1 = 0.
336     w01 = 0.
337 guez 3
338 guez 157 nday = 0
339 guez 175 call phyredem0(lmt_pas = day_step / iphysiq, itau_phy = 0)
340 guez 157
341     call nf95_inq_varid(ncid_restartphy, "trs", varid)
342     call nf95_put_var(ncid_restartphy, varid, null_array)
343    
344 guez 175 call phyredem(pctsrf, tsoil(:, 1, :), tsoil, qsolsrf, &
345     pack(qsol_2d, dyn_phy), snsrf, albe, evap, null_array, null_array, &
346     solsw, sollw, null_array, null_array, frugs, agesno, zmea, zstd, &
347     zsig, zgam, zthe, zpic, zval, t_ancien, q_ancien, rnebcon, ratqs, &
348     clwcon, null_array, sig1, w01)
349 guez 3
350     END SUBROUTINE etat0
351    
352     end module etat0_mod

  ViewVC Help
Powered by ViewVC 1.1.21