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

Contents of /trunk/dyn3d/etat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 276 - (show annotations)
Thu Jul 12 14:49:20 2018 UTC (5 years, 10 months ago) by guez
Original Path: trunk/dyn3d/etat0.f
File size: 11517 byte(s)
Move procedure read_serre from module read_serre_m to module
dynetat0_m, to avoid side effet on variables of module dynetat0_m.

Create procedure set_unit_nml to avoid side effect on variable of
module unit_nml_m.

Downgrade pctsrf from variable of module etat0_m to argument of etat0
and limit to avoid side effect on pctsrf.

Move variable zmasq from module dimphy to module phyetat0_m to avoid
side effect on zmasq.

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, pa, disvert
19 use dynetat0_m, only: day_ref, annee_ref, xprimp025, xprimm025, rlatu1, &
20 rlatu2, rlatu, rlatv, yprimu1, yprimu2, rlonu, rlonv, xprimu, xprimv
21 use dynredem0_m, only: dynredem0
22 use dynredem1_m, only: dynredem1
23 use exner_hyb_m, only: exner_hyb
24 use fxhyp_m, only: fxhyp
25 use fyhyp_m, only: fyhyp
26 use geopot_m, only: geopot
27 use grille_m_m, only: grille_m
28 use grid_change, only: init_dyn_phy, dyn_phy
29 use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra, nbsrf
30 use iniadvtrac_m, only: iniadvtrac
31 use inifilr_m, only: inifilr
32 use massdair_m, only: massdair
33 use netcdf, only: nf90_nowrite
34 use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, nf95_put_var, &
35 nf95_inq_varid, nf95_open
36 use nr_util, only: pi, assert
37 use phyetat0_m, only: rlat, rlon, itau_phy, zmasq
38 use phyredem0_m, only: phyredem0, ncid_restartphy
39 use phyredem_m, only: phyredem
40 use q_sat_m, only: q_sat
41 use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz
42 use regr_pr_o3_m, only: regr_pr_o3
43 use startdyn, only: start_init_dyn
44 USE start_init_orog_m, only: start_init_orog, mask
45 use start_init_phys_m, only: start_init_phys
46 use start_inter_3d_m, only: start_inter_3d
47 use test_disvert_m, only: test_disvert
48 use unit_nml_m, only: unit_nml
49
50 REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1)
51 ! surface geopotential, in m2 s-2
52
53 REAL, intent(out):: pctsrf(:, :) ! (klon, nbsrf)
54 ! "pctsrf(i, :)" is the composition of the surface at horizontal
55 ! position "i".
56
57 ! Local:
58
59 REAL, dimension(iim + 1, jjm + 1, llm):: ucov, t3d, teta
60 REAL vcov(iim + 1, jjm, llm)
61
62 REAL q(iim + 1, jjm + 1, llm, nqmx)
63 ! (mass fractions of trace species
64 ! "q(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)"
65 ! and pressure level "pls(i, j, l)".)
66
67 real qsat(iim + 1, jjm + 1, llm) ! mass fraction of saturating water vapor
68 REAL qsolsrf(klon, nbsrf), snsrf(klon, nbsrf)
69 REAL albe(klon, nbsrf), evap(klon, nbsrf)
70 REAL tsoil(klon, nsoilmx, nbsrf)
71 REAL null_array(klon)
72 REAL solsw(klon), sollw(klon)
73 !IM "slab" ocean
74 REAL frugs(klon, nbsrf), agesno(klon, nbsrf)
75 REAL rugmer(klon)
76 real, dimension(iim + 1, jjm + 1):: zmea_2d, zstd_2d, zsig_2d, zgam_2d
77 real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d
78 real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, ps
79 REAL zmea(klon), zstd(klon)
80 REAL zsig(klon), zgam(klon)
81 REAL zthe(klon)
82 REAL zpic(klon), zval(klon)
83 REAL t_ancien(klon, llm), q_ancien(klon, llm)
84 real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
85
86 ! D\'eclarations pour lecture glace de mer :
87 INTEGER iml_lic, jml_lic
88 INTEGER ncid, varid
89 REAL, ALLOCATABLE:: dlon_lic(:), dlat_lic(:)
90 REAL, ALLOCATABLE:: fraclic(:, :) ! fraction land ice
91 REAL flic_tmp(iim + 1, jjm + 1) ! fraction land ice temporary
92
93 INTEGER l, ji
94
95 REAL pk(iim + 1, jjm + 1, llm) ! fonction d'Exner aux milieux des couches
96 real pks(iim + 1, jjm + 1)
97 REAL masse(iim + 1, jjm + 1, llm)
98 REAL phi(iim + 1, jjm + 1, llm)
99 real sig1(klon, llm) ! section adiabatic updraft
100 real w01(klon, llm) ! vertical velocity within adiabatic updraft
101
102 real pls(iim + 1, jjm + 1, llm)
103 ! (pressure at mid-layer of LMDZ grid, in Pa)
104 ! "pls(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
105 ! for layer "l")
106
107 REAL p3d(iim + 1, jjm + 1, llm+1) ! pressure at layer interfaces, in Pa
108 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
109 ! for interface "l")
110
111 namelist /etat0_nml/ day_ref, annee_ref
112
113 !---------------------------------
114
115 print *, "Call sequence information: etat0"
116
117 print *, "Enter namelist 'etat0_nml'."
118 read(unit=*, nml=etat0_nml)
119 write(unit_nml, nml=etat0_nml)
120
121 CALL iniconst
122
123 ! Construct a grid:
124
125 pa = 5e4
126 CALL disvert
127 call test_disvert
128
129 CALL fyhyp(rlatu, rlatv, rlatu2, yprimu2, rlatu1, yprimu1)
130 CALL fxhyp(xprimm025, rlonv, xprimv, rlonu, xprimu, xprimp025)
131
132 rlatu(1) = pi / 2.
133 rlatu(jjm + 1) = -rlatu(1)
134
135 CALL inigeom
136 CALL inifilr
137
138 rlat(1) = 90.
139 rlat(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * 180. / pi
140 ! (with conversion to degrees)
141 rlat(klon) = - 90.
142
143 rlon(1) = 0.
144 rlon(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * 180. / pi
145 ! (with conversion to degrees)
146 rlon(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 null_array = 0.
213 rugmer = 0.001
214 zmea = pack(zmea_2d, dyn_phy)
215 zstd = pack(zstd_2d, dyn_phy)
216 zsig = pack(zsig_2d, dyn_phy)
217 zgam = pack(zgam_2d, dyn_phy)
218 zthe = pack(zthe_2d, dyn_phy)
219 zpic = pack(zpic_2d, dyn_phy)
220 zval = pack(zval_2d, dyn_phy)
221
222 ! On initialise les sous-surfaces.
223 ! Lecture du fichier glace de terre pour fixer la fraction de terre
224 ! et de glace de terre :
225
226 call nf95_open("landiceref.nc", nf90_nowrite, ncid)
227
228 call nf95_inq_varid(ncid, 'longitude', varid)
229 call nf95_gw_var(ncid, varid, dlon_lic)
230 iml_lic = size(dlon_lic)
231
232 call nf95_inq_varid(ncid, 'latitude', varid)
233 call nf95_gw_var(ncid, varid, dlat_lic)
234 jml_lic = size(dlat_lic)
235
236 call nf95_inq_varid(ncid, 'landice', varid)
237 ALLOCATE(fraclic(iml_lic, jml_lic))
238 call nf95_get_var(ncid, varid, fraclic)
239
240 call nf95_close(ncid)
241
242 ! Interpolation sur la grille T du mod\`ele :
243 PRINT *, 'Dimensions de "landiceref.nc"'
244 print *, "iml_lic = ", iml_lic
245 print *, "jml_lic = ", jml_lic
246
247 ! Si les coordonn\'ees sont en degr\'es, on les transforme :
248 IF (MAXVAL(dlon_lic) > pi) THEN
249 dlon_lic = dlon_lic * pi / 180.
250 ENDIF
251 IF (maxval(dlat_lic) > pi) THEN
252 dlat_lic = dlat_lic * pi/ 180.
253 ENDIF
254
255 flic_tmp(:iim, :) = grille_m(dlon_lic, dlat_lic, fraclic, rlonv(:iim), &
256 rlatu)
257 flic_tmp(iim + 1, :) = flic_tmp(1, :)
258
259 ! Passage sur la grille physique
260 pctsrf = 0.
261 pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)
262 ! Ad\'equation avec le maque terre/mer
263 WHERE (pctsrf(:, is_lic) < EPSFRA) pctsrf(:, is_lic) = 0.
264 WHERE (zmasq < EPSFRA) pctsrf(:, is_lic) = 0.
265 where (zmasq <= EPSFRA) pctsrf(:, is_ter) = zmasq
266 where (zmasq > EPSFRA)
267 where (pctsrf(:, is_lic) >= zmasq)
268 pctsrf(:, is_lic) = zmasq
269 pctsrf(:, is_ter) = 0.
270 elsewhere
271 pctsrf(:, is_ter) = zmasq - pctsrf(:, is_lic)
272 where (pctsrf(:, is_ter) < EPSFRA)
273 pctsrf(:, is_ter) = 0.
274 pctsrf(:, is_lic) = zmasq
275 end where
276 end where
277 end where
278
279 ! Sous-surface oc\'ean et glace de mer (pour d\'emarrer on met glace
280 ! de mer \`a 0) :
281 pctsrf(:, is_oce) = 1. - zmasq
282 WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.
283
284 ! V\'erification que la somme des sous-surfaces vaut 1 :
285 ji = count(abs(sum(pctsrf, dim = 2) - 1.) > EPSFRA)
286 IF (ji /= 0) then
287 PRINT *, 'Bad surface percentages for ', ji, 'points'
288 end IF
289
290 ! Calcul interm\'ediaire :
291 CALL massdair(p3d, masse)
292
293 forall (l = 1:llm)
294 masse(:, 1, l) = SUM(aire_2d(:iim, 1) * masse(:iim, 1, l)) / apoln
295 masse(:, jjm + 1, l) = &
296 SUM(aire_2d(:iim, jjm + 1) * masse(:iim, jjm + 1, l)) / apols
297 END forall
298
299 call iniadvtrac
300 CALL geopot(teta, pk , pks, phis, phi)
301 CALL caldyn0(ucov, vcov, teta, ps, pk, phis, phi)
302 CALL dynredem0(day_ref, phis)
303 CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = 0)
304
305 ! Initialisations :
306 snsrf = 0.
307 albe(:, is_ter) = 0.08
308 albe(:, is_lic) = 0.6
309 albe(:, is_oce) = 0.5
310 albe(:, is_sic) = 0.6
311 evap = 0.
312 qsolsrf = 150.
313 tsoil = spread(spread(pack(tsol_2d, dyn_phy), 2, nsoilmx), 3, nbsrf)
314 solsw = 165.
315 sollw = -53.
316 t_ancien = 273.15
317 q_ancien = 0.
318 agesno = 0.
319
320 frugs(:, is_oce) = rugmer
321 frugs(:, is_ter) = MAX(1e-5, zstd * zsig / 2)
322 frugs(:, is_lic) = MAX(1e-5, zstd * zsig / 2)
323 frugs(:, is_sic) = 0.001
324 clwcon = 0.
325 rnebcon = 0.
326 ratqs = 0.
327 sig1 = 0.
328 w01 = 0.
329
330 nday = 0
331 itau_phy = 0 ! side effect
332 call phyredem0
333
334 call nf95_inq_varid(ncid_restartphy, "trs", varid)
335 call nf95_put_var(ncid_restartphy, varid, null_array)
336
337 call phyredem(pctsrf, tsoil(:, 1, :), tsoil, qsolsrf, &
338 pack(qsol_2d, dyn_phy), snsrf, albe, evap, null_array, null_array, &
339 solsw, sollw, null_array, null_array, frugs, agesno, zmea, zstd, &
340 zsig, zgam, zthe, zpic, zval, t_ancien, q_ancien, rnebcon, ratqs, &
341 clwcon, null_array, sig1, w01)
342
343 END SUBROUTINE etat0
344
345 end module etat0_m

  ViewVC Help
Powered by ViewVC 1.1.21