source: trunk/LMDZ.MARS/libf/phymars/pbl_parameters.F @ 662

Last change on this file since 662 was 662, checked in by acolaitis, 13 years ago

minor modif to pbl_param. Now that we compute wstar taking into account diffusion, we dont need to recompute it in pbl_param.

File size: 16.8 KB
Line 
1      SUBROUTINE pbl_parameters(ngrid,nlay,ps,pplay,pz0,
2     & pg,zzlay,zzlev,pu,pv,wstar_in,hfmax,zmax,pts,ph,z_out,n_out,
3     & T_out,u_out,ustar,tstar,L_mo,vhf,vvv)
4      IMPLICIT NONE
5!=======================================================================
6!
7!   Anlysis of the PBL from input temperature, wind field and thermals outputs.
8!
9!   ------- 
10!
11!   Author: Arnaud Colaitis 09/01/12
12!   -------
13!
14!   Arguments:
15!   ----------
16!
17!   inputs:
18!   ------
19!     ngrid            size of the horizontal grid
20!     nlay             size of the vertical grid
21!     pz0(ngrid)       surface roughness length
22!     pg               gravity (m s -2)
23!     zzlay(ngrid,nlay)   height of mid-layers
24!     zzlev(ngrid,nlay+1)   height of layers interfaces
25!     pu(ngrid,nlay)   u component of the wind
26!     pv(ngrid,nlay)   v component of the wind
27!     wstar_in(ngrid)  free convection velocity in PBL
28!     hfmax(ngrid)     maximum vertical turbulent heat flux in thermals
29!     zmax(ngrid)      height reached by the thermals (pbl height)
30!     pts(ngrid)       surface temperature
31!     ph(ngrid,nlay)   potential temperature T*(p/ps)^kappa
32!     z_out(n_out)     heights of interpolation
33!     n_out            number of points for interpolation
34!
35!   outputs:
36!   ------
37!
38!     Teta_out(ngrid,n_out)  interpolated teta
39!     u_out(ngrid,n_out)     interpolated u
40!     ustar(ngrid)     friction velocity
41!     tstar(ngrid)     friction temperature
42!     L_mo(ngrid)      monin_obukhov length
43!
44!
45!=======================================================================
46!
47!-----------------------------------------------------------------------
48!   Declarations:
49!   -------------
50
51#include "comcstfi.h"
52#include "callkeys.h"
53
54!   Arguments:
55!   ----------
56
57      INTEGER, INTENT(IN) :: ngrid,nlay,n_out
58      REAL, INTENT(IN) :: pz0(ngrid),ps(ngrid),pplay(ngrid,nlay)
59      REAL, INTENT(IN) :: pg,zzlay(ngrid,nlay),zzlev(ngrid,nlay)
60      REAL, INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay)
61      REAL, INTENT(IN) :: wstar_in(ngrid),hfmax(ngrid),zmax(ngrid)
62      REAL, INTENT(IN) :: pts(ngrid),ph(ngrid,nlay)
63      REAL, INTENT(IN) :: z_out(n_out)
64
65!    Outputs:
66!    --------
67
68      REAL, INTENT(OUT) :: T_out(ngrid,n_out),u_out(ngrid,n_out)
69      REAL Teta_out(ngrid,n_out)
70      REAL, INTENT(OUT) :: ustar(ngrid), tstar(ngrid)
71      REAL, INTENT(OUT) :: L_mo(ngrid)
72
73!   Local:
74!   ------
75
76      INTEGER ig,k,n
77      REAL karman,nu
78      DATA karman,nu/.41,0.001/
79      SAVE karman,nu
80
81!    Local(2):
82!    ---------
83     
84      REAL zout
85      REAL rib(ngrid)  ! Bulk Richardson number
86      REAL fm(ngrid) ! stability function for momentum
87      REAL fh(ngrid) ! stability function for heat
88      REAL z1z0,z1z0t ! ratios z1/z0 and z1/z0T
89          ! phim = 1+betam*zeta   or   (1-bm*zeta)**am
90          ! phih = alphah + betah*zeta    or   alphah(1.-bh*zeta)**ah
91      REAL betam, betah, alphah, bm, bh, lambda
92          ! ah and am are assumed to be -0.25 and -0.5 respectively
93      REAL cdn(ngrid),chn(ngrid)  ! neutral momentum and heat drag coefficient
94      REAL pz0t        ! initial thermal roughness length. (local)
95      REAL ric         ! critical richardson number
96      REAL reynolds(ngrid)    ! reynolds number for UBL
97      REAL prandtl(ngrid)     ! prandtl number for UBL
98      REAL pz0tcomp(ngrid)     ! computed z0t
99      REAL ite
100      REAL residual,zcd0,z1
101      REAL pcdv(ngrid),pcdh(ngrid)
102      REAL zu2(ngrid)                  ! Large-scale wind at first layer
103      REAL pbl_teta(ngrid)             ! mixed-layer potential temperature
104      INTEGER pbl_height_index(ngrid)  ! index of nearest vertical grid point for zmax
105      REAL dteta(ngrid,nlay),x(ngrid)  ! potential temperature gradient and z/zi
106      REAL dvhf(ngrid),dvvv(ngrid)     ! dimensionless vertical heat flux and
107                                       ! dimensionless vertical velocity variance
108      REAL vhf(ngrid),vvv(ngrid)       ! vertical heat flux and vertical velocity variance
109      INTEGER ii(1)
110! temporary
111      INTEGER dimout
112
113!------------------------------------------------------------------------
114!------------------------------------------------------------------------
115! PART I : RICHARDSON/REYNOLDS/THERMAL_ROUGHNESS/STABILITY_COEFFICIENTS
116!------------------------------------------------------------------------
117!------------------------------------------------------------------------
118
119      DO n=1,n_out
120
121c Initialisation :
122
123      L_mo(:)=0.
124      ustar(:)=0.
125      tstar(:)=0.
126      zout=z_out(n)
127      reynolds(:)=0.
128      pz0t = 0.
129      pz0tcomp(:) = 0.1*pz0(:)
130      rib(:)=0.
131      pcdv(:)=0.
132      pcdh(:)=0.
133
134      ! this formulation assumes alphah=1., implying betah=betam
135      ! We use Dyer et al. parameters, as they cover a broad range of Richardson numbers :
136
137      bm=16.            !UBL
138      bh=16.            !UBL
139      alphah=1.
140      betam=5.         !SBL
141      betah=5.         !SBL
142      lambda=(sqrt(bh/bm))/alphah
143      ric=betah/(betam**2)
144      DO ig=1,ngrid
145       ite=0.
146       residual=abs(pz0tcomp(ig)-pz0t)
147!       zu2(ig)=MAX(pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1)
148!     &                                ,(0.3*wstar_in(ig))**2)
149       zu2(ig)=pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1)
150     &     + (log(1.+0.7*wstar_in(ig) + 2.3*wstar_in(ig)**2))**2
151
152       DO WHILE((residual .gt. 0.01*pz0(ig)) .and.  (ite .lt. 10.))
153
154         pz0t=pz0tcomp(ig)
155         IF (zu2(ig) .ne. 0.) THEN
156            ! Richardson number formulation proposed by D.E. England et al. (1995)
157          rib(ig) = (pg/pts(ig))
158     &        *sqrt(zzlay(ig,1)*pz0(ig))
159     &        *(((log(zzlay(ig,1)/pz0(ig)))**2)/(log(zzlay(ig,1)/pz0t)))
160     &        *(ph(ig,1)-pts(ig))/zu2(ig)
161         ELSE
162            print*,'warning, infinite Richardson at surface'
163            print*,pu(ig,1),pv(ig,1)
164            rib(ig) = ric
165         ENDIF
166
167         z1z0=zzlay(ig,1)/pz0(ig)
168         z1z0t=zzlay(ig,1)/pz0t
169
170         cdn(ig)=karman/log(z1z0)
171         cdn(ig)=cdn(ig)*cdn(ig)
172         chn(ig)=cdn(ig)*log(z1z0)/log(z1z0t)
173
174         ! STABLE BOUNDARY LAYER :
175         IF (rib(ig) .gt. 0.) THEN
176            ! From D.E. England et al. (95)
177            prandtl(ig)=1.
178            if(rib(ig) .lt. ric) then
179               ! Assuming alphah=1. and bh=bm for stable conditions :
180               fm(ig)=((ric-rib(ig))/ric)**2
181               fh(ig)=((ric-rib(ig))/ric)**2
182            else
183               ! For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface
184               fm(ig)=0.
185               fh(ig)=0.
186            endif
187         ! UNSTABLE BOUNDARY LAYER :
188         ELSE
189            ! From D.E. England et al. (95)
190            fm(ig)=sqrt(1.-lambda*bm*rib(ig))
191            fh(ig)=(1./alphah)*((1.-lambda*bh*rib(ig))**0.5)*
192     &                     (1.-lambda*bm*rib(ig))**0.25
193            prandtl(ig)=alphah*((1.-lambda*bm*rib(ig))**0.25)/
194     &             ((1.-lambda*bh*rib(ig))**0.5)
195         ENDIF
196 
197        reynolds(ig)=karman*sqrt(fm(ig))
198     &              *sqrt(zu2(ig))
199     &              *pz0(ig)/(log(z1z0)*nu)
200        pz0tcomp(ig)=pz0(ig)*exp(-karman*7.3*
201     &              (reynolds(ig)**0.25)*(prandtl(ig)**0.5))
202        residual = abs(pz0t-pz0tcomp(ig))
203        ite = ite+1
204
205       ENDDO  ! of while
206       pz0t=0.
207
208! Drag computation:
209
210         pcdv(ig)=cdn(ig)*fm(ig)
211         pcdh(ig)=chn(ig)*fh(ig)
212       
213      ENDDO ! of ngrid
214
215!------------------------------------------------------------------------
216!------------------------------------------------------------------------
217! PART II : USTAR/TSTAR/U_OUT/TETA_OUT COMPUTATIONS
218!------------------------------------------------------------------------
219!------------------------------------------------------------------------
220
221
222! u* theta* computation
223! and Monin Obukhov length:
224
225      DO ig=1,ngrid
226         IF (rib(ig) .ge. ric) THEN
227           ustar(ig)=0.
228           tstar(ig)=0.
229           L_mo(ig)=0.
230         ELSE
231           ustar(ig)=sqrt(pcdv(ig))
232     &        *sqrt(pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1))
233           tstar(ig)=-pcdh(ig)*(pts(ig)-ph(ig,1))
234     &        /sqrt(pcdv(ig))
235           L_mo(ig)=pts(ig)*(ustar(ig)**2)/(pg*karman*tstar(ig))  !as defined here, L is positive for SBL, negative for UBL
236         ENDIF
237      ENDDO
238
239! Monin Obukhov length:
240
241!      DO ig=1,ngrid
242!         IF (rib(ig) .ge. ric) THEN
243!            L_mo(ig)=0.
244!         ELSE
245!            L_mo(ig)=pts(ig)*(ustar(ig)**2)/(pg*karman*tstar(ig))  !as defined here, L is positive for SBL, negative for UBL
246!         ENDIF
247!      ENDDO
248
249! Interpolation:
250
251      DO ig=1,ngrid
252        IF(zout .lt. pz0tcomp(ig)) THEN
253           u_out(ig,n)=0.
254           Teta_out(ig,n)=pts(ig)
255        ELSEIF (L_mo(ig) .gt. 0.) THEN
256!           u_out(ig,n)=(ustar(ig)/karman)*log(zout/pz0(ig)) +
257!     &            5.*(ustar(ig)/(karman*L_mo(ig)))*(zout-pz0(ig))
258!           Teta_out(ig,n)=pts(ig)+(tstar(ig)/(prandtl(ig)*karman))
259!     &                            *log(zout/pz0tcomp(ig)) +
260!     &            5.*(tstar(ig)/(prandtl(ig)*karman*L_mo(ig)))
261!     &                            *(zout-pz0tcomp(ig))
262          IF ((zout/L_mo(ig).lt.ric).and.(pz0(ig).lt.ric)) THEN
263!     &  ((zout/L_mo(ig).gt.ric).and.(pz0(ig).gt.ric))  ) THEN
264
265            u_out(ig,n)=(ustar(ig)/karman)*(
266     &  log((ric-pz0(ig)/L_mo(ig))/(ric-zout/L_mo(ig)))
267     & +log(zout/pz0(ig))
268     & )
269          ELSEIF ((zout/L_mo(ig).gt.ric).and.(pz0(ig).gt.ric)) THEN
270
271            u_out(ig,n)=(ustar(ig)/karman)*(
272     &  log((zout-ric*L_mo(ig))/(pz0(ig)-ric*L_mo(ig)))
273     & +log(pz0(ig)/zout)
274     & )
275
276          ELSEIF ((zout/L_mo(ig).gt.ric).and.(pz0(ig).lt.ric)) THEN
277          !integral of the stability function does not converge
278
279
280           u_out(ig,n)=sqrt(pu(ig,1)**2 + pv(ig,1)**2)
281
282
283          ENDIF
284          IF((zout/L_mo(ig).lt.ric).and.(pz0tcomp(ig).lt.ric)) THEN
285!     &  ((zout/L_mo(ig).gt.ric).and.(pz0tcomp(ig).gt.ric))  ) THEN
286
287           Teta_out(ig,n)=pts(ig)+(tstar(ig)/(prandtl(ig)*karman))*(
288     &  log((ric-pz0tcomp(ig)/L_mo(ig))/(ric-zout/L_mo(ig)))
289     & +log(zout/pz0tcomp(ig))
290     & )
291          ELSEIF ((zout/L_mo(ig).gt.ric).and.(pz0tcomp(ig).gt.ric)) THEN
292
293           Teta_out(ig,n)=pts(ig)+(tstar(ig)/(prandtl(ig)*karman))*(
294     &  log((zout-ric*L_mo(ig))/(pz0tcomp(ig)-ric*L_mo(ig)))
295     & +log(pz0tcomp(ig)/zout)
296     & )
297
298          ELSEIF ((zout/L_mo(ig).gt.ric).and.(pz0tcomp(ig).lt.ric)) THEN
299          !integral of the stability function does not converge
300
301           Teta_out(ig,n)=ph(ig,1)
302
303          ENDIF
304
305        ELSEIF (L_mo(ig) .lt. 0.) THEN
306
307          IF(L_mo(ig) .gt. -1000.) THEN
308           
309           u_out(ig,n)=(ustar(ig)/karman)*(
310     &  2.*atan((1.-16.*zout/L_mo(ig))**0.25)
311     & -2.*atan((1.-16.*pz0(ig)/L_mo(ig))**0.25)
312     & -2.*log(1.+(1.-16.*zout/L_mo(ig))**0.25)
313     & +2.*log(1.+(1.-16.*pz0(ig)/L_mo(ig))**0.25)
314     & -   log(1.+sqrt(1.-16.*zout/L_mo(ig)))
315     & +   log(1.+sqrt(1.-16.*pz0(ig)/L_mo(ig)))
316     & +   log(zout/pz0(ig))
317     &                                  )
318
319           Teta_out(ig,n)=MAX(pts(ig)+(tstar(ig)/(prandtl(ig)*karman))*(
320     &  2.*log(1.+sqrt(1.-16.*pz0tcomp(ig)/L_mo(ig)))
321     & -2.*log(1.+sqrt(1.-16.*zout/L_mo(ig)))
322     & +   log(zout/pz0tcomp(ig))
323     &                             ),ph(ig,1))
324
325          ELSE
326
327! We have to treat the case where L is very large and negative,
328! corresponding to a very nearly stable atmosphere (but not quite !)
329! i.e. teta* <0 and almost zero. We choose the cutoff
330! corresponding to L_mo < -1000 and use a 3rd order taylor expansion. The difference
331! between the two expression for z/L = -1/1000 is -1.54324*10^-9
332! (we do that to avoid using r*4 precision, otherwise, we get -inf values)         
333
334           u_out(ig,n)=(ustar(ig)/karman)*(
335     &     (4./L_mo(ig))*(zout-pz0(ig))
336     &   + (20./(L_mo(ig))**2)*(zout**2-pz0(ig)**2)
337     &   + (160./(L_mo(ig))**3)*(zout**3-pz0(ig)**3)
338     &   + log(zout/pz0(ig))
339     &                                   )
340
341           Teta_out(ig,n)=MAX(pts(ig)+(tstar(ig)/(prandtl(ig)*karman))*(
342     &     (8./L_mo(ig))*(zout-pz0tcomp(ig))
343     &   + (48./(L_mo(ig))**2)*(zout**2-pz0tcomp(ig)**2)
344     &   + (1280./(3.*(L_mo(ig))**3))*(zout**3-pz0tcomp(ig)**3)
345     &   + log(zout/pz0tcomp(ig))
346     &                             ),ph(ig,1))
347
348          ENDIF
349        ELSE
350           u_out(ig,n)=sqrt(pu(ig,1)**2 + pv(ig,1)**2)
351           Teta_out(ig,n)=ph(ig,1)
352        ENDIF
353        IF(zout .lt. pz0(ig)) THEN
354           u_out(ig,n)=0.
355        ENDIF
356
357! Final check
358        IF (L_mo(ig) .gt. 0) THEN
359           IF (Teta_out(ig,n) .gt. ph(ig,1)) THEN
360              Teta_out(ig,n)=ph(ig,1)
361           ELSEIF (Teta_out(ig,n) .lt. pts(ig)) THEN
362              Teta_out(ig,n)=pts(ig)
363           ENDIF
364        ELSEIF (L_mo(ig) .lt. 0) THEN
365           IF (Teta_out(ig,n) .lt. ph(ig,1)) THEN
366              Teta_out(ig,n)=ph(ig,1)
367           ELSEIF (Teta_out(ig,n) .gt. pts(ig)) THEN
368              Teta_out(ig,n)=pts(ig)
369           ENDIF
370        ENDIF
371
372        IF (u_out(ig,n) .gt. sqrt(pu(ig,1)**2 + pv(ig,1)**2)) THEN
373           u_out(ig,n)=sqrt(pu(ig,1)**2 + pv(ig,1)**2)
374        ENDIF
375
376      ENDDO
377
378! when using convective adjustment without thermals, a vertical potential temperature
379! profile is assumed up to the thermal roughness length. Hence, theoretically, theta
380! interpolated at any height in the surface layer is theta at the first level.
381
382      IF ((.not.calltherm).and.(calladj)) THEN
383       Teta_out(:,n)=ph(:,1)
384       u_out(:,n)=(sqrt(cdn(:))*sqrt(pu(:,1)*pu(:,1)+pv(:,1)*pv(:,1))
385     &                                /karman)*log(zout/pz0(:))
386      ENDIF
387              T_out(:,n) = Teta_out(:,n)*(exp(
388     &   (zout/zzlay(:,1))*(log(pplay(:,1)/ps))
389     &                             )
390     &                         )**rcp
391
392      ENDDO   !of n=1,n_out
393
394!------------------------------------------------------------------------
395!------------------------------------------------------------------------
396! PART III : WSTAR COMPUTATION
397!------------------------------------------------------------------------
398!------------------------------------------------------------------------
399
400! Detection of the mixed-layer potential temperature
401! ------------
402
403! Nearest index for the pbl height
404
405      IF (calltherm) THEN
406
407      pbl_height_index(:)=1
408
409      DO k=1,nlay-1
410         DO ig=1, ngrid
411            IF (abs(zmax(ig)-zzlay(ig,k)) .lt.
412     &              abs(zmax(ig)-zzlay(ig,pbl_height_index(ig)))) THEN
413               pbl_height_index(ig)=k
414            ENDIF
415         ENDDO
416      ENDDO
417
418! Potential temperature gradient
419
420      dteta(:,nlay)=0.
421      DO k=1,nlay-1
422         DO ig=1, ngrid
423         dteta(ig,k) = (ph(ig,k+1)-ph(ig,k))/(zzlay(ig,k+1)-zzlay(ig,k))
424         ENDDO
425      ENDDO
426
427! Computation of the pbl mixed layer temperature
428
429      DO ig=1, ngrid
430         ii=MINLOC(abs(dteta(ig,1:pbl_height_index(ig))))
431         pbl_teta(ig) = ph(ig,ii(1))
432      ENDDO
433
434!------------------------------------------------------------------------
435!------------------------------------------------------------------------
436! PART IV : VERTICAL_VELOCITY_VARIANCE/VERTICAL_TURBULENT_FLUX PROFILES
437!------------------------------------------------------------------------
438!------------------------------------------------------------------------
439
440! We follow Spiga et. al 2010 (QJRMS)
441! ------------
442
443      DO ig=1, ngrid
444         IF (zmax(ig) .gt. 0.) THEN
445            x(ig) = zout/zmax(ig)
446         ELSE
447            x(ig) = 999.
448         ENDIF
449      ENDDO
450
451      DO ig=1, ngrid
452         ! dimensionless vertical heat flux
453         IF (x(ig) .le. 0.3) THEN
454            dvhf(ig) = ((-3.85/log(x(ig)))+0.07*log(x(ig)))
455     &                                       *exp(-4.61*x(ig))
456         ELSEIF (x(ig) .le. 1.) THEN
457            dvhf(ig) = -1.52*x(ig) + 1.24
458         ELSE
459            dvhf(ig) = 0.
460         ENDIF
461         ! dimensionless vertical velocity variance
462         IF (x(ig) .le. 1.) THEN
463            dvvv(ig) = 2.05*(x(ig)**(2./3.))*(1.-0.64*x(ig))**2
464         ELSE
465            dvvv(ig) = 0.
466         ENDIF
467      ENDDO
468
469      vhf(:) = dvhf(:)*hfmax(:)
470      vvv(:) = dvvv(:)*(wstar_in(:))**2
471
472      ENDIF ! of if calltherm
473
474            call WRITEDIAGFI(ngrid,'rib_pbl',
475     &   'Richardson in pbl parameter','m/s',2,rib)
476            call WRITEDIAGFI(ngrid,'cdn_pbl',
477     &   'neutral momentum coef','m/s',2,cdn)
478            call WRITEDIAGFI(ngrid,'fm_pbl',
479     &   'momentum stab function','m/s',2,fm)
480            call WRITEDIAGFI(ngrid,'uv',
481     &   'wind norm first layer','m/s',2,sqrt(zu2))
482            call WRITEDIAGFI(ngrid,'uvtrue',
483     &   'wind norm first layer','m/s',2,sqrt(pu(:,1)**2+pv(:,1)**2))
484            call WRITEDIAGFI(ngrid,'chn_pbl',
485     &   'neutral momentum coef','m/s',2,chn)
486            call WRITEDIAGFI(ngrid,'fh_pbl',
487     &   'momentum stab function','m/s',2,fh)
488            call WRITEDIAGFI(ngrid,'B_pbl',
489     &   'flottabilité','m/',2,(ph(:,1)-pts(:))/pts(:))
490            call WRITEDIAGFI(ngrid,'zot_pbl',
491     &   'flottabilité','ms',2,pz0tcomp)
492       call WRITEDIAGFI(ngrid,'zz1','flottabilité','m',2,zzlay(:,1))
493
494      RETURN
495      END
Note: See TracBrowser for help on using the repository browser.