GFS Operational Physics Documentation  gsm/branches/DTC/phys-doc-all phys-doc-all R82971
Orographic Gravity Wave Drag and Mountain Blocking

This subroutine includes orographic gravity wave drag and mountain blocking. More...

Detailed Description

This subroutine includes orographic gravity wave drag and mountain blocking.

The time tendencies of zonal and meridional wind are altered to include the effect of mountain induced gravity wave drag from subgrid scale orography including convective breaking, shear breaking and the presence of critical levels.

Collaboration diagram for Orographic Gravity Wave Drag and Mountain Blocking:
subroutine gwdps (IM, IX, IY, KM, A, B, C, U1, V1, T1, Q1, KPBL, PRSI, DEL, PRSL, PRSLK, PHII, PHIL, DELTIM, KDT, HPRIME, OC, OA4, CLX4, THETA, SIGMA, GAMMA, ELVMAX, DUSFC, DVSFC, G, CP, RD, RV, IMX, nmtvr, cdmbgwd, me, lprnt, ipr)
 

Function/Subroutine Documentation

subroutine gwdps ( integer  IM,
integer  IX,
integer  IY,
integer  KM,
real(kind=kind_phys), dimension(iy,km)  A,
real(kind=kind_phys), dimension(iy,km)  B,
real(kind=kind_phys), dimension(iy,km)  C,
real(kind=kind_phys), dimension(ix,km)  U1,
real(kind=kind_phys), dimension(ix,km)  V1,
real(kind=kind_phys), dimension(ix,km)  T1,
real(kind=kind_phys), dimension(ix,km)  Q1,
integer, dimension(im)  KPBL,
real(kind=kind_phys), dimension(ix,km+1)  PRSI,
real(kind=kind_phys), dimension(ix,km)  DEL,
real(kind=kind_phys), dimension(ix,km)  PRSL,
real(kind=kind_phys), dimension(ix,km)  PRSLK,
real(kind=kind_phys), dimension(ix,km+1)  PHII,
real(kind=kind_phys), dimension(ix,km)  PHIL,
real(kind=kind_phys)  DELTIM,
integer  KDT,
real(kind=kind_phys), dimension(im)  HPRIME,
real(kind=kind_phys), dimension(im)  OC,
real(kind=kind_phys), dimension(iy,4)  OA4,
real(kind=kind_phys), dimension(iy,4)  CLX4,
real(kind=kind_phys), dimension(im)  THETA,
real(kind=kind_phys), dimension(im)  SIGMA,
real(kind=kind_phys), dimension(im)  GAMMA,
real(kind=kind_phys), dimension(im)  ELVMAX,
real(kind=kind_phys), dimension(im)  DUSFC,
real(kind=kind_phys), dimension(im)  DVSFC,
real(kind=kind_phys)  G,
real(kind=kind_phys)  CP,
real(kind=kind_phys)  RD,
real(kind=kind_phys)  RV,
integer  IMX,
integer  nmtvr,
real(kind=kind_phys), dimension(2)  cdmbgwd,
integer  me,
logical  lprnt,
integer  ipr 
)
Parameters
[in]IMhorizontal number of used pts
[in]IXhorizontal dimension
[in]IYhorizontal number of used pts
[in]KMvertical layer dimension
[in,out]Anon-linear tendency for v wind component
[in,out]Bnon-linear tendency for u wind component
[in,out]Cnon-linear tendency for temperature (not used)
[in]U1zonal wind component of model layer wind (m/s)
[in]V1meridional wind component of model layer wind (m/s)
[in]T1model layer mean temperature (K)
[in]Q1model layer mean specific humidity
[in]KPBLindex for the PBL top layer
[in]PRSIpressure at layer interfaces
[in]DELpositive increment of p/psfc across layer
[in]PRSLmean layer pressure
[in]PRSLKExner function at layer
[in]PHIIinterface geopotential ( \(m^2/s^2\))
[in]PHILlayer geopotential ( \(m^2/s^2\))
[in]DELTIMphysics time step in seconds
[in]KDTnumber of the current time step
[in]HPRIMEorographic standard deviation (m) (mtnvar(:,1))
[in]OCorographic Convexity (mtnvar(:,2))
[in]OA4orographic Asymmetry (mtnvar(:,3:6))
[in]CLX4Lx, the fractional area covered by the subgrid-scale orography higher than a critical height for a grid box with the interval \( \triangle x \) (mtnvar(:,7:10))
[in]THETAthe angle of the mtn with that to the east (x) axis (mtnvar(:,11))
[in]SIGMAorographic slope (mtnvar(:,13))
[in]GAMMAorographic anisotropy (mtnvar(:,12))
[in]ELVMAXorographic maximum (mtnvar(:,14))
[out]DUSFCu component of surface stress
[out]DVSFCv component of surface stress
[in]Gsee physcons::con_g
[in]CPsee physcons::con_cp
[in]RDsee physcons::con_tird
[in]RVsee physcons::con_rv
[in]IMXnumber of longitude points
[in]NMTVRnumber of topographic variables such as variance etc used in the GWD parameterization,current operational, nmtvr=14
[in]CDMBGWDmultiplication factors for cdmb and gwd
[in]MEpe number - used for debug prints
[in]LPRNTlogical print flag
[in]IPRcheck print point for debugging

General Algorithm

— Subgrid Mountain Blocking Section

  • Compute Brunt-Vaisala Frequency \(N\).
  • Find the dividing streamline height starting from the level above the maximum mountain height and processing downward.
  • Compute wind speed UDS

    \[ UDS=\max(\sqrt{U1^2+V1^2},minwnd) \]

    where \( minwnd=0.1 \), \(U1\) and \(V1\) are zonal and meridional wind components of model layer wind.
  • The dividing streamline height (idxzb), of a subgrid scale obstable, is found by comparing the potential (PE) and kinetic energies (EK) of the upstream large scale wind and subgrid scale air parcel movements. the dividing streamline is found when \(PE\geq EK\). Mountain-blocked flow is defined to exist between the surface and the dividing streamline height ( \(h_d\)), which can be found by solving an integral equation for \(h_d\):

    \[ \frac{U^{2}(h_{d})}{2}=\int_{h_{d}}^{H} N^{2}(z)(H-z)dz \]

    where \(H\) is the maximum subgrid scale elevation within the grid box of actual orography, \(h\), obtained from the GTOPO30 dataset from the U.S. Geological Survey.
  • Calculate \(ZLEN\), which sums up a number of contributions of elliptic obstables.

    \[ ZLEN=\sqrt{[\frac{h_{d}-z}{z+h'}]} \]

    where \(z\) is the height, \(h'\) is the orographic standard deviation (HPRIME).
  • Calculate the drag coefficient to vary with the aspect ratio of the obstable as seen by the incident flow (see eq.14 in Lott and Miller (1997) [42])

    \[ R=\frac{\cos^{2}\psi+\gamma\sin^{2}\psi}{\gamma\cos^{2}\psi+\sin^{2}\psi} \]

    where \(\psi\), which is derived from THETA, is the angle between the incident flow direction and the normal ridge direcion. \(\gamma\) is the orographic anisotropy (GAMMA).
  • In each model layer below the dividing streamlines, a drag from the blocked flow is exerted by the obstacle on the large scale flow. The drag per unit area and per unit height is written (eq.15 in Lott and Miller (1997) [42]):

    \[ D_{b}(z)=-C_{d}\max(2-\frac{1}{R},0)\rho\frac{\sigma}{2h'}ZLEN\max(\cos\psi,\gamma\sin\psi)\frac{UDS}{2} \]

    where \(C_{d}\) is a specified constant, \(\sigma\) is the orographic slope.

— Orographic Gravity Wave Drag Section

  • Calculate the reference level index: kref=max(2,KPBL+1). where KPBL is the index for the PBL top layer.
  • Calculate low-level horizontal wind direction, the derived orographic asymmetry parameter (OA), and the derived Lx (CLX).
  • Calculate enhancement factor (E),number of mountans (m') and aspect ratio constant.
    As in eq.(4.9),(4.10),(4.11) in Kim and Arakawa (1995) [36], we define m' and E in such a way that they depend on the geometry and location of the subgrid-scale orography through OA and the nonlinearity of flow above the orography through Fr. OC, which is the orographic convexity, and statistically determine how protruded (sharp) the subgrid-scale orography is, is included in the saturation flux G' in such a way that G' is proportional to OC. The forms of E,m' and G' are:

    \[ E(OA,F_{r_{0}})=(OA+2)^{\delta} \]

    \[ \delta=C_{E}F_{r_{0}}/F_{r_{c}} \]

    \[ m'(OA,CLX)=C_{m}\triangle x(1+CLX)^{OA+1} \]

    \[ G'(OC,F_{r_{0}})=\frac{F_{r_{0}}^2}{F_{r_{0}}^2+a^{2}} \]

    \[ a^{2}=C_{G}OC^{-1} \]

    where \(F_{r_{c}}(=1)\) is the critical Froude number, \(F_{r_{0}}\) is the Froude number. \(C_{E}\), \(C_{m}\), \(C_{G}\) are constants.
  • Calculate the reference-level drag \(\tau_{0}\) (eq.(4.8) in Kim and Arakawa (1995) [36]):

    \[ \tau_0=E\frac{m'}{\triangle x}\frac{\rho_{0}U_0^3}{N_{0}}G' \]

    where \(E\), \(m'\), and \(G'\) are the enhancement factor, "the number of mountains", and the flux function defined above, respectively.
  • Compute the drag above the reference level ( \(k\geq kref\)):
    • Calculate the ratio of the Scorer parameter ( \(R_{scor}\)).
      From a series of experiments, Kim and Arakawa (1995) [36] found that the magnitude of drag divergence tends to be underestimated by the revised scheme in low-level downstream regions with wave breaking. Therefore, at low levels when OA > 0 (i.e., in the "downstream" region) the saturation hypothesis is replaced by the following formula based on the ratio of the the Scorer parameter:

      \[ R_{scor}=\min \left[\frac{\tau_i}{\tau_{i+1}},1\right] \]

    • The drag above the reference level is expressed as:

      \[ \tau=\frac{m'}{\triangle x}\rho NUh_d^2 \]

      where \(h_{d}\) is the displacement wave amplitude. In the absence of wave breaking, the displacement amplitude for the \(i^{th}\) layer can be expressed using the drag for the layer immediately below. Thus, assuming \(\tau_i=\tau_{i+1}\), we can get:

      \[ h_{d_i}^2=\frac{\triangle x}{m'}\frac{\tau_{i+1}}{\rho_{i}N_{i}U_{i}} \]

  • The minimum Richardson number ( \(Ri_{m}\)) or local wave-modified Richardson number, which determines the onset of wave breaking, is expressed in terms of \(R_{i}\) and \(F_{r_{d}}=Nh_{d}/U\):

    \[ Ri_{m}=\frac{Ri(1-Fr_{d})}{(1+\sqrt{Ri}\cdot Fr_{d})^{2}} \]

    see eq.(4.6) in Kim and Arakawa (1995) [36].
    • Check stability to employ the 'saturation hypothesis' of Lindzen (1981) [39] except at tropospheric downstream regions.
      Wave breaking occurs when \(Ri_{m}<Ri_{c}=0.25\). Then Lindzen's wave saturation hypothesis resets the displacement amplitude \(h_{d}\) to that corresponding to \(Ri_{m}=0.25\), we obtain the critical \(h_{d}\)(or \(h_{c}\)) expressed in terms of the mean values of \(U\), \(N\), and \(Ri\) ( eq.(4.7) in Kim and Arakawa (1995) [36]):

      \[ h_{c}=\frac{U}{N}\left\{2(2+\frac{1}{\sqrt{Ri}})^{1/2}-(2+\frac{1}{\sqrt{Ri}})\right\} \]

      if \(Ri_{m}\leq Ri_{c}\), obtain \(\tau\) from the drag above the reference level by using \(h_{c}\) computed above; otherwise \(\tau\) is unchanged (note: scaled by the ratio of the Scorer paramter).
  • Calculate outputs: A, B, DUSFC, DVSFC (see parameter description).
    • Below the dividing streamline height (k < idxzb), mountain blocking( \(D_{b}\)) is applied.
    • Otherwise (k>= idxzb), orographic GWD ( \(\tau\)) is applied.

Definition at line 188 of file gwdps.f.

Referenced by gbphys().

Here is the caller graph for this function: