GammaLib  2.1.0.dev
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GCTAEdisp2D Class Reference

CTA 2D energy dispersion class. More...

#include <GCTAEdisp2D.hpp>

Inheritance diagram for GCTAEdisp2D:
GCTAEdisp GBase

Classes

class  edisp_ereco_kern
 

Public Member Functions

 GCTAEdisp2D (void)
 Void constructor. More...
 
 GCTAEdisp2D (const GFilename &filename)
 File constructor. More...
 
 GCTAEdisp2D (const GCTAEdisp2D &edisp)
 Copy constructor. More...
 
virtual ~GCTAEdisp2D (void)
 Destructor. More...
 
GCTAEdisp2Doperator= (const GCTAEdisp2D &edisp)
 Assignment operator. More...
 
double operator() (const GEnergy &ereco, const GEnergy &etrue, const double &theta=0.0, const double &phi=0.0, const double &zenith=0.0, const double &azimuth=0.0) const
 Return energy dispersion in units of MeV \(^{-1}\). More...
 
void clear (void)
 Clear energy dispersion. More...
 
GCTAEdisp2Dclone (void) const
 Clone energy dispersion. More...
 
std::string classname (void) const
 Return class name. More...
 
void load (const GFilename &filename)
 Load energy dispersion from FITS file. More...
 
GFilename filename (void) const
 Return filename. More...
 
GEnergy mc (GRan &ran, const GEnergy &etrue, const double &theta=0.0, const double &phi=0.0, const double &zenith=0.0, const double &azimuth=0.0) const
 Simulate energy dispersion. More...
 
GEbounds ereco_bounds (const GEnergy &etrue, const double &theta=0.0, const double &phi=0.0, const double &zenith=0.0, const double &azimuth=0.0) const
 Return observed energy interval that contains the energy dispersion. More...
 
GEbounds etrue_bounds (const GEnergy &ereco, const double &theta=0.0, const double &phi=0.0, const double &zenith=0.0, const double &azimuth=0.0) const
 Return true energy interval that contains the energy dispersion. More...
 
double prob_erecobin (const GEnergy &ereco_min, const GEnergy &ereco_max, const GEnergy &etrue, const double &theta) const
 Return energy dispersion probability for reconstructed energy interval. More...
 
std::string print (const GChatter &chatter=NORMAL) const
 Print energy dispersion information. More...
 
void fetch (void) const
 Fetch energy dispersion. More...
 
const GCTAResponseTabletable (void) const
 Return response table. More...
 
void table (const GCTAResponseTable &table)
 Assign response table. More...
 
void read (const GFitsTable &table)
 Read energy dispersion from FITS table. More...
 
void write (GFitsBinTable &table) const
 Write energy dispersion into FITS binary table. More...
 
void save (const GFilename &filename, const bool &clobber=false) const
 Save energy dispersion table into FITS file. More...
 
- Public Member Functions inherited from GCTAEdisp
 GCTAEdisp (void)
 Void constructor. More...
 
 GCTAEdisp (const GCTAEdisp &edisp)
 Copy constructor. More...
 
virtual ~GCTAEdisp (void)
 Destructor. More...
 
GCTAEdispoperator= (const GCTAEdisp &edisp)
 Assignment operator. More...
 
- Public Member Functions inherited from GBase
virtual ~GBase (void)
 Destructor. More...
 

Private Member Functions

void init_members (void)
 Initialise class members. More...
 
void copy_members (const GCTAEdisp2D &edisp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
GEnergy etrue (const int &ietrue) const
 Get true energy. More...
 
double migra (const int &imigra) const
 Get migration. More...
 
double theta (const int &itheta) const
 Get offset angle in radiaus. More...
 
void compute_ereco_bounds (const double &theta=0.0, const double &phi=0.0, const double &zenith=0.0, const double &azimuth=0.0) const
 Compute vector of reconstructed energy bounds. More...
 
void compute_etrue_bounds (const double &theta=0.0, const double &phi=0.0, const double &zenith=0.0, const double &azimuth=0.0) const
 Compute vector of true energy boundaries. More...
 
void set_table (void)
 Set table. More...
 
void set_boundaries (void)
 Set energy dispersion boundaries. More...
 
void set_max_edisp (void)
 Set array of maximum energy dispersion values. More...
 
double get_max_edisp (const GEnergy &etrue, const double &theta) const
 Get maximum energy dispersion value. More...
 
void normalize_table (void)
 Normalize energy dispersion table. More...
 
int table_index (const int &ietrue, const int &imigra, const int &itheta) const
 Return index of response table element. More...
 
int table_stride (const int &axis) const
 Return stride of response table axis. More...
 
double table_value (const int &base_ll, const int &base_lr, const int &base_rl, const int &base_rr, const double &wgt_el, const double &wgt_er, const double &wgt_tl, const double &wgt_tr, const int &offset) const
 Return bi-linearly interpolate table value for given migration bin. More...
 
double table_value (const int &base_ll, const int &base_lr, const int &base_rl, const int &base_rr, const double &wgt_el, const double &wgt_er, const double &wgt_tl, const double &wgt_tr, const double &migra) const
 Return bi-linearly interpolate table value for given migration value. More...
 
void smooth_table (void)
 Smoothed energy dispersion table. More...
 
GNdarray smooth_array (const GNdarray &array, const int &nstart, const int &nstop, const double &limit) const
 Smoothed array. More...
 
double smoothed_array_value (const int &inx, const GNdarray &array, const int &nodes, const double &limit) const
 Get smoothed array value. More...
 
void get_moments (const int &itheta, GNdarray *mean, GNdarray *rms, GNdarray *total) const
 Compute moments. More...
 
void get_mean_rms (const GNdarray &array, double *mean, double *rms) const
 Compute mean and root mean square of migration array. More...
 
GNdarray gaussian_array (const double &mean, const double &rms, const double &total) const
 Return Gaussian approximation of energy dispersion array. More...
 

Private Attributes

GFilename m_filename
 Name of Edisp response file. More...
 
GCTAResponseTable m_edisp
 Edisp response table. More...
 
bool m_fetched
 Signals that Edisp has been fetched. More...
 
int m_inx_etrue
 True energy index. More...
 
int m_inx_migra
 Migration index. More...
 
int m_inx_theta
 Theta index. More...
 
int m_inx_matrix
 Matrix. More...
 
double m_logEsrc_min
 Minimum logE (log10(E/TeV)) More...
 
double m_logEsrc_max
 Maximum logE (log10(E/TeV)) More...
 
double m_migra_min
 Minimum migration. More...
 
double m_migra_max
 Maximum migration. More...
 
double m_theta_min
 Minimum theta (radians) More...
 
double m_theta_max
 Maximum theta (radians) More...
 
std::vector< double > m_max_edisp
 Maximum Edisp. More...
 
bool m_ereco_bounds_computed
 
bool m_etrue_bounds_computed
 
double m_last_theta_ereco
 
double m_last_theta_etrue
 
GEnergy m_last_etrue
 
GEnergy m_last_ereco
 
int m_index_ereco
 
int m_index_etrue
 
std::vector< GEboundsm_ereco_bounds
 
std::vector< GEboundsm_etrue_bounds
 

Additional Inherited Members

- Protected Member Functions inherited from GCTAEdisp
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GCTAEdisp &edisp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 

Detailed Description

CTA 2D energy dispersion class.

This class implements the energy dispersion for the CTA 2D response. The CTA 2D energy dispersion is in fact a 3-dimensional function, where the energy dispersion

\[ E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta) \]

is given as function of true energy \(E_{\rm true}\), reconstructed energy \(E_{\rm reco}\) and offset angle \(\theta\). The GCTAEdisp2D::operator() returns the value of the function \(E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta)\) in units of MeV \(^{-1}\) with

\[ \int_{E_{\rm reco}^{\rm min}}^{E_{\rm reco}^{\rm max}} E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta) \, dE_{\rm reco} = 1 \]

The energy dispersion is stored internally as a 3-dimensional response table, spanned by true energy \(E_{\rm true}\), migration \(m=E_{\rm reco}/E_{\rm true}\) and offset angle \(\theta\), with

\[ E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta) = \frac{E_{\rm disp}(m | E_{\rm true}, \theta)}{E_{\rm true}} \]

and

\[ \int_{m_{\rm min}}^{m_{\rm max}} E_{\rm disp}(m | E_{\rm true}, \theta) \, dm = 1 \]

Definition at line 90 of file GCTAEdisp2D.hpp.

Constructor & Destructor Documentation

GCTAEdisp2D::GCTAEdisp2D ( void  )

Void constructor.

Constructs empty energy dispersion.

Definition at line 76 of file GCTAEdisp2D.cpp.

References init_members().

Referenced by clone().

GCTAEdisp2D::GCTAEdisp2D ( const GFilename filename)
explicit

File constructor.

Parameters
[in]filenameFITS file name.

Constructs energy dispersion from a FITS file.

Definition at line 93 of file GCTAEdisp2D.cpp.

References init_members(), and load().

GCTAEdisp2D::GCTAEdisp2D ( const GCTAEdisp2D edisp)

Copy constructor.

Parameters
[in]edispEnergy dispersion.

Constructs energy dispersion by copying from another energy dispersion.

Definition at line 114 of file GCTAEdisp2D.cpp.

References copy_members(), and init_members().

GCTAEdisp2D::~GCTAEdisp2D ( void  )
virtual

Destructor.

Destructs energy dispersion.

Definition at line 132 of file GCTAEdisp2D.cpp.

References free_members().

Member Function Documentation

std::string GCTAEdisp2D::classname ( void  ) const
inlinevirtual

Return class name.

Returns
String containing the class name ("GCTAEdisp2D").

Implements GCTAEdisp.

Definition at line 260 of file GCTAEdisp2D.hpp.

void GCTAEdisp2D::clear ( void  )
virtual

Clear energy dispersion.

Clears energy dispersion.

Implements GCTAEdisp.

Definition at line 270 of file GCTAEdisp2D.cpp.

References GCTAEdisp::free_members(), free_members(), GCTAEdisp::init_members(), and init_members().

Referenced by load().

GCTAEdisp2D * GCTAEdisp2D::clone ( void  ) const
virtual

Clone energy dispersion.

Returns
Deep copy of energy dispersion.

Returns a pointer to a deep copy of the point spread function.

Implements GCTAEdisp.

Definition at line 292 of file GCTAEdisp2D.cpp.

References GCTAEdisp2D().

void GCTAEdisp2D::compute_ereco_bounds ( const double &  theta = 0.0,
const double &  phi = 0.0,
const double &  zenith = 0.0,
const double &  azimuth = 0.0 
) const
private

Compute vector of reconstructed energy bounds.

Parameters
[in]thetaOffset angle (radians).
[in]phiAzimuth angle (radians).
[in]zenithZenith angle (radians).
[in]azimuthAzimuth angle (radians).

Computes for all true energies the energy boundaries of the reconstructed energies covered by valid migration matrix elements. Only matrix elements with values >= 1.0e-12 are considered as valid elements. In case that no matrix elements are found of a given true energy, the interval of observed energies will be set to 0,0.

Definition at line 1186 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), etrue(), GEnergy::log10TeV(), m_edisp, m_ereco_bounds, m_inx_etrue, m_inx_matrix, m_inx_migra, m_inx_theta, migra(), GEnergy::TeV(), and theta().

Referenced by ereco_bounds().

void GCTAEdisp2D::compute_etrue_bounds ( const double &  theta = 0.0,
const double &  phi = 0.0,
const double &  zenith = 0.0,
const double &  azimuth = 0.0 
) const
private

Compute vector of true energy boundaries.

Parameters
[in]thetaOffset angle (radians).
[in]phiAzimuth angle (radians).
[in]zenithZenith angle (radians).
[in]azimuthAzimuth angle (radians).

Computes for all observed energies the energy boundaries of the true energies covered by valid migration matrix elements. Only matrix elements with values >= 1.0e-12 are considered as valid elements. In case that no matrix elements are found of a given observed energy, the interval of true energies will be set to 0,0.

Definition at line 1294 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), GEnergy::clear(), etrue(), m_edisp, m_etrue_bounds, m_inx_etrue, m_inx_migra, migra(), and operator()().

Referenced by etrue_bounds().

GEbounds GCTAEdisp2D::ereco_bounds ( const GEnergy etrue,
const double &  theta = 0.0,
const double &  phi = 0.0,
const double &  zenith = 0.0,
const double &  azimuth = 0.0 
) const
virtual

Return observed energy interval that contains the energy dispersion.

Parameters
[in]etrueTrue photon energy.
[in]thetaOffset angle in camera system (radians).
[in]phiAzimuth angle in camera system (radians). Not used.
[in]zenithZenith angle in Earth system (radians). Not used.
[in]azimuthAzimuth angle in Earth system (radians). Not used.
Returns
Observed energy boundaries.

Returns the band of observed energies outside of which the energy dispersion becomes negligible for a given true energy etrue and offset angle theta.

Implements GCTAEdisp.

Definition at line 588 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_lo(), compute_ereco_bounds(), etrue(), fetch(), m_edisp, m_ereco_bounds, m_ereco_bounds_computed, m_index_ereco, m_inx_etrue, m_last_etrue, m_last_theta_ereco, GEnergy::TeV(), and theta().

Referenced by mc(), and normalize_table().

GEnergy GCTAEdisp2D::etrue ( const int &  ietrue) const
private

Get true energy.

Parameters
[in]ietrueTrue energy index.
Returns
True energy.

Definition at line 1122 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), m_edisp, m_inx_etrue, sqrt(), and GEnergy::TeV().

Referenced by compute_ereco_bounds(), compute_etrue_bounds(), ereco_bounds(), mc(), normalize_table(), prob_erecobin(), and set_max_edisp().

GEbounds GCTAEdisp2D::etrue_bounds ( const GEnergy ereco,
const double &  theta = 0.0,
const double &  phi = 0.0,
const double &  zenith = 0.0,
const double &  azimuth = 0.0 
) const
virtual

Return true energy interval that contains the energy dispersion.

Parameters
[in]erecoReconstructed event energy.
[in]thetaOffset angle in camera system (radians).
[in]phiAzimuth angle in camera system (radians). Not used.
[in]zenithZenith angle in Earth system (radians). Not used.
[in]azimuthAzimuth angle in Earth system (radians). Not used.
Returns
True energy boundaries.

Returns the band of true photon energies outside of which the energy dispersion becomes negligible for a given observed energy ereco and offset angle theta.

Implements GCTAEdisp.

Definition at line 657 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_lo(), compute_etrue_bounds(), fetch(), m_edisp, m_etrue_bounds, m_etrue_bounds_computed, m_index_etrue, m_inx_etrue, m_last_ereco, m_last_theta_etrue, GEnergy::TeV(), and theta().

