source: branches/jerome/src_flexwrf_v3.1/tke_partition_my.f90 @ 16

Last change on this file since 16 was 16, checked in by jebri, 10 years ago

sources for flexwrf v3.1

File size: 9.1 KB
Line 
1!***********************************************************************
2!* Copyright 2012,2013                                                *
3!* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine,       *
4!* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,*
5!* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso,          *
6!* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann,   *
7!*                                                                     *
8!* This file is part of FLEXPART WRF                                   *
9!*                                                                     *
10!* FLEXPART is free software: you can redistribute it and/or modify    *
11!* it under the terms of the GNU General Public License as published by*
12!* the Free Software Foundation, either version 3 of the License, or   *
13!* (at your option) any later version.                                 *
14!*                                                                     *
15!* FLEXPART is distributed in the hope that it will be useful,         *
16!* but WITHOUT ANY WARRANTY; without even the implied warranty of      *
17!* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
18!* GNU General Public License for more details.                        *
19!*                                                                     *
20!* You should have received a copy of the GNU General Public License   *
21!* along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
22!***********************************************************************
23
24      subroutine tke_partition_my(z, &
25   ust,wst,ol,h,zeta,sigu,sigv,tlu,tlv,tlw, &
26   sigw,dsigwdz,dsigw2dz,uprof,vprof,tkeprof,pttprof,indz,indzp)
27
28!                      i
29!*******************************************************************************
30!                                                                              *
31!     Computation of \sigma_u,v,w,dsigwdz, and dsigw2dz based on TKE from WRF  *
32!                                                                              *
33!     Author: A. Stohl                                                         *
34!                                                                              *
35!     4 December 1997                                                          *
36!                                                                              *
37!*******************************************************************************
38!                                                                              *
39! Variables:                                                                   *
40! dsigwdz [1/s]     vertical gradient of sigw                                  *
41! ol [m]            Obukhov length                                             *
42! sigu, sigv, sigw  standard deviations of turbulent velocity fluctuations     *
43! tlu [s]           Lagrangian time scale for the along wind component.        *
44! tlv [s]           Lagrangian time scale for the cross wind component.        *
45! tlw [s]           Lagrangian time scale for the vertical wind component.     *
46! ust, ustar [m/s]  friction velocity                                          *
47! wst, wstar [m/s]  convective velocity scale                                  *
48!                                                                              *
49!*******************************************************************************
50! 12 JUNE 2007  compute sigu,sigv,sigw,dsigwdz, and dsigw2dz from TKE          *
51! fu2               fraction for u2
52! fv2               fraction for v2
53! fw2               fraction for w2
54!
55! 25 JUNE 2007 merged from hanna_tke and tke_partition.f
56! 11 Sep 2007 implement different formulations for growing and decaying turbulences
57 
58!      include 'includepar'
59!      include 'includecom'
60!      include 'includehanna'
61!      include 'includeinterpol'
62  use par_mod
63  use com_mod
64! use hanna_mod
65!  use interpol_mod
66
67      implicit none
68      real :: z,zzz,zzz1,dz,dz1,dz2,tke_z,dpttdz
69      real :: fu2,fv2,fw2,ftotal,ygu,ygv,ygm,ygh,yl,ylmax,ysm,ysh
70      real :: siguprof(2),sigvprof(2),sigwprof(2)
71      real :: ya1
72      real :: ya2
73      real :: yb1
74      real :: yb2
75      real :: yc1
76      real :: yse
77      real :: yap
78      real :: e1,e2,e3,e4,e5,rf,rfc,rf1,rf2,ch,cm
79      real :: er,smr,shr,dudz,dvdz,tke_mid
80      integer :: k,indz,indzp
81     real :: uprof(nzmax),vprof(nzmax),tkeprof(nzmax),pttprof(nzmax)
82    real :: ust,wst,ol,h,zeta,sigu,sigv,tlu,tlv,tlw
83    real :: sigw,dsigwdz,dsigw2dz
84
85!      data ya1,ya2,yb1,yb2,yc1,yse,yap/
86!     +   0.92,0.74,16.6,10.1,0.08,0.20,0.17/
87       ya1=0.92
88       ya2=0.74
89       yb1=16.6
90       yb2=10.1
91       yc1=0.08
92       yse=0.20
93       yap=0.17
94
95!- interpolate tke
96          dz=1./(height(indzp)-height(indz))
97          dz1=(z - height(indz))*dz
98          dz2=(height(indzp)-z)*dz
99 
100          tke_z=max(0.01,dz1*tkeprof(indzp)+dz2*tkeprof(indz))
101
102!- compute turbulent length scale
103!  yl -- turbulence length scale
104!  ylmax -- max yl
105          yl=0.01
106          ylmax=0.01
107          do k=1,nz-1
108            zzz=0.5*(height(k)+height(k+1))
109            zzz1=height(k+1)-height(k)
110            yl=yl+zzz*sqrt(0.5*(tkeprof(k)+tkeprof(k+1)))*zzz1
111            ylmax=ylmax+zzz1*sqrt(0.5*(tkeprof(k)+tkeprof(k+1)))
112          enddo
113            ylmax=0.1*yl/ylmax
114            yl=0.4*(z+0.1)/(1.0+0.4*(z+0.1)/ylmax)
115            if (pttprof(indzp) .gt. pttprof(indz)) then
116
117              dpttdz = (pttprof(indzp)-pttprof(indz))*dz
118
119              yl=min(yl,0.75*sqrt(2.0*tke_z/(9.8/300.0)/dpttdz))
120!            endif
121!- compute nondimensional vertical gradients
122          dudz=(uprof(indzp)-uprof(indz))*dz
123          dvdz=(vprof(indzp)-vprof(indz))*dz
124         ygu=yl/sqrt(2.0*tke_z)*dudz
125         ygv=yl/sqrt(2.0*tke_z)*dvdz
126         ygm=ygu*ygu+ygv*ygv+0.001       ! in case of zero
127         ygh=-yl*yl/2.0/tke_z*dpttdz*9.8/300.0
128
129!--  compute SM,SH
130!-   sm,sh nondimensional eddy diffusivities
131
132         ysm=1.0-3.0*ya2*(7.0*ya1+yb2)*ygh+ &
133              27.0*ya1*ya2*ya2*(4*ya1+yb2)*ygh*ygh+ &
134              6.0*ya1*ya1*(1.0-3.0*ya2*(yb2-3.0*ya2)*ygh)*ygm
135         ysm=ya1*(1-3*yc1-3*ya2*(yb2*(1-3*yc1)-12*ya1*yc1-3*ya2) &
136                 *ygh)/ysm
137
138         ysh=ya2*(1-6*ya1*ysm*ygm)/(1-3*ya2*ygh*(4*ya1+yb2))   
139
140
141! -- compute the equlibrium TKE, er
142        e1=yb1-6.0*ya1
143        e2=12*ya1+yb1+3.0*yb2
144        e3=yb1*(1.0-3.0*yc1)-6.0*ya1
145        e4=yb1*(1.0-3.0*yc1)+12*ya1+9.0*ya2
146        e5=yb1+3.0*ya1+3.0*yb2
147        ch=ya2*e2/yb1
148        cm=ya1*e4/ya2/e5
149        rf1=e3/e4
150        rf2=e1/e5
151        rfc=e1/e2
152        rf=-(ysh/ysm)*( ygh/ygm )
153!          write(*,*)'ysh,ysm',ysh,ysm
154!          write(*,*)'ygh,ygm,dpttdz, ptt',ygh,ygm,dpttdz,pttprof(indzp)
155!          if (rf .le. 0) write(*,*)'aaaaaaaaaaaaaaa'       
156        shr = 0.0
157        if (rf .ge. rfc) then
158            shr=ch*1.0e-4            ! Ri > critical value, then =
159        elseif ((1.0-rf) .ne. 0.0) then
160            shr=ch*(rfc-rf)/(1.0-rf)
161        endif
162
163        smr = 0.0
164        if (rf .ge. rf1) then
165           smr=cm*shr*1.0e-4
166        elseif((rf2-rf) .ne. 0.0 .and. rf .lt. rf2) then
167           smr=cm*(rf1-rf)/(rf2-rf)*shr
168        endif
169       
170        er=0.5*yb1*yl*yl*(smr*dudz*dudz*dvdz*dvdz- &
171                          shr*dpttdz*9.8/300.0)
172        er=max(er,0.01)
173        if (tke_z .ge. er ) then                    ! decaying turbuelnce 
174          tke_mid=tke_z
175        else
176          tke_mid=er
177          ysm=sqrt(tke_z/er)*smr
178          ysh=sqrt(tke_z/er)*shr
179        endif
180 
181!- fractions of u2,v2,w2
182!        print*,'in tke',tke_mid,ya1,ysm,ygm,ysh,ygh,ygu
183         fu2=0.333-2*ya1*(ysm*ygm+ysh*ygh)+6*ya1*ysm*ygu*ygu
184         fv2=0.333-2*ya1*(ysm*ygm+ysh*ygh)+6*ya1*ysm*ygv*ygv
185         fw2=0.333-2*ya1*(ysm*ygm+ysh*ygh)+6*ya1*ysh*ygh
186!        fw2=0.333-2*ya1*(ysm*ygm+ysh*ygh)+6*ya1*ysm*ygh*ygh
187!        fw2=fw2*0.5
188         if(fu2*fv2*fw2 .le. 0.0) then
189           fu2=0.333
190           fv2=0.333
191           fw2=0.333
192         endif
193
194        else ! test on pttprof
195
196           fu2=0.333
197           fv2=0.333
198           fw2=0.333
199          tke_mid=tke_z
200        endif
201          !ftotal=(fu2+fv2+fw2)
202          !fu2=fu2/ftotal
203          !fv2=fv2/ftotal
204          !fw2=fw2/ftotal
205
206! JB: make the partitionning equal
207!          fu2=0.333
208!          fv2=0.333
209!          fw2=0.333
210!
211         sigu=sqrt(2.0*tke_mid*fu2)
212         sigv=sqrt(2.0*tke_mid*fv2)
213         sigw=sqrt(2.0*tke_mid*fw2)
214
215         tlu=2.0*yl/sigu
216         tlv=2.0*yl/sigv
217         tlw=2.0*yl/sigw
218
219         tlu=max(10.,tlu)
220         tlv=max(10.,tlv)
221         tlw=max(30.,tlw)
222
223!- dsigw/dz,dsigw2/dz ; assuming fraction is not changed
224
225         dsigwdz=(sqrt(2.0*tkeprof(indzp)*fw2)- &
226                  sqrt(2.0*tkeprof(indz)*fw2))*dz
227         dsigwdz=max(1.0e-10,dsigwdz)
228
229         dsigw2dz=(2.0*tkeprof(indzp)*fw2-2.0*tkeprof(indz)*fw2)*dz
230         dsigw2dz=max(1.0e-10,dsigw2dz)
231
232!         write(*,*)'tke=',tke_z, 'er=',er
233!         write(*,*)'rfc,rf1,rf2,rf=',rfc,rf1,rf2,rf
234!          write(*,*)'smr=',smr,'ysm=',ysm
235!         write(*,*)'shr=',shr,'ysh=',ysh
236!         write(*,*)'fu2,fv2,fw2=',fu2,fv2,fw2
237!         write(*,*)'ftotal=',ftotal
238!         if (fu2*fv2*fw2 .le. 0) stop
239end subroutine tke_partition_my
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG