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

Diff of /trunk/libf/dyn3d/etat0.f90

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

revision 12 by guez, Mon Jul 21 16:05:07 2008 UTC revision 66 by guez, Thu Sep 20 13:00:41 2012 UTC
# Line 6  module etat0_mod Line 6  module etat0_mod
6    IMPLICIT NONE    IMPLICIT NONE
7    
8    REAL pctsrf(klon, nbsrf)    REAL pctsrf(klon, nbsrf)
9      ! ("pctsrf(i, :)" is the composition of the surface at horizontal
10      !  position "i")
11    
12    private nbsrf, klon    private nbsrf, klon
13    
# Line 15  contains Line 17  contains
17    
18      ! From "etat0_netcdf.F", version 1.3 2005/05/25 13:10:09      ! From "etat0_netcdf.F", version 1.3 2005/05/25 13:10:09
19    
20      ! This subroutine creates "masque".      use caldyn0_m, only: caldyn0
21        use comconst, only: dtvr, daysec, cpp, kappa
     USE ioipsl, only: flinget, flinclo, flinopen_nozoom, flininfo, histclo  
   
     USE start_init_orog_m, only: start_init_orog, masque, phis  
     use start_init_phys_m, only: qsol_2d  
     use startdyn, only: start_inter_3d, start_init_dyn  
     use dimens_m, only: iim, jjm, llm, nqmx  
     use paramet_m, only: ip1jm, ip1jmp1  
     use comconst, only: dtvr, daysec, cpp, kappa, pi  
     use comdissnew, only: lstardis, nitergdiv, nitergrot, niterh, &  
          tetagdiv, tetagrot, tetatemp  
     use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra  
     use comvert, only: ap, bp, preff, pa  
     use dimphy, only: zmasq  
     use conf_gcm_m, only: day_step, iphysiq, dayref, anneeref  
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
24      use serre, only: alphax      use disvert_m, only: ap, bp, preff, pa
25        use conf_gcm_m, only: day_step, iphysiq, dayref, anneeref
26        use dimens_m, only: iim, jjm, llm, nqmx
27        use dimphy, only: zmasq
28      use dimsoil, only: nsoilmx      use dimsoil, only: nsoilmx
29      use temps, only: itau_dyn, itau_phy, annee_ref, day_ref, dt      use dynredem0_m, only: dynredem0
30      use clesphys2, only: ok_orodr, nbapp_rad      use dynredem1_m, only: dynredem1
31        use exner_hyb_m, only: exner_hyb
32        USE flincom, only: flinclo, flinopen_nozoom, flininfo
33        use geopot_m, only: geopot
34      use grid_atob, only: grille_m      use grid_atob, only: grille_m
35      use grid_change, only: init_dyn_phy, dyn_phy      use grid_change, only: init_dyn_phy, dyn_phy
36      use q_sat_m, only: q_sat      use histclo_m, only: histclo
37      use exner_hyb_m, only: exner_hyb      use indicesol, only: is_oce, is_sic, is_ter, is_lic, epsfra
38      use advtrac_m, only: iniadvtrac      use iniadvtrac_m, only: iniadvtrac
39        use inidissip_m, only: inidissip
40        use inifilr_m, only: inifilr
41        use inigeom_m, only: inigeom
42        use netcdf, only: nf90_nowrite
43        use netcdf95, only: nf95_open, nf95_close, nf95_get_var, nf95_inq_varid
44        use nr_util, only: pi
45        use paramet_m, only: ip1jm, ip1jmp1
46        use phyredem_m, only: phyredem
47      use pressure_var, only: pls, p3d      use pressure_var, only: pls, p3d
48      use dynredem0_m, only: dynredem0      use q_sat_m, only: q_sat
49      use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz      use regr_lat_time_coefoz_m, only: regr_lat_time_coefoz
50      use regr_pr_o3_m, only: regr_pr_o3      use regr_pr_o3_m, only: regr_pr_o3
51        use serre, only: alphax
52        USE start_init_orog_m, only: start_init_orog, mask, phis
53        use start_init_phys_m, only: start_init_phys
54        use startdyn, only: start_init_dyn
55        use start_inter_3d_m, only: start_inter_3d
56        use temps, only: itau_phy, annee_ref, day_ref
57    
58      ! Variables local to the procedure:      ! Variables local to the procedure:
59    
# Line 53  contains Line 61  contains
61      ! (latitude and longitude of a point of the scalar grid identified      ! (latitude and longitude of a point of the scalar grid identified
62      ! by a simple index, in °)      ! by a simple index, in °)
63    
64      REAL, dimension(iim + 1, jjm + 1, llm):: uvent, t3d, tpot      REAL, dimension(iim + 1, jjm + 1, llm):: ucov, t3d, tpot
65      REAL vvent(iim + 1, jjm, llm)      REAL vcov(iim + 1, jjm, llm)
66    
67      REAL q3d(iim + 1, jjm + 1, llm, nqmx)      REAL q3d(iim + 1, jjm + 1, llm, nqmx)
68      ! (mass fractions of trace species      ! (mass fractions of trace species
# Line 76  contains Line 84  contains
84      REAL rugmer(klon)      REAL rugmer(klon)
85      real, dimension(iim + 1, jjm + 1):: relief, zstd_2d, zsig_2d, zgam_2d      real, dimension(iim + 1, jjm + 1):: relief, 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, psol      real, dimension(iim + 1, jjm + 1):: tsol_2d, qsol_2d, psol
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)
     REAL rugsrel(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      ! déclarations pour lecture glace de mer      ! déclarations pour lecture glace de mer
96      INTEGER iml_lic, jml_lic, llm_tmp, ttm_tmp      INTEGER iml_lic, jml_lic, llm_tmp, ttm_tmp
97      INTEGER itaul(1), fid      INTEGER itaul(1), fid, ncid, varid
98      REAL lev(1), date      REAL lev(1), date
99      REAL, ALLOCATABLE:: lon_lic(:, :), lat_lic(:, :)      REAL, ALLOCATABLE:: lon_lic(:, :), lat_lic(:, :)
100      REAL, ALLOCATABLE:: dlon_lic(:), dlat_lic(:)      REAL, ALLOCATABLE:: dlon_lic(:), dlat_lic(:)
# Line 104  contains Line 111  contains
111      REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm)      REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm)
112      REAL w(ip1jmp1, llm)      REAL w(ip1jmp1, llm)
113      REAL phystep      REAL phystep
114      INTEGER radpas      real trash
115    
116      !---------------------------------      !---------------------------------
117    
118      print *, "Call sequence information: etat0"      print *, "Call sequence information: etat0"
119    
     ! Construct a grid:  
   
120      dtvr = daysec / real(day_step)      dtvr = daysec / real(day_step)
121      print *, 'dtvr = ', dtvr      print *, 'dtvr = ', dtvr
122    
123        ! Construct a grid:
124    
125      pa = 5e4      pa = 5e4
126      CALL iniconst      CALL iniconst
127      CALL inigeom      CALL inigeom
# Line 131  contains Line 138  contains
138      lonfi(klon) = 0.      lonfi(klon) = 0.
139    
140      call start_init_orog(relief, zstd_2d, zsig_2d, zgam_2d, zthe_2d, zpic_2d, &      call start_init_orog(relief, zstd_2d, zsig_2d, zgam_2d, zthe_2d, zpic_2d, &
141           zval_2d) ! also compute "masque" and "phis"           zval_2d) ! also compute "mask" and "phis"
142      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
143      zmasq = pack(masque, dyn_phy)      zmasq = pack(mask, dyn_phy)
144      PRINT *, 'Masque construit'      PRINT *, 'Masque construit'
145    
146      CALL start_init_dyn(tsol_2d, psol) ! also compute "qsol_2d"      call start_init_phys(tsol_2d, qsol_2d)
147        CALL start_init_dyn(tsol_2d, psol)
148    
149      ! Compute pressure on intermediate levels:      ! Compute pressure on intermediate levels:
150      forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * psol(:, :)      forall(l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * psol
151      CALL exner_hyb(psol, p3d, pks, pk)      CALL exner_hyb(psol, p3d, pks, pk)
152      IF (MINVAL(pk) == MAXVAL(pk)) stop '"pk" should not be constant'      IF (MINVAL(pk) == MAXVAL(pk)) then
153           print *, '"pk" should not be constant'
154      pls(:, :, :) = preff * (pk(:, :, :) / cpp)**(1. / kappa)         stop 1
155      PRINT *, "minval(pls(:, :, :)) = ", minval(pls(:, :, :))      end IF
156      print *, "maxval(pls(:, :, :)) = ", maxval(pls(:, :, :))  
157        pls = preff * (pk / cpp)**(1. / kappa)
158      uvent(:, :, :) = start_inter_3d('U', rlonv, rlatv, pls)      PRINT *, "minval(pls) = ", minval(pls)
159      forall (l = 1: llm) uvent(:iim, :, l) = uvent(:iim, :, l) * cu_2d(:iim, :)      print *, "maxval(pls) = ", maxval(pls)
160      uvent(iim+1, :, :) = uvent(1, :, :)  
161        call start_inter_3d('U', rlonv, rlatv, pls, ucov)
162      vvent(:, :, :) = start_inter_3d('V', rlonu, rlatu(:jjm), pls(:, :jjm, :))      forall (l = 1: llm) ucov(:iim, :, l) = ucov(:iim, :, l) * cu_2d(:iim, :)
163      forall (l = 1: llm) vvent(:iim, :, l) = vvent(:iim, :, l) * cv_2d(:iim, :)      ucov(iim+1, :, :) = ucov(1, :, :)
164      vvent(iim + 1, :, :) = vvent(1, :, :)  
165        call start_inter_3d('V', rlonu, rlatu(:jjm), pls(:, :jjm, :), vcov)
166      t3d(:, :, :) = start_inter_3d('TEMP', rlonu, rlatv, pls)      forall (l = 1: llm) vcov(:iim, :, l) = vcov(:iim, :, l) * cv_2d(:iim, :)
167      PRINT *,  'minval(t3d(:, :, :)) = ', minval(t3d(:, :, :))      vcov(iim + 1, :, :) = vcov(1, :, :)
168      print *, "maxval(t3d(:, :, :)) = ", maxval(t3d(:, :, :))  
169        call start_inter_3d('TEMP', rlonu, rlatv, pls, t3d)
170        PRINT *,  'minval(t3d) = ', minval(t3d)
171        print *, "maxval(t3d) = ", maxval(t3d)
172    
173      tpot(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)      tpot(:iim, :, :) = t3d(:iim, :, :) * cpp / pk(:iim, :, :)
174      tpot(iim + 1, :, :) = tpot(1, :, :)      tpot(iim + 1, :, :) = tpot(1, :, :)
# Line 168  contains Line 179  contains
179      ENDDO      ENDDO
180    
181      ! Calcul de l'humidité à saturation :      ! Calcul de l'humidité à saturation :
182      qsat(:, :, :) = q_sat(t3d, pls)      qsat = q_sat(t3d, pls)
183      PRINT *, "minval(qsat(:, :, :)) = ", minval(qsat(:, :, :))      PRINT *, "minval(qsat) = ", minval(qsat)
184      print *, "maxval(qsat(:, :, :)) = ", maxval(qsat(:, :, :))      print *, "maxval(qsat) = ", maxval(qsat)
185      IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'      IF (MINVAL(qsat) == MAXVAL(qsat)) stop '"qsat" should not be constant'
186    
187      ! Water vapor:      ! Water vapor:
188      q3d(:, :, :, 1) = 0.01 * start_inter_3d('R', rlonu, rlatv, pls) * qsat      call start_inter_3d('R', rlonu, rlatv, pls, q3d(:, :, :, 1))
189        q3d(:, :, :, 1) = 0.01 * q3d(:, :, :, 1) * qsat
190      WHERE (q3d(:, :, :, 1) < 0.) q3d(:, :, :, 1) = 1E-10      WHERE (q3d(:, :, :, 1) < 0.) q3d(:, :, :, 1) = 1E-10
191      DO l = 1, llm      DO l = 1, llm
192         q3d(:, 1, l, 1) = SUM(aire_2d(:, 1) * q3d(:, 1, l, 1)) / apoln         q3d(:, 1, l, 1) = SUM(aire_2d(:, 1) * q3d(:, 1, l, 1)) / apoln
# Line 192  contains Line 204  contains
204         q3d(:, :, :, 5) = q3d(:, :, :, 5)  * 48. / 29.         q3d(:, :, :, 5) = q3d(:, :, :, 5)  * 48. / 29.
205      end if      end if
206    
207      tsol(:) = pack(tsol_2d, dyn_phy)      tsol = pack(tsol_2d, dyn_phy)
208      qsol(:) = pack(qsol_2d, dyn_phy)      qsol = pack(qsol_2d, dyn_phy)
209      sn(:) = 0. ! snow      sn = 0. ! snow
210      radsol(:) = 0.      radsol = 0.
211      tslab(:) = 0. ! IM "slab" ocean      tslab = 0. ! IM "slab" ocean
212      seaice(:) = 0.      seaice = 0.
213      rugmer(:) = 0.001      rugmer = 0.001
214      zmea(:) = pack(relief, dyn_phy)      zmea = pack(relief, dyn_phy)
215      zstd(:) = pack(zstd_2d, dyn_phy)      zstd = pack(zstd_2d, dyn_phy)
216      zsig(:) = pack(zsig_2d, dyn_phy)      zsig = pack(zsig_2d, dyn_phy)
217      zgam(:) = pack(zgam_2d, dyn_phy)      zgam = pack(zgam_2d, dyn_phy)
218      zthe(:) = pack(zthe_2d, dyn_phy)      zthe = pack(zthe_2d, dyn_phy)
219      zpic(:) = pack(zpic_2d, dyn_phy)      zpic = pack(zpic_2d, dyn_phy)
220      zval(:) = pack(zval_2d, dyn_phy)      zval = pack(zval_2d, dyn_phy)
   
     rugsrel(:) = 0.  
     IF (ok_orodr) rugsrel(:) = MAX(1.e-05, zstd(:) * zsig(:) / 2)  
221    
222      ! On initialise les sous-surfaces:      ! On initialise les sous-surfaces.
223      ! Lecture du fichier glace de terre pour fixer la fraction de terre      ! Lecture du fichier glace de terre pour fixer la fraction de terre
224      ! et de glace de terre:      ! et de glace de terre :
225      CALL flininfo("landiceref.nc", iml_lic, jml_lic, llm_tmp, &      CALL flininfo("landiceref.nc", iml_lic, jml_lic, llm_tmp, &
226           ttm_tmp, fid)           ttm_tmp, fid)
227      ALLOCATE(lat_lic(iml_lic, jml_lic))      ALLOCATE(lat_lic(iml_lic, jml_lic))
# Line 220  contains Line 229  contains
229      ALLOCATE(dlon_lic(iml_lic))      ALLOCATE(dlon_lic(iml_lic))
230      ALLOCATE(dlat_lic(jml_lic))      ALLOCATE(dlat_lic(jml_lic))
231      ALLOCATE(fraclic(iml_lic, jml_lic))      ALLOCATE(fraclic(iml_lic, jml_lic))
232      CALL flinopen_nozoom("landiceref.nc", iml_lic, jml_lic, &      CALL flinopen_nozoom(iml_lic, jml_lic, &
233           llm_tmp, lon_lic, lat_lic, lev, ttm_tmp, itaul, date, dt,  &           llm_tmp, lon_lic, lat_lic, lev, ttm_tmp, itaul, date, trash,  &
234           fid)           fid)
     CALL flinget(fid, 'landice', iml_lic, jml_lic, llm_tmp, ttm_tmp &  
          , 1, 1, fraclic)  
235      CALL flinclo(fid)      CALL flinclo(fid)
236        call nf95_open("landiceref.nc", nf90_nowrite, ncid)
237        call nf95_inq_varid(ncid, 'landice', varid)
238        call nf95_get_var(ncid, varid, fraclic)
239        call nf95_close(ncid)
240    
241      ! Interpolation sur la grille T du modèle :      ! Interpolation sur la grille T du modèle :
242      PRINT *, 'Dimensions de "landice"'      PRINT *, 'Dimensions de "landice"'
# Line 233  contains Line 244  contains
244      print *, "jml_lic = ", jml_lic      print *, "jml_lic = ", jml_lic
245    
246      ! Si les coordonnées sont en degrés, on les transforme :      ! Si les coordonnées sont en degrés, on les transforme :
247      IF (MAXVAL( lon_lic(:, :) ) > pi)  THEN      IF (MAXVAL( lon_lic ) > pi)  THEN
248         lon_lic(:, :) = lon_lic(:, :) * pi / 180.         lon_lic = lon_lic * pi / 180.
249      ENDIF      ENDIF
250      IF (maxval( lat_lic(:, :) ) > pi) THEN      IF (maxval( lat_lic ) > pi) THEN
251         lat_lic(:, :) = lat_lic(:, :) * pi/ 180.         lat_lic = lat_lic * pi/ 180.
252      ENDIF      ENDIF
253    
254      dlon_lic(:) = lon_lic(:, 1)      dlon_lic = lon_lic(:, 1)
255      dlat_lic(:) = lat_lic(1, :)      dlat_lic = lat_lic(1, :)
256    
257      flic_tmp(:iim, :) = grille_m(dlon_lic, dlat_lic, fraclic, rlonv(:iim), &      flic_tmp(:iim, :) = grille_m(dlon_lic, dlat_lic, fraclic, rlonv(:iim), &
258           rlatu)           rlatu)
259      flic_tmp(iim + 1, :) = flic_tmp(1, :)      flic_tmp(iim + 1, :) = flic_tmp(1, :)
260    
261      ! Passage sur la grille physique      ! Passage sur la grille physique
262      pctsrf(:, :)=0.      pctsrf = 0.
263      pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)      pctsrf(:, is_lic) = pack(flic_tmp, dyn_phy)
264      ! Adéquation avec le maque terre/mer      ! Adéquation avec le maque terre/mer
265      WHERE (pctsrf(:, is_lic) < EPSFRA ) pctsrf(:, is_lic) = 0.      WHERE (pctsrf(:, is_lic) < EPSFRA ) pctsrf(:, is_lic) = 0.
266      WHERE (zmasq(:) < EPSFRA) pctsrf(:, is_lic) = 0.      WHERE (zmasq < EPSFRA) pctsrf(:, is_lic) = 0.
267      pctsrf(:, is_ter) = zmasq(:)      pctsrf(:, is_ter) = zmasq
268      where (zmasq(:) > EPSFRA)      where (zmasq > EPSFRA)
269         where (pctsrf(:, is_lic) >= zmasq(:))         where (pctsrf(:, is_lic) >= zmasq)
270            pctsrf(:, is_lic) = zmasq(:)            pctsrf(:, is_lic) = zmasq
271            pctsrf(:, is_ter) = 0.            pctsrf(:, is_ter) = 0.
272         elsewhere         elsewhere
273            pctsrf(:, is_ter) = zmasq(:) - pctsrf(:, is_lic)            pctsrf(:, is_ter) = zmasq - pctsrf(:, is_lic)
274            where (pctsrf(:, is_ter) < EPSFRA)            where (pctsrf(:, is_ter) < EPSFRA)
275               pctsrf(:, is_ter) = 0.               pctsrf(:, is_ter) = 0.
276               pctsrf(:, is_lic) = zmasq(:)               pctsrf(:, is_lic) = zmasq
277            end where            end where
278         end where         end where
279      end where      end where
280    
281      ! Sous-surface océan et glace de mer (pour démarrer on met glace      ! Sous-surface océan et glace de mer (pour démarrer on met glace
282      ! de mer à 0) :      ! de mer à 0) :
283      pctsrf(:, is_oce) = 1. - zmasq(:)      pctsrf(:, is_oce) = 1. - zmasq
284      WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.      WHERE (pctsrf(:, is_oce) < EPSFRA) pctsrf(:, is_oce) = 0.
285    
286      ! Vérification que somme des sous-surfaces vaut 1:      ! Vérification que somme des sous-surfaces vaut 1:
287      ji = count(abs(sum(pctsrf(:, :), dim = 2) - 1. ) > EPSFRA)      ji = count(abs(sum(pctsrf, dim = 2) - 1.) > EPSFRA)
288      IF (ji /= 0) PRINT *, 'Problème répartition sous maille pour', ji, 'points'      IF (ji /= 0) then
289           PRINT *, 'Problème répartition sous maille pour ', ji, 'points'
290        end IF
291    
292      ! Calcul intermédiaire:      ! Calcul intermédiaire:
293      CALL massdair(p3d, masse)      CALL massdair(p3d, masse)
# Line 289  contains Line 302  contains
302    
303      ! Initialisation pour traceurs:      ! Initialisation pour traceurs:
304      call iniadvtrac      call iniadvtrac
305      ! Ecriture:      CALL inidissip
     CALL inidissip(lstardis, nitergdiv, nitergrot, niterh, tetagdiv, &  
          tetagrot, tetatemp)  
     itau_dyn = 0  
306      itau_phy = 0      itau_phy = 0
307      day_ref = dayref      day_ref = dayref
308      annee_ref = anneeref      annee_ref = anneeref
309    
310      CALL geopot(ip1jmp1, tpot, pk , pks,  phis  , phi)      CALL geopot(ip1jmp1, tpot, pk , pks,  phis, phi)
311      CALL caldyn0(0, uvent, vvent, tpot, psol, masse, pk, phis, phi, w, &      CALL caldyn0(ucov, vcov, tpot, psol, masse, pk, phis, phi, w, pbaru, &
312           pbaru, pbarv, 0)           pbarv)
313      CALL dynredem0("start.nc", dayref, phis)      CALL dynredem0("start.nc", dayref, phis)
314      CALL dynredem1("start.nc", 0., vvent, uvent, tpot, q3d, masse, psol)      CALL dynredem1("start.nc", vcov, ucov, tpot, q3d, masse, psol, itau=0)
315    
316      ! Ecriture état initial physique:      ! Ecriture état initial physique:
     print *, 'dtvr = ', dtvr  
317      print *, "iphysiq = ", iphysiq      print *, "iphysiq = ", iphysiq
     print *, "nbapp_rad = ", nbapp_rad  
318      phystep   = dtvr * REAL(iphysiq)      phystep   = dtvr * REAL(iphysiq)
     radpas    = NINT (86400./phystep/ nbapp_rad)  
