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

INTEGRAL/SPI instrument response function class. More...

#include <GSPIResponse.hpp>

Inheritance diagram for GSPIResponse:
GResponse GBase

Public Member Functions

 GSPIResponse (void)
 Void constructor. More...
 
 GSPIResponse (const GSPIResponse &rsp)
 Copy constructor. More...
 
 GSPIResponse (const GFilename &rspname)
 Response constructor. More...
 
virtual ~GSPIResponse (void)
 Destructor. More...
 
virtual GSPIResponseoperator= (const GSPIResponse &rsp)
 Assignment operator. More...
 
virtual void clear (void)
 Clear instance. More...
 
virtual GSPIResponseclone (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 INTEGRAL/SPI instrument response for a photon. 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 INTEGRAL/SPI response information. More...
 
void rspname (const GFilename &rspname)
 Set response name. More...
 
const GFilenamerspname (void) const
 Get response group file name. More...
 
bool is_precomputed (void) const
 Signals if response is precomputed. More...
 
const double & energy_keV (void) const
 Return line IRF energy in keV. More...
 
const double & dlogE (void) const
 Return logarithmic step size for continuum IRFs. More...
 
const double & gamma (void) const
 Return power-law index for continuum IRFs. More...
 
void set (const GSPIObservation &obs, const GEnergy &energy=GEnergy())
 Set response for a specific observation. More...
 
double irf_value (const GSkyDir &srcDir, const GSPIEventBin &bin, const int &ireg) const
 Return value of INTEGRAL/SPI instrument response for sky direction and event bin. More...
 
double zenith (const int &ipt, const GSkyDir &dir) const
 Return zenith angle of sky direction for pointing in radians. More...
 
double azimuth (const int &ipt, const GSkyDir &dir) const
 Return azimuth angle of sky direction for pointing in radians. More...
 
void read (const GFits &fits)
 Read SPI response from FITS object. More...
 
void write (GFits &fits) const
 Write SPI response into FITS object. More...
 
void load (const GFilename &filename)
 Load SPI response from file. More...
 
void save (const GFilename &filename, const bool &clobber=false) const
 Save SPI response into file. 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 GSPIResponse &rsp)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void read_detids (const GFits &fits)
 Read detector identifiers from FITS object. More...
 
void read_energies (const GFits &fits)
 Read energies from FITS object. More...
 
void write_detids (GFits &fits) const
 Write detector identifiers into FITS object. More...
 
void write_energies (GFits &fits) const
 Write energies into FITS object. More...
 
void load_irfs (const int &region)
 Load Instrument Response Functions. More...
 
GSkyMap load_irf (const GFilename &rspname) const
 Load IRF as sky map. More...
 
GSkyMap compute_irf (const double &emin, const double &emax) const
 Compute as sky map. More...
 
void set_wcs (const GFitsImage *image)
 Set IRF image limits. More...
 
void set_detids (const GSPIEventCube *cube)
 Set vector of detector identifiers used by the observation. More...
 
void set_cache (const GSPIEventCube *cube)
 Set computation cache. More...
 
int irf_detid (const int &detid) const
 Convert detector identifier into IRF detector identifier. More...
 
double irf_weight (const double &beta, const double &emin, const double &emax) const
 Compute weight of logarithmic energy bin. More...
 

Private Attributes

GFilename m_rspname
 File name of response group. More...
 
std::vector< int > m_detids
 Vector of detector IDs. More...
 
GNodeArray m_energies
 Node array of IRF energies. More...
 
GEbounds m_ebounds
 Energy bounaries of IRF. More...
 
GSkyMap m_irfs
 IRFs stored as sky map. More...
 
double m_energy_keV
 IRF line energy (optional) More...
 
double m_dlogE
 Logarithmic energy step for IRF band. More...
 
double m_gamma
 Power-law spectral index for IRF band. More...
 
std::vector< GSkyDirm_spix
 SPI pointing direction. More...
 
