source: trunk/LMDZ.GENERIC/libf/phystd/evap_generic.F90 @ 2725

Last change on this file since 2725 was 2725, checked in by aslmd, 2 years ago

switch to epsi_generic & RLVTT_generic

File size: 2.4 KB
Line 
1subroutine evap_generic(ngrid,nlayer,nq,dtime,pt,pq,pdq,pdt,igcm_generic_vap, &
2                        igcm_generic_ice,dqevap,dtevap,qevap,tevap)
3    Use generic_cloud_common_h
4    Use tracer_h
5    implicit none
6!==================================================================
7!     
8!     Purpose
9!     -------
10!     Evaporation of all generic tracers in the atmopshere.
11!     
12!     Authors
13!     -------
14!     Adapted from evap.F by Lucas Teinturier (2022)
15!     evap.F is adapted from the LMDTERRE code by B. Charnay (2010)
16!     Original author Z. X. Li (1993)
17!     
18!==================================================================
19    !Inputs
20    integer, intent(in) :: ngrid,nlayer,nq
21    real, intent(in) :: pt(ngrid,nlayer) !Temperature (K)
22    real, intent(in) :: pq(ngrid,nlayer,nq) ! Tracers (kg/kg)
23    real, intent(in) :: pdt(ngrid,nlayer) !Tendency on temperature (K/s)
24    real, intent(in) :: pdq(ngrid,nlayer,nq) ! Tendency on tracer (kg/kg/s)
25    real, intent(in) :: dtime ! physical timestep
26    integer,intent(in) :: igcm_generic_vap ! index of the vapor generic tracer
27    integer,intent(in) :: igcm_generic_ice !index of the ice generic tracer
28
29    ! Ouputs
30    real,intent(out) :: qevap(ngrid,nlayer,nq)
31    real,intent(out) :: tevap(ngrid,nlayer)
32    real,intent(out) :: dqevap(ngrid,nlayer)
33    real,intent(out) :: dtevap(ngrid,nlayer)
34
35    ! Local variables
36    real zlvdcp
37    integer l,ig ! Used for loops on the layers and grid and the tracers respectively
38
39    !   Evaporate all the ice from the tracer
40    zlvdcp = RLVTT_generic/cpp ! RLVTT_generic is the latent heat of vaporization (comes from generic_cloud_common_h attention)
41 
42    DO l=1,nlayer
43        DO ig=1,ngrid
44            qevap(ig,l,igcm_generic_vap) = pq(ig,l,igcm_generic_vap) &
45                    +pdq(ig,l,igcm_generic_vap)*dtime
46            qevap(ig,l,igcm_generic_ice) = pq(ig,l,igcm_generic_ice) &
47                    +pdq(ig,l,igcm_generic_ice)*dtime
48            tevap(ig,l) = pt(ig,l)+pdt(ig,l)*dtime
49            dqevap(ig,l) = max(0.,qevap(ig,l,igcm_generic_ice))/dtime
50            dtevap(ig,l) = -dqevap(ig,l)*zlvdcp
51            qevap(ig,l,igcm_generic_vap) = qevap(ig,l,igcm_generic_vap) &
52                    +dqevap(ig,l)*dtime
53            qevap(ig,l,igcm_generic_ice) = 0.
54            tevap(ig,l) = tevap(ig,l) + dtevap(ig,l)*dtime
55        ENDDO ! ngrid
56    ENDDO ! nlayer
57
58end subroutine evap_generic
59
Note: See TracBrowser for help on using the repository browser.