def read_system(self,star_path='demo_star.txt',planet_path='demo_planet.txt',obs_path='demo_observations.txt'): """Reads in the stellar, planet and observation parameters from file; performing tests on the input and lifting the read variables to the class object. Parameters ---------- star_path : str Path to parameter file defining the star. planet_path: str Path to the parameter file defining the planet and its orbit. obs_path: str Path to the parameter file defining the timestamps of the observations. """ planetparams = open(planet_path,'r').read().splitlines() starparams = open(star_path,'r').read().splitlines() obsparams = open(obs_path,'r').read().splitlines() self.velStar = float(starparams[0].split()[0]) self.stelinc = float(starparams[1].split()[0]) self.drr = float(starparams[2].split()[0]) self.T = float(starparams[3].split()[0]) self.Z = float(starparams[4].split()[0]) self.logg = float(starparams[5].split()[0]) self.u1 = float(starparams[6].split()[0]) self.u2 = float(starparams[7].split()[0]) self.mus = float(starparams[8].split()[0]) self.R = float(starparams[9].split()[0]) self.sma_Rs = float(planetparams[0].split()[0]) self.ecc = float(planetparams[1].split()[0]) self.omega = float(planetparams[2].split()[0]) self.orbinc = float(planetparams[3].split()[0]) self.pob = float(planetparams[4].split()[0])#Obliquity. self.Rp_Rs = float(planetparams[5].split()[0]) self.orb_p = float(planetparams[6].split()[0]) self.transitC = float(planetparams[7].split()[0]) self.mode = planetparams[8].split()[0]#This is a string. times = [] #These are in JD-24000000.0 or in orbital phase. self.exptimes = [] for i in obsparams: times.append(float(i.split()[0])) self.times = np.array(times) self.Nexp = len(self.times)#Number of exposures. if self.mode == 'times': for i in obsparams: self.exptimes.append(float(i.split()[1])) self.exptimes = np.array(self.exptimes) try: test.typetest(self.wave_start,float,varname='wave_start in input') test.nantest(self.wave_start,varname='wave_start in input') test.notnegativetest(self.wave_start,varname='wave_start in input') test.notnegativetest(self.velStar,varname='velStar in input') test.notnegativetest(self.stelinc,varname='stelinc in input') #add all the other input parameters except ValueError as err: print("Parser: ",err.args) if self.mus != 0: self.mus = np.linspace(0.0,1.0,self.mus)#Uncomment this to run in CLV mode with SPECTRUM.
def clip_spectrum(wl,fx,wlmin,wlmax,pad=0): """This crops a spectrum wl,fx to wlmin and wlmax. If pad !=0, it should be set to a positive velocity by which the output spectrum will be padded around wlmin,wlmax, which is returned as well as the narrower cropped spectrum. Parameters ---------- wl : np.ndarray() The stellar model wavelength(s) in nm. fx : np.ndarray() The stellar model flux. wlmin: float The minimum cropping wavelength. wlmax: float The maximum cropping wavelength. pad: float, m/s, optional A velocity corresponding to a shift around which the cropped array will be padded. Returns ------- wl,fx: np.array(), np.array() The wavelength and flux of the integrated spectrum. """ import lib.test as test import numpy as np #First run standard tests on input test.typetest(wl,np.ndarray,varname='wl in clip_spectrum') test.typetest(fx,np.ndarray,varname='fx in clip_spectrum') test.typetest(wlmin,[int,float],varname='wlmin in clip_spectrum') test.typetest(wlmax,[int,float],varname='wlmax in clip_spectrum') test.typetest(pad,[int,float],varname='pad in clip_spectrum') test.notnegativetest(pad,varname='pad in clip_spectrum') test.nantest(wlmin,varname='wlmin in clip_spectrum') test.nantest(wlmax,varname='wlmax in clip_spectrum') test.nantest(pad,varname='pad in clip_spectrum') if wlmin >= wlmax: raise ValueError('wlmin in clip_spectrum should be smaller than wlmax.') wlc = wl[(wl >= wlmin) & (wl <= wlmax)]#This is the wavelength grid onto which we will interpolate the final result. fxc = fx[(wl >= wlmin) & (wl <= wlmax)] if pad > 0: wlmin_wide = wlmin/doppler(pad) wlmax_wide = wlmax*doppler(pad) wlc_wide = wl[(wl >= wlmin_wide) & (wl <= wlmax_wide)] fxc_wide = fx[(wl >= wlmin_wide) & (wl <= wlmax_wide)] return(wlc,fxc,wlc_wide,fxc_wide) else: return(wlc,fxc)
def build_spectrum_limb_resolved(wl, fx_list, mu_list, wlmin, wlmax, x, y, vel_grid): """WRITE THIS. Parameters ---------- """ #I copy paste as much as possible from above. The roles of wlc and wlc_wide have #changed because wl,fx is already narrow by construction. Therefore in order to #crop the spectrum with margin, I actually need to crop the wl axis inwards. #To do this, I needed to convert clip_spectrum to crop_spectrum. import numpy as np import lib.operations as ops import lib.stellar_spectrum as spectrum import sys import lib.test as test wlc_wide = wl #Standard tests on input test.typetest(mu_list, np.ndarray, varname='mu_list in build_spectrum_limb_resolved') test.dimtest(mu_list, [len(fx_list)], varname='mu_list in build_spectrum_limb_resolved') input_tests_global(wlc_wide, fx_list[0], wlmin, wlmax, x, y, vel_grid, vel_grid, fname='build_spectrum_limb_resolved') wlc, fxc = ops.crop_spectrum( wl, fx_list[0], 1.0 * np.nanmax(np.abs(vel_grid)) ) #I do this for only one spectrum because I only care about wlc F = 0 #output # start = time.time() for i in range(len(x)): for j in range(len(y)): if np.isnan(vel_grid[j, i]) == False: mu = 1 - np.sqrt(x[i]**2 + y[j]**2) diff = abs(mu_list - mu) index = np.argmin(diff) if mu_list[index] > 0: fxc_wide = fx_list[index] # print(i,j,x[i],y[j],mu,mu_list[index]) F += ops.shift(wlc, wlc_wide, fxc_wide, vel_grid[j, i]) #*flux_grid[j,i] statusbar(i, len(x)) # print(time.time()-start) return (wlc, F)
def vactoair(wlnm): """Convert vaccuum to air wavelengths. Parameters ---------- wlnm : float, np.array() The wavelength(s) in nm. Returns ------- wl: float, np.array() The wavelengths. """ import lib.test as test import numpy #First run standard tests on the input test.typetest(wlnm,[int,float,numpy.ndarray],varname='wlnm in vactoair') test.notnegativetest(wlnm,varname='wlnm in vactoair') test.nantest(wlnm,varname='wlnm in vactoair') wlA = wlnm*10.0 s = 1e4/wlA f = 1.0 + 5.792105e-2/(238.0185e0 - s**2) + 1.67917e-3/( 57.362e0 - s**2) return(wlA/f/10.0)
def airtovac(wlnm): """Convert air to vaccuum wavelengths. Parameters ---------- wlnm : float, np.array() The wavelength(s) in nm. Returns ------- wl: float, np.array() The wavelengths. """ import lib.test as test import numpy #First run standard tests on the input test.typetest(wlnm,[int,float,numpy.ndarray],varname='wlnm in airtovac') test.notnegativetest(wlnm,varname='wlnm in airtovac') test.nantest(wlnm,varname='wlnm in airtovac') #Would still need to test that wl is in a physical range. wlA=wlnm*10.0 s = 1e4 / wlA n = 1 + 0.00008336624212083 + 0.02408926869968 / (130.1065924522 - s**2) + 0.0001599740894897 / (38.92568793293 - s**2) return(wlA*n/10.0)
def input_tests_local(xp, yp, RpRs): """This wraps input tests for the local integrator functions, as these are the same in both cases.""" import lib.test as test #xp and yp should be ints or floats, and may not be NaN. test.typetest(xp, [int, float], varname='xp in build_local_spectrum_fast') test.typetest(yp, [int, float], varname='xp in build_local_spectrum_fast') test.nantest(xp, varname='xp in build_local_spectrum_fast') test.nantest(yp, varname='yp in build_local_spectrum_fast') #RpRs should be a non-negative float, and may not be NaN: test.typetest(RpRs, float, varname='RpRs in build_local_spectrum_fast') test.nantest(RpRs, varname='RpRs in build_local_spectrum_fast') test.notnegativetest(RpRs, varname='RpRs in build_local_spectrum_fast')
def smooth(fx,w,mode='gaussian',edge_degree=1): """This function takes a spectrum, and blurs it using either a Gaussian kernel or a box kernel, which have a FWHM width of w px everywhere. Meaning that the width changes dynamically on a constant d-lambda grid. """ import numpy as np import lib.test as test import pdb test.typetest(w,float,'w') test.typetest(fx,np.ndarray,'fx') test.typetest(mode,str,'mode') test.typetest(edge_degree,int,'edge_degree') truncsize=8.0#The gaussian is truncated at 8 sigma. shape=np.shape(fx) sig_w = w / 2*np.sqrt(2.0*np.log(2)) #Transform FWHM to Gaussian sigma. In km/s. trunc_dist=np.round(sig_w*truncsize).astype(int) #First define the kernel. kw=int(np.round(truncsize*sig_w*2.0)) if kw % 2.0 != 1.0:#This is to make sure that the kernel has an odd number of #elements, and that it is symmetric around zero. kw+=1 kx=findgen(kw) kx-=np.mean(kx)#This must be centered around zero. Doing a hardcoded check: if (-1.0)*kx[-1] != kx[0]: print(kx) raise Exception("ERROR in box_smooth: Kernel could not be made symmetric somehow. Attempted kernel grid is printed above. Kernel width is %s pixels." % kw) if mode == 'gaussian': k=gaussian(kx,1.0,0.0,sig_w) k/=np.sum(k) return(convolve(fx,k,edge_degree))
def calc_flux_stellar(x,y,u1,u2): import numpy as np import lib.operations as ops import lib.test as test test.typetest(x,np.ndarray,varname='x in calc_flux_stellar') test.typetest(y,np.ndarray,varname='x in calc_flux_stellar') test.typetest(u1,[int,float],varname='i_stellar in calc_flux_stellar') test.typetest(u2,[int,float],varname='vel_eq in calc_flux_stellar') test.nantest(x,varname='x in calc_flux_stellar') test.nantest(y,varname='y in calc_flux_stellar') test.nantest(u1,varname='x in calc_flux_stellar') test.nantest(u2,varname='y in calc_flux_stellar') z,x_full,y_full = calc_z(x,y) flux_grid = z*0.0 for i in range(len(x)): for j in range(len(y)): if np.sqrt(x[i]**2+y[j]**2) <= 1.0: flux_grid[j,i]=ops.limb_darkening((np.sqrt(x[i]**2+y[j]**2)),u1,u2)*(z[j,i]*0.0+1.0)#Last multiplication is to put NaNs back into place. # plotting.plot_star_2D(x,y,mu_grid,cmap="hot",quantities=['','',''],units=['','',''],noshow=False) flux_grid /= np.nansum(flux_grid) return(flux_grid)
def crop_spectrum(wl,fx,pad): """This crops a spectrum wl,fx to wlmin and wlmax. If pad !=0, it should be set to a positive velocity by which the output spectrum will be padded within wlmin,wlmax, to allow for velocity shifts. The difference with clip_spectrum above is that this routine pads towards the inside (only returning the narrow spectrum), while clip_spectrum pads towards the outside, returning both the cropped and the padded cropped spectra. Parameters ---------- wl : np.ndarray() The stellar model wavelength(s) in nm. fx : np.ndarray() The stellar model flux. pad: float, m/s A velocity corresponding to a shift around which the cropped array will be padded. Returns ------- wl,fx: np.array(), np.array() The wavelength and flux of the integrated spectrum. """ import lib.test as test import numpy as np wlmin = min(wl) wlmax = max(wl) #First run standard tests on input test.typetest(wl,np.ndarray,varname= 'wl in crop_spectrum') test.typetest(fx,np.ndarray,varname= 'fx in crop_spectrum') test.typetest(pad,[int,float],varname= 'pad in crop_spectrum') test.notnegativetest(pad,varname= 'pad in crop_spectrum') test.nantest(pad,varname= 'pad in crop_spectrum') if wlmin >= wlmax: raise ValueError('wlmin in crop_spectrum should be smaller than wlmax.') wlmin_wide = wlmin*doppler(pad)#Crop towards the inside wlmax_wide = wlmax/doppler(pad)#Crop towards the inside wlc_narrow = wl[(wl >= wlmin_wide) & (wl <= wlmax_wide)] fxc_narrow = fx[(wl >= wlmin_wide) & (wl <= wlmax_wide)] return(wlc_narrow,fxc_narrow)
def limb_darkening(z,a1,a2): """Evaluate quadratic limb darkening, taken straight from wikipedia, eq 1. (https://en.wikipedia.org/wiki/Limb_darkening). Provide the z-coordinate of grid cell i,j and multiply the resulting weight against spectrum i,j. z = sin(psi), with psi the angle of incidence. psi = arcsin(z) I(z) = I(0) * (a0 + a1*cos(psi) + a2*cos(psi)**2 + ....) with a0 = 1 - (a1+a2+...). Parameters ---------- z : float, np.array() The projected distance from the stellar center in units of stellar radii, bounded between 0 and 1. u1,u2: float, Linear and quadratic limb-darkening coefficients. Returns ------- w: float The weight, <= 1. """ import lib.test as test #First run standard tests on the input test.typetest(z,[int,float],varname='z in limb_darkening') test.typetest(a1,[int,float],varname='a1 in limb_darkening') test.typetest(a2,[int,float],varname='a2 in limb_darkening') test.nantest(z,varname='z in limb_darkening') test.nantest(a1,varname='a1 in limb_darkening') test.nantest(a2,varname='a2 in limb_darkening') import numpy as np if z > 1 or z < 0: raise ValueError("z coordinate should be in [0,1].") psi = np.arcsin(z) a0 = 1-a1-a2 I = a0+a1*np.cos(psi)+a2*np.cos(psi)**2 return(I)
def get_spectrum(T, logg, Z, a): """Querying a PHOENIX photosphere model, either from disk or from the PHOENIX website if the spectrum hasn't been downloaded before, and returning the names of the files in the data subfolder. These may then be read by the wrapper function read_spectrum() below. However, if limbs is set to a positive value, the code is switched to a mode in which the CLV effect is taken into account. In this event, PHOENIX spectra cannot be used, and the SPECTRUM package will be used to calculate spectra. For this functionality to work, the user needs to have built the SPECTRUM source code (included in this distribution). See the documentation for instructions. Parameters ---------- T : int, float The model photosphere temperature. Acceptable values are: 2300 - 7000 in steps of 100, and 7200 - 12000 in steps of 200. logg : int, float The model log(g) value. Acceptable values are: 0.0 - 6.0 in steps of 0.5. Z : int, float The model metallicity [Fe/H] value. Acceptable values are: -4.0, -3.0, -2.0 and -1.5 to +1.0 in steps of 0.5. a : int, float The model alpha element enhancement [alpha/M]. Acceptable values are: -0.2 to 1.2 in steps of 0.2, but only for Fe/H of -3.0 to 0.0. Returns ------- wlname, specname: str, str The names of the files containing the wavelength and flux axes of the requested spectrum. """ import requests import shutil import urllib.request as request from contextlib import closing import sys import os.path import lib.test as test #First run standard tests on the input test.typetest(T, [int, float], varname='T in get_spectrum') test.typetest(logg, [int, float], varname='logg in get_spectrum') test.typetest(Z, [int, float], varname='metallicity in get_spectrum') test.typetest(a, [int, float], varname='alpha in get_spectrum') test.nantest(T, varname='T in get_spectrum') test.nantest(logg, varname='logg in get_spectrum') test.nantest(Z, varname='Z in get_spectrum') test.nantest(a, varname='a in get_spectrum') test.notnegativetest(T, varname='T in get_spectrum') #This is where phoenix spectra are located. root = 'ftp://phoenix.astro.physik.uni-goettingen.de/HiResFITS/' #We assemble a combination of strings to parse the user input into the URL, z_string = '{:.1f}'.format(float(Z)) if Z > 0: z_string = '+' + z_string else: z_string = '-' + z_string a_string = '' if a > 0: a_string = '.Alpha=+' + '{:.2f}'.format(float(a)) if a < 0: a_string = '.Alpha=' + '{:.2f}'.format(float(a)) t_string = str(int(T)) if T < 10000: t_string = '0' + t_string g_string = '-' + '{:.2f}'.format(float(logg)) #These are URLS for the input files. waveurl = root + 'WAVE_PHOENIX-ACES-AGSS-COND-2011.fits' specurl = root + 'PHOENIX-ACES-AGSS-COND-2011/Z' + z_string + a_string + '/lte' + t_string + g_string + z_string + a_string + '.PHOENIX-ACES-AGSS-COND-2011-HiRes.fits' #These are the output filenames, they will also be returned so that the wrapper #of this function can take them in. wavename = 'data/PHOENIX/WAVE.fits' if os.path.isdir('data/PHOENIX/') == False: os.makedirs('data/PHOENIX/') specname = 'data/PHOENIX/lte' + t_string + g_string + z_string + a_string + '.PHOENIX-ACES-AGSS-COND-2011-HiRes.fits' #If it doesn't already exists, we download them, otherwise, we just pass them on: if os.path.exists(wavename) == False: with closing(request.urlopen(waveurl)) as r: with open(wavename, 'wb') as f: shutil.copyfileobj(r, f) if os.path.exists(specname) == False: print(specurl) with closing(request.urlopen(specurl)) as r: with open(specname, 'wb') as f: shutil.copyfileobj(r, f) return (wavename, specname)
def blur_spec(wl,spec,dv,truncsize = 20.0): """This function takes a spectrum, and blurs it using either a Gaussian kernel or a box kernel, which have a FWHM width of dv km/s everywhere. Meaning that the width changes dynamically on a constant d-lambda grid. Because the kernel needs to be recomputed on each element of the wavelength axis individually, this operation is much slower than convolution with a constant kernel, in which a simple shifting of the array, rather than a recomputation of the kernel is sufficient.""" # print('I MAY NOT WANT TO USE BLUR-SPEC BECAUSE IT IS SLOW, AT LEAST IN BOX MODE.') # print('AND I HAVE NOT THOROUGHLY BENCHMARKED IT.') import numpy as np import pdb from matplotlib import pyplot as plt import time import astropy.constants as const import numpy as np import lib.test as test # Do not perform tests because this thing is in a double forloop. # test.typetest(dv,float,varname='dv in doppler(dv)') c = const.c.value#Comes out in m/s. test.typetest(dv,float,varname='dv in blur_spec') test.typetest(wl,np.ndarray,varname='w in blur_specl') test.typetest(spec,np.ndarray,varname='spec in blur_spec') test.typetest(truncsize,float,varname='truncsize in blur_spec') #truncsize=8.0#The gaussian is truncated at 8 sigma. sig_dv = dv / 2*np.sqrt(2.0*np.log(2)) #Transform FWHM to Gaussian sigma. In km/s. d_kernel=np.array([-1,0,1])/2.0 deriv = convolve(wl,d_kernel) #l*dv/c=dl dwl=wl*dv/const.c.value*1000.0 sig_wl=wl*sig_dv/const.c.value*1000.0#in nm sig_px=sig_wl/deriv trunc_dist=np.round(sig_px*truncsize).astype(int) spec_blurred=spec*0.0 summm = [] # pdb.set_trace() for i in range(0,len(wl)): #Im going to select wl in a bin so that I dont need to evaluate a gaussian over millions of points that are all zero binstart=max([0,i-trunc_dist[i]]) binend=i+trunc_dist[i] k = gaussian(wl[binstart:binend],1.0,wl[i],sig_wl[i]) k_n=k/np.sum(k) summm.append(np.sum(k)) # try: # print(len(k_n),len(spec[binstart:binend]),len(wl[binstart:binend])) try: spec_blurred[i]=np.sum(k_n*spec[binstart:binend]) except: print(len(wl),len(spec)) pdb.set_trace() # except: # pdb.set_trace() #To speed up, need to select wl and then append with zeroes. <= what does that mean? Jens 03 mar 18 # plt.plot(wl,summm) # plt.show() return(spec_blurred)
def read_spectrum(T, logg, metallicity=0.0, alpha=0.0): """Wrapper for the function get_spectrum() above, that checks that the input T, log(g), metallicity and alpha are in accordance with what is provided by PHOENIX (as of November 1st, 2019), and passes them on to get_spectrum(). Parameters ---------- T : int,float The model photosphere temperature. Acceptable values are: 2300 - 7000 in steps of 100, and 7200 - 12000 in steps of 200. logg : int,float The model log(g) value. Acceptable values are: 0.0 - 6.0 in steps of 0.5. metallicity : int,float (optional, default = 0) The model metallicity [Fe/H] value. Acceptable values are: -4.0, -3.0, -2.0 and -1.5 to +1.0 in steps of 0.5. If no location is given, the ``location`` attribute of the Time object is used. alpha : int,float (optional, default = 0) The model alpha element enhancement [alpha/M]. Acceptable values are: -0.2 to 1.2 in steps of 0.2, but only for Fe/H of -3.0 to 0.0. Returns ------- wl,f : np.array(),np.array() The wavelength (nm) and flux (erg/s/cm^2/cm) axes of the requested spectrum. """ import numpy as np import sys import astropy.io.fits as fits import lib.test as test phoenix_factor = np.pi #This is a factor to correct the PHOENIX spectra to produce the correctly calibrated output flux. #If you compare PHOENIX to SPECTRUM, a factor of 3 seems to be missing. Brett Morris encountered this problem when trying #to use PHOENIX to predict the received flux of real sources, to construct an ETC. He also needed a factor of 3, which he #interpreted as a missing factor of pi. So pi it is, until something else is deemed better (or until PHOENIX updates #their grid). #First run standard tests on the input: test.typetest(T, [int, float], varname='T in read_spectrum') test.typetest(logg, [int, float], varname='logg in read_spectrum') test.typetest(metallicity, [int, float], varname='metallicity in read_spectrum') test.typetest(alpha, [int, float], varname='alpha in read_spectrum') test.nantest(T, varname='T in get_spectrum') test.nantest(logg, varname='logg in get_spectrum') test.nantest(metallicity, varname='Z in get_spectrum') test.nantest(alpha, varname='a in get_spectrum') test.notnegativetest(T, varname='T in get_spectrum') #These contain the acceptable values. T_a = np.concatenate((np.arange(2300, 7100, 100), np.arange(7200, 12200, 200))) logg_a = np.arange(0, 6.5, 0.5) FeH_a = np.concatenate((np.arange(-4, -1, 1), np.arange(-1.5, 1.5, 0.5))) alpha_a = np.arange(0, 1.6, 0.2) - 0.2 #Check that the input is contained in them: if T in T_a and metallicity in FeH_a and alpha in alpha_a and logg in logg_a: if (metallicity < -3 or metallicity > 0.0) and alpha != 0: print( 'Error: Alpha element enhancement != 0 only available for -3.0<[Fe/H]<0.0' ) sys.exit() #If so, retrieve the spectra: wavename, specname = get_spectrum(T, logg, metallicity, alpha) f = fits.getdata(specname) w = fits.getdata(wavename) return (w / 10.0, f / phoenix_factor) #Implicit unit conversion here. else: print( 'Error: Provided combination of T, log(g), Fe/H and a/M is out of bounds.' ) print('T = %s, log(g) = %s, Fe/H = %s, a/M = %s' % (T, logg, metallicity, alpha)) print('The following values are accepted:') print('T:') print(T_a) print('Log(g):') print(logg_a) print('Fe/H:') print(FeH_a) print('a/M:') print(alpha_a) print( 'Alpha element enhancement != 0 only available for -3.0<[Fe/H]<0.0' ) sys.exit()
def input_tests_global(wl, fx, wlmin, wlmax, x, y, vel_grid, flux_grid, fname=''): """This wraps the input tests for all the integrator functions, as these are all the same.""" import lib.test as test import numpy as np import lib.operations as ops #wl and fx should be numpy arrays with the same length. test.typetest(wl, np.ndarray, varname='wl in ' + fname) test.typetest(fx, np.ndarray, varname='fx in ' + fname) test.dimtest(fx, [len(wl)], varname='fx in ' + fname) #wlmin and wlmax should be floats or ints, and wlmin < wlmax: test.typetest(wlmin, [int, float], varname='wlmin in ' + fname) test.typetest(wlmax, [int, float], varname='wlmax in ' + fname) if wlmax <= wlmin: raise ValueError('Wlmax (%s) in ' + fname + ' should be smaller than wlmin (%s)' % (wlmax, wlmin)) #x and y should be numpy arrays, contain no NaNs and their dimensions should match #the velocity and flux grids, which should be 2D numpy arrays. test.typetest(x, np.ndarray, varname='x in ' + fname) test.typetest(y, np.ndarray, varname='y in ' + fname) test.nantest(x, varname='x in ' + fname) test.nantest(y, varname='x in ' + fname) test.dimtest(vel_grid, [len(y), len(x)], varname='y in ' + fname) test.dimtest(flux_grid, [len(y), len(x)], varname='y in ' + fname) test.typetest(vel_grid, np.ndarray, varname='vel_grid in ' + fname) test.typetest(flux_grid, np.ndarray, varname='flux_grid in ' + fname) dvm = np.nanmax(np.abs(vel_grid)) if wlmin / ops.doppler(dvm) <= np.nanmin(wl): raise Exception( "Value error in " + fname + ": wlmin (%s) is smaller than min(wl) (%s) after accounting for the maximum velocity shift (%s)" % (wlmin / ops.doppler(dvm), min(wl), dvm)) if wlmax * ops.doppler(dvm) >= np.nanmax(wl): raise Exception( "Value error in " + fname + ": wlmax (%s) is greater than max(wl) (%s) after accounting for the maximum velocity shift (%s)" % (wlmax * ops.doppler(dvm), max(wl), dvm)) #fname should be a string. test.typetest(fname, str, varname='fname in input_tests_global')
def build_local_spectrum_limb_resolved(xp, yp, zp, RpRs, wl, fx_list, mu_list, wlmin, wlmax, x, y, vel_grid): """This is used when SPECTRUM was called to provide spectra wth mu angles, in which case the flux map doesn't exist and the computation is different from the disk-averaged PHOENIX spectra. WRITE THIS MORE Parameters ---------- xp : float The x position of the planet in units of x (see below). yp : float The y position of the planet in units of y (see below). RpRs : float The radius of the planet in units of stellar radii. wl : np.array() The stellar model wavelength(s) in nm. fx : np.array() The stellar model flux. wlmin: float The minimum wavelength to be considered, in units of wl. wlmax: float The maximum wavelength to be considered, in units of wl. x,y: The x and y axis of the velocity grid, typically in units of stellar radius. vel_grid: 2D np.array() The velocity grid of the stellar disk. Values outside of the disk are to be set to NaN. Returns ------- wl,fx: np.array(), np.array() The wavelength and flux of the integrated spectrum, flux: float The lightcurve of the remainder of the stellar disk. mask: The inverse mask that shows the stellar disk with the planet masked out. """ import numpy as np import lib.operations as ops import lib.test as test import warnings import sys import pdb wlc_wide = wl #Standard tests on input: input_tests_local(xp, yp, RpRs) input_tests_global(wl, fx_list[0], wlmin, wlmax, x, y, vel_grid, vel_grid, fname='build_local_spectrum_limb_resolved') test.typetest(mu_list, np.ndarray, varname='mu_list in build_local_spectrum_limb_resolved') test.dimtest(mu_list, [len(fx_list)], varname='mu_list in build_local_spectrum_limb_resolved') #We start by creating a mask that selects just the area of the star that is #covered by the planet, as well as its inverse which we like to return for #plotting purposes. x_full = np.tile(x, (len(y), 1)) - xp y_full = np.tile(y, (len(x), 1)).T - yp d = np.sqrt(x_full**2 + y_full**2) if zp < 0: d += np.inf di = d * 1.0 d[d > RpRs] = np.nan #Mask out everything but the location of the planet. with warnings.catch_warnings(): warnings.simplefilter("ignore", category=RuntimeWarning) if np.nanmin(d) <= RpRs: di[d <= RpRs] = np.nan #out the location of the planet. Inverse of d. if zp < 0.0: d *= np.nan #Force the planet to not be in front of the disk if the z-coordinate is such that it is behind the star. di[np.isnan(di)] = 1.0 #started as 1.0, set back to 1.0 mask = (0.0 * vel_grid + 1.0) * ( d * 0.0 + 1.0 ) #Set that to 1.0 and multiply with flux grid. Nansum coming! mask_i = (0.0 * vel_grid + 1.0) * (di * 0.0 + 1.0) #Inverse of mask. #MOVE ALL OF THIS INTO A WRAPPER? I THINK THIS IS NOW COPYPASTED 3 TIMES? #ACTUALLY, NOW I AM USING VEL GRID INS TEAD OF FLUX GRID. THAT MEANS THAT THE #FLUX IS NO LONGER CALCULATED PROPERLY. NO LIGHTCURVE! # wlc,fxc,wlc_wide,fxc_wide = ops.clip_spectrum(wl,fx,wlmin,wlmax,pad=2.0*np.nanmax(np.abs(vel_grid))) wlc, fxc = ops.crop_spectrum( wl, fx_list[0], 1.0 * np.nanmax(np.abs(vel_grid)) ) #I do this for only one spectrum because I only care about wlc F = 0 #output # start = time.time() # pdb.set_trace() for i in range(len(x)): for j in range(len(y)): if np.isnan(mask[j, i]) == False: mu = 1 - np.sqrt(x[i]**2 + y[j]**2) diff = np.abs(mu_list - mu) index = np.argmin(diff) if mu_list[index] > 0: fxc_wide = fx_list[index] # print(i,j,x[i],y[j],mu,mu_list[index]) F += ops.shift(wlc, wlc_wide, fxc_wide, vel_grid[j, i]) #*flux_grid[j,i] # print(time.time()-start) return (wlc, F, np.nansum(mask_i), (di * 0.0) + 1.0)
def calc_vel_stellar(x,y,i_stellar, vel_eq, diff_rot_rate, proj_obliquity): """ based on Cegla+2016. See Fig. 3 and equations 2 - 8 https://arxiv.org/pdf/1602.00322.pdf It takes the stellar parameters and then calculates the stellar velocities in all bins for one quarter of the stellar disk input: x, y: 1D numpy arrays to create the stellar grid in units of stellar radius i_stellar: inclination in degrees vel_eq: equatorial stellar velocity diff_rot_rate: differential rotation rate proj_obliquity: projected obliquity output: vel_stellar_grid: 2D numpy array of stellar velocities over one quarter of the stellar disk """ import numpy as np import pdb import lib.test as test #Standard tests on input test.typetest(x,np.ndarray,varname='x in calc_vel_stellar') test.typetest(y,np.ndarray,varname='x in calc_vel_stellar') test.typetest(i_stellar,float,varname='i_stellar in calc_vel_stellar') test.typetest(vel_eq,float,varname='vel_eq in calc_vel_stellar') test.typetest(diff_rot_rate,float,varname='diff_rot_rate in calc_vel_stellar') test.typetest(proj_obliquity,float,varname='proj_obliquity in calc_vel_stellar') test.nantest(x,varname='x in calc_vel_stellar') test.nantest(y,varname='y in calc_vel_stellar') test.nantest(i_stellar,varname='i_stellar in calc_vel_stellar') test.nantest(vel_eq,varname='vel_eq in calc_vel_stellar') test.nantest(diff_rot_rate,varname='diff_rot_rate in calc_vel_stellar') test.nantest(proj_obliquity,varname='proj_obliquity in calc_vel_stellar') test.notnegativetest(vel_eq,varname='vel_eq in calc_vel_stellar') #Careful! all angles have to be in radians! #Think carefully about which ones of these have to be transposed #convert angles to rad alpha = np.radians(proj_obliquity) beta = np.radians(i_stellar) #pre calculate matrices z,x_full,y_full = calc_z(x,y) #equation 8 from Cegla+2016 # vel_stellar_grid = (x_full*np.cos(alpha)-y_full*np.sin(alpha))*vel_eq*np.sin(beta)*(1.-diff_rot_rate*(z*np.sin(np.pi/2.-beta)+np.cos(np.pi/2.-beta)*(x_full*np.sin(alpha)-y_full*np.cos(alpha)))) #this has the star aligned to the coordinate system. Like this the projected obliquity is not incorporated here. The tilt of the star compared to the planet has to be incorporated in the planet's path. vel_stellar_grid = x_full*vel_eq*np.sin(beta)*(1.-diff_rot_rate*(z*np.sin(np.pi/2.-beta)+np.cos(np.pi/2.-beta)*y_full)) return(vel_stellar_grid)
def compute_spectrum(T, logg, Z, mu, wlmin, wlmax, macroturbulence=0.0, mode='an', loud=False): """If the CLV effect needs to be taken into account, PHOENIX spectra cannot be used, and the SPECTRUM package will be used to calculate spectra. For this functionality to work, the user needs to have built the SPECTRUM source code (included in this distribution). See the documentation for instructions on how to do this. This code then acts as a wrapper for SPECTRUM. SPECTRUM needs to be made in the lib/SPECTRUM folder, in the same place where it is provided along with this distribution. This distribution also provides the necessary linelist files, as well as grids of stellar atmosphere files obtained from Kurucz. After spectrum is made, this should work out of the box. Provide a temperature, metallicity and logg, as well as a mu value or an np.array of mu values for the wrapper to loop over. The accepted values are different from those when using the disk-integrated PHOENIX models. T_a = np.concatenate((np.arange(3500,13250,250),np.arange(14000,1000,51000))) logg_a = np.arange(0,5.5,0.5) Z_a = np.round(np.concatenate((np.arange(-5,0,0.5),np.arange(-0.3,0.4,0.1),np.array([0.5,1.0]))),decimals=1)#I need to do a rounding here because np.concatenate() makes a numerical error on the middle array... Weird! Parameters ---------- T : int, float The model photosphere temperature. Acceptable values are: 3500 - 13,000 in steps of 250, and 14,000 - 50,000 in steps of 1,000. logg : int, float The model log(g) value. Acceptable values are: 0.0 - 5.0 in steps of 0.5. Z : int, float The model metallicity [Fe/H] value. Acceptable values are: -5.0 - -0.5 in steps of 0.5; -0.3 - +0.3 in steps of 0.1; 0.5; 1.0. wlmin: float The minimum wavelength to be considered, in units of wl. wlmax: float The maximum wavelength to be considered, in units of wl. mu: float,np.array Limb angle(s) (cos(theta)) macroturbulence: int,float Macroturbulent velocity (which broadens the lines) in km/s mode: str a: Prompt user for custom statoms.dat i: Isotope mode. For use with lukeall2.iso.lst linelist. n: Silent running. m/M: Specific intensities, normalised (m) or true (M). f: Gives output in specific intensities. Implied if m or M is set. Returns ------- wl,fx The wavelengths (nm) and fluxes of the computed spectrum. If mu has more than 1 element, then fx is a list of flux axes, each of which matches to wl, and each of which matches to its respective mu value. """ import subprocess import numpy as np import lib.test as test import os import pdb import astropy.io.ascii as ascii from lib.integrate import statusbar import lib.operations as ops from os import path import sys kurucz_root = 'data/KURUCZ/' #THIS IS THE SECOND TIME THIS IS DEFINED. SHOULD MATCH THE DEFINITION IN test_KURUCZ()! spath = './lib/SPECTRUM/selectmod' #This points to the supermod function. SPECTRUM = './lib/SPECTRUM/spectrum' #This points to spectrum. test.dir_exists('./lib/SPECTRUM/', varname='SPECTRUM root folder') test.dir_exists(kurucz_root, varname='Kurucz model folder') test.file_exists(spath, varname=spath + 'in compute_spectrum()') test.file_exists( SPECTRUM, varname=spath + 'in compute_spectrum()' ) #These four tests need to be moved into some dedicated test routine for SPECTRUM? #Hardcoded not to use isotope mode because itneeds the massive lukeall2.iso.lst that I dont want to upload to github. #mode = 'ainMf' #Mode of calling SPECTRUM. if 'i' in mode: isotope = True else: isotope = False if (('m' in mode) or ('M' in mode)) == False: mu = 0.0 if 'a' in mode == False: print("WARNING: Compute_spectrum() mode does not include 'a'.") print("This will break SPECTRUM as it expects the stdatom file to be") print( "in the root folder. 'a' points it to where it is supposed to be.") print("Adding 'a' to mode variable.") mode += 'a' # if c_norm == True: # mode = mode.lower()#Lower the case of the M. # if isotope == False: # mode = mode.replace('i','') # if disk_integrated == True: # mode = mode.replace('m','').replace('M','') # mu = 0.0 #Standard tests on input. test.typetest(T, [int, float], varname='T in compute_spectrum') test.typetest(logg, [int, float], varname='log(g) in compute_spectrum') test.typetest(Z, [int, float], varname='Metallicity in compute_spectrum') test.typetest(mu, [float, np.ndarray], varname='mu angle(s) in compute_spectrum') test.typetest(wlmin, [int, float], varname='wlmin in compute_spectrum') test.typetest(wlmax, [int, float], varname='wlmax in compute_spectrum') test.typetest(macroturbulence, [int, float], varname='macroturbulence in compute_spectrum') test.nantest(mu, 'mu angle(s) in compute_spectrum') if wlmax <= wlmin: raise ValueError( 'Wlmax (%s) in compute_spectrum() should be smaller than wlmin (%s)' % (wlmax, wlmin)) if wlmin < 90: raise ValueError( 'Wlmin (%s) in compute_spectrum() should be larger than 90 nm.' % wlmin) if wlmax > 4000: raise ValueError( 'Wlmax (%s) in compute_spectrum() should be less than than 4,000 nm.' % wlmin) #These hard-code the acceptable values. T_a, logg_a, Z_a = available_kurucz_models() if T in T_a and Z in Z_a and logg in logg_a: #This builds the Kurucz model atmosphere file, extracted from a 'supermodel' #grid. This follows the manual of SPECTRUM, section 4.10. #We first construct the filename of the supermodel, using the provided metallicity. mpath = construct_kurucz_path(Z, kurucz_root) test.file_exists(mpath, varname='Kurucz supermodel') #Check that this supermodel actually exists. It should, if test_KURUCZ() was run (which tests all allowed values of Z). opath = kurucz_root + 'star.mod' #This will contain the temporary .mod file that is the model extracted from the Kurucz supermodel file. rpath = 'data/SPECTRUM.res' #This will contain the temporary response file used to < into the bash call of SPECTRUM. if isotope == True: lpath = 'lib/SPECTRUM/lukeall2.iso.lst' else: lpath = 'lib/SPECTRUM/merged_linelist.lst' #I made this linelist by manually mergin the other linelists. Will be part of our SPECTRUM fork. test.file_exists(lpath, varname=lpath + ' (linelist file) in compute_spectrum()') sopath = 'data/SPECTRUM.spc' #Spectrum output path, temporary. try: res = subprocess.check_output( [spath, mpath, opath, str(T), str(logg)] ) #This executes the command line as >>> selectmod supermodel.dat output.mod teff logg #So now we have made an atmosphere model file that SPECTRUM can use. except: print('Unexpected file input error.') print('Apparently the provided combination of T, log(g), Fe/H is') print('not included in this (%s) Kurucz model grid.' % mpath) print('T = %s, log(g) = %s, Fe/H = %s' % (T, logg, Z)) if isinstance(mu, float): #Convert into a list so that we can loop. mu = [mu] fx = [] #This will contain the flux axes, even if there is only one. bashcmd = SPECTRUM + ' ' + mode + ' < ' + rpath if loud == False: bashcmd += ' > temp' #Write to temp to suppress output. i = 0 #Just a counter for the statusbar. print('------ Executing bash cmd: ' + bashcmd) for m in mu: #Loop over all mu angles, keeping the atmosphere model constant. statusbar(i, mu) lines = [] lines.append(opath) lines.append(lpath) lines.append('lib/SPECTRUM/stdatom.dat') if isotope == True: lines.append('lib/SPECTRUM/isotope.iso') lines.append(sopath) lines.append(str(macroturbulence)) if ('m' in mode) or ('M' in mode): lines.append(str(m)) lines.append(str(wlmin * 10.0) + ',' + str(wlmax * 10.0)) #Convert to angstroms lines.append( '0.01') #Hardcoded to 0.01A, i.e. R = 500,000 at 500nm. # line1 = opath#Input .mod file. # line3 = 'lib/SPECTRUM/stdatom.dat' # line5 = 'data/SPECTRUM.spc'#The output file. # line6 = str(macroturbulence)#Macroturbulence, default 2kms. # line7 = str(m) # line8 = str(wlmin*10.0)+','+str(wlmax*10.0)#Convert to angstroms # line9 = '0.01' # if isotope == True:#This is hardcoded to be off, at present. # line2 = 'lib/SPECTRUM/lukeall2.iso.lst'#The linelist. # line4 = 'lib/SPECTRUM/isotope.iso'#The isotope file. # lines = [line1,line2,line3,line4,line5,line6,line7,line8,line9] # else: # line2 = 'lib/SPECTRUM/merged_linelist.lst' # lines = [line1,line2,line3,line5,line6,line7,line8,line9]#Load non-isotope linelist and ingore isotope file on line 4. with open(rpath, 'w') as f: for line in lines: f.write( line + '\n' ) #This works because the last line of the response #file must end in a carriage return. void = os.system( bashcmd ) #I use os.system() against the advice of teh internet because subprocesses.run() hangs forever. # print('Finished computing SPECTRUM') # with open('temp','w') as f: # process=subprocess.Popen(bashcmd.split(),stdout=subprocess.PIPE,shell=False) # stdout = process.communicate() test.file_exists(sopath, varname='SPECTRUM output (' + sopath + ')') d = ascii.read(sopath) wl = d.columns[ 0].data #This is overwritten each time, but that doesn't matter because each time its the same. fx.append(d.columns[1].data) #Now we remove all the scratch files we made in this loop. # os.remove('temp') os.remove(sopath) os.remove(rpath) i += 1 os.remove(opath) if len(mu) == 1: return (ops.airtovac(wl / 10.0), fx[0] ) #Convert back to nm and do airtovac. else: return (ops.airtovac(wl / 10.0), fx ) #Convert back to nm and do airtovac. else: print( 'Error: Provided combination of T, log(g), Fe/H is out of bounds.') print('T = %s, log(g) = %s, Fe/H = %s' % (T, logg, Z)) print('The following values are accepted:') print('T:') print(T_a) print('Log(g):') print(logg_a) print('Fe/H:') print(list(Z_a)) sys.exit()