source: trunk/LMDZ.MARS/libf/dyn3d/sponge.F @ 937

Last change on this file since 937 was 575, checked in by emillour, 13 years ago

Mars GCM:

correction in nirco2abs: io and ico2 should be declared as integers
Update sponge: remove posibility of specifying 'hsponge', all modes

apply to the last upper "nsponge" layers (default nsponge=3) and
sponge time scale doubles with every level, starting from top.

EM

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