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

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

--- AC 03/08/2011 ---
M 267 libf/phymars/physiq.F
<> Minor modification to pass Ch from vdifc to meso_inc_les

M 267 libf/phymars/surflayer_interpol.F
<> Major modification to the formulation of integrals

Now stable for most cases. Some cases with highly negative Monin Obukhov length
remain to be explored.

M 267 libf/phymars/vdif_cd.F
<> Added gustiness to the Richardson computation. Gustiness factor is for now of beta=1., after

several comparisons with LES aerodynamic conductances. May be subject to a minor change (+/- 0.1)
in the near future. (almost transparent for the user)

M 267 libf/phymars/vdifc.F
<> Minor modifications relative to variables.

M 267 libf/phymars/calltherm_mars.F90
<> Added a comment on a "sensitive" parameter that should not be changed without knowing the consequence !

M 267 libf/phymars/meso_inc/meso_inc_les.F
<> Changed the definition for HFX computation in the LES (to be discussed with Aymeric). New definition yields

very similar results to old one and follows a strict definition of what HFX should be.


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