GLAMERDOC++
Gravitational Lensing Code Library
Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Friends | List of all members
LensHaloParticles< PType > Class Template Reference

A class that represents the lensing by a collection of simulation particles. More...

#include <particle_halo.h>

Inheritance diagram for LensHaloParticles< PType >:
Inheritance graph
[legend]
Collaboration diagram for LensHaloParticles< PType >:
Collaboration graph
[legend]

Public Member Functions

 LensHaloParticles (const std::string &simulation_filename, SimFileFormat format, PosType redshift, int Nsmooth, const COSMOLOGY &cosmo, Point_2d theta_rotate, bool recenter, bool my_multimass, double inv_area, PosType MinPSize=0, PosType rescale_mass=1.0, bool verbose=false)
 
 LensHaloParticles (std::vector< PType > &pvector, float redshift, const COSMOLOGY &cosmo, Point_2d theta_rotate, bool recenter, double my_inv_area, float MinPSize=0, double max_range=-1, bool verbose=false)
 
 LensHaloParticles (PType *begin, size_t n, float redshift, const COSMOLOGY &cosmo, Point_2d theta_rotate, bool recenter, double my_inv_area, float MinPSize=0, double max_range=-1, bool verbose=false)
 this constructor does not take possession of the particles More...
 
 LensHaloParticles (LensHaloParticles &&h)
 
LensHaloParticles< PType > & operator= (LensHaloParticles< PType > &&h)
 
void force_halo (double *alpha, KappaType *kappa, KappaType *gamma, KappaType *phi, double const *xcm, bool subtract_point=false, PosType screening=1.0)
 
size_t getN () const
 
void rotate (Point_2d theta)
 rotate the simulation around the origin of the simulation coordinates, (radians) More...
 
Point_3d CenterOfMass ()
 get current center of mass in input coordinates
 
Point_3d DensestPoint ()
 get the densistt point in input coordinates
 
void readPositionFileASCII (const std::string &filename)
 Reads number of particle and particle positons into Npoint and xp from a ASCII file. More...
 
- Public Member Functions inherited from LensHalo
 LensHalo ()
 Shell constructor.
 
 LensHalo (PosType z, const COSMOLOGY &cosmo)
 
 LensHalo (const LensHalo &h)
 
 LensHalo (LensHalo &&h)
 
LensHalooperator= (const LensHalo &h)
 
LensHalooperator= (LensHalo &&h)
 
float get_Rmax () const
 this can be used to tag types of LensHalos More...
 
float getRsize () const
 get the Rsize which is the size of the halo in Mpc
 
float get_mass () const
 get the mass solar units
 
float get_rscale () const
 get the scale radius in Mpc
 
PosType getZlens () const
 get the redshift
 
void getX (PosType *MyPosHalo) const
 get the position of the Halo in physical Mpc on the lens plane
 
PosType operator[] (int i) const
 returns position of the Halo in physical Mpc on the lens plane
 
void setTheta (PosType PosX, PosType PosY)
 set the position of the Halo in radians
 
void setTheta (PosType *PosXY)
 set the position of the Halo in radians
 
void setTheta (const Point_2d &p)
 set the position of the Halo in radians
 
void getTheta (PosType *MyPosHalo) const
 get the position of the Halo in radians
 
void setDist (COSMOLOGY &co)
 Set the angular size distance to the halo. This should be the distance to the lens plane.
 
PosType getDist () const
 
void displayPos ()
 
virtual void initFromFile (float my_mass, long *seed, float vmax, float r_halfmass)
 initialize from a simulation file
 
virtual void initFromMassFunc (float my_mass, float my_Rsize, float my_rscale, PosType my_slope, long *seed)
 initialize from a mass function
 
virtual void set_RsizeRmax (float my_Rsize)
 set Rsize (in Mpc) and reset Rmax
 
void set_mass (float my_mass)
 set mass (in solar masses)
 
virtual void set_rscale (float my_rscale)
 set scale radius (in Mpc)
 
void setZlens (PosType my_zlens, const COSMOLOGY &cosmo)
 set redshift
 
void setRsize (PosType R)
 
void setZlensDist (PosType my_zlens, const COSMOLOGY &cos)
 
void setMass (PosType m)
 
virtual void set_slope (PosType my_slope)
 set slope
 
virtual PosType get_slope ()
 get slope
 
