def make_sky_box(self, padding=5.): maxx = max(np.max(self.ran.x), np.max(self.cat.x)) minx = min(np.min(self.ran.x), np.min(self.cat.x)) maxy = max(np.max(self.ran.y), np.max(self.cat.y)) miny = min(np.min(self.ran.y), np.min(self.cat.y)) maxz = max(np.max(self.ran.z), np.max(self.cat.z)) minz = min(np.min(self.ran.z), np.min(self.cat.z)) dx = maxx - minx dy = maxy - miny dz = maxz - minz x0 = 0.5 * (maxx + minx) y0 = 0.5 * (maxy + miny) z0 = 0.5 * (maxz + minz) box = max([dx, dy, dz ]) + 2 * padding # a bit bigger than strictly necessary xmin = x0 - box / 2 ymin = y0 - box / 2 zmin = z0 - box / 2 self.xmin = xmin self.ymin = ymin self.zmin = zmin self.box_length = box print('Box size [Mpc/h]: %0.3f' % self.box_length) # this is clearly a major underestimate mean_dens = np.sum(self.cat.weight) / box**3. # starting estimate for bin size self.nbins = int( np.floor(box / (0.5 * (4 * np.pi * mean_dens / 3.)**(-1. / 3)))) self.binsize = self.box_length / self.nbins if self.verbose: print('Initial bin size [Mpc/h]: %0.2f, nbins = %d' % (self.binsize, self.nbins)) # now approximately check true survey volume ran = self.ran rhor = np.zeros((self.nbins, self.nbins, self.nbins), dtype='float64') fastmodules.allocate_gal_cic(rhor, ran.x, ran.y, ran.z, ran.weight, ran.size, self.xmin, self.ymin, self.zmin, self.box_length, self.nbins, 1.) # re-estimate the mean density: this will still be a slight underestimate filled_cells = np.sum(rhor.flatten() > 0) mean_dens = np.sum(self.cat.weight) / (filled_cells * self.binsize**3.) # thus get better choice of bin size (this is sufficient for current purposes) self.nbins = int( np.floor(box / (0.5 * (4 * np.pi * mean_dens / 3.)**(-1. / 3)))) self.binsize = self.box_length / self.nbins print('Final bin size [Mpc/h]: %0.2f, nbins = %d' % (self.binsize, self.nbins)) # choose an appropriate smoothing scale self.smooth = mean_dens**(-1. / 3) print('Smoothing scale [Mpc/h]: %0.2f' % self.smooth) sys.stdout.flush() return mean_dens
def run_voidfinder(self): # create the folder in which to store various raw outputs raw_dir = self.output_folder + "rawVoxelInfo/" if not os.access(raw_dir, os.F_OK): os.makedirs(raw_dir) # get the path to where the C executables are stored binpath = os.path.dirname(__file__).replace('python_tools', 'bin/') if self.is_box: # measure the galaxy density field if self.verbose: print('Allocating galaxies in cells...') sys.stdout.flush() rhog = np.zeros((self.nbins, self.nbins, self.nbins), dtype='float64') fastmodules.allocate_gal_cic(rhog, self.cat.x, self.cat.y, self.cat.z, self.cat.weight, self.cat.size, self.xmin, self.ymin, self.zmin, self.box_length, self.nbins, 1) # smooth with pre-determined smoothing scale if self.verbose: print('Smoothing galaxy density field ...') sys.stdout.flush() rhog = gaussian_filter(rhog, self.smooth / self.binsize, mode='wrap') # then normalize number counts to get density in units of mean (i.e. 1 + delta) fastmodules.normalize_rho_box(rhog, self.cat.size) self.rhoflat = rhog.flatten() self.mask_cut = np.zeros( self.nbins**3, dtype='int') # we don't mask any voxels in a box else: # measure the galaxy density field if self.verbose: print('Allocating galaxies in cells...') sys.stdout.flush() rhog = np.zeros((self.nbins, self.nbins, self.nbins), dtype='float64') fastmodules.allocate_gal_cic(rhog, self.cat.x, self.cat.y, self.cat.z, self.cat.weight, self.cat.size, self.xmin, self.ymin, self.zmin, self.box_length, self.nbins, 1) if self.verbose: print('Allocating randoms in cells...') sys.stdout.flush() rhor = np.zeros((self.nbins, self.nbins, self.nbins), dtype='float64') fastmodules.allocate_gal_cic(rhor, self.ran.x, self.ran.y, self.ran.z, self.ran.weight, self.ran.size, self.xmin, self.ymin, self.zmin, self.box_length, self.nbins, 1) # identify "empty" cells for later cuts on void catalogue mask_cut = np.zeros(self.nbins**3, dtype='int') fastmodules.survey_mask(mask_cut, rhor, self.ran_min) self.mask_cut = mask_cut # smooth both galaxy and randoms with pre-determined smoothing scale if self.verbose: print('Smoothing density fields ...') sys.stdout.flush() rhog = gaussian_filter(rhog, self.smooth / self.binsize, mode='nearest') rhor = gaussian_filter(rhor, self.smooth / self.binsize, mode='nearest') rho = np.empty((self.nbins, self.nbins, self.nbins), dtype='float64') fastmodules.normalize_rho_survey(rho, rhog, rhor, self.alpha, self.ran_min) self.rhoflat = rho.flatten() # write this to file for jozov-grid to read rhogflat = np.array(self.rhoflat, dtype=np.float32) with open(raw_dir + 'density_n%d.dat' % self.nbins, 'w') as F: rhogflat.tofile(F, format='%f') # now call jozov-grid cmd = [ binpath + "jozov-grid", "v", raw_dir + "density_n%d.dat" % self.nbins, raw_dir + self.handle, str(self.nbins) ] subprocess.call(cmd) # postprocess void data self.postprocess_voids() # if reqd, find superclusters if self.find_clusters: print( "\n ==== bonus: overdensity-finding with voxel-based method ==== " ) sys.stdout.flush() cmd = [ binpath + "jozov-grid", "c", raw_dir + "density_n%d.dat" % self.nbins, raw_dir + self.handle, str(self.nbins) ] subprocess.call(cmd) self.postprocess_clusters() print(" ==== Finished with voxel-based method ==== ") sys.stdout.flush()
def iterate(self, iloop, save_wisdom=1, debug=False): cat = self.cat ran = self.ran binsize = self.binsize beta = self.beta bias = self.bias f = self.f nbins = self.nbins print("Loop %d" % iloop) # -- Creating arrays for FFTW if iloop == 0: # first iteration requires initialization delta = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') deltak = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') rho = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') rhok = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') psi_x = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') psi_y = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') psi_z = pyfftw.empty_aligned((nbins, nbins, nbins), dtype='complex128') # -- Initialize FFT objects and load wisdom if available wisdom_file = "wisdom." + str(nbins) + "." + str(self.nthreads) if os.path.isfile(wisdom_file): print('Reading wisdom from ', wisdom_file) g = open(wisdom_file, 'r') wisd = json.load(g) for i in range(len(wisd)): wisd[i] = wisd[i].encode('utf-8') wisd = tuple(wisd) pyfftw.import_wisdom(wisd) g.close() print('Creating FFTW objects...') sys.stdout.flush() fft_obj = pyfftw.FFTW(delta, deltak, axes=[0, 1, 2], threads=self.nthreads) fft_obj_inplace = pyfftw.FFTW(delta, delta, axes=[0, 1, 2], threads=self.nthreads) ifft_obj = pyfftw.FFTW(deltak, psi_x, axes=[0, 1, 2], threads=self.nthreads, direction='FFTW_BACKWARD') kr = fftfreq(nbins, d=binsize) * 2 * np.pi * self.smooth norm = np.exp(-0.5 * (kr[:, None, None]**2 + kr[None, :, None]**2 + kr[None, None, :]**2)) if self.is_box: deltar = 0 else: if self.verbose: print('Allocating randoms in cells...') sys.stdout.flush() deltar = np.zeros((nbins, nbins, nbins), dtype='float64') fastmodules.allocate_gal_cic(deltar, ran.x, ran.y, ran.z, ran.weight, ran.size, self.xmin, self.ymin, self.zmin, self.box, nbins, 1.) if self.verbose: print('Smoothing...') sys.stdout.flush() # NOTE - we do the smoothing via FFTs rather than scipy's gaussian_filter because if using several # threads for pyfftw it is much faster this way (if only using 1 thread gains are negligible) rho = deltar + 0.0j fft_obj(input_array=rho, output_array=rhok) fastmodules.mult_norm(rhok, rhok, norm) ifft_obj(input_array=rhok, output_array=rho) deltar = rho.real else: delta = self.delta deltak = self.deltak deltar = self.deltar rho = self.rho rhok = self.rhok psi_x = self.psi_x psi_y = self.psi_y psi_z = self.psi_z fft_obj_inplace = self.fft_obj_inplace fft_obj = self.fft_obj ifft_obj = self.ifft_obj norm = self.norm # -- Allocate galaxies and randoms to grid with CIC method # -- using new positions if self.verbose: print('Allocating galaxies in cells...') sys.stdout.flush() deltag = np.zeros((nbins, nbins, nbins), dtype='float64') fastmodules.allocate_gal_cic(deltag, cat.newx, cat.newy, cat.newz, cat.weight, cat.size, self.xmin, self.ymin, self.zmin, self.box, nbins, 1.) if self.verbose: print('Smoothing galaxy density field ...') sys.stdout.flush() # NOTE - smoothing via FFTs rho = deltag + 0.0j fft_obj(input_array=rho, output_array=rhok) fastmodules.mult_norm(rhok, rhok, norm) ifft_obj(input_array=rhok, output_array=rho) deltag = rho.real if self.verbose: print('Computing density fluctuations, delta...') sys.stdout.flush() if self.is_box: # simply normalize based on (constant) mean galaxy number density fastmodules.normalize_delta_box(delta, deltag, cat.size) else: # normalize using the randoms, avoiding possible divide-by-zero errors fastmodules.normalize_delta_survey(delta, deltag, deltar, self.alpha, self.ran_min) del deltag # deltag no longer required anywhere if self.verbose: print('Fourier transforming delta field...') sys.stdout.flush() fft_obj_inplace(input_array=delta, output_array=delta) # -- delta/k**2 k = fftfreq(self.nbins, d=binsize) * 2 * np.pi fastmodules.divide_k2(delta, delta, k) # now solve the basic building block: IFFT[-i k delta(k)/(b k^2)] if self.verbose: print('Inverse Fourier transforming to get psi...') sys.stdout.flush() fastmodules.mult_kx(deltak, delta, k, bias) ifft_obj(input_array=deltak, output_array=psi_x) fastmodules.mult_ky(deltak, delta, k, bias) ifft_obj(input_array=deltak, output_array=psi_y) fastmodules.mult_kz(deltak, delta, k, bias) ifft_obj(input_array=deltak, output_array=psi_z) # from grid values of Psi_est = IFFT[-i k delta(k)/(b k^2)], compute the values at the galaxy positions if self.verbose: print('Calculating shifts...') sys.stdout.flush() shift_x, shift_y, shift_z = self.get_shift(cat, psi_x.real, psi_y.real, psi_z.real, use_newpos=True) # now we update estimates of the Psi field in the following way: if iloop == 0: # starting estimate chosen according to Eq. 12 of Burden et al 2015, in order to improve convergence if self.is_box: # line-of-sight direction is along the z-axis (hard-coded) psi_dot_rhat = shift_z shift_z -= beta / (1 + beta) * psi_dot_rhat else: # line-of-sight direction determined by galaxy coordinates psi_dot_rhat = (shift_x * cat.x + shift_y * cat.y + shift_z * cat.z) / cat.dist shift_x -= beta / (1 + beta) * psi_dot_rhat * cat.x / cat.dist shift_y -= beta / (1 + beta) * psi_dot_rhat * cat.y / cat.dist shift_z -= beta / (1 + beta) * psi_dot_rhat * cat.z / cat.dist # given estimate of Psi, subtract approximate RSD to get estimate of real-space galaxy positions if self.is_box: # line-of-sight direction is along the z-axis (hard-coded) cat.newz = cat.z + f * shift_z # check PBC cat.newz[cat.newz >= cat.box_length] -= cat.box_length cat.newz[cat.newz < 0] += cat.box_length else: psi_dot_rhat = (shift_x * cat.x + shift_y * cat.y + shift_z * cat.z) / cat.dist cat.newx = cat.x + f * psi_dot_rhat * cat.x / cat.dist cat.newy = cat.y + f * psi_dot_rhat * cat.y / cat.dist cat.newz = cat.z + f * psi_dot_rhat * cat.z / cat.dist # for debugging: if self.verbose and debug: if self.is_box: print( 'Debug: first 10 x,y,z shifts and old and new z positions') for i in range(10): print('%0.3f %0.3f %0.3f %0.3f %0.3f' % (shift_x[i], shift_y[i], shift_z[i], cat.z[i], cat.newz[i])) else: print( 'Debug: first 10 x,y,z shifts and old and new observer distances' ) for i in range(10): oldr = np.sqrt(cat.x[i]**2 + cat.y[i]**2 + cat.z[i]**2) newr = np.sqrt(cat.newx[i]**2 + cat.newy[i]**2 + cat.newz[i]**2) print('%0.3f %0.3f %0.3f %0.3f %0.3f' % (shift_x[i], shift_y[i], shift_z[i], oldr, newr)) # in the next loop of iteration, these new positions are used to compute next approximation of # the (real-space) galaxy density, and then this is used to get new estimate of Psi, etc. # at the end of the iterations, newx, newy, newz should be the real-space galaxy positions (or best # estimate thereof) self.deltar = deltar self.delta = delta self.deltak = deltak self.rho = rho self.rhok = rhok self.psi_x = psi_x self.psi_y = psi_y self.psi_z = psi_z self.fft_obj_inplace = fft_obj_inplace self.fft_obj = fft_obj self.ifft_obj = ifft_obj self.norm = norm # -- save wisdom wisdom_file = "wisdom." + str(nbins) + "." + str(self.nthreads) if iloop == 0 and save_wisdom and not os.path.isfile(wisdom_file): wisd = pyfftw.export_wisdom() wisd = list(wisd) for i in range(len(wisd)): wisd[i] = wisd[i].decode('utf-8') f = open(wisdom_file, 'w') json.dump(wisd, f) f.close() print('Wisdom saved at', wisdom_file)