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

Contents of /trunk/dyn3d/etat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 313 - (show annotations)
Mon Dec 10 15:54:30 2018 UTC (5 years, 5 months ago) by guez
File size: 10637 byte(s)
Remove module temps. Move variable itau_dyn from module temps to
module dynetat0_m, where it is defined.

Split module dynetat0_m into dynetat0_m and dynetat0_chosen_m. The
motivation is to create smaller modules. Procedures principal_cshift
and invert_zoomx had to stay in dynetat0_m because of circular
dependency. Now we will be able to move them away. Module variables
which are chosen by the user, not computed, in program ce0l go to
dynetat0_chosen_m: day_ref, annee_ref, clon, clat, grossismx,
grossismy, dzoomx, dzoomy, taux, tauy.

Move variable "pa" from module disvert_m to module
dynetat0_chosen_m. Define "pa" in dynetat0_chosen rather than etat0.

Define day_ref and annee_ref in procedure read_serre rather than
etat0.

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

  ViewVC Help
Powered by ViewVC 1.1.21