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

Annotation of /trunk/dyn3d/etat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 113 - (hide annotations)
Thu Sep 18 19:56:46 2014 UTC (9 years, 8 months ago) by guez
File size: 11689 byte(s)
Moved the call to read_serre out of conf_gcm so that it can be called
only in the program ce0l, not in gcm. In gcm, variables of module
serre are read from start file. Added reading of dzoomx, dzoomy, taux,
tauy from start file, in dynetat0. Those variables were written by
dynredem0 but not read.

Removed possibility fxyhypb = false, because the geometric part of the
program is such a mess. Could then remove variables transx, transy,
alphax, alphay, pxo, pyo of module serre.

Bug fix in tau2alpha: missing save attributes. The first call to
tau2alpha needs to compute dxdyu and dxdyv regardless of value of
argument type, because they will be needed for subsequent calls to
tau2alpha with various values of argument type.

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

  ViewVC Help
Powered by ViewVC 1.1.21