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

Annotation of /trunk/dyn3d/etat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 328 - (hide annotations)
Thu Jun 13 14:40:06 2019 UTC (4 years, 11 months ago) by guez
File size: 10656 byte(s)
Change all `.f` suffixes to `.f90`. (The opposite was done in revision
82.)  Because of change of philosopy in GNUmakefile: we already had a
rewritten rule for `.f`, so it does not make the makefile longer to
replace it by a rule for `.f90`. And it spares us options of
makedepf90 and of the compiler. Also we prepare the way for a simpler
`CMakeLists.txt`.

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

  ViewVC Help
Powered by ViewVC 1.1.21