/[lmdze]/trunk/Sources/phylmd/Interface_surf/calcul_fluxs.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/Interface_surf/calcul_fluxs.f

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

trunk/phylmd/Interface_surf/calcul_fluxs.f revision 104 by guez, Thu Sep 4 10:05:52 2014 UTC trunk/Sources/phylmd/Interface_surf/calcul_fluxs.f revision 206 by guez, Tue Aug 30 12:52:46 2016 UTC
# Line 4  module calcul_fluxs_m Line 4  module calcul_fluxs_m
4    
5  contains  contains
6    
7    SUBROUTINE calcul_fluxs(nisurf, dtime, tsurf, p1lay, cal, beta, coef1lay, &    SUBROUTINE calcul_fluxs(dtime, tsurf, p1lay, cal, beta, coef1lay, ps, &
8         ps, qsurf, radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, petAcoef, &         qsurf, radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, petAcoef, &
9         peqAcoef, petBcoef, peqBcoef, tsurf_new, evap, fluxlat, fluxsens, &         peqAcoef, petBcoef, peqBcoef, tsurf_new, evap, fluxlat, flux_t, &
10         dflux_s, dflux_l)         dflux_s, dflux_l)
11    
12      ! Cette routine calcule les fluxs en h et q à l'interface et une      ! Cette routine calcule les fluxs en h et q à l'interface et une
# Line 14  contains Line 14  contains
14    
15      ! L. Fairhead April 2000      ! L. Fairhead April 2000
16    
     USE abort_gcm_m, ONLY: abort_gcm  
     USE indicesol, ONLY: is_ter  
17      USE fcttre, ONLY: dqsatl, dqsats, foede, foeew, qsatl, qsats, thermcep      USE fcttre, ONLY: dqsatl, dqsats, foede, foeew, qsatl, qsats, thermcep
     USE interface_surf, ONLY: run_off  
18      use nr_util, only: assert_eq      use nr_util, only: assert_eq
19      USE suphec_m, ONLY: rcpd, rd, retv, rkappa, rlstt, rlvtt, rtt      USE suphec_m, ONLY: rcpd, rd, retv, rlstt, rlvtt, rtt
20      USE yoethf_m, ONLY: r2es, r5ies, r5les, rvtmp2      USE yoethf_m, ONLY: r2es, r5ies, r5les, rvtmp2
21    
     integer, intent(IN):: nisurf ! surface a traiter  
22      real, intent(IN):: dtime      real, intent(IN):: dtime
23      real, intent(IN):: tsurf(:) ! (knon) temperature de surface      real, intent(IN):: tsurf(:) ! (knon) température de surface
24      real, intent(IN):: p1lay(:) ! (knon) pression 1er niveau (milieu de couche)  
25        real, intent(IN):: p1lay(:) ! (knon)
26        ! pression première couche (milieu de couche)
27    
28      real, intent(IN):: cal(:) ! (knon) capacité calorifique du sol      real, intent(IN):: cal(:) ! (knon) capacité calorifique du sol
29      real, intent(IN):: beta(:) ! (knon) evap reelle      real, intent(IN):: beta(:) ! (knon) évaporation réelle
30      real, intent(IN):: coef1lay(:) ! (knon) coefficient d'échange      real, intent(IN):: coef1lay(:) ! (knon) coefficient d'échange
31      real, intent(IN):: ps(:) ! (knon) pression au sol      real, intent(IN):: ps(:) ! (knon) pression au sol
32      real, intent(OUT):: qsurf(:) ! (knon) humidite de l'air au dessus du sol      real, intent(OUT):: qsurf(:) ! (knon) humidité de l'air au-dessus du sol
33      real, intent(IN):: radsol(:) ! (knon) rayonnement net au sol (LW + SW)  
34        real, intent(IN):: radsol(:) ! (knon)
35        ! rayonnement net au sol (longwave + shortwave)
36    
37      real, intent(IN):: dif_grnd(:) ! (knon)      real, intent(IN):: dif_grnd(:) ! (knon)
38      ! coefficient diffusion vers le sol profond      ! coefficient de diffusion vers le sol profond
39    
40      real, intent(IN):: t1lay(:), q1lay(:), u1lay(:), v1lay(:) ! (knon)      real, intent(IN):: t1lay(:), q1lay(:), u1lay(:), v1lay(:) ! (knon)
41    
42      real, intent(IN):: petAcoef(:), peqAcoef(:) ! (knon)      real, intent(IN):: petAcoef(:), peqAcoef(:) ! (knon)
43      ! coefficients A de la résolution de la couche limite pour t et q      ! coefficients A de la résolution de la couche limite pour T et q
44    
45      real, intent(IN):: petBcoef(:), peqBcoef(:) ! (knon)      real, intent(IN):: petBcoef(:), peqBcoef(:) ! (knon)
46      ! petBcoef coeff. B de la resolution de la CL pour t      ! coefficients B de la résolution de la couche limite pour t et q
     ! peqBcoef coeff. B de la resolution de la CL pour q  
47    
48      real, intent(OUT):: tsurf_new(:) ! (knon) température au sol      real, intent(OUT):: tsurf_new(:) ! (knon) température au sol
49      real, intent(OUT):: evap(:), fluxlat(:), fluxsens(:) ! (knon)      real, intent(OUT):: evap(:) ! (knon)
50      ! fluxlat flux de chaleur latente  
51      ! fluxsens flux de chaleur sensible      real, intent(OUT):: fluxlat(:), flux_t(:) ! (knon)
52        ! flux de chaleur latente et sensible
53    
54      real, intent(OUT):: dflux_s(:), dflux_l(:) ! (knon)      real, intent(OUT):: dflux_s(:), dflux_l(:) ! (knon)
55      ! Dérivées des flux dF/dTs (W m-2 K-1)      ! dérivées des flux de chaleurs sensible et latente par rapport à
56      ! dflux_s derivee du flux de chaleur sensible / Ts      ! Ts (W m-2 K-1)
     ! dflux_l derivee du flux de chaleur latente / Ts  
57    
58      ! Local:      ! Local:
59      integer i      integer i
     real, dimension(size(ps)) :: zx_mh, zx_nh, zx_oh  
     real, dimension(size(ps)) :: zx_mq, zx_nq, zx_oq  
     real, dimension(size(ps)) :: zx_pkh, zx_dq_s_dt, zx_qsat, zx_coef  
     real, dimension(size(ps)) :: zx_sl, zx_k1  
     real, dimension(size(ps)) :: zx_q_0 , d_ts  
     logical zdelta  
     real zcvm5, zx_qs, zcor, zx_dq_s_dh  
     real :: bilan_f, fq_fonte  
     REAL :: subli, fsno  
     REAL :: qsat_new, q1_new  
60      integer knon ! nombre de points a traiter      integer knon ! nombre de points a traiter
61        real, dimension(size(ps)):: mh, oh, mq, nq, oq, dq_s_dt, coef ! (knon)
62        real qsat(size(ps)) ! mass fraction
63        real sl(size(ps)) ! chaleur latente d'évaporation ou de sublimation
64        logical delta
65        real zcor
66      real, parameter:: t_grnd = 271.35, t_coup = 273.15      real, parameter:: t_grnd = 271.35, t_coup = 273.15
67    
68      !---------------------------------------------------------------------      !---------------------------------------------------------------------
# Line 75  contains Line 71  contains
71           size(coef1lay), size(ps), size(qsurf), size(radsol), size(dif_grnd), &           size(coef1lay), size(ps), size(qsurf), size(radsol), size(dif_grnd), &
72           size(t1lay), size(q1lay), size(u1lay), size(v1lay), size(petAcoef), &           size(t1lay), size(q1lay), size(u1lay), size(v1lay), size(petAcoef), &
73           size(peqAcoef), size(petBcoef), size(peqBcoef), size(tsurf_new), &           size(peqAcoef), size(petBcoef), size(peqBcoef), size(tsurf_new), &
74           size(evap), size(fluxlat), size(fluxsens), size(dflux_s), &           size(evap), size(fluxlat), size(flux_t), size(dflux_s), &
75           size(dflux_l)/), "calcul_fluxs knon")           size(dflux_l)/), "calcul_fluxs knon")
76    
77      if (size(run_off) /= knon .AND. nisurf == is_ter) then      ! Traitement de l'humidité du sol
78         print *, 'Bizarre, le nombre de points continentaux'  
79         print *, 'a change entre deux appels. J''arrete.'      IF (thermcep) THEN
80         call abort_gcm('calcul_fluxs', 'Pb run_off', 1)         DO i = 1, knon
81      endif            delta = rtt >= tsurf(i)
82              qsat(i) = MIN(0.5, r2es * FOEEW(tsurf(i), delta) / ps(i))
83      ! Traitement humidite du sol            zcor = 1. / (1. - retv * qsat(i))
84              qsat(i) = qsat(i) * zcor
85      evap = 0.            dq_s_dt(i) = RCPD * FOEDE(tsurf(i), delta, merge(R5IES * RLSTT, &
86      fluxsens=0.                 R5LES * RLVTT, delta) / RCPD / (1. + RVTMP2 * q1lay(i)), &
87      fluxlat=0.                 qsat(i), zcor) / RLVTT
88      dflux_s = 0.         ENDDO
89      dflux_l = 0.      ELSE
90           DO i = 1, knon
91      ! zx_qs = qsat en kg/kg            IF (tsurf(i) < t_coup) THEN
92                 qsat(i) = qsats(tsurf(i)) / ps(i)
93      DO i = 1, knon               dq_s_dt(i) = RCPD * dqsats(tsurf(i), qsat(i)) / RLVTT
        zx_pkh(i) = (ps(i)/ps(i))**RKAPPA  
        IF (thermcep) THEN  
           zdelta= rtt >= tsurf(i)  
           zcvm5 = merge(R5IES*RLSTT, R5LES*RLVTT, zdelta)  
           zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q1lay(i))  
           zx_qs= r2es * FOEEW(tsurf(i), zdelta)/ps(i)  
           zx_qs=MIN(0.5, zx_qs)  
           zcor=1./(1.-retv*zx_qs)  
           zx_qs=zx_qs*zcor  
           zx_dq_s_dh = FOEDE(tsurf(i), zdelta, zcvm5, zx_qs, zcor) &  
                /RLVTT / zx_pkh(i)  
        ELSE  
           IF (tsurf(i).LT.t_coup) THEN  
              zx_qs = qsats(tsurf(i)) / ps(i)  
              zx_dq_s_dh = dqsats(tsurf(i), zx_qs)/RLVTT &  
                   / zx_pkh(i)  
94            ELSE            ELSE
95               zx_qs = qsatl(tsurf(i)) / ps(i)               qsat(i) = qsatl(tsurf(i)) / ps(i)
96               zx_dq_s_dh = dqsatl(tsurf(i), zx_qs)/RLVTT &               dq_s_dt(i) = RCPD * dqsatl(tsurf(i), qsat(i)) / RLVTT
                   / zx_pkh(i)  
