Fit Radial Velocity

You can use Octofitter to fit radial velocity data, either alone or in combination with other kinds of data. Multiple instruments (up to 10) are supported, as are gaussian processes to model stellar activity..

Note

Radial velocity modelling is supported in Octofitter via the extension package OctofitterRadialVelocity. To install it, run pkg> add OctofitterRadialVelocity

For this example, we will fit the orbit of the planet K2-131 to perform the same fit as in the RadVel Gaussian Process Fitting tutorial.

We will use the following packages:

using Octofitter
using OctofitterRadialVelocity
using PlanetOrbits
using CairoMakie
using PairPlots
using CSV
using DataFrames
using Distributions

We will start by downloading and preparing a table of radial velocity measurements, and create a StarAbsoluteRVLikelihood object to hold them.

The following functions allow you to directly load data from various public RV databases:

  • HARPS_DR1_rvs("star-name")
  • HARPS_RVBank_observations("star-name")
  • Lick_rvs("star-name")
  • HIRES_rvs("star-name")

Basic Fit

Start by fitting a 1-planet Keplerian model.

Calling those functions with the name of a star will return a StarAbsoluteRVLikelihood table.

If you would like to manually specify RV data, use the following format:

rv_like = StarAbsoluteRVLikelihood(
    # epoch is in units of MJD. `jd` is a helper function to convert.
    # you can also put `years2mjd(2016.1231)`.
    # rv and σ_rv are in units of meters/second
    (;inst_idx=1, epoch=jd(2455110.97985),  rv=−6.54, σ_rv=1.30),
    (;inst_idx=1, epoch=jd(2455171.90825),  rv=−3.33, σ_rv=1.09),
    # ... etc.
)

For this example, to replicate the results of RadVel, we will download their example data for K2-131. and format it for Octofitter:

rv_file = download("https://raw.githubusercontent.com/California-Planet-Search/radvel/master/example_data/k2-131.txt")
rv_dat_raw = CSV.read(rv_file, DataFrame, delim=' ')
rv_dat = DataFrame();
rv_dat.epoch = jd2mjd.(rv_dat_raw.time)
rv_dat.rv = rv_dat_raw.mnvel
rv_dat.σ_rv = rv_dat_raw.errvel
tels = sort(unique(rv_dat_raw.tel))
rv_dat.inst_idx = map(rv_dat_raw.tel) do tel
    findfirst(==(tel), tels)
end
rvlike = StarAbsoluteRVLikelihood(
    rv_dat,
    instrument_names=["harps-n", "psf"],
)
StarAbsoluteRVLikelihood Table with 4 columns and 70×1 rows:
        epoch    rv        σ_rv   inst_idx
      ┌───────────────────────────────────
 1,1  │ 57782.2  -6682.78  3.71   1
 2,1  │ 57783.1  -6710.55  5.95   1
 3,1  │ 57783.2  -6698.96  8.76   1
 4,1  │ 57812.1  -6672.32  4.0    1
 5,1  │ 57812.2  -6672.0   4.22   1
 6,1  │ 57812.2  -6690.31  4.63   1
 7,1  │ 57813.0  -6701.27  4.12   1
 8,1  │ 57813.1  -6700.71  5.34   1
 9,1  │ 57813.1  -6695.97  3.67   1
 10,1 │ 57813.1  -6705.87  4.23   1
 11,1 │ 57813.1  -6694.27  5.12   1
 12,1 │ 57813.2  -6699.99  5.0    1
 13,1 │ 57813.2  -6696.17  10.43  1
 14,1 │ 57836.0  -6675.4   7.03   1
 15,1 │ 57836.0  -6665.92  7.48   1
 16,1 │ 57836.0  -6661.9   6.03   1
 17,1 │ 57836.1  -6657.92  5.08   1
  ⋮   │    ⋮        ⋮        ⋮       ⋮

Now, create a planet. We can use the RadialVelocityOrbit type from PlanetOrbits.jl that requires fewer parameters (eg no inclination or longitude of ascending node). We could instead use a Visual{KepOrbit} or similar if we wanted to include these parameters and visualize the orbit in the plane of the sky.

