Welcome to the transmission spectrum model/retrieval tutorial!
For this particular setup, the atmosphere is parameterized within the classic "free retrieval" framework. In the free retrieval, typically, molecular/atomic abundances are assumed to be constant with altitude and chemically un-related (e.g., any abundance combinations are allowed). In principle, we can retrieve anything with an opacity (cross-sections, CK-coefficients). This tuotiral/code has these opacities for H2O, CH4, CO, CO2, NH3, HCN, H2S, PH3 C2H2, C2H6, Na, K, TiO, VO, FeH, H, H2, He, e-, H- either through "line" or "continuum" opacities (e.g., H- free-free dependes on both the e- and H abundances). Here we will assume any remaining gas is H2+He (at solar He/H2 ratio). This may not be a good assumption for ultra hot Jupiters, where H becomes dominent, or for very high metallicity objects where other species like H2O, CO2, N2, or CO become major background gases. Free retrievals tend to work best in the "run-of-the-mill-hot-but-not-too-hot-Jupiter" realm. Be weary, however, as free retrievals can trick one into a good fit with seemingly unphysical solutions. None-the-less, they do provide some insight.
For the temperature profile, here, we assume the 3-parameter temperature Guillot 2010/Parmentier et al. 2014 analytic formulism (see Line et al. 2013a for implementation details).
The transmission spectrum routine closely follows the equations (and figure) in Tinetti et al. 2012 (as described in the tutorial text). Instead of using line-by-line, or "sampled" cross-sections, this implementation uses the "correlated-K" method (see Lacis & Oinas 1990, or more recently Amundsen et al. 2017). Correlated-K is advantageous as it preserves the wavelength bin"integrated" precision as line-by-line but with far less demanding computation. See the "opacity" tutorial for more details on how correlated-K works.
There are two different cloud prescriptions built in. The first is the Ackerman & Marley 2001 "eddy-sed" approach that self-consistently computes the vertical particle size distribution given a sedimentation factor, $f_{sed}$ and an eddy mixing factor (K$_{zz}$) from some cloud base pressure and intrinsic condensate mixing ratio. The classic "power-law haze" and "grey cloud" prescripton is also included.
This specific notebook goes through the steps to generate the forward model, and illustrate how to actually perform the free retrieval. However, the retrievals are bust run on a compute cluster or a node with more than 4 cores. We will use the "benchmark" system, WASP-43b as our example utilizing the HST WFC3 data presented in Kreidberg et al. 2014b.
Software Requirements: This runs in the python 3 anaconda environment. It is also crucial that anaconda numba is installed as many of the routines are optimized using numba's "@jit" decorator (http://numba.pydata.org/).
This first segment loads in the routines from fm.py and the correlated-K coefficients. There are two sets of correlated-K coefficients (which I've called "xsecs" here). There are ones taylored for HST WFC3+STIS (xsects_HST function in fm.py) and JWST (xsects_JWST in fm.py). The WFC3+STIS correlated-K coefficients are generated at an R=200 longwards of 1 $\mu$m (up to 5$\mu$m and R=500 from 0.3 - 1 $\mu$m. Note...these are not sampled cross-sections so each resolution element at that R is correctly computed and matches line-by-line when binned to that same R.
Note that the "core" set of routines are all in fm.py. If you want to know more about what is in the sausage, look into fm.py.
#import all of the functions in fm, namely, the CK-coefficients (may take a minute)
from chimera import *
%matplotlib notebook
/Users/michaelline/anaconda3/lib/python3.6/site-packages/h5py/__init__.py:36: FutureWarning: Conversion of the second argument of issubdtype from `float` to `np.floating` is deprecated. In future, it will be treated as `np.float64 == np.dtype(float).type`. from ._conv import register_converters as _register_converters /Users/michaelline/anaconda3/lib/python3.6/site-packages/numpy/core/fromnumeric.py:2957: RuntimeWarning: Mean of empty slice. out=out, **kwargs) /Users/michaelline/anaconda3/lib/python3.6/site-packages/numpy/core/_methods.py:80: RuntimeWarning: invalid value encountered in double_scalars ret = ret.dtype.type(ret / rcount)
Cross-sections Loaded
stellar_file = 'sum_star.h5'
temp = 5000
logmh = 0
logg = 4.0
stellar_db = 'phoenix'
make_stellar(temp,logmh,logg,stellar_db,stellar_file)
#preload CK-coeffs--a giant array/variable to be passed--inputs are lower wavenumber, upper wavenumber
#between 2000 and 30000 cm-1 for HST--R=200 > 1 um, then R=500 < 1 um
#to convert between microns and wavenumbers-- wavelength [um] = 10,000/wavenumber [cm-1]
#make sure xsec wavenumber/wavelength range is *larger* than data wavelength range
wnomin = 2000
wnomax = 30000
observatory='HST'
directory = os.path.join(os.getcwd(),'..','..','ABSCOEFF_CK')
xsecs=xsects(wnomin, wnomax, observatory, directory,stellar_file=stellar_file)
This segement defines the various atmospheric quantities and assignes them values for the generation of a simple transmission spectrum. A description of each parameter along with a reasonable range of values is given as a comment following the assigned value. All of the parameters are then put into the parameter "state-vector" array, x.
#setup "input" parameters. We are defining our 1D atmosphere with these
#the parameters
#planet/star system params--xRp is the "Rp" free parameter, M right now is fixed, but could be free param
Rp=1.036 #Planet radius in Jupiter Radii--this will be forced to be 10 bar radius--arbitrary (scaling to this is free par)
Rstar=0.667 # #Stellar Radius in Solar Radii
M =2.034 #Mass in Jupiter Masses
D=0.01526 #semimajor axis in AU
#TP profile params (3--Guillot 2010, Parmentier & Guillot 2013--see Line et al. 2013a for implementation)
Tirr=1400 #Irradiation temperature as defined in Guillot 2010
logKir=-1.5 #TP profile IR opacity (log there-of) controlls the "vertical" location of the gradient
logg1=-0.7 #single channel Vis/IR (log) opacity. Controls the delta T between deep T and TOA T
Tint=200 #interior temperature...this would be the "effective temperature" if object were not irradiated
#log Gas abundances
H2O=-4.
CH4=-10.
CO=-4.
CO2=-10.
NH3=-10.
N2=-4.
HCN=-10.
H2S=-5.
PH3=-6.
C2H2=-10.
C2H6=-10.
Na=-6.
K=-7.
TiO=-10.
VO=-10.
FeH=-10.
H=-10.
em=-10.
hm=-10.
#Ackerman & Marley 2001 Cloud parameters--physically motivated with Mie particles
logKzz=7 #log Kzz (cm2/s)--valid range: 2 - 11 -- higher values make larger particles
fsed=2.0 #sediminetation efficiency--valid range: 0.5 - 5--lower values make "puffier" more extended cloud
logPbase=-1.0 #cloud base pressure--valid range: -6.0 - 1.5
logCldVMR=-5.5 #cloud condensate base mixing ratio (e.g, see Fortney 2005)--valid range: -15 - -2.0
#simple 'grey+rayleigh' parameters just in case you don't want to use a physically motivated cloud
#(most are just made up anyway since we don't really understand all of the micro-physics.....)
logKcld = -40 #uniform in altitude and in wavelength "grey" opacity (it's a cross-section)--valid range: -50 - -10
logRayAmp = -30 #power-law haze amplitude (log) as defined in des Etangs 2008 "0" would be like H2/He scat--valid range: -30 - 3
RaySlope = 0 #power law index 4 for Rayleigh, 0 for "gray". Valid range: 0 - 6
#10 bar radiuss scaling param (only used in transmission)
xRp=0.991
#stuffing all variables into state vector array
x=np.array([Tirr, logKir,logg1,Tint, 0, 0, 0,0, Rp*xRp, Rstar, M, logKzz, fsed,logPbase,logCldVMR, logKcld, logRayAmp, RaySlope])
#gas scaling factors to mess with turning on various species
#set to "0" to turn off a gas. Otherwise keep set at 1
#thermochemical gas profile scaling factors
# 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21
#H2O CH4 CO CO2 NH3 N2 HCN H2S PH3 C2H2 C2H6 Na K TiO VO FeH H H2 He e- h- mmw
gas_scale=np.array([H2O,CH4,CO,CO2,NH3,N2,HCN,H2S,PH3,C2H2,C2H6,Na,K,TiO,VO ,FeH,H,-50.,-50.,em, hm,-50.]) #
#WASP43b transmission spectrum from kreidberg et al. 2014...a 3 column ascii file
wlgrid, y_meas, err=np.loadtxt('w43b_trans.txt').T
Here we call the forward model routine "fx" (think F(x)) from fm.py. fx controls the input values and calls the relevent functions to compute the transmission spectrum. The inputs into fx are the parameter state vector, "x", the data wavelength grid, "wlgrid", the gas scaling factors (for turning off particular gases), "gas_scale", and the correlated-K tables, "xsects". Fx then returns the simulated model spectrum ($(R_p/R_{\star})^2$) at the native CK-table resolution, "y_mod", the native wavenumber grid, "wno", the data wavelength grid binned model spectrum, "y_binned". The "atm" array contains the generated temperature-pressure profile and gas mixing ratio profiles generated under the chemically consistent assumption.
#calling forward model, fx. This will produce the (Rp/Rstar)^2 spectrum....
y_binned,y_mod,wno,atm=fx_trans_free(x,wlgrid,gas_scale, xsecs) #returns binned model spectrum, higher res model spectrum, wavenumber grid, and vertical abundance profiles from chemistry
print('DONE')
DONE
from matplotlib.pyplot import *
from matplotlib.ticker import FormatStrFormatter
%matplotlib notebook
#unpacking variables
#P is in bars
#T is in K
#H2O, CH4,CO,CO2,NH3,Na,K,TiO,VO,C2H2,HCN,H2S,FeH,H2,He are gas mixing ratio profiles
#qc is the condensate abundance profile given an "f_sed" value and cloud base pressure
#r_eff is the effective cloud droplet radius given (see A&M 2001 or Charnay et al. 2017)
#f_r is the mixing ratio array for each of the cloud droplet sizes.
P,T, H2O, CH4,CO,CO2,NH3,Na,K,TiO,VO,C2H2,HCN,H2S,FeH,H2,He,H,e, Hm,qc,r_eff,f_r=atm
fig2, ax1=subplots()
#feel free to plot whatever you want here....
ax1.semilogx(H2O,P,'b',ls='--',lw=2,label='H2O')
ax1.semilogx(CH4,P,'black',ls='--',lw=2,label='CH4')
ax1.semilogx(CO,P,'g',ls='--',lw=2,label='CO')
ax1.semilogx(CO2,P,'orange',ls='--',lw=2,label='CO2')
ax1.semilogx(NH3,P,'darkblue',ls='--',lw=2,label='NH3')
ax1.semilogx(Na,P,'b',lw=2,label='Na')
ax1.semilogx(K,P,'g',lw=2,label='K')
ax1.semilogx(TiO,P,'k',lw=2,label='TiO')
ax1.semilogx(VO,P,'orange',lw=2,label='VO')
ax1.semilogx(qc,P,'gray',lw=1,ls='--',label='Cond. VMR.') #<---- A&M Cloud Condensate VMR profile (not droplets)
ax1.set_xlabel('Mixing Ratio',fontsize=20)
ax1.set_ylabel('Pressure [bar]',fontsize=20)
ax1.semilogy()
ax1.legend(loc=4,frameon=False)
ax1.axis([1E-9,1,100,1E-7])
#plotting TP profile on other x-axis
ax2=ax1.twiny()
ax2.semilogy(T,P,'r-',lw='4',label='TP')
ax2.set_xlabel('Temperature [K]',color='r',fontsize=20)
ax2.axis([0.8*T.min(),1.2*T.max(),100,1E-6])
for tl in ax2.get_xticklabels(): tl.set_color('r')
ax2.legend(loc=1,frameon=False)
savefig('./plots/atmosphere_transmission_WFC3_FREE.pdf',fmt='pdf')
show()
#close()
#finally doing some plotting
#and the usual matplotlib shenanigans
ymin=np.min(y_binned)*1E2*0.99
ymax=np.max(y_binned)*1E2*1.01
fig1, ax=subplots()
xlabel('$\lambda$ ($\mu$m)',fontsize=18)
ylabel('(R$_{p}$/R$_{*}$)$^{2} \%$',fontsize=18)
minorticks_on()
errorbar(wlgrid, y_meas*100, yerr=err*100, xerr=None, fmt='Dk', label='Data')
plot(wlgrid, y_binned*1E2,'ob',label='Binned Model')
plot(1E4/wno, y_mod*1E2, label='Model')
ax.set_xscale('log')
ax.set_xticks([0.3, 0.5,0.8,1,1.4, 2, 3, 4, 5])
ax.axis([0.3,5,ymin,ymax])
ax.get_xaxis().set_major_formatter(matplotlib.ticker.ScalarFormatter())
ax.tick_params(length=10,width=1,labelsize='large',which='major')
legend(frameon=False)
savefig('./plots/transmission_spectrum_WFC3_FREE.pdf',fmt='pdf')
show()
#close()
from matplotlib.pyplot import *
from matplotlib.ticker import FormatStrFormatter
#log Gas abundances
H2O=-4.
CH4=-10.
CO=-4.
CO2=-10.
NH3=-10.
N2=-4.
HCN=-10.
H2S=-5.
PH3=-6.
C2H2=-10.
C2H6=-10.
Na=-6.
K=-7.
TiO=-15.
VO=-15.
FeH=-15.
H=-10.
em=-10.
hm=-10.
#Ackerman & Marley 2001 Cloud parameters--physically motivated with Mie particles
logKzz=7 #log Kzz (cm2/s)--valid range: 2 - 11 -- higher values make larger particles
fsed=2.0 #sediminetation efficiency--valid range: 0.5 - 5--lower values make "puffier" more extended cloud
logPbase=-1.0 #cloud base pressure--valid range: -6.0 - 1.5
logCldVMR=-15.5 #cloud condensate base mixing ratio (e.g, see Fortney 2005)--valid range: -15 - -2.0
#simple 'grey+rayleigh' parameters just in case you don't want to use a physically motivated cloud
#(most are just made up anyway since we don't really understand all of the micro-physics.....)
logKcld = -30 #uniform in altitude and in wavelength "grey" opacity (it's a cross-section)--valid range: -50 - -10
logRayAmp =2 #power-law haze amplitude (log) as defined in des Etangs 2008 "0" would be like H2/He scat--valid range: -30 - 3
RaySlope = 4 #power law index 4 for Rayleigh, 0 for "gray". Valid range: 0 - 6
#10 bar radiuss scaling param (only used in transmission)
xRp=0.991
#stuffing all variables into state vector array
x=np.array([Tirr, logKir,logg1,Tint, 0, 0, 0,0, Rp*xRp, Rstar, M, logKzz, fsed,logPbase,logCldVMR, logKcld, logRayAmp, RaySlope])
#gas scaling factors to mess with turning on various species
#set to "0" to turn off a gas. Otherwise keep set at 1
#thermochemical gas profile scaling factors
# 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21
#H2O CH4 CO CO2 NH3 N2 HCN H2S PH3 C2H2 C2H6 Na K TiO VO FeH H H2 He e- h- mmw
gas_scale=np.array([H2O,CH4,CO,CO2,NH3,N2,HCN,H2S,PH3,C2H2,C2H6,Na,K,TiO,VO ,FeH,H,-50.,-50.,em, hm,-50.]) #
y_binned,y_mod,wno,atm=fx_trans_free(x,wlgrid,gas_scale, xsecs) #returns binned model spectrum, higher res model spectrum, wavenumber grid, and vertical abundance profiles from chemistry
gas_scale2=np.ones(len(gas_scale))*gas_scale
gas_scale2[0:-1]=-50
y_binned2,y_mod2,wno2,atm2=fx_trans_free(x,wlgrid,gas_scale2, xsecs) #returns binned model spectrum, higher res model spectrum, wavenumber grid, and vertical abundance profiles from chemistry
ymin=np.min(y_binned)*1E2*0.99
ymax=np.max(y_binned)*1E2*1.01
fig1, ax=subplots()
xlabel('$\lambda$ ($\mu$m)',fontsize=18)
ylabel('(R$_{p}$/R$_{*}$)$^{2} \%$',fontsize=18)
minorticks_on()
errorbar(wlgrid, y_meas*100, yerr=err*100, xerr=None, fmt='Dk', label='Data')
plot(wlgrid, y_binned*1E2,'ob',label='Binned Model')
plot(1E4/wno, y_mod*1E2, label='Model')
plot(1E4/wno, y_mod2*1E2,label='Cloud+H2/He Continuum')
ax.set_xscale('log')
ax.set_xticks([0.3, 0.5,0.8,1,1.4, 2, 3, 4, 5])
ax.axis([0.3,5,ymin,ymax])
ax.get_xaxis().set_major_formatter(matplotlib.ticker.ScalarFormatter())
ax.tick_params(length=10,width=1,labelsize='large',which='major')
legend(frameon=False)
show()
#close()
The "retrieval" is performed using the DYNESTY (https://dynesty.readthedocs.io/en/latest/index.html) nested sampling suite. It's basically like all of the others (e.g., multinest, pymultnest, nestle etc.), though it's more flexible in terms of sampling methods optimized for certain numbers of parameters. This example just uses "generic" settings.
#set up a "dynest" nested sampling run--see https://dynesty.readthedocs.io/en/latest/index.html
#a super cool useful comparision of all MCMC/Multinests out there...
#http://mattpitkin.github.io/samplers-demo/pages/samplers-samplers-everywhere/
#(feel free to mix and match samplers--probably not a bad idea...)
#for safty, just reloading everything again
import numpy as np
from matplotlib import pyplot as plt
import dynesty
from multiprocessing import Pool
from fm import *
import pickle
#load crosssections between wnomin and wnomax
xsecs=xsects_HST(6000,9400) #make sure this range is *larger* than the data wavelength grid (but not by too much)
#6000 cm-1 (1.68 um) to 9400 cm-1 (1.06 um)
/Users/michaelline/anaconda3/lib/python3.6/site-packages/numpy/core/fromnumeric.py:2957: RuntimeWarning: Mean of empty slice. out=out, **kwargs) /Users/michaelline/anaconda3/lib/python3.6/site-packages/numpy/core/_methods.py:80: RuntimeWarning: invalid value encountered in double_scalars ret = ret.dtype.type(ret / rcount)
Cross-sections Loaded
This computes "chi-square" log-likelihood function. The input value is "theta" which is the same as "x"--the parameter state vector, though just for the parameters we care about. These can be a sub-set of the full "x" defined above. The first block defines a bunch of parameters--same as what gets passed into fx--to a generic default value. Thes parameter values get overridden with the values in the "theta (parameters to retrieve) vector. In this particular example only Tirr, H2O,CO,NH3,CH4, logKcld, xRp are retrieved. All other values passed into fx are assigned the other values.
#defining log-likelihood function
# log-likelihood
def loglike(theta):
#setting default parameters---will be fixed to these values unless replaced with 'theta'
#planet/star system params--xRp is the "Rp" free parameter, M right now is fixed, but could be free param
Rp= 1.036#0.930#*x[4]# Planet radius in Jupiter Radii--this will be forced to be 10 bar radius--arbitrary (scaling to this is free par)
Rstar=0.667#0.598 #Stellar Radius in Solar Radii
M =2.034#1.78 #Mass in Jupiter Masses
#TP profile params (3--Guillot 2010, Parmentier & Guillot 2013--see Line et al. 2013a for implementation)
Tirr=1400#1500#x[0]#544.54 #terminator **isothermal** temperature--if full redistribution this is equilibrium temp
logKir=-5 #TP profile IR opacity controlls the "vertical" location of the gradient
logg1=0.0 #single channel Vis/IR opacity. Controls the delta T between deep T and TOA T
Tint=0.
#A&M Cloud parameters--includes full multiple scattering (for realzz) in both reflected and emitted light
logKzz=7 #log Rayleigh Haze Amplitude (relative to H2)
fsed=3.0 #haze slope--4 is Rayeigh, 0 is "gray" or flat.
logPbase=-1.0 #gray "large particle" cloud opacity (-35 - -25)
logCldVMR=-25.0 #cloud fraction
#simple 'grey+rayleigh' parameters--non scattering--just pure extinction
logKcld = -40
logRayAmp = -30
RaySlope = 0
H2O=-15.
CH4=-15.
CO=-15.
CO2=-15.
NH3=-15.
N2=-15.
HCN=-15.
H2S=-15.
PH3=-15.
C2H2=-15.
C2H6=-15.
Na=-15.
K=-15.
TiO=-15.
VO=-15.
FeH=-15.
H=-15.
em=-15.
hm=-15.
#unpacking parameters to retrieve (these override the fixed values above)
Tirr, H2O,CO,NH3,CH4, logKcld, xRp=theta
##all values required by forward model go here--even if they are fixed
x=np.array([Tirr, logKir,logg1, Tint,0,0,0,0, Rp*xRp, Rstar, M, logKzz, fsed,logPbase,logCldVMR, logKcld, logRayAmp, RaySlope])
# 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21
#H2O CH4 CO CO2 NH3 N2 HCN H2S PH3 C2H2 C2H6 Na K TiO VO FeH H H2 He e- h- mmw
gas_scale=np.array([H2O,CH4,CO,CO2,NH3,N2,HCN,H2S,PH3,C2H2,C2H6,Na,K,TiO,VO ,FeH,H,-50.,-50.,em, hm,-50.]) #
foo=fx_trans_free(x,wlgrid,gas_scale,xsecs)
y_binned=foo[0]
loglikelihood=-0.5*np.sum((y_meas-y_binned)**2/err**2) #your typical "quadratic" or "chi-square"
return loglikelihood
Nested samplers use this "cube" concept whereby all of the parameter values are transformed on the interval [0,1]. This makes the sampling more-or-less scale indepenent. The prior_transform function maps those back onto the actual parameter value ranges we want. In this case, the prior ranges are all uniform over the range as defined here. Live points will be drawn from this hypercube uniformly. These "live point" parameter values are then assesed by the log-likelihood function.
#defining prior cube (cube is a standard multi-nest way of doing priors)
def prior_transform(utheta):
Tirr, H2O,CO,CO2,CH4, logKcld, xRp=utheta
#uniform prior ranges--each "variable", say Tirr is sampled over the interval [0,1]--the numbers here transform that
Tirr = 1900 * Tirr + 100 #Tirr uniform from 400 - 3000K (add the lower value to the multiplier to get upper bound)
H2O=12*H2O-12
CO=12*CO-12
CO2=12*CO2-12
CH4=12*CH4-12
logKcld=20*logKcld-45
xRp=1*xRp+0.5
return Tirr, H2O,CO,CO2,CH4, logKcld, xRp
This segment first loads in the "3 column" data file...Here it is just a pickle, but this could be replaced with a 3-column ascii file read in--wavelength grid, data values, error bars. The other knobs are self-explanatory. The number of live points is "problem dependent". Some problems can get away with less. I wouldn't go below 100. It's safer to use more (e.g., just like walkers in emcee). For "real science" I prefer 1000+ live points to make sure the posterior is well sampled and no modes are missed. Of course, this takes longer to run, but no sense in getting the the wrong answer faster!
#setting up other dynesty run params and loading in the data
#WASP43b transmission spectrum...a 3 column ascii file
wlgrid, y_meas, err=np.loadtxt('w43b_trans.txt').T
outname='./OUTPUT/dyn_output_trans_WFC3_FREE.pic' #dynesty output file name (saved as a pickle)
Nparam=7 #number of parameters--make sure it is the same as what is in prior and loglike
Nproc=4 #number of processors for multi processing--best if you can run on a 12 core+ node or something
Nlive=500 #number of nested sampling live points
This calls the Nested Sampler function to compute the posterior. This may take a few hours depending on the number of parameters and number of live points.
#running the "standard" nested sampler. Again, see https://dynesty.readthedocs.io/en/latest/index.html for details
#depending on the number of live points (I like lots, usually 1000+, but I'm paranoid), number of params, and number of
#processors on your computer, this could take some time. -- your computer will make loud noises..1-2 hours with 4 cores
#in real life, you should probably run this on a multi-core machine (8+)
#in really real life, you should probably use pymultinest...it tends to "converge" faster with fewer likelihood
#evaluations than does dynest, even with the same number of live points (see: https://mattpitkin.github.io/samplers-demo/pages/samplers-samplers-everywhere/)
#e.g., when I ran this with 56 cores (2x28 core nodes) with multinest it took ~20 min.
import time
pool = Pool(processes=Nproc)
dsampler = dynesty.NestedSampler(loglike, prior_transform, ndim=Nparam,
bound='multi', sample='auto', nlive=Nlive,
update_interval=3., pool=pool, queue_size=Nproc)
#this executes and runs it
t1=time.time()
dsampler.run_nested()
t2=time.time()
print("Run Time:", t2-t1)
#extracting results from sampler object
dres = dsampler.results
#dumping as a pickle
pickle.dump(dres,open(outname,'wb'))
#some real time dynesty output/status will pop up down here (this one took 51931 lnL calls over 4309s = 1.2 hrs)
Plots the corner plots and spectral fit plot. Note: All of the plotting below can be used for a "previously" generated run. Just start here and load in the sampler output.
#plotting dynest runs corner plot
from matplotlib import pyplot as plt
from dynesty import plotting as dyplot
import pickle
labels=['Tirr', 'H2O', 'CO', 'CO2', 'CH4' ,'logKcld','xRp']
#import past run
#samples=pickle.load(open('dyn_output_100LP.pic','rb')) #an example 100 live point run
samples=pickle.load(open('./OUTPUT/dyn_output_trans_WFC3_FREE.pic','rb'))
#samples=pickle.load(open('dyn_output_1000LP.pic','rb')) #an example 1000 live point run
#printing evidence:
print('ln(Z)= ', samples.logz[-1])
# corner plot
#NOTE...DEFAULT CONFIDENCE INTERVAL IS THE 95%!!!
fig, axes = dyplot.cornerplot(samples,smooth=0.05, color='blue',show_titles=True, labels=labels,title_kwargs={'y': 1.04}, fig=plt.subplots(7, 7, figsize=(12, 12)))
plt.savefig('./plots/Dynesty_WFC3_free_stair_pairs.pdf',fmt='pdf')
plt.show()
#plt.close()
ln(Z)= -26.070611115262167
Generating spectra from parameters of a subset of samples drawn from the posterior. These spectra, as always, are then summarized with their median, 1-, 2-, sigma confidence intervals.
import numpy as np
xsecs=xsects_HST(2000, 30000)
Nspectra=200
#loading in data again just to be safe
wlgrid, y_meas, err=np.loadtxt('w43b_trans.txt').T
#setting up default parameter values--SET THESE TO SAME VALUES AS IN LOG-LIKE FUNCTION
#planet/star system params--xRp is the "Rp" free parameter, M right now is fixed, but could be free param
#setting default parameters---will be fixed to these values unless replaced with 'theta'
#planet/star system params--xRp is the "Rp" free parameter, M right now is fixed, but could be free param
Rp= 1.036#0.930#*x[4]# Planet radius in Jupiter Radii--this will be forced to be 10 bar radius--arbitrary (scaling to this is free par)
Rstar=0.667#0.598 #Stellar Radius in Solar Radii
M =2.034#1.78 #Mass in Jupiter Masses
#TP profile params (3--Guillot 2010, Parmentier & Guillot 2013--see Line et al. 2013a for implementation)
Tirr=1400#1500#x[0]#544.54 #terminator **isothermal** temperature--if full redistribution this is equilibrium temp
logKir=-1.5 #TP profile IR opacity controlls the "vertical" location of the gradient
logg1=-0.7 #single channel Vis/IR opacity. Controls the delta T between deep T and TOA T
Tint=200.
#A&M Cloud parameters--includes full multiple scattering (for realzz) in both reflected and emitted light
logKzz=7 #log Rayleigh Haze Amplitude (relative to H2)
fsed=3.0 #haze slope--4 is Rayeigh, 0 is "gray" or flat.
logPbase=-1.0 #gray "large particle" cloud opacity (-35 - -25)
logCldVMR=-25.0 #cloud fraction
#simple 'grey+rayleigh' parameters--non scattering--just pure extinction
logKcld = -40
logRayAmp = -30
RaySlope = 0
H2O=-15.
CH4=-15.
CO=-15.
CO2=-15.
NH3=-15.
N2=-15.
HCN=-15.
H2S=-15.
PH3=-15.
C2H2=-15.
C2H6=-15.
Na=-15.
K=-15.
TiO=-15.
VO=-15.
FeH=-15.
H=-15.
em=-15.
hm=-15.
#weighting the posterior samples for appropriate random drawing
from dynesty import utils as dyfunc
samp, wts = samples.samples, np.exp(samples.logwt - samples.logz[-1])
samples2 = dyfunc.resample_equal(samp, wts)
#choosing random indicies to draw from properly weighted posterior samples
draws=np.random.randint(0, samples2.shape[0], Nspectra)
Nwno_bins=xsecs[2].shape[0]
y_mod_array=np.zeros((Nwno_bins, Nspectra))
y_binned_array=np.zeros((len(wlgrid), Nspectra))
for i in range(Nspectra):
print(i)
#make sure this is the same as in log-Like
Tirr, H2O,CO,CO2,CH4, logKcld, xRp=samples2[draws[i],:]
x=np.array([Tirr, logKir,logg1, Tint,0,0,0,0, Rp*xRp, Rstar, M, logKzz, fsed,logPbase,logCldVMR, logKcld, logRayAmp, RaySlope])
print(samples.samples[draws[i],:])
gas_scale=np.array([H2O,CH4,CO,CO2,NH3,N2,HCN,H2S,PH3,C2H2,C2H6,Na,K,TiO,VO ,FeH,H,-50.,-50.,em, hm,-50.]) #
y_binned,y_mod,wno,atm=fx_trans_free(x,wlgrid,gas_scale,xsecs)
y_mod_array[:,i]=y_mod
y_binned_array[:,i]=y_binned
#saving these arrays since it takes a few minutes to generate
pickle.dump([wlgrid, y_meas, err, y_binned_array, wno, y_mod_array],open('./OUTPUT/spectral_samples_trans_dyn_wfc3_free.pic','wb'))
/Users/michaelline/anaconda3/lib/python3.6/site-packages/numpy/core/fromnumeric.py:2957: RuntimeWarning: Mean of empty slice. out=out, **kwargs) /Users/michaelline/anaconda3/lib/python3.6/site-packages/numpy/core/_methods.py:80: RuntimeWarning: invalid value encountered in double_scalars ret = ret.dtype.type(ret / rcount)
Cross-sections Loaded 0 [ 8.25002782e+02 -1.00816618e+01 -6.86196367e+00 -1.19918977e-02 -5.99301772e+00 -2.72575041e+01 9.82173666e-01] 1 [880.09792159 -5.39522694 -4.7082032 -11.13418576 -8.38708626 -43.52710534 0.99787328] 2 [434.82712759 -7.94462418 -9.48849615 -8.29330645 -5.99966581 -41.45502567 1.00316492] 3 [724.66139926 -5.21535775 -10.03713856 -11.35614525 -6.99120368 -41.16845652 0.99852339] 4 [ 1.95154519e+03 -6.28940622e+00 -4.88271386e+00 -3.38630466e+00 -1.12913493e+01 -3.21845776e+01 1.37253093e+00] 5 [137.1537112 -1.82701421 -5.83466432 -2.95016278 -9.18382428 -30.10621449 0.99769601] 6 [457.5507551 -1.68527893 -9.22804235 -10.15536262 -5.84058898 -32.55449285 0.99725451] 7 [478.00916569 -3.82114963 -2.69572006 -10.17370597 -8.47896676 -39.03022875 0.99819003] 8 [906.2856652 -3.76567629 -10.08050133 -10.34276094 -9.18432202 -29.02793739 1.07781098] 9 [549.72131916 -1.54692799 -5.58958661 -4.39836165 -3.76129032 -29.24902483 0.99509267] 10 [ 9.65250764e+02 -1.11430024e+01 -9.94171909e+00 -6.31174802e-02 -1.19835345e+01 -3.32773480e+01 9.98004392e-01] 11 [502.17946755 -0.55438818 -8.14609091 -10.28931802 -7.89367528 -25.11732604 0.99553915] 12 [651.63768406 -4.64913296 -4.9389888 -11.79426326 -5.88916501 -34.82000289 0.99811212] 13 [685.43903906 -4.61964462 -9.43052842 -1.21829102 -8.9825166 -27.16325746 0.99390621] 14 [ 1.52085075e+03 -3.54090650e+00 -1.20033933e+00 -6.48583234e-02 -9.21460214e+00 -3.69450557e+01 8.95901430e-01] 15 [495.36129839 -8.41203262 -4.43071567 -6.14355547 -0.93276652 -41.87709749 1.27810956] 16 [ 3.29160771e+02 -3.15110597e-01 -9.32968245e-01 -5.00586937e+00 -1.17097368e+01 -3.66010787e+01 9.95695176e-01] 17 [ 1.76723288e+03 -1.07947149e+01 -1.33533425e+00 -1.00728594e+01 -3.97894741e+00 -2.68186003e+01 9.58576817e-01] 18 [540.25702042 -3.65640493 -8.86451427 -10.69824913 -6.81637519 -44.85137825 0.99784847] 19 [228.74271801 -2.5721904 -4.42090749 -5.1068006 -3.6762154 -40.29668352 0.99731317] 20 [466.90810044 -0.72687764 -4.4495029 -3.68862478 -6.78565416 -39.13166266 0.97517185] 21 [222.20656292 -3.43459583 -2.21469499 -3.72751441 -3.77207874 -43.17933115 0.99747436] 22 [ 1.58313583e+03 -2.92085031e-03 -9.18961840e-01 -6.04802191e+00 -4.20434284e+00 -2.80591902e+01 9.91519732e-01] 23 [949.92886606 -0.9681169 -2.18650432 -10.74789598 -6.88103006 -30.64838464 0.99352189] 24 [ 1.63271600e+03 -2.87500630e+00 -6.56451572e-02 -8.98339473e+00 -6.57975070e+00 -3.84976110e+01 9.95434926e-01] 25 [381.78542267 -4.49284167 -7.62491127 -11.17945402 -3.75115704 -31.44843439 0.9976027 ] 26 [ 1.26589022e+03 -5.93643912e+00 -7.69563448e+00 -5.52971144e+00 -4.33080529e+00 -2.90886801e+01 9.75723919e-01] 27 [663.70495024 -5.25313867 -5.20179488 -9.59309675 -8.75552049 -35.86159852 0.99836777] 28 [221.51746298 -6.71443904 -9.12120285 -6.2239867 -11.93973898 -30.10884221 0.84309063] 29 [ 1.13917596e+03 -1.97572495e+00 -2.10500279e+00 -1.08737284e+01 -5.47517682e+00 -4.44499805e+01 9.91819013e-01] 30 [ 1.79816948e+03 -7.70793967e+00 -3.17428075e+00 -9.13183096e+00 -5.48395547e+00 -3.46593860e+01 9.90253150e-01] 31 [368.79091677 -4.26331043 -7.92269461 -11.80761346 -5.05736453 -37.30144933 0.99773302] 32 [466.93966866 -3.38697984 -2.83187885 -11.09605638 -10.15023418 -37.43519128 0.99830053] 33 [500.73001522 -3.19142442 -3.31023701 -11.56580653 -7.88231181 -32.74401747 0.99794867] 34 [ 1.32815566e+03 -4.66755385e+00 -1.37972482e+00 -1.79751858e+00 -9.59906129e+00 -3.31054615e+01 1.30075289e+00] 35 [418.31270699 -1.53224507 -4.13662303 -8.96218024 -10.46185288 -33.24438337 0.99738563] 36 [400.86183453 -5.94069152 -9.24525009 -2.77196317 -4.81504403 -30.4858743 0.98136976] 37 [ 3.68405284e+02 -5.88185264e-01 -2.87841072e-01 -1.09762362e+01 -4.28888077e+00 -4.26463426e+01 9.99389635e-01] 38 [568.96557522 -9.51347525 -7.23271647 -11.63414669 -4.43150487 -31.2085586 0.99704492] 39 [407.74436642 -10.89515771 -8.30374649 -2.70426877 -9.55171395 -27.94615357 1.0100423 ] 40 [525.26686421 -3.81169968 -9.79316962 -7.62655222 -5.82463693 -37.99068361 0.99756675] 41 [726.70459286 -4.1724994 -5.70680997 -10.64715231 -6.62447999 -36.8348248 0.99773466] 42 [ 1.66300444e+03 -3.50672134e+00 -1.89142095e+00 -3.75666626e+00 -6.57791759e-01 -3.64121702e+01 9.92212192e-01] 43 [ 5.81121795e+02 -7.11864381e+00 -1.05507915e+01 -8.96685883e+00 -2.45250202e-01 -2.89982710e+01 9.96132610e-01] 44 [ 9.79157653e+02 -6.90136331e+00 -1.02294343e+01 -4.18603822e-01 -8.08050911e+00 -4.46746288e+01 8.61769644e-01] 45 [559.65468389 -1.29111894 -2.81862437 -7.02241686 -11.6729928 -34.47022526 0.99661 ] 46 [546.00647275 -4.04296672 -11.30114346 -5.57724465 -3.69582534 -36.90654844 0.99657374] 47 [879.94185332 -10.2230063 -5.94368156 -7.7795128 -3.03937443 -40.72550428 0.99523499] 48 [834.58154328 -6.55449357 -11.85122459 -1.06037309 -4.82629698 -44.50078727 0.98980483] 49 [ 1.85503113e+03 -6.39153747e+00 -1.52273396e+00 -5.89643647e+00 -6.34363120e+00 -3.42535143e+01 9.97570656e-01] 50 [257.72102239 -2.26133117 -4.22131909 -4.38533891 -11.26426564 -43.70439686 0.9973994 ] 51 [494.056412 -4.16897899 -8.02536933 -7.60999677 -8.20479611 -33.57558614 0.99848012] 52 [537.18005664 -4.41005849 -4.89591017 -6.69353197 -7.13946703 -34.76581239 0.99828547] 53 [ 1.66271014e+03 -4.73072311e-01 -8.91067114e+00 -6.67787722e+00 -9.67203445e+00 -3.30533031e+01 9.93476759e-01] 54 [ 1.31960802e+03 -3.54237853e-01 -6.75486714e+00 -1.12855420e+01 -8.65388443e+00 -3.67868515e+01 9.95041183e-01] 55 [ 1.39302157e+03 -6.64075455e+00 -5.48982081e+00 -7.38528931e+00 -5.10380537e+00 -3.63888663e+01 9.84071571e-01] 56 [ 1.38821250e+03 -9.39545407e-01 -3.12069762e+00 -9.02424972e+00 -5.72199918e+00 -2.54222775e+01 7.54735285e-01] 57 [485.41070327 -1.25219083 -11.12015275 -4.52161907 -7.43760402 -34.45681332 0.99690593] 58 [ 1.05654764e+03 -3.11448751e+00 -8.50565316e+00 -1.13897917e+01 -7.43114557e-01 -3.48272527e+01 9.93237290e-01] 59 [514.52555519 -1.41452694 -9.39854865 -8.76190503 -11.11520153 -42.35052563 0.99699106] 60 [293.29564354 -11.04012579 -6.03816488 -10.31006377 -1.74348998 -38.76442134 0.9970706 ] 61 [853.00110737 -11.27189245 -9.43078731 -11.11213748 -1.27735011 -25.6646967 0.9894778 ] 62 [ 1.29692635e+03 -7.56604668e+00 -1.17575520e+01 -8.66893623e+00 -4.67402633e+00 -2.99871664e+01 1.15473167e+00] 63 [303.74413861 -4.23533608 -4.45480776 -5.21222385 -2.93765006 -41.81687034 0.99730976] 64 [276.75991519 -2.37455417 -5.59225439 -7.43288331 -8.93469874 -42.52032367 0.99789533] 65 [ 1.20578635e+03 -8.42438496e+00 -6.11009160e+00 -1.87492851e+00 -4.82644394e+00 -3.80293913e+01 9.40863105e-01] 66 [ 1.95811088e+03 -6.58080760e+00 -5.93798949e+00 -6.23153380e+00 -6.10423836e+00 -4.49769502e+01 7.09458886e-01] 67 [427.39347126 -2.3294651 -3.9381686 -8.4677459 -5.29419003 -33.85436662 0.99735851] 68 [484.66758954 -10.41008962 -0.59311705 -5.38620596 -8.81589565 -25.97387153 1.06345094] 69 [879.26828875 -8.38648311 -7.32403462 -7.73312246 -1.64501567 -39.98891484 0.9892201 ] 70 [646.61713526 -5.44867209 -4.3473526 -6.86329271 -10.27826194 -40.72940919 0.99852224] 71 [617.71486676 -11.62343155 -5.73323912 -4.84043613 -8.0685746 -33.10092726 0.99901173] 72 [454.38420908 -9.98761274 -5.50922561 -3.08684799 -11.21679679 -28.79353507 1.00140576] 73 [ 1.32649269e+03 -5.93738449e+00 -8.22718959e+00 -6.47485520e+00 -1.04787733e+01 -3.50503765e+01 9.97305430e-01] 74 [463.19349173 -1.08833797 -9.60342889 -8.1576795 -7.90719 -41.16849806 0.99731812] 75 [ 1.87875993e+03 -2.39448916e+00 -3.24363732e-01 -9.60768815e+00 -6.44223394e+00 -3.20998439e+01 9.95632488e-01] 76 [452.27854407 -0.64430644 -10.6855115 -10.94823858 -7.78335663 -37.34066235 0.99763333] 77 [575.7883574 -8.33048485 -9.86850467 -7.72388594 -4.61322149 -37.42816832 0.99843516] 78 [ 1.27992348e+03 -4.85402995e+00 -2.18600968e+00 -1.76192366e+00 -2.15636657e+00 -4.14578718e+01 8.83284751e-01] 79 [549.36345191 -2.00643568 -0.89843479 -6.42642194 -9.03248322 -40.9011983 0.99789862] 80 [707.06934244 -1.0927203 -3.1429159 -10.25158047 -7.10734989 -31.91068171 1.1206625 ] 81 [ 1.93099444e+03 -7.54508761e+00 -2.29304190e+00 -2.91887465e+00 -9.49879835e+00 -3.89991736e+01 9.86830024e-01] 82 [ 1.00006836e+03 -5.30850008e-01 -1.16325179e+01 -2.73393891e+00 -1.14231153e-01 -2.58727988e+01 7.67083378e-01] 83 [ 1.27466339e+03 -4.52476492e+00 -7.16538598e+00 -4.50146181e+00 -7.51899937e+00 -3.98783528e+01 9.95409380e-01] 84 [ 1.17495339e+03 -8.67964117e-01 -2.15621221e+00 -3.37681511e+00 -4.88659558e+00 -2.85248133e+01 9.88307929e-01] 85 [ 9.44160729e+02 -2.93622147e+00 -2.95117949e+00 -2.83663110e+00 -4.11166859e-01 -3.98863303e+01 9.95497175e-01] 86 [472.09647758 -2.82882197 -6.57769123 -10.77711351 -5.36812743 -34.6519138 0.99757951] 87 [ 1.68236250e+03 -3.67979781e+00 -4.59147596e+00 -4.94401031e+00 -9.47472124e+00 -2.52797020e+01 9.82666578e-01] 88 [507.99157278 -2.50125097 -1.74889638 -9.14064063 -7.00947424 -38.13952529 0.99744258] 89 [663.9753842 -3.69460511 -11.12795466 -11.29099284 -5.76646563 -32.57089053 0.99748999] 90 [ 1.82916307e+03 -3.04551919e+00 -2.30273181e-01 -8.51045259e+00 -3.08469960e+00 -2.59871150e+01 1.05649768e+00] 91 [714.38153447 -4.41513219 -4.20117479 -11.72258717 -3.96243134 -30.15590113 0.99570508] 92 [589.87053658 -3.9668154 -2.29559607 -5.85600973 -8.47098078 -31.54555671 0.99811899] 93 [549.66110667 -10.29066507 -0.63523998 -10.84891131 -11.29640057 -39.85535035 0.99960348] 94 [480.00797662 -1.8320824 -11.33331189 -8.20110613 -3.63501587 -28.60168764 0.98413626] 95 [324.4299068 -1.56683328 -5.07007751 -6.1818635 -8.78515644 -43.3283287 0.99752478] 96 [ 1.13819136e+03 -1.93826622e+00 -6.82162647e-01 -1.06598552e+01 -5.43040134e+00 -3.75085003e+01 9.96466605e-01] 97 [ 7.32531088e+02 -4.69805102e+00 -1.01457765e+01 -2.77653537e+00 -1.80825029e-01 -2.97244617e+01 1.03238211e+00] 98 [ 1.21490270e+03 -4.49396630e-01 -5.45025103e+00 -3.37258925e+00 -7.89757464e+00 -2.82436996e+01 9.92260672e-01] 99 [ 1.93099444e+03 -7.54508761e+00 -2.29304190e+00 -2.91887465e+00 -9.49879835e+00 -3.89991736e+01 9.86830024e-01] 100 [723.92878378 -5.98012126 -2.07539992 -3.77962021 -1.04864221 -36.23495548 0.98352848] 101 [384.36243191 -7.99010859 -9.76891402 -5.05715357 -7.16979065 -29.47070163 0.99649515] 102 [517.47221116 -5.32191458 -7.40063459 -5.91349076 -3.56896507 -33.20769601 0.99721735] 103 [277.03685255 -2.40951139 -8.02255345 -5.53766102 -8.08268564 -33.95896793 1.01893731] 104 [1290.61615945 -4.49567194 -3.49830287 -3.12245421 -1.32252428 -32.50893578 1.30338009] 105 [252.69800595 -9.80386678 -4.95840951 -1.01781134 -4.90434706 -29.49983596 0.83657504] 106 [ 1.75315051e+03 -4.84839784e+00 -5.96278257e+00 -1.08977754e+01 -4.68227153e+00 -4.17680431e+01 1.19074811e+00] 107 [505.79132392 -1.0669323 -11.21923471 -10.39816383 -7.52204872 -27.54931004 0.99566936] 108 [467.0641429 -10.0537644 -7.29668874 -4.68311203 -2.35997241 -35.34403246 0.99157526] 109 [845.35764594 -4.19447531 -2.51437659 -6.8269608 -5.09246797 -43.26223268 0.99381253] 110 [592.62063078 -4.480328 -7.55286249 -10.45813654 -8.50315458 -37.66973254 0.99830344] 111 [543.82530916 -7.91670364 -11.83161518 -8.27917836 -3.92859644 -32.85864792 0.99777205] 112 [401.66697035 -7.48142048 -2.06306457 -9.20689621 -8.40694992 -37.06506474 0.92607075] 113 [742.44166584 -7.174359 -2.86784376 -8.17422657 -1.20031942 -44.78628615 0.9944324 ] 114 [824.21830433 -2.25738175 -3.62676317 -5.93357336 -10.57570069 -28.05177362 0.93159811] 115 [329.97491653 -7.85454269 -2.51301687 -11.45865891 -2.80360703 -32.58808704 0.99866553] 116 [504.18823467 -2.30925313 -1.86159405 -10.39386036 -6.05530795 -34.67673753 0.99709908] 117 [ 1.04001359e+03 -7.02133044e+00 -8.66028655e+00 -6.99990436e+00 -3.78237954e+00 -3.10380549e+01 9.95209945e-01] 118 [274.878811 -5.49159608 -7.80260211 -7.94095614 -10.57848947 -44.42611228 0.9986575 ] 119 [ 1.85434286e+03 -5.80675876e+00 -3.49302829e+00 -1.00579096e+01 -1.10456588e+01 -2.72133778e+01 9.73993112e-01] 120 [ 2.12700933e+02 -1.15533454e+01 -3.32408429e+00 -1.08704809e+01 -2.29570048e-02 -3.88561504e+01 1.10571150e+00] 121 [616.73722133 -1.44988366 -9.5651106 -10.76334898 -7.36762056 -29.69487856 0.99540598] 122 [ 1.26190747e+03 -8.52495874e+00 -1.12508084e+01 -7.17102475e+00 -5.41888754e+00 -3.43799867e+01 9.96546131e-01] 123 [ 3.42674339e+02 -4.84605938e+00 -2.66552667e+00 -2.49613782e-02 -5.92266913e+00 -4.00549147e+01 1.01540558e+00] 124 [ 1.32815566e+03 -4.66755385e+00 -1.37972482e+00 -1.79751858e+00 -9.59906129e+00 -3.31054615e+01 1.30075289e+00] 125 [ 1.66878724e+03 -5.59516444e-01 -6.15049605e-01 -4.89473671e+00 -8.51801299e+00 -3.40176924e+01 9.95821571e-01] 126 [427.39347126 -2.3294651 -3.9381686 -8.4677459 -5.29419003 -33.85436662 0.99735851] 127 [ 1.66164519e+03 -1.76196792e+00 -8.91818567e+00 -3.71366915e+00 -1.02904728e+00 -4.11845807e+01 9.89147493e-01] 128 [309.51215206 -2.24928932 -3.74562776 -4.82626016 -11.61390816 -37.57516435 0.99803551] 129 [929.91996264 -9.11127767 -3.07292524 -7.36056115 -11.93592989 -34.71967812 0.99684782] 130 [225.93614774 -0.63072793 -3.84294233 -6.73678721 -7.25901266 -39.52477391 0.74954997] 131 [372.07400853 -1.58603928 -2.96290459 -6.42451177 -4.15041532 -41.05077774 0.99712273] 132 [ 1.25353932e+03 -4.66298734e+00 -1.33731736e+00 -9.11609645e+00 -2.62163026e+00 -2.64048398e+01 9.88904946e-01] 133 [714.38153447 -4.41513219 -4.20117479 -11.72258717 -3.96243134 -30.15590113 0.99570508] 134 [376.14173279 -6.41562482 -10.25712206 -11.47015286 -5.3705378 -40.45794546 0.99851952] 135 [723.16918069 -3.96420957 -11.40224639 -10.62313758 -8.31869718 -34.27349212 0.99758021] 136 [ 1.45758546e+03 -9.75473685e-01 -3.02985706e+00 -1.00558835e+01 -1.00818046e+01 -3.54342169e+01 9.90375211e-01] 137 [ 1.17673139e+03 -2.07530384e+00 -3.33357967e+00 -1.14190825e+01 -1.21615138e+00 -4.33384637e+01 7.41224494e-01] 138 [399.65247707 -8.98731321 -10.74693306 -10.8509333 -4.74247822 -33.28493137 1.06825078] 139 [187.86790309 -10.41210584 -0.26156714 -9.23466313 -5.89815998 -28.95824694 0.99450936] 140 [ 1.29643527e+03 -1.03079705e+01 -3.61834104e+00 -5.84271883e+00 -6.52747696e+00 -3.34376776e+01 9.85700498e-01] 141 [684.21684412 -1.07125899 -6.66860201 -7.85214555 -5.32023804 -44.51424445 0.99607943] 142 [525.88275627 -4.10360344 -4.77289523 -11.81855138 -6.1955703 -28.29449619 1.23478294] 143 [806.34636603 -1.45533478 -6.43637648 -3.99610681 -10.36121109 -34.10420271 0.99485848] 144 [ 1.75763486e+03 -4.11762355e+00 -2.27644858e+00 -7.71936279e+00 -9.19368951e+00 -3.36658403e+01 9.88904484e-01] 145 [227.71073492 -4.82802508 -11.13927199 -3.01113616 -1.56963637 -32.52118676 0.99729119] 146 [611.76391405 -1.72985108 -2.70737756 -10.66774176 -9.45149266 -34.14743368 0.9964256 ] 147 [370.08608013 -7.20166122 -7.93366202 -10.64905342 -2.99573377 -30.38903391 0.99694476] 148 [435.91910423 -2.39268838 -8.99748425 -8.87949741 -9.50953487 -37.70272537 0.99747989] 149 [708.96854429 -3.84675351 -6.86046453 -7.03304311 -5.83184043 -35.08541469 0.99742446] 150 [731.22182962 -3.88932756 -4.12401321 -7.15806356 -6.77697136 -44.19648686 0.97911952] 151 [893.68786595 -5.2718812 -3.33047809 -8.57982631 -10.62998136 -39.34578019 0.99656875] 152 [ 1.41967056e+03 -7.42934481e+00 -7.23939415e+00 -5.18135151e+00 -3.41759650e+00 -3.82060691e+01 1.02861520e+00] 153 [285.30013164 -11.27565209 -9.13752201 -10.26820128 -11.72039643 -36.71720219 0.81985405] 154 [421.81338841 -1.0619625 -6.78755098 -10.09408492 -2.95011562 -34.89256195 0.99692304] 155 [424.80029856 -2.68803432 -4.64497489 -2.66613278 -4.69502418 -37.05522863 0.99407248] 156 [942.00438725 -8.04897622 -11.50460599 -6.30197288 -10.58005365 -31.07354105 0.99157664] 157 [462.60419746 -1.59340837 -3.53029685 -10.46682569 -11.98269177 -41.16842059 0.9970955 ] 158 [ 1.94870187e+03 -8.16197793e+00 -6.06529058e+00 -1.68342715e+00 -1.09352305e+01 -3.78017564e+01 1.04427576e+00] 159 [ 1.98161044e+03 -6.91662409e+00 -3.45100668e+00 -9.31970988e+00 -8.64622986e+00 -4.19671825e+01 9.91449217e-01] 160 [ 1.69845595e+03 -7.91949999e+00 -4.84282814e+00 -4.04841874e+00 -7.74757199e+00 -2.97043659e+01 9.77528806e-01] 161 [211.97627765 -1.91602787 -4.51324698 -11.70288115 -4.08820425 -41.48895278 0.99662221] 162 [ 1.82675975e+03 -1.03078464e+01 -2.23210695e+00 -4.77465431e+00 -9.26556372e+00 -3.61236145e+01 9.72294720e-01] 163 [ 1.18529767e+03 -1.26722416e+00 -5.62218433e+00 -1.00275960e+01 -1.69966765e+00 -3.89390338e+01 8.93473380e-01] 164 [ 1.18882828e+03 -7.94166727e+00 -2.55494079e+00 -1.11501401e+01 -1.30809997e+00 -2.80436040e+01 9.71377561e-01] 165 [ 1.31871785e+03 -6.08562283e+00 -1.77578095e+00 -1.02338700e+01 -2.44803440e+00 -4.33609447e+01 9.97453166e-01] 166 [613.56140416 -4.93404312 -6.20917777 -9.24214775 -7.03452405 -35.93451122 0.99830598] 167 [288.54167429 -10.68054012 -3.50574757 -8.27268432 -2.17571133 -29.95294376 0.99681462] 168 [793.92502936 -6.57775288 -5.62401133 -5.0614586 -2.83704419 -42.3973281 0.99199869] 169 [363.47945667 -2.56046613 -6.08008776 -6.22548494 -9.51685365 -36.50356106 0.99767906] 170 [374.10139399 -1.13101031 -1.5804855 -0.88146341 -6.88947447 -28.4715479 0.95864294] 171 [ 1.09969300e+03 -9.26688672e+00 -2.99470021e+00 -1.17842654e+01 -8.96271005e+00 -3.91346895e+01 9.98754837e-01] 172 [ 1.39860499e+03 -3.62462276e+00 -5.83543886e+00 -2.23094503e+00 -7.34364797e+00 -3.65364153e+01 1.01373488e+00] 173 [628.65754666 -5.66469177 -0.85848302 -10.36137032 -1.2223761 -38.54255188 0.99766491] 174 [ 1.46707718e+03 -1.94763162e+00 -7.88911857e+00 -4.93750948e+00 -6.77167502e-01 -3.09012225e+01 9.92229153e-01] 175 [ 6.76016827e+02 -1.41928832e+00 -4.98074204e+00 -4.63350641e+00 -2.30735501e-01 -2.91181576e+01 1.01021924e+00] 176 [ 1.88058533e+03 -1.16786451e+01 -3.68000189e+00 -9.55362697e+00 -1.01007117e+01 -3.03384063e+01 9.91660565e-01] 177 [ 1.67685074e+03 -1.15296050e-01 -1.12737827e+01 -1.17268370e+01 -1.00596580e+01 -2.95280446e+01 7.00795185e-01] 178 [555.88645591 -3.67816515 -4.34458858 -11.40369981 -9.88502078 -37.93316465 0.99744136] 179 [ 6.76016827e+02 -1.41928832e+00 -4.98074204e+00 -4.63350641e+00 -2.30735501e-01 -2.91181576e+01 1.01021924e+00] 180 [897.59054181 -11.74761962 -9.44859712 -4.82885772 -7.65295798 -35.48933105 1.00224107] 181 [498.32403828 -3.37535702 -4.69801518 -6.67444412 -8.98915273 -41.77407497 0.99784509] 182 [ 7.17442967e+02 -2.55740664e+00 -9.62730237e-02 -8.57541850e+00 -1.06912927e+01 -3.85018320e+01 9.98425877e-01] 183 [457.77682979 -10.35724017 -2.1962905 -10.17253292 -7.80377486 -33.941042 1.04731343] 184 [ 1.30565612e+03 -5.52893661e+00 -7.32190554e+00 -7.83423391e+00 -3.07916670e+00 -4.24369009e+01 9.69756547e-01] 185 [ 9.37038428e+02 -5.50288865e+00 -7.85261128e-01 -1.00506982e+01 -8.30701635e+00 -3.25236743e+01 9.98501722e-01] 186 [598.78906797 -3.36880698 -10.02287946 -11.78832783 -7.15314861 -42.73291791 0.99736335] 187 [576.74862978 -3.75544466 -10.63632668 -5.51763178 -4.53693764 -38.08275684 0.99690741] 188 [370.08608013 -7.20166122 -7.93366202 -10.64905342 -2.99573377 -30.38903391 0.99694476] 189 [529.02479988 -10.52489351 -9.88387358 -10.89304005 -5.08714164 -26.20205904 0.99179548] 190 [578.19741972 -4.3410244 -11.1255267 -7.52035404 -9.36443683 -34.12610018 0.99801632] 191 [466.55228705 -2.4845256 -7.82034787 -7.96345104 -8.6663657 -31.37724787 0.99715191] 192 [642.83285645 -3.90070551 -9.9038303 -6.47488648 -10.79976194 -31.82267344 0.99727807] 193 [ 1.98220506e+03 -4.52220737e+00 -6.57061964e+00 -4.97733016e+00 -1.15888641e+01 -4.45588211e+01 9.64265651e-01] 194 [537.80417963 -4.87984415 -1.93054569 -9.63912128 -10.25131243 -32.59256711 0.99835109] 195 [514.73967013 -2.7387374 -3.64050537 -10.30749826 -5.49753187 -30.49135963 0.99663202] 196 [ 9.44454675e+02 -8.62149107e+00 -1.94906430e+00 -5.61073942e+00 -3.51581562e+00 -3.66610493e+01 5.14860524e-01] 197 [659.07316934 -3.44443807 -6.5008343 -9.37039988 -11.10547428 -43.44929409 0.99816589] 198 [889.70720513 -7.08347555 -9.21791399 -2.94479604 -7.89546088 -36.66933409 0.98530948] 199 [463.19349173 -1.08833797 -9.60342889 -8.1576795 -7.90719 -41.16849806 0.99731812]
wlgrid, y_meas, err, y_binned_array, wno, y_mod_array=pickle.load(open('./OUTPUT/spectral_samples_trans_dyn_wfc3_free.pic','rb'))
y_median=np.zeros(wno.shape[0])
y_high_1sig=np.zeros(wno.shape[0])
y_high_2sig=np.zeros(wno.shape[0])
y_low_1sig=np.zeros(wno.shape[0])
y_low_2sig=np.zeros(wno.shape[0])
for i in range(wno.shape[0]):
percentiles=np.percentile(y_mod_array[i,:],[4.55, 15.9, 50, 84.1, 95.45])
y_low_2sig[i]=percentiles[0]
y_low_1sig[i]=percentiles[1]
y_median[i]=percentiles[2]
y_high_1sig[i]=percentiles[3]
y_high_2sig[i]=percentiles[4]
from matplotlib.pyplot import *
from matplotlib.ticker import FormatStrFormatter
ymin=np.min(y_meas)*1E2*0.995
ymax=np.max(y_meas)*1E2*1.005
fig1, ax=subplots()
xlabel('$\lambda$ ($\mu$m)',fontsize=14)
ylabel('(R$_{p}$/R$_{*}$)$^{2} \%$',fontsize=14)
minorticks_on()
#for i in range(20): plot(wlgrid, y_binned_array[:,i]*100.,alpha=0.5,color='red')
#for i in range(20): plot(1E4/wno, y_mod_array[:,i]*100.,alpha=0.5,color='red')
fill_between(1E4/wno[::-1],y_low_2sig[::-1]*100,y_high_2sig[::-1]*100,facecolor='r',alpha=0.5,edgecolor='None')
fill_between(1E4/wno[::-1],y_low_1sig[::-1]*100,y_high_1sig[::-1]*100,facecolor='r',alpha=1.,edgecolor='None')
errorbar(wlgrid, y_meas*100, yerr=err*100, xerr=None, fmt='Dk')
plot(1E4/wno, y_median*1E2)
ax.set_xscale('log')
ax.set_xticks([0.3, 0.5,0.8,1,1.4, 2, 3, 4, 5])
ax.axis([0.3,5.0,ymin,ymax])
ax.get_xaxis().set_major_formatter(matplotlib.ticker.ScalarFormatter())
ax.tick_params(length=5,width=1,labelsize='small',which='major')
savefig('./plots/dyn_transmission_spectrum_fits_WFC3_FREE.pdf',fmt='pdf')
show()
#close()
PyMultiNest (a python wrapper for MultiNest): Buchner et al. 2014
Paper: https://ui.adsabs.harvard.edu/abs/2014A%26A...564A.125B/abstract)
GitHub: https://johannesbuchner.github.io/PyMultiNest/
MultiNest: Feroz & Hobson 2008
Paper: https://ui.adsabs.harvard.edu/abs/2009MNRAS.398.1601F/abstract
GitHub:https://github.com/farhanferoz/MultiNest
It is highly advised that you read the relevant sections of these papers and online documentation
While dynesty is easy to install and use, it's not quite as an efficient as a sampler as MultiNest/PyMultiNest (at least that I have found). Here we will re-do the above retrieval, but with the PyMultiNest package. The MultiNest package is readily parallelizeble on any compute cluster or multicore desktop/laptop (using mpi). Be forewarned, installing PyMultiNest can be a challenge as it is a combination of python, c, and fortran codes. Most compute cluster admins are able to install this package (both within slurm and torque). Below I outline the steps that I took on my 2016 MacBook Pro with High Sierra (10.13.6) based on a combination of https://www.astrobetter.com/wiki/MultiNest+Installation+Notes and https://exoai.github.io/software/taurex/installation.
Assuming you have some version of MacPorts, execute the following in this order from a terminal command line (note this can take several hours):
OK, time to run. Sorry, you will have to leave the comfort and safety of this notebook, but I'll walk you through it here. We will be using "mpirun" to run on your multiple processors on your laptop which I have not figured out how to do from a jupyter notebook (maybe you can use spawn or whatever, if your ambitious).
Open the routine "call_pymultinest_transmission.py". Have a look, digest it, read the comments, etc. It should look farily similar to the above dynesty fundtions (e.g., a loglikelihood, prior cube, etc.). No need to modify anything for this. However, you would modify this script if you wanted to add/remove particular parameters or change the wavelength/wavenumber ranch of the run.
Go to the folder (in terminal) where you have downloaded the CHIMERA code (mainly, where "call_pymultinest_transmission.py" lives). To run using one CPU, simpley type the following:
python call_pymultinest_transmission_wfc3_FREE.py
Wait a minute. You should see some useless unimportant error messages pop out followed by "Cross-sections Loaded" followed shortely there after by
$*************************$
MultiNest v3.10
Copyright Farhan Feroz & Mike Hobson
Release Jul 2015
no. of live points = 500
dimensionality = 8
$*************************$
Starting MultiNest
generating live points
$------------------$
Eventually (a few minutes) you whould see more output:
live points generated, starting sampling
Do not be alarmed by the $*****$. It just means that there are too many digits to print out for the inital ln(Z) values. Eventually normalish looking numbers (usually negative) print out here. Ok, feel free to kill this now, because we will use mpi to make it faster.
To do mpi, if everything installed correctly, you just need to type the following into the same command line:
mpirun -np 4 python call_pymultinest_transmission_wfc3_FREE.py
Note, you can use more than 4 cpu's/cores if your computer has them. Use as many as you have. Now, just wait. You should see similar output as for the single cpu, ubt at a faster rate. For the default setup, with 4 cpu's, this takes $\sim$1hr 10 min (15452 calls). When it is done it will print out this:
ln(ev)= -25.249086782326863 +/- 0.12562433089471531
Total Likelihood Evaluations: 15452
Sampling finished. Exiting MultiNest
analysing data from ./pmn$_$transmission/template_.txt
analysing data from ./pmn$_$transmission/template_.txt
analysing data from ./pmn$_$transmission/template_.txt
analysing data from ./pmn$_$transmission/template_.txt
When this is done, it is time to plot. Open the routine, "plot_PMN_transmission.py". This will make the usual obnoxious corner plots, a reconstructed "TP" profile (which is kind of meaningless), and will generate sample spectra from random parameter vectors drawn from the posterior. This last part takes some time. When this is complete a nice 1-,2-sigma spectral spread plot will pop up. Feel free to manipulate this to make plots to your liking. Note, the spectral draws will be saved as a python pickle so that you can twiddle with plot adjustments without haveing to regenerate that each time.
That's it! You are done! So after having worked throught his notebook you should feel comfortable playing around with the forward model to gain an understanding of how each parameter influences the spectrum, comfortable running both the dynesty and pymultinest samplers, and plotting up standard output/results. Good luck. Note, there is no warrenty, if you get an unpysical answer, that's on you =)
Feel free to move onto the emission tutorial!