source: trunk/LMDZ.COMMON/libf/dyn3dpar/sponge_mod_p.F90 @ 3553

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