source: trunk/LMDZ.MARS/libf/phymars/vdif_cd.F @ 1009

Last change on this file since 1009 was 765, checked in by acolaitis, 12 years ago

Added a low threshold subgrid gustiness in richardson computations for the LES, after observations of unreallistically high HFX values for very low winds conditions.

File size: 7.7 KB
RevLine 
[256]1      SUBROUTINE vdif_cd(ngrid,nlay,pz0,
[765]2     &           pg,pz,pu,pv,wstar,pts,ph,pcdv,pcdh
3#ifdef MESOSCALE
4     &                ,flag_LES
5#endif
6     &                  )
[38]7      IMPLICIT NONE
8c=======================================================================
9c
10c   Subject: computation of the surface drag coefficient using the
11c   -------  approch developed by Loui for ECMWF.
12c
13c   Author: Frederic Hourdin  15 /10 /93
[256]14c   Modified by : Arnaud Colaitis 03/08/11
[38]15c   -------
16c
17c   Arguments:
18c   ----------
19c
20c   inputs:
21c   ------
22c     ngrid            size of the horizontal grid
23c     pg               gravity (m s -2)
[256]24c     pz(ngrid,nlay)   height of layers
25c     pu(ngrid,nlay)   u component of the wind
26c     pv(ngrid,nlay)   v component of the wind
27c     pts(ngrid)       surface temperature
[38]28c     ph(ngrid)        potential temperature T*(p/ps)^kappa
29c
30c   outputs:
31c   --------
32c     pcdv(ngrid)      Cd for the wind
33c     pcdh(ngrid)      Cd for potential temperature
34c
35c=======================================================================
36c
37c-----------------------------------------------------------------------
38c   Declarations:
39c   -------------
40
[256]41#include "comcstfi.h"
[284]42#include "callkeys.h"
[256]43
[38]44c   Arguments:
45c   ----------
46
[256]47      INTEGER, INTENT(IN) :: ngrid,nlay
48      REAL, INTENT(IN) :: pz0(ngrid)
49      REAL, INTENT(IN) :: pg,pz(ngrid,nlay)
50      REAL, INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay)
51      REAL, INTENT(IN) :: pts(ngrid),ph(ngrid,nlay)
[499]52      REAL, INTENT(IN) :: wstar(ngrid)
[256]53      REAL, INTENT(OUT) :: pcdv(ngrid),pcdh(ngrid) ! momentum and heat drag coefficient
[38]54
[765]55#ifdef MESOSCALE
56      LOGICAL, INTENT(IN) :: flag_LES     !! pour LES avec isfflx!=0
57#endif
58
[38]59c   Local:
60c   ------
61
62      INTEGER ig
63
[267]64      REAL karman,nu    ! Von Karman constant and fluid kinematic viscosity
[38]65      LOGICAL firstcal
[267]66      DATA karman,nu/.41,0.001/
[38]67      DATA firstcal/.true./
[267]68      SAVE karman,nu
[38]69
[256]70c    Local(2):
71c    ---------
[276]72      REAL z1,zcd0
[256]73
74      REAL rib(ngrid)  ! Bulk Richardson number
[268]75      REAL rig(ngrid)  ! Gradient Richardson number
[256]76      REAL fm(ngrid) ! stability function for momentum
77      REAL fh(ngrid) ! stability function for heat
78      REAL z1z0,z1z0t ! ratios z1/z0 and z1/z0T
79
80c phim = 1+betam*zeta   or   (1-bm*zeta)**am
81c phih = alphah + betah*zeta    or   alphah(1.-bh*zeta)**ah
82      REAL betam, betah, alphah, bm, bh, lambda
83c ah and am are assumed to be -0.25 and -0.5 respectively
84
85      REAL cdn(ngrid),chn(ngrid)  ! neutral momentum and heat drag coefficient
86      REAL pz0t        ! initial thermal roughness length. (local)
87      REAL ric         ! critical richardson number
[265]88      REAL reynolds(ngrid)    ! reynolds number for UBL
[256]89      REAL prandtl(ngrid)     ! prandtl number for UBL
90      REAL pz0tcomp(ngrid)     ! computed z0t
91      REAL ite
92      REAL residual
[276]93      REAL zu2(ngrid)
[38]94c-----------------------------------------------------------------------
95c   couche de surface:
96c   ------------------
97
[256]98c Original formulation :
[38]99
[329]100      if(.not.callrichsl) then
[284]101
102      DO ig=1,ngrid
103         z1=1.E+0 + pz(ig,1)/pz0(ig)
104         zcd0=karman/log(z1)
105         zcd0=zcd0*zcd0
106         pcdv(ig)=zcd0
107         pcdh(ig)=zcd0
108      ENDDO
[276]109     
[256]110!      print*,'old : cd,ch; ',pcdv,pcdh
[284]111      else
[256]112
[284]113      reynolds(:)=0.
114
[256]115c New formulation (AC) :
116
117c phim = 1+betam*zeta   or   (1-bm*zeta)**am
118c phih = alphah + betah*zeta    or   alphah(1.-bh*zeta)**ah
[265]119c am=-0.25, ah=-0.5
[256]120
121      pz0t = 0.     ! for the sake of simplicity
122      pz0tcomp(:) = 0.1*pz0(:)
123      rib(:)=0.
[276]124
[256]125      pcdv(:)=0.
126      pcdh(:)=0.
127
128c this formulation assumes alphah=1., implying betah=betam
129c We use Dyer et al. parameters, as they cover a broad range of Richardson numbers :
130      bm=16.            !UBL
131      bh=16.            !UBL
132      alphah=1.
133      betam=5.         !SBL
134      betah=5.         !SBL
135      lambda=(sqrt(bh/bm))/alphah
136      ric=betah/(betam**2)
137
[38]138      DO ig=1,ngrid
139
[256]140         ite=0.
141         residual=abs(pz0tcomp(ig)-pz0t)
142
143         do while((residual .gt. 0.01*pz0(ig)) .and.  (ite .lt. 10.))
144
145         pz0t=pz0tcomp(ig)
146
147         if ((pu(ig,1) .ne. 0.) .or. (pv(ig,1) .ne. 0.)) then
148
149c Classical Richardson number formulation
150
151c         rib(ig) = (pg/ph(ig,1))*pz(ig,1)*(ph(ig,1)-pts(ig))
152c     &           /(pu(ig,1)*pu(ig,1) + pv(ig,1)*pv(ig,1))
153
154c Richardson number formulation proposed by D.E. England et al. (1995)
155
[499]156!         zu2=MAX(pu(ig,1)*pu(ig,1) + pv(ig,1)*pv(ig,1),0.25*wstar(ig)**2)
[268]157!         zu2=pu(ig,1)*pu(ig,1) + pv(ig,1)*pv(ig,1)
[636]158!         zu2(ig)=MAX(pu(ig,1)*pu(ig,1) + pv(ig,1)*pv(ig,1),             &
159!     &      (0.3*wstar(ig))**2)
[648]160          zu2(ig)=pu(ig,1)*pu(ig,1) + pv(ig,1)*pv(ig,1)
[636]161     &     + (log(1.+0.7*wstar(ig) + 2.3*wstar(ig)**2))**2
[765]162#ifdef MESOSCALE
163          if(flag_LES) then
164             zu2(ig)=MAX(zu2(ig),1.)
165          endif
166#endif
[499]167!       zu2(ig)=pu(ig,1)*pu(ig,1) + pv(ig,1)*pv(ig,1) + (0.5*wstar(ig))**2
[268]168
[499]169               ! we add the wstar to simulate
[268]170               ! bulk Ri changes due to subgrid wind feeding the thermals
171
172!          rig(ig) = (pg/ph(ig,1))*((pz(ig,1)-pz0(ig))**2
173!     &         /(pz(ig,1)-pz0t))*(ph(ig,1)-pts(ig))
174!     &         /zu2
175
[604]176          rib(ig) = (pg/pts(ig))
[268]177!     &      *pz(ig,1)*pz0(ig)/sqrt(pz(ig,1)*pz0t)
178     &      *sqrt(pz(ig,1)*pz0(ig))
[256]179     &      *(((log(pz(ig,1)/pz0(ig)))**2)/(log(pz(ig,1)/pz0t)))
180     &      *(ph(ig,1)-pts(ig))
[276]181     &  /zu2(ig)
[256]182
183         else
184         print*,'warning, infinite Richardson at surface'
185         print*,pu(ig,1),pv(ig,1)
186         rib(ig) = ric          ! traiter ce cas !
187         endif
188
189         z1z0=pz(ig,1)/pz0(ig)
190         z1z0t=pz(ig,1)/pz0t
191
192         cdn(ig)=karman/log(z1z0)
193         cdn(ig)=cdn(ig)*cdn(ig)
194         chn(ig)=cdn(ig)*log(z1z0)/log(z1z0t)
195
196c Stable case :
197      if (rib(ig) .gt. 0.) then
198c From D.E. England et al. (95)
199      prandtl(ig)=1.
200         if(rib(ig) .lt. ric) then
201c Assuming alphah=1. and bh=bm for stable conditions :
202            fm(ig)=((ric-rib(ig))/ric)**2
203            fh(ig)=((ric-rib(ig))/ric)**2
204         else
[260]205c For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface
[256]206            fm(ig)=0.
207            fh(ig)=0.
208         endif
209c Unstable case :
210      else
211c From D.E. England et al. (95)
212         fm(ig)=sqrt(1.-lambda*bm*rib(ig))
213         fh(ig)=(1./alphah)*((1.-lambda*bh*rib(ig))**0.5)*
214     &                     (1.-lambda*bm*rib(ig))**0.25
215         prandtl(ig)=alphah*((1.-lambda*bm*rib(ig))**0.25)/
216     &             ((1.-lambda*bh*rib(ig))**0.5)
217      endif
218
[276]219       reynolds(ig)=karman*sqrt(fm(ig))
220     &      *sqrt(zu2(ig))
221c     &      *sqrt(pu(ig,1)**2 + pv(ig,1)**2)
[267]222     &       *pz0(ig)/(log(z1z0)*nu)
[256]223       pz0tcomp(ig)=pz0(ig)*exp(-karman*7.3*
224     &              (reynolds(ig)**0.25)*(prandtl(ig)**0.5))
225
226         
227      residual = abs(pz0t-pz0tcomp(ig))
228      ite = ite+1
229!      print*, "iteration nnumber, residual",ite,residual
230
231      enddo  ! of while
232
233          pz0t=0.
234
235c Drag computation :
236
237         pcdv(ig)=cdn(ig)*fm(ig)
238         pcdh(ig)=chn(ig)*fh(ig)
[38]239       
[256]240      ENDDO
[265]241!
[256]242!      print*,'new : cd,ch; ',pcdv,pcdh
[38]243
[256]244! Some useful diagnostics :
[38]245
[284]246!       call WRITEDIAGFI(ngrid,'RiB',
[268]247!     &              'Bulk Richardson nb','no units',
[256]248!     &                         2,rib)
[268]249!                call WRITEDIAGFI(ngrid,'RiG',
250!     &              'Grad Richardson nb','no units',
251!     &                         2,rig)
[256]252!        call WRITEDIAGFI(ngrid,'Pr',
253!     &              'Prandtl nb','no units',
254!     &                         0,prandtl)
[276]255!       call WRITEDIAGFI(ngrid,'Re',
[256]256!     &              'Reynolds nb','no units',
257!     &                         0,reynolds)
258!        call WRITEDIAGFI(ngrid,'z0tcomp',
259!     &              'computed z0t','m',
260!     &                         2,pz0tcomp)
261
[276]262
[284]263      endif !of if call richardson surface layer
264
[38]265      RETURN
266      END
Note: See TracBrowser for help on using the repository browser.