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

Contents of /trunk/dyn3d/etat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 277 - (show annotations)
Thu Jul 12 15:56:17 2018 UTC (5 years, 10 months ago) by guez
Original Path: trunk/dyn3d/etat0.f
File size: 11244 byte(s)
Move fxhyp and fyhyp to module dynetat0_m to avoid side effect on
variables of module dynetat0_m. A downside is that we need to link
heavyside, coefpoly and tanh_cautious into the gcm and test_fxhyp
executables.

We must move invert_zoom_x and principal_cshift to module dynetat0_m
to avoid circular dependency.

Move definition of rlatu(1) and rlatu(jjm + 1) inside fyhyp to avoid
side effect on rlatu.

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

  ViewVC Help
Powered by ViewVC 1.1.21