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

CTA instrument response function class. More...

#include <GCTAResponseIrf.hpp>

Inheritance diagram for GCTAResponseIrf:
GCTAResponse GResponse GBase

Public Member Functions

 GCTAResponseIrf (void)
 Void constructor. More...
 
 GCTAResponseIrf (const GCTAResponseIrf &rsp)
 Copy constructor. More...
 
 GCTAResponseIrf (const GXmlElement &xml)
 XML constructor. More...
 
 GCTAResponseIrf (const std::string &rspname, const GCaldb &caldb)
 Response constructor. More...
 
virtual ~GCTAResponseIrf (void)
 Destructor. More...
 
virtual GCTAResponseIrfoperator= (const GCTAResponseIrf &rsp)
 Assignment operator. More...
 
virtual void clear (void)
 Clear instance. More...
 
virtual GCTAResponseIrfclone (void) const
 Clone instance. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual bool is_valid (void) const
 Signal if response is valid. More...
 
virtual bool use_edisp (void) const
 Signal if response uses energy dispersion. More...
 
virtual bool use_tdisp (void) const
 Signal if time dispersion will be used. More...
 
virtual bool apply_edisp (void) const
 Signal if energy dispersion should be applied. More...
 
virtual void apply_edisp (const bool &apply_edisp) const
 Signal if energy dispersion should be applied. More...
 
virtual double irf (const GEvent &event, const GPhoton &photon, const GObservation &obs) const
 Return value of instrument response function. More...
 