@planet b RadialVelocityOrbit begin
    e = 0
    ω = 0.0

    # To match RadVel, we set a prior on Period and calculate semi-major axis from it
    P ~ truncated(Normal(0.3693038/365.256360417, 0.0000091/365.256360417),lower=0.00001)
    a = cbrt(system.M * b.P^2) # note the equals sign.

    τ ~ UniformCircular(1.0)
    tp = b.τ*b.P*365.256360417 + 57782 # reference epoch for τ. Choose an MJD date near your data.

    # minimum planet mass [jupiter masses]. really m*sin(i)
    mass ~ LogUniform(0.001, 10)
end


@system k2_132 begin
    # total mass [solar masses]
    M ~ truncated(Normal(0.82, 0.02),lower=0) # (Baines & Armstrong 2011).

    # HARPS-N
    rv0_1 ~ Normal(0,10000) # m/s
    jitter_1 ~ LogUniform(0.01,100) # m/s

    # FPS
    rv0_2 ~ Normal(0,10000) # m/s
    jitter_2 ~ LogUniform(0.01,100) # m/s
end rvlike b
System model k2_132
Derived:
  
Priors:
  M	Truncated(Distributions.Normal{Float64}(μ=0.82, σ=0.02); lower=0.0)
  rv0_1	Distributions.Normal{Float64}(μ=0.0, σ=10000.0)
  jitter_1	Distributions.LogUniform{Float64}(a=0.01, b=100.0)
  rv0_2	Distributions.Normal{Float64}(μ=0.0, σ=10000.0)
  jitter_2	Distributions.LogUniform{Float64}(a=0.01, b=100.0)
Planet b
Derived:
  e, ω, a, τ, tp, 
Priors:
  P	Truncated(Distributions.Normal{Float64}(μ=0.001011081092683449, σ=2.4914008313533153e-8); lower=1.0e-5)
  τy	Distributions.Normal{Float64}(μ=0.0, σ=1.0)
  τx	Distributions.Normal{Float64}(μ=0.0, σ=1.0)
  mass	Distributions.LogUniform{Float64}(a=0.001, b=10.0)
Octofitter.UnitLengthPrior{:τy, :τx}: √(τy^2+τx^2) ~ LogNormal(log(1), 0.02)


StarAbsoluteRVLikelihood Table with 4 columns and 70×1 rows:
        epoch    rv        σ_rv   inst_idx
      ┌───────────────────────────────────
 1,1  │ 57782.2  -6682.78  3.71   1
 2,1  │ 57783.1  -6710.55  5.95   1
 3,1  │ 57783.2  -6698.96  8.76   1
 4,1  │ 57812.1  -6672.32  4.0    1
 5,1  │ 57812.2  -6672.0   4.22   1
 6,1  │ 57812.2  -6690.31  4.63   1
 7,1  │ 57813.0  -6701.27  4.12   1
 8,1  │ 57813.1  -6700.71  5.34   1
 9,1  │ 57813.1  -6695.97  3.67   1
 10,1 │ 57813.1  -6705.87  4.23   1
 11,1 │ 57813.1  -6694.27  5.12   1
 12,1 │ 57813.2  -6699.99  5.0    1
 13,1 │ 57813.2  -6696.17  10.43  1
 14,1 │ 57836.0  -6675.4   7.03   1
 15,1 │ 57836.0  -6665.92  7.48   1
 16,1 │ 57836.0  -6661.9   6.03   1
 17,1 │ 57836.1  -6657.92  5.08   1
  ⋮   │    ⋮        ⋮        ⋮       ⋮

Note how the rvlike object was attached to the k2_132 system instead of the planet. This is because the observed radial velocity is of the star, and is caused by any/all orbiting planets.

The rv0 and jitter parameters specify priors for the instrument-specific offset and white noise jitter standard deviation. The _i index matches the inst_idx used to create the observation table.

Note also here that the mass variable is really msini, or the minimum mass of the planet.

We can now prepare our model for sampling.

model = Octofitter.LogDensityModel(k2_132)
LogDensityModel for System k2_132 of dimension 9 with fields .ℓπcallback and .∇ℓπcallback

Sample:

using Random
rng = Random.Xoshiro(0)

chain = octofit(rng, model)
Chains MCMC chain (1000×27×1 Array{Float64, 3}):

Iterations        = 1:1:1000
Number of chains  = 1
Samples per chain = 1000
Wall duration     = 28.07 seconds
Compute duration  = 28.07 seconds
parameters        = M, rv0_1, jitter_1, rv0_2, jitter_2, b_P, b_τy, b_τx, b_mass, b_e, b_ω, b_a, b_τ, b_tp
internals         = n_steps, is_accept, acceptance_rate, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size, is_adapt, loglike, logpost, tree_depth, numerical_error

