Module: NumRu::GPhys::EP_Flux

Extended by:
Misc::EMath
Includes:
Misc::EMath
Defined in:
lib/numru/gphys/ep_flux.rb

Constant Summary collapse

Deriv_methods =

C_p = 1004.0 # specific heat at constant pressure of the earth’s atmosphere [J.K-1.kg-1]

R   = 287.0  # gas constant per unit mass for dry air of the earth [J.K-1.kg-1]
[ 'cderiv', 'threepoint_O2nd_deriv' ]
@@scale_height =

<<< module variable >>>

UNumeric.new(7000,  "m")
@@radius =

radius of the planet

UNumeric.new(6.37E6,"m")
@@rot_period =

rotation period of the planet

UNumeric.new(8.64E4,"s")
@@g_forces =
UNumeric.new(9.81,"m.s-2")
@@p00 =

gravitational acceleration in the surface

UNumeric.new(1.0E5,"Pa")
@@cp =
UNumeric.new(1004.0, "J.K-1.kg-1")
@@gas_const =

specific heat at constant pressure

UNumeric.new(287.0, "J.K-1.kg-1")
@@deriv_method =

gas constant per molecular mass

Proc.new{|*args|
  GPhys::Derivative::threepoint_O2nd_deriv(*args)
}

Class Method Summary collapse

Class Method Details

.cpObject



497
498
499
# File 'lib/numru/gphys/ep_flux.rb', line 497

def cp
	@@cp
end

.cp=(cp) ⇒ Object



500
501
502
503
# File 'lib/numru/gphys/ep_flux.rb', line 500

def cp=(cp)
	@@cp = cp
	return nil
end

.deriv(*args) ⇒ Object



542
543
544
# File 'lib/numru/gphys/ep_flux.rb', line 542

def deriv(*args)
  @@deriv_method.call(*args)
end

.div_sphere(gp_fy, gp_fz, yzdims = [0,1]) ⇒ Object

Raises:

  • (ArgumentError)


628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
# File 'lib/numru/gphys/ep_flux.rb', line 628

def div_sphere(gp_fy, gp_fz, yzdims=[0,1])
	raise ArgumentError,"yzdims's size (#{yzdims.size}) must be 2." if yzdims.size != 2
	## get axis and name
	ax_lat = gp_fy.axis(yzdims[0])    # Axis of latitude
	ax_z   = gp_fy.axis(yzdims[1])    # Axis of vertical
  lat_nm, z_nm = ax_lat.pos.name, ax_z.pos.name
	gp_lat, gp_z = make_gphys(ax_lat, ax_z)
	## convert
	gp_z =   to_z_if_pressure(gp_z)   # P => z=-H*log(P/P00) (units-based)
	gp_lat = to_rad_if_deg(gp_lat)    # deg => rad (unit convesion)

	## replace grid (without duplicating data)
	grid = gp_fy.grid_copy
	cp_grid = gp_fy.grid_copy         # saved to use in outputs
	grid.axis(lat_nm).pos = gp_lat.data
	grid.axis(z_nm).pos = gp_z.data
	gp_fy = GPhys.new(grid, gp_fy.data)
	gp_fz = GPhys.new(grid, gp_fz.data)

	## d_F_phi_dz
	a_cos_lat = @@radius * cos(gp_lat)
  remove_0_at_poles(a_cos_lat)	
	d_gp_fy_d_phi = deriv(gp_fy * cos(gp_lat), lat_nm)
	## d_F_z_dz
	d_gp_fz_d_z =   deriv(gp_fz, z_nm)
	f_div = ( d_gp_fy_d_phi / a_cos_lat )  + d_gp_fz_d_z

	f_div.data.name = "epflx_div"
	f_div.data.set_att("long_name", "EP Flux divergence")
	## convert with past grid
	return GPhys.new(cp_grid, f_div.data)
end

.eddy_products(gp_u, gp_v, gp_w, gp_t, dimname) ⇒ Object



