source: flexpart.git/src/cbl.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

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