97            ENDIF            ENDIF
98         ENDIF         ENDDO
99         zx_dq_s_dt(i) = RCPD * zx_pkh(i) * zx_dq_s_dh      ENDIF
100         zx_qsat(i) = zx_qs  
101         zx_coef(i) = coef1lay(i) &      coef = coef1lay * (1. + SQRT(u1lay**2 + v1lay**2)) * p1lay / (RD * t1lay)
102              * (1.0+SQRT(u1lay(i)**2+v1lay(i)**2)) &      sl = merge(RLSTT, RLVTT, tsurf < RTT)
103              * p1lay(i)/(RD*t1lay(i))  
104        ! Q
105      ENDDO      oq = 1. - (beta * coef * peqBcoef * dtime)
106        mq = beta * coef * (peqAcoef - qsat + dq_s_dt * tsurf) / oq
107      ! === Calcul de la temperature de surface ===      nq = beta * coef * (- 1. * dq_s_dt) / oq
108    
109      ! zx_sl = chaleur latente d'evaporation ou de sublimation      ! H
110        oh = 1. - (coef * petBcoef * dtime)
111      do i = 1, knon      mh = coef * petAcoef / oh
112         zx_sl(i) = RLVTT      dflux_s = - (coef * RCPD)/ oh
113         if (tsurf(i) .LT. RTT) zx_sl(i) = RLSTT  
114         zx_k1(i) = zx_coef(i)      ! Tsurface
115      enddo      tsurf_new = (tsurf + cal / RCPD * dtime * (radsol + mh + sl * mq) &
116             + dif_grnd * t_grnd * dtime) / (1. - dtime * cal / RCPD * (dflux_s &
117      do i = 1, knon           + sl * nq) + dtime * dif_grnd)
118         ! Q  
119         zx_oq(i) = 1. - (beta(i) * zx_k1(i) * peqBcoef(i) * dtime)      evap = - mq - nq * tsurf_new
120         zx_mq(i) = beta(i) * zx_k1(i) * &      fluxlat = - evap * sl
121              (peqAcoef(i) - zx_qsat(i) &      flux_t = mh + dflux_s * tsurf_new
122              + zx_dq_s_dt(i) * tsurf(i)) &      dflux_l = sl * nq
123              / zx_oq(i)  
124         zx_nq(i) = beta(i) * zx_k1(i) * (-1. * zx_dq_s_dt(i)) &      ! Nouvelle valeur de l'humidité au dessus du sol :
125              / zx_oq(i)      qsurf = (peqAcoef - peqBcoef * evap * dtime) * (1. - beta) + beta * (qsat &
126             + dq_s_dt * (tsurf_new - tsurf))
        ! H  
        zx_oh(i) = 1. - (zx_k1(i) * petBcoef(i) * dtime)  
        zx_mh(i) = zx_k1(i) * petAcoef(i) / zx_oh(i)  
        zx_nh(i) = - (zx_k1(i) * RCPD * zx_pkh(i))/ zx_oh(i)  
   
        ! Tsurface  
        tsurf_new(i) = (tsurf(i) + cal(i)/(RCPD * zx_pkh(i)) * dtime * &  
             (radsol(i) + zx_mh(i) + zx_sl(i) * zx_mq(i)) &  
             + dif_grnd(i) * t_grnd * dtime)/ &  
             ( 1. - dtime * cal(i)/(RCPD * zx_pkh(i)) * ( &  
             zx_nh(i) + zx_sl(i) * zx_nq(i)) &  
             + dtime * dif_grnd(i))  
   
   
        ! Y'a-t-il fonte de neige?  
   
        ! fonte_neige = (nisurf /= is_oce) .AND. &  
        ! & (snow(i) > epsfra .OR. nisurf == is_sic .OR. nisurf == is_lic) &  
        ! & .AND. (tsurf_new(i) >= RTT)  
        ! if (fonte_neige) tsurf_new(i) = RTT  
        d_ts(i) = tsurf_new(i) - tsurf(i)  
        ! zx_h_ts(i) = tsurf_new(i) * RCPD * zx_pkh(i)  
        ! zx_q_0(i) = zx_qsat(i) + zx_dq_s_dt(i) * d_ts(i)  
        !== flux_q est le flux de vapeur d'eau: kg/(m**2 s) positive vers bas  
        !== flux_t est le flux de cpt (energie sensible): j/(m**2 s)  
        evap(i) = - zx_mq(i) - zx_nq(i) * tsurf_new(i)  
        fluxlat(i) = - evap(i) * zx_sl(i)  
        fluxsens(i) = zx_mh(i) + zx_nh(i) * tsurf_new(i)  
        ! Derives des flux dF/dTs (W m-2 K-1):  
        dflux_s(i) = zx_nh(i)  
        dflux_l(i) = (zx_sl(i) * zx_nq(i))  
        ! Nouvelle valeure de l'humidite au dessus du sol  
        qsat_new=zx_qsat(i) + zx_dq_s_dt(i) * d_ts(i)  
        q1_new = peqAcoef(i) - peqBcoef(i)*evap(i)*dtime  
        qsurf(i)=q1_new*(1.-beta(i)) + beta(i)*qsat_new  
     ENDDO  
127    
128    END SUBROUTINE calcul_fluxs    END SUBROUTINE calcul_fluxs
129    

Legend:
Removed from v.104  
changed lines
  Added in v.206

  ViewVC Help
Powered by ViewVC 1.1.21