.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/end_to_end/plot_spotted_transit.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_examples_end_to_end_plot_spotted_transit.py: Observe a transit of a spotted star. ======================================== This example demonstrates stellar contamination of a transit spectrum. A transit can change the spot coverage of a star and produce a signal that is difficult to distinguish from atmospheric absorption. We aim to simulate data from :cite:t:`2023ApJ...948L..11M`. .. GENERATED FROM PYTHON SOURCE LINES 11-28 .. code-block:: Python from pathlib import Path import numpy as np import matplotlib.pyplot as plt from astropy import units as u import libpypsg from VSPEC import ObservationModel,PhaseAnalyzer from VSPEC import params from VSPEC.config import MSH SEED = 10 PRE_TRANSIT = 8 IN_TRANSIT = 13 libpypsg.docker.set_url_and_run() .. rst-class:: sphx-glr-script-out .. code-block:: none Saved settings to /home/runner/.libpypsg/settings.json Reloading settings... .. GENERATED FROM PYTHON SOURCE LINES 29-34 Create the needed configurations -------------------------------- :cite:t:`2023ApJ...948L..11M` observed super-Earth GJ 486b with JWST NIRSpec/G395H .. GENERATED FROM PYTHON SOURCE LINES 34-252 .. code-block:: Python # Instrument inst = params.InstrumentParameters( telescope=params.SingleDishParameters.jwst(), bandpass=params.BandpassParameters( wl_blue=2.87*u.um, wl_red=5.14*u.um, resolving_power=200, wavelength_unit=u.um, flux_unit=u.Unit('W m-2 um-1') ), detector=params.DetectorParameters( beam_width=5*u.arcsec, integration_time=0.5*u.s, ccd=params.ccdParameters( pixel_sampling=8, read_noise=16.8*u.electron, dark_current=0.005*u.electron/u.s, throughput=0.3, emissivity=0.1, temperature=50*u.K ) ) ) # Observation observation = params.ObservationParameters( observation_time=3.53*u.hr, integration_time=8*u.min ) # PSG psg_kwargs = dict( gcm_binning=200, phase_binning=1, nmax=0, lmax=0, use_continuum_stellar=True, continuum=['Rayleigh', 'Refraction'], ) psg_params = params.psgParameters( use_molecular_signatures=True, **psg_kwargs ) psg_no_atm = params.psgParameters( use_molecular_signatures=False, **psg_kwargs ) # Star and Planet star_teff = 3343*u.K star_rad = 0.339*u.R_sun a_rstar = 11.229 # Moran+23 Table 1 rp_rstar = 0.03709 # Moran+23 Table 1 planet_rad = star_rad * rp_rstar orbit_rad = star_rad * a_rstar orbit_period = 1.467119*u.day planet_rot_period = orbit_period star_rot_period = np.inf*u.s # assume the star does not change. planet_mass = 2.82*u.M_earth star_mass = 0.323*u.M_sun inclination = 89.06*u.deg planet_norm_factor = { 'rock_quiet': 1, 'rock_spotted': 0.97, 'water_quiet': 0.98, } # We will multiply the radius by these scalars rather than fit # later (because I have run this code before and know how much # each model is off by). print(f'Planet radius: {planet_rad.to_value(u.R_earth)} R_earth') print(f'Semimajor axis: {orbit_rad.to_value(u.AU)} AU') observation_angle = (2*np.pi*u.rad * observation.observation_time/orbit_period).to(u.deg) initial_phase = 180*u.deg - 0.5*observation_angle planet_kwargs = dict( name='GJ486b', gravity=params.GravityParameters('kg',planet_mass), semimajor_axis=orbit_rad, orbit_period=orbit_period, rotation_period=planet_rot_period, eccentricity=0, obliquity=0*u.deg, obliquity_direction=0*u.deg, init_phase=initial_phase, init_substellar_lon=0*u.deg ) quiet_rock_planet = params.PlanetParameters( radius=planet_rad*planet_norm_factor['rock_quiet'], **planet_kwargs ) quiet_water_planet = params.PlanetParameters( radius=planet_rad*planet_norm_factor['water_quiet'], **planet_kwargs ) spotted_rock_planet = params.PlanetParameters( radius=planet_rad*planet_norm_factor['rock_spotted'], **planet_kwargs ) system_params = params.SystemParameters( distance=8.07*u.pc, inclination=inclination, phase_of_periastron=0*u.deg ) star_dict = { 'teff': star_teff, 'radius': star_rad } planet_dict = {'semimajor_axis': orbit_rad} gcm_dict = { 'nlayer': 40, 'nlon': 60, 'nlat': 30, 'epsilon': 100, 'albedo': 0.3, 'emissivity': 1.0, 'lat_redistribution': 0.0, 'gamma': 1.0, 'psurf': 1*u.bar, 'ptop': 1e-10*u.bar, 'wind': {'U': '0 m/s','V':'0 m/s'}, 'molecules':{'H2O':0.99} } # Create two sets of GCM Parameters h2o_atm = {'molecules':{'H2O':0.99}} gcm_h2o = params.gcmParameters.from_dict({ 'star':star_dict, 'planet':planet_dict, 'gcm':{'vspec':dict(gcm_dict,**h2o_atm),'mean_molec_weight':18} }) star_kwargs = dict( psg_star_template='M', teff=star_teff, mass=star_mass, radius=star_rad, period=star_rot_period, misalignment=0*u.deg, misalignment_dir=0*u.deg, ld=params.LimbDarkeningParameters.lambertian(), faculae=params.FaculaParameters.none(), flares=params.FlareParameters.none(), granulation=params.GranulationParameters.none(), grid_params=100000, ) quiet_star = params.StarParameters( spots=params.SpotParameters.none(), **star_kwargs ) spotted_star = params.StarParameters( spots=params.SpotParameters( distribution='iso', initial_coverage=0.11, area_mean=500*MSH, area_logsigma=0.2, teff_umbra=2700*u.K, teff_penumbra=2700*u.K, equillibrium_coverage=0.0, burn_in=0*u.s, growth_rate=0.0/u.day, decay_rate=0*MSH/u.day, initial_area=10*MSH ), **star_kwargs ) # Set parameters for simulation header_kwargs = dict( seed = SEED, spec_grid = params.VSPECGridParameters( max_teff=3400*u.K,min_teff=2300*u.K, impl_bin='rust',impl_interp='scipy',fail_on_missing=False ) ) internal_params_kwargs = dict( system=system_params, obs=observation, gcm=gcm_h2o, inst=inst ) # Make the three cases params_rock_quiet = params.InternalParameters( header=params.Header(data_path=Path('.vspec/rock_quiet'),**header_kwargs), star = quiet_star, psg=psg_no_atm, planet=quiet_rock_planet, **internal_params_kwargs ) params_h2o_quiet = params.InternalParameters( header=params.Header(data_path=Path('.vspec/h2o_quiet'),**header_kwargs), star = quiet_star, psg=psg_params, planet=quiet_water_planet, **internal_params_kwargs ) params_rock_spotted = params.InternalParameters( header=params.Header(data_path=Path('.vspec/rock_spotted'),**header_kwargs), star = spotted_star, psg=psg_no_atm, planet=spotted_rock_planet, **internal_params_kwargs ) .. rst-class:: sphx-glr-script-out .. code-block:: none Planet radius: 1.371472837835719 R_earth Semimajor axis: 0.017702612840063636 AU .. GENERATED FROM PYTHON SOURCE LINES 253-257 Run VSPEC for the simplest case ------------------------------- We read in the config file and run the model. .. GENERATED FROM PYTHON SOURCE LINES 257-262 .. code-block:: Python model_rock_quiet = ObservationModel(params_rock_quiet) model_rock_quiet.build_planet() model_rock_quiet.build_spectra() .. rst-class:: sphx-glr-script-out .. code-block:: none Starting at phase 161.9544290544939 deg, observe for 3.53 h in 26 steps Phases = [161.95 163.32 164.68 166.04 167.41 168.77 170.13 171.5 172.86 174.22 175.59 176.95 178.31 179.68 181.04 182.4 183.77 185.13 186.49 187.86 189.22 190.58 191.95 193.31 194.67 196.03 197.4 ] deg Build Planet: 0%| | 0/27 [00:00 .. GENERATED FROM PYTHON SOURCE LINES 402-406 Make the figure again +++++++++++++++++++++ This time without the bare rock .. GENERATED FROM PYTHON SOURCE LINES 406-429 .. code-block:: Python fig, ax = plt.subplots(1,1,figsize=(5,3.5),tight_layout=True) for data,label,color in zip( [data_h2o_quiet,data_rock_spotted], ['H2O', 'Rock+Spots'], ['xkcd:azure','xkcd:golden yellow'] ): pre_transit = PRE_TRANSIT in_transit = IN_TRANSIT unocculted_spectrum = data.spectrum('total',pre_transit) occulted_spectrum = data.spectrum('total',in_transit) lost_to_transit = unocculted_spectrum-occulted_spectrum transit_depth = (lost_to_transit/unocculted_spectrum).to_value(u.dimensionless_unscaled) ax.plot(data.wavelength,transit_depth*1e6,label=label,color=color) ax.errorbar(moran_x,moran_y,yerr=moran_yerr, fmt='o',color='k',label='Moran+23 Eureka!') ax.set_xlabel('Wavelength (um)') ax.set_ylabel('Transit depth (ppm)') ax.legend() .. image-sg:: /auto_examples/end_to_end/images/sphx_glr_plot_spotted_transit_005.png :alt: plot spotted transit :srcset: /auto_examples/end_to_end/images/sphx_glr_plot_spotted_transit_005.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. rst-class:: sphx-glr-timing **Total running time of the script:** (2 minutes 40.511 seconds) .. _sphx_glr_download_auto_examples_end_to_end_plot_spotted_transit.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_spotted_transit.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_spotted_transit.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: plot_spotted_transit.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_