mirror of
https://bitbucket.org/cosmicvoids/vide_public.git
synced 2025-07-04 15:21:11 +00:00
Added double-checking of flagged galaxies to ensure survey boundaries are respected. Python components of launchers now add their status to log files. Rearranged some warnings for broken tessellations due to guard point encounters, since it's not always a bad thing.
This commit is contained in:
parent
5d93a8a737
commit
ad2d3722cc
5 changed files with 122 additions and 36 deletions
|
@ -131,7 +131,7 @@ def figureOutMask(galFile, nside, outMaskFile):
|
|||
# figures out which galaxies live on a mask or redshift edge
|
||||
def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
||||
zmin, zmax, omegaM, useComoving, boundaryWidth,
|
||||
meanPartSep, outputDir):
|
||||
meanPartSep, outputDir, log):
|
||||
|
||||
if useComoving:
|
||||
zmin = comovingDistance(zmin, Om=omegaM)*LIGHT_SPEED
|
||||
|
@ -140,12 +140,13 @@ def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
|||
zmin *= LIGHT_SPEED
|
||||
zmax *= LIGHT_SPEED
|
||||
|
||||
log.write(" Reading contour map...\n")
|
||||
contourMap = healpy.read_map(contourFile)
|
||||
nside = healpy.get_nside(contourMap)
|
||||
npix = healpy.nside2npix(nside)
|
||||
|
||||
# load in galaxies
|
||||
#galPos = np.genfromtxt(galFile)
|
||||
log.write(" Loading galaxies...\n")
|
||||
# TODO - WHY IS THIS FASTER THAN np.column_stack???
|
||||
galPos = np.genfromtxt(outputDir+"/galaxies.txt")
|
||||
with open(galFile, 'rb') as File:
|
||||
|
@ -181,6 +182,7 @@ def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
|||
#galPos = np.column_stack((x,y,z))
|
||||
#print(galPos.shape)
|
||||
|
||||
log.write(" Building k-d tree...\n")
|
||||
flagList = np.zeros(len(galPos[:,0]))
|
||||
galTree = scipy.spatial.cKDTree(galPos)
|
||||
|
||||
|
@ -188,10 +190,14 @@ def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
|||
# using the "ray marching" algorithm: follow rays along lines of sight
|
||||
# of all mask edges, flagging nearest neighbor galaxies as we go
|
||||
|
||||
raySteps = np.arange(zmin, zmax, meanPartSep)
|
||||
log.write(" Begin ray-marching...\n")
|
||||
rayDist = meanPartSep
|
||||
raySteps = np.arange(zmin, zmax, rayDist)
|
||||
log.write(" Will take %d steps per ray with %.2e distance between steps\n" % (len(raySteps), rayDist))
|
||||
#print(meanPartSep, len(raySteps))
|
||||
|
||||
contourPixels = np.nonzero(contourMap)[0]
|
||||
log.write(" We have %d rays to work with\n" % (len(contourPixels)))
|
||||
#print(contourPixels)
|
||||
for pixel in contourPixels:
|
||||
#print("Working with pixel %d" % pixel)
|
||||
|
@ -207,6 +213,7 @@ def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
|||
# flag galaxies near redsfhit boundaries
|
||||
# TODO - save time by only covering portion of sphere that has data
|
||||
|
||||
log.write(" Flagging galaxies near redshift bounds...\n")
|
||||
sphereIndices = np.arange(len(contourMap))
|
||||
vec = healpy.pix2vec(nside, sphereIndices)
|
||||
vec = np.asarray(vec)
|
||||
|
@ -221,10 +228,12 @@ def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
|||
flagList[nearest] = 3
|
||||
|
||||
# output flag information
|
||||
log.write(" Saving galaxy flags to file...\n")
|
||||
np.savetxt(edgeGalFile, flagList, fmt="%d")
|
||||
|
||||
# paint galaxy flags onto healpix map for diagnostics
|
||||
# TODO - drop this when done testing
|
||||
log.write(" Saving diagnostic map to file...\n")
|
||||
flagMap = np.zeros(len(contourMap))
|
||||
justEdgeRA = RA[flagList == 1]
|
||||
justEdgeDec = Dec[flagList == 1]
|
||||
|
@ -234,7 +243,7 @@ def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
|||
ipix = healpy.ang2pix(nside, theta, phi)
|
||||
np.put(flagMap, ipix, 1)
|
||||
|
||||
healpy.write_map("./flagged_galaxies.fits", flagMap, overwrite=True,
|
||||
healpy.write_map(outputDir+"/flagged_galaxies.fits", flagMap, overwrite=True,
|
||||
dtype=np.dtype('float64'))
|
||||
|
||||
|
||||
|
|
Loading…
Add table
Add a link
Reference in a new issue