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

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

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 76 by guez, Fri Nov 15 18:45:49 2013 UTC revision 78 by guez, Wed Feb 5 17:51:07 2014 UTC
# Line 7  module etat0_mod Line 7  module etat0_mod
7    
8    REAL pctsrf(klon, nbsrf)    REAL pctsrf(klon, nbsrf)
9    ! ("pctsrf(i, :)" is the composition of the surface at horizontal    ! ("pctsrf(i, :)" is the composition of the surface at horizontal
10    !  position "i")    ! position "i")
11    
12    private nbsrf, klon    private nbsrf, klon
13    
# Line 20  contains Line 20  contains
20      use caldyn0_m, only: caldyn0      use caldyn0_m, only: caldyn0
21      use comconst, only: dtvr, daysec, cpp, kappa      use comconst, only: dtvr, daysec, cpp, kappa
22      use comgeom, only: rlatu, rlonv, rlonu, rlatv, aire_2d, apoln, apols, &      use comgeom, only: rlatu, rlonv, rlonu, rlatv, aire_2d, apoln, apols, &
23           cu_2d, cv_2d           cu_2d, cv_2d, inigeom
24      use conf_gcm_m, only: day_step, iphysiq, dayref, anneeref      use conf_gcm_m, only: day_step, iphysiq, dayref, anneeref
25      use dimens_m, only: iim, jjm, llm, nqmx      use dimens_m, only: iim, jjm, llm, nqmx
26      use dimphy, only: zmasq      use dimphy, only: zmasq
# Line 36  contains Line 36  contains
36      use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra      use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra
37      use iniadvtrac_m, only: iniadvtrac      use iniadvtrac_m, only: iniadvtrac
38      use inifilr_m, only: inifilr      use inifilr_m, only: inifilr
     use inigeom_m, only: inigeom  
39      use massdair_m, only: massdair      use massdair_m, only: massdair
40      use netcdf, only: nf90_nowrite      use netcdf, only: nf90_nowrite
41      use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, &      use netcdf95, only: nf95_close, nf95_get_var, nf95_gw_var, &
# Line 50  contains Line 49  contains
49      use regr_pr_o3_m, only: regr_pr_o3      use regr_pr_o3_m, only: regr_pr_o3
50      use serre, only: alphax      use serre, only: alphax
51      use startdyn, only: start_init_dyn      use startdyn, only: start_init_dyn
52      USE start_init_orog_m, only: start_init_orog, mask, phis      USE start_init_orog_m, only: start_init_orog, mask
53      use start_init_phys_m, only: start_init_phys      use start_init_phys_m, only: start_init_phys
54      use start_inter_3d_m, only: start_inter_3d      use start_inter_3d_m, only: start_inter_3d
55      use temps, only: itau_phy, annee_ref, day_ref      use temps, only: itau_phy, annee_ref, day_ref
# Line 82  contains Line 81  contains
81      real seaice(klon) ! kg m-2      real seaice(klon) ! kg m-2
82      REAL frugs(klon, nbsrf), agesno(klon, nbsrf)      REAL frugs(klon, nbsrf), agesno(klon, nbsrf)
83      REAL rugmer(klon)      REAL rugmer(klon)
84      real, dimension(iim + 1, jjm + 1):: relief, zstd_2d, zsig_2d, zgam_2d      REAL phis(iim + 1, jjm + 1) ! surface geopotential, in m2 s-2
85        real, dimension(iim + 1, jjm + 1):: zmea_2d, zstd_2d, zsig_2d, zgam_2d
86      real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d      real, dimension(iim + 1, jjm + 1):: zthe_2d, zpic_2d, zval_2d
87      real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, ps      real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, ps
88      REAL zmea(klon), zstd(klon)      REAL zmea(klon), zstd(klon)
89      REAL zsig(klon), zgam(klon)      REAL zsig(klon), zgam(klon)
90      REAL zthe(klon)      REAL zthe(klon)
91      REAL zpic(klon), zval(klon)      REAL zpic(klon), zval(klon)
92      REAL t_ancien(klon, llm), q_ancien(klon, llm)      !      REAL t_ancien(klon, llm), q_ancien(klon, llm)
93      REAL run_off_lic_0(klon)      REAL run_off_lic_0(klon)
94      real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)      real clwcon(klon, llm), rnebcon(klon, llm), ratqs(klon, llm)
95    
# Line 138  contains Line 138  contains
138      ! (with conversion to degrees)      ! (with conversion to degrees)
139      lonfi(klon) = 0.      lonfi(klon) = 0.
140    
141      call start_init_orog(relief, zstd_2d, zsig_2d, zgam_2d, zthe_2d, zpic_2d, &      call start_init_orog(phis, zmea_2d, zstd_2d, zsig_2d, zgam_2d, zthe_2d, &
142           zval_2d) ! also compute "mask" and "phis"           zpic_2d, zval_2d) ! also compute "mask"
143      call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points      call init_dyn_phy ! define the mask "dyn_phy" for distinct grid points
144      zmasq = pack(mask, dyn_phy)      zmasq = pack(mask, dyn_phy)
145      PRINT *, 'Masque construit'      PRINT *, 'Masque construit'
146    
147      call start_init_phys(tsol_2d, qsol_2d)      call start_init_phys(tsol_2d, qsol_2d)
148      CALL start_init_dyn(tsol_2d, ps)      CALL start_init_dyn(tsol_2d, phis, ps)
149    
150      ! Compute pressure on intermediate levels:      ! Compute pressure on intermediate levels:
151      forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps      forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
# Line 168  contains Line 168  contains
168      vcov(iim + 1, :, :) = vcov(1, :, :)      vcov(iim + 1, :, :) = vcov(1, :, :)
169    
170      call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)      call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
171      PRINT *,  'minval(t3d) = ', minval(t3d)      PRINT *, 'minval(t3d) = ', minval(t3d)
172      print *, "maxval(t3d) = ", maxval(t3d)      print *, "maxval(t3d) = ", maxval(t3d)
173    
174      teta(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)      teta(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)
# Line 202  contains Line 202  contains
202         call regr_lat_time_coefoz         call regr_lat_time_coefoz
203         call regr_pr_o3(q(:, :, :, 5))         call regr_pr_o3(q(:, :, :, 5))
204         ! Convert from mole fraction to mass fraction:         ! Convert from mole fraction to mass fraction:
205         q(:, :, :, 5) = q(:, :, :, 5)  * 48. / 29.         q(:, :, :, 5) = q(:, :, :, 5) * 48. / 29.
206      end if      end if
207    
208      tsol = pack(tsol_2d, dyn_phy)      tsol = pack(tsol_2d, dyn_phy)
# Line 212  contains Line 212  contains
212      tslab = 0. ! IM "slab" ocean      tslab = 0. ! IM "slab" ocean
213      seaice = 0.      seaice = 0.
214      rugmer = 0.001      rugmer = 0.001
215      zmea = pack(relief, dyn_phy)      zmea = pack(zmea_2d, dyn_phy)
216      zstd = pack(zstd_2d, dyn_phy)      zstd = pack(zstd_2d, dyn_phy)
217      zsig = pack(zsig_2d, dyn_phy)      zsig = pack(zsig_2d, dyn_phy)
218      zgam = pack(zgam_2d, dyn_phy)      zgam = pack(zgam_2d, dyn_phy)
# Line 246  contains Line 246  contains
246      print *, "jml_lic = ", jml_lic      print *, "jml_lic = ", jml_lic
247    
248      ! Si les coordonnées sont en degrés, on les transforme :      ! Si les coordonnées sont en degrés, on les transforme :
249      IF (MAXVAL(dlon_lic) > pi)  THEN      IF (MAXVAL(dlon_lic) > pi) THEN
250         dlon_lic = dlon_lic * pi / 180.         dlon_lic = dlon_lic * pi / 180.
251      ENDIF      ENDIF
252      IF (maxval(dlat_lic) > pi) THEN      IF (maxval(dlat_lic) > pi) THEN
# Line 295  contains Line 295  contains
295    
296      print *, 'ALPHAX = ', alphax      print *, 'ALPHAX = ', alphax
297    
298      forall  (l = 1:llm)      forall (l = 1:llm)
299         masse(:, 1, l) = SUM(aire_2d(:iim, 1) * masse(:iim, 1, l)) / apoln         masse(:, 1, l) = SUM(aire_2d(:iim, 1) * masse(:iim, 1, l)) / apoln
300         masse(:, jjm + 1, l) = &         masse(:, jjm + 1, l) = &
301              SUM(aire_2d(:iim, jjm + 1) * masse(:iim, jjm + 1, l)) / apols              SUM(aire_2d(:iim, jjm + 1) * masse(:iim, jjm + 1, l)) / apols
# Line 307  contains Line 307  contains
307      day_ref = dayref      day_ref = dayref
308      annee_ref = anneeref      annee_ref = anneeref
309    
310      CALL geopot(teta, pk , pks,  phis, phi)      CALL geopot(teta, pk , pks, phis, phi)
311      CALL caldyn0(ucov, vcov, teta, ps, masse, pk, phis, phi, w, pbaru, &      CALL caldyn0(ucov, vcov, teta, ps, masse, pk, phis, phi, w, pbaru, &
312           pbarv)           pbarv)
313      CALL dynredem0("start.nc", dayref, phis)      CALL dynredem0("start.nc", dayref, phis)
# Line 315  contains Line 315  contains
315    
316      ! Ecriture état initial physique:      ! Ecriture état initial physique:
317      print *, "iphysiq = ", iphysiq      print *, "iphysiq = ", iphysiq
318      phystep   = dtvr * REAL(iphysiq)      phystep = dtvr * REAL(iphysiq)
319      print *, 'phystep = ', phystep      print *, 'phystep = ', phystep
320    
321      ! Initialisations :      ! Initialisations :
# Line 350  contains Line 350  contains
350      seaice = 0.      seaice = 0.
351    
352      frugs(:, is_oce) = rugmer      frugs(:, is_oce) = rugmer
353      frugs(:, is_ter) = MAX(1.e-05, zstd * zsig / 2)      frugs(:, is_ter) = MAX(1e-5, zstd * zsig / 2)
354      frugs(:, is_lic) = MAX(1.e-05, zstd * zsig / 2)      frugs(:, is_lic) = MAX(1e-5, zstd * zsig / 2)
355      frugs(:, is_sic) = 0.001      frugs(:, is_sic) = 0.001
356      fder = 0.      fder = 0.
357      clwcon = 0.      clwcon = 0.

Legend:
Removed from v.76  
changed lines
  Added in v.78

  ViewVC Help
Powered by ViewVC 1.1.21