std::vector< double > m_posang
 Position angle of Y axis (CEL, radians) More...
 
bool m_has_wcs
 Has WCS information. More...
 
double m_wcs_xmin
 Minimum X value (radians) More...
 
double m_wcs_ymin
 Minimum Y value (radians) More...
 
double m_wcs_xmax
 Maximum X value (radians) More...
 
double m_wcs_ymax
 Maximum Y value (radians) More...
 
double m_wcs_xbin
 X value bin size (radians) More...
 
double m_wcs_ybin
 Y value bin size (radians) More...
 
double m_wcs_xpix_max
 Maximum X pixel index. More...
 
double m_wcs_ypix_max
 Maximum Y pixel index. More...
 
double m_max_zenith
 Maximum zenith angle (radians) More...
 

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_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

INTEGRAL/SPI instrument response function class.

The INTEGRAL/SPI instrument response function class defines the function that translates from physical quantities to measured quantities.

Todo:
Complete the class description.

Definition at line 63 of file GSPIResponse.hpp.

Constructor & Destructor Documentation

GSPIResponse::GSPIResponse ( void  )

Void constructor.

Creates an empty INTEGRAL/SPI response.

Definition at line 81 of file GSPIResponse.cpp.

References init_members().

Referenced by clone().

GSPIResponse::GSPIResponse ( const GSPIResponse rsp)

Copy constructor.

Parameters
[in]rspINTEGRAL/SPI response.

Definition at line 118 of file GSPIResponse.cpp.

References copy_members(), and init_members().

GSPIResponse::GSPIResponse ( const GFilename rspname)
explicit

Response constructor.

Constructs an INTEGRAL/SPI response for an observation using a response group filename.

This constructor simply stores the file name of a response group, the actual loading will be done using the set() method.

Definition at line 100 of file GSPIResponse.cpp.

References init_members(), m_rspname, and rspname().

GSPIResponse::~GSPIResponse ( void  )
virtual

Destructor.

Definition at line 134 of file GSPIResponse.cpp.

References free_members().

Member Function Documentation

double GSPIResponse::azimuth ( const int &  ipt,
const GSkyDir dir 
) const
inline

Return azimuth angle of sky direction for pointing in radians.

Parameters
[in]iptPointing index.
[in]dirSky direction.
Returns
Azimuth angle (radians).

Returns azimuth angle of sky direction for pointing in radians.

Definition at line 309 of file GSPIResponse.hpp.

References m_posang, m_spix, and gammalib::twopi.

Referenced by irf_value().

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

Return class name.

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

Implements GResponse.

Definition at line 163 of file GSPIResponse.hpp.

void GSPIResponse::clear ( void  )
virtual

Clear instance.

Clears INTEGRAL/SPI response by resetting all class members to an initial state. Any information that was present before will be lost.

Implements GResponse.

Definition at line 197 of file GSPIResponse.cpp.

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

Referenced by GSPIObservation::init_members().

GSPIResponse * GSPIResponse::clone ( void  ) const
virtual

Clone instance.

Returns
Pointer to deep copy of INTEGRAL/SPI response.

Implements GResponse.

Definition at line 217 of file GSPIResponse.cpp.

References GSPIResponse().

GSkyMap GSPIResponse::compute_irf ( const double &  emin,
const double &  emax 
) const
private

Compute as sky map.

Parameters
[in]eminMinimum energy (keV).
[in]emaxMaximum energy (keV).

Compute IRF for an energy band specified by an energy interval. If the width of the energy interval is zero a line IRF will be computed.

Definition at line 1216 of file GSPIResponse.cpp.

References exp(), GNodeArray::inx_left(), GNodeArray::inx_right(), irf(), irf_weight(), gammalib::ln10, log10(), m_dlogE, m_energies, m_gamma, m_irfs, GSkyMap::nmaps(), GSkyMap::npix(), GNodeArray::set_value(), GSkyMap::shape(), GNodeArray::wgt_left(), and GNodeArray::wgt_right().