bool get_flag_elliptical ()
 flag=True if halo elliptical
 
void set_flag_elliptical (bool ell)
 
bool get_switch_flag ()
 
void set_switch_flag (bool swt)
 flag permits case distinction in force_halo_asym for elliptical NFWs only (get_switch_flag==true), in latter case the mass_norm_factor^2 is used instead of mass_norm_factor.
 
virtual void setCosmology (const COSMOLOGY &cosmo)
 used for elliptical NFWs only, in that case get_switch_flag==true More...
 
bool compareZ (PosType z)
 force tree calculation for stars More...
 
EllipMethod getEllipMethod () const
 stars More...
 
std::vector< double > get_mod ()
 get vector of Fourier modes, which are calculated in the constructors of the LensHaloes when main_ellip_method is set to 'Fourier'
 
virtual std::size_t Nparams () const
 get the number of halo parameters
 
virtual PosType getParam (std::size_t p) const
 get the value of a scaled halo parameter by index
 
virtual PosType setParam (std::size_t p, PosType value)
 set the value of a scaled halo parameter by index
 
virtual void printCSV (std::ostream &, bool header=false) const
 print the halo parameters in CSV format
 
PosType MassBy2DIntegation (PosType R)
 Prints star parameters; if show_stars is true, prints data for single stars. More...
 
PosType MassBy1DIntegation (PosType R)
 calculates the mass within radius R by integating alpha on a ring and using Gauss' law, used only for testing
 
PosType test_average_gt (PosType R)
 calculates the average gamma_t for LensHalo::test()
 
PosType test_average_kappa (PosType R)
 
void set_norm_factor ()
 
void set_rsize (float my_rsize)
 set radius rsize beyond which interpolation values between alpha_ellip and alpha_iso are computed
 
float get_rsize ()
 
bool test ()
 perform some basic consistancy checks for halo More...
 
size_t getID () const
 
void setID (size_t id)
 
PosType renormalization (PosType r_max)
 
PixelMap map_variables (LensingVariable lensvar, size_t Nx, size_t Ny, double res)
 Map a PixelMap of the surface, density, potential and potential gradient centred on (0,0) in LensHalo coordinates. More...
 

Static Public Member Functions

static void makeSIE (std::string new_filename, PosType redshift, double particle_mass, double total_mass, double sigma, double q, Utilities::RandomNumbers_NR &ran)
 This is a test static function that makes a truncated SIE out of particles and puts it into a file in the right format for constructing a LensHaloParticles. More...
 
static LensHaloParticles< ParticleTypeSimpleSIE (PosType redshift, double particle_mass, double total_mass, double sigma, double q, int Nneighbors, COSMOLOGY &cosmo, Utilities::RandomNumbers_NR &ran)
 This is a test static function that makes a truncated SIE out of particles and puts it into a file in the right format for constructing a LensHaloParticles. More...
 
static void calculate_smoothing (int Nsmooth, PType *pp, size_t Npoints, bool verbose=false)
 
static void writeSizes (const std::string &filename, int Nsmooth, const PType *pp, size_t Npoints)
 
static bool readSizesFile (const std::string &filename, PType *pp, size_t Npoints, int Nsmooth, PosType min_size)
 
- Static Public Member Functions inherited from LensHalo
static const int get_Nmod ()
 get length of mod array, which is Nmod. Not to be confused with getNmodes in the class LensHaloFit
 

Protected Member Functions

 LensHaloParticles (float redshift, const COSMOLOGY &cosmo)
 
 LensHaloParticles (PType *pdata, size_t Nparticles, float redshift, const COSMOLOGY &cosmo, Point_2d theta_rotate, bool recenter, double my_inv_area, bool verbose=false)
 
void rotate_particles (PosType theta_x, PosType theta_y)
 
void assignParams (InputParams &params)
 
void set_up (float redshift, const COSMOLOGY &cosmo, Point_2d theta_rotate, double max_range, bool recenter, bool verbose)
 
- Protected Member Functions inherited from LensHalo
PosType alpha_int (PosType x) const
 Calculates potential (phi_int) from alpha_h. If flag is_alphah_a_table is True it takes and integrates directly the gfunction instead of alpha_h. The gfunction is used for the InterpolationTable used in alpha_h. Setting the flag to False speeds up the calculation of phi_h.
 
PosType norm_int (PosType r_max)
 
