source: trunk/LMDZ.COMMON/libf/dyn3d/top_bound.F @ 1441

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

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

File size: 6.1 KB
Line 
1!
2! $Id: top_bound.F 1793 2013-07-18 07:13:18Z emillour $
3!
4      SUBROUTINE top_bound(vcov,ucov,teta,masse,dt,ducov)
5
6      USE comvert_mod, ONLY: presnivs,scaleheight,preff
7      USE comconst_mod, ONLY: iflag_top_bound,tau_top_bound,
8     .                  mode_top_bound
9
10      IMPLICIT NONE
11c
12#include "dimensions.h"
13#include "paramet.h"
14#include "comgeom2.h"
15
16
17c ..  DISSIPATION LINEAIRE A HAUT NIVEAU, RUN MESO,
18C     F. LOTT DEC. 2006
19c                                 (  10/12/06  )
20
21c=======================================================================
22c
23c   Auteur:  F. LOTT 
24c   -------
25c
26c   Objet:
27c   ------
28c
29c   Dissipation linéaire (ex top_bound de la physique)
30c
31c=======================================================================
32
33! top_bound sponge layer model:
34! Quenching is modeled as: A(t)=Am+A0*exp(-lambda*t)
35! where Am is the zonal average of the field (or zero), and lambda the inverse
36! of the characteristic quenching/relaxation time scale
37! Thus, assuming Am to be time-independent, field at time t+dt is given by:
38! A(t+dt)=A(t)-(A(t)-Am)*(1-exp(-lambda*dt))
39! Moreover lambda can be a function of model level (see below), and relaxation
40! can be toward the average zonal field or just zero (see below).
41
42! NB: top_bound sponge is only called from leapfrog if ok_strato=.true.
43
44! sponge parameters: (loaded/set in conf_gcm.F ; stored in comconst.h)
45!    iflag_top_bound=0 for no sponge
46!    iflag_top_bound=1 for sponge over 4 topmost layers
47!    iflag_top_bound=2 for sponge from top to ~1% of top layer pressure
48!    mode_top_bound=0: no relaxation
49!    mode_top_bound=1: u and v relax towards 0
50!    mode_top_bound=2: u and v relax towards their zonal mean
51!    mode_top_bound=3: u,v and pot. temp. relax towards their zonal mean
52!    tau_top_bound : inverse of charactericstic relaxation time scale at
53!                       the topmost layer (Hz)
54
55
56#include "comdissipn.h"
57#include "iniprint.h"
58
59c   Arguments:
60c   ----------
61
62      real,intent(inout) :: ucov(iip1,jjp1,llm) ! covariant zonal wind
63      real,intent(inout) :: vcov(iip1,jjm,llm) ! covariant meridional wind
64      real,intent(inout) :: teta(iip1,jjp1,llm) ! potential temperature
65      real,intent(in) :: masse(iip1,jjp1,llm) ! mass of atmosphere
66      real,intent(in) :: dt ! time step (s) of sponge model
67      real,intent(out) :: ducov(iip1,jjp1,llm) ! increment on ucov due to sponge
68
69c   Local:
70c   ------
71
72      REAL massebx(iip1,jjp1,llm),masseby(iip1,jjm,llm),zm
73      REAL uzon(jjp1,llm),vzon(jjm,llm),tzon(jjp1,llm)
74     
75      integer i
76      REAL,SAVE :: rdamp(llm) ! quenching coefficient
77      real,save :: lambda(llm) ! inverse or quenching time scale (Hz)
78
79      LOGICAL,SAVE :: first=.true.
80
81      INTEGER j,l
82     
83      if (first) then
84         if (iflag_top_bound.eq.1) then
85! sponge quenching over the topmost 4 atmospheric layers
86             lambda(:)=0.
87             lambda(llm)=tau_top_bound
88             lambda(llm-1)=tau_top_bound/2.
89             lambda(llm-2)=tau_top_bound/4.
90             lambda(llm-3)=tau_top_bound/8.
91         else if (iflag_top_bound.eq.2) then
92! sponge quenching over topmost layers down to pressures which are
93! higher than 100 times the topmost layer pressure
94             lambda(:)=tau_top_bound
95     s       *max(presnivs(llm)/presnivs(:)-0.01,0.)
96         endif
97
98! quenching coefficient rdamp(:)
99!         rdamp(:)=dt*lambda(:) ! Explicit Euler approx.
100         rdamp(:)=1.-exp(-lambda(:)*dt)
101
102         write(lunout,*)'TOP_BOUND mode',mode_top_bound
103         write(lunout,*)'Sponge layer coefficients'
104         write(lunout,*)'p (Pa)  z(km)  tau(s)   1./tau (Hz)'
105         do l=1,llm
106           if (rdamp(l).ne.0.) then
107             write(lunout,'(6(1pe12.4,1x))')
108     &        presnivs(l),log(preff/presnivs(l))*scaleheight,
109     &           1./lambda(l),lambda(l)
110           endif
111         enddo
112         first=.false.
113      endif ! of if (first)
114
115      CALL massbar(masse,massebx,masseby)
116
117      ! compute zonal average of vcov and u
118      if (mode_top_bound.ge.2) then
119       do l=1,llm
120        do j=1,jjm
121          vzon(j,l)=0.
122          zm=0.
123          do i=1,iim
124! NB: we can work using vcov zonal mean rather than v since the
125! cv coefficient (which relates the two) only varies with latitudes
126            vzon(j,l)=vzon(j,l)+vcov(i,j,l)*masseby(i,j,l)
127            zm=zm+masseby(i,j,l)
128          enddo
129          vzon(j,l)=vzon(j,l)/zm
130        enddo
131       enddo
132
133       do l=1,llm
134        do j=2,jjm ! excluding poles
135          uzon(j,l)=0.
136          zm=0.
137          do i=1,iim
138            uzon(j,l)=uzon(j,l)+massebx(i,j,l)*ucov(i,j,l)/cu(i,j)
139            zm=zm+massebx(i,j,l)
140          enddo
141          uzon(j,l)=uzon(j,l)/zm
142        enddo
143       enddo
144      else ! ucov and vcov will relax towards 0
145        vzon(:,:)=0.
146        uzon(:,:)=0.
147      endif ! of if (mode_top_bound.ge.2)
148
149      ! compute zonal average of potential temperature, if necessary
150      if (mode_top_bound.ge.3) then
151       do l=1,llm
152        do j=2,jjm ! excluding poles
153          zm=0.
154          tzon(j,l)=0.
155          do i=1,iim
156            tzon(j,l)=tzon(j,l)+teta(i,j,l)*masse(i,j,l)
157            zm=zm+masse(i,j,l)
158          enddo
159          tzon(j,l)=tzon(j,l)/zm
160        enddo
161       enddo
162      endif ! of if (mode_top_bound.ge.3)
163
164      if (mode_top_bound.ge.1) then
165       ! Apply sponge quenching on vcov:
166       do l=1,llm
167        do i=1,iip1
168          do j=1,jjm
169            vcov(i,j,l)=vcov(i,j,l)
170     &                  -rdamp(l)*(vcov(i,j,l)-vzon(j,l))
171          enddo
172        enddo
173       enddo
174
175       ! Apply sponge quenching on ucov:
176       do l=1,llm
177        do i=1,iip1
178          do j=2,jjm ! excluding poles
179            ducov(i,j,l)=-rdamp(l)*(ucov(i,j,l)-cu(i,j)*uzon(j,l))
180            ucov(i,j,l)=ucov(i,j,l)
181     &                  +ducov(i,j,l)
182          enddo
183        enddo
184       enddo
185      endif ! of if (mode_top_bound.ge.1)
186
187      if (mode_top_bound.ge.3) then
188       ! Apply sponge quenching on teta:
189       do l=1,llm
190        do i=1,iip1
191          do j=2,jjm ! excluding poles
192            teta(i,j,l)=teta(i,j,l)
193     &                  -rdamp(l)*(teta(i,j,l)-tzon(j,l))
194          enddo
195        enddo
196       enddo
197      endif ! of if (mode_top_bound.ge.3)
198   
199      END
Note: See TracBrowser for help on using the repository browser.