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

Contents of /trunk/dyn3d/etat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 350 - (show annotations)
Mon Dec 23 16:07:24 2019 UTC (4 years, 4 months ago) by guez
File size: 10728 byte(s)
Add argument `itau_phy_redem` to phyredem0

Add argument `itau_phy_redem` to phyredem0. Motivation: being able to
call phyredem0 with `itau_phy_redem` = 0 in program ce0l, avoiding the
side effect on nday. We can now add attribute protected to variable
nday of module `conf_gcm`. And we do not need to set `itau_phy` to 0
in procedure `phyetat0_new`.

Remove variable `prt_level` of module `conf_gcm`, which was only used
in procedure advnx. Motivation: cluttering the output is not viable,
even for debug. Better use a debugger.

Forgot to add `dyn3d/ADVN/CMakeLists.txt` in previous revision.

1 module etat0_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE etat0(phis, pctsrf)
8
9 ! From "etat0_netcdf.F", version 1.3, 2005/05/25 13:10:09
10
11 use caldyn0_m, only: caldyn0
12 use comconst, only: cpp, kappa, iniconst
13 use comgeom, only: aire_2d, apoln, apols, cu_2d, cv_2d, inigeom
14 use conf_gcm_m, only: nday
15 use dimensions, only: iim, jjm, llm, nqmx
16 use dimphy, only: klon
17 use dimsoil, only: nsoilmx
18 use disvert_m, only: ap, bp, preff, disvert
19 use dynetat0_m, only: rlatu, rlatv, rlonu, rlonv, fyhyp, fxhyp
20 use dynetat0_chosen_m, only: day_ref
21 use dynredem0_m, only: dynredem0
22 use dynredem1_m, only: dynredem1
23 use exner_hyb_m, only: exner_hyb
24 use geopot_m, only: geopot
25 use grille_m_m, only: grille_m
26 use grid_change, only: init_dyn_phy, dyn_phy
27 use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra, nbsrf
28 use infotrac_init_m, only: infotrac_init
29 use inifilr_m, only: inifilr
30 use massdair_m, only: massdair
31 use netcdf, only: nf90_nowrite
32 use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, nf95_put_var, &
33 nf95_inq_varid, nf95_open
34 use nr_util, only: pi, assert
35 use phyetat0_m, only: masque, phyetat0_new
36 use phyredem0_m, only: phyredem0, ncid_restartphy
37 use phyredem_m, only: phyredem
38 use q_sat_m, only: q_sat
39 use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz
40 use regr_pr_o3_m, only: regr_pr_o3
41 use startdyn, only: start_init_dyn
42 USE start_init_orog_m, only: start_init_orog
43 use start_init_phys_m, only: start_init_phys
44 use start_inter_3d_m, only: start_inter_3d
45 use test_disvert_m, only: test_disvert
46
47 REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1)
48 ! surface geopotential, in m2 s-2
49
50 REAL, intent(out):: pctsrf(:, :) ! (klon, nbsrf)
51 ! "pctsrf(i, :)" is the composition of the surface at horizontal
52 ! position "i".
53
54 ! Local:
55
56 REAL, dimension(iim + 1, jjm + 1, llm):: ucov, t3d, teta
57 REAL vcov(iim + 1, jjm, llm)
58
59 REAL q(iim + 1, jjm + 1, llm, nqmx)
60 ! (mass fractions of trace species
61 ! "q(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)"
62 ! and pressure level "pls(i, j, l)".)
63
64 real qsat(iim + 1, jjm + 1, llm) ! mass fraction of saturating water vapor
65 REAL qsolsrf(klon, nbsrf), fsnow(klon, nbsrf)
66 REAL falbe(klon, nbsrf)
67 REAL ftsoil(klon, nsoilmx, nbsrf)
68 REAL null_array(klon)
69 REAL solsw(klon), sollw(klon)
70 !IM "slab" ocean
71 REAL frugs(klon, nbsrf), agesno(klon, nbsrf)
72 REAL rugmer(klon)
73 real, dimension(iim + 1, jjm + 1):: zmea_2d, zstd_2d, zsig_2d, zgam_2d
74 real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d
75
76 real tsol_2d(iim + 1, jjm + 1)
77 ! both soil temperature and surface temperature, in K
78
79 real, dimension(iim + 1, jjm + 1):: qsol_2d, ps
80 REAL zmea(klon), zstd(klon)
81 REAL zsig(klon), zgam(klon)
82 REAL zthe(klon)
83 REAL zpic(klon), zval(klon)
84 REAL t_ancien(klon, llm), q_ancien(klon, llm)
85 real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
86
87 ! D\'eclarations pour lecture glace de mer :
88 INTEGER iml_lic, jml_lic
89 INTEGER ncid, varid
90 REAL, ALLOCATABLE:: dlon_lic(:), dlat_lic(:)
91 REAL, ALLOCATABLE:: landice(:, :) ! fraction land ice
92 REAL flic_tmp(iim + 1, jjm + 1) ! fraction land ice temporary
93
94 INTEGER l, ji
95
96 REAL pk(iim + 1, jjm + 1, llm) ! fonction d'Exner aux milieux des couches
97 real pks(iim + 1, jjm + 1)
98 REAL masse(iim + 1, jjm + 1, llm)
99 REAL phi(iim + 1, jjm + 1, llm)
100 real sig1(klon, llm) ! section adiabatic updraft
101 real w01(klon, llm) ! vertical velocity within adiabatic updraft
102
103 real pls(iim + 1, jjm + 1, llm)
104 ! (pressure at mid-layer of LMDZ grid, in Pa)
105 ! "pls(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
106 ! for layer "l")
107
108 REAL p3d(iim + 1, jjm + 1, llm+1) ! pressure at layer interfaces, in Pa
109 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
110 ! for interface "l")
111
112 !---------------------------------
113
114 print *, "Call sequence information: etat0"
115
116 CALL iniconst
117
118 ! Construct a grid:
119
120 CALL disvert
121 call test_disvert
122 CALL fyhyp
123 CALL fxhyp
124 CALL inigeom
125 CALL inifilr
126 call start_init_orog(phis, zmea_2d, zstd_2d, zsig_2d, zgam_2d, zthe_2d, &
127 zpic_2d, zval_2d) ! also compute "mask"
128 call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points
129 call phyetat0_new
130
131 call start_init_phys(tsol_2d, qsol_2d)
132 CALL start_init_dyn(tsol_2d, phis, ps)
133
134 ! Compute pressure on intermediate levels:
135 forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
136 CALL exner_hyb(ps, p3d, pks, pk)
137 call assert(MINVAL(pk) /= MAXVAL(pk), '"pk" should not be constant')
138
139 pls = preff * (pk / cpp)**(1. / kappa)
140 PRINT *, "minval(pls) = ", minval(pls)
141 print *, "maxval(pls) = ", maxval(pls)
142
143 call start_inter_3d('U', rlonv, rlatv, pls, ucov)
144 forall (l = 1: llm) ucov(:iim, :, l) = ucov(:iim, :, l) * cu_2d(:iim, :)
145 ucov(iim+1, :, :) = ucov(1, :, :)
146
147 call start_inter_3d('V', rlonu, rlatu(:jjm), pls(:, :jjm, :), vcov)
148 forall (l = 1: llm) vcov(:iim, :, l) = vcov(:iim, :, l) * cv_2d(:iim, :)
149 vcov(iim + 1, :, :) = vcov(1, :, :)
150
151 call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
152 PRINT *, 'minval(t3d) = ', minval(t3d)
153 print *, "maxval(t3d) = ", maxval(t3d)
154
155 teta(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)
156 teta(iim + 1, :, :) = teta(1, :, :)
157 DO l = 1, llm
158 teta(:, 1, l) = SUM(aire_2d(:, 1) * teta(:, 1, l)) / apoln
159 teta(:, jjm + 1, l) = SUM(aire_2d(:, jjm + 1) * teta(:, jjm + 1, l)) &
160 / apols
161 ENDDO
162
163 ! Calcul de l'humidit\'e \`a saturation :
164 qsat = q_sat(t3d, pls)
165 PRINT *, "minval(qsat) = ", minval(qsat)
166 print *, "maxval(qsat) = ", maxval(qsat)
167 IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'
168
169 ! Water vapor:
170 call start_inter_3d('R', rlonu, rlatv, pls, q(:, :, :, 1))
171 q(:, :, :, 1) = 0.01 * q(:, :, :, 1) * qsat
172 WHERE (q(:, :, :, 1) < 0.) q(:, :, :, 1) = 1E-10
173 DO l = 1, llm
174 q(:, 1, l, 1) = SUM(aire_2d(:, 1) * q(:, 1, l, 1)) / apoln
175 q(:, jjm + 1, l, 1) &
176 = SUM(aire_2d(:, jjm + 1) * q(:, jjm + 1, l, 1)) / apols
177 ENDDO
178
179 q(:, :, :, 2:4) = 0. ! liquid water, radon and lead
180
181 if (nqmx >= 5) then
182 ! Ozone:
183 call regr_lat_time_coefoz
184 call regr_pr_o3(p3d, q(:, :, :, 5))
185 ! Convert from mole fraction to mass fraction:
186 q(:, :, :, 5) = q(:, :, :, 5) * 48. / 29.
187 end if
188
189 null_array = 0.
190 rugmer = 0.001
191 zmea = pack(zmea_2d, dyn_phy)
192 zstd = pack(zstd_2d, dyn_phy)
193 zsig = pack(zsig_2d, dyn_phy)
194 zgam = pack(zgam_2d, dyn_phy)
195 zthe = pack(zthe_2d, dyn_phy)
196 zpic = pack(zpic_2d, dyn_phy)
197 zval = pack(zval_2d, dyn_phy)
198
199 ! On initialise les sous-surfaces.
200 ! Lecture du fichier glace de terre pour fixer la fraction de terre
201 ! et de glace de terre :
202
203 call nf95_open("landiceref.nc", nf90_nowrite, ncid)
204
205 call nf95_inq_varid(ncid, 'longitude', varid)
206 call nf95_gw_var(ncid, varid, dlon_lic)
207 iml_lic = size(dlon_lic)
208
209 call nf95_inq_varid(ncid, 'latitude', varid)
210 call nf95_gw_var(ncid, varid, dlat_lic)
211 jml_lic = size(dlat_lic)
212
213 call nf95_inq_varid(ncid, 'landice', varid)
214 ALLOCATE(landice(iml_lic, jml_lic))
215 call nf95_get_var(ncid, varid, landice)
216
217 call nf95_close(ncid)
218
219 ! Interpolation sur la grille T du mod\`ele :
220 PRINT *, 'Dimensions de "landiceref.nc"'
221 print *, "iml_lic = ", iml_lic
222 print *, "jml_lic = ", jml_lic
223
224 ! Si les coordonn\'ees sont en degr\'es, on les transforme :
225 IF (MAXVAL(dlon_lic) > pi) THEN
226 dlon_lic = dlon_lic * pi / 180.
227 ENDIF
228 IF (maxval(dlat_lic) > pi) THEN
229 dlat_lic = dlat_lic * pi/ 180.
230 ENDIF
231
232 flic_tmp(:iim, :) = grille_m(dlon_lic, dlat_lic, landice, rlonv(:iim), &
233 rlatu)
234 flic_tmp(iim + 1, :) = flic_tmp(1, :)
235
236 ! Passage sur la grille physique :
237 pctsrf = 0.
238 pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)
239
240 ! Ad\'equation avec le maque terre/mer :
241 WHERE (pctsrf(:, is_lic) < EPSFRA) pctsrf(:, is_lic) = 0.
242 WHERE (masque < EPSFRA) pctsrf(:, is_lic) = 0.
243 where (masque <= EPSFRA) pctsrf(:, is_ter) = masque
244 where (masque > EPSFRA)
245 where (pctsrf(:, is_lic) >= masque)
246 pctsrf(:, is_lic) = masque
247 pctsrf(:, is_ter) = 0.
248 elsewhere
249 pctsrf(:, is_ter) = masque - pctsrf(:, is_lic)
250 where (pctsrf(:, is_ter) < EPSFRA)
251 pctsrf(:, is_ter) = 0.
252 pctsrf(:, is_lic) = masque
253 end where
254 end where
255 end where
256
257 ! Sous-surface oc\'ean et glace de mer (pour d\'emarrer on met glace
258 ! de mer \`a 0) :
259 pctsrf(:, is_oce) = 1. - masque
260 WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.
261
262 ! V\'erification que la somme des sous-surfaces vaut 1 :
263 ji = count(abs(sum(pctsrf, dim = 2) - 1.) > EPSFRA)
264 IF (ji /= 0) then
265 PRINT *, 'Bad surface percentages for ', ji, 'points'
266 end IF
267
268 ! Calcul interm\'ediaire :
269 CALL massdair(p3d, masse)
270
271 forall (l = 1:llm)
272 masse(:, 1, l) = SUM(aire_2d(:iim, 1) * masse(:iim, 1, l)) / apoln
273 masse(:, jjm + 1, l) = &
274 SUM(aire_2d(:iim, jjm + 1) * masse(:iim, jjm + 1, l)) / apols
275 END forall
276
277 call infotrac_init
278 CALL geopot(teta, pk , pks, phis, phi)
279 CALL caldyn0(ucov, vcov, teta, ps, pk, phis, phi)
280 CALL dynredem0(day_ref, phis)
281 CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = 0)
282
283 ! Initialisations :
284 fsnow = 0.
285 falbe(:, is_ter) = 0.08
286 falbe(:, is_lic) = 0.6
287 falbe(:, is_oce) = 0.5
288 falbe(:, is_sic) = 0.6
289 qsolsrf = 150.
290 ftsoil = spread(spread(pack(tsol_2d, dyn_phy), 2, nsoilmx), 3, nbsrf)
291 solsw = 165.
292 sollw = -53.
293 t_ancien = 273.15
294 q_ancien = 0.
295 agesno = 0.
296
297 frugs(:, is_oce) = rugmer
298 frugs(:, is_ter) = MAX(1e-5, zstd * zsig / 2)
299 frugs(:, is_lic) = MAX(1e-5, zstd * zsig / 2)
300 frugs(:, is_sic) = 0.001
301 clwcon = 0.
302 rnebcon = 0.
303 ratqs = 0.
304 sig1 = 0.
305 w01 = 0.
306
307 call phyredem0(0)
308
309 call nf95_inq_varid(ncid_restartphy, "trs", varid)
310 call nf95_put_var(ncid_restartphy, varid, null_array)
311
312 call phyredem(pctsrf, ftsoil(:, 1, :), ftsoil, qsolsrf, &
313 pack(qsol_2d, dyn_phy), fsnow, falbe, null_array, null_array, solsw, &
314 sollw, null_array, null_array, frugs, agesno, zmea, zstd, zsig, zgam, &
315 zthe, zpic, zval, t_ancien, q_ancien, rnebcon, ratqs, clwcon, &
316 null_array, sig1, w01)
317
318 END SUBROUTINE etat0
319
320 end module etat0_m

  ViewVC Help
Powered by ViewVC 1.1.21