In [33]:
import numpy as np
import matplotlib.pyplot as plt
%matplotlib inline
from scipy.optimize import fmin_cg, minimize
import h5py
c = 2.99792458e8   # m/s

In [34]:
def doppler(v):
    frac = (1. - v/c) / (1. + v/c)
    return np.sqrt(frac)

def gamma(v):
    return 1. / np.sqrt(1. - (v/c) ** 2)

def dlndopplerdv(v):
    dv = doppler(v)
    return -1. * gamma(v) / (2. * c) * (1. / dv  + dv)

def state(v, xs, xps):
    '''
    outputs: (M, Mp, v, xs, ms, mps, ehs, bees, seas)
    M and Mp are the lengths of data and model wavelength grids
    v is the RV
    xs is the wavelength values of the data grid
    ms is the data index m at which there is an interpolated model value
    mps is the model index m' from which we interpolate to ms
    ehs, bees, and seas go into the coefficients of interpolation
    '''
    # every input must be 1-d
    M = len(xs)
    Mp = len(xps)
    xps_shifted = xps + np.log(doppler(v))
    ms = np.arange(M)
    mps = np.searchsorted(xps_shifted, xs, side='left')
    good = (mps > 0) * (mps < Mp)
    ms = ms[good]
    mps = mps[good]
    ehs = xps_shifted[mps] - xs[ms]
    bees = xs[ms] - xps_shifted[mps - 1]
    seas = ehs + bees
    return (M, Mp, v, xs, ms, mps, ehs, bees, seas)

def Pdot(state, vec):
    # takes state and model flux vector, returns (shifted) model interpolated into data space
    # unpack state
    M, Mp, v, xs, ms, mps, ehs, bees, seas = state
    # do shit
    result = np.zeros(M)
    result[ms] = vec[mps - 1] * ehs / seas + vec[mps] * bees / seas
    return result

def dotP(state, vec):
    # takes state and data flux vector, returns data interpolated into (shifted) model space
    # unpack state
    M, Mp, v, xs, ms, mps, ehs, bees, seas = state
    # do shit
    result = np.zeros(Mp)
    result[mps - 1] += vec[ms] * ehs / seas
    result[mps] += vec[ms] * bees / seas
    return result

def dotdPdv(state, vec):
    # unpack state
    M, Mp, v, xs, ms, mps, ehs, bees, seas = state
    # do shit
    result = np.zeros(Mp)
    foos = vec[ms] / seas * dlndopplerdv(v) # * xs[ms] ??
    result[mps - 1] += foos
    result[mps] -= foos
    return result

def dPdotdv(state, vec):
    # unpack state
    M, Mp, v, xs, ms, mps, ehs, bees, seas = state
    # do shit
    result = np.zeros(M)
    result[ms] = (vec[mps - 1] - vec[mps]) * dlndopplerdv(v) / seas
    return result

In [35]:
with h5py.File('../data/hip30037.hdf5', 'r') as f:
    N = 16
    inds = (f['xs'][:] > 5900.) & (f['xs'][:] < 6000.)
    data = np.copy(f['data'])[:N,inds]
    data_xs = np.log(np.copy(f['xs'][inds]))
    ivars = np.copy(f['ivars'])[:N,inds]
    true_rvs = np.copy(f['true_rvs'])[:N]
    bervs = np.copy(f['berv'])[:N] * -1.e3

for i in xrange(len(data)):
    data[i] /= np.median(data[i])
    
data = np.log(data)

In [36]:
plt.scatter(np.arange(N),true_rvs - bervs)
print true_rvs


[ 50341.31574258  50361.43008207  51079.54982724  50977.90950712
  50827.16467053  47564.37201857  47214.5807874   46409.68404177
  48131.20385621  47776.03953678  47391.24726139  46978.18009382
  40937.59474406  40744.2769875   41153.4526148   41892.27624778]

In [37]:
plt.plot(data_xs, data[0,:])
plt.plot(data_xs, data[10,:])
plt.ylim([-0.05, 0.05])


Out[37]:
(-0.05, 0.05)

In [38]:
def make_template(all_data, rvs, xs, dx):
    """
    `all_data`: `[N, M]` array of pixels
    `rvs`: `[N]` array of RVs
    `xs`: `[M]` array of wavelength values
    `dx`: linear spacing desired for template wavelength grid (A)
    """
    (N,M) = np.shape(all_data)
    all_xs = np.empty_like(all_data)
    for i in range(N):
        all_xs[i,:] = xs - np.log(doppler(rvs[i])) # shift to rest frame
    all_data, all_xs = np.ravel(all_data), np.ravel(all_xs)
    tiny = 10.
    template_xs = np.arange(min(all_xs)-tiny*dx, max(all_xs)+tiny*dx, dx)
    template_ys = np.nan + np.zeros_like(template_xs)
    for i,t in enumerate(template_xs):
        ind = (all_xs >= t-dx/2.) & (all_xs < t+dx/2.)
        if np.sum(ind) > 0:
            template_ys[i] = np.nanmedian(all_data[ind])
    ind_nan = np.isnan(template_ys)
    template_ys.flat[ind_nan] = np.interp(template_xs[ind_nan], template_xs[~ind_nan], template_ys[~ind_nan]) #np.interp(template_xs[ind_nan], template_xs[~ind_nan], template_ys[~ind_nan])
    return template_xs, template_ys

