source: trunk/LMDZ.MARS/libf/phymars/thermcell.F90 @ 164

Last change on this file since 164 was 161, checked in by acolaitis, 14 years ago

================================================
======== IMPLEMENTATION OF THERMALS ============
================================================

Author: A. Colaitis (2011-06-16)

The main goal of this revision is to start including the thermals into the model
for development purposes. Users should not use the thermals yet, as
several major configuration changes still need to be done.

This version includes :

  • updraft and downdraft parametrizations
  • velocity in the thermal, including drag
  • plume height analysis
  • closure equation
  • updraft transport of heat, tracers and momentum
  • downdraft transport of heat

This model should not be used without upcoming developments, namely :

  • downdraft transport of tracers and momentum
  • updraft & downdraft transport of q2 (tke)
  • revision of vdif_kc to compute q2 for non-stratified cases

Thermals could also include in a later revision :

  • momentum loss during transport (horizontal drag)

Compilation of the thermals has been successfully tested on ifort, gfortran and pgf90

================================================
================================================

M libf/phymars/callkeys.h
M libf/phymars/inifis.F

Added new control flags to call the thermals :

  • calltherm (false by default) <- to call thermals
  • outptherm (false by default) <- to output thermal-related diagnostics (for dev purposes)

================================================

M libf/phymars/vdifc.F
------> added a temporary output for thermal-related diagnostics

M libf/phymars/testphys1d.F
------> added treatment for a initialization from a profile of neutral gas (ar)

-> will be transformed in a decaying tracer for thermal diagnostics

M libf/phymars/physiq.F
------> added a section to call the thermals

-> changed the call to convadj
-> added thermal-related outputs for diagnostics

M libf/phymars/convadj.F
------> takes now into account the height of thermals to execute convective adjustment

=> note : convective adjustment needs to be activated when using thermals, in case of a

second instable layer above the thermals

================================================

A libf/phymars/calltherm_interface.F90
------> Interface between physiq.F and the thermals

A libf/phymars/calltherm_mars.F90
------> Routine running the sub-timestep of the thermals

A libf/phymars/thermcell_main_mars.F90
------> Main thermals routine specific to Martian physics

A libf/phymars/thermcell_dqupdown.F90
------> Thermals subroutine computing transport of quantities by updrafts and downdrafts

A libf/phymars/thermcell.F90
------> Module including parameters from the Earth to Mars importation. Will disappear in future dev

================================================
================================================

  • Property svn:executable set to *
File size: 4.7 KB
Line 
1MODULE thermcell
2!      USE ioipsl_getincom
3      IMPLICIT NONE
4
5      INTEGER :: iflag_thermals
6      REAL :: r_aspect_thermals,l_mix_thermals,tau_thermals
7      INTEGER :: w2di_thermals
8      INTEGER :: nsplit_thermals
9      INTEGER :: isplit
10      INTEGER :: iflag_coupl,iflag_clos,iflag_wake
11      INTEGER :: iflag_thermals_ed,iflag_thermals_optflux
12      REAL :: RG,RD
13      REAL :: rlvtt,rcpd,rtt,r2es
14      REAL :: retv,r5les,r5ies,rkappa
15      REAL :: R4LES,R4IES,R3IES,R3LES
16      INTEGER :: klon,klev
17      INTEGER :: prt_level,lunout
18      real,allocatable :: rlatd(:)
19      real,allocatable :: rlond(:)
20
21! Thermodynamic constants. The [OK] means that the constant has been made
22! compatible with the Martian model
23! ----------------------------------
24
25! RG : mars mean gravity field                          [OK]
26! RD : dry air constant = 1000 R/Mair                   [OK]
27! rlvtt : value of Lv(Tt) (vapo. latent heat at Tt)     [dep. on air ?]
28! rcpd : Cp of dry air = 7/2 RD for a perfect gas       [OK,CHECK VALUE with gcm one]
29!rcpd = r/(r/CP), avec r/CP = 1/4.4 pour LES, 0.256793 pour gcm (1./3.89419), et r=1000.R/Mair=RD
30! rtt : triple point temperature Tt                     [OK]
31! r2es :                                                [Probably OK]
32! retv = restt*RD/RV :                                  [-]
33! restt : saturation pressure at Tt es(Tt) = 611.14 Pa  [dep. on air ?]
34! RV : vapor constant = 1000 R/Mh2o                     [OK]
35! r3les :                                               [Probably OK]
36! r3ies :                                               [Probably OK]
37! r4les :                                               [Probably OK]
38! r4ies :                                               [Probably OK]
39! r5ies = r3ies*(rtt-r4les)                             [Probably OK]
40! rkappa = RD/RCPD = r/cp = 0.256793                    [OK]
41! FOEEW : vapeur d'eau saturante                        [OK]
42! FOEDE : derive par rapport a la temperature           [OK]
43      PARAMETER (RG=3.72,RD=191.182)
44      PARAMETER (rlvtt=2.5008E+06,rcpd=RD*3.89419)
45      PARAMETER (rtt=273.16,r2es=253.156)
46      PARAMETER (retv=1.41409,r5les=4097.93)
47      PARAMETER (r5ies=5807.81,rkappa=1./3.89419)
48      PARAMETER (R4IES=7.66,R4LES=35.86)
49      PARAMETER (R3IES=21.875,R3LES=17.269)
50
51! r_aspect_thermals : rapport d'aspect : parameter defined in FH's HDR
52!   length over height ratio of the thermals in the alimentation layer (check)
53!   it's value is advised to be set to 2 in HDR           [OK]
54! l_mix_thermals : mixing length => lambda : parameter for plume penetration in
55! the atm (the plum size will decrease with sqrt(lambda) : USELESS
56! w2di_thermals : ?
57! tau_thermals : relaxation length of the thermals : the model implements the thermals
58!   tendancies with a typical relaxation time. Typical value is said to be 1800 s in the
59!   code however, it is found at 720 s in the .def files  [OK]
60! nsplit_thermals : subdivisions of the timestep in the calculations of the plume.
61!   Increase the calculation precision greatly, but also computation time
62!   typical value is found at 1 in .def (no split)        [OK]
63! prt_level : print level for the gcm printed outputs     [OK]
64! lunout : output channel for prt level                   [OK]
65! iflag_thermals : choice of the thermals version. 15 and 16 are the newest ones
66!   and 16 activates "bidouille stratocu" which is required for now
67!                                                         [OK]
68! iflag_thermals_ed : choice of the plume version
69!   8 is the normal version, 9 is the working version of AJ, 10 is CR
70!   also used for cases .ge.1 in thermcell_height         [OK]
71! iflag_coupl : USELESS ?                                 [useless ?]
72! iflag_clos : USELESS ?                                  [useless ?]
73! iflag_wake : USELESS ?                                  [useless ?]
74      PARAMETER (r_aspect_thermals = 3.)
75      PARAMETER (l_mix_thermals = 30.)
76      PARAMETER (w2di_thermals = 1)
77!      PARAMETER (tau_thermals = 720.)
78!      PARAMETER (nsplit_thermals = 1)
79      PARAMETER (prt_level=0,lunout=6)
80      PARAMETER (iflag_thermals=20,iflag_thermals_ed=10)
81      PARAMETER (iflag_thermals_optflux=1,iflag_coupl=1)
82      PARAMETER (iflag_clos=2)
83
84      CONTAINS
85        FUNCTION foede(PTARG,PDELARG,P5ARG,PQSARG,PCOARG)
86          IMPLICIT NONE
87          REAL :: foede
88          REAL, INTENT(IN) :: PTARG,PDELARG,P5ARG,PQSARG,PCOARG
89          foede = PQSARG*PCOARG*P5ARG &
90       & / (PTARG-(R4LES*(1.-PDELARG)+R4IES*PDELARG))**2
91        END FUNCTION foede
92
93        FUNCTION foeew(PTARG,PDELARG)
94          IMPLICIT NONE
95          REAL :: foeew
96          REAL, INTENT(IN) :: PTARG,PDELARG
97          foeew = EXP (                                   &
98       &          (R3LES*(1.-PDELARG)+R3IES*PDELARG) * (PTARG-RTT)        &
99       & / (PTARG-(R4LES*(1.-PDELARG)+R4IES*PDELARG)) )
100
101        END FUNCTION foeew
102
103! Variables which have been moved to .def and called independantly :
104! nsplit_thermals, tau
105
106END MODULE thermcell
107
108
Note: See TracBrowser for help on using the repository browser.