1 | module cpdet_mod |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | ! ADAPTATION OF GCM TO CP(T) |
---|
6 | !====================================================================== |
---|
7 | ! S. Lebonnois, 10/2010 |
---|
8 | ! |
---|
9 | ! Cp must be computed using cpdet(t) to be valid |
---|
10 | ! |
---|
11 | ! The Exner function is still pk = RCPD*(play/pref)**RKAPPA |
---|
12 | ! (RCPD=cpp, RKAPPA=kappa) |
---|
13 | ! |
---|
14 | ! One goes from T to teta (potential temperature) using t2tpot(t,teta,pk) |
---|
15 | ! One goes from teta to T using tpot2t(teta,t,pk) |
---|
16 | ! |
---|
17 | !====================================================================== |
---|
18 | |
---|
19 | contains |
---|
20 | |
---|
21 | SUBROUTINE ini_cpdet |
---|
22 | |
---|
23 | USE control_mod, ONLY: planet_type |
---|
24 | IMPLICIT none |
---|
25 | !====================================================================== |
---|
26 | ! Initialization of nu_venus and t0_venus |
---|
27 | !====================================================================== |
---|
28 | |
---|
29 | ! for cpp, nu_venus and t0_venus: |
---|
30 | #include "comconst.h" |
---|
31 | |
---|
32 | if (planet_type.eq."venus") then |
---|
33 | nu_venus=0.35 |
---|
34 | t0_venus=460. |
---|
35 | else |
---|
36 | nu_venus=0. |
---|
37 | t0_venus=0. |
---|
38 | endif |
---|
39 | |
---|
40 | return |
---|
41 | end subroutine ini_cpdet |
---|
42 | |
---|
43 | !====================================================================== |
---|
44 | !====================================================================== |
---|
45 | |
---|
46 | FUNCTION cpdet(t) |
---|
47 | |
---|
48 | USE control_mod, ONLY: planet_type |
---|
49 | IMPLICIT none |
---|
50 | |
---|
51 | ! for cpp, nu_venus and t0_venus: |
---|
52 | #include "comconst.h" |
---|
53 | |
---|
54 | real,intent(in) :: t |
---|
55 | real cpdet |
---|
56 | |
---|
57 | if (planet_type.eq."venus") then |
---|
58 | cpdet = cpp*(t/t0_venus)**nu_venus |
---|
59 | else |
---|
60 | cpdet = cpp |
---|
61 | endif |
---|
62 | |
---|
63 | return |
---|
64 | end function cpdet |
---|
65 | |
---|
66 | !====================================================================== |
---|
67 | !====================================================================== |
---|
68 | |
---|
69 | SUBROUTINE t2tpot(npoints, yt, yteta, ypk) |
---|
70 | !====================================================================== |
---|
71 | ! Arguments: |
---|
72 | ! |
---|
73 | ! yt --------input-R- Temperature |
---|
74 | ! yteta-------output-R- Temperature potentielle |
---|
75 | ! ypk --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA |
---|
76 | ! |
---|
77 | !====================================================================== |
---|
78 | |
---|
79 | USE control_mod, ONLY: planet_type |
---|
80 | IMPLICIT NONE |
---|
81 | |
---|
82 | ! for cpp, nu_venus and t0_venus: |
---|
83 | #include "comconst.h" |
---|
84 | |
---|
85 | integer,intent(in) :: npoints |
---|
86 | REAL,intent(in) :: yt(npoints), ypk(npoints) |
---|
87 | REAL,intent(out) :: yteta(npoints) |
---|
88 | |
---|
89 | if (planet_type.eq."venus") then |
---|
90 | yteta = yt**nu_venus & |
---|
91 | & - nu_venus * t0_venus**nu_venus * log(ypk/cpp) |
---|
92 | yteta = yteta**(1./nu_venus) |
---|
93 | else |
---|
94 | yteta = yt * cpp/ypk |
---|
95 | endif |
---|
96 | |
---|
97 | return |
---|
98 | end subroutine t2tpot |
---|
99 | |
---|
100 | !====================================================================== |
---|
101 | !====================================================================== |
---|
102 | |
---|
103 | SUBROUTINE tpot2t(npoints,yteta, yt, ypk) |
---|
104 | !====================================================================== |
---|
105 | ! Arguments: |
---|
106 | ! |
---|
107 | ! yteta--------input-R- Temperature potentielle |
---|
108 | ! yt -------output-R- Temperature |
---|
109 | ! ypk --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA |
---|
110 | ! |
---|
111 | !====================================================================== |
---|
112 | |
---|
113 | USE control_mod, ONLY: planet_type |
---|
114 | IMPLICIT NONE |
---|
115 | |
---|
116 | ! for cpp, nu_venus and t0_venus: |
---|
117 | #include "comconst.h" |
---|
118 | |
---|
119 | integer,intent(in) :: npoints |
---|
120 | REAL,intent(in) :: yteta(npoints), ypk(npoints) |
---|
121 | REAL,intent(out) :: yt(npoints) |
---|
122 | |
---|
123 | if (planet_type.eq."venus") then |
---|
124 | yt = yteta**nu_venus & |
---|
125 | & + nu_venus * t0_venus**nu_venus * log(ypk/cpp) |
---|
126 | yt = yt**(1./nu_venus) |
---|
127 | else |
---|
128 | yt = yteta * ypk/cpp |
---|
129 | endif |
---|
130 | |
---|
131 | return |
---|
132 | end subroutine tpot2t |
---|
133 | |
---|
134 | !====================================================================== |
---|
135 | !====================================================================== |
---|
136 | ! Routines pour les calculs paralleles |
---|
137 | !====================================================================== |
---|
138 | #ifdef CPP_PARA |
---|
139 | !====================================================================== |
---|
140 | |
---|
141 | SUBROUTINE t2tpot_p(nlon,nlev, yt, yteta, ypk) |
---|
142 | ! Parallel version of t2tpot, for an arbitrary number of columns |
---|
143 | USE control_mod, only : planet_type |
---|
144 | USE parallel_lmdz, only : OMP_CHUNK |
---|
145 | IMPLICIT none |
---|
146 | |
---|
147 | ! for cpp, nu_venus and t0_venus: |
---|
148 | #include "comconst.h" |
---|
149 | |
---|
150 | integer,intent(in) :: nlon,nlev |
---|
151 | real,intent(in) :: yt(nlon,nlev) |
---|
152 | real,intent(out) :: yteta(nlon,nlev) |
---|
153 | real,intent(in) :: ypk(nlon,nlev) |
---|
154 | ! local variable: |
---|
155 | integer :: l |
---|
156 | |
---|
157 | if (planet_type.eq."venus") then |
---|
158 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
159 | do l=1,nlev |
---|
160 | yteta(:,l)=yt(:,l)**nu_venus & |
---|
161 | & -nu_venus*t0_venus**nu_venus* & |
---|
162 | & log(ypk(:,l)/cpp) |
---|
163 | yteta(:,l)=yteta(:,l)**(1./nu_venus) |
---|
164 | enddo |
---|
165 | !$OMP END DO |
---|
166 | else |
---|
167 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
168 | do l=1,nlev |
---|
169 | yteta(:,l)=yt(:,l)*cpp/ypk(:,l) |
---|
170 | enddo |
---|
171 | !$OMP END DO |
---|
172 | endif ! of if (planet_type.eq."venus") |
---|
173 | |
---|
174 | end subroutine t2tpot_p |
---|
175 | |
---|
176 | !====================================================================== |
---|
177 | !====================================================================== |
---|
178 | |
---|
179 | SUBROUTINE t2tpot_glo_p(yt, yteta, ypk) |
---|
180 | ! Parallel version of t2tpot, over the full dynamics (scalar) grid |
---|
181 | ! (more efficient than multiple calls to t2tpot_p() with slices of data) |
---|
182 | USE parallel_lmdz, only : jj_begin,jj_end,OMP_CHUNK |
---|
183 | USE control_mod, only : planet_type |
---|
184 | IMPLICIT none |
---|
185 | ! for iip1, jjp1 and llm |
---|
186 | #include "dimensions.h" |
---|
187 | #include "paramet.h" |
---|
188 | ! for cpp, nu_venus and t0_venus: |
---|
189 | #include "comconst.h" |
---|
190 | |
---|
191 | real,intent(in) :: yt(iip1,jjp1,llm) |
---|
192 | real,intent(out) :: yteta(iip1,jjp1,llm) |
---|
193 | real,intent(in) :: ypk(iip1,jjp1,llm) |
---|
194 | ! local variable: |
---|
195 | integer :: j,l |
---|
196 | integer :: jjb,jje |
---|
197 | |
---|
198 | jjb=jj_begin |
---|
199 | jje=jj_end |
---|
200 | |
---|
201 | if (planet_type.eq."venus") then |
---|
202 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
203 | do l=1,llm |
---|
204 | yteta(:,jjb:jje,l)=yt(:,jjb:jje,l)**nu_venus & |
---|
205 | & -nu_venus*t0_venus**nu_venus* & |
---|
206 | & log(ypk(:,jjb:jje,l)/cpp) |
---|
207 | yteta(:,jjb:jje,l)=yteta(:,jjb:jje,l)**(1./nu_venus) |
---|
208 | enddo |
---|
209 | !$OMP END DO |
---|
210 | else |
---|
211 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
212 | do l=1,llm |
---|
213 | yteta(:,jjb:jje,l)=yt(:,jjb:jje,l)*cpp/ypk(:,jjb:jje,l) |
---|
214 | enddo |
---|
215 | !$OMP END DO |
---|
216 | endif ! of if (planet_type.eq."venus") |
---|
217 | |
---|
218 | end subroutine t2tpot_glo_p |
---|
219 | |
---|
220 | !====================================================================== |
---|
221 | !====================================================================== |
---|
222 | |
---|
223 | SUBROUTINE tpot2t_p(nlon,nlev,yteta,yt,ypk) |
---|
224 | ! Parallel version of tpot2t, for an arbitrary number of columns |
---|
225 | USE control_mod, only : planet_type |
---|
226 | USE parallel_lmdz, only : OMP_CHUNK |
---|
227 | IMPLICIT none |
---|
228 | ! for cpp, nu_venus and t0_venus: |
---|
229 | #include "comconst.h" |
---|
230 | |
---|
231 | integer,intent(in) :: nlon,nlev |
---|
232 | real,intent(out) :: yt(nlon,nlev) |
---|
233 | real,intent(in) :: yteta(nlon,nlev) |
---|
234 | real,intent(in) :: ypk(nlon,nlev) |
---|
235 | |
---|
236 | ! local variable: |
---|
237 | integer :: l |
---|
238 | |
---|
239 | if (planet_type.eq."venus") then |
---|
240 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
241 | do l=1,nlev |
---|
242 | yt(:,l)=yteta(:,l)**nu_venus & |
---|
243 | & +nu_venus*t0_venus**nu_venus* & |
---|
244 | & log(ypk(:,l)/cpp) |
---|
245 | yt(:,l)=yt(:,l)**(1./nu_venus) |
---|
246 | enddo |
---|
247 | !$OMP END DO |
---|
248 | else |
---|
249 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
250 | do l=1,nlev |
---|
251 | yt(:,l)=yteta(:,l)*ypk(:,l)/cpp |
---|
252 | enddo |
---|
253 | !$OMP END DO |
---|
254 | endif ! of if (planet_type.eq."venus") |
---|
255 | end subroutine tpot2t_p |
---|
256 | |
---|
257 | !====================================================================== |
---|
258 | !====================================================================== |
---|
259 | |
---|
260 | SUBROUTINE tpot2t_glo_p(yteta,yt,ypk) |
---|
261 | ! Parallel version of tpot2t, over the full dynamics (scalar) grid |
---|
262 | ! (more efficient than multiple calls to tpot2t_p() with slices of data) |
---|
263 | USE parallel_lmdz, only : jj_begin,jj_end,OMP_CHUNK |
---|
264 | USE control_mod, only : planet_type |
---|
265 | IMPLICIT none |
---|
266 | ! for iip1, jjp1 and llm |
---|
267 | #include "dimensions.h" |
---|
268 | #include "paramet.h" |
---|
269 | ! for cpp, nu_venus and t0_venus: |
---|
270 | #include "comconst.h" |
---|
271 | |
---|
272 | real,intent(out) :: yt(iip1,jjp1,llm) |
---|
273 | real,intent(in) :: yteta(iip1,jjp1,llm) |
---|
274 | real,intent(in) :: ypk(iip1,jjp1,llm) |
---|
275 | ! local variable: |
---|
276 | integer :: j,l |
---|
277 | integer :: jjb,jje |
---|
278 | |
---|
279 | jjb=jj_begin |
---|
280 | jje=jj_end |
---|
281 | |
---|
282 | if (planet_type.eq."venus") then |
---|
283 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
284 | do l=1,llm |
---|
285 | yt(:,jjb:jje,l)=yteta(:,jjb:jje,l)**nu_venus & |
---|
286 | & +nu_venus*t0_venus**nu_venus* & |
---|
287 | & log(ypk(:,jjb:jje,l)/cpp) |
---|
288 | yt(:,jjb:jje,l)=yt(:,jjb:jje,l)**(1./nu_venus) |
---|
289 | enddo |
---|
290 | !$OMP END DO |
---|
291 | else |
---|
292 | !$OMP DO SCHEDULE(STATIC,OMP_CHUNK) |
---|
293 | do l=1,llm |
---|
294 | yt(:,jjb:jje,l)=yteta(:,jjb:jje,l)*ypk(:,jjb:jje,l)/cpp |
---|
295 | enddo |
---|
296 | !$OMP END DO |
---|
297 | endif ! of if (planet_type.eq."venus") |
---|
298 | end subroutine tpot2t_glo_p |
---|
299 | |
---|
300 | !====================================================================== |
---|
301 | #endif |
---|
302 | !====================================================================== |
---|
303 | ! Fin routines specifiques parallele |
---|
304 | !====================================================================== |
---|
305 | !====================================================================== |
---|
306 | ! |
---|
307 | ! ATTENTION |
---|
308 | ! |
---|
309 | ! Si un jour on a besoin, il faudra coder les routines |
---|
310 | ! dt2dtpot / dtpto2dt |
---|
311 | ! |
---|
312 | !====================================================================== |
---|
313 | !====================================================================== |
---|
314 | end module cpdet_mod |
---|