GLAMERDOC++
Gravitational Lensing Code Library
Public Member Functions | Public Attributes | List of all members
MakeParticleLenses Class Reference

A class for constructing LensHalos from particles in a data file. More...

#include <particle_halo.h>

Public Member Functions

std::vector< size_t > getnp ()
 returns number of particles of each type
 
std::vector< float > getmp ()
 returns mass of particles of each type, if 0 they can have different masses
 
double sim_redshift ()
 returns original redshift of snapshot, redshifts of the halos can be changed
 
 MakeParticleLenses (const std::string &filename, SimFileFormat format, int Nsmooth, bool recenter, bool ignore_type_in_smoothing=false)
 
 MakeParticleLenses (const std::string &filename, bool recenter)
 
void Recenter (Point_3d<> x)
 recenter the particles to 3d point in physical Mpc/h units If the halos have already been created they will be destroyed.
 
void CreateHalos (const COSMOLOGY &cosmo, double redshift, double inv_area)
 
void radialCut (Point_3d<> center, double radius)
 remove particles that are beyond radius (Mpc/h) from center More...
 
void cylindricalCut (Point_2d center, double radius)
 remove particles that are beyond cylindrical radius (Mpc/h) of center More...
 
Point_3d getCenterOfMass () const
 returns the original center of mass of all the particles
 
Point_3d densest_particle () const
 returns the location of the densest particle in (Mpc/h)
 
void getBoundingBox (Point_3d<> &Xmin, Point_3d<> &Xmax) const
 return the maximum and minimum coordinates of the particles in each dimension in for the original simulation in Mpc/h
 
double getZoriginal ()
 

Public Attributes

std::vector< LensHaloParticles< ParticleType< float > > * > halos
 vector of LensHalos, one for each type of particle type
 
std::vector< ParticleType< float > > data
 

Detailed Description

A class for constructing LensHalos from particles in a data file.

The particle data is stored in this structure so the LensHaloParticles should not be copied and then this object allowed to be destroyed. The halos will be destroyed when this structure is destroyed.

A separate LensHaloParticles is made for each type of particle that is present in the gadget file. The nearest N neighbour smoothing is done in 3D on construction separately for each type of particle. The smoothing sizes are automatically saved to files and used again if the class is constructed again with the same file and smoothing number.

On construction the LensHalos are not constructed. You nead to run the CreateHalos() method the them to be created.

The data file formats are :

gadget2 - standard Gadget-2 output. A different LensHalo is made for each type of particle. The nearest neighbour smoothing scales are calculated within particle type. The sph density of gas particles are not used.

csv3,csv4,csv5,csv6 - CSV ascii format without header. The first three columns are the positions. Next columns are used for the other formats being and interpreted as (column 4) masses are in Msun/h, (column 5) the paricle smoothing size in Mpc/h and (column 6) an integer for type of particle. There can be more columns in the file than are uesed. In the case of csv6, when there are more then one type of halo each type will be in a differeent LensHaloParticles with differnt smoothing.

glmb - This is a binary format internal to GLAMER used to store the positions, masses and sizes of the particles. If GLAMER has generated one, it should be all that is needed to recreate the LensHaloParticles.

ascii2 - This is the original ascii GLAMER format. Three floats for positions in comoving Mpc (no h factor). The lines "# nparticles ...." and "# mass ...." must be in header at the top of the file. # is otherwise a comment character. Only one type of particle in a single input file.

example of use:

COSMOLOGY cosmo(Planck);

double zs = 2,zl = 0.5; double Dl = cosmo.coorDist(zl);

std::string filename = "DataFiles/snap_058_centered.txt";

MakeParticleLenses halomaker(filename,csv4,30,false);

Point_3d Xmax,Xmin;

halomaker.getBoundingBox(Xmin, Xmax);

Point_3d c_mass = halomaker.getCenterOfMass(); Point_2d center;

center[0] = c_mass[0]; center[1] = c_mass[1];

cut out a cylinder, could also do a ball halomaker.cylindricalCut(center,(Xmax[0]-Xmin[0]/2));

long seed = 88277394; Lens lens(&seed,zs);

double range = (Xmax[0]-Xmin[0])*1.05*cosmo.gethubble()/Dl; // angular range of simulation

center *= cosmo.gethubble()/Dl; // convert to angular coordinates

halomaker.CreateHalos(cosmo,zl);

for(auto h : halomaker.halos){ lens.insertMainHalo(h,zl, true); }

GridMap gridmap(&lens, 2049,center.x,range);

PixelMap pmap = gridmap.writePixelMapUniform(KAPPA); pmap.printFITS("!" + filename + ".kappa.fits");

pmap = gridmap.writePixelMapUniform(ALPHA1); pmap.printFITS("!" + filename + ".alpha1.fits");

pmap = gridmap.writePixelMapUniform(ALPHA2); pmap.printFITS("!" + filename + ".alpha2.fits");

< >

Constructor & Destructor Documentation

◆ MakeParticleLenses() [1/2]

MakeParticleLenses::MakeParticleLenses ( const std::string &  filename,
SimFileFormat  format,
int  Nsmooth,
bool  recenter,
bool  ignore_type_in_smoothing = false 
)
Parameters
filenamepath / root name of gadget-2 snapshot
formatformat of data file
Nsmoothnumber of nearest neighbors used for smoothing
recenterrecenter so that the LensHalos are centered on the center of mass
ignore_type_in_smoothingused only when format == gadget2, nearest neighbour smoothing is done amongst particles by type if set to false

◆ MakeParticleLenses() [2/2]

MakeParticleLenses::MakeParticleLenses ( const std::string &  filename,
bool  recenter 
)
Parameters
filenamepath / name of glmb file
recenterrecenter so that the LensHalos are centered on the center of mass

Member Function Documentation

◆ CreateHalos()

void MakeParticleLenses::CreateHalos ( const COSMOLOGY &  cosmo,
double  redshift,
double  inv_area 
)

how do you set up inv_area in general ???

◆ cylindricalCut()

void MakeParticleLenses::cylindricalCut ( Point_2d  center,
double  radius 
)

remove particles that are beyond cylindrical radius (Mpc/h) of center

resort by type

◆ radialCut()

void MakeParticleLenses::radialCut ( Point_3d<>  center,
double  radius 
)

remove particles that are beyond radius (Mpc/h) from center

resort by type


The documentation for this class was generated from the following files: