Module: NumRu::GAnalysis::QG_sphere_div

Extended by:
QG_common, QG_sphere_common
Defined in:
lib/numru/ganalysis/qg.rb

Overview

QG on sphere with divergence in the geostrophic wind

The geostrophic wind is defined as

[ug, vg] = 1/f curl gpd,

where gpd is the geostrophic geopotential, which is the deviation of geopotential from some globally uniform background. A background can be defined as the global mean geopotential, as is done in the method gph2gpd_gpref.

Class Method Summary collapse

Methods included from QG_common

cut_bottom, div_waf, extend_bottom, gp2gpref, gpd2qzz, gph2gpd_gpref, gph2gpref, gpref2n2, waf_tn2001

Methods included from QG_sphere_common

f, f_mask0

Class Method Details

.gpd2q(gpd, n2) ⇒ Object



858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
# File 'lib/numru/ganalysis/qg.rb', line 858

def gpd2q(gpd, n2)
  ug, vg = gpd2ug_vg(gpd)

  avor = Planet::absvor_s(ug,vg)
  avor.name = "qgavor"
  avor.long_name = "QG abs vor"

  f = f_mask0(gpd)
  qzz = gpd2qzz(gpd, n2) * f

  q = avor + qzz
  q.name = "q"
  q.long_name = "QG PV"

  [q, avor, qzz, ug, vg]
end

.gpd2qb(gpd, n2) ⇒ Object

geopotential height (gpd: deviation from the reference profie) to the extended QGPV

qb: q including the contribution from the bottom boundary



878
879
880
881
882
883
# File 'lib/numru/ganalysis/qg.rb', line 878

def gpd2qb(gpd, n2)
  gpde = extend_bottom(gpd, nil)
  n2e = extend_bottom(n2, nil)
  results = gpd2q(gpde, n2e)
  results.collect{|z| cut_bottom(z)}
end

.gpd2ug_vg(gpd) ⇒ Object



885
886
887
888
889
890
891
892
893
894
895
# File 'lib/numru/ganalysis/qg.rb', line 885

def gpd2ug_vg(gpd)
  f = f_mask0(gpd)
  gpx, gpy = Planet::grad_s(gpd)
  vg = gpx/f
  ug = gpy* (-1/f)
  ug.name = "ug"
  vg.name = "vg"
  ug.long_name = "ug"
  vg.long_name = "vg"
  [ug, vg]
end

.gph2q(gph) ⇒ Object

geopotential height to quasi-geostrophic potential vorticity (QGPV)



837
838
839
840
841
# File 'lib/numru/ganalysis/qg.rb', line 837

def gph2q(gph)
  gpd, gpref = gph2gpd_gpref(gph)
  n2 = gpref2n2(gpref)
  gpd2q(gpd, n2)
end

.gph2qb(gph) ⇒ Object

same as gph2q, but the QGPV is extended to reflect the lowest-level temperature anomalies



844
845
846
847
848
# File 'lib/numru/ganalysis/qg.rb', line 844

def gph2qb(gph)
  gpd, gpref = gph2gpd_gpref(gph)
  n2 = gpref2n2(gpref)
  gpd2qb(gpd, n2)
end

.gph2ug_vg(gph) ⇒ Object

geopotential height -> geostrophic winds



851
852
853
854
855
# File 'lib/numru/ganalysis/qg.rb', line 851

def gph2ug_vg(gph)
  gpx, gpy = Planet::grad_s(gph * Met::g)
  f = f_mask0(gph)
  [ -gpy/f, gpx/f]
end