source: trunk/LMDZ.VENUS/libf/phyvenus/cpdet_phy_mod.F90 @ 3137

Last change on this file since 3137 was 2135, checked in by slebonnois, 6 years ago

SL, Venus: new keys for flexibility cp0/cp(T) and Held-Suarez type physics

File size: 4.4 KB
Line 
1module cpdet_phy_mod
2
3! Cpdet module to use in physics (same as Cpdet module in dynamics)
4! Warning: Changing hard coded formulaes in this module should also be
5!          done in dynamics cpdet_mod module  (and vice-versa) !
6
7implicit none
8
9real,save :: cpp ! reference Cp
10real,save :: nu_venus
11real,save :: t0_venus
12
13contains
14
15      SUBROUTINE init_cpdet_phy(cpp_,nu_venus_,t0_venus_)
16      ! initialize module variables
17      REAL,INTENT(IN) :: cpp_ ! cpp from dynamics
18      REAL,INTENT(IN) :: nu_venus_ ! nu_venus from dynamics
19      REAL,INTENT(IN) :: t0_venus_ ! t0_venus from dynamics
20     
21      cpp=cpp_
22      nu_venus=nu_venus_
23      t0_venus=t0_venus_
24     
25      END SUBROUTINE init_cpdet_phy
26
27!======================================================================
28
29      FUNCTION cpdet(t)
30
31      IMPLICIT none
32
33! for cpp, nu_venus and t0_venus:
34
35      real,intent(in) :: t
36      real cpdet
37
38      if (t0_venus.ne.0.) then
39       cpdet = cpp*(t/t0_venus)**nu_venus
40      else
41       cpdet = cpp
42      endif
43
44      end function cpdet
45     
46!======================================================================
47
48      SUBROUTINE t2tpot(npoints, yt, yteta, ypk)
49!======================================================================
50! Arguments:
51!
52! yt   --------input-R- Temperature
53! yteta-------output-R- Temperature potentielle
54! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
55!
56!======================================================================
57
58      IMPLICIT NONE
59
60      integer,intent(in) :: npoints
61      REAL,intent(in) :: yt(npoints), ypk(npoints)
62      REAL,intent(out) :: yteta(npoints)
63     
64!----------------------
65! ATMOSPHERE PROFONDE
66      real ypklim,ypklim2,mmm0
67      real ratio_mod(npoints)
68      integer k
69! theta_mod = theta * ratio_mod => ratio_mod = mmm0/mmm(p)
70
71      ypklim  = cpp*(  6e6/9.2e6)**0.1914
72      ypklim2 = cpp*(8.9e6/9.2e6)**0.1914
73      mmm0 = 43.44
74     
75      DO k = 1, npoints
76        ratio_mod(k) = 1.
77! ATM PROFONDE DESACTIVEE !
78        if (   (1 .EQ. 0)  .AND.(ypk(k).gt.ypklim)) then
79           ratio_mod(k) = mmm0 /                                        &
80     &  (mmm0+0.56*(log(ypklim/cpp)-log(ypk(k)/cpp))                    &
81     &            /(log(ypklim/cpp)-log(ypklim2/cpp)))
82           ratio_mod(k) = max(ratio_mod(k),mmm0/(mmm0+0.56))
83        endif
84      ENDDO
85!----------------------
86
87      if (t0_venus.ne.0.) then
88       yteta = yt**nu_venus                                          &
89     &            - nu_venus * t0_venus**nu_venus * log(ypk/cpp)
90       yteta = yteta**(1./nu_venus)
91      else
92       yteta = yt * cpp/ypk
93      endif
94
95!----------------------
96! ATMOSPHERE PROFONDE
97      yteta = yteta*ratio_mod
98!----------------------
99
100      end subroutine t2tpot
101
102!======================================================================
103
104      SUBROUTINE tpot2t(npoints,yteta, yt, ypk)
105!======================================================================
106! Arguments:
107!
108! yteta--------input-R- Temperature potentielle
109! yt   -------output-R- Temperature
110! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
111!
112!======================================================================
113
114      IMPLICIT NONE
115
116      integer,intent(in) :: npoints
117      REAL,intent(in) :: yteta(npoints), ypk(npoints)
118      REAL,intent(out) :: yt(npoints)
119     
120!----------------------
121! ATMOSPHERE PROFONDE
122      real ypklim,ypklim2,mmm0
123      real ratio_mod(npoints)
124      integer k
125! theta_mod = theta * ratio_mod => ratio_mod = mmm0/mmm(p)
126
127      ypklim  = cpp*(  6e6/9.2e6)**0.1914
128      ypklim2 = cpp*(8.9e6/9.2e6)**0.1914
129      mmm0 = 43.44
130
131      DO k = 1, npoints
132        ratio_mod(k) = 1.
133! ATM PROFONDE DESACTIVEE !
134        if (   (1 .EQ. 0)  .AND.(ypk(k).gt.ypklim)) then
135           ratio_mod(k) = mmm0 /                                        &
136     &  (mmm0+0.56*(log(ypklim/cpp)-log(ypk(k)/cpp))                    &
137     &            /(log(ypklim/cpp)-log(ypklim2/cpp)))
138           ratio_mod(k) = max(ratio_mod(k),mmm0/(mmm0+0.56))
139        endif
140      ENDDO
141!----------------------
142
143!----------------------
144! ATMOSPHERE PROFONDE
145!----------------------
146      if (t0_venus.ne.0.) then
147     yt = (yteta/ratio_mod)**nu_venus                               &
148!----------------------
149         + nu_venus * t0_venus**nu_venus * log(ypk/cpp)
150     yt = yt**(1./nu_venus)
151      else
152     yt = (yteta/ratio_mod) * ypk/cpp
153      endif
154 
155 
156      end subroutine tpot2t
157
158end module cpdet_phy_mod
Note: See TracBrowser for help on using the repository browser.