Module: NumRu::GAnalysis::QG
- Extended by:
- QG_common
- Defined in:
- lib/numru/ganalysis/qg.rb
Overview
module QG: quasi-geostrophic analysis module for Cartesian coordinates
Defined Under Namespace
Classes: Uninitialized
Constant Summary collapse
- @@bp =
a BetaPlane to be initialized by set_lat0
Uninitialized.new
Class Method Summary collapse
-
.bp ⇒ Object
returns the BetaPlane object created by initialization (((<set_lat0>))).
-
.div_h(gphys_u, gphys_v) ⇒ Object
horizontal divergence (Cartesian).
-
.f0 ⇒ Object
Returns the current f0 (the Coriolis parameter at the reference latitude).
-
.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 (Cartesian).
- .grad_hx(gphys) ⇒ Object
- .grad_hy(gphys) ⇒ Object
-
.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.
-
.psi2Qvector(psi) ⇒ Object
QG stream function -> the Q-vector.
-
.psi2ug_vg(psi) ⇒ Object
QG stream function -> geostrophic winds.
-
.psi_T2Qvector(psi, temp, p = nil) ⇒ Object
same as psi2Qvector, but temperature is given independently.
-
.set_lat0(lat0_or_latary) ⇒ Object
Initialize the QG module by setting a reference latitude.
Methods included from QG_common
cut_bottom, div_waf, extend_bottom, gp2gpref, gpd2qzz, gph2gpd_gpref, gph2gpref, gpref2n2, waf_tn2001
Class Method Details
.bp ⇒ Object
returns the BetaPlane object created by initialization (((<set_lat0>)))
476 477 478 |
# File 'lib/numru/ganalysis/qg.rb', line 476 def bp @@bp end |
.div_h(gphys_u, gphys_v) ⇒ Object
horizontal divergence (Cartesian)
630 631 632 |
# File 'lib/numru/ganalysis/qg.rb', line 630 def div_h(gphys_u, gphys_v) @@bp.div_h(gphys_u, gphys_v) end |
.f0 ⇒ Object
Returns the current f0 (the Coriolis parameter at the reference latitude)
481 |
# File 'lib/numru/ganalysis/qg.rb', line 481 def f0; @@bp.f0; end |
.gph2psi(gph, gpref) ⇒ Object
geopotential height -> QG stream function
515 516 517 518 519 520 521 |
# File 'lib/numru/ganalysis/qg.rb', line 515 def gph2psi(gph, gpref) gpd = gph * Met::g - gpref psi = gpd / @@bp.f0 psi.name = "psi" psi.long_name = "QG stream function" psi end |
.gph2psi_gpref(gph) ⇒ Object
geopotential height -> QG stream function and the reference geopotential
506 507 508 509 510 511 512 |
# File 'lib/numru/ganalysis/qg.rb', line 506 def gph2psi_gpref(gph) gpd, gpref = gph2gpd_gpref(gph) psi = gpd / @@bp.f0 psi.name = "psi" psi.long_name = "QG stream function" [psi, gpref] end |
.gph2q(gph) ⇒ Object
geopotential height to quasi-geostrophic potential vorticity (QGPV)
486 487 488 489 490 |
# File 'lib/numru/ganalysis/qg.rb', line 486 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
493 494 495 496 497 |
# File 'lib/numru/ganalysis/qg.rb', line 493 def gph2qb(gph) psi, gpref = gph2psi_gpref(gph) n2 = gpref2n2(gpref) psi2qb(psi, n2) end |
.gph2ug_vg(gph) ⇒ Object
geopotential height -> geostrophic winds
500 501 502 503 |
# File 'lib/numru/ganalysis/qg.rb', line 500 def gph2ug_vg(gph) psi, gpref = gph2psi_gpref(gph) psi2ug_vg(psi) end |
.grad_h(gphys) ⇒ Object
horizontal gradient (Cartesian)
619 620 621 |
# File 'lib/numru/ganalysis/qg.rb', line 619 def grad_h(gphys) @@bp.grad_h(gphys) end |
.grad_hx(gphys) ⇒ Object
625 626 627 |
# File 'lib/numru/ganalysis/qg.rb', line 625 def grad_hx(gphys) @@bp.grad_x(gphys) end |
.grad_hy(gphys) ⇒ Object
622 623 624 |
# File 'lib/numru/ganalysis/qg.rb', line 622 def grad_hy(gphys) @@bp.grad_y(gphys) end |
.psi2q(psi, n2, perturbation = false) ⇒ Object
QG stream function -> QGPV
524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 |
# File 'lib/numru/ganalysis/qg.rb', line 524 def psi2q(psi, n2, perturbation=false) x, y = @@bp.get_x_y(psi) bc = GPhys::Derivative::CYCLIC_OR_LINEAR f0 = @@bp.f0 vor = psi.deriv2nd(0,bc,x) + psi.deriv2nd(1,bc,y) if !perturbation avor = vor + (f0 + @@bp.beta*y) avor.name = "qgavor" avor.long_name = "QG abs vor" else vor.name = "qgvor" vor.long_name = "QG vorticity" avor = vor end qzz = gpd2qzz(psi, n2) * (f0*f0) q = avor + qzz q.name = "q" q.long_name = "QG PV" [q, avor, qzz] end |
.psi2qb(psi, n2, perturbation = false) ⇒ Object
same as psi2q, but the QGPV is extended to reflect the lowest-level temperature anomalies
547 548 549 550 551 552 |
# File 'lib/numru/ganalysis/qg.rb', line 547 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 |
.psi2Qvector(psi) ⇒ Object
QG stream function -> the Q-vector
568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 |
# File 'lib/numru/ganalysis/qg.rb', line 568 def psi2Qvector(psi) bc = GPhys::Derivative::CYCLIC_OR_LINEAR f0 = @@bp.f0 x, y = @@bp.get_x_y(psi) p = Met.get_prs(psi).convert_units("Pa") psi_x = psi.threepoint_O2nd_deriv(0,bc,x) psi_y = psi.threepoint_O2nd_deriv(1,bc,y) psi_xp = psi_x.threepoint_O2nd_deriv(2,bc,p) psi_yp = psi_y.threepoint_O2nd_deriv(2,bc,p) psi_xy = psi_x.threepoint_O2nd_deriv(1,bc,y) psi_xx = psi.deriv2nd(0,bc,x) psi_yy = psi.deriv2nd(1,bc,y) q1 = (-psi_xy*psi_xp + psi_xx*psi_yp) * f0 q2 = (-psi_yy*psi_xp + psi_xy*psi_yp) * f0 q1.name = q1.long_name = "Q1" q2.name = q2.long_name = "Q2" [q1,q2] end |
.psi2ug_vg(psi) ⇒ Object
QG stream function -> geostrophic winds
555 556 557 558 559 560 561 562 563 564 565 |
# File 'lib/numru/ganalysis/qg.rb', line 555 def psi2ug_vg(psi) bc = GPhys::Derivative::CYCLIC_OR_LINEAR x, y = @@bp.get_x_y(psi) vg = psi.cderiv(0,bc,x) ug = -psi.threepoint_O2nd_deriv(1,bc,y) ug.name = "ug" vg.name = "vg" ug.long_name = "ug" vg.long_name = "vg" [ug, vg] end |
.psi_T2Qvector(psi, temp, p = nil) ⇒ Object
same as psi2Qvector, but temperature is given independently
p (nil or UNumeric or VArray or..) : specify pressure if the input data does not have a pressure axis
591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 |
# File 'lib/numru/ganalysis/qg.rb', line 591 def psi_T2Qvector(psi, temp, p=nil) bc = GPhys::Derivative::CYCLIC_OR_LINEAR f0 = @@bp.f0 x, y = @@bp.get_x_y(psi) if p if p.respond_to?(:convert_units) p = p.convert_units("Pa") else # UNumeric p = p.convert2("Pa") end else p = Met.get_prs(psi).convert_units("Pa") end psi_xy = psi.threepoint_O2nd_deriv(0,bc,x).threepoint_O2nd_deriv(1,bc,y) psi_xx = psi.deriv2nd(0,bc,x) psi_yy = psi.deriv2nd(1,bc,y) t_x = temp.threepoint_O2nd_deriv(0,bc,x) t_y = temp.threepoint_O2nd_deriv(1,bc,y) q1 = (psi_xy*t_x - psi_xx*t_y) * (Met::R / p) q2 = (psi_yy*t_x - psi_xy*t_y) * (Met::R / p) #puts "@@@ psi_T2Qvector @@@", psi.units, psi_xx.units, t_x.units, q1.units q1.name = q1.long_name = "Q1" q2.name = q2.long_name = "Q2" [q1,q2] end |