Skip to content

Commit

Permalink
Add in check that hierarchical parent in simbad still passes DLR<5 cut
Browse files Browse the repository at this point in the history
  • Loading branch information
alexandergagliano committed Feb 10, 2024
1 parent 6165e67 commit 1bf4814
Showing 1 changed file with 19 additions and 7 deletions.
26 changes: 19 additions & 7 deletions astro_ghost/ghostHelperFunctions.py
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@
from astro_ghost.starSeparation import separateStars_STRM, separateStars_South
from astro_ghost.sourceCleaning import clean_dict, removePS1Duplicates, getColors, makeCuts
from astro_ghost.stellarLocus import calc_7DCD
from astro_ghost.DLR import chooseByDLR, chooseByGladeDLR
from astro_ghost.DLR import chooseByDLR, chooseByGladeDLR, calc_DLR
import importlib_resources
import requests
import pickle
Expand Down Expand Up @@ -151,6 +151,8 @@ def checkSimbadHierarchy(df, verbose=False):
host_DF = df.copy()
HierarchicalHostedSNe = []
parents = []

transientCoord = SkyCoord(host_DF['TransientRA'].values[0], host_DF['TransientDEC'].values[0], unit=u.deg)
for idx, row in host_DF.iterrows():
# cone search of the best-fit host in SIMBAD - if it gets it right,
#replace the info with the parent information!
Expand All @@ -173,9 +175,6 @@ def checkSimbadHierarchy(df, verbose=False):
if (tap_pandas['main_id'].values[0].startswith("VIRTUAL PARENT")) or (tap_pandas['otype'].values[0] == 'GrG'):
continue

if verbose:
print("Warning! Host of %s is the hierarchical child of another object in Simbad, choosing parent as host instead..." % row.TransientName)

# query PS1 for correct host
a = ps1cone(tap_pandas.loc[0, 'ra'], tap_pandas.loc[0, 'dec'], 10./3600)
if a:
Expand All @@ -187,8 +186,21 @@ def checkSimbadHierarchy(df, verbose=False):
parent['TransientRA'] = row.TransientRA
parent['TransientDEC'] = row.TransientDEC
parent = getNEDInfo(parent)
HierarchicalHostedSNe.append(row.TransientName)
parents.append(parent)

#only choose the new source if the DLR < 5.
r_a = np.nanmax([parent['%sKronRad'%i] for i in 'grizy'])
bands = 'grizy'
best_band = bands[np.argmax([row['%sKronRad'%i] for i in bands])]
dist, DLR = calc_DLR(Angle(row['TransientRA'], unit=u.deg), Angle(row['TransientDEC'], unit=u.deg), parent['raMean'].values[0], parent['decMean'].values[0], r_a, r_a, parent, best_band)
if DLR < 5:
if verbose:
print("Warning! Host of %s is the hierarchical child of another object in Simbad with dist/DLR<5, choosing parent as host instead..." % row.TransientName)
HierarchicalHostedSNe.append(row.TransientName)
parents.append(parent)
else:
if verbose:
print("Host of %s is the hierarchical child of another object in Simbad but dist/DLR>5, keeping existing host." % row.TransientName)

if len(parents)>0:
parentDF = pd.concat(parents)
else:
Expand Down Expand Up @@ -693,7 +705,7 @@ def findNewHosts(transientName, snCoord, snClass, verbose=False, starcut='gentle

#new low-z method (beta) - before we do anything else, find and associate with GLADE
if GLADE:
# snag the GWGC, which we'll use to get radii
# snag the GWGC, which we'll use to get radii
stream = importlib_resources.files(__name__).joinpath('gwgc_good.csv')
GWGC = pd.read_csv(stream)
if verbose:
Expand Down

0 comments on commit 1bf4814

Please sign in to comment.