void GCTAEdisp2D::fetch ( void  ) const

Fetch energy dispersion.

Exceptions
GException::file_errorFile not found. Unable to load energy dispersion.

Fetches the energy dispersion by reading it from a FITS file.

If the filename contains no extension name the method scans the HDUCLASS keywords of all extensions and loads the energy dispersion from the first extension for which HDUCLAS4=EDISP_2D.

Otherwise, the background will be loaded from the ENERGY DISPERSION extension.

This method does nothing if the energy dispersion is already loaded, if there is nothing to fetch, or if the m_filename member is empty.

The method is thread save. The method checks whether the file from which the energy dispersion should be loaded actually exists.

Definition at line 734 of file GCTAEdisp2D.cpp.

References GFits::close(), GFilename::exists(), GFilename::extname(), gammalib::extname_cta_edisp2d, G_FETCH, gammalib::gadf_hduclas4(), GFilename::is_empty(), m_fetched, m_filename, read(), GFits::table(), and table().

Referenced by ereco_bounds(), etrue_bounds(), mc(), operator()(), print(), prob_erecobin(), table(), and write().

GFilename GCTAEdisp2D::filename ( void  ) const
inlinevirtual

Return filename.

Returns
Name of FITS file from which energy dispersion was loaded.

Implements GCTAEdisp.

Definition at line 272 of file GCTAEdisp2D.hpp.

References m_filename.

Referenced by load().

void GCTAEdisp2D::free_members ( void  )
private

Delete class members.

Definition at line 1109 of file GCTAEdisp2D.cpp.

Referenced by clear(), operator=(), and ~GCTAEdisp2D().

GNdarray GCTAEdisp2D::gaussian_array ( const double &  mean,
const double &  rms,
const double &  total 
) const
private

Return Gaussian approximation of energy dispersion array.

Parameters
[in]meanGaussian mean.
[in]rmsGaussian rms.
[in]totalGaussian total.
Returns
Gaussian approximation of energy dispersion array.

Returns a Gaussian approximation of the energy dispersion array by computing the mean migration value and its root mean square and by using these values as the centre and the width of a Gaussian function. The Gaussian function is normalized so that the sum of the output array is unity.

Definition at line 2196 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_nodes(), exp(), m_edisp, m_inx_migra, and GNodeArray::size().

Referenced by smooth_table().

double GCTAEdisp2D::get_max_edisp ( const GEnergy etrue,
const double &  theta 
) const
private

Get maximum energy dispersion value.

Parameters
[in]etrueTrue photon energy.
[in]thetaOffset angle (radians).
Returns
Maximum energy dispersion value.

For a given true energy etrue and offset angle theta, returns the maximum energy dispersion value that may occur. This method makes use of the internal array pre-computed by set_max_edisp().

If set_max_edisp() was not called, the method returns 0.

Definition at line 1556 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_nodes(), GNodeArray::inx_left(), GNodeArray::inx_right(), GEnergy::log10TeV(), m_edisp, m_inx_etrue, m_inx_theta, m_max_edisp, and GNodeArray::set_value().

Referenced by mc().

void GCTAEdisp2D::get_mean_rms ( const GNdarray array,
double *  mean,
double *  rms 
) const
private

Compute mean and root mean square of migration array.

Parameters
[in]arrayEnergy dispersion array.
[out]meanPointer to mean migration value.
[out]rmsPointer to root mean square value.

Computes the mean and the root mean square of the migration array. If the migration array is empty the mean and root mean square values are set to zero.

The method does not check whether the mean and rms pointers are valid.

Definition at line 2147 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_nodes(), m_edisp, m_inx_migra, migra(), GNodeArray::size(), sqrt(), and sum().

Referenced by get_moments().

void GCTAEdisp2D::get_moments ( const int &  itheta,
GNdarray mean,
GNdarray rms,
GNdarray total 
) const
private

Compute moments.

Parameters
[in]ithetaOffset angle index.
[out]meanPointer to mean array.
[out]rmsPointer to rms array.
[out]totalPointer to total array.

Computes the mean and root mean square values as function of true energy for a given offset angle.

Definition at line 2092 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), get_mean_rms(), m_edisp, m_inx_etrue, m_inx_matrix, m_inx_migra, and sum().

Referenced by smooth_table().

void GCTAEdisp2D::load ( const GFilename filename)
virtual

Load energy dispersion from FITS file.

Parameters
[in]filenameFITS file name.

Loads the energy dispersion from a FITS file.

The method does not actually load the FITS file, which will only be loaded on request. Only the filename is stored. See the fetch() method for the actually loading of the energy dispersion.

Implements GCTAEdisp.

Definition at line 401 of file GCTAEdisp2D.cpp.

References clear(), filename(), and m_filename.

Referenced by GCTAEdisp2D().

GEnergy GCTAEdisp2D::mc ( GRan ran,
const GEnergy etrue,
const double &  theta = 0.0,
const double &  phi = 0.0,
const double &  zenith = 0.0,
const double &  azimuth = 0.0 
) const
virtual