void force_halo_sym (PosType *alpha, KappaType *kappa, KappaType *gamma, KappaType *phi, PosType const *xcm, bool subtract_point=false, PosType screening=1.0)
 returns the lensing quantities of a ray in center of mass coordinates for a symmetric halo More...
 
void force_halo_asym (PosType *alpha, KappaType *kappa, KappaType *gamma, KappaType *phi, PosType const *xcm, bool subtract_point=false, PosType screening=1.0)
 
bool force_point (PosType *alpha, KappaType *kappa, KappaType *gamma, KappaType *phi, PosType const *xcm, PosType rcm2, bool subtract_point, PosType screening)
 
void assignParams (InputParams &params, bool needRsize)
 read in parameters from a parameterfile in InputParams params
 
void error_message1 (std::string name, std::string filename)
 read in star parameters. This is valid for all halos and not overloaded. More...
 
virtual PosType alpha_h (PosType x) const
 
virtual KappaType kappa_h (PosType x) const
 
virtual KappaType gamma_h (PosType x) const
 
virtual KappaType phi_h (PosType x) const
 
virtual KappaType phi_int (PosType x) const
 
virtual PosType ffunction (PosType x) const
 
virtual PosType gfunction (PosType x) const
 
virtual PosType dgfunctiondx (PosType x)
 
virtual PosType bfunction (PosType x)
 
virtual PosType dhfunction (PosType x) const
 
virtual PosType ddhfunction (PosType x, bool numerical)
 
virtual PosType dddhfunction (PosType x, bool numerical)
 
virtual PosType bnumfunction (PosType x)
 
virtual PosType dbfunction (PosType x)
 
virtual PosType ddbfunction (PosType x)
 
virtual PosType dmoddb (int whichmod, PosType q, PosType b)
 
virtual PosType ddmoddb (int whichmod, PosType q, PosType b)
 
virtual PosType dmoddq (int whichmod, PosType q, PosType b)
 
virtual PosType ddmoddq (int whichmod, PosType q, PosType b)
 
void faxial (PosType x, PosType theta, PosType f[])
 If set to true the correct normalization is applied for asymmetric NFW profiles, the mass_norm_factor is different for the other halos. More...
 
void faxial0 (PosType theta, PosType f0[])
 
void faxial1 (PosType theta, PosType f1[])
 
void faxial2 (PosType theta, PosType f2[])
 
void gradial (PosType r, PosType g[])
 Derivatives of the potential damping factor with respect to r ... TODO: come up with a better damping faction.
 
void gradial2 (PosType r, PosType mu, PosType sigma, PosType g[])
 
void felliptical (PosType x, PosType q, PosType theta, PosType f[], PosType g[])
 Calculate the derivatives of the G function = r*sqrt(cos(theta)^2 + q(r)^2 sin(theta)) More...
 
virtual void gamma_asym (PosType x, PosType theta, PosType gamma[])
 
virtual PosType kappa_asym (PosType x, PosType theta)
 
virtual void alphakappagamma_asym (PosType x, PosType theta, PosType alpha[], PosType *kappa, PosType gamma[], PosType *phi)
 Pseudo-elliptical profiles by Phi(G)-Ansatz. More...
 
virtual void alphakappagamma1asym (PosType x, PosType theta, PosType alpha[2], PosType *kappa, PosType gamma[], PosType *phi)
 Elliptical profiles by Fourier-Ansatz. More...
 
virtual void alphakappagamma2asym (PosType x, PosType theta, PosType alpha[2], PosType *kappa, PosType gamma[], PosType *phi)
 
virtual void alphakappagamma3asym (PosType x, PosType theta, PosType alpha[2], PosType *kappa, PosType gamma[], PosType *phi)
 
virtual PosType alpha_ell (PosType x, PosType theta)
 
double fourier_coeff (double n, double q, double beta)
 Calculates fourier-coefficients for power law halo.
 
double IDAXDM (double lambda, double a2, double b2, double x[], double rmax, double mo)
 
double IDAYDM (double lambda, double a2, double b2, double x[], double rmax, double mo)
 
double SCHRAMMKN (double n, double x[], double rmax)
 
double SCHRAMMJN (double n, double x[], double rmax)
 
double SCHRAMMI (double x[], double rmax)
 
void calcModes (double q, double beta, double rottheta, PosType newmod[])
 Calculates the modes for fourier expansion of power law halo. All the modes are relative to the zero mode to conserve mass throughout the calculation of kappa etc.
 
