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

Last change on this file since 657 was 648, checked in by emillour, 13 years ago

Mars GCM:

  • some syntax corrections in thermcall_main_mars, vdif_cd, pbl_parameters which cause problems when compiling with some strict compilers (g95, gfortran)
  • added an initialisation of 'varian' in initracer for cases when using conrath dust; because that value can be is used elsewhere (e.g. surfacearea)

JYC+EM

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