Simulate energy dispersion.

Parameters
[in]ranRandom number generator.
[in]etrueTrue photon energy.
[in]thetaOffset angle in camera system (radians).
[in]phiAzimuth angle in camera system (radians). Not used.
[in]zenithZenith angle in Earth system (radians). Not used.
[in]azimuthAzimuth angle in Earth system (radians). Not used.
Returns
Reconstructed energy.
Exceptions
GException::invalid_return_valueNo energy dispersion information found for parameters or energy dispersion matrix is empty.

Draws reconstructed energy value given a true energy etrue and offset angle theta. If no energy dispersion information is available the method will return the true photon energy.

Implements GCTAEdisp.

Definition at line 483 of file GCTAEdisp2D.cpp.

References GEbounds::emax(), GEbounds::emin(), ereco_bounds(), etrue(), fetch(), G_MC, get_max_edisp(), operator()(), GEnergy::print(), gammalib::rad2deg, gammalib::str(), GEnergy::TeV(), and GRan::uniform().

double GCTAEdisp2D::migra ( const int &  imigra) const
private

Get migration.

Parameters
[in]imigraMigration index.
Returns
Migration.

Definition at line 1143 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), m_edisp, and m_inx_migra.

Referenced by compute_ereco_bounds(), compute_etrue_bounds(), get_mean_rms(), and operator()().

void GCTAEdisp2D::normalize_table ( void  )
private

Normalize energy dispersion table.

Normalize the energy dispersion table using

\[ \int_{E_{\rm reco}^{\rm min}}^{E_{\rm reco}^{\rm max}} E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta) \, dE_{\rm reco} = 1 \]

where \(E_{\rm reco}_{\rm min}\) and \(E_{\rm reco}_{\rm max}\) is the minimum and maximum migration value for a given true energy as returned by the ebounds_obs() method, and \(E_{\rm disp}(E_{\rm true}, E_{\rm reco}, \theta)\) is the energy dispersion returned by the operator(), given in units of \(MeV^{-1}\).

The normalisation is performed for each \(E_{\rm true}\) and \(\theta\) bin using a Romberg integration method. Since the normalisation may affect the energy boundaries \(E_{\rm reco}^{\rm min}\) and \(E_{\rm reco}^{\rm max}\), two normalisation passes are performed to stabilize the result.

Definition at line 1629 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), GEbounds::emax(), GEbounds::emin(), GIntegral::eps(), ereco_bounds(), etrue(), GEnergy::log10MeV(), m_edisp, m_inx_etrue, m_inx_matrix, m_inx_migra, m_inx_theta, GIntegral::romberg(), GEbounds::size(), sum(), table_index(), table_stride(), and theta().

Referenced by set_table().

double GCTAEdisp2D::operator() ( const GEnergy ereco,
const GEnergy etrue,
const double &  theta = 0.0,
const double &  phi = 0.0,
const double &  zenith = 0.0,
const double &  azimuth = 0.0 
) const
virtual

Return energy dispersion in units of MeV \(^{-1}\).

Parameters
[in]erecoReconstructed photon energy.
[in]etrueTrue photon energy.
[in]thetaOffset angle in camera system (radians).
[in]phiAzimuth angle in camera system (radians). Not used.
[in]zenithZenith angle in Earth system (radians). Not used.
[in]azimuthAzimuth angle in Earth system (radians). Not used.
Returns
Energy dispersion (MeV \(^{-1}\))

Returns the energy dispersion

\[ E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta) = \frac{E_{\rm disp}(m | E_{\rm true}, \theta)}{E_{\rm true}} \]

in units of MeV \(^{-1}\) where \(E_{\rm reco}\) is the reconstructed energy, \(E_{\rm true}\) is the true energy, and \(\theta\) is the offset angle.

Implements GCTAEdisp.

Definition at line 203 of file GCTAEdisp2D.cpp.

References fetch(), GEnergy::log10TeV(), m_edisp, m_inx_etrue, m_inx_matrix, m_inx_migra, m_inx_theta, m_logEsrc_max, m_logEsrc_min, m_migra_max, m_migra_min, m_theta_max, m_theta_min, GEnergy::MeV(), migra(), and theta().

Referenced by compute_etrue_bounds(), and mc().

GCTAEdisp2D & GCTAEdisp2D::operator= ( const GCTAEdisp2D edisp)

Assignment operator.

Parameters
[in]edispEnergy dispersion.
Returns
Energy dispersion.

Assigns energy dispersion.

Definition at line 156 of file GCTAEdisp2D.cpp.

References copy_members(), free_members(), init_members(), and GCTAEdisp::operator=().

std::string GCTAEdisp2D::print ( const GChatter chatter = NORMAL) const
virtual

Print energy dispersion information.

Parameters
[in]chatterChattiness.
Returns
String containing energy dispersion information.

Implements GCTAEdisp.

Definition at line 962 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axes(), GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), fetch(), m_edisp, m_filename, m_inx_etrue, m_inx_migra, m_inx_theta, gammalib::parformat(), SILENT, and gammalib::str().

double GCTAEdisp2D::prob_erecobin ( const GEnergy ereco_min,
const GEnergy ereco_max,
const GEnergy etrue,
const double &  theta 
) const
virtual

Return energy dispersion probability for reconstructed energy interval.

Parameters
[in]ereco_minMinimum of reconstructed energy interval.
[in]ereco_maxMaximum of reconstructed energy interval.
[in]etrueTrue energy.
[in]thetaOffset angle (radians).
Returns
Integrated energy dispersion probability.

Computes

\[ \int_{E_{\rm reco}^{\rm min}}^{E_{\rm reco}^{\rm max}} E_{\rm disp}(E_{\rm reco} | E_{\rm true}, \theta) \, dE_{\rm reco} \]

where \(E_{\rm reco}\) is the reconstructed energy, \(E_{\rm true}\) is the true energy and \(\theta\) is the offset angle.

The method takes into account that the energy dispersion data are stored in a matrix and uses the trapezoidal rule for integration of the tabulated data. This assures the fastest possible integration.

Implements GCTAEdisp.

Definition at line 834 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_nodes(), etrue(), fetch(), GNodeArray::inx_left(), GNodeArray::inx_right(), GEnergy::log10TeV(), m_edisp, m_inx_etrue, m_inx_migra, m_inx_theta, m_logEsrc_max, m_logEsrc_min, m_migra_max, m_migra_min, m_theta_max, m_theta_min, GNodeArray::set_value(), GNodeArray::size(), table_index(), table_stride(), table_value(), GNodeArray::wgt_left(), and GNodeArray::wgt_right().

void GCTAEdisp2D::read ( const GFitsTable table)

Read energy dispersion from FITS table.

Parameters
[in]tableFITS table.
Exceptions
GException::invalid_valueResponse table is not three-dimensional.

Reads the energy dispersion form the FITS table. The following column names are mandatory:

ENERG_LO - True energy lower bin boundaries (alternative name: ETRUE_LO)
ENERG_HI - True energy upper bin boundaries (alternative name: ETRUE_HI)
MIGRA_LO - Migration lower bin boundaries
MIGRA_HI - Migration upper bin boundaries
THETA_LO - Offset angle lower bin boundaries
THETA_HI - Offset angle upper bin boundaries
MATRIX   - Migration matrix

The data are stored in the m_edisp member. The energy axis will be set to \(log_{10} E_{\rm true}\), the offset angle axis to radians.

The method assures that the energy dispersion is properly normalised.

Definition at line 342 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axes(), GCTAResponseTable::clear(), G_READ, m_edisp, GCTAResponseTable::read(), set_table(), and gammalib::str().

Referenced by fetch().

void GCTAEdisp2D::save ( const GFilename filename,
const bool &  clobber = false 
) const

Save energy dispersion table into FITS file.

Parameters
[in]filenameFITS file name.
[in]clobberOverwrite existing file?

Saves energy dispersion into a FITS file. If a file with the given filename does not yet exist it will be created, otherwise the method opens the existing file. The method will create a (or replace an existing) energy dispersion extension. The extension name can be specified as part of the filename, or if no extension name is given, is assumed to be ENERGY DISPERSION.

An existing file will only be modified if the clobber flag is set to true.

Definition at line 430 of file GCTAEdisp2D.cpp.

References GFitsHDU::extname(), GFilename::extname(), gammalib::extname_cta_edisp2d, GFits::remove(), table(), GFilename::url(), and write().

Referenced by set_table().

void GCTAEdisp2D::set_boundaries ( void  )
private

Set energy dispersion boundaries.

Sets the data members m_logEsrc_min, m_logEsrc_max, m_migra_min, m_migra_max, m_theta_min and m_theta_max that define the validity range of the energy dispersion.

Definition at line 1468 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), gammalib::deg2rad, log10(), m_edisp, m_inx_etrue, m_inx_migra, m_inx_theta, m_logEsrc_max, m_logEsrc_min, m_migra_max, m_migra_min, m_theta_max, and m_theta_min.

Referenced by set_table().

void GCTAEdisp2D::set_max_edisp ( void  )
private

Set array of maximum energy dispersion values.

Sets an internal array of maximum energy dispersion values for each given true photon energy and offset angle in the response table. This array will be used for Monte Carlo simulations.

Definition at line 1499 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), etrue(), m_edisp, m_inx_etrue, m_inx_matrix, m_inx_migra, m_inx_theta, m_max_edisp, GEnergy::MeV(), table_index(), and table_stride().

Referenced by set_table().

void GCTAEdisp2D::set_table ( void  )
private

Set table.

After assigning or loading a response table, performs all necessary steps to set the table.

For CTA, the method applies a kludge that smoothes the energy dispersion table to get rid of numerical fluctuations.

Sets the data members m_inx_etrue m_inx_migra m_inx_theta m_inx_matrix m_logEsrc_min m_logEsrc_max m_migra_min m_migra_max m_theta_min m_theta_max m_max_edisp