def subtract_template(data_xs, data, model_xs_t, model_ys_t, rvs_t):
    (N,M) = np.shape(data)
    data_sub = np.copy(data)
    for n,v in enumerate(rvs_t):
        s = state(v, data_xs, model_xs_t)
        model_ys_t_shifted = Pdot(s, model_ys_t)
        data_sub[n,:] -= np.ravel(model_ys_t_shifted)
        if n == 0:
            plt.plot(data_xs, data[n,:], color='k')
            plt.plot(data_xs, data_sub[n,:], color='blue')
            plt.plot(data_xs, np.ravel(model_ys_t_shifted), color='red')
    return data_sub

In [39]:
#x0_star = -np.copy(bervs)
x0_star = -np.copy(true_rvs)
x0_star -= np.mean(x0_star)
x0_t = np.zeros(N)
model_xs_star, model_ys_star = make_template(data, x0_star, data_xs, np.log(6000.01) - np.log(6000.))
model_xs_t, model_ys_t = make_template(data, x0_t, data_xs, np.log(6000.01) - np.log(6000.))

In [40]:
plt.plot(model_xs_star, model_ys_star, color='k', alpha=0.5)
plt.plot(model_xs_t, model_ys_t, color='red', alpha=0.5)
plt.plot(data_xs, data[0,:], color='blue', alpha=0.5)
plt.ylim([-0.1, 0.1])


Out[40]:
(-0.1, 0.1)

In [41]:
def rv_lnprior(rvs):
    return -0.5 * np.mean(rvs)**2/1.**2

def drv_lnprior_dv(rvs):
    return np.zeros_like(rvs) - np.mean(rvs)/1.**2/len(rvs)

def lnlike_star(rvs_star, rvs_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, model_ys_t):
    try:
        N = len(rvs_star)
    except:
        N = 1  
    lnlike = 0.
    dlnlike_dv = np.zeros(N)
    for n in range(N):
        state_star = state(rvs_star[n], data_xs, model_xs_star)
        pd_star = Pdot(state_star, model_ys_star)
        state_t = state(rvs_t[n], data_xs, model_xs_t)
        pd_t = Pdot(state_t, model_ys_t)
        pd = pd_star + pd_t
        lnlike += -0.5 * np.sum((data[n,:] - pd)**2 * ivars[n,:])
        dpd_dv = dPdotdv(state_star, model_ys_star)
        dlnlike_dv[n] = np.sum((data[n,:] - pd) * ivars[n,:] * dpd_dv)
    lnpost = lnlike + rv_lnprior(rvs_star)
    dlnpost_dv = dlnlike_dv + drv_lnprior_dv(rvs_star)
    return -1 * lnpost, -1 * dlnpost_dv

def lnlike_t(rvs_t, rvs_star, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, model_ys_t):
    try:
        N = len(rvs_t)
    except:
        N = 1  
    lnlike = 0.
    dlnlike_dv = np.zeros(N)
    for n in range(N):
        state_star = state(rvs_star[n], data_xs, model_xs_star)
        pd_star = Pdot(state_star, model_ys_star)
        state_t = state(rvs_t[n], data_xs, model_xs_t)
        pd_t = Pdot(state_t, model_ys_t)
        pd = pd_star + pd_t
        lnlike += -0.5 * np.sum((data[n,:] - pd)**2 * ivars[n,:])
        dpd_dv = dPdotdv(state_t, model_ys_t)
        dlnlike_dv[n] = np.sum((data[n,:] - pd) * ivars[n,:] * dpd_dv)
    lnpost = lnlike + rv_lnprior(rvs_t)
    dlnpost_dv = dlnlike_dv + drv_lnprior_dv(rvs_t)
    return -1 * lnpost, -1 * dlnpost_dv

In [42]:
def model_ys_lnprior(w):
    return -0.5 * np.sum(w**2)/100.**2


def dmodel_ys_lnprior_dw(w):
    return -1.*w / 100.**2


