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

Last change on this file since 1769 was 1238, checked in by emillour, 11 years ago

Mars GCM and common dynamics:

Common dynamics:

  • correction in inidissip (only matters in Martian case)
  • added correction in addfi on theta to account for surface pressure change.

Mars GCM:
Some fixes and updates:

  • addfi (dyn3d): Add correction on theta when surface pressure changes
  • vdif_cd (phymars): Correction for coefficients in stable nighttime case
  • jthermcalc (aeronomars): Fix for some pathological cases (further investigations on the origin of these is needed)

EM

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