Definition at line 1416 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis(), GCTAResponseTable::axis_log10(), GCTAResponseTable::axis_radians(), GCTAResponseTable::has_axis(), m_edisp, m_inx_etrue, m_inx_matrix, m_inx_migra, m_inx_theta, normalize_table(), save(), set_boundaries(), set_max_edisp(), smooth_table(), gammalib::strip_whitespace(), GCTAResponseTable::table(), and GCTAResponseTable::telescope().

Referenced by read(), and table().

GNdarray GCTAEdisp2D::smooth_array ( const GNdarray array,
const int &  nstart,
const int &  nstop,
const double &  limit 
) const
private

Smoothed array.

Parameters
[in]arrayArray.
[in]nstartNumber of nodes used for regression at first array value.
[in]nstopNumber of nodes used for regression at last array value.
[in]limitLimit for array element exclusion.
Returns
Smoothed array.

Returns a smoothed array that is computed by locally adjusting a linear regression law to the array elements.

Definition at line 1949 of file GCTAEdisp2D.cpp.

References GNdarray::size(), and smoothed_array_value().

Referenced by smooth_table().

void GCTAEdisp2D::smooth_table ( void  )
private

Smoothed energy dispersion table.

This method implements the kludge for CTA that reduces the statistical noise in the energy dispersion matrix.

Definition at line 1867 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), gaussian_array(), get_moments(), m_edisp, m_inx_etrue, m_inx_matrix, m_inx_migra, m_inx_theta, and smooth_array().

Referenced by set_table().

double GCTAEdisp2D::smoothed_array_value ( const int &  inx,
const GNdarray array,
const int &  nodes,
const double &  limit 
) const
private

Get smoothed array value.

Parameters
[in]inxIndex of array element.
[in]arrayArray.
[in]nodesNumber of nodes used for regression.
[in]limitLimit for array element exclusion.
Returns
Smoothed array value.

Returns a smoothed array value that is derived from adjusting a linear law using regression to nodes adjacent array elements. Array elements with values below limit are excluded from the regression.

Definition at line 1984 of file GCTAEdisp2D.cpp.

References GNdarray::size().

Referenced by smooth_array().

const GCTAResponseTable & GCTAEdisp2D::table ( void  ) const
inline

Return response table.

Returns
Reference to response table.

Definition at line 285 of file GCTAEdisp2D.hpp.

References fetch(), and m_edisp.

Referenced by fetch(), save(), and table().

void GCTAEdisp2D::table ( const GCTAResponseTable table)

Assign response table.

Parameters
[in]tableResponse table.

Assigns the response table for an energy dispersion.

Definition at line 305 of file GCTAEdisp2D.cpp.

References m_edisp, set_table(), and table().

int GCTAEdisp2D::table_index ( const int &  ietrue,
const int &  imigra,
const int &  itheta 
) const
private

Return index of response table element.

Parameters
[in]ietrueTrue energy index.
[in]imigraMigration index.
[in]ithetaOffset index.
Returns
Index of response table element.

Definition at line 1720 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), m_edisp, m_inx_etrue, m_inx_migra, and m_inx_theta.

Referenced by normalize_table(), prob_erecobin(), and set_max_edisp().

int GCTAEdisp2D::table_stride ( const int &  axis) const
private

Return stride of response table axis.

Parameters
[in]axisResponse table axis.
Returns
Stride of response table axis.

Definition at line 1745 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_bins(), and m_edisp.

Referenced by normalize_table(), prob_erecobin(), set_max_edisp(), and table_value().

double GCTAEdisp2D::table_value ( const int &  base_ll,
const int &  base_lr,
const int &  base_rl,
const int &  base_rr,
const double &  wgt_el,
const double &  wgt_er,
const double &  wgt_tl,
const double &  wgt_tr,
const int &  offset 
) const
private

Return bi-linearly interpolate table value for given migration bin.

Parameters
[in]base_llBase index for left true energy and left offset angle.
[in]base_lrBase index for left true energy and right offset angle.
[in]base_rlBase index for right true energy and left offset angle.
[in]base_rrBase index for right true energy and right offset angle.
[in]wgt_elWeighting for left true energy.
[in]wgt_erWeighting for right true energy.
[in]wgt_tlWeighting for left offset angle.
[in]wgt_trWeighting for right offset angle.
[in]offsetOffset of migration bin with respect to base indices.
Returns
Bi-linearly interpolate table value.

Definition at line 1777 of file GCTAEdisp2D.cpp.

References m_edisp, and m_inx_matrix.

Referenced by prob_erecobin(), and table_value().

double GCTAEdisp2D::table_value ( const int &  base_ll,
const int &  base_lr,
const int &  base_rl,
const int &  base_rr,
const double &  wgt_el,
const double &  wgt_er,
const double &  wgt_tl,
const double &  wgt_tr,
const double &  migra 
) const
private

Return bi-linearly interpolate table value for given migration value.

