source: trunk/LMDZ.GENERIC/libf/dyn3d/sponge_mod.F90 @ 1422

Last change on this file since 1422 was 1422, checked in by milmd, 10 years ago

In GENERIC, MARS and COMMON models replace some include files by modules (usefull for decoupling physics with dynamics).

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