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

Last change on this file since 1723 was 1658, checked in by slebonnois, 8 years ago

SL: deep atmosphere commented until further work

File size: 4.1 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      cpdet = cpp*(t/t0_venus)**nu_venus
39
40      end function cpdet
41     
42!======================================================================
43
44      SUBROUTINE t2tpot(npoints, yt, yteta, ypk)
45!======================================================================
46! Arguments:
47!
48! yt   --------input-R- Temperature
49! yteta-------output-R- Temperature potentielle
50! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
51!
52!======================================================================
53
54      IMPLICIT NONE
55
56      integer,intent(in) :: npoints
57      REAL,intent(in) :: yt(npoints), ypk(npoints)
58      REAL,intent(out) :: yteta(npoints)
59     
60!----------------------
61! ATMOSPHERE PROFONDE
62      real ypklim,ypklim2,mmm0
63      real ratio_mod(npoints)
64      integer k
65! theta_mod = theta * ratio_mod => ratio_mod = mmm0/mmm(p)
66
67      ypklim  = cpp*(  6e6/9.2e6)**0.1914
68      ypklim2 = cpp*(8.9e6/9.2e6)**0.1914
69      mmm0 = 43.44
70     
71      DO k = 1, npoints
72        ratio_mod(k) = 1.
73! ATM PROFONDE DESACTIVEE !
74        if (   (1 .EQ. 0)  .AND.(ypk(k).gt.ypklim)) then
75           ratio_mod(k) = mmm0 /                                        &
76     &  (mmm0+0.56*(log(ypklim/cpp)-log(ypk(k)/cpp))                    &
77     &            /(log(ypklim/cpp)-log(ypklim2/cpp)))
78           ratio_mod(k) = max(ratio_mod(k),mmm0/(mmm0+0.56))
79        endif
80      ENDDO
81!----------------------
82
83      yteta = yt**nu_venus                                          &
84     &            - nu_venus * t0_venus**nu_venus * log(ypk/cpp)
85      yteta = yteta**(1./nu_venus)
86!----------------------
87! ATMOSPHERE PROFONDE
88      yteta = yteta*ratio_mod
89!----------------------
90
91      end subroutine t2tpot
92
93!======================================================================
94
95      SUBROUTINE tpot2t(npoints,yteta, yt, ypk)
96!======================================================================
97! Arguments:
98!
99! yteta--------input-R- Temperature potentielle
100! yt   -------output-R- Temperature
101! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
102!
103!======================================================================
104
105      IMPLICIT NONE
106
107      integer,intent(in) :: npoints
108      REAL,intent(in) :: yteta(npoints), ypk(npoints)
109      REAL,intent(out) :: yt(npoints)
110     
111!----------------------
112! ATMOSPHERE PROFONDE
113      real ypklim,ypklim2,mmm0
114      real ratio_mod(npoints)
115      integer k
116! theta_mod = theta * ratio_mod => ratio_mod = mmm0/mmm(p)
117
118      ypklim  = cpp*(  6e6/9.2e6)**0.1914
119      ypklim2 = cpp*(8.9e6/9.2e6)**0.1914
120      mmm0 = 43.44
121
122      DO k = 1, npoints
123        ratio_mod(k) = 1.
124! ATM PROFONDE DESACTIVEE !
125        if (   (1 .EQ. 0)  .AND.(ypk(k).gt.ypklim)) then
126           ratio_mod(k) = mmm0 /                                        &
127     &  (mmm0+0.56*(log(ypklim/cpp)-log(ypk(k)/cpp))                    &
128     &            /(log(ypklim/cpp)-log(ypklim2/cpp)))
129           ratio_mod(k) = max(ratio_mod(k),mmm0/(mmm0+0.56))
130        endif
131      ENDDO
132!----------------------
133
134!----------------------
135! ATMOSPHERE PROFONDE
136!----------------------
137     yt = (yteta/ratio_mod)**nu_venus                               &
138!----------------------
139         + nu_venus * t0_venus**nu_venus * log(ypk/cpp)
140     yt = yt**(1./nu_venus)
141 
142      end subroutine tpot2t
143
144end module cpdet_phy_mod
Note: See TracBrowser for help on using the repository browser.