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

@ -18,9 +18,9 @@
+*/
#include <healpix_map.h>
#include <healpix_map_fitsio.h>
#include <pointing.h>
#include <boost/format.hpp>
#include <iostream>
#include <fstream>
@ -30,7 +30,6 @@
#include "contour_pixels.hpp"
#include <netcdf>
#include <CosmoTool/fortran.hpp>
#include <gsl/gsl_interp.h>
#include <gsl/gsl_integration.h>
#define LIGHT_SPEED 299792.458
@ -65,10 +64,10 @@ struct ParticleData
vector<double> redshift;
vector<double> catalogID;
vector<long> uniqueID;
int id_mask;
// PMS
int mask_index;
// END PMS
int edgeFlag = 0;
vector<Position> pos;
double box[3][2];
double Lmax;
@ -122,15 +121,13 @@ void loadData(const string& fname, NYU_VData & data)
}
void placeGalaxiesInCube(NYU_VData& data, ParticleData& output_data,
bool useComoving, double omegaM)
{
bool useComoving, double omegaM) {
double d2r = M_PI/180;
gsl_function expanF;
expanF.function = &expanFun;
struct my_expan_params expanParams;
double maxZ = 2.0, z, result, error, *dL, *redshifts;
int numZ = 1000, iZ;
double result, error;
size_t nEval;
expanParams.Om = omegaM;
@ -138,20 +135,6 @@ void placeGalaxiesInCube(NYU_VData& data, ParticleData& output_data,
expanParams.wa = 0.0;
expanF.params = &expanParams;
dL = (double *) malloc(numZ * sizeof(double));
redshifts = (double *) malloc(numZ * sizeof(double));
for (iZ = 0; iZ < numZ; iZ++) {
z = iZ * maxZ/numZ;
//gsl_integration_qng(&expanF, 0.0, z, 1.e-6, 1.e-6, &result, &error, &nEval);
dL[iZ] = result;
redshifts[iZ] = z;
}
gsl_interp *interp = gsl_interp_alloc(gsl_interp_linear, numZ);
gsl_interp_init(interp, redshifts, dL, numZ);
gsl_interp_accel *acc = gsl_interp_accel_alloc();
output_data.pos.resize(data.size());
output_data.id_gal.resize(data.size());
output_data.ra.resize(data.size());
@ -168,43 +151,36 @@ void placeGalaxiesInCube(NYU_VData& data, ParticleData& output_data,
for (int i = 0; i < data.size(); i++)
{
double ra = data[i].ra*d2r, dec = data[i].dec*d2r;
double Dc = data[i].cz;
Position& p = output_data.pos[i];
if (useComoving) {
//double pos = gsl_interp_eval(interp, redshifts, dL, data[i].cz, acc);
// Maubert - Lower boundary in redshift set to 0 to be consistent with pruneVoids (was 1.e-6 before).
gsl_integration_qng(&expanF, 0.0, data[i].cz/LIGHT_SPEED,
1.e-6,
1.e-6, &result, &error, &nEval);
double Dc = result*LIGHT_SPEED;
p.xyz[0] = Dc*cos(ra)*cos(dec);
p.xyz[1] = Dc*sin(ra)*cos(dec);
p.xyz[2] = Dc*sin(dec);
} else {
p.xyz[0] = data[i].cz*cos(ra)*cos(dec);
p.xyz[1] = data[i].cz*sin(ra)*cos(dec);
p.xyz[2] = data[i].cz*sin(dec);
1.e-6, 1.e-6, &result, &error, &nEval);
Dc = result*LIGHT_SPEED;
}
//printf("CREATE %e %e\n", data[i].cz, sqrt(p.xyz[0]*p.xyz[0] + p.xyz[1]*p.xyz[1] + p.xyz[2]*p.xyz[2]));
p.xyz[0] = Dc*cos(ra)*cos(dec);
p.xyz[1] = Dc*sin(ra)*cos(dec);
p.xyz[2] = Dc*sin(dec);
output_data.id_gal[i] = data[i].index;
output_data.ra[i] = ra;
output_data.dec[i] = dec;
output_data.redshift[i] = data[i].cz;
output_data.uniqueID[i] = data[i].uniqueID;
for (int j = 0; j < 3; j++)
{
if (p.xyz[j] > output_data.box[j][0])
output_data.box[j][0] = p.xyz[j];
if (p.xyz[j] < output_data.box[j][1])
output_data.box[j][1] = p.xyz[j];
}
//printf("INSERT GAL %d %e %e %e\n", output_data.id_gal[i], p.xyz[0], p.xyz[1], p.xyz[2]);
for (int j = 0; j < 3; j++) {
if (p.xyz[j] > output_data.box[j][0])
output_data.box[j][0] = p.xyz[j];
if (p.xyz[j] < output_data.box[j][1])
output_data.box[j][1] = p.xyz[j];
}
}
// normalize box
float left = 1.e99;
float right = -1.e99;
// normalize the box volume
float left = INFINITY;
float right = -INFINITY;
for (int j = 0; j < 3; j++) {
if (output_data.box[j][1] < left) left = output_data.box[j][1];
if (output_data.box[j][0] > right) right = output_data.box[j][0];
@ -214,24 +190,27 @@ void placeGalaxiesInCube(NYU_VData& data, ParticleData& output_data,
output_data.box[j][0] = right;
}
double Rmax = -1;
for (int j = 0; j < 3; j++) {
Rmax = max(Rmax, max(output_data.box[j][0], -output_data.box[j][1]));
}
output_data.Lmax = Rmax;
cout << format("Galaxy position generated: %d galaxies") % output_data.pos.size() << endl;
cout << format("box is %g < x < %g; %g < y < %g; %g < z < %g")
% (1e-2*output_data.box[0][1]) % (1e-2*output_data.box[0][0])
% (1e-2*output_data.box[1][1]) % (1e-2*output_data.box[1][0])
% (1e-2*output_data.box[2][1]) % (1e-2*output_data.box[2][0]) << endl;
gsl_interp_free(interp);
}
void generateSurfaceMask(prepObservation_info& args ,
void flagEdgeGalaxies(prepObservation_info& args ,
Healpix_Map<float>& mask,
vector<int>& pixel_list,
vector<int>& full_mask_list,
vector<int>& contourPixels,
NYU_VData& data,
ParticleData& output_data,
bool useComoving,
double omegaM)
{
double omegaM) {
//Maubert - Needed for comobile distance in mock_sphere
gsl_function expanF;
@ -242,38 +221,14 @@ void generateSurfaceMask(prepObservation_info& args ,
expanParams.wa = 0.0;
expanF.params = &expanParams;
double result, error ;
double result, error;
size_t nEval;
//End Maubert - Needed for comobile distance in mock_sphere
// Find the first free index
int idx = -1;
int insertion = 0;
double volume = pixel_list.size()*1.0/mask.Npix()*4*M_PI;
int numToInsert;
for (int i = 0; i < output_data.id_gal.size(); i++)
{
if (idx < output_data.id_gal[i])
idx = output_data.id_gal[i]+1;
}
output_data.id_mask = idx;
// PMS
// TODO - REMOVE THIS
output_data.mask_index = output_data.id_gal.size();
// END PMS
cout << "Generate surface mask..." << endl;
double Rmax = -1;
for (int j = 0; j < 3; j++)
{
Rmax = max(Rmax, max(output_data.box[j][0], -output_data.box[j][1]));
}
output_data.Lmax = Rmax;
// PMS - write a small text file with galaxy position (for diagnostic purposes)
// write a small text file with galaxy position (for diagnostic purposes)
FILE *fp;
fp = fopen("galaxies.txt", "w");
for (int i = 0; i < data.size(); i++) {
@ -284,208 +239,82 @@ void generateSurfaceMask(prepObservation_info& args ,
(p.xyz[2]));
}
fclose(fp);
// END PMS
cout << format("Rmax is %g, surface volume is %g") % (Rmax/100) % (volume/(4*M_PI)) << endl;
volume *= Rmax*Rmax*Rmax/3/1e6;
numToInsert = (int)floor(volume*args.density_fake_arg);
// TEST NOT USING MOCK PARTICLES
numToInsert = 0;
// END TEST
cout << format("3d volume to fill: %g (Mpc/h)^3") % volume << endl;
/* NOTE: temporarily moved to python for quick debugging. Will move back to
here once it's all sorted
cout << format("Will insert %d particles") % numToInsert << endl;
fp = fopen("mock_galaxies.txt", "w");
double pct = 0;
for (int i = 0; i < numToInsert; i++) {
double new_pct = i*100./numToInsert;
if (new_pct-pct > 5.) {
pct = new_pct;
cout << format(" .. %3.0f %%") % pct << endl;
}
Position p;
bool stop_here;
do {
int p0 = (int)floor(drand48()*pixel_list.size());
vec3 v = mask.pix2vec(pixel_list[p0]);
double r = Rmax*pow(drand48(),1./3);
p.xyz[0] = v.x * r;
p.xyz[1] = v.y * r;
p.xyz[2] = v.z * r;
stop_here = true;
for (int j = 0; j < 3; j++) {
if (p.xyz[j] > output_data.box[j][0] ||
p.xyz[j] < output_data.box[j][1])
stop_here = false;
}
}
while (!stop_here);
// PMS : write mock galaxies to a small file for diagnostic purposes
fprintf(fp, "%e %e %e\n",
(p.xyz[0]),
(p.xyz[1]),
(p.xyz[2]));
// END PMS
output_data.pos.push_back(p);
output_data.id_gal.push_back(idx);
output_data.ra.push_back(-1);
output_data.dec.push_back(-1);
output_data.redshift.push_back(-1);
output_data.uniqueID.push_back(-1);
//printf("INSERT MOCK %d %e %e %e\n", idx, p.xyz[0], p.xyz[1], p.xyz[2]);
insertion++;
}
fclose(fp);
// PMS
// TEST - insert mock galaxies along box edge - this is for tesselation safety
fp = fopen("mock_boundary.txt", "w");
double dx[3];
dx[0] = output_data.box[0][1] - output_data.box[0][0];
dx[1] = output_data.box[1][1] - output_data.box[1][0];
dx[2] = output_data.box[2][1] - output_data.box[2][0];
int nPart = 100;
// TEST
for (int iDir = 0; iDir < 0; iDir++) {
for (int iFace = 0; iFace < 0; iFace++) {
//for (int iDir = 0; iDir < 3; iDir++) {
//for (int iFace = 0; iFace < 2; iFace++) {
int iy = (iDir + 1) % 3;
int iz = (iDir + 2) % 3;
for (int i = 0; i < nPart; i++) {
for (int j = 0; j < nPart; j++) {
Position p;
p.xyz[iDir] = output_data.box[iDir][iFace];
p.xyz[iy] = i * dx[iy]/nPart + output_data.box[iy][0];
p.xyz[iz] = j * dx[iz]/nPart + output_data.box[iz][0];
output_data.pos.push_back(p);
output_data.id_gal.push_back(idx);
output_data.ra.push_back(-1);
output_data.dec.push_back(-1);
output_data.redshift.push_back(-1);
output_data.uniqueID.push_back(-1);
insertion++;
fprintf(fp, "%e %e %e\n",
(p.xyz[0]),
(p.xyz[1]),
(p.xyz[2]));
}
}
}
}
fclose(fp);
// END PMS
// PMS
// TEST - insert mock galaxies along spheres of survey redshift boundaries
fp = fopen("mock_sphere.txt", "w");
//Maubert - insert mock galaxies according to useComoving specification
// Added & compute rmin & rmax out of loop
double rmin ;
double rmax ;
// convert redshift boundaries to covmoving if necessary
double rmin;
double rmax;
if (useComoving) {
gsl_integration_qng(&expanF, 0.0, args.zMin_arg,
1.e-6,
1.e-6, &result, &error, &nEval);
rmin = result* LIGHT_SPEED;
gsl_integration_qng(&expanF, 0.0, args.zMin_arg, 1.e-6, 1.e-6, &result,
&error, &nEval);
rmin = result*LIGHT_SPEED;
gsl_integration_qng(&expanF, 0.0, args.zMax_arg,
1.e-6,
1.e-6, &result, &error, &nEval);
rmax = result* LIGHT_SPEED;
gsl_integration_qng(&expanF, 0.0, args.zMax_arg, 1.e-6, 1.e-6, &result,
&error, &nEval);
rmax = result*LIGHT_SPEED;
} else {
rmin = args.zMin_arg * LIGHT_SPEED;
rmax = args.zMax_arg * LIGHT_SPEED;
rmin = args.zMin_arg * LIGHT_SPEED;
rmax = args.zMax_arg * LIGHT_SPEED;
}
// TEST NOT USING BOUNDARY PARTICLES
for (int q = 0; q < 0; q++) {
//for (int q = 0; q < full_mask_list.size(); q++) {
vec3 v = mask.pix2vec(full_mask_list[q]);
Position p;
if (rmin > 0.) {
p.xyz[0] = v.x * rmin;
p.xyz[1] = v.y * rmin;
p.xyz[2] = v.z * rmin;
output_data.pos.push_back(p);
output_data.id_gal.push_back(idx);
output_data.ra.push_back(-1);
output_data.dec.push_back(-1);
output_data.redshift.push_back(-1);
output_data.uniqueID.push_back(-1);
insertion++;
fprintf(fp, "%e %e %e\n",
(p.xyz[0]),
(p.xyz[1]),
(p.xyz[2]));
}
p.xyz[0] = v.x * rmax;
p.xyz[1] = v.y * rmax;
p.xyz[2] = v.z * rmax;
output_data.pos.push_back(p);
output_data.id_gal.push_back(idx);
output_data.ra.push_back(-1);
output_data.dec.push_back(-1);
output_data.redshift.push_back(-1);
output_data.uniqueID.push_back(-1);
insertion++;
fprintf(fp, "%e %e %e\n",
(p.xyz[0]),
(p.xyz[1]),
(p.xyz[2]));
}
fclose(fp);
// END PMS
cout << format("Done. Inserted %d particles.") % insertion << endl;
}
void saveData(ParticleData& pdata)
{
NcFile f("particles.nc", NcFile::replace);
NcDim d = f.addDim("space", 3);
NcDim p = f.addDim("Np", pdata.pos.size());
NcVar v = f.addVar("particles", ncDouble, {d, p});
double *x = new double[pdata.pos.size()];
for (int j = 0; j < 3; j++)
{
for (int i = 0; i < pdata.pos.size(); i++)
x[i] = pdata.pos[i].xyz[j];
v.putVar({size_t(j), 0}, {1, pdata.pos.size()}, x);
}
v = f.addVar("id_gal", ncInt, std::vector<NcDim>({p}));
v.putVar(&pdata.id_gal[0]);
delete[] x;
}
double dx = args.meanPartSep_arg;
int nSteps = floor( (rmax - rmin) / dx);
cout << "Assumed resolution element: " << dx << " " << nSteps << endl;
void saveForZobov(ParticleData& pdata, const string& fname, const string& paramname)
// 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
// TODO - replace this with faster kd-tree search
cout << "Flagging galaxies on edges of survey..." << endl;
// explore rays along mask contours
for (int pixel : contourPixels) {
cout << "Working with pixel " << pixel << endl;
vec3 v = mask.pix2vec(pixel);
//cout << v*rmin << " " << v*rmax << endl;
// march along single ray and find nearest neighbors
for (int n = 0; n <= nSteps; n++) {
double r = n*dx + rmin;
vec3 rayPos = v*r;
double x = rayPos.x;
double y = rayPos.y;
double z = rayPos.z;
//cout << "Step " << n << " " << rayPos << endl;
// scan all galaxies
double minDist = INFINITY;
double dist = 0;
int closestGal = -1;
for (int i = 0; i < data.size(); i++) {
Position& galPos = output_data.pos[i];
dist = pow(galPos.xyz[0] - x, 2) +
pow(galPos.xyz[1] - y, 2) +
pow(galPos.xyz[2] - z, 2);
if (dist < minDist) closestGal = i;
} // galaxy search
} // marching along one ray
} // all contours
// flag galaxies near redshift boundaries
cout << "Flagging galaxies at redshift boundaries..." << endl;
*/
} // end flagEdgeGalaxies
void saveForZobov(ParticleData& pdata, const string& fname,
const string& paramname)
{
UnformattedWrite f(fname);
static const char axis[] = { 'X', 'Y', 'Z' };
@ -562,21 +391,6 @@ void saveForZobov(ParticleData& pdata, const string& fname, const string& paramn
v.putVar({0}, {size_t(nOutputPart)}, &pdata.id_gal[0]);
//v2->put(expansion_fac, pdata.pos.size());
//delete[] expansion_fac;
/*
FILE *infoFile = fopen("sample_info.txt", "w");
fprintf(infoFile, "x_min = %f\n", -Lmax/100.);
fprintf(infoFile, "x_max = %f\n", Lmax/100.);
fprintf(infoFile, "y_min = %f\n", -Lmax/100.);
fprintf(infoFile, "y_max = %f\n", Lmax/100.);
fprintf(infoFile, "z_min = %f\n", -Lmax/100.);
fprintf(infoFile, "z_max = %f\n", Lmax/100.);
fprintf(infoFile, "mask_index = %d\n", pdata.mask_index);
fprintf(infoFile, "total_particles = %d\n", pdata.pos.size());
fclose(infoFile);
*/
}
int main(int argc, char **argv)
@ -608,40 +422,40 @@ int main(int argc, char **argv)
prepObservation_conf_print_version();
cout << "Loading data " << args_info.catalog_arg << "..." << endl;
cout << "Loading galaxy data " << args_info.catalog_arg << "..." << endl;
vector<NYU_Data> data;
Healpix_Map<float> o_mask;
vector<int> pixel_list;
vector<int> full_mask_list;
vector<int> contourPixels;
ParticleData output_data;
loadData(args_info.catalog_arg, data);
cout << "Loading mask..." << endl;
Healpix_Map<float> mask;
Healpix_Map<float> o_mask;
int newNside = args_info.nsideForContour_arg;
read_Healpix_map_from_fits(args_info.mask_arg, o_mask);
Healpix_Map<float> mask;
if (newNside == -1) newNside = o_mask.Nside();
mask.SetNside(128, RING);
mask.SetNside(newNside, RING);
mask.Import(o_mask);
computeContourPixels(mask,pixel_list);
computeMaskPixels(mask,full_mask_list);
// We compute a cube holding all the galaxies + the survey surface mask
computeContourPixels(mask, contourPixels);
cout << "Placing galaxies..." << endl;
placeGalaxiesInCube(data, output_data, args_info.useComoving_flag,
args_info.omegaM_arg);
generateSurfaceMask(args_info, mask, pixel_list, full_mask_list,
data, output_data,args_info.useComoving_flag,
args_info.omegaM_arg);
args_info.omegaM_arg);
//cout << "Flagging edge galaxies..." << endl;
flagEdgeGalaxies(args_info, mask, contourPixels,
data, output_data,args_info.useComoving_flag,
args_info.omegaM_arg);
saveForZobov(output_data, args_info.output_arg, args_info.params_arg);
// saveData(output_data);
// PMS
// PMS - TODO REMOVE THIS
FILE *fp = fopen("mask_index.txt", "w");
fprintf(fp, "%d", output_data.mask_index);
fclose(fp);
@ -650,6 +464,7 @@ int main(int argc, char **argv)
fprintf(fp, "%d", output_data.pos.size());
fclose(fp);
printf("Done!\n");
// END PMS
// END PMS - TODO REMOVE THIS
return 0;
}

View file

@ -17,3 +17,7 @@ option "params" - "Output parameters of the datacube" string required
option "useComoving" - "Convert to real space using LCDM cosmology" flag off
option "omegaM" - "Omega Matter for fiducial cosmology" double optional default="0.27"
option "nsideForContour" - "HEALPix NSIDE resolution for figuring out mask contours" int optional default="-1"
option "meanPartSep" - "Estimated mean tracer seperation in h^3 / Mpc^3" double optional default="1"