def dlnlike_star_dw_star(model_ys_star, rvs_star, rvs_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t):
    try:
        N = len(rvs_star)
    except:
        N = 1  
    lnlike = 0.
    Mp = len(model_xs_star)
    dlnlike_dw = np.zeros(Mp)
    for n in range(N):
        state_star = state(rvs_star[n], data_xs, model_xs_star)
        pd_star = Pdot(state_star, model_ys_star)
        state_t = state(rvs_t[n], data_xs, model_xs_t)
        pd_t = Pdot(state_t, model_ys_t)
        pd = pd_star + pd_t
        dp_star = dotP(state_star, (data[n,:] - pd)*ivars[n,:]) 
        lnlike += -0.5 * np.sum((data[n,:] - pd)**2 * ivars[n,:])
        dlnlike_dw += dp_star
    lnprior = model_ys_lnprior(model_ys_star)
    dlnprior = dmodel_ys_lnprior_dw(model_ys_star)
    return -lnlike - lnprior, -dlnlike_dw - dlnprior

def dlnlike_t_dw_t(model_ys_t, rvs_star, rvs_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t):
    try:
        N = len(rvs_t)
    except:
        N = 1  
    lnlike = 0.
    Mp = len(model_xs_t)
    dlnlike_dw = np.zeros(Mp)
    for n in range(N):
        state_star = state(rvs_star[n], data_xs, model_xs_star)
        pd_star = Pdot(state_star, model_ys_star)
        state_t = state(rvs_t[n], data_xs, model_xs_t)
        pd_t = Pdot(state_t, model_ys_t)
        pd = pd_star + pd_t
        dp_t = dotP(state_t, (data[n,:] - pd)*ivars[n,:]) 
        lnlike += -0.5 * np.sum((data[n,:] - pd)**2 * ivars[n,:])
        dlnlike_dw += dp_t
    lnprior = model_ys_lnprior(model_ys_t)
    dlnprior = dmodel_ys_lnprior_dw(model_ys_t)
    return -lnlike - lnprior, -dlnlike_dw - dlnprior

In [43]:
def improve_telluric_model(model_ys_t, x0_star, x0_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, step_scale=5e-7):
    w = np.copy(model_ys_t)
    for i in range(50):
        lnlike, dlnlike_dw = dlnlike_t_dw_t(w, x0_star, x0_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t)
        w -= step_scale * dlnlike_dw    
    return w

def improve_star_model(model_ys_star, x0_star, x0_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t, step_scale=5e-7):
    w = np.copy(model_ys_star)
    for i in range(50):
        lnlike, dlnlike_dw = dlnlike_star_dw_star(w, x0_star, x0_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t)
        w -= step_scale * dlnlike_dw 
    return w

In [44]:
for n in range(5):
    print "Fitting stellar RVs..."
    soln_star =  minimize(lnlike_star, x0_star, args=(x0_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, model_ys_t),
                     method='BFGS', jac=True, options={'disp':True, 'gtol':1.e-2, 'eps':1.5e-5})['x']
    
    model_ys_t = improve_telluric_model(model_ys_t, soln_star, x0_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t)
    model_ys_star = improve_star_model(model_ys_star, soln_star, x0_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t)
    print "Star model improved. Fitting telluric RVs..."
    soln_t =  minimize(lnlike_t, x0_t, args=(soln_star, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, model_ys_t),
                     method='BFGS', jac=True, options={'disp':True, 'gtol':1.e-2, 'eps':1.5e-5})['x']
    
    model_ys_t = improve_telluric_model(model_ys_t, soln_star, soln_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t)
    model_ys_star = improve_star_model(model_ys_star, soln_star, soln_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t)
    
    x0_star = soln_star
    x0_t = soln_t

    print "iter {0}: star std = {1:.2f}, telluric std = {2:.2f}".format(n, np.std(soln_star + true_rvs), np.std(soln_t))
    plt.plot(np.arange(N), soln_star + true_rvs - np.mean(soln_star + true_rvs), color='k')
    plt.plot(np.arange(N), soln_t - np.mean(soln_t), color='red')
    plt.show()


Fitting stellar RVs...
Warning: Desired error not necessarily achieved due to precision loss.
         Current function value: 1572098.899122
         Iterations: 79
         Function evaluations: 129
         Gradient evaluations: 124
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 191044.601894
         Iterations: 56
         Function evaluations: 68
         Gradient evaluations: 68
iter 0: star std = 1162.91, telluric std = 154.63
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 124613.203765
         Iterations: 42
         Function evaluations: 44
         Gradient evaluations: 44
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 95621.874463
         Iterations: 40
         Function evaluations: 43
         Gradient evaluations: 43
iter 1: star std = 508.23, telluric std = 136.63
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 83565.492810
         Iterations: 36
         Function evaluations: 38
         Gradient evaluations: 38
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 76058.326615
         Iterations: 49
         Function evaluations: 54
         Gradient evaluations: 54