754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
# File 'lib/numru/gphys/ep_flux.rb', line 754

def eddy_products(gp_u, gp_v, gp_w, gp_t, dimname)
	# get zonal_eddy 
	u_dash = gp_u - gp_u.mean(dimname)
	v_dash = gp_v - gp_v.mean(dimname)
	w_dash = gp_w - gp_w.mean(dimname)
	t_dash = gp_t - gp_t.mean(dimname)

	# get eddy_product 
	uv_dash = u_dash*v_dash  # u'v'
	vt_dash = v_dash*t_dash  # v't'
	uw_dash = u_dash*w_dash  # u'w'

	# set attribute
	uv_dash.data.set_att("long_name", "U'V'")
	vt_dash.data.set_att("long_name", "V'T'")
	uw_dash.data.set_att("long_name", "U'W'")
	uv_dash.data.rename!("uv_dash")
	vt_dash.data.rename!("vt_dash")
	uw_dash.data.rename!("uw_dash")	

	return uv_dash.mean(dimname), vt_dash.mean(dimname), uw_dash.mean(dimname)
end

.ep_full_sphere(gp_u, gp_v, gp_w, gp_t, flag_temp_or_theta = true, xyzdims = [0,1,2]) ⇒ Object

<<< calculation method >>> ———————————————

Raises:

  • (ArgumentError)


547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
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
617
618
619
620
621
622
623
624
625
626
# File 'lib/numru/gphys/ep_flux.rb', line 547

def ep_full_sphere(gp_u, gp_v, gp_w, gp_t, 
                   flag_temp_or_theta=true, xyzdims=[0,1,2]) ## get axis and name
	raise ArgumentError,"xyzdims's size (#{xyzdims.size}) must be 3." if xyzdims.size != 3 
	ax_lon = gp_u.axis(xyzdims[0]) # Axis of longitude 
	ax_lat = gp_u.axis(xyzdims[1]) # Axis of latitude 
	ax_z =   gp_u.axis(xyzdims[2]) # Axis of vertical 
	lon_nm, lat_nm, z_nm = ax_lon.pos.name, ax_lat.pos.name, ax_z.pos.name
	gp_lon, gp_lat, gp_z = make_gphys(ax_lon, ax_lat, ax_z)
	
	## convert axes
	gp_z = to_z_if_pressure(gp_z)     # P => z=-H*log(P/P00) (units-based)
	gp_lon = to_rad_if_deg(gp_lon)    # deg => rad (unit convesion)
	gp_lat = to_rad_if_deg(gp_lat)    # deg => rad (unit convesion)
	gp_w = to_w_if_omega(gp_w, gp_z)  # dP/dt => dz/dt (units-based)
	gp_t = to_theta_if_temperature(gp_t, gp_z, flag_temp_or_theta)	
                # temperature => potential temperature (if flag is true)

	## replace grid (without duplicating data)
	grid = gp_u.grid_copy
	old_grid = gp_u.grid_copy                 # saved to use in outputs
	grid.axis(lon_nm).pos = gp_lon.data       # in radian
	grid.axis(lat_nm).pos = gp_lat.data       # in radian
	grid.axis(z_nm).pos = gp_z.data           # log-p height
	gp_u = GPhys.new(grid, gp_u.data)
	gp_v = GPhys.new(grid, gp_v.data)
	gp_w = GPhys.new(grid, gp_w.data)
	gp_t = GPhys.new(grid, gp_t.data)
	## get each term
  #  needed in F_y and F_z
	uv_dash, vt_dash, uw_dash = eddy_products(gp_u, gp_v, gp_w, gp_t, lon_nm)
	theta_mean = gp_t.mean(lon_nm)
	dtheta_dz = deriv(theta_mean, z_nm)
	cos_lat = cos(gp_lat)
  a_cos_lat = @@radius * cos_lat
  	a_cos_lat.data.rename!('a_cos_lat')
	a_cos_lat.data.set_att('long_name', 'radius * cos_lat')
  remove_0_at_poles(a_cos_lat)
  #  needed in F_y only
	u_mean = gp_u.mean(lon_nm)
	du_dz  = deriv(u_mean, z_nm)
  #  needed in F_z only
	f_cor = 2 * (2 * PI / @@rot_period) * sin(gp_lat) 
  	f_cor.data.rename!('f_cor')
	f_cor.data.set_att('long_name', 'Coriolis parameter')
	ducos_dphi = deriv( u_mean * cos_lat, lat_nm)
	avort = (-ducos_dphi/a_cos_lat) + f_cor        # -- absolute vorticity
	avort.data.units = "s-1"
	avort.data.rename!('avort')
	avort.data.set_att('long_name', 'zonal mean absolute vorticity')

	## F_y, F_z
	sigma = exp(-gp_z/@@scale_height)
	epflx_y = ( - uv_dash + du_dz*vt_dash/dtheta_dz ) * cos_lat * sigma
	epflx_z = ( - uw_dash + avort*vt_dash/dtheta_dz ) * cos_lat * sigma
	epflx_y.data.name = "epflx_y"; epflx_z.data.name = "epflx_z"
	epflx_y.data.set_att("long_name", "EP flux y component")
	epflx_z.data.set_att("long_name", "EP flux z component")

	## v_rmean, w_rmean
	z_nm = gp_z.data.name    # change z_nm from pressure to z
	v_mean = gp_v.mean(lon_nm); w_mean = gp_w.mean(lon_nm)
	v_rmean = ( v_mean - deriv( (vt_dash/dtheta_dz*sigma), z_nm )/sigma )
	w_rmean = ( w_mean + deriv( (vt_dash/dtheta_dz*cos_lat), lat_nm )/a_cos_lat )
	v_rmean.data.name = "v_rmean"; w_rmean.data.name = "w_rmean"
	v_rmean.data.set_att("long_name", "residual zonal mean V")
	w_rmean.data.set_att("long_name", "residual zonal mean W")

	## convert with past grid
	gp_ary = [] # grid convertes gphyss into 
	grid_xmean = old_grid.delete_axes(lon_nm)
	[epflx_y, epflx_z, v_rmean, w_rmean, gp_lat, gp_z, u_mean, theta_mean, 
   uv_dash, vt_dash, uw_dash, dtheta_dz].each {|gp|  
	  if grid_xmean.shape.size != gp.shape.size
	    gp_ary << gp
	  else
	    gp_ary << GPhys.new(grid_xmean, gp.data) #back to the original grid
	  end
	}
	return gp_ary
end

.g_forcesObject



483
484
485
# File 'lib/numru/gphys/ep_flux.rb', line 483

def g_forces
	@@g_forces
end

.g_forces=(g) ⇒ Object



479
480
481
482
# File 'lib/numru/gphys/ep_flux.rb', line 479

def g_forces=(g)
	@@g_forces = g
	return nil
end

.gas_constObject



504
505
506
# File 'lib/numru/gphys/ep_flux.rb', line 504

def gas_const
	@@gas_const
end

.gas_const=(r) ⇒ Object



507
508
509
510
# File 'lib/numru/gphys/ep_flux.rb', line 507

def gas_const=(r)
	@@gas_const = r
	return nil
end

.get_constantsObject



521
522
523
524
# File 'lib/numru/gphys/ep_flux.rb', line 521

def get_constants
	return @@scale_height, @@radius, @@rot_period, @@g_forces, 
                                                 @@p00, @@cp, @@gas_const
end

.make_gphys(*ax_ary) ⇒ Object



661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
# File 'lib/numru/gphys/ep_flux.rb', line 661

def make_gphys(*ax_ary) 
           # it will be lost when new grid.rb released : to use delete_ax
	gp_ary = []           
	ax_ary.each{|ax|
	  if ax.is_a?(Axis)
	    ax_data = ax.pos
	    ax_grid = Grid.new(ax) 
	  elsif ax.is_a?(VArray)
	    ax_data = ax
	    ax_grid = Grid.new(Axis.new().set_pos(ax))
	  end
	  gp = GPhys.new(ax_grid, ax_data)
	  gp_ary << gp
	}
	return gp_ary
end

.p00Object



490
491
492
# File 'lib/numru/gphys/ep_flux.rb', line 490

def p00
	@@p00
end

.p00=(p00) ⇒ Object



493
494
495
496
# File 'lib/numru/gphys/ep_flux.rb', line 493

def p00=(p00)
	@@p00 = p00
	return nil
end

.preparate_for_vector_on_merdional_section(xax, zax) ⇒ Object



792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
# File 'lib/numru/gphys/ep_flux.rb', line 792

def preparate_for_vector_on_merdional_section(xax, zax)
	gp_x, gp_z = make_gphys(xax, zax) # make gphys from axis	
	gp_phi = to_rad_if_deg(gp_x)      # deg => rad (unit convesion)
	gp_aphi = @@radius * gp_phi       # radius * phi
	# check zax units, if proportional to z or p
	if ( gp_z.data.units =~ Units.new('Pa') )
	  was_proportional_to_p = true 
	elsif ( gp_z.data.units =~ Units.new('m') )
	  was_proportional_to_p = false
	else
	  raise ArgumentError,'unit of zax #{gp_z.data.units} must be 
                         compatible to length or pressure.'
	end
  gp_z = to_z_if_pressure(gp_z)   # convert to z if gp_z is pressure
	gp_z[0] = +1E-6 if gp_z.data.val[0] == -0.0
	return gp_aphi.data, gp_z.data, was_proportional_to_p
end

.radiusObject



469
470
471
# File 'lib/numru/gphys/ep_flux.rb', line 469

def radius
	@@radius
end

.radius=(a) ⇒ Object



472
473
474
475
# File 'lib/numru/gphys/ep_flux.rb', line 472

def radius=(a)
	@@radius = a
	return nil
end

.remove_0_at_poles(a_cos_lat) ⇒ Object



777
778
779
780
781
782
783
784
785
786
787
788
789
790
# File 'lib/numru/gphys/ep_flux.rb', line 777

def remove_0_at_poles(a_cos_lat)
	eps = 1e-6
	if ( (a_cos_lat.val[0]/@@radius).abs.val < eps )
	  a_cos_lat[0] = (a_cos_lat.val[0] + a_cos_lat.val[1])/2
	end
	if ( (a_cos_lat.val[-1]/@@radius).abs.val < eps )
	  a_cos_lat[-1] = (a_cos_lat.val[-1] + a_cos_lat.val[-2])/2
	end
	if a_cos_lat.min.val <= 0
	  raise "Illegal cos(phi) data. phi must between -pi/2 and +pi/2 " +
          "and aligned in increasing or decreasing order."
	end
	nil
end

.rot_periodObject



476
477
478
# File 'lib/numru/gphys/ep_flux.rb', line 476

def rot_period
	@@rot_period
end

.rot_period=(rp) ⇒ Object



486
487
488
489
# File 'lib/numru/gphys/ep_flux.rb', line 486

def rot_period=(rp)
	@@rot_period = rp
	return nil
end

.scale_heightObject

<<< access to constants method >>> ————————————-



462
463
464
# File 'lib/numru/gphys/ep_flux.rb', line 462

def scale_height 
	@@scale_height 
end

.scale_height=(h) ⇒ Object



465
466
467
468
# File 'lib/numru/gphys/ep_flux.rb', line 465

def scale_height=(h)
	@@scale_height = h 
	return nil
end

.set_constants(scale_height, radius, rot_period, g_forces, p00, cp, gas_const) ⇒ Object



511
512
513
514
515
516
517
518
519
520
# File 'lib/numru/gphys/ep_flux.rb', line 511

def set_constants(scale_height, radius, rot_period, g_forces, p00, cp, gas_const)
	@@scale_height = scale_height
	@@radius       = radius
	@@rot_period   = rot_period
	@@g_forces     = g_forces
	@@p00          = p00
	@@cp           = cp
	@@gas_const    = gas_const
	return nil
end

.set_deriv_method(method_name) ⇒ Object

<<< derivation method >>> ———————————————



528
529
530
531
532
533
534
535
536
537
538
539
540
# File 'lib/numru/gphys/ep_flux.rb', line 528

def set_deriv_method( method_name )
  if Deriv_methods.include?( method_name )
    @@deriv_method = eval <<-EOS
	  Proc.new{|*args|
	    GPhys::Derivative::#{method_name}(*args)
	  }
    EOS
  else
    raise ArgumentError, "Unsupported method: #{method_name}. " +
	    "(Supported are #{Deriv_methods.inspect}.)"
  end
	nil
end

.strm_rmean(gp_v, yzdims = [0,1]) ⇒ Object

Raises:

  • (ArgumentError)


810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
# File 'lib/numru/gphys/ep_flux.rb', line 810

def strm_rmean(gp_v, yzdims=[0,1])

	raise ArgumentError,"yzdims's size (#{yzdims.size}) must be 2." if yzdims.size != 2
	## get axis and name
	lat_dim, z_dim = yzdims             # Index of dims
	ax_lat = gp_v.axis(lat_dim)      # Axis of latitude
	ax_z   = gp_v.axis(z_dim)        # Axis of vertical
  lat_nm, z_nm = ax_lat.pos.name, ax_z.pos.name
	gp_lat, gp_z = make_gphys(ax_lat, ax_z)
	## convert
	gp_lat =   to_rad_if_deg(gp_lat)      # deg. to rad.
	gp_p   =   to_p_if_altitude(gp_z)     # z => p=p00exp(-z/H) (units-based) and "Pa"

	## copy grid 
	grid =    gp_v.grid_copy

	## calculate stream function
	na_v  = gp_v.data.val                  # for integration box 
	int_v = gp_v.data.val.dup.fill!(0.0)   # for integration box 
	pres  = gp_p.data.val.dup
	if pres[0] < pres[-1]
	  int_v[*([true]*z_dim+[0, false])] = 0.5*(na_v[*([true] + [0, false])])*pres[0]
	  1.upto( pres.size-1 ) do |idx|
	    dp = (pres[idx] - pres[idx-1])
	    int_v[*([true]*z_dim+[idx, false])] = \
   0.5 * (na_v[*([true] + [idx-1, false])] + na_v[*([true] + [idx, false])]) * dp \
 + int_v[*([true] + [idx-1, false])]
	  end
	else
	  int_v[*([true]*z_dim+[-1, false])] = 0.5*(na_v[*([true] + [-1, false])])*pres[-1]
	  ( pres.size-2 ).downto(0) do |idx|
	    dp = (pres[idx] - pres[idx+1])
	    int_v[*([true]*z_dim+[idx, false])] = \
 0.5 * (na_v[*([true] + [idx+1, false])] + na_v[*([true] + [idx, false])]) * dp \
 + int_v[*([true] + [idx+1, false])]
	  end
	end
	int_v = VArray.new( int_v, gp_v.data, gp_v.data.name )
	int_v.units = Units.new("Pa.m.s-1")
	gp_int_v   = GPhys.new(grid, int_v)
  strm_rmean = gp_int_v * cos(gp_lat) * 2 * PI * @@radius / @@g_forces

	## change attribute
	strm_rmean.name = "strm_rmean"
	strm_rmean.set_att("long_name", "Residual mean mass stream function")

	## convert with past grid
	return strm_rmean
end

.to_p_if_altitude(gp_z) ⇒ Object



724
725
726
727
728
729
730
731
732
733
734
735
736
737
# File 'lib/numru/gphys/ep_flux.rb', line 724

def to_p_if_altitude(gp_z) 
                      # number in units is not considerd operater as log.
	if ( gp_z.data.units =~ Units.new('m') )
	  h = @@scale_height.convert(gp_z.units)
	  gp_z = @@p00*exp(-gp_z/h)
	  gp_z.data.set_att('long_name', "p").rename!("p")
	elsif ( gp_z.data.units =~ Units.new('Pa') )
	  gp_z = gp_z.convert_units(Units.new("Pa"))
	else
	  raise ArgumentError,"units of gp_z (#{gp_z.data.units}) 
                         must be dimention of pressure or length."
	end
	return gp_z
end

.to_rad_if_deg(gp) ⇒ Object



739
740
741
742
743
744
745
746
747
748
749
750
751
752
# File 'lib/numru/gphys/ep_flux.rb', line 739

def to_rad_if_deg(gp)
	if gp.data.units =~ Units.new("degrees")
	  gp = gp.convert_units(Units.new('rad'))
	  gp.units = Units[""]	  
	  gp
	elsif gp.data.units =~ Units.new('rad') 
	  gp.data = gp.data.copy
	  gp.data.units = Units[""]	  
	  gp
	else
	  raise ArgumentError,"units of gp #{gp.data.units} must be equal to deg or radian."
	end
	return gp
end

.to_theta_if_temperature(gp_t, z, flag_temp_or_theta = true) ⇒ Object



697
698
699
700
701
702
703
704
705
706
# File 'lib/numru/gphys/ep_flux.rb', line 697

def to_theta_if_temperature(gp_t, z, flag_temp_or_theta=true) 
                                         # it is only for z coordinate!!!
	if flag_temp_or_theta
	  gp_un = gp_t.data.units
	  gp_t = gp_t.convert_units(Units.new("K"))
	  gp_t = gp_t*exp((@@gas_const/@@cp)*z/@@scale_height)
	  gp_t.data.set_att('long_name', "Potential Temperature")
	end
	return gp_t
end

.to_w_if_omega(gp, z) ⇒ Object

it is only for z coordinate!!!



678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
# File 'lib/numru/gphys/ep_flux.rb', line 678

def to_w_if_omega(gp, z) # it is only for z coordinate!!!
  gp_units = gp.data.units
	if gp_units =~ Units.new("Pa/s")
	  pr = @@p00*exp(-z/@@scale_height)
	  gp_un = gp_units
	  pr = pr.convert_units(gp_un*Units.new('s'))
	  gp = gp*(-@@scale_height/pr)
	  gp.data.rename!("wwnd")
	  gp.data.set_att('long_name', "log-P vertical wind")
	elsif gp_units =~ Units.new("m/s") 
	  gp = gp.convert_units(Units.new('m/s'))
	else
	  raise ArgumentError,"units of gp.data (#{gp.data.units}) 
                         must be dimention of pressure/time 
                                           or length/time."
	end
	return gp
end

.to_z_if_pressure(gp_z) ⇒ Object



708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
# File 'lib/numru/gphys/ep_flux.rb', line 708

def to_z_if_pressure(gp_z) 
                      # number in units is not considerd operater as log.
	if ( gp_z.data.units =~ Units.new('Pa') )
	  p00 = @@p00.convert(gp_z.units)
    gp_z = -@@scale_height*log(gp_z.to_type(NArray::DFLOAT)/p00)
      # it will be change if GPhys is modified for scalor production
	  gp_z.data.set_att('long_name', "z").rename!("z")
	elsif ( gp_z.data.units =~ Units.new('m') )
	  gp_z = gp_z.convert_units(Units.new("m"))
	else
	  raise ArgumentError,"units of gp_z (#{gp_z.data.units}) 
                         must be dimention of pressure or length."
	end
	return gp_z
end