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

Last change on this file was 3826, checked in by emillour, 5 months ago

Venus PCM:

Some preliminary tidying with respect to OpenMP when initializing the
physics. Turned suphec into a module in the process.
EM

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