iter 2: star std = 270.43, telluric std = 129.04
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 72303.204879
         Iterations: 35
         Function evaluations: 37
         Gradient evaluations: 37
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 69811.222788
         Iterations: 30
         Function evaluations: 37
         Gradient evaluations: 37
iter 3: star std = 153.16, telluric std = 137.53
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 68477.637779
         Iterations: 33
         Function evaluations: 35
         Gradient evaluations: 35
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 67546.589287
         Iterations: 40
         Function evaluations: 45
         Gradient evaluations: 45
iter 4: star std = 93.11, telluric std = 155.17

In [45]:
plt.scatter(np.arange(N), soln_star+true_rvs)


Out[45]:
<matplotlib.collections.PathCollection at 0x1197925d0>

In [46]:
plt.plot(model_xs_star, model_ys_star, 'k', alpha=0.5)
plt.plot(model_xs_t, model_ys_t, 'r', alpha=0.5)
plt.ylim([-0.05, 0.05])


Out[46]:
(-0.05, 0.05)

In [47]:
print len(data_xs)
print len(model_xs_star)
print len(model_xs_t)


10000
10124
10104

In [48]:
for n in range(5):
    print "Fitting stellar RVs..."
    soln_star =  minimize(lnlike_star, x0_star, args=(x0_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, model_ys_t),
                     method='BFGS', jac=True, options={'disp':True, 'gtol':1.e-2, 'eps':1.5e-5})['x']
    
    model_ys_t = improve_telluric_model(model_ys_t, soln_star, x0_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t)
    model_ys_star = improve_star_model(model_ys_star, soln_star, x0_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t)
    print "Star model improved. Fitting telluric RVs..."
    soln_t =  minimize(lnlike_t, x0_t, args=(soln_star, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t, model_ys_t),
                     method='BFGS', jac=True, options={'disp':True, 'gtol':1.e-2, 'eps':1.5e-5})['x']
    
    model_ys_t = improve_telluric_model(model_ys_t, soln_star, soln_t, data_xs, data, ivars, model_xs_star, model_ys_star, model_xs_t)
    model_ys_star = improve_star_model(model_ys_star, soln_star, soln_t, data_xs, data, ivars, model_xs_star, model_xs_t, model_ys_t)
    
    x0_star = soln_star
    x0_t = soln_t

    print "iter {0}: star residual std = {1:.2f}, telluric std = {2:.2f}".format(n, np.std(soln_star + true_rvs), np.std(soln_t))
    plt.plot(np.arange(N), soln_star + true_rvs - np.mean(soln_star + true_rvs), color='k')
    plt.plot(np.arange(N), soln_t - np.mean(soln_t), color='red')
    plt.show()


Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 67003.163672
         Iterations: 29
         Function evaluations: 31
         Gradient evaluations: 31
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 66608.387046
         Iterations: 41
         Function evaluations: 48
         Gradient evaluations: 48
iter 0: star residual std = 62.78, telluric std = 173.01
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 66350.890302
         Iterations: 24
         Function evaluations: 26
         Gradient evaluations: 26
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 66157.012417
         Iterations: 37
         Function evaluations: 49
         Gradient evaluations: 49
iter 1: star residual std = 47.14, telluric std = 187.51
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 66019.922305
         Iterations: 25
         Function evaluations: 27
         Gradient evaluations: 27
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 65911.112708
         Iterations: 33
         Function evaluations: 40
         Gradient evaluations: 40
iter 2: star residual std = 38.61, telluric std = 198.98
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 65827.172719
         Iterations: 19
         Function evaluations: 21
         Gradient evaluations: 21
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 65756.677613
         Iterations: 33
         Function evaluations: 39
         Gradient evaluations: 39
iter 3: star residual std = 33.23, telluric std = 209.31
Fitting stellar RVs...
Optimization terminated successfully.
         Current function value: 65698.754549
         Iterations: 16
         Function evaluations: 19
         Gradient evaluations: 19
Star model improved. Fitting telluric RVs...
Optimization terminated successfully.
         Current function value: 65649.812565
         Iterations: 16
         Function evaluations: 19
         Gradient evaluations: 19
iter 4: star residual std = 29.67, telluric std = 216.29

In [50]:
plt.scatter(np.arange(N), soln_star+true_rvs - np.mean(soln_star+true_rvs))


Out[50]:
<matplotlib.collections.PathCollection at 0x11484c850>

In [63]:
plt.plot(np.arange(N), soln_star - np.mean(soln_star) + bervs, 'ko')
plt.plot(np.arange(N), -true_rvs + np.mean(true_rvs) + bervs, 'r.')


Out[63]:
[<matplotlib.lines.Line2D at 0x114a68d50>]

In [ ]: