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

Annotation of /trunk/Sources/dyn3d/etat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 107 - (hide annotations)
Thu Sep 11 15:09:15 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/dyn3d/etat0.f
File size: 11751 byte(s)
Imported procedure grilles_gcm_sub from LMDZ. Had then to transform
local variable phis of etat to argument.

Replaced calls to lnblnk by calls to trim.

Removed arguments nlat, klevel and griscal of filtreg. Replaced
integer arguments ifiltre and iaire by logical arguments direct and
intensive.

Changed default values of guide_t and guide_q to false.

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

  ViewVC Help
Powered by ViewVC 1.1.21