source: trunk/LMDZ.GENERIC/libf/phystd/evap.F @ 2176

Last change on this file since 2176 was 1384, checked in by emillour, 10 years ago

Generic GCM:

  • Some code cleanup: turning comcstfi.h into module comcstfi_mod.F90

EM

File size: 2.2 KB
Line 
1      subroutine evap(ngrid,nlayer,nq,dtime,pt, pq, pdq, pdt, 
2     $     dqevap,dtevap, qevap, tevap)
3       
4      use watercommon_h
5      USE tracer_h
6
7      implicit none
8
9!==================================================================
10!     
11!     Purpose
12!     -------
13!     Evaporation of all water in the atmopshere.
14!     
15!     Authors
16!     -------
17!     Adapted from the LMDTERRE code by B. Charnay (2010)
18!     Original author Z. X. Li (1993)
19!     
20!==================================================================
21
22      INTEGER ngrid,nlayer,nq
23
24! Arguments:
25      REAL pt(ngrid,nlayer)
26      REAL pq(ngrid,nlayer,nq)
27      REAL pdt(ngrid,nlayer)
28      REAL pdq(ngrid,nlayer,nq)
29      REAL dqevap(ngrid,nlayer)
30      REAL dtevap(ngrid,nlayer)
31      REAL qevap(ngrid,nlayer,nq)
32      REAL dtime
33 
34! Local:
35      REAL tevap(ngrid,nlayer)
36      REAL zlvdcp
37      REAL zlsdcp
38      REAL zdelta
39      INTEGER l,ig
40
41!
42! Re-evaporer l'eau liquide nuageuse
43!
44
45      DO l=1,nlayer
46        DO ig=1,ngrid
47         qevap(ig,l,igcm_h2o_vap)=pq(ig,l,igcm_h2o_vap)
48     s                            +pdq(ig,l,igcm_h2o_vap)*dtime
49         qevap(ig,l,igcm_h2o_ice)=pq(ig,l,igcm_h2o_ice)
50     s                            +pdq(ig,l,igcm_h2o_ice)*dtime
51         tevap(ig,l)=pt(ig,l)+pdt(ig,l)*dtime
52        ENDDO
53      ENDDO
54
55      DO l = 1, nlayer
56        DO ig = 1, ngrid
57         zlvdcp=RLVTT/RCPD!/(1.0+RVTMP2*qevap(ig,l,igcm_h2o_vap))
58         zlsdcp=RLSTT/RCPD!/(1.0+RVTMP2*qevap(ig,l,igcm_h2o_vap))
59         ! ignoring qevap term creates huge difference when qevap large!!!
60
61         zdelta = MAX(0.,SIGN(1.,T_h2O_ice_liq-tevap(ig,l))) ! what is this?
62                                                  ! for division between water / liquid
63         dqevap(ig,l) = MAX(0.0,qevap(ig,l,igcm_h2o_ice))/dtime
64         dtevap(ig,l) = - dqevap(ig,l)*RLVTT/RCPD ! exactly as in largescale.F
65!         dtevap(ig,l) = - dqevap(ig,l)
66!     s                       * (zlvdcp*(1.-zdelta)+zlsdcp*zdelta)
67         qevap(ig,l,igcm_h2o_vap) = qevap(ig,l,igcm_h2o_vap)   
68     s                              +dqevap(ig,l)*dtime
69         qevap(ig,l,igcm_h2o_ice) = 0.0
70         tevap(ig,l) = tevap(ig,l)+dtevap(ig,l)*dtime
71
72        ENDDO
73      ENDDO
74
75      END
Note: See TracBrowser for help on using the repository browser.