/[lmdze]/trunk/phylmd/Interface_surf/calcul_fluxs.f90
ViewVC logotype

Diff of /trunk/phylmd/Interface_surf/calcul_fluxs.f90

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

revision 278 by guez, Thu Apr 19 17:54:55 2018 UTC revision 279 by guez, Fri Jul 20 14:30:23 2018 UTC
# Line 4  module calcul_fluxs_m Line 4  module calcul_fluxs_m
4    
5  contains  contains
6    
7    SUBROUTINE calcul_fluxs(dtime, tsurf, p1lay, cal, beta, coef1lay, ps, &    SUBROUTINE calcul_fluxs(dtime, tsurf, p1lay, cal, beta, coef1lay, ps, qsurf, &
8         qsurf, radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, petAcoef, &         radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, petAcoef, peqAcoef, &
9         peqAcoef, petBcoef, peqBcoef, tsurf_new, evap, fluxlat, flux_t, &         petBcoef, peqBcoef, tsurf_new, evap, fluxlat, flux_t, dflux_s, dflux_l)
        dflux_s, dflux_l)  
10    
11      ! Cette routine calcule les flux en h et q à l'interface et une      ! Cette routine calcule les flux en h et q à l'interface et une
12      ! température de surface.      ! température de surface.
13    
14      ! L. Fairhead, April 2000      ! L. Fairhead, April 2000
15    
16      USE fcttre, ONLY: foede, foeew      ! Note that, if cal = 0, beta = 1 and dif_grnd = 0, then tsurf_new
17        ! = tsurf and qsurf = qsat.
18    
19      use nr_util, only: assert_eq      use nr_util, only: assert_eq
20    
21        USE fcttre, ONLY: foede, foeew
22      USE suphec_m, ONLY: rcpd, rd, retv, rlstt, rlvtt, rtt      USE suphec_m, ONLY: rcpd, rd, retv, rlstt, rlvtt, rtt
23      USE yoethf_m, ONLY: r2es, r5ies, r5les, rvtmp2      USE yoethf_m, ONLY: r2es, r5ies, r5les, rvtmp2
24    
# Line 57  contains Line 60  contains
60    
61      ! Local:      ! Local:
62      integer i      integer i
63      integer knon ! nombre de points a traiter      integer knon ! nombre de points \`a traiter
64      real, dimension(size(ps)):: mh, oh, mq, nq, oq, dq_s_dt, coef ! (knon)      real, dimension(size(ps)):: mh, oh, mq, nq, oq, dq_s_dt, coef ! (knon)
65      real qsat(size(ps)) ! (knon) mass fraction      real qsat(size(ps)) ! (knon) mass fraction
66      real sl(size(ps)) ! (knon) chaleur latente d'évaporation ou de sublimation      real sl(size(ps)) ! (knon) chaleur latente d'évaporation ou de sublimation
# Line 67  contains Line 70  contains
70    
71      !---------------------------------------------------------------------      !---------------------------------------------------------------------
72    
73      knon = assert_eq((/size(tsurf), size(p1lay), size(cal), size(beta), &      knon = assert_eq([size(tsurf), size(p1lay), size(cal), size(beta), &
74           size(coef1lay), size(ps), size(qsurf), size(radsol), size(dif_grnd), &           size(coef1lay), size(ps), size(qsurf), size(radsol), size(dif_grnd), &
75           size(t1lay), size(q1lay), size(u1lay), size(v1lay), size(petAcoef), &           size(t1lay), size(q1lay), size(u1lay), size(v1lay), size(petAcoef), &
76           size(peqAcoef), size(petBcoef), size(peqBcoef), size(tsurf_new), &           size(peqAcoef), size(petBcoef), size(peqBcoef), size(tsurf_new), &
77           size(evap), size(fluxlat), size(flux_t), size(dflux_s), &           size(evap), size(fluxlat), size(flux_t), size(dflux_s), &
78           size(dflux_l)/), "calcul_fluxs knon")           size(dflux_l)], "calcul_fluxs knon")
79    
80      ! Traitement de l'humidité du sol      ! Traitement de l'humidité du sol
81    
# Line 92  contains Line 95  contains
95      ! Q      ! Q
96      oq = 1. - beta * coef * peqBcoef * dtime      oq = 1. - beta * coef * peqBcoef * dtime
97      mq = beta * coef * (peqAcoef - qsat + dq_s_dt * tsurf) / oq      mq = beta * coef * (peqAcoef - qsat + dq_s_dt * tsurf) / oq
98      nq = beta * coef * (- 1. * dq_s_dt) / oq      nq = - beta * coef * dq_s_dt / oq
99    
100      ! H      ! H
101      oh = 1. - (coef * petBcoef * dtime)      oh = 1. - coef * petBcoef * dtime
102      mh = coef * petAcoef / oh      mh = coef * petAcoef / oh
103      dflux_s = - (coef * RCPD)/ oh      dflux_s = - coef * RCPD / oh
104    
     ! Tsurface  
105      tsurf_new = (tsurf + cal / RCPD * dtime * (radsol + mh + sl * mq) &      tsurf_new = (tsurf + cal / RCPD * dtime * (radsol + mh + sl * mq) &
106           + dif_grnd * t_grnd * dtime) / (1. - dtime * cal / RCPD * (dflux_s &           + dif_grnd * t_grnd * dtime) / (1. - dtime * cal / RCPD * (dflux_s &
107           + sl * nq) + dtime * dif_grnd)           + sl * nq) + dtime * dif_grnd)
   
108      evap = - mq - nq * tsurf_new      evap = - mq - nq * tsurf_new
109      fluxlat = - evap * sl      fluxlat = - evap * sl
110      flux_t = mh + dflux_s * tsurf_new      flux_t = mh + dflux_s * tsurf_new
111      dflux_l = sl * nq      dflux_l = sl * nq
   
     ! Nouvelle valeur de l'humidité au dessus du sol :  
112      qsurf = (peqAcoef - peqBcoef * evap * dtime) * (1. - beta) + beta * (qsat &      qsurf = (peqAcoef - peqBcoef * evap * dtime) * (1. - beta) + beta * (qsat &
113           + dq_s_dt * (tsurf_new - tsurf))           + dq_s_dt * (tsurf_new - tsurf))
114    

Legend:
Removed from v.278  
changed lines
  Added in v.279

  ViewVC Help
Powered by ViewVC 1.1.21