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

Last change on this file since 374 was 253, checked in by emillour, 13 years ago

Generic GCM

  • Massive update to version 0.7

EM+RW

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