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

Last change on this file since 2800 was 2616, checked in by romain.vande, 3 years ago

LMDZ_MARS RV : Open_MP;
Put all the "save" variables as "!$OMP THREADPRIVATE" in phymars.
The code can now be tested, see README for more info

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