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:
Paul M. Sutter 2025-01-07 20:04:29 +08:00
parent 62dd66be79
commit 3dce2593d9
9 changed files with 348 additions and 454 deletions

View file

@ -69,7 +69,7 @@ class Sample:
nickName = "dim"
outputDir = ""
maskFile = "rast_window_512.fits"
nsideForMask = 128
nsideForContour = 128
selFunFile = "czselfunc.all.dr72dim.dat"
zBoundary = (0.0, 0.1)
zBoundaryMpc = (0., 300)
@ -78,7 +78,8 @@ class Sample:
zRange = (0.0, 0.1)
omegaM = 0.27
minVoidRadius = -1
fakeDensity = 0.01
meanPartSep = 1 # calculated mean particle separation
fakeDensity = 0.01 # TODO - remove
hasWeightedVolumes = False
profileBinSize = 2 # Mpc
autoNumInStack = -1 # set to >0 to automatically generate stacks of size N
@ -101,9 +102,9 @@ class Sample:
stacks = []
def __init__(self, dataFile="", fullName="", dataUnit=1,
nickName="", maskFile="", nsideForMask=128, selFunFile="",
nickName="", maskFile="", nsideForContour=128, selFunFile="",
zBoundary=(), zRange=(), zBoundaryMpc=(), boundaryWidth=0.1,
shiftSimZ=False,
shiftSimZ=False, meanPartSep = 1,
minVoidRadius=-1, fakeDensity=0.01, volumeLimited=True,
numAPSlices=1, hasWeightedVolumes=False,
includeInHubble=True, partOfCombo=False, isCombo=False,
@ -118,7 +119,7 @@ class Sample:
self.fullName = fullName
self.nickName = nickName
self.maskFile = maskFile
self.nsideForMask = nsideForMask
self.nsideForContour = nsideForContour
self.selFunFile = selFunFile
self.zBoundary = zBoundary
self.zBoundaryMpc = zBoundaryMpc
@ -126,6 +127,7 @@ class Sample:
self.shiftSimZ = shiftSimZ
self.zRange = zRange
self.minVoidRadius = minVoidRadius
self.meanPartSep = meanPartSep
self.fakeDensity = fakeDensity
self.hasWeightedVolumes = hasWeightedVolumes
self.volumeLimited = volumeLimited

View file

@ -39,17 +39,20 @@ from backend.cosmologyTools import *
from backend.surveyTools import *
import pickle
import scipy.interpolate as interpolate
import time
NetCDFFile = Dataset
ncFloat = 'f8' # Double precision
LIGHT_SPEED = 299792.458
#LIGHT_SPEED = 299792.458
# -----------------------------------------------------------------------------
def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
outputDir=None, figDir=None, logFile=None, useComoving=False,
continueRun=None, regenerate=False):
startTime = time.time()
if sample.dataType == "observation":
sampleName = sample.fullName
@ -67,14 +70,29 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
datafile = inputDataDir+"/"+sample.dataFile
if sample.maskFile == "":
sample.maskFile = outputDir + "/constructed_mask.fits"
figureOutMask(datafile, sample.nsideForMask, sample.maskFile)
if sample.nsideForContour == -1:
sample.nsideForContour = 128
sample.maskFile = outputDir + "/constructed_mask.fits"
figureOutMask(datafile, sample.nsideForContour, sample.maskFile)
# compute mean particle separation
(boxVol, nbar) = getSurveyProps(sample.maskFile, sample.zRange[0],
sample.zRange[1], sample.zRange[0], sample.zRange[1], "all",
sample.omegaM, useComoving=useComoving)
numTracers = int(open(outputDir+"/mask_index.txt", "r").read())
sample.meanPartSep = (1.*numTracers/boxVol/nbar)**(-1/3.)
# flag edge galaxies
galFile = outputDir + "galaxies.txt"
edgeGalFile = outputDir + "/galaxy_edge_flags.txt"
edgeMaskFile = outputDir + "/mask_edge_map.fits"
findEdgeGalaxies(datafile, sample.maskFile, edgeGalFile, edgeMaskFile,
#edgeMaskFile = outputDir + "/mask_edge_map.fits"
contourFile = outputDir + "/contour_map.fits"
findEdgeGalaxies(galFile, sample.maskFile, edgeGalFile, contourFile,
sample.zBoundary[0], sample.zBoundary[1], sample.omegaM,
useComoving, sample.boundaryWidth)
useComoving, sample.boundaryWidth, sample.meanPartSep)
if useComoving:
useComovingFlag = "useComoving"
@ -92,12 +110,15 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
%s
%s
omegaM %g
nsideForContour %g
meanPartSep %g
""" % (datafile, sample.maskFile, outputFile,
outputDir+"/zobov_slice_"+sampleName+".par",
sample.zBoundary[0], sample.zBoundary[1], sample.fakeDensity,
useComovingFlag, inputParameterFlag, sample.omegaM)
useComovingFlag, inputParameterFlag, sample.omegaM,
sample.nsideForContour, sample.meanPartSep)
parmFile = os.getcwd()+"/generate_"+sample.fullName+".par"
parmFile = os.getcwd()+"/prep_"+sample.fullName+".par"
if regenerate or not (continueRun and jobSuccessful(logFile, "Done!\n")):
with open(parmFile, mode="wt") as f:
@ -106,9 +127,11 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
with open(logFile, 'wt') as log:
subprocess.call([binPath, arg1], stdout=log, stderr=log)
if jobSuccessful(logFile, "Done!\n"):
print("done")
endTime = time.time()
walltime = endTime - startTime
print("done (%.2fs elapsed)" % walltime)
else:
print("FAILED!")
print("FAILED! See log file for details.")
exit(-1)
else:
@ -118,7 +141,6 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
if os.access("contour_map.fits", os.F_OK):
os.system("mv %s %s" % ("contour_map.fits", outputDir))
os.system("mv %s %s" % ("mask_map.fits", outputDir))
if os.access("comoving_distance.txt", os.F_OK):
os.system("mv %s %s" % ("comoving_distance.txt", outputDir))
@ -129,15 +151,26 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
if os.access("galaxies.txt", os.F_OK):
os.system("mv %s %s" % ("galaxies.txt", outputDir))
os.system("mv %s %s" % ("mock_galaxies.txt", outputDir))
os.system("mv %s %s" % ("mock_boundary.txt", outputDir))
os.system("mv %s %s" % ("mock_sphere.txt", outputDir))
#os.system("mv %s %s" % ("galaxy_edge_flags.txt", outputDir))
else: # simulation
sampleName = sample.fullName
datafile = inputDataDir+"/"+sample.dataFile
# compute mean particle separation
iX = float(sample.mySubvolume[0])
iY = float(sample.mySubvolume[1])
xMin = iX/sample.numSubvolumes * sample.boxLen
yMin = iY/sample.numSubvolumes * sample.boxLen
xMax = (iX+1)/sample.numSubvolumes * sample.boxLen
yMax = (iY+1)/sample.numSubvolumes * sample.boxLen
zMin = sample.zBoundaryMpc[0]
zMax = sample.zBoundaryMpc[1]
boxVol = (xMax-xMin)*(yMax-yMin)*(zMax-zMin)
sample.meanPartSep = (1.*numTracers/boxVol)**(-1/3.)
# check if the final subsampling is done
lastSample = sample.subsample.split(', ')[-1]
doneLine = "Done! %5.2e\n" % float(lastSample)
@ -245,7 +278,7 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
cmd = "%s --configFile=%s" % (binPath,parmFile)
log = open(logFile, 'a')
arg1 = "--configFile=%s" % parmFile
subprocess.call(cmd, stdout=log, stderr=log, shell=True)
log.close()
@ -257,7 +290,7 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
doneLine = "Done! %5.2e\n" % keepFraction
if not jobSuccessful(logFile, doneLine):
print("FAILED!") ### dies here for now
print("FAILED! See log file for details.") ### dies here for now
exit(-1)
prevSubSample = thisSubSample
@ -281,29 +314,6 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
os.system("mv %s %s" % ("total_particles.txt", outputDir))
#os.system("mv %s %s" % ("sample_info.txt", outputDir))
# add to sample info file
if sample.dataType == "observation":
(boxVol, nbar) = getSurveyProps(sample.maskFile, sample.zRange[0],
sample.zRange[1], sample.zRange[0], sample.zRange[1], "all",
sample.omegaM, useComoving=useComoving)
else:
iX = float(sample.mySubvolume[0])
iY = float(sample.mySubvolume[1])
xMin = iX/sample.numSubvolumes * sample.boxLen
yMin = iY/sample.numSubvolumes * sample.boxLen
xMax = (iX+1)/sample.numSubvolumes * sample.boxLen
yMax = (iY+1)/sample.numSubvolumes * sample.boxLen
zMin = sample.zBoundaryMpc[0]
zMax = sample.zBoundaryMpc[1]
boxVol = (xMax-xMin)*(yMax-yMin)*(zMax-zMin)
nbar = 1.0
numTracers = int(open(outputDir+"/mask_index.txt", "r").read())
numTotal = int(open(outputDir+"/total_particles.txt", "r").read())
meanSep = (1.*numTracers/boxVol/nbar)**(-1/3.)
# save this sample's information
with open(outputDir+"/sample_info.dat", mode='wb') as output:
pickle.dump(sample, output, pickle.HIGHEST_PROTOCOL)
@ -325,9 +335,8 @@ def launchPrep(sample, binPath, workDir=None, inputDataDir=None,
fp.write("Number of simulation subvolumes: %s\n" % sample.numSubvolumes)
fp.write("My subvolume index: %s\n" % sample.mySubvolume)
fp.write("Estimated volume (cubic Mpc/h): %g\n" % boxVol)
fp.write("Number of real (non-boundary) tracers: %d\n" % numTracers)
fp.write("Total number of tracers: %d\n" % numTotal)
fp.write("Estimated mean tracer separation (Mpc/h): %g\n" % meanSep)
fp.write("Total number of tracers: %d\n" % numTracers)
fp.write("Estimated mean tracer separation (Mpc/h): %g\n" % sample.meanPartSep)
fp.write("Minimum void size actually used (Mpc/h): %g\n" % sample.minVoidRadius)
fp.close()
@ -336,6 +345,8 @@ def launchZobov(sample, binPath, outputDir=None, logDir=None, continueRun=None,
numZobovDivisions=None, numZobovThreads=None,
mergingThreshold=0.2):
startTime = time.time()
sampleName = sample.fullName
datafile = outputDir+"zobov_slice_"+sampleName
@ -490,9 +501,11 @@ def launchZobov(sample, binPath, outputDir=None, logDir=None, continueRun=None,
os.unlink(fileName)
if jobSuccessful(logFile, "Done!\n"):
print("done")
endTime = time.time()
walltime = endTime - startTime
print("done (%.2fs elapsed)" % walltime)
else:
print("FAILED!")
print("FAILED! See log file for details.")
exit(-1)
else:
@ -507,6 +520,8 @@ def launchPrune(sample, binPath,
continueRun=None, useComoving=False, mergingThreshold=0.2,
boundaryTolerance=1.0):
startTime = time.time()
sampleName = sample.fullName
numVoids = sum(1 for line in \
@ -580,9 +595,11 @@ def launchPrune(sample, binPath,
if jobSuccessful(logFile, "NetCDF: Not a valid ID\n") or \
jobSuccessful(logFile, "Done!\n"):
print("done")
endTime = time.time()
walltime = endTime - startTime
print("done (%.2fs elapsed)" % walltime)
else:
print("FAILED!")
print("FAILED! See log file for details.")
#exit(-1)
else:
@ -597,6 +614,8 @@ def launchVoidOverlap(sample1, sample2, sample1Dir, sample2Dir,
overlapFrac=0.25,
matchMethod=None, strictMatch=False):
startTime = time.time()
sampleName1 = sample1.fullName
sampleName2 = sample2.fullName
@ -663,7 +682,9 @@ def launchVoidOverlap(sample1, sample2, sample1Dir, sample2Dir,
log.close()
#if jobSuccessful(logFile, "Done!\n"):
print("done")
endTime = time.time()
walltime = endTime - startTime
print("done (%.2fs elapsed)" % walltime)
#else:
# print "FAILED!"
# exit(-1)
@ -707,7 +728,7 @@ def launchVelocityStack(sample, stack, binPath,
if jobSuccessful(logFile, "Done!\n"):
print("done")
else:
print("FAILED!")
print("FAILED! See log file for details.")
exit(-1)
else:

View file

@ -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