Parameters
[in]base_llBase index for left true energy and left offset angle.
[in]base_lrBase index for left true energy and right offset angle.
[in]base_rlBase index for right true energy and left offset angle.
[in]base_rrBase index for right true energy and right offset angle.
[in]wgt_elWeighting for left true energy.
[in]wgt_erWeighting for right true energy.
[in]wgt_tlWeighting for left offset angle.
[in]wgt_trWeighting for right offset angle.
[in]migraMigration value.
Returns
Bi-linearly interpolate table value.

Definition at line 1818 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_nodes(), GNodeArray::inx_left(), GNodeArray::inx_right(), m_edisp, m_inx_migra, GNodeArray::set_value(), table_stride(), table_value(), GNodeArray::wgt_left(), and GNodeArray::wgt_right().

double GCTAEdisp2D::theta ( const int &  itheta) const
private

Get offset angle in radiaus.

Parameters
[in]ithetaOffset angle index.
Returns
Offset angle (radians).

Definition at line 1160 of file GCTAEdisp2D.cpp.

References GCTAResponseTable::axis_hi(), GCTAResponseTable::axis_lo(), gammalib::deg2rad, m_edisp, and m_inx_theta.

Referenced by compute_ereco_bounds(), ereco_bounds(), etrue_bounds(), normalize_table(), and operator()().

void GCTAEdisp2D::write ( GFitsBinTable table) const

Write energy dispersion into FITS binary table.

Parameters
[in]tableFITS binary table.

Writes the energy dispersion into the FITS binary table.

Todo:
Add keywords.

Definition at line 377 of file GCTAEdisp2D.cpp.

References fetch(), m_edisp, and GCTAResponseTable::write().

Referenced by save().

Member Data Documentation

std::vector<GEbounds> GCTAEdisp2D::m_ereco_bounds
mutableprivate

Definition at line 249 of file GCTAEdisp2D.hpp.

Referenced by compute_ereco_bounds(), copy_members(), ereco_bounds(), and init_members().

bool GCTAEdisp2D::m_ereco_bounds_computed
mutableprivate

Definition at line 241 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), ereco_bounds(), and init_members().

std::vector<GEbounds> GCTAEdisp2D::m_etrue_bounds
mutableprivate

Definition at line 250 of file GCTAEdisp2D.hpp.

Referenced by compute_etrue_bounds(), copy_members(), etrue_bounds(), and init_members().

bool GCTAEdisp2D::m_etrue_bounds_computed
mutableprivate

Definition at line 242 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), etrue_bounds(), and init_members().

bool GCTAEdisp2D::m_fetched
mutableprivate

Signals that Edisp has been fetched.

Definition at line 227 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), fetch(), and init_members().

GFilename GCTAEdisp2D::m_filename
mutableprivate

Name of Edisp response file.

Definition at line 225 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), fetch(), filename(), init_members(), load(), and print().

int GCTAEdisp2D::m_index_ereco
mutableprivate

Definition at line 247 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), ereco_bounds(), and init_members().

int GCTAEdisp2D::m_index_etrue
mutableprivate

Definition at line 248 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), etrue_bounds(), and init_members().

int GCTAEdisp2D::m_inx_matrix
private
GEnergy GCTAEdisp2D::m_last_ereco
mutableprivate

Definition at line 246 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), etrue_bounds(), and init_members().

GEnergy GCTAEdisp2D::m_last_etrue
mutableprivate

Definition at line 245 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), ereco_bounds(), and init_members().

double GCTAEdisp2D::m_last_theta_ereco
mutableprivate

Definition at line 243 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), ereco_bounds(), and init_members().

double GCTAEdisp2D::m_last_theta_etrue
mutableprivate

Definition at line 244 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), etrue_bounds(), and init_members().

double GCTAEdisp2D::m_logEsrc_max
private

Maximum logE (log10(E/TeV))

Definition at line 233 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), init_members(), operator()(), prob_erecobin(), and set_boundaries().

double GCTAEdisp2D::m_logEsrc_min
private

Minimum logE (log10(E/TeV))

Definition at line 232 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), init_members(), operator()(), prob_erecobin(), and set_boundaries().

std::vector<double> GCTAEdisp2D::m_max_edisp
private

Maximum Edisp.

Definition at line 238 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), get_max_edisp(), init_members(), and set_max_edisp().

double GCTAEdisp2D::m_migra_max
private

Maximum migration.

Definition at line 235 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), init_members(), operator()(), prob_erecobin(), and set_boundaries().

double GCTAEdisp2D::m_migra_min
private

Minimum migration.

Definition at line 234 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), init_members(), operator()(), prob_erecobin(), and set_boundaries().

double GCTAEdisp2D::m_theta_max
private

Maximum theta (radians)

Definition at line 237 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), init_members(), operator()(), prob_erecobin(), and set_boundaries().

double GCTAEdisp2D::m_theta_min
private

Minimum theta (radians)

Definition at line 236 of file GCTAEdisp2D.hpp.

Referenced by copy_members(), init_members(), operator()(), prob_erecobin(), and set_boundaries().


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