Referenced by set().

void GSPIResponse::copy_members ( const GSPIResponse rsp)
private

Copy class members.

Parameters
[in]rspINTEGRAL/SPI response function.

Definition at line 761 of file GSPIResponse.cpp.

References m_detids, m_dlogE, m_ebounds, m_energies, m_energy_keV, m_gamma, m_has_wcs, m_irfs, m_max_zenith, m_posang, m_rspname, m_spix, m_wcs_xbin, m_wcs_xmax, m_wcs_xmin, m_wcs_xpix_max, m_wcs_ybin, m_wcs_ymax, m_wcs_ymin, and m_wcs_ypix_max.

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

const double & GSPIResponse::dlogE ( void  ) const
inline

Return logarithmic step size for continuum IRFs.

Returns
Logarithmic step size for continuum IRFs.

Returns the logarithmic step size for the computation of continuum IRFs.

Definition at line 263 of file GSPIResponse.hpp.

References m_dlogE.

GEbounds GSPIResponse::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.
Todo:
Implement this method if you need energy dispersion.

Implements GResponse.

Definition at line 442 of file GSPIResponse.cpp.

References G_EBOUNDS.

const double & GSPIResponse::energy_keV ( void  ) const
inline

Return line IRF energy in keV.

Returns
Line IRF energy (keV).

Returns the energy in keV for a line IRF. If the IRF is a continuum IRF the method returns 0.

Definition at line 249 of file GSPIResponse.hpp.

References m_energy_keV.

Referenced by GSPIObservation::write().

void GSPIResponse::free_members ( void  )
private

Delete class members.

Definition at line 795 of file GSPIResponse.cpp.

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

const double & GSPIResponse::gamma ( void  ) const
inline

Return power-law index for continuum IRFs.

Returns
Power-law index for continuum IRFs.

Returns the power-law index for the computation of continuum IRFs.

Definition at line 277 of file GSPIResponse.hpp.

References m_gamma.

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

Return value of INTEGRAL/SPI instrument response for a photon.

Parameters
[in]eventObserved event.
[in]photonIncident photon.
[in]obsObservation.
Returns
Instrument response \((cm^2 sr^{-1})\)
Exceptions
GException::invalid_argumentObservation is not a INTEGRAL/SPI observation.

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

\[ R(p'|p) = \]

Todo:

Write down formula

Describe in detail how the response is computed.

Implements GResponse.

Definition at line 245 of file GSPIResponse.cpp.

References GPhoton::dir(), GEvents::ebounds(), GPhoton::energy(), GObservation::events(), G_IRF, GSPIEventBin::iebin(), irf_value(), gammalib::is_infinite(), gammalib::is_notanumber(), and GSPIEventBin::livetime().

Referenced by compute_irf(), irf_value(), load_irf(), load_irfs(), and set().

int GSPIResponse::irf_detid ( const int &  detid) const
private

Convert detector identifier into IRF detector identifier.

Parameters
[in]detidSPI event detector identifier.
Returns
IRF detector identifier

Converts a SPI event detector identifier into an IRF detector identifier. TODO: Describe how this is done and why

Definition at line 1485 of file GSPIResponse.cpp.

Referenced by irf_value(), and set_detids().

double GSPIResponse::irf_value ( const GSkyDir srcDir,
const GSPIEventBin bin,
const int &  ireg 
) const

Return value of INTEGRAL/SPI instrument response for sky direction and event bin.

Parameters
[in]srcDirSky direction.
[in]binINTEGRAL/SPI event bin.
[in]iregIRF region (0: photo peak).
Returns
Instrument response \((cm^2 sr^{-1})\)

Returns the instrument response function for a given sky direction and event bin. The value of the IRF is bilinearly interpolated from the pre-computed IRFs cube that is stored in the class.

Definition at line 325 of file GSPIResponse.cpp.

References azimuth(), cos(), GSPIInstDir::detid(), GSPIEventBin::dir(), GSPIEventBin::iebin(), GSPIEventBin::ipt(), irf(), irf_detid(), m_irfs, m_max_zenith, m_wcs_xbin, m_wcs_xmin, m_wcs_xpix_max, m_wcs_ybin, m_wcs_ymin, m_wcs_ypix_max, GSkyMap::nx(), GSkyMap::shape(), sin(), and zenith().

Referenced by irf().

double GSPIResponse::irf_weight ( const double &  beta,
const double &  emin,
const double &  emax 
) const
private

Compute weight of logarithmic energy bin.

Parameters
[in]beta1-gamma.
[in]eminLogarithmic minimum energy.
[in]emaxLogarithmic maximum energy.

Definition at line 1513 of file GSPIResponse.cpp.

References abs(), exp(), and gammalib::ln10.

Referenced by compute_irf().

bool GSPIResponse::is_precomputed ( void  ) const
inline

Signals if response is precomputed.

Returns
True if response is precomputed.

Signals if the response was precomputed.

Definition at line 234 of file GSPIResponse.hpp.

References GEbounds::is_empty(), and m_ebounds.

void GSPIResponse::load ( const GFilename filename)

Load SPI response from file.

Parameters
[in]filenameResponse file name.

Loads SPI response from response file.

Definition at line 610 of file GSPIResponse.cpp.

References GFits::close(), and read().

Referenced by GSPIObservation::read().

GSkyMap GSPIResponse::load_irf ( const GFilename irfname) const
private

Load IRF as sky map.

Parameters
[in]irfnameIRF file name.

Loads an IRF FITS file as a sky map. The sky map is returned in ARC projection that is a zenith equidistant projection, which is the projection that is natively used to store the IRFs.

Definition at line 1003 of file GSPIResponse.cpp.

References abs(), GFits::close(), gammalib::deg2rad, G_LOAD_IRF, GFits::image(), GFitsHDU::integer(), irf(), m_has_wcs, m_max_zenith, m_wcs_xbin, m_wcs_xmax, m_wcs_xmin, m_wcs_xpix_max, m_wcs_ybin, m_wcs_ymax, m_wcs_ymin, m_wcs_ypix_max, GFitsImage::pixel(), GFitsHDU::real(), GSkyMap::shape(), and GFilename::url().

Referenced by load_irfs().

void GSPIResponse::load_irfs ( const int &  region)
private

Load Instrument Response Functions.

Parameters
[in]regionIRF region (-1 = load all regions)

The method requires that the required detector identifiers were previously setup using the set_detids() method.

Definition at line 1110 of file GSPIResponse.cpp.

References GNodeArray::append(), GNodeArray::clear(), GSkyMap::clear(), GFits::close(), G_LOAD_IRFS, irf(), load_irf(), m_detids, m_energies, m_irfs, m_rspname, GSkyMap::nmaps(), GSkyMap::nx(), GSkyMap::ny(), GFilename::path(), GSkyMap::shape(), gammalib::spi_hdu(), gammalib::spi_num_hdus(), and GFilename::url().

Referenced by set().

double GSPIResponse::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 416 of file GSPIResponse.cpp.

References G_NROI.

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

Assignment operator.

Parameters
[in]rspINTEGRAL/SPI response.
Returns
INTEGRAL/SPI response.

Assign INTEGRAL/SPI response to this object. The assignment performs a deep copy of all information, hence the original object from which the assignment was made can be destroyed after this operation without any loss of information.

Definition at line 161 of file GSPIResponse.cpp.

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

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

Print INTEGRAL/SPI response information.

Parameters
[in]chatterChattiness.
Returns
String containing INTEGRAL/SPI response information.

Implements GResponse.

Definition at line 659 of file GSPIResponse.cpp.

References m_dlogE, m_energy_keV, m_gamma, m_irfs, m_max_zenith, m_rspname, m_wcs_xmax, m_wcs_xmin, m_wcs_ymax, m_wcs_ymin, GSkyMap::nx(), GSkyMap::ny(), gammalib::parformat(), gammalib::rad2deg, GSkyMap::shape(), SILENT, gammalib::str(), and GFilename::url().

Referenced by GSPIObservation::print().

void GSPIResponse::read ( const GFits fits)

Read SPI response from FITS object.

Parameters
[in]fitsFITS object.

Reads the SPI response from FITS object.

Definition at line 558 of file GSPIResponse.cpp.

References GFits::image(), m_irfs, GSkyMap::read(), read_detids(), read_energies(), and set_wcs().

Referenced by load().

void GSPIResponse::read_detids ( const GFits fits)
private

Read detector identifiers from FITS object.

Parameters
[in]fitsFITS object.

Read the detector identifiers from a FITS file into the m_detids member.

Definition at line 809 of file GSPIResponse.cpp.

References GFitsTableCol::integer(), m_detids, GFitsTable::nrows(), and GFits::table().

Referenced by read().

void GSPIResponse::read_energies ( const GFits fits)
private

Read energies from FITS object.

Parameters
[in]fitsFITS object.

Read the energies from a FITS file. If the FITS file contains energy boundaries then read them, otherwise read the energy vector.

Definition at line 842 of file GSPIResponse.cpp.

References GNodeArray::append(), GNodeArray::clear(), GEbounds::clear(), GFits::contains(), gammalib::extname_ebounds, gammalib::extname_energies, GFitsHDU::has_card(), GEnergy::keV(), m_dlogE, m_ebounds, m_energies, m_energy_keV, m_gamma, GFitsTable::nrows(), GEbounds::read(), GFitsTableCol::real(), GFitsHDU::real(), GNodeArray::reserve(), GFitsHDU::string(), and GFits::table().

Referenced by read().

void GSPIResponse::rspname ( const GFilename rspname)
inline

Set response name.

Parameters
[in]rspnameResponse group file name.

Sets the response group file name.

Definition at line 205 of file GSPIResponse.hpp.

References m_rspname, and rspname().

Referenced by GSPIObservation::read().

const GFilename & GSPIResponse::rspname ( void  ) const
inline

Get response group file name.

Returns
Response group file name.

Returns the response group file name.

Definition at line 220 of file GSPIResponse.hpp.

References m_rspname.

Referenced by GSPIResponse(), and rspname().

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

Save SPI response into file.

Parameters
[in]filenameResponse file name.
[in]clobberOverwrite existing FITS file?

Saves SPI response into response file.

Definition at line 634 of file GSPIResponse.cpp.

References GFits::close(), GFits::saveto(), and write().

void GSPIResponse::set ( const GSPIObservation obs,
const GEnergy energy = GEnergy() 
)

Set response for a specific observation.

Parameters
[in]obsINTEGRAL/SPI observation.
[in]energyLine energy

Set response for a specific INTEGRAL/SPI observation.

If the energy argument is set to a positive value, the IRF will be computed for the specified line energy.

Definition at line 467 of file GSPIResponse.cpp.

References GEbounds::append(), GNodeArray::clear(), GEbounds::clear(), GSkyMap::clear(), compute_irf(), GEvents::ebounds(), GObservation::events(), irf(), GEnergy::keV(), load_irfs(), m_detids, m_ebounds, m_energies, m_energy_keV, m_has_wcs, m_irfs, GSPIEventCube::naxis(), GSkyMap::nmaps(), GSkyMap::npix(), set_cache(), set_detids(), and GSkyMap::shape().

Referenced by GSPIObservation::read().

void GSPIResponse::set_cache ( const GSPIEventCube cube)
private

Set computation cache.

Parameters
[in]cubeINTEGRAL/SPI event cube.

Setup of two vectors for fast coordinate transformation into the instrument system. The first vector m_spix stores the SPI pointing direction (the X direction) while the second vector stores the position angle in celestial coordinates of the SPI Y direction.

Definition at line 1448 of file GSPIResponse.cpp.

References m_posang, m_spix, GSPIEventCube::naxis(), gammalib::pihalf, GSkyDir::posang(), GSPIEventCube::spi_x(), and GSPIEventCube::spi_z().

Referenced by set().

void GSPIResponse::set_detids ( const GSPIEventCube cube)
private

Set vector of detector identifiers used by the observation.

Parameters
[in]cubeINTEGRAL/SPI event cube.

Sets the vector of detector identifiers that is used by the observation. The method scans the detector identifiers for all pointings and builds a vector of unique detector identifiers in the order they were encountered in the event cube.

Definition at line 1408 of file GSPIResponse.cpp.

References GSPIInstDir::detid(), GSPIEventCube::dir(), irf_detid(), m_detids, and GSPIEventCube::naxis().

Referenced by set().

void GSPIResponse::set_wcs ( const GFitsImage image)
private

Set IRF image limits.

Parameters
[in]imageIRF FITS image.

Computes the IRF image limits.

Definition at line 1362 of file GSPIResponse.cpp.

References abs(), gammalib::deg2rad, GFitsHDU::integer(), m_has_wcs, m_max_zenith, m_wcs_xbin, m_wcs_xmax, m_wcs_xmin, m_wcs_xpix_max, m_wcs_ybin, m_wcs_ymax, m_wcs_ymin, m_wcs_ypix_max, and GFitsHDU::real().

Referenced by read().

bool GSPIResponse::use_edisp ( void  ) const
inlinevirtual

Signal if energy dispersion will be used.

Returns
False.
Todo:
Implement method as needed.

Implements GResponse.

Definition at line 177 of file GSPIResponse.hpp.

bool GSPIResponse::use_tdisp ( void  ) const
inlinevirtual

Signal if time dispersion will be used.

Returns
False.
Todo:
Implement method as needed.

Implements GResponse.

Definition at line 191 of file GSPIResponse.hpp.

void GSPIResponse::write ( GFits fits) const

Write SPI response into FITS object.

Parameters
[in,out]fitsFITS object.

Writes the SPI response into FITS object.

Definition at line 587 of file GSPIResponse.cpp.

References m_irfs, GSkyMap::write(), write_detids(), and write_energies().

Referenced by save().

void GSPIResponse::write_detids ( GFits fits) const
private

Write detector identifiers into FITS object.

Parameters
[in,out]fitsFITS object.

Writes the detector identifiers stored into the m_detids member into a FITS object.

Definition at line 911 of file GSPIResponse.cpp.

References GFitsTable::append(), GFits::append(), GFitsHDU::extname(), and m_detids.

Referenced by write().

void GSPIResponse::write_energies ( GFits fits) const
private

Write energies into FITS object.

Parameters
[in,out]fitsFITS object.

Writes the energy nodes stored into the m_energies member into a FITS object. The energies are written in units of keV.

Definition at line 945 of file GSPIResponse.cpp.

References GFitsTable::append(), GFits::append(), GFitsHDU::card(), GFitsHDU::extname(), gammalib::extname_ebounds, gammalib::extname_energies, GEbounds::is_empty(), m_dlogE, m_ebounds, m_energies, m_energy_keV, m_gamma, GNodeArray::size(), GFits::table(), GFitsTableCol::unit(), and GEbounds::write().

Referenced by write().

double GSPIResponse::zenith ( const int &  ipt,
const GSkyDir dir 
) const
inline

Return zenith angle of sky direction for pointing in radians.

Parameters
[in]iptPointing index.
[in]dirSky direction.
Returns
Zenith angle (radians).

Returns zenith angle of sky direction for pointing in radians.

Definition at line 293 of file GSPIResponse.hpp.

References m_spix.

Referenced by irf_value().

Member Data Documentation

std::vector<int> GSPIResponse::m_detids
private

Vector of detector IDs.

Definition at line 133 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), load_irfs(), read_detids(), set(), set_detids(), and write_detids().

double GSPIResponse::m_dlogE
private

Logarithmic energy step for IRF band.

Definition at line 138 of file GSPIResponse.hpp.

Referenced by compute_irf(), copy_members(), dlogE(), init_members(), print(), read_energies(), and write_energies().

GEbounds GSPIResponse::m_ebounds
private

Energy bounaries of IRF.

Definition at line 135 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), is_precomputed(), read_energies(), set(), and write_energies().

