source: LMDZ5/trunk/libf/phymar/Mod_PHY_AT_kkl.f90 @ 2089

Last change on this file since 2089 was 2089, checked in by Laurent Fairhead, 10 years ago

Inclusion de la physique de MAR


Integration of MAR physics

File size: 8.3 KB
RevLine 
[2089]1      module Mod_PHY_AT_kkl
2
3!--------------------------------------------------------------------------+
4!                                                     Sun 12-May-2013  MAR |
5!     module Mod_PHY_AT_kkl contains the main (prognostic) variables of    |
6!                Turbulent Vertical Diffusion Scheme                       |
7!                                                                          |
8!     version 3.p.4.1 created by H. Gallee,           Tue 12-Mar-2013      |
9!           Last Modification by H. Gallee,           Sun 12-May-2013      |
10!                                                                          |
11!--------------------------------------------------------------------------+
12
13
14      use Mod_Real
15
16
17      IMPLICIT NONE
18
19
20
21! E-e and K-l models parameters
22! -----------------------------
23
24      real(kind=real8), SAVE   ::   TKEmin       = 0.0001   ! Minimum SBL    turbulent kinetic   energy
25
26      real(kind=real8), SAVE   ::      cmub      = 0.0900   ! Ee Model Parameter (Bintanja , 2000, BLM (95),       mid p.355)
27      real(kind=real8), SAVE   ::   sqrcmub      = 3.333    !
28      real(kind=real8), SAVE   ::     c1epb      = 1.46     !
29      real(kind=real8), SAVE   ::     c2epb      = 1.92     !
30      real(kind=real8), SAVE   ::     sigeb      = 0.862    !
31      real(kind=real8), SAVE   ::     sigkb      = 1.000    !
32
33      real(kind=real8), SAVE   ::      cmud      = 0.0330   ! Ee Model Parameter (Duynkerke, 1988, JAS (45), (19), top p.868)
34      real(kind=real8), SAVE   ::   sqrcmud      = 5.500    !                    (c_mu)^1/2=(0.033)^1/2=5.50           p.869
35      real(kind=real8), SAVE   ::     c1epd      = 1.46     !                                                          p.868
36      real(kind=real8), SAVE   ::     c2epd      = 1.83     !                                                          p.868
37      real(kind=real8), SAVE   ::     siged      = 0.420    !
38      real(kind=real8), SAVE   ::     sigkd      = 1.000    !
39
40      real(kind=real8), SAVE   ::      cmuk      = 0.0900   ! Ee Model Parameter (Kitada   , 1987, BLM (41),       top p.220)
41      real(kind=real8), SAVE   ::   sqrcmuk      = 3.333    !                    (c_mu)^1/2=(0.090)^1/2=3.333
42      real(kind=real8), SAVE   ::     c1epk      = 1.44     !                                                      top p.220
43      real(kind=real8), SAVE   ::     c2epk      = 1.92     !
44      real(kind=real8), SAVE   ::     sigek      = 0.769    !
45      real(kind=real8), SAVE   ::     sigkk      = 1.000    !
46
47      real(kind=real8), SAVE   ::   sqrcmut      = 4.0000   ! Kl Model Parameter (Schayes & Thunis, 1990, Contrib. 60 Inst.Astr.Geoph. p.8)
48      real(kind=real8), SAVE   ::     siget      = 0.420    !
49      real(kind=real8), SAVE   ::     sigkt      = 1.200    !
50
51      real(kind=real8), SAVE   ::    betahr      = 2.0000   ! Ee Model Parameter (Huang and Raman,  1991, BLM (55), p.386 and (A22) p.405)
52
53      real(kind=real8), SAVE   ::      cmu                  !
54      real(kind=real8), SAVE   ::   sqrcmu                  !
55      real(kind=real8), SAVE   ::     c1ep                  !
56      real(kind=real8), SAVE   ::     c2ep                  !
57      real(kind=real8), SAVE   ::     sige                  !
58      real(kind=real8), SAVE   ::     sigk                  !
59      real(kind=real8), SAVE   ::   vK_inv                  ! Inverse      of  Von-Karman Constant
60
61! #KA real(kind=real8), SAVE   ::   zz__KA       = 5.0000   ! Height below which use a vertical weighted box filter of TKE, e  [m.agl]
62! #KA integer, SAVE            ::   mz__KA                  ! Level  below which use a vertical weighted box filter of TKE, e 
63! #KA logical            ::   log_KA       = .true.   ! Swich  deciding use of a vertical weighted box filter of TKE, e 
64
65     
66
67! PHY_AT INPUT        Variables
68! -----------------------------
69
70      integer, SAVE         ,ALLOCATABLE ,dimension(:)    ::  ii__AT  ! WORK   point   i Coordinate
71      integer, SAVE         ,ALLOCATABLE ,dimension(:)    ::  jj__AT  ! WORK   point   j Coordinate
72      integer, SAVE         ,ALLOCATABLE ,dimension(:,:)  ::  ikl_AT  ! WORK   point vec Coordinate
73
74      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  var_AT  ! Dummy    to Diffuse                                        [x]
75      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  Ac0_AT  ! Tridiagonal Matrix Coefficient A: Common Factor        [m2/s3]
76      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  Cc0_AT  ! Tridiagonal Matrix Coefficient C: Common Factor        [m2/s3]
77      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  Ac__AT  ! Tridiagonal Matrix Coefficient A: Common Factor (t)     [s/m2]
78      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  Cc__AT  ! Tridiagonal Matrix Coefficient C: Common Factor (t)     [s/m2]
79      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  Kz0_AT  ! Vertical  Turbulent Diffusion Coefficient (MINIMUM)     [m2/s]
80      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  Kz__AT  ! Vertical  Turbulent Diffusion Coefficient               [m2/s]
81      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  Kzm_AT  ! Vertical  Turbulent Diffusion Coefficient (Momentum)    [m2/s]
82      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  Kzh_AT  ! Vertical  Turbulent Diffusion Coefficient (Scalars)     [m2/s]
83      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  Kzh0AT  ! Vertical  Turbulent Diffusion Coefficient (Scalars)     [m2/s]
84      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  A___AT  ! Tridiagonal Matrix Coefficient A                           [-]
85      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  B___AT  ! Tridiagonal Matrix Coefficient B                           [-]
86      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  C___AT  ! Tridiagonal Matrix Coefficient C                           [-]
87      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  D___AT  ! Independant Term               D                           [x]
88      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  P___AT  ! Auxiliary   Term               P                           [-]
89      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  Q___AT  ! Auxiliary   Term               Q                           [-]
90      real*16         ,ALLOCATABLE ,dimension(:)    ::  X___AT  ! Auxiliary   Unknown            X                           [x]
91
92      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  LMO_AT  ! Monin-Obukhov     Length         (Grid Cell Average)       [m]
93      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:)    ::  zi__AT  ! Inversion         Height         (Grid Cell Average)[m a.g.l.]
94
95
96
97! PHY_AT INPUT/OUTPUT Variables
98! -----------------------------
99
100      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  TKE_AT  ! Turbulent Kinetic Energy                               [m2/s2]
101      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  TrT_AT  ! Turbulent Kinetic Energy Transport                     [m2/s2]
102      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  eps_AT  ! Turbulent Kinetic Energy Dissipation                   [m2/s3]
103
104
105
106! PHY_AT OUTPUT       Variables
107! -----------------------------
108
109      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dua_AT  ! Wind Speed (x-direc.) Tendency                          [m/s2]
110      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dva_AT  ! Wind Speed (y-direc.) Tendency                          [m/s2]
111      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dpktAT  ! Potential Temperature Tendency, divided by p0**(R/Cp)   [x/s]
112      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dqv_AT  ! Specific  Humidity    Tendency                      [kg/kg/s]
113      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dqw_AT  ! Cloud Droplets Concen.Tendency                      [kg/kg/s]
114      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dqi_AT  ! Cloud Crystals Concen.Tendency                      [kg/kg/s]
115      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dCi_AT  ! CCNi           Concen.Tendency                          [1/s]
116      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dqs_AT  ! Snow  Particls Concen.Tendency                      [kg/kg/s]
117      real(kind=real8), SAVE,ALLOCATABLE ,dimension(:,:)  ::  dqr_AT  ! Rain  Drops    Concen.Tendency                      [kg/kg/s]
118
119
120
121      end module Mod_PHY_AT_kkl
Note: See TracBrowser for help on using the repository browser.