319      print *, 'phystep = ', phystep      print *, 'phystep = ', phystep
     print *, "radpas = ", radpas  
320    
321      ! Initialisations :      ! Initialisations :
322      tsolsrf(:, is_ter) = tsol      tsolsrf(:, is_ter) = tsol
# Line 326  contains Line 332  contains
332      albe(:, is_oce) = 0.5      albe(:, is_oce) = 0.5
333      albe(:, is_sic) = 0.6      albe(:, is_sic) = 0.6
334      alblw = albe      alblw = albe
335      evap(:, :) = 0.      evap = 0.
336      qsolsrf(:, is_ter) = 150.      qsolsrf(:, is_ter) = 150.
337      qsolsrf(:, is_lic) = 150.      qsolsrf(:, is_lic) = 150.
338      qsolsrf(:, is_oce) = 150.      qsolsrf(:, is_oce) = 150.
# Line 340  contains Line 346  contains
346      q_ancien = 0.      q_ancien = 0.
347      agesno = 0.      agesno = 0.
348      !IM "slab" ocean      !IM "slab" ocean
349      tslab(:) = tsolsrf(:, is_oce)      tslab = tsolsrf(:, is_oce)
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(1.e-05, zstd * zsig / 2)
354      frugs(:, is_lic) = MAX(1.e-05, zstd(:) * zsig(:) / 2)      frugs(:, is_lic) = MAX(1.e-05, zstd * zsig / 2)
355      frugs(:, is_sic) = 0.001      frugs(:, is_sic) = 0.001
356      fder = 0.      fder = 0.
357      clwcon = 0.      clwcon = 0.
# Line 353  contains Line 359  contains
359      ratqs = 0.      ratqs = 0.
360      run_off_lic_0 = 0.      run_off_lic_0 = 0.
361    
362      call phyredem("startphy.nc", radpas, latfi, lonfi, pctsrf, &      call phyredem("startphy.nc", latfi, lonfi, pctsrf, &
363           tsolsrf, tsoil, tslab, seaice, qsolsrf, qsol, snsrf, albe, alblw, &           tsolsrf, tsoil, tslab, seaice, qsolsrf, qsol, snsrf, albe, alblw, &
364           evap, rain_fall, snow_fall, solsw, sollw, fder, radsol, frugs, &           evap, rain_fall, snow_fall, solsw, sollw, fder, radsol, frugs, &
365           agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, rugsrel, &           agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, &
366           t_ancien, q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0)           t_ancien, q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0)
367      CALL histclo      CALL histclo
368    

Legend:
Removed from v.12  
changed lines
  Added in v.66

  ViewVC Help
Powered by ViewVC 1.1.21