GNodeArray GSPIResponse::m_energies
private

Node array of IRF energies.

Definition at line 134 of file GSPIResponse.hpp.

Referenced by compute_irf(), copy_members(), init_members(), load_irfs(), read_energies(), set(), and write_energies().

double GSPIResponse::m_energy_keV
private

IRF line energy (optional)

Definition at line 137 of file GSPIResponse.hpp.

Referenced by copy_members(), energy_keV(), init_members(), print(), read_energies(), set(), and write_energies().

double GSPIResponse::m_gamma
private

Power-law spectral index for IRF band.

Definition at line 139 of file GSPIResponse.hpp.

Referenced by compute_irf(), copy_members(), gamma(), init_members(), print(), read_energies(), and write_energies().

bool GSPIResponse::m_has_wcs
mutableprivate

Has WCS information.

Definition at line 144 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), load_irf(), set(), and set_wcs().

GSkyMap GSPIResponse::m_irfs
private

IRFs stored as sky map.

Definition at line 136 of file GSPIResponse.hpp.

Referenced by compute_irf(), copy_members(), init_members(), irf_value(), load_irfs(), print(), read(), set(), and write().

double GSPIResponse::m_max_zenith
mutableprivate

Maximum zenith angle (radians)

Definition at line 153 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), print(), and set_wcs().

std::vector<double> GSPIResponse::m_posang
private

Position angle of Y axis (CEL, radians)

Definition at line 143 of file GSPIResponse.hpp.

Referenced by azimuth(), copy_members(), init_members(), and set_cache().

GFilename GSPIResponse::m_rspname
private

File name of response group.

Definition at line 132 of file GSPIResponse.hpp.

Referenced by copy_members(), GSPIResponse(), init_members(), load_irfs(), print(), and rspname().

std::vector<GSkyDir> GSPIResponse::m_spix
private

SPI pointing direction.

Definition at line 142 of file GSPIResponse.hpp.

Referenced by azimuth(), copy_members(), init_members(), set_cache(), and zenith().

double GSPIResponse::m_wcs_xbin
mutableprivate

X value bin size (radians)

Definition at line 149 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), and set_wcs().

double GSPIResponse::m_wcs_xmax
mutableprivate

Maximum X value (radians)

Definition at line 147 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), load_irf(), print(), and set_wcs().

double GSPIResponse::m_wcs_xmin
mutableprivate

Minimum X value (radians)

Definition at line 145 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), print(), and set_wcs().

double GSPIResponse::m_wcs_xpix_max
mutableprivate

Maximum X pixel index.

Definition at line 151 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), and set_wcs().

double GSPIResponse::m_wcs_ybin
mutableprivate

Y value bin size (radians)

Definition at line 150 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), and set_wcs().

double GSPIResponse::m_wcs_ymax
mutableprivate

Maximum Y value (radians)

Definition at line 148 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), load_irf(), print(), and set_wcs().

double GSPIResponse::m_wcs_ymin
mutableprivate

Minimum Y value (radians)

Definition at line 146 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), print(), and set_wcs().

double GSPIResponse::m_wcs_ypix_max
mutableprivate

Maximum Y pixel index.

Definition at line 152 of file GSPIResponse.hpp.

Referenced by copy_members(), init_members(), irf_value(), load_irf(), and set_wcs().


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