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
-
.cos_phi(gphys) ⇒ Object
cosine of latitude.
-
.div_h(fx, fy) ⇒ Object
horizontal divergence (spherical).
-
.div_waf(*args) ⇒ Object
divergence of wave activity flux (redirected to ((<div_waf>))).
-
.gph2psi(gph, gpref) ⇒ Object
geopotential height -> QG stream function.
-
.gph2psi_gpref(gph) ⇒ Object
geopotential height -> QG stream function and the reference geopotential.
-
.gph2q(gph) ⇒ Object
geopotential height to quasi-geostrophic potential vorticity (QGPV).
-
.gph2qb(gph) ⇒ Object
same as gph2q, but the QGPV is extended to reflect the lowest-level temperature anomalies.
-
.gph2ug_vg(gph) ⇒ Object
geopotential height -> geostrophic winds.
-
.grad_h(gphys) ⇒ Object
horizontal gradient (spherical).
-
.psi2q(psi, n2, perturbation = false) ⇒ Object
QG stream function -> QGPV.
-
.psi2qb(psi, n2, perturbation = false) ⇒ Object
same as psi2q, but the QGPV is extended to reflect the lowest-level temperature anomalies.
-
.psi2ug_vg(psi) ⇒ Object
QG stream function -> geostrophic winds.
-
.waf_plumb1986_B1(psi, n2) ⇒ Object
Flux of the pseudo-momentum in y direction by Plumb (1986).
-
.waf_plumb1986_B2(psi, n2) ⇒ Object
Flux of the pseudo-momentum in x direction by Plumb (1986).
Methods included from QG_common
cut_bottom, div_waf, extend_bottom, gp2gpref, gpd2qzz, gph2gpd_gpref, gph2gpref, gpref2n2
Methods included from QG_sphere_common
Class Method Details
.cos_phi(gphys) ⇒ Object
cosine of latitude
496 497 498 499 |
# File 'lib/numru/ganalysis/qg.rb', line 496 def cos_phi(gphys) lam, phi, = Planet::get_lambda_phi(gphys) phi.cos end |
.div_h(fx, fy) ⇒ Object
horizontal divergence (spherical)
529 530 531 |
# File 'lib/numru/ganalysis/qg.rb', line 529 def div_h(fx, fy) Planet::div_s(fx, fy) end |
.div_waf(*args) ⇒ Object
divergence of wave activity flux (redirected to ((<div_waf>)))
538 539 540 |
# File 'lib/numru/ganalysis/qg.rb', line 538 def self.div_waf(*args) super(*args) # defined in QG_common end |
.gph2psi(gph, gpref) ⇒ Object
geopotential height -> QG stream function
461 462 463 464 465 466 467 468 |
# File 'lib/numru/ganalysis/qg.rb', line 461 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
451 452 453 454 455 456 457 458 |
# File 'lib/numru/ganalysis/qg.rb', line 451 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)
430 431 432 433 434 |
# File 'lib/numru/ganalysis/qg.rb', line 430 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
437 438 439 440 441 |
# File 'lib/numru/ganalysis/qg.rb', line 437 def gph2qb(gph) psi, gpref = gph2psi_gpref(gph) n2 = gpref2n2(gpref) psi2qb(psi, n2) end |
.gph2ug_vg(gph) ⇒ Object
geopotential height -> geostrophic winds
444 445 446 447 448 |
# File 'lib/numru/ganalysis/qg.rb', line 444 def gph2ug_vg(gph) gpx, gpy = Planet::grad_s(gph * Met::g) f = f_mask0(gph) [ -gpy/f, gpx/f] end |
.grad_h(gphys) ⇒ Object
horizontal gradient (spherical)
524 525 526 |
# File 'lib/numru/ganalysis/qg.rb', line 524 def grad_h(gphys) Planet::grad_s(gphys) end |
.psi2q(psi, n2, perturbation = false) ⇒ Object
QG stream function -> QGPV
471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 |
# File 'lib/numru/ganalysis/qg.rb', line 471 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
503 504 505 506 507 508 |
# File 'lib/numru/ganalysis/qg.rb', line 503 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
511 512 513 514 515 516 517 518 519 520 521 |
# File 'lib/numru/ganalysis/qg.rb', line 511 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).
574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 |
# File 'lib/numru/ganalysis/qg.rb', line 574 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).
547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 |
# File 'lib/numru/ganalysis/qg.rb', line 547 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 |