source: flexpart.git/src/cbl.f90 @ 4c64400

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 4c64400 was 5f9d14a, checked in by Espen Sollum ATMOS <eso@…>, 9 years ago

Updated wet depo scheme

  • Property mode set to 100644
File size: 9.2 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
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
22subroutine 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  i    o
24     
25    use par_mod, only:pi
26    use com_mod, only:ldirect
27   
28    implicit none
29!=======================================================================================================================================================
30!=============== CBL skewed vertical profiles and formulation of LHH 1996 with profile of w3 from LHB 2000                                      ========
31!=============== LHH formulation has been modified to account for variable density profiles and backward in time or forward in time simulations ========
32!=============== see Cassiani et al. BLM 2014 doi  for explanations and references                                                              ========
33!=======================================================================================================================================================
34
35    real :: usurad2,usurad2p,C0,costluar4,eps
36    parameter  (usurad2=0.7071067812,usurad2p=0.3989422804,C0=3,costluar4=0.66667,eps=0.000001)
37
38       
39    integer :: flagrein
40    real :: wp,zp,ust,wst,h,dens,ddens,sigmaw,dsigmawdz,tlw,rhoa,rhograd
41    real ::fluarw,fluarw2
42    real ::w3,w2
43    real ::dw3,dw2
44    real ::wb,wa
45    real ::deltawa,deltawb
46    real ::wold,wold2,wold_z
47    real ::pa,pb,alfa
48    real ::Phi,Q,ptot 
49    real :: timedir
50    real ::cuberoot
51    real ::z0,ol,transition   
52   
53
54    real :: &
55    erf, &
56    aperfa, &
57    aperfb, &
58    ath, &
59    bth
60
61    real ::  &
62    pow, &
63    z, &   
64    skew, &
65    skew2, &
66    radw2, &
67    rluarw, &
68    xluarw, &
69    aluarw, &
70    bluarw, &
71    sigmawa, &
72    sigmawb, &
73    dskew, &
74    dradw2, &
75    dfluarw, &
76    drluarw, &
77    dxluarw, &
78    daluarw, &
79    dbluarw, &
80    dsigmawa, &
81    dsigmawb, &
82    dwa, &
83    dwb, &
84    sigmawa2, &
85    sigmawb2
86   
87   
88    dens=rhoa      !put to 1 for test constant density simulation
89    ddens=rhograd  !put to 0 for test constant density simulation
90             
91   
92    timedir=ldirect !ldirect contains direction of time forward (1) or backward(-1)
93    !========================= assegnazione z ===========================================================
94    z=(zp/h)
95   
96    !================== stability transition function see Cassiani et al(2015) BLM ======================
97    transition=1.
98    !if (ol.lt.-50) transition=((sin(((ol+100.)/100.)*pi))-1.)/2.
99    if (-h/ol.lt.15) transition=((sin((((-h/ol)+10.)/10.)*pi)))/2.+0.5
100
101    !========================= momento secondo ==========================================================
102   
103    w2=(sigmaw*sigmaw)
104    dw2=(2.*sigmaw*dsigmawdz)
105   
106
107    !=================== dissipazione =======================================
108   
109    alfa=2.*w2/(C0*tlw)
110
111    !========================================================================
112   
113    wold=timedir*wp
114   
115    !=========================================================================
116    !------------------------------ momento terzo ============================
117                               
118    w3=((1.2*z*((1.-z)**(3./2.)))+eps)*(wst**3)*transition
119    dw3=(1.2*(((1.-z)**(3./2.))+z*1.5*((1.-z)**(1./2.))*(-1.)))*(wst**3)*(1./h)*transition
120   
121    !===========================================================================0
122   
123    skew=w3/(w2**1.5)
124    skew2=skew*skew
125    dskew=(dw3*w2**(1.5)-w3*1.5*w2**0.5*dw2)/w2**3
126    radw2=w2**0.5
127    dradw2=0.5*w2**(-0.5)*dw2
128    !costluar4=0.66667  ! costante da LHH
129    fluarw=costluar4*(cuberoot(skew))   !skew**(1./3.)
130    fluarw2=fluarw*fluarw
131   
132    if (skew.ne.0) then
133     
134        dfluarw=costluar4*(1./3.)*cuberoot(skew**(-2.))*dskew
135   
136        rluarw=(1.+fluarw2)**3.*skew2/((3.+fluarw2)**2.*fluarw2)  !-> r
137        xluarw=(1.+fluarw2)**1.5*skew/((3.+fluarw2)*fluarw)    !----> r^1/2
138       
139        drluarw=( ((3.*(1.+fluarw2)**2*(2.*fluarw*dfluarw)*skew2)+ &
140            (1.+fluarw2)**3*2.*skew*dskew) *(3.+fluarw2)**2.*fluarw2 - &
141            (1.+fluarw2)**3*skew2* &
142            ( (2.*(3.+fluarw2)*(2.*fluarw*dfluarw)*fluarw2) + &
143            (3.+fluarw2)**2*2.*fluarw*dfluarw) )/ &
144            (((3.+fluarw2)**2.*fluarw2)**2)
145                                     
146        dxluarw=( ((1.5*(1.+fluarw2)**0.5*(2.*fluarw*dfluarw)*skew)+ &
147            (1.+fluarw2)**1.5*dskew) *(3.+fluarw2)*fluarw - &
148            (1.+fluarw2)**1.5*skew* &
149            (3.*dfluarw+3*fluarw2*dfluarw) )/ &
150            (((3.+fluarw2)*fluarw)**2)
151   
152    else       
153        dfluarw=0.
154        rluarw=0.
155        drluarw=0.
156        xluarw=0.
157        dxluarw=0.
158    end if
159
160
161
162   aluarw=0.5*(1.-xluarw/(4.+rluarw)**0.5)
163   bluarw=1.-aluarw
164
165   daluarw=-0.5*(  (dxluarw*(4.+rluarw)**0.5) - &
166           (0.5*xluarw*(4.+rluarw)**(-0.5)*drluarw) ) &
167           /(4.+rluarw)
168   dbluarw=-daluarw
169   
170   sigmawa=radw2*(bluarw/(aluarw*(1.+fluarw2)))**0.5
171   sigmawb=radw2*(aluarw/(bluarw*(1.+fluarw2)))**0.5
172
173   dsigmawa=dradw2*(bluarw/(aluarw*(1.+fluarw2)))**0.5+ &
174            radw2*( &
175            (0.5*(bluarw/(aluarw*(1.+fluarw2)))**(-0.5))      * &
176            ( &
177            (dbluarw*(aluarw*(1.+fluarw2))- &
178            bluarw*(daluarw*(1.+fluarw2)+aluarw*2.*fluarw*dfluarw)) &
179            /((aluarw*(1.+fluarw2))**2) &
180            ) &
181            )
182  dsigmawb=dradw2*(aluarw/(bluarw*(1.+fluarw2)))**0.5+ &
183           radw2*( &
184           (0.5*(aluarw/(bluarw*(1.+fluarw2)))**(-0.5))      * &
185           ( &
186           (daluarw*(bluarw*(1.+fluarw2))- &
187           aluarw*(dbluarw*(1.+fluarw2)+bluarw*2.*fluarw*dfluarw)) &
188           /((bluarw*(1.+fluarw2))**2) &
189           ) &
190           )
191   
192           wa=(fluarw*sigmawa)
193           wb=(fluarw*sigmawb)
194           dwa=dfluarw*sigmawa+fluarw*dsigmawa
195           dwb=dfluarw*sigmawb+fluarw*dsigmawb
196
197           deltawa=wold-wa
198           deltawb=wold+wb
199           wold2=wold*wold
200           sigmawa2=sigmawa*sigmawa
201           sigmawb2=sigmawb*sigmawb
202                           
203           if (abs(deltawa).gt.6.*sigmawa.and.abs(deltawb).gt.6.*sigmawb) flagrein=1
204
205           pa=(usurad2p*(1./sigmawa))*(exp(-(0.5*((deltawa/sigmawa)**2.))))
206           pb=(usurad2p*(1./sigmawb))*(exp(-(0.5*((deltawb/sigmawb)**2.))))
207                                   
208           ptot=dens*aluarw*pa+dens*bluarw*pb
209                               
210           aperfa=deltawa*usurad2/sigmawa
211           aperfb=deltawb*usurad2/sigmawb
212
213           Phi=-0.5* &
214               (aluarw*dens*dwa+dens*wa*daluarw+aluarw*wa*ddens)*erf(aperfa) &
215               +sigmawa*(aluarw*dens*dsigmawa*(wold2/sigmawa2+1.)+ &
216               sigmawa*dens*daluarw+sigmawa*ddens*aluarw+ &
217               aluarw*wold*dens/sigmawa2*(sigmawa*dwa-wa*dsigmawa))*pa &
218               +0.5* &
219               (bluarw*dens*dwb+wb*dens*dbluarw+wb*bluarw*ddens)*erf(aperfb) &
220               +sigmawb*(bluarw*dens*dsigmawb*(wold2/sigmawb2+1.)+ &
221               sigmawb*dens*dbluarw+sigmawb*ddens*bluarw+ &
222               bluarw*wold*dens/sigmawb2*(-sigmawb*dwb+wb*dsigmawb))*pb
223
224               Q=timedir*((aluarw*dens*deltawa/sigmawa2)*pa+ &
225               (bluarw*dens*deltawb/sigmawb2)*pb)
226
227               ath=(1./ptot)*(-(C0/2.)*alfa*Q+phi)
228               bth=sqrt(C0*alfa)
229              !bth=sngl(sigmaw*sqrt(2.*tlw))
230
231    return
232    end
233                           
234
235
236
237
238    FUNCTION CUBEROOT (X) RESULT (Y)
239
240    IMPLICIT NONE
241
242    real, INTENT(IN) :: X
243    real:: Y
244
245    real, PARAMETER :: THIRD = 0.333333333
246
247
248    Y = SIGN((ABS(X))**THIRD, X)
249
250    RETURN
251
252    END FUNCTION CUBEROOT
253   
254   
255   
256
257    FUNCTION CUBEROOTD (X) RESULT (Y)
258
259    IMPLICIT NONE
260
261    DOUBLE PRECISION, INTENT(IN) :: X
262    DOUBLE PRECISION :: Y
263
264    DOUBLE PRECISION, PARAMETER :: THIRD = 0.33333333333333333333333333333333333333333333333333333333333333333333333333333333333D0
265
266
267    Y = SIGN((ABS(X))**THIRD, X)
268
269    RETURN
270
271    END FUNCTION CUBEROOTD
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG