mirror of
https://bitbucket.org/cosmicvoids/vide_public.git
synced 2025-07-05 07:41:11 +00:00
renamed src to source
This commit is contained in:
parent
4dcaf3959b
commit
4d9c5ab2c1
71 changed files with 7 additions and 3 deletions
5
c_source/util/CMakeLists.txt
Normal file
5
c_source/util/CMakeLists.txt
Normal file
|
@ -0,0 +1,5 @@
|
|||
|
||||
add_library(zobovTool loadZobov.cpp particleInfo.cpp contour_pixels.cpp)
|
||||
|
||||
set_target_properties(zobovTool PROPERTIES COMPILE_FLAGS ${OpenMP_CXX_FLAGS} LINK_FLAGS ${OpenMP_CXX_FLAGS})
|
||||
target_link_libraries(zobovTool ${COSMOTOOL_LIBRARY} ${GSL_LIBRARIES} ${NETCDF_LIBRARIES} ${HDF5_CXX_LIBRARIES} ${HDF5_HL_LIBRARIES} ${HDF5_LIBRARIES} ${DL_LIBRARY})
|
107
c_source/util/contour_pixels.cpp
Normal file
107
c_source/util/contour_pixels.cpp
Normal file
|
@ -0,0 +1,107 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/contour_pixels.cpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#include <vector>
|
||||
#include <healpix_map.h>
|
||||
#include <healpix_map_fitsio.h>
|
||||
#include "contour_pixels.hpp"
|
||||
|
||||
using namespace std;
|
||||
|
||||
static const bool DEBUG = true;
|
||||
|
||||
void computeFilledPixels(Healpix_Map<float>& m, vector<int>& filled)
|
||||
{
|
||||
filled.clear();
|
||||
for (int p = 0; p < m.Npix(); p++)
|
||||
if (m[p] > 0)
|
||||
filled.push_back(p);
|
||||
}
|
||||
|
||||
void computeContourPixels(Healpix_Map<float>& m, vector<int>& contour)
|
||||
{
|
||||
contour.clear();
|
||||
for (int p = 0; p < m.Npix(); p++)
|
||||
{
|
||||
fix_arr<int, 8> result;
|
||||
|
||||
m.neighbors(p, result);
|
||||
for (int q = 0; q < 8; q++)
|
||||
{
|
||||
if (result[q] < 0)
|
||||
continue;
|
||||
|
||||
float delta = (m[p]-0.5)*(m[result[q]]-0.5);
|
||||
if (delta < 0)
|
||||
{
|
||||
contour.push_back(p);
|
||||
// This is boundary go to next pixel
|
||||
break;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
if (DEBUG)
|
||||
{
|
||||
Healpix_Map<int> contour_map;
|
||||
|
||||
contour_map.SetNside(m.Nside(), RING);
|
||||
contour_map.fill(0);
|
||||
for (int p = 0; p < contour.size(); p++)
|
||||
{
|
||||
contour_map[contour[p]]=1;
|
||||
}
|
||||
|
||||
fitshandle h;
|
||||
h.create("!contour_map.fits");
|
||||
write_Healpix_map_to_fits(h, contour_map, planckType<int>());
|
||||
}
|
||||
}
|
||||
|
||||
void computeMaskPixels(Healpix_Map<float>& m, vector<int>& contour)
|
||||
{
|
||||
for (int p = 0; p < m.Npix(); p++)
|
||||
{
|
||||
|
||||
if (m[p]>0)
|
||||
{
|
||||
contour.push_back(p);
|
||||
// This is boundary go to next pixel
|
||||
}
|
||||
}
|
||||
|
||||
if (DEBUG)
|
||||
{
|
||||
Healpix_Map<int> contour_map;
|
||||
|
||||
contour_map.SetNside(m.Nside(), RING);
|
||||
contour_map.fill(0);
|
||||
for (int p = 0; p < contour.size(); p++)
|
||||
{
|
||||
contour_map[contour[p]]=1;
|
||||
}
|
||||
|
||||
fitshandle h;
|
||||
h.create("!mask_map.fits");
|
||||
write_Healpix_map_to_fits(h, contour_map, planckType<int>());
|
||||
}
|
||||
}
|
||||
|
33
c_source/util/contour_pixels.hpp
Normal file
33
c_source/util/contour_pixels.hpp
Normal file
|
@ -0,0 +1,33 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/contour_pixels.hpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#ifndef __CONTOUR_PIXELS_HPP
|
||||
#define __CONTOUR_PIXELS_HPP
|
||||
|
||||
#include <vector>
|
||||
#include <healpix_map.h>
|
||||
|
||||
|
||||
void computeContourPixels(Healpix_Map<float>& map, std::vector<int>& contour);
|
||||
void computeFilledPixels(Healpix_Map<float>& map, std::vector<int>& contour);
|
||||
void computeMaskPixels(Healpix_Map<float>& map, std::vector<int>& contour);
|
||||
|
||||
#endif
|
54
c_source/util/gslIntegrate.hpp
Normal file
54
c_source/util/gslIntegrate.hpp
Normal file
|
@ -0,0 +1,54 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/gslIntegrate.hpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#ifndef __MYGSL_INTEGRATE_HPP
|
||||
#define __MYGSL_INTEGRATE_HPP
|
||||
|
||||
#include <gsl/gsl_integration.h>
|
||||
|
||||
template<typename FunT>
|
||||
double gslSpecialFunction(double x, void *param)
|
||||
{
|
||||
FunT *f = (FunT *)param;
|
||||
|
||||
return (*f)(x);
|
||||
}
|
||||
|
||||
template<typename FunT>
|
||||
double gslIntegrate(FunT& v, double a, double b, double prec, int NPTS = 1024)
|
||||
{
|
||||
gsl_integration_workspace *w = gsl_integration_workspace_alloc(NPTS);
|
||||
gsl_function f;
|
||||
double result;
|
||||
double abserr;
|
||||
|
||||
f.function = &gslSpecialFunction<FunT>;
|
||||
f.params = &v;
|
||||
|
||||
gsl_integration_qag(&f, a, b, prec, 0, NPTS, GSL_INTEG_GAUSS61,
|
||||
w, &result, &abserr);
|
||||
|
||||
gsl_integration_workspace_free(w);
|
||||
|
||||
return result;
|
||||
}
|
||||
|
||||
#endif
|
245
c_source/util/loadZobov.cpp
Normal file
245
c_source/util/loadZobov.cpp
Normal file
|
@ -0,0 +1,245 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/loadZobov.cpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#include <cassert>
|
||||
#include <string>
|
||||
#include <fstream>
|
||||
#include <iostream>
|
||||
#include <vector>
|
||||
#include <cstdlib>
|
||||
#include <sstream>
|
||||
#include <algorithm>
|
||||
#include "loadZobov.hpp"
|
||||
#include <CosmoTool/fortran.hpp>
|
||||
|
||||
using namespace std;
|
||||
using namespace CosmoTool;
|
||||
|
||||
bool loadZobov(const char *descName, const char *adjName, const char *voidsName,
|
||||
const char *volName, ZobovRep& z)
|
||||
{
|
||||
ifstream descFile(descName);
|
||||
ifstream adjFile(adjName);
|
||||
ifstream volFile(voidsName);
|
||||
int32_t numParticles, numZones, numPinZone;
|
||||
int32_t totalParticles;
|
||||
int32_t numVoids;
|
||||
int32_t minParticlesInZone, maxParticlesInZone;
|
||||
|
||||
adjFile.read((char *)&numParticles, sizeof(numParticles));
|
||||
adjFile.read((char *)&numZones, sizeof(numZones));
|
||||
if (!adjFile)
|
||||
return false;
|
||||
|
||||
cout << "Number of particles = " << numParticles << endl;
|
||||
cout << "Number of zones = " << numZones << endl;
|
||||
|
||||
totalParticles = 0;
|
||||
|
||||
minParticlesInZone = -1;
|
||||
maxParticlesInZone = -1;
|
||||
|
||||
z.allZones.resize(numZones);
|
||||
for (int zone = 0; zone < numZones; zone++)
|
||||
{
|
||||
adjFile.read((char *)&numPinZone, sizeof(numPinZone));
|
||||
if (!adjFile)
|
||||
{
|
||||
cout << "Problem on the zone " << zone << " / " << numZones << endl;
|
||||
return false;
|
||||
}
|
||||
|
||||
z.allZones[zone].pId.resize(numPinZone);
|
||||
adjFile.read((char *)&z.allZones[zone].pId[0], sizeof(int)*numPinZone);
|
||||
|
||||
if (maxParticlesInZone < 0 || numPinZone > maxParticlesInZone)
|
||||
maxParticlesInZone = numPinZone;
|
||||
|
||||
if (minParticlesInZone < 0 || numPinZone < minParticlesInZone)
|
||||
minParticlesInZone = numPinZone;
|
||||
|
||||
totalParticles += numPinZone;
|
||||
}
|
||||
cout << "Zoned " << totalParticles << endl;
|
||||
|
||||
cout << "Minimum number of particles in zone = " << minParticlesInZone << endl;
|
||||
cout << "Maximum number of particles in zone = " << maxParticlesInZone << endl;
|
||||
|
||||
if (totalParticles != numParticles)
|
||||
{
|
||||
cerr << "The numbers of particles are inconsistent ! (" << totalParticles << " vs " << numParticles << ")"<< endl;
|
||||
abort();
|
||||
}
|
||||
|
||||
volFile.read((char *)&numVoids, sizeof(numVoids));
|
||||
if (!volFile)
|
||||
return false;
|
||||
|
||||
cout << "Number of voids = " << numVoids << endl;
|
||||
|
||||
z.allVoids.resize(numVoids);
|
||||
for (int v = 0; v < numVoids; v++)
|
||||
{
|
||||
int32_t numZinV;
|
||||
|
||||
volFile.read((char *)&numZinV, sizeof(numZinV));
|
||||
if (!volFile)
|
||||
return false;
|
||||
|
||||
z.allVoids[v].zId.resize(numZinV);
|
||||
|
||||
int *zId = new int[numZinV];
|
||||
|
||||
volFile.read((char *)zId, sizeof(int)*numZinV);
|
||||
for (int k = 0; k < numZinV; k++)
|
||||
z.allVoids[v].zId[k] = zId[k];
|
||||
std::sort(&z.allVoids[v].zId[0], &z.allVoids[v].zId[numZinV]);
|
||||
|
||||
delete[] zId;
|
||||
}
|
||||
|
||||
if (volName != 0)
|
||||
{
|
||||
cout << "Loading particle volumes (requested)" << endl;
|
||||
ifstream f(volName);
|
||||
int numParticles;
|
||||
|
||||
if (!f)
|
||||
{
|
||||
cerr << "No such file " << volName << endl;
|
||||
abort();
|
||||
}
|
||||
|
||||
f.read((char *)&numParticles, sizeof(int));
|
||||
z.particleVolume.resize(numParticles);
|
||||
f.read((char *)&z.particleVolume[0], sizeof(float)*numParticles);
|
||||
}
|
||||
|
||||
cout << "Loading description" << endl;
|
||||
int lineid = 1;
|
||||
|
||||
string line;
|
||||
getline(descFile, line);
|
||||
lineid++;
|
||||
getline(descFile, line);
|
||||
lineid++;
|
||||
getline(descFile, line);
|
||||
lineid++;
|
||||
while (!descFile.eof())
|
||||
{
|
||||
istringstream lineStream(line.c_str());
|
||||
int orderId, volId, coreParticle, numParticlesInZone, numZonesInVoid, numInVoid;
|
||||
float coreDensity, volumeZone, volumeVoid, densityContrast;
|
||||
float probability;
|
||||
|
||||
lineStream
|
||||
>> orderId
|
||||
>> volId
|
||||
>> coreParticle
|
||||
>> coreDensity
|
||||
>> volumeZone
|
||||
>> numParticlesInZone
|
||||
>> numZonesInVoid
|
||||
>> volumeVoid
|
||||
>> numInVoid
|
||||
>> densityContrast
|
||||
>> probability;
|
||||
if (!lineStream)
|
||||
{
|
||||
cerr << "Error in text stream at line " << lineid << endl;
|
||||
abort();
|
||||
}
|
||||
|
||||
/*
|
||||
sscanf(line.c_str(), "%d %d %d %f %f %d %d %f %d %f %f\n", &orderId, &volId,
|
||||
&coreParticle, &coreDensity, &volumeZone, &numParticlesInZone,
|
||||
&numZonesInVoid,
|
||||
&volumeVoid, &numInVoid, &densityContrast, &probability);
|
||||
*/
|
||||
|
||||
z.allVoids[volId].proba = probability;
|
||||
z.allVoids[volId].volume = volumeVoid;
|
||||
z.allVoids[volId].numParticles = numInVoid;
|
||||
z.allVoids[volId].coreParticle = coreParticle;
|
||||
|
||||
// Sanity check
|
||||
int actualNumber = 0;
|
||||
for (int j = 0; j < z.allVoids[volId].zId.size(); j++)
|
||||
{
|
||||
int zzid = z.allVoids[volId].zId[j];
|
||||
assert(zzid < z.allZones.size());
|
||||
actualNumber += z.allZones[zzid].pId.size();
|
||||
}
|
||||
|
||||
if (actualNumber != numInVoid)
|
||||
{
|
||||
cerr << "Sanity check failed."
|
||||
<< " The number of particles in the description ("
|
||||
<< numInVoid
|
||||
<< ") is different from the one in the file ("
|
||||
<< actualNumber << ")" << endl;
|
||||
}
|
||||
getline(descFile, line);
|
||||
lineid++;
|
||||
}
|
||||
|
||||
cout << "Done loading" << endl;
|
||||
|
||||
|
||||
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
bool loadZobovParticles(const char *fname, std::vector<ZobovParticle>& particles)
|
||||
{
|
||||
UnformattedRead f(fname);
|
||||
int N;
|
||||
|
||||
f.beginCheckpoint();
|
||||
N = f.readInt32();
|
||||
f.endCheckpoint();
|
||||
|
||||
particles.resize(N);
|
||||
|
||||
f.beginCheckpoint();
|
||||
for (int i = 0; i < N; i++)
|
||||
{
|
||||
particles[i].x = f.readReal32();
|
||||
}
|
||||
f.endCheckpoint();
|
||||
|
||||
f.beginCheckpoint();
|
||||
for (int i = 0; i < N; i++)
|
||||
{
|
||||
particles[i].y = f.readReal32();
|
||||
}
|
||||
f.endCheckpoint();
|
||||
|
||||
f.beginCheckpoint();
|
||||
for (int i = 0; i < N; i++)
|
||||
{
|
||||
particles[i].z = f.readReal32();
|
||||
}
|
||||
f.endCheckpoint();
|
||||
|
||||
return true;
|
||||
}
|
60
c_source/util/loadZobov.hpp
Normal file
60
c_source/util/loadZobov.hpp
Normal file
|
@ -0,0 +1,60 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/loadZobov.hpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#ifndef __LOAD_ZOBOV_HPP
|
||||
#define __LOAD_ZOBOV_HPP
|
||||
|
||||
#include <vector>
|
||||
|
||||
struct ZobovZone
|
||||
{
|
||||
std::vector<int> pId;
|
||||
};
|
||||
|
||||
struct ZobovVoid
|
||||
{
|
||||
std::vector<int> zId;
|
||||
float proba;
|
||||
int numParticles, coreParticle;
|
||||
float volume;
|
||||
float barycenter[3];
|
||||
float nearestBoundary;
|
||||
};
|
||||
|
||||
struct ZobovRep
|
||||
{
|
||||
std::vector<ZobovZone> allZones;
|
||||
std::vector<ZobovVoid> allVoids;
|
||||
std::vector<float> particleVolume;
|
||||
};
|
||||
|
||||
struct ZobovParticle
|
||||
{
|
||||
float x, y, z;
|
||||
};
|
||||
|
||||
bool loadZobov(const char *descName,
|
||||
const char *adjName, const char *voidName,
|
||||
const char *volName, ZobovRep& z);
|
||||
|
||||
bool loadZobovParticles(const char *fname, std::vector<ZobovParticle>& particles);
|
||||
|
||||
#endif
|
137
c_source/util/particleInfo.cpp
Normal file
137
c_source/util/particleInfo.cpp
Normal file
|
@ -0,0 +1,137 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/particleInfo.cpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
#include "particleInfo.hpp"
|
||||
|
||||
#include <CosmoTool/fortran.hpp>
|
||||
#include <cstdlib>
|
||||
#include <netcdf>
|
||||
|
||||
using namespace std;
|
||||
using namespace CosmoTool;
|
||||
using namespace netCDF;
|
||||
|
||||
template <bool failure>
|
||||
double retrieve_attr_safe_double(NcFile& f, const char* name, double defval) {
|
||||
NcGroupAtt a = f.getAtt(name);
|
||||
if (a.isNull()) {
|
||||
if (failure) abort();
|
||||
return defval;
|
||||
}
|
||||
if (a.getAttLength() != 1) {
|
||||
abort();
|
||||
}
|
||||
double x;
|
||||
a.getValues(&x);
|
||||
return x;
|
||||
}
|
||||
|
||||
template <bool failure>
|
||||
int retrieve_attr_safe_int(NcFile& f, const char* name, int defval) {
|
||||
NcGroupAtt a = f.getAtt(name);
|
||||
if (a.isNull()) {
|
||||
if (failure) abort();
|
||||
return defval;
|
||||
}
|
||||
if (a.getAttLength() != 1) {
|
||||
abort();
|
||||
}
|
||||
int x;
|
||||
a.getValues(&x);
|
||||
return x;
|
||||
}
|
||||
|
||||
bool loadParticleInfo(ParticleInfo& info, const std::string& particles,
|
||||
const std::string& extra_info) {
|
||||
int numpart;
|
||||
int isObservation;
|
||||
|
||||
try {
|
||||
NcFile f_info(extra_info, NcFile::read);
|
||||
|
||||
info.ranges[0][0] =
|
||||
retrieve_attr_safe_double<true>(f_info, "range_x_min", 0);
|
||||
info.ranges[0][1] =
|
||||
retrieve_attr_safe_double<true>(f_info, "range_x_max", 0);
|
||||
info.ranges[1][0] =
|
||||
retrieve_attr_safe_double<true>(f_info, "range_y_min", 0);
|
||||
info.ranges[1][1] =
|
||||
retrieve_attr_safe_double<true>(f_info, "range_y_max", 0);
|
||||
info.ranges[2][0] =
|
||||
retrieve_attr_safe_double<true>(f_info, "range_z_min", 0);
|
||||
info.ranges[2][1] =
|
||||
retrieve_attr_safe_double<true>(f_info, "range_z_max", 0);
|
||||
info.mask_index = retrieve_attr_safe_int<true>(f_info, "mask_index", 0);
|
||||
isObservation = retrieve_attr_safe_int<false>(f_info, "is_observation", 0);
|
||||
|
||||
for (int i = 0; i < 3; i++)
|
||||
info.length[i] = info.ranges[i][1] - info.ranges[i][0];
|
||||
|
||||
try {
|
||||
UnformattedRead f(particles);
|
||||
|
||||
float mul, offset;
|
||||
|
||||
f.beginCheckpoint();
|
||||
numpart = f.readInt32();
|
||||
f.endCheckpoint();
|
||||
|
||||
info.particles.resize(numpart);
|
||||
|
||||
offset = info.ranges[0][0];
|
||||
// TEST PMS NON-COBIC BOXES
|
||||
// mul = 1.0;
|
||||
mul = info.ranges[0][1] - info.ranges[0][0];
|
||||
f.beginCheckpoint();
|
||||
for (int i = 0; i < numpart; i++)
|
||||
info.particles[i].x = mul * f.readReal32();
|
||||
f.endCheckpoint();
|
||||
|
||||
offset = info.ranges[1][0];
|
||||
// mul = 1.0;
|
||||
mul = info.ranges[1][1] - info.ranges[1][0];
|
||||
f.beginCheckpoint();
|
||||
for (int i = 0; i < numpart; i++)
|
||||
info.particles[i].y = mul * f.readReal32();
|
||||
f.endCheckpoint();
|
||||
|
||||
offset = info.ranges[2][0];
|
||||
// mul = 1.0;
|
||||
mul = info.ranges[2][1] - info.ranges[2][0];
|
||||
f.beginCheckpoint();
|
||||
for (int i = 0; i < numpart; i++)
|
||||
info.particles[i].z = mul * f.readReal32();
|
||||
f.endCheckpoint();
|
||||
|
||||
if (!isObservation) {
|
||||
for (int i = 0; i < numpart; i++) {
|
||||
info.particles[i].x += info.ranges[0][0];
|
||||
info.particles[i].y += info.ranges[1][0];
|
||||
info.particles[i].z += info.ranges[2][0];
|
||||
}
|
||||
}
|
||||
} catch (NoSuchFileException const& e) {
|
||||
return false;
|
||||
}
|
||||
|
||||
} catch (exceptions::NcCantRead const&) {
|
||||
return false;
|
||||
}
|
||||
return true;
|
||||
}
|
46
c_source/util/particleInfo.hpp
Normal file
46
c_source/util/particleInfo.hpp
Normal file
|
@ -0,0 +1,46 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/particleInfo.hpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#ifndef _PARTICLE_INFO_HEADER_HPP
|
||||
#define _PARTICLE_INFO_HEADER_HPP
|
||||
|
||||
#include <vector>
|
||||
#include <string>
|
||||
|
||||
struct ParticleData {
|
||||
float x, y, z;
|
||||
};
|
||||
|
||||
typedef std::vector<ParticleData> ParticleVector;
|
||||
|
||||
struct ParticleInfo
|
||||
{
|
||||
ParticleVector particles;
|
||||
float ranges[3][2];
|
||||
float length[3];
|
||||
int mask_index; // PMS
|
||||
};
|
||||
|
||||
bool loadParticleInfo(ParticleInfo& info,
|
||||
const std::string& particles,
|
||||
const std::string& extra_info);
|
||||
|
||||
#endif
|
367
c_source/util/voidTree.hpp
Normal file
367
c_source/util/voidTree.hpp
Normal file
|
@ -0,0 +1,367 @@
|
|||
/*+
|
||||
VIDE -- Void IDentification and Examination -- ./c_tools/libzobov/voidTree.hpp
|
||||
Copyright (C) 2010-2014 Guilhem Lavaux
|
||||
Copyright (C) 2011-2014 P. M. Sutter
|
||||
|
||||
This program is free software; you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation; version 2 of the License.
|
||||
|
||||
This program is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License along
|
||||
with this program; if not, write to the Free Software Foundation, Inc.,
|
||||
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
|
||||
+*/
|
||||
|
||||
|
||||
|
||||
#ifndef _VOID_TREE_HPP
|
||||
#define _VOID_TREE_HPP
|
||||
|
||||
#include <iostream>
|
||||
#include <stdint.h>
|
||||
#include "loadZobov.hpp"
|
||||
#include <list>
|
||||
#include <set>
|
||||
#include <vector>
|
||||
|
||||
struct VoidNode
|
||||
{
|
||||
int vid;
|
||||
VoidNode *parent;
|
||||
std::list<VoidNode *> children;
|
||||
};
|
||||
|
||||
struct VoidNodeOnDisk
|
||||
{
|
||||
int vid;
|
||||
int parent;
|
||||
};
|
||||
|
||||
class VoidTree
|
||||
{
|
||||
protected:
|
||||
uint32_t totalNumNodes, activeNodes;
|
||||
VoidNode *nodes;
|
||||
VoidNode *rootNode;
|
||||
ZobovRep& zobov;
|
||||
public:
|
||||
typedef std::list<VoidNode *> VoidList;
|
||||
|
||||
void dumpTree(std::ostream& o)
|
||||
{
|
||||
VoidNodeOnDisk data;
|
||||
|
||||
o.write((char*)&activeNodes, sizeof(uint32_t));
|
||||
for (uint32_t i = 0; i < activeNodes; i++)
|
||||
{
|
||||
data.vid = nodes[i].vid;
|
||||
if (nodes[i].parent == 0)
|
||||
data.parent = -1;
|
||||
else
|
||||
data.parent = nodes[i].parent - nodes;
|
||||
o.write((char *)&data, sizeof(VoidNodeOnDisk));
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
int lookupParent(int voidId, const std::vector<std::list<int> >& voids_for_zones)
|
||||
{
|
||||
int lastSize = 0x7fffffff;
|
||||
int goodParent = -1;
|
||||
ZobovVoid &ref_void = zobov.allVoids[voidId];
|
||||
const std::list<int>& candidateList = voids_for_zones[ref_void.zId.front()];
|
||||
std::list<int>::const_iterator iter_candidate = candidateList.begin();
|
||||
|
||||
// std::cout << "candidate list size is " << candidateList.size() << std::endl;
|
||||
|
||||
while (iter_candidate != candidateList.end())
|
||||
{
|
||||
int vid_candidate = *iter_candidate;
|
||||
|
||||
if (vid_candidate == voidId)
|
||||
break;
|
||||
++iter_candidate;
|
||||
}
|
||||
if (iter_candidate == candidateList.end())
|
||||
{
|
||||
// std::cout << "Failure to lookup parent" << std::endl;
|
||||
return -1;
|
||||
}
|
||||
|
||||
// voidId must be in the list.
|
||||
// assert(iter_candidate != candidateList.end());
|
||||
// Go back
|
||||
|
||||
iter_candidate = candidateList.end();
|
||||
|
||||
int vid_good_candidate = -1;
|
||||
int old_good_candidate_size = zobov.allZones.size()+1;
|
||||
|
||||
do
|
||||
{
|
||||
int vid_candidate;
|
||||
|
||||
--iter_candidate;
|
||||
|
||||
vid_candidate = *iter_candidate;
|
||||
std::vector<int>& candidate_zIds = zobov.allVoids[vid_candidate].zId;
|
||||
|
||||
if (voidId == vid_candidate)
|
||||
continue;
|
||||
|
||||
if (candidate_zIds.size() < ref_void.zId.size())
|
||||
{
|
||||
continue;
|
||||
}
|
||||
|
||||
int counter = 0;
|
||||
// All zones id are sorted in each void. So we just have parse the
|
||||
// vector of zones and check whether all the zones in ref_void.zId is
|
||||
// in iter_candidate->zId, the list is analyzed only once.
|
||||
// THOUGHT: candidateList may contain directly the information. It would suffice to have the void ids sorted according to volume. Then we just have to jump to the indice just smaller than voidId.
|
||||
|
||||
int k = 0;
|
||||
for (int j = 0; j < candidate_zIds.size() && k < ref_void.zId.size(); j++)
|
||||
{
|
||||
if (candidate_zIds[j] == ref_void.zId[k])
|
||||
k++;
|
||||
else if (candidate_zIds[j] > ref_void.zId[k])
|
||||
break;
|
||||
}
|
||||
if (k==ref_void.zId.size())
|
||||
{
|
||||
if (candidate_zIds.size() < old_good_candidate_size)
|
||||
{
|
||||
vid_good_candidate = vid_candidate;
|
||||
old_good_candidate_size = candidate_zIds.size();
|
||||
}
|
||||
// std::cout << "Found parent " << vid_candidate << std::endl;
|
||||
// return vid_candidate;
|
||||
}
|
||||
|
||||
// Go bigger, though I would say we should not to.
|
||||
}
|
||||
while (iter_candidate != candidateList.begin()) ;
|
||||
//if (vid_good_candidate < 0)
|
||||
// std::cout << "Failure to lookup parent (2)" << std::endl;
|
||||
return vid_good_candidate;
|
||||
}
|
||||
|
||||
VoidTree(ZobovRep& rep, std::istream& disk)
|
||||
: zobov(rep)
|
||||
{
|
||||
totalNumNodes = rep.allVoids.size();
|
||||
|
||||
disk.read((char *)&activeNodes, sizeof(uint32_t));
|
||||
nodes = new VoidNode[activeNodes];
|
||||
rootNode = 0;
|
||||
for (uint32_t i = 0; i < activeNodes; i++)
|
||||
{
|
||||
VoidNodeOnDisk data;
|
||||
|
||||
disk.read((char *)&data, sizeof(data));
|
||||
nodes[i].vid = data.vid;
|
||||
if (data.parent < 0)
|
||||
{
|
||||
if (rootNode != 0)
|
||||
{
|
||||
std::cerr << "Multiple root to the tree !!!" << std::endl;
|
||||
abort();
|
||||
}
|
||||
nodes[i].parent = 0;
|
||||
rootNode = &nodes[i];
|
||||
}
|
||||
else
|
||||
{
|
||||
nodes[i].parent = nodes + data.parent;
|
||||
nodes[i].parent->children.push_back(&nodes[i]);
|
||||
}
|
||||
}
|
||||
|
||||
computeMaxDepth();
|
||||
computeChildrenByNode();
|
||||
}
|
||||
|
||||
VoidTree(ZobovRep& rep)
|
||||
: zobov(rep)
|
||||
{
|
||||
totalNumNodes = rep.allVoids.size();
|
||||
|
||||
std::vector<std::list<int> > voids_for_zones;
|
||||
|
||||
voids_for_zones.resize(rep.allZones.size());
|
||||
for (int i = 0; i < rep.allVoids.size(); i++)
|
||||
{
|
||||
ZobovVoid& v = rep.allVoids[i];
|
||||
for (int j = 0; j < v.zId.size(); j++)
|
||||
voids_for_zones[v.zId[j]].push_back(i);
|
||||
}
|
||||
|
||||
// One additional for the mega-root
|
||||
nodes = new VoidNode[totalNumNodes+1];
|
||||
|
||||
for (int i = 0; i <= totalNumNodes; i++)
|
||||
{
|
||||
nodes[i].vid = i;
|
||||
nodes[i].parent = 0;
|
||||
}
|
||||
|
||||
std::cout << "Linking voids together..." << std::endl;
|
||||
double volMin = 0;// 4*M_PI/3*pow(4.*512/500.,3);
|
||||
int inserted = 0;
|
||||
for (int i = 0; i < rep.allVoids.size(); i++)
|
||||
{
|
||||
if (rep.allVoids[i].volume < volMin) continue;
|
||||
|
||||
int p = lookupParent(i, voids_for_zones);
|
||||
if ((i % 1000) == 0) std::cout << i << std::endl;
|
||||
|
||||
if (p >= 0)
|
||||
{
|
||||
nodes[p].children.push_back(&nodes[i]);
|
||||
nodes[i].parent = &nodes[p];
|
||||
}
|
||||
|
||||
inserted++;
|
||||
}
|
||||
|
||||
assert(inserted <= totalNumNodes);
|
||||
rootNode = &nodes[inserted];
|
||||
rootNode->vid = -1;
|
||||
rootNode->parent = 0;
|
||||
|
||||
for (int i = 0; i < inserted; i++)
|
||||
if (nodes[i].parent == 0)
|
||||
{
|
||||
nodes[i].parent = rootNode;
|
||||
rootNode->children.push_back(&nodes[i]);
|
||||
}
|
||||
activeNodes = inserted+1;
|
||||
|
||||
computeMaxDepth();
|
||||
computeChildrenByNode();
|
||||
}
|
||||
|
||||
~VoidTree()
|
||||
{
|
||||
delete[] nodes;
|
||||
}
|
||||
|
||||
int _depth_computer(VoidNode *node)
|
||||
{
|
||||
VoidList::iterator i = node->children.begin();
|
||||
int d = 0;
|
||||
|
||||
while (i != node->children.end())
|
||||
{
|
||||
d = std::max(d,_depth_computer(*i));
|
||||
++i;
|
||||
}
|
||||
|
||||
return d+1;
|
||||
}
|
||||
|
||||
void computeMaxDepth()
|
||||
{
|
||||
std::cout << "maximum depth is " << _depth_computer(rootNode) << std::endl;
|
||||
}
|
||||
|
||||
struct _children_stat {
|
||||
int num, min_num, max_num, num_zero,num_one, num_multiple;
|
||||
};
|
||||
|
||||
void _children_computer(VoidNode *node, _children_stat& s)
|
||||
{
|
||||
VoidList::iterator i = node->children.begin();
|
||||
int d = 0, j = 0;
|
||||
|
||||
while (i != node->children.end())
|
||||
{
|
||||
_children_computer(*i, s);
|
||||
++i;
|
||||
++j;
|
||||
}
|
||||
s.num += j;
|
||||
if (j!= 0)
|
||||
s.min_num = std::min(s.min_num, j);
|
||||
else s.num_zero ++;
|
||||
if (j==1) s.num_one++;
|
||||
if (j>1) s.num_multiple++;
|
||||
s.max_num = std::max(s.max_num, j);
|
||||
}
|
||||
|
||||
void computeChildrenByNode()
|
||||
{
|
||||
_children_stat s;
|
||||
s.num = 0;
|
||||
s.min_num = activeNodes+1;
|
||||
s.max_num = s.num_zero = s.num_one =s.num_multiple= 0;
|
||||
_children_computer(rootNode, s);
|
||||
std::cout << "Average children by node " << s.num*1.0/activeNodes << " , " << s.min_num << " " << s.max_num << " " << s.num_zero << " " << s.num_one << " " << s.num_multiple << std::endl;
|
||||
}
|
||||
|
||||
int getParent(int vid) const
|
||||
{
|
||||
assert(nodes[vid].parent != 0);
|
||||
return nodes[vid].parent->vid;
|
||||
}
|
||||
|
||||
const VoidList& getChildren(int vid) const
|
||||
{
|
||||
return nodes[vid].children;
|
||||
}
|
||||
|
||||
|
||||
VoidNode *getRoot() { return rootNode; }
|
||||
|
||||
template<typename T>
|
||||
void walkNode(VoidNode *node, T& traverse)
|
||||
{
|
||||
if (!traverse(node))
|
||||
return;
|
||||
|
||||
VoidList::iterator i = node->children.begin();
|
||||
|
||||
while (i != node->children.end())
|
||||
{
|
||||
walkNode(*i, traverse);
|
||||
++i;
|
||||
}
|
||||
}
|
||||
|
||||
template<typename T>
|
||||
void walk(T& traverse)
|
||||
{
|
||||
walkNode(rootNode, traverse);
|
||||
}
|
||||
|
||||
template<typename T, typename T2>
|
||||
void walkNodeWithMark(VoidNode *node, T& traverse, const T2& mark)
|
||||
{
|
||||
T2 new_mark = mark;
|
||||
|
||||
if (!traverse(node, new_mark))
|
||||
return;
|
||||
|
||||
VoidList::iterator i = node->children.begin();
|
||||
|
||||
while (i != node->children.end())
|
||||
{
|
||||
walkNodeWithMark(*i, traverse, new_mark);
|
||||
++i;
|
||||
}
|
||||
}
|
||||
|
||||
template<typename T,typename T2>
|
||||
void walkWithMark(T& traverse, T2 mark)
|
||||
{
|
||||
walkNodeWithMark(rootNode, traverse, mark);
|
||||
}
|
||||
};
|
||||
|
||||
#endif
|
Loading…
Add table
Add a link
Reference in a new issue