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

Methods included from QG_common

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

Class Method Details

.bpObject

returns the BetaPlane object created by initialization (((<set_lat0>)))



266
267
268
# File 'lib/numru/ganalysis/qg.rb', line 266

def bp
  @@bp
end

.div_h(gphys_u, gphys_v) ⇒ Object

horizontal divergence (Cartesian)



414
415
416
# File 'lib/numru/ganalysis/qg.rb', line 414

def div_h(gphys_u, gphys_v)
  @@bp.div_h(gphys_u, gphys_v)
end

.f0Object

Returns the current f0 (the Coriolis parameter at the reference latitude)



271
# File 'lib/numru/ganalysis/qg.rb', line 271

def f0; @@bp.f0; end

.gph2psi(gph, gpref) ⇒ Object

geopotential height -> QG stream function



305
306
307
308
309
310
311
# File 'lib/numru/ganalysis/qg.rb', line 305

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



296
297
298
299
300
301
302
# File 'lib/numru/ganalysis/qg.rb', line 296

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)



276
277
278
279
280
# File 'lib/numru/ganalysis/qg.rb', line 276

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



283
284
285
286
287
# File 'lib/numru/ganalysis/qg.rb', line 283

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

.gph2ug_vg(gph) ⇒ Object

geopotential height -> geostrophic winds



290
291
292
293
# File 'lib/numru/ganalysis/qg.rb', line 290

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

.grad_h(gphys) ⇒ Object

horizontal gradient (Cartesian)



409
410
411
# File 'lib/numru/ganalysis/qg.rb', line 409

def grad_h(gphys)
  @@bp.grad_h(gphys)
end

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

QG stream function -> QGPV



314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
# File 'lib/numru/ganalysis/qg.rb', line 314

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



337
338
339
340
341
342
# File 'lib/numru/ganalysis/qg.rb', line 337

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



358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
# File 'lib/numru/ganalysis/qg.rb', line 358

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



345
346
347
348
349
350
351
352
353
354
355
# File 'lib/numru/ganalysis/qg.rb', line 345

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



381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
# File 'lib/numru/ganalysis/qg.rb', line 381

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

.set_lat0(lat0_or_latary) ⇒ Object

Initialize the QG module by setting a reference latitude.



261
262
263
# File 'lib/numru/ganalysis/qg.rb', line 261

def set_lat0(lat0_or_latary)
  @@bp = BetaPlane.new(lat0_or_latary)
end