source: trunk/LMDZ.COMMON/libf/dyn3dpar/cpdet_mod.F90 @ 1017

Last change on this file since 1017 was 1017, checked in by emillour, 11 years ago

Common dynamics: (and collateral adaptations in Venus physics)
Improved cpdet routines in and additional sponge mode:

  • Additionnal sponge mode (trigered with "callsponge" flag), in line with the one used in the Generic and Martian GCM. This sponge is called whenever there is a dissipation step.
  • Improvement of the cpdet routines : created routines tpot2t_glo_p and t2tpot_glo_p which handle fields on the whole dynamics (scaler) grid, which are more efficient than calling tpot2t_p or t2tpot_p with slabs of data (generated use of intermediate copies of these chunks of data at every call)
  • Turned cpdet.F into a module cpdet_mod.F90 (and correspondingly adapted all routines in the Venus physics).

EM

File size: 9.0 KB
Line 
1module cpdet_mod
2
3implicit 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
19contains
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 t2tpot_p(nlon,nlev, yt, yteta, ypk)
104! Parallel version of t2tpot, for an arbitrary number of columns
105      USE control_mod, only : planet_type
106      IMPLICIT none
107
108! for cpp, nu_venus and t0_venus:
109#include "comconst.h"
110
111      integer,intent(in) :: nlon,nlev
112      real,intent(in) :: yt(nlon,nlev)
113      real,intent(out) :: yteta(nlon,nlev)
114      real,intent(in) :: ypk(nlon,nlev)
115! local variable:
116      integer :: l
117     
118      if (planet_type.eq."venus") then
119!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
120        do l=1,nlev
121          yteta(:,l)=yt(:,l)**nu_venus                                  &
122     &                     -nu_venus*t0_venus**nu_venus*                &
123     &                          log(ypk(:,l)/cpp)
124          yteta(:,l)=yteta(:,l)**(1./nu_venus)
125        enddo
126!$OMP END DO
127      else
128!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
129        do l=1,nlev
130          yteta(:,l)=yt(:,l)*cpp/ypk(:,l)
131        enddo
132!$OMP END DO
133      endif ! of if (planet_type.eq."venus")
134
135      end subroutine t2tpot_p
136
137!======================================================================
138!======================================================================
139
140      SUBROUTINE t2tpot_glo_p(yt, yteta, ypk)
141! Parallel version of t2tpot, over the full dynamics (scalar) grid
142! (more efficient than multiple calls to t2tpot_p() with slices of data)
143      USE parallel, only : jj_begin,jj_end
144      USE control_mod, only : planet_type
145      IMPLICIT none
146! for iip1, jjp1 and llm
147#include "dimensions.h"
148#include "paramet.h"
149! for cpp, nu_venus and t0_venus:
150#include "comconst.h"
151
152      real,intent(in) :: yt(iip1,jjp1,llm)
153      real,intent(out) :: yteta(iip1,jjp1,llm)
154      real,intent(in) :: ypk(iip1,jjp1,llm)
155! local variable:
156      integer :: j,l
157      integer :: jjb,jje
158     
159      jjb=jj_begin
160      jje=jj_end
161
162      if (planet_type.eq."venus") then
163!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
164        do l=1,llm
165          yteta(:,jjb:jje,l)=yt(:,jjb:jje,l)**nu_venus                  &
166     &                     -nu_venus*t0_venus**nu_venus*                &
167     &                          log(ypk(:,jjb:jje,l)/cpp)
168          yteta(:,jjb:jje,l)=yteta(:,jjb:jje,l)**(1./nu_venus)
169        enddo
170!$OMP END DO
171      else
172!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
173        do l=1,llm
174          yteta(:,jjb:jje,l)=yt(:,jjb:jje,l)*cpp/ypk(:,jjb:jje,l)
175        enddo
176!$OMP END DO
177      endif ! of if (planet_type.eq."venus")
178
179      end subroutine t2tpot_glo_p
180
181!======================================================================
182
183      SUBROUTINE tpot2t(npoints,yteta, yt, ypk)
184!======================================================================
185! Arguments:
186!
187! yteta--------input-R- Temperature potentielle
188! yt   -------output-R- Temperature
189! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
190!
191!======================================================================
192
193      USE control_mod, ONLY: planet_type
194      IMPLICIT NONE
195
196! for cpp, nu_venus and t0_venus:
197#include "comconst.h"
198
199      integer,intent(in) :: npoints
200      REAL,intent(in) :: yteta(npoints), ypk(npoints)
201      REAL,intent(out) :: yt(npoints)
202     
203      if (planet_type.eq."venus") then
204         yt = yteta**nu_venus                                           &
205     &       + nu_venus * t0_venus**nu_venus * log(ypk/cpp)
206         yt = yt**(1./nu_venus)
207      else
208          yt = yteta * ypk/cpp
209      endif
210 
211      return
212      end subroutine tpot2t
213
214!======================================================================
215!======================================================================
216
217      SUBROUTINE tpot2t_p(nlon,nlev,yteta,yt,ypk)
218! Parallel version of tpot2t, for an arbitrary number of columns
219      USE control_mod, only : planet_type
220      IMPLICIT none
221! for cpp, nu_venus and t0_venus:
222#include "comconst.h"
223
224      integer,intent(in) :: nlon,nlev
225      real,intent(out) :: yt(nlon,nlev)
226      real,intent(in) :: yteta(nlon,nlev)
227      real,intent(in) :: ypk(nlon,nlev)
228
229! local variable:
230      integer :: l
231
232      if (planet_type.eq."venus") then
233!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
234        do l=1,nlev
235          yt(:,l)=yteta(:,l)**nu_venus                                  &
236     &                  +nu_venus*t0_venus**nu_venus*                   &
237     &                       log(ypk(:,l)/cpp)
238          yt(:,l)=yt(:,l)**(1./nu_venus)
239        enddo
240!$OMP END DO
241      else
242!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
243        do l=1,nlev
244          yt(:,l)=yteta(:,l)*ypk(:,l)/cpp
245        enddo
246!$OMP END DO
247      endif ! of if (planet_type.eq."venus")
248      end subroutine tpot2t_p
249
250!======================================================================
251!======================================================================
252
253      SUBROUTINE tpot2t_glo_p(yteta,yt,ypk)
254! Parallel version of tpot2t, over the full dynamics (scalar) grid
255! (more efficient than multiple calls to tpot2t_p() with slices of data)
256      USE parallel, only : jj_begin,jj_end
257      USE control_mod, only : planet_type
258      IMPLICIT none
259! for iip1, jjp1 and llm
260#include "dimensions.h"
261#include "paramet.h"
262! for cpp, nu_venus and t0_venus:
263#include "comconst.h"
264
265      real,intent(out) :: yt(iip1,jjp1,llm)
266      real,intent(in) :: yteta(iip1,jjp1,llm)
267      real,intent(in) :: ypk(iip1,jjp1,llm)
268! local variable:
269      integer :: j,l
270      integer :: jjb,jje
271     
272      jjb=jj_begin
273      jje=jj_end
274
275      if (planet_type.eq."venus") then
276!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
277        do l=1,llm
278          yt(:,jjb:jje,l)=yteta(:,jjb:jje,l)**nu_venus                  &
279     &                  +nu_venus*t0_venus**nu_venus*                   &
280     &                       log(ypk(:,jjb:jje,l)/cpp)
281          yt(:,jjb:jje,l)=yt(:,jjb:jje,l)**(1./nu_venus)
282        enddo
283!$OMP END DO
284      else
285!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
286        do l=1,llm
287          yt(:,jjb:jje,l)=yteta(:,jjb:jje,l)*ypk(:,jjb:jje,l)/cpp
288        enddo
289!$OMP END DO
290      endif ! of if (planet_type.eq."venus")
291      end subroutine tpot2t_glo_p
292
293!======================================================================
294!======================================================================
295!
296! ATTENTION
297!
298! Si un jour on a besoin, il faudra coder les routines
299!    dt2dtpot / dtpto2dt
300!
301!======================================================================
302!======================================================================
303end module cpdet_mod
Note: See TracBrowser for help on using the repository browser.