Welcome to the emission/secondary-eclipse spectrum model/retrieval tutorial!
For this particular setup, the atmosphere is parameterized within the "chemically-consistent" framework as described in Kreidberg et al. 2015, but for emission. The atmospheric composition is parameterized with only a metalicity and carbon-to-oxygen ratio assumed to be in thermochemical equilibrium along the temperature-pressure profile. Originally this code would compute the gas and condensate phase mixing ratios by calling the NASA CEA routine. However, in order to remove fortran dependencies, a finely sampled pre-computed, interpolateable chemistry grid was instead produced with CEA as a function of temperature (T from 400K - 3400K in 100K increments), pressure ($log_{10}(P)$ from -7.0 (0.1$\mu$bar) - 2.4 (316 bar) in 0.1 increments), metallicity ($[M/H]$ from -2.0 (0.01$\times$) to 3.0 (1000$\times$)), and C/O ($log_{10}(C/O)$ from -2.0 (0.01) to 0.3 (2) awkwardly spaced to better sample the transition about C/O=1). All elemental abundances are scaled with respect to the Lodders 2009 solar abundance pattern. A pseudo-hack rainout approximation is made to the gas phase abundances of TiO, VO, Na, K, and FeH. In this hack, these species are set to 0 abundance at levels above where they first fall below some critical value ($10^{-8}$). This is to mimic the loss of these species from the gas phase into the condensate phase. In no case are we accounting for the loss of elemental abundances.
The 3-parameter temperature profile parameterization utilizes the Guillot 2010/Parmentier et al. 2014 analytic formulism (see Line et al. 2013a for implementation details).
This "emission" routine accounts for both stellar reflected and planetary thermal emission treated within the Toon et al. 1989 two stream + two stream source function technique (toonpy.py and toonpy_solar.py). It accounts for multiple scattering given basic cloud optical properties (extinction efficiencies, single scatter albedo, and asymmetry parameter). 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). Specifically, in multi-gas, multiple scattering atmospheres in "emission", the "resort-rebin" gas-mixing procedure has to be implemented (Molliere et al. 2015, 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. We include as correlated-K line opacites H2O, CH4, CO, CO2, NH3, HCN, H2S, C2H2, Na, K, TiO, VO, FeH and as continuum gas opacities H2-H2, H2-He CIA, and the H- bound free and free free (e.g., Arcangeli et al. 2018). See the "opacity" tutorial for more details on correlated-K.
To handle the effects of disequilibrium chemistry due to vertical mixing, we apply the "quench-pressure" approximation. We include a quench pressure parameter for the carbon-system and one for the nitrogen system (as in Morley et al. 2017 for GJ436b, and Kreidberg et al. 2018 for WASP-107b). The carbon quench pressure fixes the H2O, CO, and CH4 abundances above the quench pressure level to their abundances at the quench pressure level. Similarly, the nitrogen quench pressure fixes the N2, NH3, and HCN abundances above the quench pressure to their values at the quench pressure level. This is indeed a kludge, and a better implementation would be to use the timescale/eddy mixing prescription described in Zahnle & Marley 2015. Regardless, any non-full kinetics approach is a kludge anyway (not to mention the 1D nature of the problem...).
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.
Finally, if one doesn't like the "chemically-consistent" concept, they can use the "gas_scale" array to switch off or scale the abundances each opacity source.
This specific notebook goes through the steps to generate the forward model, and illustrate how to actually perform the 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. The WFC3+STIS correlated-K coefficients (xsects_HST) 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 that the "core" set of routines are all in fm.py. The thermal emission radiative transfer solver is toonpy.py and the incident stellar flux solver is toonpy_solar.py If you want to know more about what is in the sausage, look into these routines.
#import all of the functions in fm, namely, the CK-coefficients (may take a minute)
from chimera import *
%matplotlib notebook
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 = 25000
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 emission 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--for reflected light component
#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
#Composition parameters---assumes "chemically consistent model" described in Kreidberg et al. 2015
logMet=0.0 #. #Metallicity relative to solar log--solar is 0, 10x=1, 0.1x = -1: valid range is -1.5 - 3.0
logCtoO=-0.26 #log C-to-O ratio: log solar is -0.26: valid range is -1.0 - 0.3
logPQCarbon=-5.5 #CH4, CO, H2O Qunech pressure--forces CH4, CO, and H2O to constant value at quench pressure value: valid range -6.0 - 1.5
logPQNitrogen=-5.5 #N2, NH3 Quench pressure--forces N2 and NH3 to ""
#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, logMet, logCtoO, logPQCarbon,logPQNitrogen, Rp, Rstar, M, D, 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([1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1., 1., 1.]) #can be made free params if desired (won't affect mmw)#can be made free params if desired (won't affect mmw)
#WASP43b emission spectrum from kreidberg et al. 2014...a 3 column ascii file
wlgrid, y_meas, err=np.loadtxt('w43b_emis.txt').T
Here we call the forward model routine "fx_emis" (think F(x)) from fm.py. fx_emis controls the input values and calls the relevent functions to compute the secondary eclipse spectrum. The inputs into fx_emis 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 ($F_p/F_{\star}$) 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,wnocrop,atm,Ftoa,Fstar,Fstar_TOA,Fup_therm,Fup_ref=fx_emis(x,wlgrid,gas_scale, xsecs)
print('SPECTRUM GENERATE')
SPECTRUM GENERATE
Spaghetti plot of the model atmosphere.
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
Pavg=0.5*(P[1:]+P[:-1])
fig2, ax1=subplots()
#feel free to plot whatever you want here....
ax1.semilogx(H2O,Pavg,'b',ls='--',lw=2,label='H2O')
ax1.semilogx(CH4,Pavg,'black',ls='--',lw=2,label='CH4')
ax1.semilogx(CO,Pavg,'g',ls='--',lw=2,label='CO')
ax1.semilogx(CO2,Pavg,'orange',ls='--',lw=2,label='CO2')
ax1.semilogx(NH3,Pavg,'darkblue',ls='--',lw=2,label='NH3')
ax1.semilogx(Na,Pavg,'b',lw=2,label='Na')
ax1.semilogx(K,Pavg,'g',lw=2,label='K')
ax1.semilogx(TiO,Pavg,'k',lw=2,label='TiO')
ax1.semilogx(VO,Pavg,'orange',lw=2,label='VO')
ax1.semilogx(qc,Pavg,'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_emission_WFC3_CC.pdf',fmt='pdf')
show()
#close()
%matplotlib notebook
ymin=0
ymax=np.max(y_mod)*1E3*1.2
fig1, ax=subplots()
xlabel('$\lambda$ ($\mu$m)',fontsize=14)
ylabel('F$_p$/F$_{star}$ [$\\times 10^{-3}$]',fontsize=14)
minorticks_on()
errorbar(wlgrid, y_meas*1E3, yerr=err*1E3, xerr=None, fmt='Dk')
plot(wlgrid, y_binned*1E3,'ob')
plot(1E4/wnocrop, y_mod*1E3,color='black',label='Total')
#reflected component
plot(1E4/wnocrop, Fup_ref/Fstar*1E3*(Rp/Rstar*0.10279)**2 ,color='blue',label='Reflected Stellar')
#emission component
plot(1E4/wnocrop, Fup_therm/Fstar*1E3*(Rp/Rstar*0.10279)**2 ,color='red',label='Thermal Emission ')
legend(frameon=False)
ax.set_xscale('log')
#ax.set_yscale('log')
ax.set_xticks([0.3, 0.5,0.8,1, 2, 3, 5])
ax.axis([0.3,5,ymin,ymax])
ax.get_xaxis().set_major_formatter(matplotlib.ticker.ScalarFormatter())
ax.tick_params(length=5,width=1,labelsize='large',which='major')
savefig('./plots/emission_spectrum_WFC3_CC.pdf',fmt='pdf')
show()
from matplotlib.pyplot import *
from matplotlib.ticker import FormatStrFormatter
#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--for reflected light component
#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
#Composition parameters---assumes "chemically consistent model" described in Kreidberg et al. 2015
logMet=0.0 #. #Metallicity relative to solar log--solar is 0, 10x=1, 0.1x = -1: valid range is -1.5 - 3.0
logCtoO=-0.26 #log C-to-O ratio: log solar is -0.26: valid range is -1.0 - 0.3
logPQCarbon=-5.5 #CH4, CO, H2O Qunech pressure--forces CH4, CO, and H2O to constant value at quench pressure value: valid range -6.0 - 1.5
logPQNitrogen=-5.5 #N2, NH3 Quench pressure--forces N2 and NH3 to ""
#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=-25.5 #cloud condensate base mixing ratio (e.g, see Fortney 2005)--valid range: -15 - -2.0
#for second model
logKzz2=7 #log Kzz (cm2/s)--valid range: 2 - 11 -- higher values make larger particles
fsed2=0.5 #sediminetation efficiency--valid range: 0.5 - 5--lower values make "puffier" more extended cloud
logPbase2=1.0 #cloud base pressure--valid range: -6.0 - 1.5
logCldVMR2=-3. #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, logMet, logCtoO, logPQCarbon,logPQNitrogen, Rp, Rstar, M, D, logKzz, fsed,logPbase,logCldVMR, logKcld, logRayAmp, RaySlope])
x2=np.array([Tirr, logKir,logg1,Tint, logMet, logCtoO, logPQCarbon,logPQNitrogen, Rp, Rstar, M, D, logKzz2, fsed2,logPbase2,logCldVMR2, 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([1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1., 1., 1.]) #can be made free params if desired (won't affect mmw)#can be made free params if desired (won't affect mmw)
y_binned,y_mod,wnocrop,atm,Ftoa,Fstar,Fstar_TOA,Fup_therm,Fup_ref=fx_emis(x,wlgrid,gas_scale, xsecs)
print('SPECTRUM GENERATE1')
y_binned2,y_mod2,wnocrop2,atm2,Ftoa2,Fstar2,Fstar_TOA2,Fup_therm2,Fup_ref2=fx_emis(x2,wlgrid,gas_scale, xsecs)
print('SPECTRUM GENERATE2')
ymin=0
ymax=np.max(y_mod)*1E3*1.2
fig1, ax=subplots()
xlabel('$\lambda$ ($\mu$m)',fontsize=14)
ylabel('F$_p$/F$_{star}$ [$\\times 10^{-3}$]',fontsize=14)
minorticks_on()
errorbar(wlgrid, y_meas*1E3, yerr=err*1E3, xerr=None, fmt='Dk')
plot(wlgrid, y_binned*1E3,'ob')
plot(1E4/wnocrop, y_mod*1E3,color='black',label='Total')
#reflected component
plot(1E4/wnocrop, Fup_ref/Fstar*1E3*(Rp/Rstar*0.10279)**2 ,color='blue',label='Reflected Stellar')
#emission component
plot(1E4/wnocrop, Fup_therm/Fstar*1E3*(Rp/Rstar*0.10279)**2 ,color='red',label='Thermal Emission ')
plot(1E4/wnocrop, y_mod2*1E3,color='black',label='Total 2',ls='--')
#reflected component
plot(1E4/wnocrop, Fup_ref2/Fstar2*1E3*(Rp/Rstar*0.10279)**2 ,color='blue',label='Reflected Stellar 2',ls='--')
#emission component
plot(1E4/wnocrop, Fup_therm2/Fstar2*1E3*(Rp/Rstar*0.10279)**2 ,color='red',label='Thermal Emission 2 ',ls='--')
legend(frameon=False)
ax.set_xscale('log')
#ax.set_yscale('log')
ax.set_xticks([0.3, 0.5,0.8,1, 2, 3, 5])
ax.axis([0.3,5,ymin,ymax])
ax.get_xaxis().set_major_formatter(matplotlib.ticker.ScalarFormatter())
ax.tick_params(length=5,width=1,labelsize='large',which='major')
show()
SPECTRUM GENERATE1 SPECTRUM GENERATE2
Here, I'm assuming you've done the transmission tutorial and have installed the requisite packages for pymultinest. This is very similar.
Open the routine "call_pymultinest_emission.py". Have a look, digest it, read the comments, etc. It should look farily similar to the call_pymultinest_emission.py (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 you "call_pymultinest_emission.py" lives). To run using one CPU, simpley type the following:
python call_pymultinest_emission_wfc3_CC.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 = 9
$*************************$
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_emission_wfc3_CC.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$4.4 hrs (38068 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$_$emission/template_.txt
analysing data from ./pmn$_$emission/template_.txt
analysing data from ./pmn$_$emission/template_.txt
analysing data from ./pmn$_$emission/template_.txt
When this is done, it is time to plot. Open the routine, "plot_PMN_emission.py". This will make the usual obnoxious corner plots, a reconstructed "TP" profile, 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 emission spectrum, comfortable running both the dynesty and pymultinest samplers, and plotting up standard output/results. Good luck. Note, there is no warranty, if you get an unpysical answer, that's on you =)
Happy Secondary Eclipse Retrieving!