Summary Statistics
  parameters         mean       std      mcse    ess_bulk   ess_tail      rhat ⋯
      Symbol      Float64   Float64   Float64     Float64    Float64   Float64 ⋯

           M       0.8192    0.0202    0.0006   1125.1426   489.3575    1.0021 ⋯
       rv0_1   -6693.4505    2.2889    0.0812    805.2573   651.2369    0.9995 ⋯
    jitter_1      13.7535    1.9467    0.0683    881.4580   567.3051    0.9992 ⋯
       rv0_2       1.4364    3.4984    0.1319    705.3459   472.7397    1.0139 ⋯
    jitter_2      17.8337    2.9363    0.1122    759.9206   520.5404    1.0019 ⋯
         b_P       0.0010    0.0000    0.0000   1084.4514   691.7446    1.0012 ⋯
        b_τy       0.7520    0.1615    0.0066    732.9331   376.7255    1.0001 ⋯
        b_τx      -0.6232    0.1906    0.0066    866.5706   623.8066    0.9996 ⋯
      b_mass       0.0419    0.0087    0.0003    683.6193   362.7727    0.9991 ⋯
         b_e       0.0000    0.0000       NaN         NaN        NaN       NaN ⋯
         b_ω       0.0000    0.0000       NaN         NaN        NaN       NaN ⋯
         b_a       0.0094    0.0001    0.0000   1126.2965   489.3575    1.0022 ⋯
         b_τ      -0.1099    0.0373    0.0015    632.9820   423.8635    0.9992 ⋯
        b_tp   57781.9594    0.0138    0.0006    633.0129   423.8635    0.9992 ⋯
                                                                1 column omitted

Quantiles
  parameters         2.5%        25.0%        50.0%        75.0%        97.5%
      Symbol      Float64      Float64      Float64      Float64      Float64

           M       0.7812       0.8047       0.8190       0.8334       0.8560
       rv0_1   -6698.1343   -6694.9449   -6693.3785   -6691.9918   -6689.1428
    jitter_1      10.6321      12.4279      13.5143      14.9226      18.2011
       rv0_2      -5.7378      -0.8855       1.3597       3.7745       8.1148
    jitter_2      13.0378      15.7607      17.5097      19.6121      24.0961
         b_P       0.0010       0.0010       0.0010       0.0010       0.0010
        b_τy       0.4018       0.6615       0.7614       0.8646       1.0262
        b_τx      -0.9751      -0.7507      -0.6261      -0.5025      -0.2279
      b_mass       0.0236       0.0363       0.0423       0.0478       0.0580
         b_e       0.0000       0.0000       0.0000       0.0000       0.0000
         b_ω       0.0000       0.0000       0.0000       0.0000       0.0000
         b_a       0.0093       0.0094       0.0094       0.0095       0.0096
         b_τ      -0.1838      -0.1338      -0.1101      -0.0860      -0.0365
        b_tp   57781.9321   57781.9506   57781.9594   57781.9683   57781.9865

Excellent! Let's plot the maximum likelihood orbit:

using CairoMakie: Makie
fig = OctofitterRadialVelocity.rvpostplot(model, chain) # saved to "k2_132-rvpostplot.png"
Example block output

Create a corner plot:

using PairPlots
using CairoMakie: Makie
octocorner(model, chain, small=true)
Example block output

Gaussian Process Fit

Let us now add a Gaussian process to model stellar activity. This should improve the fit.

We start by writing a functino that creates a Gaussian process kernel from a set of system parameters. We will create a quasi-periodic kernel.

using AbstractGPs

function gauss_proc_quasi_periodic(θ_system)
    # θ_system is a named tuple of variable values for this posterior draw.
    η₁ = θ_system.gp_η₁ # h
    η₂ = θ_system.gp_η₂ # λ
    η₃ = θ_system.gp_η₃ # θ
    η₄ = θ_system.gp_η₄ # ω
    kernel = η₁^2 *
        (SqExponentialKernel() ∘ ScaleTransform(1/(η₂))) *
        (PeriodicKernel(r=[η₄]) ∘ ScaleTransform(1/(η₃)))
    gp = GP(kernel)
    return gp
end
gauss_proc_quasi_periodic (generic function with 1 method)

Now provide this function when we construct our likelihood object:

rvlike = StarAbsoluteRVLikelihood(rv_dat,
    instrument_names=["harps-n", "psf"],
    gaussian_process = gauss_proc_quasi_periodic
)
StarAbsoluteRVLikelihood Table with 4 columns and 70×1 rows:
        epoch    rv        σ_rv   inst_idx
      ┌───────────────────────────────────
 1,1  │ 57782.2  -6682.78  3.71   1
 2,1  │ 57783.1  -6710.55  5.95   1
 3,1  │ 57783.2  -6698.96  8.76   1
 4,1  │ 57812.1  -6672.32  4.0    1
 5,1  │ 57812.2  -6672.0   4.22   1
 6,1  │ 57812.2  -6690.31  4.63   1
 7,1  │ 57813.0  -6701.27  4.12   1
 8,1  │ 57813.1  -6700.71  5.34   1
 9,1  │ 57813.1  -6695.97  3.67   1
 10,1 │ 57813.1  -6705.87  4.23   1
 11,1 │ 57813.1  -6694.27  5.12   1
 12,1 │ 57813.2  -6699.99  5.0    1
 13,1 │ 57813.2  -6696.17  10.43  1
 14,1 │ 57836.0  -6675.4   7.03   1
 15,1 │ 57836.0  -6665.92  7.48   1
 16,1 │ 57836.0  -6661.9   6.03   1
 17,1 │ 57836.1  -6657.92  5.08   1
  ⋮   │    ⋮        ⋮        ⋮       ⋮

We create a new system model with added parameters gp_η₁ etc for the Gaussian process fit.

gp_explength_mean = 9.5*sqrt(2.) # sqrt(2)*tau in Dai+ 2017 [days]
gp_explength_unc = 1.0*sqrt(2.)
gp_perlength_mean = sqrt(1. /(2. *3.32)) # sqrt(1/(2*gamma)) in Dai+ 2017
gp_perlength_unc = 0.019
gp_per_mean = 9.64 # T_bar in Dai+ 2017 [days]
gp_per_unc = 0.12

# No change to planet model
@planet b RadialVelocityOrbit begin
    e = 0
    ω = 0.0
    P ~ truncated(Normal(0.3693038/365.25, 0.0000091/365.25),lower=0.00001)
    a = cbrt(system.M * b.P^2)
    mass ~ LogUniform(0.001, 10)

    τ ~ UniformCircular(1.0)
    tp =  b.τ*b.P*365.25 + 57782 # reference epoch for τ. Choose an MJD date near your data.
end


@system k2_132 begin
    M ~ truncated(Normal(0.82, 0.02),lower=0) # (Baines & Armstrong 2011).

    # HARPS-N
    rv0_1 ~ Normal(0,10000)
    jitter_1 ~ LogUniform(0.01,10)
    # FPS
    rv0_2 ~ Normal(0,10000)
    jitter_2 ~ LogUniform(0.01,10)

    # Add priors on GP kernel hyper-parameters.
    gp_η₁ ~ truncated(Normal(25,10),lower=0)
    gp_η₂ ~ truncated(Normal(gp_explength_mean,gp_explength_unc),lower=0)
    gp_η₃ ~ truncated(Normal(gp_per_mean,1),lower=0)
    gp_η₄ ~ truncated(Normal(gp_perlength_mean,gp_perlength_unc),lower=0)

end rvlike b

model = Octofitter.LogDensityModel(k2_132; autodiff=:ForwardDiff)
LogDensityModel for System k2_132 of dimension 13 with fields .ℓπcallback and .∇ℓπcallback

Sample from the model as before:

using Random
rng = Random.Xoshiro(0)

chain = octofit(
    rng, model,
    adaptation = 100,
    iterations = 100,
)
Chains MCMC chain (100×31×1 Array{Float64, 3}):

Iterations        = 1:1:100
Number of chains  = 1
Samples per chain = 100
Wall duration     = 1907.53 seconds
Compute duration  = 1907.53 seconds
parameters        = M, rv0_1, jitter_1, rv0_2, jitter_2, gp_η₁, gp_η₂, gp_η₃, gp_η₄, b_P, b_mass, b_τy, b_τx, b_e, b_ω, b_a, b_τ, b_tp
internals         = n_steps, is_accept, acceptance_rate, hamiltonian_energy, hamiltonian_energy_error, max_hamiltonian_energy_error, tree_depth, numerical_error, step_size, nom_step_size, is_adapt, loglike, logpost, tree_depth, numerical_error

Summary Statistics
  parameters         mean       std      mcse   ess_bulk   ess_tail      rhat  ⋯
      Symbol      Float64   Float64   Float64    Float64    Float64   Float64  ⋯

           M       0.8208    0.0194    0.0022    74.9873    34.3628    1.0306  ⋯
       rv0_1   -6563.4520    0.4181    0.3630     1.4059    18.3939    2.1313  ⋯
    jitter_1       0.0687    0.0136    0.0091     2.2833    17.2397    1.4087  ⋯
       rv0_2    -116.9768    0.2423    0.1871     1.5713    12.8043    1.8356  ⋯
    jitter_2       1.7523    0.9953    0.9110     1.4188    21.5943    2.1313  ⋯
       gp_η₁      65.6337    5.5222    1.0836    28.4422    66.9638    1.0285  ⋯
       gp_η₂      14.6994    1.1476    0.1764    43.6140    51.5563    0.9985  ⋯
       gp_η₃       9.5627    1.0615    0.4806     5.5638    22.4187    1.1405  ⋯
       gp_η₄       0.4083    0.0200    0.0019   123.2051    78.3393    0.9968  ⋯
         b_P       0.0010    0.0000    0.0000   107.2863    51.5338    1.0019  ⋯
      b_mass       0.0143    0.0044    0.0040     1.4444    40.5549    2.0821  ⋯
        b_τy       0.6837    0.1415    0.0387    14.2879    21.0246    1.1262  ⋯
        b_τx      -0.6515    0.1372    0.0340    17.4281    20.6246    1.0849  ⋯
         b_e       0.0000    0.0000       NaN        NaN        NaN       NaN  ⋯
         b_ω       0.0000    0.0000       NaN        NaN        NaN       NaN  ⋯
         b_a       0.0094    0.0001    0.0000    74.2250    34.3628    1.0308  ⋯
         b_τ      -0.1214    0.0297    0.0075    16.7500    14.7365    1.2178  ⋯
        b_tp   57781.9552    0.0110    0.0028    16.7500    14.7365    1.2262  ⋯
                                                                1 column omitted

Quantiles
  parameters         2.5%        25.0%        50.0%        75.0%        97.5%
      Symbol      Float64      Float64      Float64      Float64      Float64

           M       0.7897       0.8065       0.8215       0.8354       0.8584
       rv0_1   -6563.9381   -6563.8013   -6563.6235   -6563.1571   -6562.6608
    jitter_1       0.0511       0.0578       0.0667       0.0769       0.1009
       rv0_2    -117.4854    -117.1010    -116.9954    -116.8309    -116.4809
    jitter_2       0.5133       0.7950       1.6374       2.6769       3.4244
       gp_η₁      57.3208      61.4333      65.6576      69.0235      76.9867
       gp_η₂      12.4763      13.9981      14.7889      15.3610      16.7696
       gp_η₃       7.5877       8.6108       9.8695      10.1704      11.2740
       gp_η₄       0.3722       0.3934       0.4080       0.4232       0.4466
         b_P       0.0010       0.0010       0.0010       0.0010       0.0010
      b_mass       0.0085       0.0100       0.0135       0.0188       0.0215
        b_τy       0.3666       0.6108       0.6922       0.7811       0.9181
        b_τx      -0.8858      -0.7315      -0.6618      -0.5714      -0.3875
         b_e       0.0000       0.0000       0.0000       0.0000       0.0000
         b_ω       0.0000       0.0000       0.0000       0.0000       0.0000
         b_a       0.0093       0.0094       0.0094       0.0095       0.0096
         b_τ      -0.1878      -0.1378      -0.1226      -0.1062      -0.0640
        b_tp   57781.9307   57781.9491   57781.9547   57781.9608   57781.9764

For real data, we would want to increase the adaptation and iterations to about 1000 each.

Plot the results:

fig = OctofitterRadialVelocity.rvpostplot(model, chain) # saved to "k2_132-rvpostplot.png"
Example block output

Corner plot:

octocorner(model, chain, small=true) # saved to "k2_132-pairplot-small.png"
Example block output

It is expensive (per iteration) to fit Guassian processes. If you have many cores (or a cluster) available you might also try starting julia with multiple threads (julia --threads=auto) and using octofit_pigeons. Pigeons is less efficient with single-core, but has better parallelization scaling:

using Pigeons
chain, pt = octofit_pigeons(
    model,
    n_rounds=9
)
display(chain)