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

A base class for all types of lensing "halos" which are any mass distribution that cause lensing. More...

#include <lens_halos.h>

Inheritance diagram for LensHalo:
Inheritance graph
[legend]

Public Member Functions

 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...
 
virtual void force_halo (PosType *alpha, KappaType *kappa, KappaType *gamma, KappaType *phi, PosType const *xcm, bool subtract_point=false, PosType screening=1.0)
 
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 const int get_Nmod ()
 get length of mod array, which is Nmod. Not to be confused with getNmodes in the class LensHaloFit
 

Public Attributes

int tag =0
 

Protected Member Functions

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])
 

Protected Attributes

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
 

Static Protected Attributes

static const int Nmod = 32
 

Detailed Description

A base class for all types of lensing "halos" which are any mass distribution that cause lensing.

It contains the mass, maximum radius (Rsize), and scale radius (rscale) of the halo, as well as the redshift (zlens).

It has get and set functions for the members as well as virtual functions like: force_halo that compute the lensing properties – deflection, convergence, and shear.

Along with the simple set function, there are two general initialization functions, that calculate the rest of the properties based on some input lens halo parameter (e.g. mass).

initFromFile is intended to be used when the halo data is read in from a simulation file. In this case the general halo is assumed to be an NFW halo and therefore the maximum velocity and the half-mass radius need to be set. This function is overridden in derived classes and in cases where it is not applicable only the mass is taken into initializing the lensing halo.

initFromMassFunc is intended for the cases where the simulation is populated by lensing halos from a mass function. Then one needs all parameters of the halo – mass, Rsize, and rscale.

Constructor & Destructor Documentation

◆ LensHalo()

LensHalo::LensHalo ( const LensHalo h)

Identification number of halo. It is not always used.

Member Function Documentation

◆ alpha_h()

virtual PosType LensHalo::alpha_h ( PosType  x) const
inlineprotectedvirtual

point mass case r |alpha(r)| pi Sigma_crit / Mass

Reimplemented in LensHaloNFW.

◆ alphakappagamma1asym()

void LensHalo::alphakappagamma1asym ( PosType  x,
PosType  theta,
PosType  alpha[2],
PosType *  kappa,
PosType  gamma[],
PosType *  phi 
)
protectedvirtual

Elliptical profiles by Fourier-Ansatz.

with damping

◆ alphakappagamma_asym()

void LensHalo::alphakappagamma_asym ( PosType  r,
PosType  theta,
PosType  alpha[],
PosType *  kappa,
PosType  gamma[],
PosType *  phi 
)
protectedvirtual

Pseudo-elliptical profiles by Phi(G)-Ansatz.

Parameters
rRadius in Mpc (not scale lengths)
thetaangle of ray
alphaoutput deflection
kappaoutput kappa
gammaoutpot gamma

◆ analModes()

void LensHalo::analModes ( int  modnumber,
PosType  my_beta,
PosType  q,
PosType  amod[3] 
)
protected

Created on: March 28, 2014 Author: D. Leier

◆ compareZ()

bool LensHalo::compareZ ( PosType  z)
inline

force tree calculation for stars

internal compare redshift function

◆ displayPos()

void LensHalo::displayPos ( )
inline

get the position of the Halo in physical Mpc display the position of the halo

◆ error_message1()

void LensHalo::error_message1 ( std::string  name,
std::string  filename 
)
protected

read in star parameters. This is valid for all halos and not overloaded.

error message printout

◆ faxial()

void LensHalo::faxial ( PosType  x,
PosType  theta,
PosType  f[] 
)
protected

If set to true the correct normalization is applied for asymmetric NFW profiles, the mass_norm_factor is different for the other halos.

This function returns the lensing quantities for an asymmetric version of the symmetric baseclass halo.

This function should only be used by the second generation of classes derived from LensHalo.

The math needs to be PosType checked and the sign convention checked. The method used to make the lenses asymmetric is laid out in http://metcalf1.bo.astro.it/~bmetcalf/ExtraNotes/notes_elliptical.pdf Derivatives of the potential factor with respect to theta

◆ felliptical()

void LensHalo::felliptical ( PosType  x,
PosType  q,
PosType  theta,
PosType  f[],
PosType  g[] 
)
protected

Calculate the derivatives of the G function = r*sqrt(cos(theta)^2 + q(r)^2 sin(theta))

output: f[0] = G, f[1] = dG/dtheta, f[2] = ddg/dtheta^2, g[0] = R/r, g[1] = dG/dr , g[2] = ddG/dr^2, g[3] = ddG/dr/dtheta

◆ force_halo()

void LensHalo::force_halo ( PosType *  alpha,
KappaType *  kappa,
KappaType *  gamma,
KappaType *  phi,
PosType const *  xcm,
bool  subtract_point = false,
PosType  screening = 1.0 
)
virtual
Parameters
alphadeflection solar mass/Mpc
kappasurface density in units of Sigma crit
gammathree components of shear
phipotential in solar masses
xcmposition relative to center (in physical Mpc)
subtract_pointif true contribution from a point mass is subtracted
screeningthe factor by which to scale the mass for screening of the point mass subtraction

Reimplemented in LensHaloDummy, LensHaloTEPL, LensHaloTNSIE, LensHaloRealNSIE, LensHaloUniform, and LensHaloMassMap.

◆ force_halo_asym()

void LensHalo::force_halo_asym ( PosType *  alpha,
KappaType *  kappa,
KappaType *  gamma,
KappaType *  phi,
PosType const *  xcm,
bool  subtract_point = false,
PosType  screening = 1.0 
)
protected

intersecting, subtract the point particle

case distinction used for elliptical NFWs (true) only (get_switch_flag==true)

add stars for microlensing

Parameters
alphamass/Mpc
phipotential solar masses
subtract_pointif true contribution from a point mass is subtracted
screeningthe factor by which to scale the mass for screening of the point mass subtraction

◆ force_halo_sym()

void LensHalo::force_halo_sym ( PosType *  alpha,
KappaType *  kappa,
KappaType *  gamma,
KappaType *  phi,
PosType const *  xcm,
bool  subtract_point = false,
PosType  screening = 1.0 
)
protected

returns the lensing quantities of a ray in center of mass coordinates for a symmetric halo

phi is defined here in such a way that it differs from alpha by a sign (as we should have alpha = \nabla_x phi) but alpha agrees with the rest of the lensing quantities (kappa and gammas). Warning : Be careful, the sign of alpha is changed in LensPlaneSingular::force !

intersecting, subtract the point particle

add stars for microlensing

Parameters
alphasolar mass/Mpc
kappaconvergence
gammathree components of shear
phipotential solar masses
xcmposition relative to center (Mpc)
subtract_pointif true contribution from a point mass is subtracted
screeningthe factor by which to scale the mass for screening of the point mass subtraction

◆ get_Rmax()

float LensHalo::get_Rmax ( ) const
inline

this can be used to tag types of LensHalos

get the Rmax which is larger than Rsize in Mpc. This is the region exterior to which the halo will be considered a point mass. Between Rsize and Rmax the deflection and shear are interpolated.

◆ getDist()

PosType LensHalo::getDist ( ) const
inline

return current angular size distance, ie conversion between angular and special coordinates. This may not agree with the getZ() value because of the projection onto the lens plane.

◆ getEllipMethod()

EllipMethod LensHalo::getEllipMethod ( ) const
inline

stars

the method used to ellipticize a halo if fratio!=1 and halo is not NSIE

◆ kappa_asym()

PosType LensHalo::kappa_asym ( PosType  x,
PosType  theta 
)
protectedvirtual

radius in Mpc (Not x = r/rscale)

◆ map_variables()

PixelMap LensHalo::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.

Units : ALPHA - mass/PhysMpc - ALPHA / Sig_crit / Dl is the deflection in radians KAPPA - surface mass density , mass / /PhysMpc/PhysMpc - KAPPA / Sig_crit is the convergence GAMMA - mass / /PhysMpc/PhysMpc - GAMMA / Sig_crit is the shear PHI - mass - PHI / Sig_crit / Dl / Dl is the lensing potential whose angular gradient is the deflection and angular Laplacian is 2 times the convergence

centred on (0,0) in LensHalo coordinates

Parameters
lensvarlensing variable - KAPPA, ALPHA1, ALPHA2, GAMMA1, GAMMA2 or PHI
resangular resolution in radians

◆ MassBy2DIntegation()

PosType LensHalo::MassBy2DIntegation ( PosType  R)

Prints star parameters; if show_stars is true, prints data for single stars.

calculates the mass within radius R by integating kappa in theta and R, used only for testing

◆ operator=() [1/2]

LensHalo & LensHalo::operator= ( const LensHalo h)

Identification number of halo. It is not always used.

◆ operator=() [2/2]

LensHalo & LensHalo::operator= ( LensHalo &&  h)

Identification number of halo. It is not always used.

◆ setCosmology()

virtual void LensHalo::setCosmology ( const COSMOLOGY &  cosmo)
inlinevirtual

used for elliptical NFWs only, in that case get_switch_flag==true

set cosmology for halo

Reimplemented in LensHaloAnaNSIE.

◆ test()

bool LensHalo::test ( )

perform some basic consistancy checks for halo

Three tests: 1st - Mass via 1D integration vs mass via 2D integration. 2nd: gamma_t=alpha/r - kappa(R) which can be used for spherical distributions. Deviations are expected for axis ratios <1. For the latter case we use the next test. 3rd: The average along a circular aperture of gamma_t should be equal to <kappa(<R)> minus the average along a circular aperture over kappa. Note that also alpha/r - kappa is checked for consistency with kappa(<R)-<kappa(R)>. For axis ratios < 1 the factor between the two is expected to be of order O(10%).


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