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

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

sources for flexwrf v3.1

File size: 12.6 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!* FLEXPART is free software: you can redistribute it and/or modify    *
9!* it under the terms of the GNU General Public License as published by*
10!* the Free Software Foundation, either version 3 of the License, or   *
11!* (at your option) any later version.                                 *
12!*                                                                     *
13!* FLEXPART is distributed in the hope that it will be useful,         *
14!* but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15!* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16!* GNU General Public License for more details.                        *
17!*                                                                     *
18!* You should have received a copy of the GNU General Public License   *
19!* along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!***********************************************************************
21   
22    subroutine cbl(wp,zp,ust,wst,h,rhoa,rhograd,sigmaw,dsigmawdz,tlw,ptot,Q,phi,ath,bth,ol,flagrein)
23!                  i  i i  i   i  i    i     i       i         i   o   o   o   o  o   o   i     i/o
24!=============== well mixed formulation of CBL skewed vertical profiles following  LHH 1996 with profile of w3 from lHB 2000                    ========
25!=============== LHH formulation has been modified to account for variable density profiles and backward in time or forward in time simulations ========
26!=============== by  Massimo Cassiani ( mc ), NILU,  2012-2013, reference to Cassiani et al. 2013 (to be submitted...)                          ========
27!=======================================================================================================================================================
28!======================================================================================================
29! wp: particle velocity
30! zp: particle position
31! ust: velcotiy scale
32! wst: convective velcotiy scale
33! h: boundary layer top
34! rhoa: air density
35! rhograd: air densiy vertical gradient
36! sigmaw: turbulent flutuation of vertical velocity standard deviation
37! dsigmawdz: derivative of above
38! tlw: local lagrangina time scale
39! ptot: pdf value for the particle velocity in drift coefficient, see Cassiani et al. 2013, not used
40! Q: part of drift coefficient, not used
41! phi: part of drift coeffcient, not used
42! ath: drift coefficient, used
43! bth: diffusion coeffcient, sued
44! ol: Obukhov lenght
45! flagrein: set accordingly to conditon below if 1 then re-initialize particle velocity
46!======================================================================================================
47    use par_mod, only:pi
48    use com_mod, only:ldirect
49!   use ieee_arithmetic   
50    implicit none
51
52
53    real :: usurad2,usurad2p,C0,costluar4,eps
54    parameter  (usurad2=0.7071067812,usurad2p=0.3989422804,C0=2,costluar4=0.66667,eps=0.000001)
55
56    integer flagrein  !re-initlization flag for the particle velocity
57   
58    real :: wp,zp,ust,wst,h,dens,ddens,sigmaw,dsigmawdz,tlw,rhoa,rhograd
59    real ::fluarw,fluarw2
60    real ::w3,w2
61    real ::dw3,dw2
62    real ::wb,wa
63    real ::deltawa,deltawb
64    real ::wold,wold2,wold_z
65    real ::pa,pb,alfa
66    real ::Phi,Q,ptot 
67    real :: timedir
68    real ::cuberoot
69    real ::z0,ol,transition !added ol & transition with respect to cbl.f90 without transition
70   
71
72    real :: &
73    erf, &
74    aperfa, &
75    aperfb, &
76    ath, &
77    bth
78
79    real ::  &
80    pow, &
81    z, &   
82    skew, &
83    skew2, &
84    radw2, &
85    rluarw, &
86    xluarw, &
87    aluarw, &
88    bluarw, &
89    sigmawa, &
90    sigmawb, &
91    dskew, &
92    dradw2, &
93    dfluarw, &
94    drluarw, &
95    dxluarw, &
96    daluarw, &
97    dbluarw, &
98    dsigmawa, &
99    dsigmawb, &
100    dwa, &
101    dwb, &
102    sigmawa2, &
103    sigmawb2
104   
105   
106    dens=rhoa
107    ddens=rhograd
108             
109   
110    timedir=ldirect !direction of time forward (1) or backward(-1)
111    !========================= assign z ==============================
112    z=(zp/h)
113   
114    transition=1.
115    !if (ol.lt.-50) transition=((sin(((ol+100.)/100.)*pi))-1.)/2.
116    if (-h/ol.lt.15) transition=((sin((((-h/ol)+10.)/10.)*pi)))/2.+0.5  !transition fucntion to smoohtly
117    !========================= secondo moment of vertical velocity  =====================
118    !!!     w2=1.4*(z**1.5*(1.-z))**(2./3.)
119    !w2=(1.7*(z*(1.-0.7*z)*(1.-z))**(2./3.))*(wst**2)
120    w2=(sigmaw*sigmaw)
121    dw2=(2.*sigmaw*dsigmawdz)
122    !dw2=(1.7*(2./3.)*(z*(1.-0.7*z)*(1.-z))**(-1./3.)* &
123    !(((1.-0.7*z)*(1.-z))+z*(-0.7)*(1.-z)+z*(1.-0.7*z)*(-1.))) *(wst**2)*1/h
124
125    !=================== dissipation fo turbulent tke  =========================
126    !alfa=0.4 !(0.75-(0.5*z*z))**(3./2.) DISSIPAZIONE ADIMENSIONALE
127    alfa=2.*w2/(C0*tlw)
128
129    !========================================================================
130    wold=timedir*wp !time direction enter here for backward calculualtions
131    !wold_z=wp
132    ! =======================================================================
133    !------------------------------ momento terzo ============================
134    !!  w3=0.8*(w2**(3./2.))
135    !!  dw3=0.8*1.5*w2**0.5*dw2
136
137    ! dw3=((1.2*z*((1.-z)**(3./2.)))+eps)*(wst**3)
138    ! dw3=(1.2*(((1.-z)**(3./2.))+z*1.5*((1.-z)**(1./2.))*(-1.)))*(wst**3)
139    ! 3=(1.2*z*((1.-z)**(3./2.)))
140    ! w3=(1.2*(((1.-z)**(3./2.))+z*1.5*((1.-z)**(1./2.))*(-1.)))
141   
142                            !w3=((1.2*z*((1.-z)**(3./2.)))*1.5+eps)*(wst**3)  !§1.5 to increase skeweness see also initalize_cbl_vel.f90
143                            !dw3=(1.2*(((1.-z)**(3./2.))+z*1.5*((1.-z)**(1./2.))*(-1.)))*(wst**3)*(1./h)*1.5
144                           
145                            w3=((1.2*z*((1.-z)**(3./2.)))+eps)*(wst**3)*transition
146                            dw3=(1.2*(((1.-z)**(3./2.))+z*1.5*((1.-z)**(1./2.))*(-1.)))*(wst**3)*(1./h)*transition
147    !============================================================================0
148   
149   
150                            skew=w3/(w2**1.5)
151                            skew2=skew*skew
152                            dskew=(dw3*w2**(1.5)-w3*1.5*w2**0.5*dw2)/w2**3
153                            radw2=w2**0.5
154                            dradw2=0.5*w2**(-0.5)*dw2
155                            !costluar4=0.66667  ! costante da LHH
156                            fluarw=costluar4*(cuberoot(skew))   !
157                            fluarw2=fluarw*fluarw
158   
159   
160    if (skew.ne.0) then
161     
162                             dfluarw=costluar4*(1./3.)*cuberoot(skew**(-2.))*dskew
163       
164                            rluarw=(1.+fluarw2)**3.*skew2/((3.+fluarw2)**2.*fluarw2)  !-> r
165                            xluarw=(1.+fluarw2)**1.5*skew/((3.+fluarw2)*fluarw)       !
166       
167                            drluarw=( ((3.*(1.+fluarw2)**2*(2.*fluarw*dfluarw)*skew2)+ &
168                            (1.+fluarw2)**3*2.*skew*dskew) *(3.+fluarw2)**2.*fluarw2 - &
169                            (1.+fluarw2)**3*skew2* &
170                            ( (2.*(3.+fluarw2)*(2.*fluarw*dfluarw)*fluarw2) + &
171                            (3.+fluarw2)**2*2.*fluarw*dfluarw) )/ &
172                            (((3.+fluarw2)**2.*fluarw2)**2)
173                                     
174                            dxluarw=( ((1.5*(1.+fluarw2)**0.5*(2.*fluarw*dfluarw)*skew)+ &
175                            (1.+fluarw2)**1.5*dskew) *(3.+fluarw2)*fluarw - &
176                            (1.+fluarw2)**1.5*skew* &
177                            (3.*dfluarw+3*fluarw2*dfluarw) )/ &
178                            (((3.+fluarw2)*fluarw)**2)
179
180       
181    else       
182        dfluarw=0.
183        rluarw=0.
184        drluarw=0.
185        xluarw=0.
186        dxluarw=0.
187    end if
188
189
190
191                            aluarw=0.5*(1.-xluarw/(4.+rluarw)**0.5)
192                            bluarw=1.-aluarw
193
194                           daluarw=-0.5*(  (dxluarw*(4.+rluarw)**0.5) - &
195                            (0.5*xluarw*(4.+rluarw)**(-0.5)*drluarw) ) &
196                            /(4.+rluarw)
197                            dbluarw=-daluarw
198   
199
200                            sigmawa=radw2*(bluarw/(aluarw*(1.+fluarw2)))**0.5
201                            sigmawb=radw2*(aluarw/(bluarw*(1.+fluarw2)))**0.5
202
203
204                            dsigmawa=dradw2*(bluarw/(aluarw*(1.+fluarw2)))**0.5+ &
205                            radw2*( &
206                            (0.5*(bluarw/(aluarw*(1.+fluarw2)))**(-0.5))      * &
207                            ( &
208                            (dbluarw*(aluarw*(1.+fluarw2))- &
209                            bluarw*(daluarw*(1.+fluarw2)+aluarw*2.*fluarw*dfluarw)) &
210                            /((aluarw*(1.+fluarw2))**2) &
211                            ) &
212                            )
213
214                            dsigmawb=dradw2*(aluarw/(bluarw*(1.+fluarw2)))**0.5+ &
215                            radw2*( &
216                            (0.5*(aluarw/(bluarw*(1.+fluarw2)))**(-0.5))      * &
217                            ( &
218                            (daluarw*(bluarw*(1.+fluarw2))- &
219                            aluarw*(dbluarw*(1.+fluarw2)+bluarw*2.*fluarw*dfluarw)) &
220                            /((bluarw*(1.+fluarw2))**2) &
221                            ) &
222                            )
223   
224                            wa=(fluarw*sigmawa)
225                            wb=(fluarw*sigmawb)
226                            dwa=dfluarw*sigmawa+fluarw*dsigmawa
227                            dwb=dfluarw*sigmawb+fluarw*dsigmawb
228
229                            deltawa=wold-wa
230                            deltawb=wold+wb
231                            wold2=wold*wold
232                            sigmawa2=sigmawa*sigmawa
233                            sigmawb2=sigmawb*sigmawb
234
235                            pa=(usurad2p*(1./sigmawa))*(exp(-(0.5*((deltawa/sigmawa)**2.))))
236                            pb=(usurad2p*(1./sigmawb))*(exp(-(0.5*((deltawb/sigmawb)**2.))))
237                           
238                            if (abs(deltawa).gt.10.*sigmawa.and.abs(deltawb).gt.10.*sigmawb) flagrein=1  !added control flag for re-initialization of velocity
239!                           if (abs(deltawa).gt.6.*sigmawa.and.abs(deltawb).gt.6.*sigmawb) flagrein=1  !added control flag for re-initialization of velocity
240                                   
241                            ptot=dens*aluarw*pa+dens*bluarw*pb
242                               
243                            aperfa=deltawa*usurad2/sigmawa
244                            aperfb=deltawb*usurad2/sigmawb
245
246!       if ((ieee_is_nan(aperfa).or.ieee_is_nan(aperfb)).and.flagrein.eq.0) &
247!          print*,'PROBLEM',deltawa,deltawb,sigmawa,sigmawb,wp,zp,ust,wst,h,rhoa,rhograd,sigmaw,dsigmawdz,tlw,ptot,Q,phi,ath,bth,ol,flagrein
248                            Phi=-0.5* &
249                            (aluarw*dens*dwa+dens*wa*daluarw+aluarw*wa*ddens)*erf(aperfa) &
250                            +sigmawa*(aluarw*dens*dsigmawa*(wold2/sigmawa2+1.)+ &
251                            sigmawa*dens*daluarw+sigmawa*ddens*aluarw+ &
252                            aluarw*wold*dens/sigmawa2*(sigmawa*dwa-wa*dsigmawa))*pa &
253                            +0.5* &
254                            (bluarw*dens*dwb+wb*dens*dbluarw+wb*bluarw*ddens)*erf(aperfb) &
255                            +sigmawb*(bluarw*dens*dsigmawb*(wold2/sigmawb2+1.)+ &
256                            sigmawb*dens*dbluarw+sigmawb*ddens*bluarw+ &
257                            bluarw*wold*dens/sigmawb2*(-sigmawb*dwb+wb*dsigmawb))*pb
258
259
260                             Q=timedir*((aluarw*dens*deltawa/sigmawa2)*pa+ &
261                            (bluarw*dens*deltawb/sigmawb2)*pb)
262
263   
264
265
266                            ath=(1./ptot)*(-(C0/2.)*alfa*Q+phi)  !drift coefficient
267                            bth=sqrt(C0*alfa)                    !diffusion coefficient
268                           
269
270   
271
272    return
273
274
275    end
276                           
277
278
279
280
281    FUNCTION CUBEROOT (X) RESULT (Y)
282
283    IMPLICIT NONE
284
285    real, INTENT(IN) :: X
286    real:: Y
287
288    real, PARAMETER :: THIRD = 0.333333333
289
290
291    Y = SIGN((ABS(X))**THIRD, X)
292
293    RETURN
294
295    END FUNCTION CUBEROOT
296   
297   
298   
299
300    FUNCTION CUBEROOTD (X) RESULT (Y)
301
302    IMPLICIT NONE
303
304    DOUBLE PRECISION, INTENT(IN) :: X
305    DOUBLE PRECISION :: Y
306
307    DOUBLE PRECISION, PARAMETER :: THIRD = 0.33333333333333333333333333333333333333333333333333333333333333333333333333333333333D0
308
309
310    Y = SIGN((ABS(X))**THIRD, X)
311
312    RETURN
313
314    END FUNCTION CUBEROOTD
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG