Module: NumRu::GAnalysis::QG_sphere

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

Overview

QG on sphere with non-divergent but inaccurate geostrophic wind

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

.cos_phi(gphys) ⇒ Object

cosine of latitude



711
712
713
714
# File 'lib/numru/ganalysis/qg.rb', line 711

def cos_phi(gphys)
  lam, phi, = Planet::get_lambda_phi(gphys)
  phi.cos
end

.div_h(fx, fy) ⇒ Object

horizontal divergence (spherical)



750
751
752
# File 'lib/numru/ganalysis/qg.rb', line 750

def div_h(fx, fy)
  Planet::div_s(fx, fy)
end

.div_waf(*args) ⇒ Object

divergence of wave activity flux (redirected to ((<div_waf>)))



759
760
761
# File 'lib/numru/ganalysis/qg.rb', line 759

def self.div_waf(*args)
  super(*args)    # defined in QG_common
end

.gph2psi(gph, gpref) ⇒ Object

geopotential height -> QG stream function



676
677
678
679
680
681
682
683
# File 'lib/numru/ganalysis/qg.rb', line 676

def gph2psi(gph, gpref)
  gpd = gph * Met::g - gpref
  f = f_mask0(gph)
  psi = gpd / f
  psi.name = "psi"
  psi.long_name = "QG stream function"
  psi
end

.gph2psi_gpref(gph) ⇒ Object

geopotential height -> QG stream function and the reference geopotential



666
667
668
669
670
671
672
673
# File 'lib/numru/ganalysis/qg.rb', line 666

def gph2psi_gpref(gph)
  gpd, gpref = gph2gpd_gpref(gph)
  f = f_mask0(gph)
  psi = gpd / f
  psi.name = "psi"
  psi.long_name = "QG stream function"
  [psi, gpref]
end

.gph2q(gph) ⇒ Object

geopotential height to quasi-geostrophic potential vorticity (QGPV)



646
647
648
649
650
# File 'lib/numru/ganalysis/qg.rb', line 646

def gph2q(gph)
  psi, gpref = gph2psi_gpref(gph)
  n2 = gpref2n2(gpref)
  psi2q(psi, n2)
end

.gph2qb(gph) ⇒ Object

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



653
654
655
656
657
# File 'lib/numru/ganalysis/qg.rb', line 653

def gph2qb(gph)
  psi, gpref = gph2psi_gpref(gph)
  n2 = gpref2n2(gpref)
  psi2qb(psi, n2)
end

.gph2ug_vg(gph) ⇒ Object

geopotential height -> geostrophic winds



660
661
662
663
# File 'lib/numru/ganalysis/qg.rb', line 660

def gph2ug_vg(gph)
  psi, gpref = gph2psi_gpref(gph)
  psi2ug_vg(psi)
end

.grad_h(gphys) ⇒ Object

horizontal gradient (spherical)



739
740
741
# File 'lib/numru/ganalysis/qg.rb', line 739

def grad_h(gphys)
  Planet::grad_s(gphys)
end

.grad_hx(gphys) ⇒ Object



742
743
744
# File 'lib/numru/ganalysis/qg.rb', line 742

def grad_hx(gphys)
  Planet::grad_sx(gphys)
end

.grad_hy(gphys) ⇒ Object



745
746
747
# File 'lib/numru/ganalysis/qg.rb', line 745

def grad_hy(gphys)
  Planet::grad_sy(gphys)
end

.psi2q(psi, n2, perturbation = false) ⇒ Object

QG stream function -> QGPV



686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
# File 'lib/numru/ganalysis/qg.rb', line 686

def psi2q(psi, n2, perturbation=false)
  ug, vg = psi2ug_vg(psi)

  if !perturbation
    avor = Planet::absvor_s(ug,vg)
    avor.name = "qgavor"
    avor.long_name = "QG abs vor"
  else
    vor = Planet::vor_s(ug,vg)
    vor.name = "qgvor"
    vor.long_name = "QG vorticity"
    avor = vor
  end

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

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

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

.psi2qb(psi, n2, perturbation = false) ⇒ Object

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



718
719
720
721
722
723
# File 'lib/numru/ganalysis/qg.rb', line 718

def psi2qb(psi, n2, perturbation=false)
  psie = extend_bottom(psi, nil)
  n2e = extend_bottom(n2, nil)
  results = psi2q(psie, n2e, perturbation)
  results.collect{|z| cut_bottom(z)}
end

.psi2ug_vg(psi) ⇒ Object

QG stream function -> geostrophic winds



726
727
728
729
730
731
732
733
734
735
736
# File 'lib/numru/ganalysis/qg.rb', line 726

def psi2ug_vg(psi)
  f = f_mask0(psi)
  gpx, gpy = Planet::grad_s(psi)
  vg = gpx
  ug = -gpy
  ug.name = "ug"
  vg.name = "vg"
  ug.long_name = "ug"
  vg.long_name = "vg"
  [ug, vg]
end

.waf_plumb1986_B1(psi, n2) ⇒ Object

Flux of the pseudo-momentum in y direction by Plumb (1986). Specifically, B_1j in Eq.(2.9), but without the factor p. This flux is relative to the mean flow. Averaged over time (if the data is 4D).



795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
# File 'lib/numru/ganalysis/qg.rb', line 795

def waf_plumb1986_B1(psi, n2)
  psi_x, psi_y = Planet::grad_s(psi)
  psi_z = LogP.pcdata_dz( psi )
  f2 = f_mask0(psi) ** 2
  cosphi = cos_phi(psi)
  fx = -psi_x * psi_y * cosphi
  fy = (psi_x**2 - psi_y**2 + psi_z**2 * f2 / n2) * cosphi / 2.0
  fz = -psi_y * psi_z * (cosphi * f2) / n2
  fx.name = "B_11"
  fy.name = "B_12"
  fz.name = "B_13"
  fx.long_name = "x-comp of Px flux Plumb86"
  fy.long_name = "y-comp of Px flux Plumb86"
  fz.long_name = "z-comp of Px flux Plumb86"
  if psi.rank >= 4
    fx = fx.mean(-1)
    fy = fy.mean(-1)
    fz = fz.mean(-1)
  end
  [fx, fy, fz]
end

.waf_plumb1986_B2(psi, n2) ⇒ Object

Flux of the pseudo-momentum in x direction by Plumb (1986). Specifically, B_2j in Eq.(2.9), but without the factor p. This flux is relative to the mean flow. Averaged over time (if the data is 4D).



768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
# File 'lib/numru/ganalysis/qg.rb', line 768

def waf_plumb1986_B2(psi, n2)
  psi_x, psi_y = Planet::grad_s(psi)
  psi_z = LogP.pcdata_dz( psi )
  f2 = f_mask0(psi) ** 2
  cosphi = cos_phi(psi)
  fx = (psi_x**2 - psi_y**2 - psi_z**2 * f2 / n2) * cosphi / 2.0
  fy = psi_x * psi_y * cosphi 
  fz = psi_x * psi_z * (cosphi * f2) / n2
  fx.name = "B21"
  fy.name = "B22"
  fz.name = "B23"
  fx.long_name = "WAF B2x"
  fy.long_name = "WAF B2y"
  fz.long_name = "WAF B2z"
  if psi.rank >= 4
    fx = fx.mean(-1)
    fy = fy.mean(-1)
    fz = fz.mean(-1)
  end
  [fx, fy, fz]
end