mirror of
https://bitbucket.org/cosmicvoids/vide_public.git
synced 2025-07-04 15:21:11 +00:00
Implemented (yet another) new boundary handling scheme, whereby we scan radially along survey edge while flagging nearest galaxies. The prepObservation routine was significantly cleaned up to accommodate this, but it was ultimately implemented in python (surveyTools.py) for ease of prototyping, with the intent to move it back into C later.
Some general housekeeping, making sure some new parameters are passed around correctly, and removing the storage of some unused files. This update is considered HIGHLY UNSTABLE. It will almost certainly break somewhere for simulations. Still under active development.
This commit is contained in:
parent
62dd66be79
commit
3dce2593d9
9 changed files with 348 additions and 454 deletions
|
@ -21,6 +21,7 @@
|
|||
# distances, and expected void stretching
|
||||
|
||||
import numpy as np
|
||||
import scipy
|
||||
import healpy as healpy
|
||||
import os
|
||||
from backend import *
|
||||
|
@ -127,63 +128,111 @@ def figureOutMask(galFile, nside, outMaskFile):
|
|||
return mask
|
||||
|
||||
# -----------------------------------------------------------------------------
|
||||
# figures out which galaxies live on a mask edge, and also writes the edge
|
||||
# map to an auxillary file
|
||||
def findEdgeGalaxies(galFile, maskFile, edgeGalFile, edgeMaskFile,
|
||||
zmin, zmax, omegaM, useComoving, boundaryWidth):
|
||||
# figures out which galaxies live on a mask or redshift edge
|
||||
def findEdgeGalaxies(galFile, maskFile, edgeGalFile, contourFile,
|
||||
zmin, zmax, omegaM, useComoving, boundaryWidth,
|
||||
meanPartSep):
|
||||
|
||||
if useComoving:
|
||||
zmin = comovingDistance(zmin, Om=omegaM)
|
||||
zmax = comovingDistance(zmax, Om=omegaM)
|
||||
#zmin = LIGHT_SPEED/100.*comovingDistance(zmin, Om=omegaM)
|
||||
#zmax = LIGHT_SPEED/100.*comovingDistance(zmax, Om=omegaM)
|
||||
#else:
|
||||
# zmin *= LIGHT_SPEED/100.
|
||||
# zmax *= LIGHT_SPEED/100.
|
||||
zmin = comovingDistance(zmin, Om=omegaM)*LIGHT_SPEED
|
||||
zmax = comovingDistance(zmax, Om=omegaM)*LIGHT_SPEED
|
||||
else:
|
||||
zmin *= LIGHT_SPEED
|
||||
zmax *= LIGHT_SPEED
|
||||
|
||||
|
||||
mask = healpy.read_map(maskFile)
|
||||
nside = healpy.get_nside(mask)
|
||||
contourMap = healpy.read_map(contourFile)
|
||||
nside = healpy.get_nside(contourMap)
|
||||
npix = healpy.nside2npix(nside)
|
||||
edgeMask = np.zeros((npix))
|
||||
|
||||
edgeFile = open(edgeGalFile, "w")
|
||||
# load in galaxies
|
||||
galPos = np.genfromtxt(galFile)
|
||||
flagList = np.zeros(len(galPos[:,0]))
|
||||
galTree = scipy.spatial.cKDTree(galPos)
|
||||
|
||||
for line in open(galFile):
|
||||
line = line.split()
|
||||
RA = float(line[3])
|
||||
Dec = float(line[4])
|
||||
z = float(line[5])
|
||||
|
||||
if useComoving:
|
||||
z = comovingDistance(z/LIGHT_SPEED, Om=omegaM)
|
||||
else:
|
||||
z *= LIGHT_SPEED/100.
|
||||
|
||||
phi, theta = convertAngle(RA, Dec)
|
||||
|
||||
# check the mask edges
|
||||
ipix = healpy.ang2pix(nside, theta, phi)
|
||||
neighbors = healpy.get_all_neighbours(nside, ipix)
|
||||
isOnMaskEdge = any(mask[p] == 0 for p in neighbors)
|
||||
# flag galaxies near mask edges
|
||||
# using the "ray marching" algorithm: follow rays along lines of sight
|
||||
# of all mask edges, flagging nearest neighbor galaxies as we go
|
||||
|
||||
# check the redshift boundaries
|
||||
zbuffer = (zmax-zmin)*boundaryWidth
|
||||
isOnHighZEdge = (z >= zmax-zbuffer)
|
||||
isOnLowZEdge = (z <= zmin+zbuffer)
|
||||
raySteps = np.arange(zmin, zmax, meanPartSep)
|
||||
|
||||
contourPixels = np.nonzero(contourMap)[0]
|
||||
#print(contourPixels)
|
||||
for pixel in contourPixels:
|
||||
#print("Working with pixel %d" % pixel)
|
||||
vec = healpy.pix2vec(nside,pixel)
|
||||
x = raySteps * vec[0]
|
||||
y = raySteps * vec[1]
|
||||
z = raySteps * vec[2]
|
||||
ray = np.array((x,y,z)).T
|
||||
#print(ray)
|
||||
|
||||
dist, nearest = galTree.query(ray)
|
||||
flagList[nearest] = 1
|
||||
#print(nearest)
|
||||
|
||||
if isOnMaskEdge:
|
||||
edgeFile.write("1\n")
|
||||
edgeMask[ipix] = 1
|
||||
elif isOnHighZEdge:
|
||||
edgeFile.write("2\n")
|
||||
elif isOnLowZEdge:
|
||||
edgeFile.write("3\n")
|
||||
else:
|
||||
edgeFile.write("0\n")
|
||||
# flag galaxies near redsfhit boundaries
|
||||
# TODO - save time by only covering portion of sphere with data
|
||||
ds = np.sqrt(healpy.nside2pixarea(nside)) / 1000.
|
||||
phi = np.arange(0, 2*np.pi, ds*2)
|
||||
theta = np.arange(0, np.pi, ds)
|
||||
vec = healpy.ang2vec(theta, phi)
|
||||
|
||||
edgeFile.close()
|
||||
healpy.write_map(edgeMaskFile, edgeMask, overwrite=True,
|
||||
dtype=np.dtype('float64'))
|
||||
maxEdge = zmax * vec
|
||||
dist, nearest = galTree.query(maxEdge)
|
||||
#print(nearest)
|
||||
#print(galPos[nearest])
|
||||
flagList[nearest] = 2
|
||||
|
||||
minEdge = zmin * vec
|
||||
dist, nearest = galTree.query(minEdge)
|
||||
#print(nearest)
|
||||
#print(galPos[nearest])
|
||||
flagList[nearest] = 3
|
||||
|
||||
# output flag information
|
||||
np.savetxt(edgeGalFile, flagList, fmt="%d")
|
||||
|
||||
|
||||
|
||||
# # output galaxy edge flags
|
||||
# edgeFile = open(edgeGalFile, "w")
|
||||
#
|
||||
# for line in open(galFile):
|
||||
# line = line.split()
|
||||
# RA = float(line[3])
|
||||
# Dec = float(line[4])
|
||||
# z = float(line[5])
|
||||
#
|
||||
# if useComoving:
|
||||
# z = comovingDistance(z/LIGHT_SPEED, Om=omegaM)
|
||||
# else:
|
||||
# z *= LIGHT_SPEED/100.
|
||||
#
|
||||
# phi, theta = convertAngle(RA, Dec)
|
||||
#
|
||||
# # check the mask edges
|
||||
# ipix = healpy.ang2pix(nside, theta, phi)
|
||||
# neighbors = healpy.get_all_neighbours(nside, ipix)
|
||||
# isOnMaskEdge = any(mask[p] == 0 for p in neighbors)
|
||||
#
|
||||
# # check the redshift boundaries
|
||||
# zbuffer = (zmax-zmin)*boundaryWidth
|
||||
# isOnHighZEdge = (z >= zmax-zbuffer)
|
||||
# isOnLowZEdge = (z <= zmin+zbuffer)
|
||||
#
|
||||
# if isOnMaskEdge:
|
||||
# edgeFile.write("1\n")
|
||||
# edgeMask[ipix] = 1
|
||||
# elif isOnHighZEdge:
|
||||
# edgeFile.write("2\n")
|
||||
# elif isOnLowZEdge:
|
||||
#
|
||||
#edgeFile.write("3\n")
|
||||
# else:
|
||||
# edgeFile.write("0\n")
|
||||
#
|
||||
# edgeFile.close()
|
||||
# healpy.write_map(edgeMaskFile, edgeMask, overwrite=True,
|
||||
# dtype=np.dtype('float64'))
|
||||
|
||||
return
|
||||
|
|
Loading…
Add table
Add a link
Reference in a new issue