virtual double nroi (const GModelSky &model, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return integral of event probability for a given sky model over ROI. More...
 
virtual GEbounds ebounds (const GEnergy &obsEnergy) const
 Return true energy boundaries for a specific observed energy. More...
 
virtual void read (const GXmlElement &xml)
 Read response from XML element. More...
 
virtual void write (GXmlElement &xml) const
 Write response information into XML element. More...
 
virtual std::string print (const GChatter &chatter=NORMAL) const
 Print CTA response information. More...
 
GCTAEventAtommc (const double &area, const GPhoton &photon, const GObservation &obs, GRan &ran) const
 Simulate event from photon. More...
 
void caldb (const GCaldb &caldb)
 Set calibration database. More...
 
const GCaldbcaldb (void) const
 Return calibration database. More...
 
const std::string & rspname (void) const
 Return response name. More...
 
void load (const std::string &rspname)
 Load CTA response. More...
 
void load_aeff (const GFilename &filename)
 Load effective area. More...
 
void load_psf (const GFilename &filename)
 Load CTA PSF vector. More...
 
void load_edisp (const GFilename &filename)
 Load energy dispersion information. More...
 
void load_background (const GFilename &filename)
 Load background model. More...
 
void offset_sigma (const double &sigma)
 Set offset angle dependence (degrees) More...
 
double offset_sigma (void) const
 Return offset angle dependence (degrees) More...
 
const GCTAAeffaeff (void) const
 Return pointer to effective area response. More...
 
void aeff (GCTAAeff *aeff)
 Set pointer to effective area response. More...
 
const GCTAPsfpsf (void) const
 Return pointer to point spread function. More...
 
void psf (GCTAPsf *psf)
 Set pointer to point spread function. More...
 
const GCTAEdispedisp (void) const
 Return pointer to energy dispersion. More...
 
void edisp (GCTAEdisp *edisp)
 Set pointer to energy dispersion. More...
 
const GCTABackgroundbackground (void) const
 Return pointer to background model. More...
 
void background (GCTABackground *background)
 Set pointer to background model. More...
 
const double & lo_safe_thres (void) const
 Return low energy threshold from IRF. More...
 
const double & hi_safe_thres (void) const
 Return high energy threshold from IRF. More...
 
double aeff (const double &theta, const double &phi, const double &zenith, const double &azimuth, const double &srcLogEng) const
 Return effective area (in units of cm2) More...
 
double psf (const double &delta, const double &theta, const double &phi, const double &zenith, const double &azimuth, const double &srcLogEng) const
 Return point spread function (in units of sr^-1) More...
 
double psf_delta_max (const double &theta, const double &phi, const double &zenith, const double &azimuth, const double &srcLogEng) const
 Return maximum angular separation (in radians) More...
 
double edisp (const GEnergy &ereco, const GEnergy &etrue, const double &theta, const double &phi, const double &zenith, const double &azimuth) const
 Return energy dispersion (in units of MeV \(^{-1}\)) More...
 
double nroi (const GModelSky &model, const GEnergy &srcEng, const GTime &srcTime, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of sky model. More...
 
double nirf (const GPhoton &photon, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of Instrument Response Function. More...
 
double npsf (const GSkyDir &srcDir, const double &srcLogEng, const GTime &srcTime, const GCTAPointing &pnt, const GCTARoi &roi) const
 Return result of PSF integration over ROI. More...
 
- Public Member Functions inherited from GCTAResponse
 GCTAResponse (void)
 Void constructor. More...
 
 GCTAResponse (const GCTAResponse &rsp)
 Copy constructor. More...
 
virtual ~GCTAResponse (void)
 Destructor. More...
 
virtual GCTAResponseoperator= (const GCTAResponse &rsp)
 Assignment operator. More...
 
- Public Member Functions inherited from GResponse
 GResponse (void)
 Void constructor. More...
 
 GResponse (const GResponse &rsp)
 Copy constructor. More...
 
virtual ~GResponse (void)
 Destructor. More...
 
virtual GResponseoperator= (const GResponse &rsp)
 Assignment operator. More...
 
virtual double convolve (const GModelSky &model, const GEvent &event, const GObservation &obs, const bool &grad=true) const
 Convolve sky model with the instrument response. More...
 
virtual GVector convolve (const GModelSky &model, const GObservation &obs, GMatrixSparse *gradients=NULL) const
 Convolve sky model with the instrument response. More...
 
virtual double irf_spatial (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return instrument response integrated over the spatial model. More...
 
virtual GVector irf_spatial (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response vector integrated over the spatial model. More...
 
virtual void remove_response_cache (const std::string &name)
 Remove response cache for model. 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 GCTAResponseIrf &rsp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
GFilename irf_filename (const std::string &filename) const
 Return filename with appropriate extension. More...
 
double irf_ptsrc (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return value of point source instrument response function. More...
 
double irf_radial (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return IRF value for radial source model. More...
 
double irf_elliptical (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return Irf value for elliptical source model. More...
 
double irf_diffuse (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return value of diffuse source instrument response function. More...
 
double nroi_ptsrc (const GModelSky &model, const GEnergy &srcEng, const GTime &srcTime, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of point source model. More...
 
double nroi_radial (const GModelSky &model, const GEnergy &srcEng, const GTime &srcTime, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of radial source model. More...
 
double nroi_elliptical (const GModelSky &model, const GEnergy &srcEng, const GTime &srcTime, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of elliptical source model. More...
 
double nroi_diffuse (const GModelSky &model, const GEnergy &srcEng, const GTime &srcTime, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of diffuse source model. More...
 
double nroi_composite (const GModelSky &model, const GEnergy &srcEng, const GTime &srcTime, const GEnergy &obsEng, const GTime &obsTime, const GObservation &obs) const
 Return spatial integral of composite source model. More...
 

Private Attributes

GCaldb m_caldb
 Calibration database. More...
 
std::string m_rspname
 Name of the instrument response. More...
 
GCTAAeffm_aeff
 Effective area. More...
 
GCTAPsfm_psf
 Point spread function. More...
 
GCTAEdispm_edisp
 Energy dispersion. More...
 
GCTABackgroundm_background
 Background. More...
 
bool m_apply_edisp
 Apply energy dispersion. More...
 
double m_lo_safe_thres
 Safe low energy threshold. More...
 
double m_hi_safe_thres
 Safe high energy threshold. More...
 
std::string m_xml_caldb
 Calibration database string in XML file. More...
 
std::string m_xml_rspname
 Response name in XML file. More...
 

Additional Inherited Members

- Protected Member Functions inherited from GCTAResponse
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GCTAResponse &rsp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
- Protected Member Functions inherited from GResponse
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GResponse &rsp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
double eval_prob (const GModelSky &model, const GEvent &event, const GEnergy &srcEng, const GTime &srcTime, const GObservation &obs, const bool &grad) const
 Convolve sky model with the instrument response. More...
 
GVector eval_probs (const GModelSky &model, const GObservation &obs, GMatrixSparse *gradients=NULL) const
 Convolve sky model with the instrument response. More...
 
int size_edisp_vector (const GModelSky &model, const GObservation &obs, const bool &grad) const
 Return size of vector for energy dispersion computation. More...
 
GEbounds ebounds_model (const GModelSky &model) const
 Return true energy intervals for sky model. More...
 
virtual double irf_composite (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return instrument response to composite source. More...
 
virtual GVector irf_ptsrc (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to point source sky model. More...
 
virtual GVector irf_radial (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to radial source sky model. More...
 
virtual GVector irf_elliptical (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to ellipitical source sky model. More...
 
virtual GVector irf_diffuse (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to diffuse source sky model. More...
 
virtual GVector irf_composite (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to composite source sky model. More...
 
- Protected Attributes inherited from GResponse
bool m_use_irf_cache
 Control usage of irf cache. More...
 
bool m_use_nroi_cache
 Control usage of nroi cache. More...
 
int m_irf_radial_iter_theta
 Radial model integration theta iterations. More...
 
int m_irf_radial_iter_phi
 Radial model integration phi iterations. More...
 
int m_irf_elliptical_iter_theta
 Elliptical model integration theta iterations. More...
 
int m_irf_elliptical_iter_phi
 Elliptical model integration phi iterations. More...
 
double m_irf_diffuse_resolution
 Angular resolution for diffuse model. More...
 
GResponseCache m_irf_cache
 
GResponseCache m_nroi_cache
 
GResponseVectorCache m_irf_vector_cache
 

Detailed Description

CTA instrument response function class.

Definition at line 62 of file GCTAResponseIrf.hpp.

Constructor & Destructor Documentation

GCTAResponseIrf::GCTAResponseIrf ( void  )

Void constructor.

Constructs void CTA response.

Definition at line 140 of file GCTAResponseIrf.cpp.

References init_members().

Referenced by clone().

GCTAResponseIrf::GCTAResponseIrf ( const GCTAResponseIrf rsp)

Copy constructor.

Parameters
[in]rspCTA response.

Constructs CTA response by making a deep copy of an existing object.

Definition at line 157 of file GCTAResponseIrf.cpp.

References copy_members(), and init_members().

GCTAResponseIrf::GCTAResponseIrf ( const GXmlElement xml)
explicit

XML constructor.

Parameters
[in]xmlXML element.

Construct CTA response from XML element.

Definition at line 177 of file GCTAResponseIrf.cpp.

References init_members(), and read().

GCTAResponseIrf::GCTAResponseIrf ( const std::string &  rspname,
const GCaldb caldb 
)

Response constructor.

Parameters
[in]rspnameResponse file name.
[in]caldbCalibration database.

Create instance of CTA response by specifying the response name and the calibration database. The response name can be either a response identifier or a filename (see GCTAResponseIrf::load for more information).

Definition at line 200 of file GCTAResponseIrf.cpp.

References caldb(), init_members(), load(), and m_caldb.

GCTAResponseIrf::~GCTAResponseIrf ( void  )
virtual

Destructor.

Destroys instance of CTA response object.

Definition at line 222 of file GCTAResponseIrf.cpp.

References free_members().

Member Function Documentation

const GCTAAeff * GCTAResponseIrf::aeff ( void  ) const
inline
void GCTAResponseIrf::aeff ( GCTAAeff aeff)
inline

Set pointer to effective area response.

Parameters
[in]aeffPointer to effective area response.

Definition at line 361 of file GCTAResponseIrf.hpp.

References GCTAAeff::clone(), and m_aeff.

double GCTAResponseIrf::aeff ( const double &  theta,
const double &  phi,
const double &  zenith,
const double &  azimuth,
const double &  srcLogEng 
) const

Return effective area (in units of cm2)

Parameters
[in]thetaRadial offset angle of photon in camera (radians).
[in]phiPolar angle of photon in camera (radians).
[in]zenithZenith angle of telescope pointing (radians).
[in]azimuthAzimuth angle of telescope pointing (radians).
[in]srcLogEngLog10 of true photon energy (E/TeV).
Returns
Effective area in units fo cm2.
Exceptions
GException::invalid_valueNo effective area information found.

Returns the effective area as function of the true photon position in the camera system and the telescope pointing direction in the Earth system.

Definition at line 1736 of file GCTAResponseIrf.cpp.

References aeff(), G_AEFF, and m_aeff.

bool GCTAResponseIrf::apply_edisp ( void  ) const
inlinevirtual

Signal if energy dispersion should be applied.

Returns
True if energy dispersion should be applied

Implements GCTAResponse.

Definition at line 274 of file GCTAResponseIrf.hpp.

References m_apply_edisp.

Referenced by apply_edisp(), GCTAOnOffObservation::compute_arf(), and print().

void GCTAResponseIrf::apply_edisp ( const bool &  apply_edisp) const
inlinevirtual

Signal if energy dispersion should be applied.

Parameters
[in]apply_edispSet true if energy dispersion should be applied

Implements GCTAResponse.

Definition at line 285 of file GCTAResponseIrf.hpp.

References apply_edisp(), and m_apply_edisp.

const GCTABackground * GCTAResponseIrf::background ( void  ) const
inline

Return pointer to background model.

Returns
Pointer to background model.

Definition at line 426 of file GCTAResponseIrf.hpp.

References m_background.

Referenced by gammalib::cta_rsp_bkg(), read(), and write().

void GCTAResponseIrf::background ( GCTABackground background)
inline

Set pointer to background model.

Parameters
[in]backgroundPointer to background model.

Definition at line 438 of file GCTAResponseIrf.hpp.

References GCTABackground::clone(), and m_background.

void GCTAResponseIrf::caldb ( const GCaldb caldb)
inline

Set calibration database.

Parameters
[in]caldbCalibration database.

Sets the calibration database for the CTA response.

Definition at line 324 of file GCTAResponseIrf.hpp.

References caldb(), and m_caldb.

Referenced by GCTAObservation::response(), and GCTAOnOffObservation::set().

const GCaldb & GCTAResponseIrf::caldb ( void  ) const
inline

Return calibration database.

Returns
Calibration database.

Definition at line 310 of file GCTAResponseIrf.hpp.

References m_caldb.

Referenced by caldb(), GCTAResponseIrf(), load(), and read().

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

Return class name.

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

Implements GCTAResponse.

Definition at line 234 of file GCTAResponseIrf.hpp.

void GCTAResponseIrf::clear ( void  )
virtual

Clear instance.

Clears CTA response object by resetting all members to an initial state. Any information that was present in the object before will be lost.

Implements GCTAResponse.

Definition at line 285 of file GCTAResponseIrf.cpp.

References GCTAResponse::free_members(), GResponse::free_members(), free_members(), GCTAResponse::init_members(), GResponse::init_members(), and init_members().

Referenced by load().

GCTAResponseIrf * GCTAResponseIrf::clone ( void  ) const
virtual

Clone instance.

Returns
Pointer to deep copy of CTA response.

Creates a clone (deep copy) of a CTA response object.

Implements GCTAResponse.

Definition at line 309 of file GCTAResponseIrf.cpp.

References GCTAResponseIrf().

void GCTAResponseIrf::copy_members ( const GCTAResponseIrf rsp)
private
GEbounds GCTAResponseIrf::ebounds ( const GEnergy obsEnergy) const
virtual

Return true energy boundaries for a specific observed energy.

Parameters
[in]obsEnergyObserved Energy.
Returns
True energy boundaries for given observed energy.

Implements GCTAResponse.

Definition at line 525 of file GCTAResponseIrf.cpp.

References edisp(), and GCTAEdisp::etrue_bounds().

Referenced by nroi().

const GCTAEdisp * GCTAResponseIrf::edisp ( void  ) const
inline

Return pointer to energy dispersion.

Returns
Pointer to energy dispersion.

Definition at line 399 of file GCTAResponseIrf.hpp.

References m_edisp.

Referenced by GCTAOnOffObservation::compute_rmf(), ebounds(), edisp(), cta_irf_radial_kern_omega::eval(), cta_irf_elliptical_kern_omega::eval(), cta_irf_diffuse_kern_phi::eval(), irf(), mc(), nirf(), nroi(), and write().

void GCTAResponseIrf::edisp ( GCTAEdisp edisp)
inline

Set pointer to energy dispersion.

Parameters
[in]edispPointer to energy dispersion.

Definition at line 411 of file GCTAResponseIrf.hpp.

References GCTAEdisp::clone(), and m_edisp.

double GCTAResponseIrf::edisp ( const GEnergy ereco,
const GEnergy etrue,
const double &  theta,
const double &  phi,
const double &  zenith,
const double &  azimuth 
) const

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

Parameters
[in]erecoReconstructed event energy.
[in]etrueTrue photon energy.
[in]thetaRadial offset angle in camera (radians).
[in]phiPolar angle in camera (radians).
[in]zenithZenith angle of telescope pointing (radians).
[in]azimuthAzimuth angle of telescope pointing (radians).
Returns
Energy dispersion (MeV \(^{-1}\)).

Definition at line 1851 of file GCTAResponseIrf.cpp.

References edisp().

void GCTAResponseIrf::free_members ( void  )
private

Delete class members.

Definition at line 2215 of file GCTAResponseIrf.cpp.

References m_aeff, m_background, m_edisp, and m_psf.

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

const double & GCTAResponseIrf::hi_safe_thres ( void  ) const
inline

Return high energy threshold from IRF.

Returns
High energy threshold from IRF.

Definition at line 461 of file GCTAResponseIrf.hpp.

References m_hi_safe_thres.

void GCTAResponseIrf::init_members ( void  )
private

Initialise class members.

< Switched off by default

Definition at line 2159 of file GCTAResponseIrf.cpp.

References GCaldb::clear(), m_aeff, m_apply_edisp, m_background, m_caldb, m_edisp, m_hi_safe_thres, m_lo_safe_thres, m_psf, m_rspname, m_xml_caldb, and m_xml_rspname.

Referenced by clear(), GCTAResponseIrf(), and operator=().

double GCTAResponseIrf::irf ( const GEvent event,
const GPhoton photon,
const GObservation obs 
) const
virtual

Return value of instrument response function.

Parameters
[in]eventObserved event.
[in]photonIncident photon.
[in]obsObservation.
Todo:
Set polar angle phi of photon in camera system

Implements GCTAResponse.

Definition at line 324 of file GCTAResponseIrf.cpp.

References aeff(), GCTAPointing::azimuth(), gammalib::cta_dir(), gammalib::cta_pnt(), GObservation::deadc(), GPhoton::dir(), GCTAPointing::dir(), GCTAInstDir::dir(), GSkyDir::dist(), edisp(), GPhoton::energy(), G_IRF, gammalib::is_infinite(), gammalib::is_notanumber(), GEnergy::log10TeV(), psf(), psf_delta_max(), GPhoton::time(), use_edisp(), and GCTAPointing::zenith().

Referenced by irf_diffuse(), irf_elliptical(), irf_ptsrc(), and irf_radial().

double GCTAResponseIrf::irf_diffuse ( const GEvent event,
const GSource source,
const GObservation obs 
) const
privatevirtual

Return value of diffuse source instrument response function.

Parameters
[in]eventObserved event.
[in]sourceSource.
[in]obsObservation.

Integrates the product of the model and the IRF over the true photon arrival direction using

\[ \int_{0}^{\theta_{\rm max}} \sin \theta \times PSF(\theta) \int_{0}^{2\pi} S_{\rm p}(\theta, \phi | E, t) \, Aeff(\theta, \phi) \, Edisp(\theta, \phi) d\phi \]

where

  • \(S_{\rm p}(\theta, \phi | E, t)\) is the diffuse model,
  • \(PSF(\theta)\) is the azimuthally symmetric Point Spread Function,
  • \(Aeff(\theta, \phi)\) is the effective area,
  • \(Edisp(\theta, \phi)\) is the energy dispersion,
  • \(\theta\) is the distance from the PSF centre, and
  • \(\phi\) is the azimuth angle.

The integration is performed in the reference of the observed arrival direction. Integration is done first over the azimuth angle \(\phi\) and then over the offset angle \(\theta\).

The integration kernels for this method are implemented by the response helper classes cta_irf_diffuse_kern_theta and cta_irf_diffuse_kern_phi.

Reimplemented from GResponse.

Definition at line 2764 of file GCTAResponseIrf.cpp.

References GCTAPointing::azimuth(), gammalib::cta_dir(), gammalib::cta_pnt(), GObservation::deadc(), GCTAPointing::dir(), GCTAInstDir::dir(), GSkyDir::dist(), GSource::energy(), GMatrix::eulery(), GMatrix::eulerz(), GIntegral::fixed_iter(), G_IRF_DIFFUSE, irf(), gammalib::is_infinite(), gammalib::is_notanumber(), gammalib::iter_rho(), GEnergy::log10TeV(), GSource::model(), psf_delta_max(), gammalib::resolution(), GIntegral::romberg(), GSource::time(), and GCTAPointing::zenith().

double GCTAResponseIrf::irf_elliptical ( const GEvent event,
const GSource source,
const GObservation obs 
) const
privatevirtual

Return Irf value for elliptical source model.

Parameters
[in]eventObserved event.
[in]sourceSource.
[in]obsObservation.
Exceptions
GException::invalid_argumentModel is not an elliptical model.

Integrates the product of the model and the IRF over the true photon arrival direction using

\[ \int_{\rho_{\rm min}}^{\rho_{\rm max}} \sin \rho \times \int_{\omega} S_{\rm p}(\rho, \omega | E, t) \, IRF(\rho, \omega) d\omega d\rho \]

where \(S_{\rm p}(\rho, \omega | E, t)\) is the elliptical model, \(IRF(\rho, \omega)\) is the instrument response function, \(\rho\) is the distance from the model centre, and \(\omega\) is the position angle with respect to the connecting line between the model centre and the observed photon arrival direction.

The source model centre is located at \(\vec{m}\), and a spherical coordinate system is defined around this location with \((\rho,\omega)\) being the zenith and azimuth angles, respectively. The azimuth angle \((\omega)\) is counted counterclockwise from the vector that runs from the model centre \(\vec{m}\) to the measured photon direction \(\vec{p'}\).

Reimplemented from GResponse.

Definition at line 2564 of file GCTAResponseIrf.cpp.

References GCTAPointing::azimuth(), gammalib::centre(), gammalib::cta_dir(), gammalib::cta_pnt(), GObservation::deadc(), gammalib::deg2rad, GCTAPointing::dir(), GCTAInstDir::dir(), GModelSpatialElliptical::dir(), GSkyDir::dist(), GSource::energy(), GIntegral::fixed_iter(), G_IRF_ELLIPTICAL, irf(), gammalib::is_infinite(), gammalib::is_notanumber(), gammalib::iter_phi(), gammalib::iter_rho(), GEnergy::log10TeV(), GSource::model(), gammalib::pihalf, GSkyDir::posang(), GModelSpatialElliptical::posangle(), psf_delta_max(), GIntegral::romberg(), GModelSpatialElliptical::semimajor(), GModelSpatialElliptical::semiminor(), GModelSpatialElliptical::theta_max(), GSource::time(), and GCTAPointing::zenith().

GFilename GCTAResponseIrf::irf_filename ( const std::string &  filename) const
private

Return filename with appropriate extension.

Parameters
[in]filenameFile name.
Returns
File name.

Checks if the specified filename exists, and if not, checks whether a file with the added suffix .dat exists. Returns the file name with the appropriate extension.

Definition at line 2244 of file GCTAResponseIrf.cpp.

References GFilename::exists().

Referenced by load().

double GCTAResponseIrf::irf_ptsrc ( const GEvent event,
const GSource source,
const GObservation obs 
) const
privatevirtual

Return value of point source instrument response function.

Parameters
[in]eventObserved event.
[in]sourceSource.
[in]obsObservation.
Returns
Value of instrument response function for a point source.

This method returns the value of the instrument response function for a point source. The method assumes that source.model() is of type GModelSpatialPointSource.

Reimplemented from GResponse.

Definition at line 2274 of file GCTAResponseIrf.cpp.

References GModelSpatialPointSource::dir(), GSource::energy(), irf(), GSource::model(), and GSource::time().

double GCTAResponseIrf::irf_radial ( const GEvent event,
const GSource source,
const GObservation obs 
) const
privatevirtual

Return IRF value for radial source model.

Parameters
[in]eventObserved event.
[in]sourceSource.
[in]obsObservation.
Exceptions
GException::invalid_argumentModel is not a radial model.

Integrates the product of the spatial model and the instrument response function over the true photon arrival direction using

\[ \int_{\rho_{\rm min}}^{\rho_{\rm max}} \sin \rho \times S_{\rm p}(\rho | E, t) \times \int_{\omega_{\rm min}}^{\omega_{\rm max}} {\rm Irf}(\rho, \omega) d\omega d\rho \]

where \(S_{\rm p}(\rho | E, t)\) is the radial spatial model, \({\rm Irf}(\rho, \omega)\) is the instrument response function, \(\rho\) is the radial distance from the model centre, and \(\omega\) is the position angle with respect to the connecting line between the model centre and the observed photon arrival direction.

The integration is performed in the coordinate system of the source model spanned by \(\rho\) and \(\omega\) which allows to benefit from the symmetry of the source model.

The source centre is located at \(\vec{m}\), and a spherical system is defined around this location with \((\omega,\rho)\) being the azimuth and zenith angles, respectively. \(\omega=0\) is defined by the direction that connects the source centre \(\vec{m}\) to the measured photon direction \(\vec{p'}\), and \(\omega\) increases counterclockwise.

Note that this method approximates the true theta angle (angle between incident photon and pointing direction) by the measured theta angle (angle between the measured photon arrival direction and the pointing direction). Given the slow variation of the Psf shape over the field of view, this approximation should be fine. It helps in fact a lot in speeding up the computations.

Reimplemented from GResponse.

Definition at line 2338 of file GCTAResponseIrf.cpp.

References gammalib::acos(), GCTAPointing::azimuth(), gammalib::centre(), cos(), gammalib::cta_dir(), gammalib::cta_pnt(), GObservation::deadc(), gammalib::deg2rad, GCTAPointing::dir(), GCTAInstDir::dir(), GModelSpatialRadial::dir(), GSkyDir::dist(), GSource::energy(), GIntegral::fixed_iter(), G_IRF_RADIAL, irf(), gammalib::is_infinite(), gammalib::is_notanumber(), gammalib::iter_phi(), gammalib::iter_rho(), GEnergy::log10TeV(), GSource::model(), psf_delta_max(), gammalib::rad2deg, GModelSpatialRadialRing::radius(), GModelSpatialRadialShell::radius(), GIntegral::romberg(), sin(), GModelSpatialRadial::theta_max(), GSource::time(), and GCTAPointing::zenith().

bool GCTAResponseIrf::is_valid ( void  ) const
inlinevirtual

Signal if response is valid.

Returns
True if response is valid

Implements GCTAResponse.

Definition at line 246 of file GCTAResponseIrf.hpp.

References m_aeff, and m_psf.

const double & GCTAResponseIrf::lo_safe_thres ( void  ) const
inline

Return low energy threshold from IRF.

Returns
Low energy threshold from IRF.

Definition at line 450 of file GCTAResponseIrf.hpp.

References m_lo_safe_thres.

void GCTAResponseIrf::load ( const std::string &  rspname)

Load CTA response.

Parameters
[in]rspnameCTA response name.

Loads the CTA response with specified name rspname. The method first searchs for an appropriate response in the calibration database. If no appropriate response is found, the method takes the database root path and response name to build the full path to the response file, and tries to load the response from these paths.

The method sets the calibration database and response names for writing of the response information into the XML file.

Definition at line 1002 of file GCTAResponseIrf.cpp.

References caldb(), clear(), GCaldb::filename(), gammalib::filepath(), GCaldb::instrument(), irf_filename(), GFilename::is_empty(), load_aeff(), load_background(), load_edisp(), load_psf(), m_aeff, m_caldb, m_rspname, m_xml_caldb, m_xml_rspname, GCTAAeffArf::remove_thetacut(), GCaldb::rootdir(), and rspname().

Referenced by GCTAResponseIrf(), read(), and GCTAObservation::response().

void GCTAResponseIrf::load_aeff ( const GFilename filename)

Load effective area.

Parameters
[in]filenameEffective area filename.
Exceptions
GException::file_errorFile or extension not found.

Loads the effective area from a response file.

If the file is a FITS file, the method will either use the extension name specified with the filename, or if no extension name is given, search for an EFFECTIVE AREA or SPECRESP extension in the file and open the corresponding FITS table. If a column named SPECRESP is found in the table, a GCTAAeffArf object will be allocated. Otherwise, a GCTAAeff2D object will be allocated. In both cases, the method will extract the optional LO_THRES and HI_THRES safe energy thresholds from the FITS file header.

If the file is not a FITS file, it will be interpreted as a GCTAAeffPerfTable performance table.

Definition at line 1100 of file GCTAResponseIrf.cpp.

References GFits::close(), GFitsTable::contains(), GFits::contains(), GFilename::exists(), GFilename::extname(), gammalib::extname_arf, gammalib::extname_cta_aeff2d, G_LOAD_AEFF, gammalib::gadf_hduclas4(), GFitsHDU::has_card(), GFilename::has_extname(), GFilename::is_fits(), m_aeff, m_hi_safe_thres, m_lo_safe_thres, GFitsHDU::real(), and GFits::table().

Referenced by load(), and read().

void GCTAResponseIrf::load_background ( const GFilename filename)

Load background model.

Parameters
[in]filenameBackground model file name.
Exceptions
GException::file_errorFile not found.

Loads the background model from a response file.

If the file is a FITS file the method will check whether the file contains a BKG_2D or BKG_3D table, and load correspondingly a GCTABackground2D or GCTABackground3D response.

If the file is not a FITS file it will be interpreted as a GCTABackgroundPerfTable performance table.

Definition at line 1472 of file GCTAResponseIrf.cpp.

References GFits::close(), GFitsTable::contains(), GFits::contains(), GFilename::exists(), GFilename::extname(), gammalib::extname_cta_background2d, gammalib::extname_cta_background3d, G_LOAD_BACKGROUND, GFilename::has_extname(), GFilename::is_fits(), m_background, and GFits::table().

Referenced by load(), and read().

void GCTAResponseIrf::load_edisp ( const GFilename filename)

Load energy dispersion information.

Parameters
[in]filenameEnergy dispersion file name.
Exceptions
GException::file_errorFile or extension not found.

Loads the energy dispersion from a response file.

If the file is a FITS file, the method will either use the extension name specified with the filename, or if no extension name is given, search for an ENERGY DISPERSION or MATRIX extension in the file and open the corresponding FITS table. If columns named "MIGRA_LO" and "MIGRA_HI" are found in the table, a GCTAEdisp2D object will be allocated. Otherwise, a GCTAEdispRmf object will be allocated.

If the file is not a FITS file, it will be interpreted as a GCTAEdispPerfTable performance table.

Definition at line 1365 of file GCTAResponseIrf.cpp.

References GFits::close(), GFitsTable::contains(), GFits::contains(), GFilename::exists(), GFilename::extname(), gammalib::extname_cta_edisp2d, gammalib::extname_rmf, G_LOAD_EDISP, gammalib::gadf_hduclas4(), GFilename::has_extname(), GFilename::is_fits(), m_edisp, and GFits::table().

Referenced by load(), and read().

void GCTAResponseIrf::load_psf ( const GFilename filename)

Load CTA PSF vector.

Parameters
[in]filenameFITS file name.
Exceptions
GException::file_errorFile or extension not found.

Loads the point spead function from a response file.

If the file is a FITS file, the method will either use the extension name specified with the filename, or if no extension name is given, search for one of the following extension names

  POINT SPREAD FUNCTION
  PSF
  PSF_2D_TABLE

in the file and open the corresponding FITS table. If columns named GAMMA and SIGMA are found in the table, a GCTAPsfKing object will be allocated.

If columns named SCALE, SIGMA_1, AMPL_2, SIGMA_2, AMPL_3 and SIGMA_3 are found in the table, a GCTAPsf2D object will be allocated.

If columns named RAD_LO, RAD_HI, and RPSF are found in the table, a GCTAPsfTable object will be allocated.

Otherwise, a CTAPsfVector object will be allocated.

If the file is not a FITS file, it will be interpreted as a GCTAPsfPerfTable performance table.

Definition at line 1229 of file GCTAResponseIrf.cpp.

References GFits::close(), GFitsTable::contains(), GFits::contains(), GFilename::exists(), GFilename::extname(), G_LOAD_PSF, gammalib::gadf_hduclas4(), GFilename::has_extname(), GFilename::is_fits(), m_psf, and GFits::table().

Referenced by load(), and read().

GCTAEventAtom * GCTAResponseIrf::mc ( const double &  area,
const GPhoton photon,
const GObservation obs,
GRan ran 
) const

Simulate event from photon.

Parameters
[in]areaSimulation surface area.
[in]photonPhoton.
[in]obsObservation.
[in]ranRandom number generator.
Returns
Simulated event.

Simulates a CTA event using the response function from an incident photon. If the event is not detected a NULL pointer is returned.

The method also applies a deadtime correction using a Monte Carlo process, taking into account temporal deadtime variations. For this purpose, the method makes use of the time dependent GObservation::deadc method.

Todo:

Set polar angle phi of photon in camera system

Implement energy dispersion

Definition at line 559 of file GCTAResponseIrf.cpp.

References aeff(), GCTAPointing::azimuth(), gammalib::cta_pnt(), GObservation::deadc(), GCTAEventAtom::dir(), GPhoton::dir(), GCTAPointing::dir(), GSkyDir::dist(), edisp(), GPhoton::energy(), G_MC, GCTAPointing::instdir(), GEnergy::log10TeV(), GCTAPsf::mc(), GCTAEdisp::mc(), psf(), gammalib::rad2deg, GSkyDir::rotate_deg(), gammalib::str(), GEnergy::TeV(), GPhoton::time(), GRan::uniform(), use_edisp(), gammalib::warning(), and GCTAPointing::zenith().

double GCTAResponseIrf::nirf ( const GPhoton photon,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const

Return spatial integral of Instrument Response Function.

Parameters
[in]photonPhoton.
[in]obsEngObserved event energy.
[in]obsTimeObserved event time.
[in]obsObservation.

Computes the integral of the instrument response function over the Region of Interest

\[ R(E',t'|p,E,t) = \int_{\rm ROI} R(p',E',t'|p,E,t) dp' \]

Definition at line 1964 of file GCTAResponseIrf.cpp.

References aeff(), GCTAPointing::azimuth(), gammalib::cta_event_list(), gammalib::cta_obs(), GObservation::deadc(), GPhoton::dir(), GCTAPointing::dir(), GSkyDir::dist(), edisp(), GPhoton::energy(), G_NIRF, gammalib::is_infinite(), gammalib::is_notanumber(), GEnergy::log10TeV(), npsf(), nroi(), GCTAObservation::pointing(), GCTAEventList::roi(), GPhoton::time(), use_edisp(), and GCTAPointing::zenith().

Referenced by cta_nroi_radial_kern_omega::eval(), cta_nroi_elliptical_kern_omega::eval(), cta_nroi_diffuse_kern_phi::eval(), and nroi_ptsrc().

double GCTAResponseIrf::npsf ( const GSkyDir srcDir,
const double &  srcLogEng,
const GTime srcTime,
const GCTAPointing pnt,
const GCTARoi roi 
) const

Return result of PSF integration over ROI.

Parameters
[in]srcDirTrue photon direction.
[in]srcLogEngLog10 of true photon energy (E/TeV).
[in]srcTimeTrue photon arrival time (not used).
[in]pntCTA pointing.
[in]roiCTA region of interest.

This method integrates the PSF over the circular region of interest (ROI). Integration is done in a polar coordinate system centred on the PSF since the PSF is assumed to be azimuthally symmetric. The polar integration is done using the method npsf_kern_rad_azsym() that computes analytically the arclength that is comprised within the ROI.

Note that the integration is only performed when the PSF is spilling out of the ROI border, otherwise the integral is simply 1. Numerical integration is done using the standard Romberg method. The integration boundaries are computed so that only the PSF section that falls in the ROI is considered.

Todo:

Enhance romberg() integration method for small integration regions (see comment about kluge below)

Implement phi dependence in camera system

Definition at line 2054 of file GCTAResponseIrf.cpp.

References GCTAPointing::azimuth(), GCTARoi::centre(), gammalib::deg2rad, GCTAPointing::dir(), GCTAInstDir::dir(), GSkyDir::dist(), GIntegral::fixed_iter(), gammalib::is_infinite(), gammalib::is_notanumber(), psf_delta_max(), gammalib::rad2deg, GCTARoi::radius(), GIntegral::romberg(), GIntegral::trapzd(), and GCTAPointing::zenith().

Referenced by nirf().

double GCTAResponseIrf::nroi ( const GModelSky model,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const
virtual

Return integral of event probability for a given sky model over ROI.

Parameters
[in]modelSky model.
[in]obsEngObserved photon energy.
[in]obsTimeObserved photon arrival time.
[in]obsObservation.
Returns
Event probability.
Exceptions
GException::feature_not_implementedMethod is not implemented.

Computes the integral

\[ N_{\rm ROI}(E',t') = \int_{\rm ROI} P(p',E',t') dp' \]

of the event probability

\[ P(p',E',t') = \int \int \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \, dE \, dt \]

for a given sky model \(S(p,E,t)\) and response function \(R(p',E',t'|p,E,t)\) over the Region of Interest (ROI).

Implements GCTAResponse.

Definition at line 434 of file GCTAResponseIrf.cpp.

References ebounds(), edisp(), GEbounds::emax(), GEbounds::emin(), GCTAEdisp::etrue_bounds(), GModelTemporal::eval(), GModelSpectral::eval(), GIntegral::fixed_iter(), GModel::has_scales(), GObservation::instrument(), gammalib::is_infinite(), gammalib::is_notanumber(), GEnergy::MeV(), GIntegral::romberg(), GModel::scale(), GEbounds::size(), GModelSky::spectral(), GModelSky::temporal(), and use_edisp().

Referenced by cta_nroi_kern::eval(), nirf(), nroi(), nroi_composite(), nroi_diffuse(), nroi_elliptical(), nroi_ptsrc(), and nroi_radial().

double GCTAResponseIrf::nroi ( const GModelSky model,
const GEnergy srcEng,
const GTime srcTime,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const

Return spatial integral of sky model.

Parameters
[in]modelSky Model.
[in]srcEngTrue photon energy.
[in]srcTimeTrue photon arrival time.
[in]obsEngObserved event energy.
[in]obsTimeObserved event arrival time.
[in]obsObservation.

Computes the integral

\[ N_{\rm ROI}(E',t'|E,t) = \int_{\rm ROI} P(p',E',t'|E,t) dp' \]

of

\[ P(p',E',t'|E,t) = \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \]

over the Region of Interest (ROI) for a sky model \(S(p,E,t)\) and the response function \(R(p',E',t'|p,E,t)\).

Definition at line 1892 of file GCTAResponseIrf.cpp.

References GModelSpatial::code(), GResponseCache::contains(), GMODEL_SPATIAL_COMPOSITE, GMODEL_SPATIAL_DIFFUSE, GMODEL_SPATIAL_ELLIPTICAL, GMODEL_SPATIAL_POINT_SOURCE, GMODEL_SPATIAL_RADIAL, GModelSpatial::has_free_pars(), GObservation::id(), GResponse::m_nroi_cache, GResponse::m_use_nroi_cache, GModel::name(), nroi(), nroi_composite(), nroi_diffuse(), nroi_elliptical(), nroi_ptsrc(), nroi_radial(), GResponseCache::set(), and GModelSky::spatial().

double GCTAResponseIrf::nroi_composite ( const GModelSky model,
const GEnergy srcEng,
const GTime srcTime,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const
private

Return spatial integral of composite source model.

Parameters
[in]modelSky Model.
[in]srcEngTrue photon energy.
[in]srcTimeTrue photon arrival time.
[in]obsEngObserved event energy.
[in]obsTimeObserved event arrival time.
[in]obsObservation.

Computes the integral

\[ N_{\rm ROI}(E',t'|E,t) = \int_{\rm ROI} P(p',E',t'|E,t) dp' \]

of

\[ P(p',E',t'|E,t) = \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \]

over the Region of Interest (ROI) for a composite source model \(S(p,E,t)\) and the response function \(R(p',E',t'|p,E,t)\).

Definition at line 3511 of file GCTAResponseIrf.cpp.

References GModelSpatialComposite::component(), GModelSpatialComposite::components(), nroi(), GModelSpatialComposite::scale(), GModelSky::spatial(), sum(), and GModelSpatialComposite::sum_of_scales().

Referenced by nroi().

double GCTAResponseIrf::nroi_diffuse ( const GModelSky model,
const GEnergy srcEng,
const GTime srcTime,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const
private

Return spatial integral of diffuse source model.

Parameters
[in]modelSky Model.
[in]srcEngTrue photon energy.
[in]srcTimeTrue photon arrival time.
[in]obsEngObserved event energy.
[in]obsTimeObserved event arrival time.
[in]obsObservation.
Returns
Spatial integral of diffuse source model.
Exceptions
GException::invalid_argumentInvalid spatial model specified.

Computes the integral

\[ N_{\rm ROI}(E',t'|E,t) = \int_{\rm ROI} P(p',E',t'|E,t) dp' \]

of

\[ P(p',E',t'|E,t) = \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \]

over the Region of Interest (ROI) for a diffuse source model \(S(p,E,t)\) and the response function \(R(p',E',t'|p,E,t)\).

Definition at line 3376 of file GCTAResponseIrf.cpp.

References GCTAPointing::azimuth(), GCTARoi::centre(), gammalib::cta_event_list(), gammalib::cta_obs(), gammalib::deg2rad, GCTAInstDir::dir(), GMatrix::eulery(), GMatrix::eulerz(), GIntegral::fixed_iter(), G_NROI_DIFFUSE, gammalib::is_infinite(), gammalib::is_notanumber(), gammalib::iter_phi(), gammalib::iter_rho(), GEnergy::log10TeV(), nroi(), GCTAObservation::pointing(), psf_delta_max(), GCTARoi::radius(), GCTAEventList::roi(), GIntegral::romberg(), GModelSky::spatial(), and GCTAPointing::zenith().

Referenced by nroi().

double GCTAResponseIrf::nroi_elliptical ( const GModelSky model,
const GEnergy srcEng,
const GTime srcTime,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const
private

Return spatial integral of elliptical source model.

Parameters
[in]modelSky Model.
[in]srcEngTrue photon energy.
[in]srcTimeTrue photon arrival time.
[in]obsEngObserved event energy.
[in]obsTimeObserved event arrival time.
[in]obsObservation.
Returns
Spatial integral of elliptical source model.
Exceptions
GException::invalid_argumentInvalid spatial model specified.

Computes the integral

\[ N_{\rm ROI}(E',t'|E,t) = \int_{\rm ROI} P(p',E',t'|E,t) dp' \]

of

\[ P(p',E',t'|E,t) = \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \]

over the Region of Interest (ROI) for an elliptical source model \(S(p,E,t)\) and the response function \(R(p',E',t'|p,E,t)\).

Definition at line 3181 of file GCTAResponseIrf.cpp.

References GCTAPointing::azimuth(), GCTARoi::centre(), gammalib::centre(), gammalib::cta_event_list(), gammalib::cta_obs(), GSkyDir::dec_deg(), gammalib::deg2rad, GCTAInstDir::dir(), GModelSpatialElliptical::dir(), GMatrix::eulery(), GMatrix::eulerz(), GIntegral::fixed_iter(), G_NROI_ELLIPTICAL, gammalib::is_infinite(), gammalib::is_notanumber(), gammalib::iter_phi(), gammalib::iter_rho(), GEnergy::log10TeV(), nroi(), gammalib::pihalf, GCTAObservation::pointing(), GSkyDir::posang(), GModelSpatialElliptical::posangle(), psf_delta_max(), GSkyDir::ra_deg(), GCTARoi::radius(), GCTAEventList::roi(), GIntegral::romberg(), GModelSpatialElliptical::semimajor(), GModelSpatialElliptical::semiminor(), GModelSky::spatial(), GModelSpatialElliptical::theta_max(), and GCTAPointing::zenith().

Referenced by nroi().

double GCTAResponseIrf::nroi_ptsrc ( const GModelSky model,
const GEnergy srcEng,
const GTime srcTime,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const
private

Return spatial integral of point source model.

Parameters
[in]modelSky Model.
[in]srcEngTrue photon energy.
[in]srcTimeTrue photon arrival time.
[in]obsEngObserved event energy.
[in]obsTimeObserved event arrival time.
[in]obsObservation.

Computes the integral

\[ N_{\rm ROI}(E',t'|E,t) = \int_{\rm ROI} P(p',E',t'|E,t) dp' \]

of

\[ P(p',E',t'|E,t) = \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \]

over the Region of Interest (ROI) for a point source model \(S(p,E,t)\) and the response function \(R(p',E',t'|p,E,t)\).

Definition at line 2926 of file GCTAResponseIrf.cpp.

References GModelSpatialPointSource::dir(), nirf(), nroi(), and GModelSky::spatial().

Referenced by nroi().

double GCTAResponseIrf::nroi_radial ( const GModelSky model,
const GEnergy srcEng,
const GTime srcTime,
const GEnergy obsEng,
const GTime obsTime,
const GObservation obs 
) const
private

Return spatial integral of radial source model.

Parameters
[in]modelSky Model.
[in]srcEngTrue photon energy.
[in]srcTimeTrue photon arrival time.
[in]obsEngObserved event energy.
[in]obsTimeObserved event arrival time.
[in]obsObservation.
Returns
Spatial integral of radial source model.
Exceptions
GException::invalid_argumentInvalid spatial model specified.

Computes the integral

\[ N_{\rm ROI}(E',t'|E,t) = \int_{\rm ROI} P(p',E',t'|E,t) dp' \]

of

\[ P(p',E',t'|E,t) = \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \]

over the Region of Interest (ROI) for a radial source model \(S(p,E,t)\) and the response function \(R(p',E',t'|p,E,t)\).

Definition at line 2978 of file GCTAResponseIrf.cpp.

References GCTAPointing::azimuth(), GCTARoi::centre(), gammalib::centre(), gammalib::cta_event_list(), gammalib::cta_obs(), GSkyDir::dec_deg(), gammalib::deg2rad, GCTAInstDir::dir(), GModelSpatialRadial::dir(), GMatrix::eulery(), GMatrix::eulerz(), GIntegral::fixed_iter(), G_NROI_RADIAL, gammalib::is_infinite(), gammalib::is_notanumber(), gammalib::iter_phi(), gammalib::iter_rho(), GEnergy::log10TeV(), nroi(), GCTAObservation::pointing(), GSkyDir::posang(), psf_delta_max(), GSkyDir::ra_deg(), GCTARoi::radius(), GModelSpatialRadialRing::radius(), GModelSpatialRadialShell::radius(), GCTAEventList::roi(), GIntegral::romberg(), GModelSky::spatial(), GModelSpatialRadial::theta_max(), and GCTAPointing::zenith().

Referenced by nroi().

void GCTAResponseIrf::offset_sigma ( const double &  sigma)

Set offset angle dependence (degrees)

Parameters
[in]sigmaOffset angle dependence value (degrees).

Set the offset angle dependence for 1D effective area functions. The method set the sigma value in case that the effective area function is of type GCTAAeffArf or GCTAAeffPerfTable. Otherwise, nothing will be done.

Definition at line 1565 of file GCTAResponseIrf.cpp.

References m_aeff, GCTAAeffPerfTable::sigma(), and GCTAAeffArf::sigma().

double GCTAResponseIrf::offset_sigma ( void  ) const

Return offset angle dependence (degrees)

Returns
Offset angle dependence value (degrees).

Return the offset angle dependence for 1D effective area functions. The method returns the sigma value in case that the effective area function is of type GCTAAeffArf or GCTAAeffPerfTable. Otherwise, 0.0 will be returned.

Definition at line 1594 of file GCTAResponseIrf.cpp.

References m_aeff, GCTAAeffPerfTable::sigma(), and GCTAAeffArf::sigma().

GCTAResponseIrf & GCTAResponseIrf::operator= ( const GCTAResponseIrf rsp)
virtual

Assignment operator.

Parameters
[in]rspCTA response.
Returns
CTA response.

Assigns CTA response object to another CTA response object. The assignment performs a deep copy of all information, hence the original object from which the assignment has been performed can be destroyed after this operation without any loss of information.

Definition at line 249 of file GCTAResponseIrf.cpp.

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

std::string GCTAResponseIrf::print ( const GChatter chatter = NORMAL) const
virtual
const GCTAPsf * GCTAResponseIrf::psf ( void  ) const
inline

Return pointer to point spread function.

Returns
Pointer to point spread function.

Definition at line 374 of file GCTAResponseIrf.hpp.

References m_psf.

Referenced by cta_npsf_kern_rad_azsym::eval(), cta_irf_radial_kern_omega::eval(), cta_irf_elliptical_kern_omega::eval(), cta_irf_diffuse_kern_theta::eval(), irf(), mc(), psf(), and write().

void GCTAResponseIrf::psf ( GCTAPsf psf)
inline

Set pointer to point spread function.

Parameters
[in]psfPointer to point spread function.

Definition at line 386 of file GCTAResponseIrf.hpp.

References GCTAPsf::clone(), and m_psf.

double GCTAResponseIrf::psf ( const double &  delta,
const double &  theta,
const double &  phi,
const double &  zenith,
const double &  azimuth,
const double &  srcLogEng 
) const

Return point spread function (in units of sr^-1)

Parameters
[in]deltaAngular separation between true and measured photon directions (radians).
[in]thetaRadial offset angle of photon in camera (radians).
[in]phiPolar angle of photon in camera (radians).
[in]zenithZenith angle of telescope pointing (radians).
[in]azimuthAzimuth angle of telescope pointing (radians).
[in]srcLogEngLog10 of true photon energy (E/TeV).
Exceptions
GException::invalid_valueNo point spread function information found.

Returns the point spread function for a given offset angle as function of the true photon position in the camera system and the telescope pointing direction in the Earth system.

Definition at line 1776 of file GCTAResponseIrf.cpp.

References G_PSF, m_psf, and psf().

double GCTAResponseIrf::psf_delta_max ( const double &  theta,
const double &  phi,
const double &  zenith,
const double &  azimuth,
const double &  srcLogEng 
) const

Return maximum angular separation (in radians)

Parameters
[in]thetaRadial offset angle in camera (radians).
[in]phiPolar angle in camera (radians).
[in]zenithZenith angle of telescope pointing (radians).
[in]azimuthAzimuth angle of telescope pointing (radians).
[in]srcLogEngLog10 of true photon energy (E/TeV).
Exceptions
GException::invalid_valueNo point spread function information found.

This method returns the maximum angular separation between true and measured photon directions for which the PSF is non zero as function of the true photon position in the camera system and the telescope pointing direction in the Earth system.

Definition at line 1817 of file GCTAResponseIrf.cpp.

References GCTAPsf::delta_max(), G_PSF_DELTA_MAX, and m_psf.

Referenced by irf(), irf_diffuse(), irf_elliptical(), irf_radial(), npsf(), nroi_diffuse(), nroi_elliptical(), and nroi_radial().

void GCTAResponseIrf::read ( const GXmlElement xml)
virtual

Read response from XML element.

Parameters
[in]xmlXML element.

Reads information for a CTA observation from an XML element. The calibration database and response name can be specified using

<observation name="..." id="..." instrument="...">
  ...
  <parameter name="Calibration" database="..." response="..."/>
</observation>

If even more control is required over the response, individual file names can be specified using

<observation name="..." id="..." instrument="...">
  ...
  <parameter name="EffectiveArea"       file="..."/>
  <parameter name="PointSpreadFunction" file="..."/>
  <parameter name="EnergyDispersion"    file="..."/>
  <parameter name="Background"          file="..."/>
</observation>

Implements GCTAResponse.

Definition at line 675 of file GCTAResponseIrf.cpp.

References aeff(), GXmlElement::attribute(), background(), caldb(), G_READ, GXmlElement::has_attribute(), load(), load_aeff(), load_background(), load_edisp(), load_psf(), GCTAAeffArf::remove_thetacut(), GCTAAeffArf::scale(), GCTAAeffPerfTable::sigma(), GCTABackgroundPerfTable::sigma(), GCTAAeffArf::sigma(), gammalib::strip_whitespace(), GCTAAeffArf::thetacut(), gammalib::todouble(), gammalib::tolower(), gammalib::xml_get_par(), and gammalib::xml_has_par().

Referenced by GCTAResponseIrf().

const std::string & GCTAResponseIrf::rspname ( void  ) const
inline

Return response name.

Returns
Response name.

Definition at line 337 of file GCTAResponseIrf.hpp.

References m_rspname.

Referenced by load(), and GCTAOnOffObservation::set().

bool GCTAResponseIrf::use_edisp ( void  ) const
inlinevirtual

Signal if response uses energy dispersion.

Returns
True if response uses energy dispersion

Signals if the response uses energy dispersion. This implies that the apply_edisp flag has been set to true and that energy dispersion response information is available.

Implements GCTAResponse.

Definition at line 262 of file GCTAResponseIrf.hpp.

References m_apply_edisp, and m_edisp.

Referenced by cta_irf_radial_kern_omega::eval(), cta_irf_elliptical_kern_omega::eval(), cta_irf_diffuse_kern_phi::eval(), irf(), mc(), nirf(), nroi(), and print().

bool GCTAResponseIrf::use_tdisp ( void  ) const
inlinevirtual

Signal if time dispersion will be used.

Returns
False.

Implements GCTAResponse.

Definition at line 298 of file GCTAResponseIrf.hpp.

void GCTAResponseIrf::write ( GXmlElement xml) const
virtual

Write response information into XML element.

Parameters
[in]xmlXML element.

Writes information for a CTA response into an XML element. If the calibration database and response name had been specified, the following output is written

<observation name="..." id="..." instrument="...">
  ...
  <parameter name="Calibration" database="..." response="..."/>
</observation>

If even more control was required over the response and individual file names were specified, the following output is written

<observation name="..." id="..." instrument="...">
  ...
  <parameter name="EffectiveArea"       file="..."/>
  <parameter name="PointSpreadFunction" file="..."/>
  <parameter name="EnergyDispersion"    file="..."/>
  <parameter name="Background"          file="..."/>
</observation>

Implements GCTAResponse.

Definition at line 889 of file GCTAResponseIrf.cpp.

References aeff(), GXmlElement::attribute(), background(), edisp(), GCTAAeff::filename(), G_WRITE, GFilename::is_empty(), m_xml_caldb, m_xml_rspname, psf(), GCTAAeffArf::scale(), GCTAAeffPerfTable::sigma(), GCTAAeffArf::sigma(), gammalib::str(), GCTAAeffArf::thetacut(), and gammalib::xml_need_par().

Member Data Documentation

GCTAAeff* GCTAResponseIrf::m_aeff
private

Effective area.

Definition at line 214 of file GCTAResponseIrf.hpp.

Referenced by aeff(), copy_members(), free_members(), init_members(), is_valid(), load(), load_aeff(), offset_sigma(), and print().

bool GCTAResponseIrf::m_apply_edisp
mutableprivate

Apply energy dispersion.

Definition at line 218 of file GCTAResponseIrf.hpp.

Referenced by apply_edisp(), copy_members(), init_members(), and use_edisp().

GCTABackground* GCTAResponseIrf::m_background
private

Background.

Definition at line 217 of file GCTAResponseIrf.hpp.

Referenced by background(), copy_members(), free_members(), init_members(), load_background(), and print().

GCaldb GCTAResponseIrf::m_caldb
private

Calibration database.

Definition at line 212 of file GCTAResponseIrf.hpp.

Referenced by caldb(), copy_members(), GCTAResponseIrf(), init_members(), load(), and print().

GCTAEdisp* GCTAResponseIrf::m_edisp
private

Energy dispersion.

Definition at line 216 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), edisp(), free_members(), init_members(), load_edisp(), print(), and use_edisp().

double GCTAResponseIrf::m_hi_safe_thres
private

Safe high energy threshold.

Definition at line 220 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), hi_safe_thres(), init_members(), load_aeff(), and print().

double GCTAResponseIrf::m_lo_safe_thres
private

Safe low energy threshold.

Definition at line 219 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), init_members(), lo_safe_thres(), load_aeff(), and print().

GCTAPsf* GCTAResponseIrf::m_psf
private

Point spread function.

Definition at line 215 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), free_members(), init_members(), is_valid(), load_psf(), print(), psf(), and psf_delta_max().

std::string GCTAResponseIrf::m_rspname
private

Name of the instrument response.

Definition at line 213 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), init_members(), load(), print(), and rspname().

std::string GCTAResponseIrf::m_xml_caldb
private

Calibration database string in XML file.

Definition at line 223 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), init_members(), load(), and write().

std::string GCTAResponseIrf::m_xml_rspname
private

Response name in XML file.

Definition at line 224 of file GCTAResponseIrf.hpp.

Referenced by copy_members(), init_members(), load(), and write().


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