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

Interface for the COMPTEL instrument response function. More...

#include <GCOMResponse.hpp>

Inheritance diagram for GCOMResponse:
GResponse GBase

Public Member Functions

 GCOMResponse (void)
 Void constructor. More...
 
 GCOMResponse (const GCOMResponse &rsp)
 Copy constructor. More...
 
 GCOMResponse (const GCaldb &caldb, const std::string &rspname)
 Response constructor. More...
 
virtual ~GCOMResponse (void)
 Destructor. More...
 
virtual GCOMResponseoperator= (const GCOMResponse &rsp)
 Assignment operator. More...
 
virtual void clear (void)
 Clear instance. More...
 
virtual GCOMResponseclone (void) const
 Clone instance. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual bool use_edisp (void) const
 Signal if energy dispersion will be used. More...
 
virtual bool use_tdisp (void) const
 Signal if time dispersion will be used. 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 std::string print (const GChatter &chatter=NORMAL) const
 Print COMPTEL response information. 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 COMPTEL response. More...
 
void read (const GFitsImage &hdu)
 Read COMPTEL response from FITS image. More...
 
void write (GFitsImageFloat &image) const
 Write COMPTEL response into FITS image. More...
 
void load_cache (const GFilename &filename)
 Load response cache. More...
 
void save_cache (const GFilename &filename) const
 Save response cache. 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 GCOMResponse &rsp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
virtual GVector irf_ptsrc (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to point source. More...
 
virtual GVector irf_radial (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to radial source. More...
 
virtual GVector irf_elliptical (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to elliptical source. More...
 
virtual GVector irf_diffuse (const GModelSky &model, const GObservation &obs, GMatrix *gradients=NULL) const
 Return instrument response to diffuse source. More...
 
GVector irf_extended (const GModelSky &model, const GObservation &obs, const GSkyDir &model_dir, const double &theta_max, GMatrix *gradients=NULL) const
 

Private Attributes

GCaldb m_caldb
 Calibration database. More...
 
std::string m_rspname
 Response name. More...
 
std::vector< double > m_iaq
 IAQ array. More...
 
int m_phigeo_bins
 Number of Phigeo bins. More...
 
int m_phibar_bins
 Number of Phibar bins. More...
 
double m_phigeo_ref_value
 Phigeo reference value (deg) More...
 
double m_phigeo_ref_pixel
 Phigeo reference pixel (starting from 1) More...
 
double m_phigeo_bin_size
 Phigeo binsize (deg) More...
 
double m_phigeo_min
 Phigeo value of first bin (deg) More...
 
double m_phibar_ref_value
 Phigeo reference value (deg) More...
 
double m_phibar_ref_pixel
 Phigeo reference pixel (starting from 1) More...
 
double m_phibar_bin_size
 Phigeo binsize (deg) More...
 
double m_phibar_min
 Phigeo value of first bin (deg) More...
 

Friends

class GCOMDri
 

Additional Inherited Members

- 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_ptsrc (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return instrument response to point source. More...
 
virtual double irf_radial (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return instrument response to radial source. More...
 
virtual double irf_elliptical (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return instrument response to elliptical source. More...
 
virtual double irf_diffuse (const GEvent &event, const GSource &source, const GObservation &obs) const
 Return instrument response to diffuse source. 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_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

Interface for the COMPTEL instrument response function.

Definition at line 56 of file GCOMResponse.hpp.

Constructor & Destructor Documentation

GCOMResponse::GCOMResponse ( void  )

Void constructor.

Creates an empty COMPTEL response.

Definition at line 98 of file GCOMResponse.cpp.

References init_members().

Referenced by clone().

GCOMResponse::GCOMResponse ( const GCOMResponse rsp)

Copy constructor.

Parameters
[in]rspCOM response.

Definition at line 113 of file GCOMResponse.cpp.

References copy_members(), and init_members().

GCOMResponse::GCOMResponse ( const GCaldb caldb,
const std::string &  iaqname 
)

Response constructor.

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

Create COMPTEL response by loading an IAQ file from a calibration database.

Definition at line 135 of file GCOMResponse.cpp.

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

GCOMResponse::~GCOMResponse ( void  )
virtual

Destructor.

Destroys instance of COMPTEL response object.

Definition at line 157 of file GCOMResponse.cpp.

References free_members().

Member Function Documentation

void GCOMResponse::caldb ( const GCaldb caldb)
inline

Set calibration database.

Parameters
[in]caldbCalibration database.

Sets the calibration database for the COMPTEL response.

Definition at line 197 of file GCOMResponse.hpp.

References caldb(), and m_caldb.

Referenced by GCOMObservation::response().

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

Return calibration database.

Returns
Calibration database.

Definition at line 183 of file GCOMResponse.hpp.

References m_caldb.

Referenced by caldb(), GCOMResponse(), and load().

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

Return class name.

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

Implements GResponse.

Definition at line 147 of file GCOMResponse.hpp.

void GCOMResponse::clear ( void  )
virtual

Clear instance.

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

Implements GResponse.

Definition at line 220 of file GCOMResponse.cpp.

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

Referenced by GCOMObservation::init_members(), load(), and GCOMObservation::response().

GCOMResponse * GCOMResponse::clone ( void  ) const
virtual

Clone instance.

Returns
Pointer to deep copy of COMPTEL response.

Implements GResponse.

Definition at line 240 of file GCOMResponse.cpp.

References GCOMResponse().

void GCOMResponse::copy_members ( const GCOMResponse rsp)
private

Copy class members.

Parameters
[in]rspCOMPTEL response.

Definition at line 819 of file GCOMResponse.cpp.

References m_caldb, m_iaq, m_phibar_bin_size, m_phibar_bins, m_phibar_min, m_phibar_ref_pixel, m_phibar_ref_value, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, m_phigeo_ref_pixel, m_phigeo_ref_value, and m_rspname.

Referenced by GCOMResponse(), and operator=().

GEbounds GCOMResponse::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.
Exceptions
GException::feature_not_implementedMethod is not implemented.

Implements GResponse.

Definition at line 462 of file GCOMResponse.cpp.

References G_EBOUNDS.

void GCOMResponse::free_members ( void  )
private

Delete class members.

Definition at line 844 of file GCOMResponse.cpp.

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

void GCOMResponse::init_members ( void  )
private
double GCOMResponse::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.
Returns
Instrument response function ( \(cm^2 sr^{-1}\))
Exceptions
GException::invalid_argumentObservation is not a COMPTEL observation. Event is not a COMPTEL event bin.
GException::invalid_valueResponse not initialised with a valid IAQ

Returns the instrument response function for a given observed photon direction as function of the assumed true photon direction. The result is given by

\[ {\tt IRF} = \frac{{\tt IAQ} \times {\tt DRG} \times {\tt DRX}} {T \times {\tt TOFCOR}} \times {\tt PHASECOR} \]

where

  • \({\tt IRF}\) is the instrument response function ( \(cm^2 sr^{-1}\)),
  • \({\tt IAQ}\) is the COMPTEL response matrix ( \(sr^{-1}\)),
  • \({\tt DRG}\) is the geometry factor (probability),
  • \({\tt DRX}\) is the exposure ( \(cm^2 s\)),
  • \(T\) is the ontime ( \(s\)), and
  • \({\tt TOFCOR}\) is a correction that accounts for the Time of Flight selection window.
  • \({\tt PHASECOR}\) is a correction that accounts for pulsar phase selection.

The observed photon direction is spanned by the three values \(\Chi\), \(\Psi\), and \(\bar{\varphi})\). \(\Chi\) and \(\Psi\) is the scatter direction of the event, given in sky coordinates. \(\bar{\varphi}\) is the Compton scatter angle, computed from the energy deposits in the two detector planes.

Implements GResponse.

Definition at line 286 of file GCOMResponse.cpp.

References GObservation::deadc(), GCOMInstDir::dir(), GCOMEventBin::dir(), GPhoton::dir(), GSkyDir::dist_deg(), GCOMEventCube::dre(), GCOMObservation::drg(), GCOMObservation::drx(), GObservation::events(), G_IRF, gammalib::is_infinite(), gammalib::is_notanumber(), m_iaq, m_phibar_bin_size, m_phibar_bins, m_phigeo_bin_size, m_phigeo_bins, GCOMDri::map(), GCOMObservation::ontime(), GCOMDri::phase_correction(), GCOMInstDir::phibar(), GPhoton::time(), and GCOMDri::tof_correction().

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

GVector GCOMResponse::irf_diffuse ( const GModelSky model,
const GObservation obs,
GMatrix gradients = NULL 
) const
privatevirtual

Return instrument response to diffuse source.

Parameters
[in]modelSky model.
[in]obsObservation.
[out]gradientsGradients matrix.
Returns
Instrument response to diffuse source for all events in observation ( \(cm^2\)).
Exceptions
GException::invalid_argumentObservation is not a COMPTEL observation.
GException::invalid_valueResponse not initialised with a valid IAQ Incompatible IAQ encountered

Returns the instrument response to a diffuse source for all events in the observations. The diffuse source may be energy dependent.

The computation is done by integrating the diffuse model for each pixel in Chi and Psi over a circular region centred on the Chi/Psi pixel with a radius equal to the maximum Phigeo value. The radial integration is done by looping over all Phigeo bins of the response. For each Phigeo value, the azimuthal integration is done by stepping with an angular step size that corresponds to the Phigeo step size (which typically is 1 deg).

gradients is an optional matrix where the number of rows corresponds to the number of events in the observation and the number of columns corresponds to the number of spatial model parameters. Since for point sources no gradients are computed, the method does not alter the content of gradients.

Reimplemented from GResponse.

Definition at line 1394 of file GCOMResponse.cpp.

References GSkyMap::contains(), GCOMObservation::deadc(), gammalib::deg2rad, GCOMInstDir::dir(), GCOMEventBin::dir(), GCOMEventCube::dre(), GCOMObservation::drg(), GCOMObservation::drx(), GCOMEventBin::energy(), GModelSpatial::eval(), GObservation::events(), G_IRF_DIFFUSE, GObservation::id(), irf(), m_iaq, m_phibar_bins, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, GCOMDri::map(), GObservation::name(), GCOMEventCube::naxis(), GCOMObservation::ontime(), GCOMDri::phase_correction(), GSkyMap::pixels(), GSkyDir::rotate(), sin(), GCOMEventCube::size(), GModelSky::spatial(), gammalib::str(), GCOMEventBin::time(), GCOMDri::tof_correction(), and gammalib::twopi.

GVector GCOMResponse::irf_elliptical ( const GModelSky model,
const GObservation obs,
GMatrix gradients = NULL 
) const
privatevirtual

Return instrument response to elliptical source.

Parameters
[in]modelSky model.
[in]obsObservation.
[out]gradientsGradients matrix.
Returns
Instrument response to elliptical source for all events in observation ( \(cm^2\)).
Exceptions
GException::invalid_argumentObservation is not a COMPTEL observation. Event is not a COMPTEL event bin.
GException::invalid_valueResponse not initialised with a valid IAQ Incompatible IAQ encountered

Returns the instrument response to an elliptical source for all events in the observations. See GCOMResponse::irf_extended for more information.

Reimplemented from GResponse.

Definition at line 1213 of file GCOMResponse.cpp.

References GCOMObservation::deadc(), GSkyDir::dec_deg(), gammalib::deg2rad, GCOMInstDir::dir(), GCOMEventBin::dir(), GModelSpatialElliptical::dir(), GSkyDir::dist(), GCOMEventCube::dre(), GCOMObservation::drg(), GCOMObservation::drx(), GMatrix::eulery(), GMatrix::eulerz(), GObservation::events(), GIntegrals::fixed_iter(), G_IRF_ELLIPTICAL, GObservation::id(), irf(), gammalib::iter_rho(), m_iaq, m_phibar_bins, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, GCOMDri::map(), GObservation::name(), GCOMEventCube::naxis(), GCOMObservation::ontime(), GCOMDri::phase_correction(), GSkyMap::pixels(), GSkyDir::ra_deg(), GIntegrals::romberg(), GCOMEventCube::size(), GModelSky::spatial(), gammalib::str(), GModelSpatialElliptical::theta_max(), and GCOMDri::tof_correction().

GVector GCOMResponse::irf_extended ( const GModelSky model,
const GObservation obs,
const GSkyDir model_dir,
const double &  theta_max,
GMatrix gradients = NULL 
) const
private
GVector GCOMResponse::irf_ptsrc ( const GModelSky model,
const GObservation obs,
GMatrix gradients = NULL 
) const
privatevirtual

Return instrument response to point source.

Parameters
[in]modelSky model.
[in]obsObservation.
[out]gradientsGradients matrix.
Returns
Instrument response to point source for all events in observation ( \(cm^2\)).
Exceptions
GException::invalid_argumentObservation is not a COMPTEL observation. Event is not a COMPTEL event bin.
GException::invalid_valueResponse not initialised with a valid IAQ Incompatible IAQ encountered

Returns the instrument response to a point source for all events in the observations.

gradients is an optional matrix where the number of rows corresponds to the number of events in the observation and the number of columns corresponds to the number of spatial model parameters. Since for point sources no gradients are computed, the method does not alter the content of gradients.

Reimplemented from GResponse.

Definition at line 876 of file GCOMResponse.cpp.

References GCOMObservation::deadc(), GCOMInstDir::dir(), GCOMEventBin::dir(), GSkyDir::dist_deg(), GCOMEventCube::dre(), GCOMObservation::drg(), GCOMObservation::drx(), GObservation::events(), G_IRF_PTSRC, GObservation::id(), irf(), m_iaq, m_phibar_bins, m_phigeo_bin_size, m_phigeo_bins, GCOMDri::map(), GObservation::name(), GCOMEventCube::naxis(), GCOMObservation::ontime(), GCOMDri::phase_correction(), GSkyMap::pixels(), GCOMEventCube::size(), GModelSky::spatial(), gammalib::str(), and GCOMDri::tof_correction().

GVector GCOMResponse::irf_radial ( const GModelSky model,
const GObservation obs,
GMatrix gradients = NULL 
) const
privatevirtual

Return instrument response to radial source.

Parameters
[in]modelSky model.
[in]obsObservation.
[out]gradientsGradients matrix.
Returns
Instrument response to radial source for all events in observation ( \(cm^2\)).
Exceptions
GException::invalid_argumentObservation is not a COMPTEL observation. Event is not a COMPTEL event bin.
GException::invalid_valueResponse not initialised with a valid IAQ Incompatible IAQ encountered

Returns the instrument response to a radial source for all events in the observations. See GCOMResponse::irf_extended for more information.

Reimplemented from GResponse.

Definition at line 1042 of file GCOMResponse.cpp.

References GCOMObservation::deadc(), GSkyDir::dec_deg(), gammalib::deg2rad, GCOMInstDir::dir(), GCOMEventBin::dir(), GModelSpatialRadial::dir(), GSkyDir::dist(), GCOMEventCube::dre(), GCOMObservation::drg(), GCOMObservation::drx(), GMatrix::eulery(), GMatrix::eulerz(), GObservation::events(), GIntegrals::fixed_iter(), G_IRF_RADIAL, GObservation::id(), irf(), gammalib::iter_rho(), m_iaq, m_phibar_bins, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, GCOMDri::map(), GObservation::name(), GCOMEventCube::naxis(), GCOMObservation::ontime(), GCOMDri::phase_correction(), GSkyMap::pixels(), GSkyDir::ra_deg(), GIntegrals::romberg(), GCOMEventCube::size(), GModelSky::spatial(), gammalib::str(), GModelSpatialRadial::theta_max(), and GCOMDri::tof_correction().

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

Load COMPTEL response.

Parameters
[in]rspnameCOMPTEL response name.

Loads the COMPTEL response with specified name rspname.

The method first attempts to interpret rspname as a file name and to load the corresponding response.

If rspname is not a FITS file the method searches for an appropriate response in the calibration database. If no appropriate response is found, the method takes the CALDB root path and response name to build the full path to the response file, and tries to load the response from these paths.

If also this fails an exception is thrown.

Definition at line 494 of file GCOMResponse.cpp.

References caldb(), clear(), GFits::close(), GFilename::exists(), GCaldb::filename(), gammalib::filepath(), GFits::image(), GFilename::is_empty(), GFilename::is_fits(), m_caldb, m_rspname, read(), GCaldb::rootdir(), and rspname().

Referenced by GCOMResponse(), and GCOMObservation::response().

void GCOMResponse::load_cache ( const GFilename filename)

Load response cache.

Parameters
[in]filenameResponse cache filename.

Loads response cache from FITS file.

Definition at line 697 of file GCOMResponse.cpp.

References GResponseVectorCache::load(), and GResponse::m_irf_vector_cache.

Referenced by GCOMObservation::read().

double GCOMResponse::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
0.0
Exceptions
GException::feature_not_implementedMethod is not implemented.

Implements GResponse.

Definition at line 438 of file GCOMResponse.cpp.

References G_NROI.

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

Assignment operator.

Parameters
[in]rspCOMPTEL response.
Returns
COMPTEL response.

Assigns COMPTEL response object to another COMPTEL 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 184 of file GCOMResponse.cpp.

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

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

Print COMPTEL response information.

Parameters
[in]chatterChattiness.
Returns
String containing COMPTEL response information.

Implements GResponse.

Definition at line 730 of file GCOMResponse.cpp.

References EXPLICIT, m_caldb, m_phibar_bin_size, m_phibar_bins, m_phibar_min, m_phibar_ref_pixel, m_phibar_ref_value, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, m_phigeo_ref_pixel, m_phigeo_ref_value, m_rspname, gammalib::parformat(), GCaldb::print(), SILENT, gammalib::str(), and VERBOSE.

void GCOMResponse::read ( const GFitsImage image)

Read COMPTEL response from FITS image.

Parameters
[in]imageFITS image.

Read the COMPTEL response from IAQ FITS file and convert the IAQ values into a probability per steradian.

The IAQ values are divided by the solid angle of a Phigeo bin which is given by

\begin{eqnarray*} \Omega & = & 2 \pi \left[ \left( 1 - \cos \left( \varphi_{\rm geo} + \frac{1}{2} \Delta \varphi_{\rm geo} \right) \right) - \left( 1 - \cos \left( \varphi_{\rm geo} - \frac{1}{2} \Delta \varphi_{\rm geo} \right) \right) \right] \\ &=& 2 \pi \left[ \cos \left( \varphi_{\rm geo} - \frac{1}{2} \Delta \varphi_{\rm geo} \right) - \cos \left( \varphi_{\rm geo} + \frac{1}{2} \Delta \varphi_{\rm geo} \right) \right] \\ &=& 4 \pi \sin \left( \varphi_{\rm geo} \right) \sin \left( \frac{1}{2} \Delta \varphi_{\rm geo} \right) \end{eqnarray*}

Definition at line 576 of file GCOMResponse.cpp.

References gammalib::deg2rad, gammalib::fourpi, m_iaq, m_phibar_bin_size, m_phibar_bins, m_phibar_min, m_phibar_ref_pixel, m_phibar_ref_value, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, m_phigeo_ref_pixel, m_phigeo_ref_value, GFitsImage::naxes(), GFitsImage::naxis(), GFitsImage::pixel(), GFitsHDU::real(), and sin().

Referenced by load().

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

Return response name.

Returns
Response name.

Definition at line 210 of file GCOMResponse.hpp.

References m_rspname.

Referenced by load(), and GCOMObservation::write().

void GCOMResponse::save_cache ( const GFilename filename) const

Save response cache.

Parameters
[in]filenameResponse cache filename.

Saves response cache from FITS file.

Definition at line 714 of file GCOMResponse.cpp.

References GResponse::m_irf_vector_cache, and GResponseVectorCache::save().

Referenced by GCOMObservation::write().

bool GCOMResponse::use_edisp ( void  ) const
inlinevirtual

Signal if energy dispersion will be used.

Returns
False.

Implements GResponse.

Definition at line 159 of file GCOMResponse.hpp.

bool GCOMResponse::use_tdisp ( void  ) const
inlinevirtual

Signal if time dispersion will be used.

Returns
False.

Implements GResponse.

Definition at line 171 of file GCOMResponse.hpp.

void GCOMResponse::write ( GFitsImageFloat image) const

Write COMPTEL response into FITS image.

Parameters
[in]imageFITS image.

Writes the COMPTEL response into an IAQ FITS file.

Definition at line 642 of file GCOMResponse.cpp.

References GFitsHDU::card(), gammalib::deg2rad, gammalib::fourpi, m_iaq, m_phibar_bin_size, m_phibar_bins, m_phibar_ref_pixel, m_phibar_ref_value, m_phigeo_bin_size, m_phigeo_bins, m_phigeo_min, m_phigeo_ref_pixel, m_phigeo_ref_value, and sin().

Friends And Related Function Documentation

friend class GCOMDri
friend

Definition at line 59 of file GCOMResponse.hpp.

Member Data Documentation

GCaldb GCOMResponse::m_caldb
private

Calibration database.

Definition at line 125 of file GCOMResponse.hpp.

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

std::vector<double> GCOMResponse::m_iaq
private

IAQ array.

Definition at line 127 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), irf(), irf_diffuse(), irf_elliptical(), irf_ptsrc(), irf_radial(), read(), and write().

double GCOMResponse::m_phibar_bin_size
private

Phigeo binsize (deg)

Definition at line 136 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), irf(), print(), read(), and write().

int GCOMResponse::m_phibar_bins
private

Number of Phibar bins.

Definition at line 129 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), irf(), irf_diffuse(), irf_elliptical(), irf_ptsrc(), irf_radial(), print(), read(), and write().

double GCOMResponse::m_phibar_min
private

Phigeo value of first bin (deg)

Definition at line 137 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), print(), and read().

double GCOMResponse::m_phibar_ref_pixel
private

Phigeo reference pixel (starting from 1)

Definition at line 135 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), print(), read(), and write().

double GCOMResponse::m_phibar_ref_value
private

Phigeo reference value (deg)

Definition at line 134 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), print(), read(), and write().

double GCOMResponse::m_phigeo_bin_size
private

Phigeo binsize (deg)

Definition at line 132 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), irf(), irf_diffuse(), irf_elliptical(), irf_ptsrc(), irf_radial(), print(), read(), and write().

int GCOMResponse::m_phigeo_bins
private

Number of Phigeo bins.

Definition at line 128 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), irf(), irf_diffuse(), irf_elliptical(), irf_ptsrc(), irf_radial(), print(), read(), and write().

double GCOMResponse::m_phigeo_min
private

Phigeo value of first bin (deg)

Definition at line 133 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), irf_diffuse(), irf_elliptical(), irf_radial(), print(), read(), and write().

double GCOMResponse::m_phigeo_ref_pixel
private

Phigeo reference pixel (starting from 1)

Definition at line 131 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), print(), read(), and write().

double GCOMResponse::m_phigeo_ref_value
private

Phigeo reference value (deg)

Definition at line 130 of file GCOMResponse.hpp.

Referenced by copy_members(), init_members(), print(), read(), and write().

std::string GCOMResponse::m_rspname
private

Response name.

Definition at line 126 of file GCOMResponse.hpp.

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


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