source: trunk/LMDZ.MARS/libf/dyn3d/sponge_mod.F90 @ 2001

Last change on this file since 2001 was 1593, checked in by emillour, 9 years ago

Common dynamics:

  • some cleanup around unused or unecessary parameters in the dynamics (ecritphy, grireg and physic) to harmonize with LMDZ5 (Earth dyn core).
  • associated changes in LMDZ.MARS and LMDZ.GENERIC

EM

File size: 4.5 KB
Line 
1module sponge_mod
2
3implicit none
4
5! sponge parameters (set/read via conf_gcm.F)
6logical,save :: callsponge  ! do we use a sponge on upper layers
7integer,save :: mode_sponge ! sponge mode
8integer,save :: nsponge ! number of sponge layers
9real,save :: tetasponge  ! sponge time scale (s) at topmost layer
10
11
12contains
13
14      subroutine sponge(ucov,vcov,h,ps,dt,mode)
15
16! Sponge routine: Quench ucov, vcov and potential temperature near the
17!                 top of the model
18! Depending on 'mode' relaxation of variables is towards:
19! mode = 0 : h -> h_mean , ucov -> 0 , vcov -> 0
20! mode = 1 : h -> h_mean , ucov -> ucov_mean , vcov -> 0
21! mode >= 2 : h -> h_mean , ucov -> ucov_mean , vcov -> vcov_mean
22! Number of layer over which sponge is applied is 'nsponge' (read from def file)
23! Time scale for quenching at top level is given by 'tetasponge' (read from
24! def file) and doubles as level indexes decrease.
25! Quenching is modeled as: A(t)=Am+A0exp(-lambda*t)
26! where Am is the zonal average of the field (or zero), and lambda the inverse
27! of the characteristic quenching/relaxation time scale
28! Thus, assuming Am to be time-independent, field at time t+dt is given by:
29! A(t+dt)=A(t)-(A(t)-Am)*(1-exp(-lambda*dt))
30
31      USE comvert_mod, ONLY: ap,bp,preff
32
33      implicit none
34#include "dimensions.h"
35#include "paramet.h"
36#include "comdissip.h"
37#include "comgeom2.h"
38
39! Arguments:
40!------------
41      real,intent(inout) :: ucov(iip1,jjp1,llm) ! covariant zonal wind
42      real,intent(inout) :: vcov(iip1,jjm,llm) ! covariant meridional wind
43      real,intent(inout) :: h(iip1,jjp1,llm) ! potential temperature
44!      real,intent(in) :: pext(iip1,jjp1) ! extensive pressure
45      real,intent(in) :: ps(iip1,jjp1) ! surface pressure (Pa)
46      real,intent(in) :: dt   ! time step
47      integer,intent(in) :: mode  ! sponge mode
48
49!   Local:
50!   ------
51
52      real,save :: sig_s(llm) !sigma au milieu des couches
53      REAL vm,um,hm,ptot(jjp1)
54      real,save :: cst(llm)
55      real :: pext(iip1,jjp1) ! extensive pressure
56
57      INTEGER l,i,j
58      integer,save :: l0 ! layer down to which sponge is applied
59
60      real ssum
61
62      real echelle,zkm
63      logical,save :: firstcall=.true.
64
65
66
67      if (firstcall) then
68
69       ! build approximative sigma levels at midlayer
70        do l=1,llm
71          sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2.
72        enddo
73
74        l0=llm-nsponge+1
75 
76        PRINT*
77        print*,'sponge mode',mode
78        print*,'nsponge tetasponge ',nsponge,tetasponge
79        print*,'Coeffs for the sponge layer'
80        print*,'Z (km)     tau      cst'
81        do l=llm,l0,-1
82          ! double time scale with every level, starting from the top
83          cst(l)=dt/(tetasponge*2**(llm-l))
84        enddo
85
86        echelle=10.
87        do l=l0,llm
88           zkm=-echelle*log(sig_s(l))
89           print*,zkm,dt/cst(l),cst(l)
90        enddo
91        PRINT*
92
93        firstcall=.false.
94      endif ! of if (firstcall)
95
96!-----------------------------------------------------------------------
97!   calcul de la dissipation:
98!   -------------------------
99
100      pext(1:iip1,1:jjp1)=ps(1:iip1,1:jjp1)*aire(1:iip1,1:jjp1)
101
102      do j=1,jjp1
103         ptot(j)=ssum(iim,pext(1,j),1)
104      enddo
105 
106!   potential temperature
107      do l=l0,llm
108         do j=1,jjp1
109            hm=0.
110            do i=1,iim
111               hm=hm+h(i,j,l)*pext(i,j)
112            enddo
113            hm=hm/ptot(j)
114            do i=1,iim
115               h(i,j,l)=h(i,j,l)-cst(l)*(h(i,j,l)-hm)
116            enddo
117            h(iip1,j,l)=h(1,j,l)
118         enddo
119      enddo
120
121!   zonal wind
122      do l=l0,llm
123         do j=2,jjm
124            um=0.
125            if(mode.ge.1) then
126               do i=1,iim
127                  um=um+0.5*ucov(i,j,l)*(pext(i,j)+pext(i+1,j)) &
128                     /cu(i,j)
129               enddo
130               um=um/ptot(j)
131            endif
132            do i=1,iim
133               ucov(i,j,l)=ucov(i,j,l)-cst(l)*(ucov(i,j,l)-um*cu(i,j))
134            enddo
135            ucov(iip1,j,l)=ucov(1,j,l)
136         enddo
137      enddo
138
139!  meridional wind
140      do l=l0,llm
141         do j=1,jjm
142            vm=0.
143            if(mode.ge.2) then
144               do i=1,iim
145                  vm=vm+vcov(i,j,l)*(pext(i,j)+pext(i,j+1)) &
146                     /cv(i,j)
147               enddo
148               vm=vm/(ptot(j)+ptot(j+1))
149            endif
150            do i=1,iim
151               vcov(i,j,l)=vcov(i,j,l)-cst(l)*(vcov(i,j,l)-vm*cv(i,j))
152            enddo
153            vcov(iip1,j,l)=vcov(1,j,l)
154         enddo
155      enddo
156
157      end subroutine sponge
158     
159end module sponge_mod
160
Note: See TracBrowser for help on using the repository browser.