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

Annotation of /trunk/dyn3d/etat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 279 - (hide annotations)
Fri Jul 20 14:30:23 2018 UTC (5 years, 10 months ago) by guez
Original Path: trunk/dyn3d/etat0.f
File size: 10860 byte(s)
fqcalving was saved in physiq and had intent inout in pbl_surface. So
we could set fqcalving to 0 only once per run. The point is fqcalving
must be defined everywhere for the computation of the average over all
surfaces, even values that get multiplied by pctsrf = 0. I find it
clearer to set fqcalving to 0 at every call of pbl_surface. This is
more expensive but allows to give intent out to fqcalving in
pbl_surface and remove the save attribute in physiq.

Add zxfqcalving output netCDF variable (following LMDZ).

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

  ViewVC Help
Powered by ViewVC 1.1.21