void calcModesB (PosType x, double q, double beta, double rottheta, PosType newmod[])
 
void calcModesC (PosType beta_r, double q, double rottheta, PosType newmod[])
 
virtual PosType InterpolateModes (int whichmod, PosType q, PosType b)
 
void analModes (int modnumber, PosType my_beta, PosType q, PosType amod[3])
 

Static Protected Member Functions

static void smooth_ (TreeSimple< PType > *tree3d, PType *xp, size_t N, int Nsmooth)
 

Protected Attributes

Point_3d mcenter
 
Point_3d densest_point
 
PType * pp
 
std::vector< PType > trash_collector
 
PosType min_size
 
bool multimass
 
Utilities::Geometry::SphericalPoint center
 
size_t Npoints
 
PosType inv_area
 
std::string simfile
 
std::string sizefile
 
TreeQuadParticles< PType > * qtree
 
- Protected Attributes inherited from LensHalo
float Rsize = 0
 
float mass
 
PosType Dist
 
PosType mnorm
 
float Rmax
 
PosType beta
 
float Rmax_to_Rsize_ratio = 1.2
 The factor by which Rmax is larger than Rsize.
 
float rscale
 scale length or core size. Different meaning in different cases. Not used in NSIE case.
 
EllipMethod main_ellip_method
 
PosType xmax
 
PosType mass_norm_factor =1
 This is Rsize/rscale !!
 
float pa
 
float fratio =1.0
 
bool elliptical_flag = false
 
bool switch_flag = false
 
PosType mod [Nmod]
 
PosType mod1 [Nmod]
 
PosType mod2 [Nmod]
 
PosType r_eps
 

Friends

class MakeParticleLenses
 

Additional Inherited Members

- Public Attributes inherited from LensHalo
int tag =0
 
- Static Protected Attributes inherited from LensHalo
static const int Nmod = 32
 

Detailed Description

template<typename PType>
class LensHaloParticles< PType >

A class that represents the lensing by a collection of simulation particles.

You can create a LensHaloParticles<> directly from a file, but it is recommended that you use the MakeParticleLenses class to create them and then move them to a Lens object.

Smoothing is done according to the density of particles in 3D. Smoothing sizes are either read in from a file (names simulation_filename + "." + Nsmooth + "sizes") or calculated if the file does not exist (in which case the file is created). This can be time and memory consuming when there are a large number of particles.

Input format: ASCII - a table of 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.

More input formats will be added in the future.

Constructor & Destructor Documentation

◆ LensHaloParticles() [1/5]

template<typename PType >
LensHaloParticles< PType >::LensHaloParticles ( const std::string &  simulation_filename,
SimFileFormat  format,
PosType  redshift,
int  Nsmooth,
const COSMOLOGY &  cosmo,
Point_2d  theta_rotate,
bool  recenter,
bool  my_multimass,
double  inv_area,
PosType  MinPSize = 0,
PosType  rescale_mass = 1.0,
bool  verbose = false 
)
Parameters
simulation_filenamename of data files
formatformat of data file
redshiftredshift of origin
Nsmoothnumber of neighbours for adaptive smoothing
cosmocosmology
theta_rotaterotation of particles around the origin
recentercenter on center of mass
my_multimassset to true is particles have different sizes
inv_areainverse area for mass compensation
MinPSizeminimum particle size in Mpc, overrides nearest neighbors aor input sizes for small particle sizes
rescale_massrescale particle masses

◆ LensHaloParticles() [2/5]

template<typename PType >
LensHaloParticles< PType >::LensHaloParticles ( std::vector< PType > &  pvector,
float  redshift,
const COSMOLOGY &  cosmo,
Point_2d  theta_rotate,
bool  recenter,
double  my_inv_area,
float  MinPSize = 0,
double  max_range = -1,
bool  verbose = false 
)
inline
Parameters
pvectorlist of particles pdata[][i] should be the position in physical Mpc, the class takes possession of the data and leaves the vector empty
redshiftredshift of origin
cosmocosmology
theta_rotaterotation of particles around the origin
recentercenter on center of mass
my_inv_areainverse area for mass compensation
MinPSizeminimum particle size
max_rangeif set this will cause the tree not be fully construct down to the bucket size outside this range

◆ LensHaloParticles() [3/5]

template<typename PType >
LensHaloParticles< PType >::LensHaloParticles ( PType *  begin,
size_t  n,
float  redshift,
const COSMOLOGY &  cosmo,
Point_2d  theta_rotate,
bool  recenter,
double  my_inv_area,
float  MinPSize = 0,
double  max_range = -1,
bool  verbose = false 
)
inline

this constructor does not take possession of the particles

Parameters
beginpointer to first particles pdata[][i] should be the position in physical Mpc, the class takes possession of the data and leaves the vector empty
redshiftredshift of origin
cosmocosmology
theta_rotaterotation of particles around the origin
recentercenter on center of mass
my_inv_areainverse area for mass compensation
MinPSizeminimum particle size
max_rangeif set this will cause the tree not be fully construct down to the bucket size outside this range

◆ LensHaloParticles() [4/5]

template<typename PType >
LensHaloParticles< PType >::LensHaloParticles ( float  redshift,
const COSMOLOGY &  cosmo 
)
inlineprotected
Parameters
redshiftredshift of origin
cosmocosmology

◆ LensHaloParticles() [5/5]

template<typename PType >
LensHaloParticles< PType >::LensHaloParticles ( PType *  pdata,
size_t  Nparticles,
float  redshift,
const COSMOLOGY &  cosmo,
Point_2d  theta_rotate,
bool  recenter,
double  my_inv_area,
bool  verbose = false 
)
inlineprotected
Parameters
pdataparticle data (all physical distances)
redshiftredshift of origin
cosmocosmology
theta_rotaterotation of particles around the origin
recentercenter on center of mass
my_inv_areainverse area for mass compensation

Member Function Documentation

◆ makeSIE()

template<typename PType >
void LensHaloParticles< PType >::makeSIE ( std::string  new_filename,
PosType  redshift,
double  particle_mass,
double  total_mass,
double  sigma,
double  q,
Utilities::RandomNumbers_NR ran 
)
static

This is a test static function that makes a truncated SIE out of particles and puts it into a file in the right format for constructing a LensHaloParticles.

This is useful for calculating the level of shot noise and finite source size. The particles are distributed in 3D according to the SIE profile with only the perpendicular coordinates (1st and 2nd) distorted into an elliptical shape. If the halo is rotated from the original orientation it will not be a oblate spheroid.

Parameters
new_filenamefile name to store the particles
redshiftredshift of particles
particle_massparticle mass
total_masstotal mass of SIE
sigmavelocity dispersion in km/s
qaxis ratio

◆ readPositionFileASCII()

template<typename PType >
void LensHaloParticles< PType >::readPositionFileASCII ( const std::string &  filename)

Reads number of particle and particle positons into Npoint and xp from a ASCII file.

Data file must have the lines "# nparticles ***" and "# mass ***" in the header. All header lines must begin with a "# "

Coordinates of particles are in physical Mpc units.

◆ rotate()

template<typename PType >
void LensHaloParticles< PType >::rotate ( Point_2d  theta)

rotate the simulation around the origin of the simulation coordinates, (radians)

rotate simulation

◆ set_up()

template<typename PType >
void LensHaloParticles< PType >::set_up ( float  redshift,
const COSMOLOGY &  cosmo,
Point_2d  theta_rotate,
double  max_range,
bool  recenter,
bool  verbose 
)
protected
Parameters
redshiftredshift of origin
cosmocosmology
theta_rotaterotation of particles around the origin
recentercenter on center of mass

◆ SIE()

template<typename PType >
LensHaloParticles< ParticleTypeSimple > LensHaloParticles< PType >::SIE ( PosType  redshift,
double  particle_mass,
double  total_mass,
double  sigma,
double  q,
int  Nneighbors,
COSMOLOGY &  cosmo,
Utilities::RandomNumbers_NR ran 
)
static

This is a test static function that makes a truncated SIE out of particles and puts it into a file in the right format for constructing a LensHaloParticles.

This is useful for calculating the level of shot noise and finite source size. The particles are distributed in 3D according to the SIE profile with only the perpendicular coordinates (1st and 2nd) distorted into an elliptical shape. If the halo is rotated from the original orientation it will not be a oblate spheroid.

Parameters
redshiftredshift of particles
particle_massparticle mass
total_masstotal mass of SIE
sigmavelocity dispersion in km/s
qaxis ratio
Nneighborsnumber of neighbor particles used in smoothing

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