Package ivs :: Package asteroseismology :: Module displacements
[hide private]
[frames] | no frames]

Module displacements

source code

Calculate the stellar surface displacements due to pulsations.

Functions [hide private]
    Helper functions
 
legendre_(l, m, x)
Legendre polynomial.
source code
 
sph_harm(theta, phi, l=2, m=1)
Spherical harmonic according to Townsend, 2002.
source code
 
dsph_harm_dtheta(theta, phi, l=2, m=1)
Derivative of spherical harmonic wrt colatitude.
source code
 
dsph_harm_dphi(theta, phi, l=2, m=1)
Derivative of spherical harmonic wrt longitude.
source code
 
norm_J(l, m)
Normalisation factor
source code
 
norm_atlp1(l, m, Omega, k)
Omega is actually spin parameter (Omega_rot/omega_freq)
source code
 
norm_atlm1(l, m, Omega, k)
Omega is actually spin parameter (Omega_rot/omega_freq)
source code
    Displacement fields
 
radial(theta, phi, l, m, t)
Radial displacement, see Zima 2006.
source code
 
colatitudinal(theta, phi, l, m, t, Omega, k) source code
 
longitudinal(theta, phi, l, m, t, Omega, k) source code
 
surface(theta, phi, l, m, t, Omega=0.1, k=1., asl=0.2, radius=1.) source code
 
observables(theta, phi, teff, logg, l, m, t, Omega=0.1, k=1., asl=0.2, radius=1., delta_T=0.01+0j, delta_g=0.01+0.5j) source code
Function Details [hide private]

legendre_(l, m, x)

source code 

Legendre polynomial.

Check equation (3) from Townsend, 2002:

>>> ls,x = [0,1,2,3,4,5],cos(linspace(0,pi,100))
>>> check = 0
>>> for l in ls:
...     for m in range(-l,l+1,1):
...         Ppos = legendre_(l,m,x)
...         Pneg = legendre_(l,-m,x)
...         mycheck = Pneg,(-1)**m * factorial(l-m)/factorial(l+m) * Ppos
...         check += sum(abs(mycheck[0]-mycheck[1])>1e-10)
>>> print check
0

sph_harm(theta, phi, l=2, m=1)

source code 

Spherical harmonic according to Townsend, 2002.

This function is memoized: once a spherical harmonic is computed, the result is stored in memory

>>> theta,phi = mgrid[0:pi:20j,0:2*pi:40j]
>>> Ylm20 = sph_harm(theta,phi,2,0)
>>> Ylm21 = sph_harm(theta,phi,2,1)
>>> Ylm22 = sph_harm(theta,phi,2,2)
>>> Ylm2_2 = sph_harm(theta,phi,2,-2)
>>> p = figure()
>>> p = gcf().canvas.set_window_title('test of function <sph_harm>')
>>> p = subplot(411);p = title('l=2,m=0');p = imshow(Ylm20.real,cmap=cm.RdBu)
>>> p = subplot(412);p = title('l=2,m=1');p = imshow(Ylm21.real,cmap=cm.RdBu)
>>> p = subplot(413);p = title('l=2,m=2');p = imshow(Ylm22.real,cmap=cm.RdBu)
>>> p = subplot(414);p = title('l=2,m=-2');p = imshow(Ylm2_2.real,cmap=cm.RdBu)

dsph_harm_dtheta(theta, phi, l=2, m=1)

source code 

Derivative of spherical harmonic wrt colatitude.

Using Y_l^m(theta,phi).

Equation:

   sin(theta)*dY/dtheta = (l*J_{l+1}^m * Y_{l+1}^m - (l+1)*J_l^m * Y_{l-1,m})

E.g.: Phd thesis of Joris De Ridder

dsph_harm_dphi(theta, phi, l=2, m=1)

source code 

Derivative of spherical harmonic wrt longitude.

Using Y_l^m(theta,phi).

Equation:

   dY/dphi = i*m*Y

radial(theta, phi, l, m, t)

source code 

Radial displacement, see Zima 2006.

t in phase units