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

COMPTEL Data Space class. More...

#include <GCOMDri.hpp>

Inheritance diagram for GCOMDri:
GBase

Public Member Functions

 GCOMDri (void)
 Void constructor. More...
 
 GCOMDri (const GFilename &filename)
 File name constructor. More...
 
 GCOMDri (const GSkyMap &map, const double &phimin=0.0, const double &phibin=0.0, const int &nphibin=0)
 Sky map constructor. More...
 
 GCOMDri (const GCOMDri &dri)
 Copy constructor. More...
 
virtual ~GCOMDri (void)
 Destructor. More...
 
GCOMDrioperator= (const GCOMDri &dri)
 Assignment operator. More...
 
double & operator[] (const int &index)
 DRI bin access operators. More...
 
const double & operator[] (const int &index) const
 DRI bin access operators (const version) More...
 
virtual void clear (void)
 Clear COMPTEL Data Space. More...
 
virtual GCOMDriclone (void) const
 Clone COMPTEL Data Space. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual std::string print (const GChatter &chatter=NORMAL) const
 Print COMPTEL Data Space. More...
 
int size (void) const
 Return number of bins. More...
 
int nchi (void) const
 Return number of Chi bins. More...
 
int npsi (void) const
 Return number of Psi bins. More...
 
int nphibar (void) const
 Return number of Phibar bins. More...
 
const GSkyMapmap (void) const
 Return DRI sky map. More...
 
const std::string & name (void) const
 Return DRI cube name. More...
 
void name (const std::string &name)
 Set DRI cube name. More...
 
const GEboundsebounds (void) const
 Return energy boundaries of DRI cube. More...
 
void ebounds (const GEbounds &ebounds)
 Set energy boundaries of DRI cube. More...
 
const GGtigti (void) const
 Return Good Time Intervals of DRI cube. More...
 
void gti (const GGti &gti)
 Set Good Time Intervals of DRI cube. More...
 
const double & phimin (void) const
 Return minimum Compton scatter angle of DRI cube. More...
 
const double & phibin (void) const
 Return Compton scatter angle bin of DRI cube. More...
 
const double & tof_correction (void) const
 Return ToF correction factor. More...
 
void tof_correction (const double &tofcor)
 Set ToF correction factor. More...
 
const double & phase_correction (void) const
 Return pulsar phase correction factor. More...
 
void phase_correction (const double &phasecor)
 Set pulsar phase correction factor. More...
 
const int & num_superpackets (void) const
 Return number of superpackets read for DRI. More...
 
void num_superpackets (const int &number)
 Set number of superpackets read for DRI. More...
 
const int & num_used_superpackets (void) const
 Return number of superpackets used for DRI. More...
 
void num_used_superpackets (const int &number)
 Set number of superpackets used for DRI. More...
 
const int & num_skipped_superpackets (void) const
 Return number of superpackets skipped for DRI. More...
 
void num_skipped_superpackets (const int &number)
 Set number of superpackets skipped for DRI. More...
 
void compute_dre (const GCOMObservation &obs, const GCOMSelection &select=GCOMSelection(), const double &zetamin=5.0)
 Compute event cube. More...
 
void compute_drg (const GCOMObservation &obs, const GCOMSelection &select=GCOMSelection(), const double &zetamin=5.0)
 Compute geometry cube. More...
 
void compute_drx (const GCOMObservation &obs, const GCOMSelection &select=GCOMSelection())
 Compute DRX exposure map. More...
 
void compute_drm (const GCOMObservation &obs, const GModel &model)
 Compute DRM model. More...
 
double cone_content (const GSkyDir &dir, const double &armmin, const double &armmax) const
 Compute content in cone. More...
 
void load (const GFilename &filename)
 Load COMPTEL Data Space from DRI FITS file. More...
 
void save (const GFilename &filename, const bool &clobber=false) const
 Save COMPTEL Data Space into DRI FITS file. More...
 
void read (const GFitsImage &image)
 Read COMPTEL Data Space from DRI FITS image. More...
 
void write (GFits &fits, const std::string &extname="") const
 Write COMPTEL Data Space into FITS image. More...
 
- Public Member Functions inherited from GBase
virtual ~GBase (void)
 Destructor. More...
 

Protected Member Functions

void init_members (void)
 Initialise class members. More...
 
void copy_members (const GCOMDri &dri)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void init_cube (void)
 Initialise DRI cube. More...
 
void init_statistics (void)
 Initialise computation statistics. More...
 
bool use_superpacket (const GCOMOad &oad, const GCOMTim &tim, const GCOMSelection &select)
 Check if superpacket should be used. More...
 
void read_attributes (const GFitsHDU *hdu)
 Read DRI attributes from FITS HDU. More...
 
void write_attributes (GFitsHDU *hdu) const
 Write DRI attributes into FITS HDU. More...
 
double compute_geometry (const int &tjd, const double &theta, const double &phi, const GCOMSelection &select, const GCOMStatus &status) const
 Compute DRG geometry factor. More...
 
double compute_surface (const double &x1, const double &y1, const double &r1, const double &x2, const double &y2, const double &r2) const
 Compute surface of overlap between two circles. More...
 
double compute_overlap (const double &x1, const double &y1, const double &r1, const double &x2, const double &y2, const double &r2, const double &x3, const double &y3, const double &r3) const
 Compute overlap between three circles. More...
 
void compute_tof_correction (void)
 Compute ToF correction. More...
 

Protected Attributes

std::string m_name
 Data cube name. More...
 
GSkyMap m_dri
 Data cube. More...
 
GEbounds m_ebounds
 Energy boundaries of data cube. More...
 
GGti m_gti
 Good Time Intervals of data cube. More...
 
double m_phimin
 Phibar minimum (deg) More...
 
double m_phibin
 Phibar binsize (deg) More...
 
double m_tofcor
 ToF correction. More...
 
double m_phasecor
 Pulsar phase correction. More...
 
GTime m_tstart
 Selection start time. More...
 
GTime m_tstop
 Selection stop time. More...
 
int m_num_superpackets
 Number of superpackets. More...
 
int m_num_used_superpackets
 Number of used superpackets. More...
 
int m_num_skipped_superpackets
 Number of skipped superpackets. More...
 
std::string m_drw_method
 DRW method. More...
 
GFitsBinTable m_drw_table
 DRW binary table to append to the FITS file. More...
 
std::string m_drw_status
 DRW fitting status. More...
 
double m_drw_fprompt
 DRW fitted fprompt parameter. More...
 
double m_drw_e_fprompt
 DRW fprompt parameter error. More...
 
int m_drw_iter
 DRW fitting iterations. More...
 
bool m_has_selection
 Signal that selection was applied. More...
 
GCOMSelection m_selection
 Selection parameters. More...
 
double m_zetamin
 Minimum zeta angle. More...
 

Friends

class GCOMDris
 

Detailed Description

COMPTEL Data Space class.

Definition at line 62 of file GCOMDri.hpp.

Constructor & Destructor Documentation

GCOMDri::GCOMDri ( void  )

Void constructor.

Definition at line 88 of file GCOMDri.cpp.

References init_members().

Referenced by clone().

GCOMDri::GCOMDri ( const GFilename filename)
explicit

File name constructor.

Parameters
[in]filenameCOMPTEL Data Space FITS file name.

Definition at line 103 of file GCOMDri.cpp.

References init_members(), and load().

GCOMDri::GCOMDri ( const GSkyMap map,
const double &  phimin = 0.0,
const double &  phibin = 0.0,
const int &  nphibin = 0 
)

Sky map constructor.

Parameters
[in]mapSky map defining the DRI cube.
[in]phiminMinimum Phibar angle (deg).
[in]phibinBin size of Phibar angle (deg).
[in]nphibinNumber of Phibar bins.

Constructs a DRI cube from a sky map and a Phibar binning definition.

Definition at line 126 of file GCOMDri.cpp.

References init_members(), m_dri, m_phibin, m_phimin, map(), GSkyMap::nmaps(), phibin(), and phimin().

GCOMDri::GCOMDri ( const GCOMDri dri)

Copy constructor.

Parameters
[in]driCOMPTEL Data Space.

Definition at line 159 of file GCOMDri.cpp.

References copy_members(), and init_members().

GCOMDri::~GCOMDri ( void  )
virtual

Destructor.

Definition at line 175 of file GCOMDri.cpp.

References free_members().

Member Function Documentation

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

Return class name.

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

Implements GBase.

Definition at line 192 of file GCOMDri.hpp.

void GCOMDri::clear ( void  )
virtual

Clear COMPTEL Data Space.

Implements GBase.

Definition at line 227 of file GCOMDri.cpp.

References free_members(), and init_members().

Referenced by GCOMModelDRM::init_members(), GCOMEventCube::init_members(), GCOMObservation::init_members(), and read().

GCOMDri * GCOMDri::clone ( void  ) const
virtual

Clone COMPTEL Data Space.

Returns
Pointer to deep copy of COMPTEL Data Space.

Implements GBase.

Definition at line 245 of file GCOMDri.cpp.

References GCOMDri().

void GCOMDri::compute_dre ( const GCOMObservation obs,
const GCOMSelection select = GCOMSelection(),
const double &  zetamin = 5.0 
)

Compute event cube.

Parameters
[in]obsCOMPTEL observation.
[in]selectSelection set.
[in]zetaminMinimum Earth horizon - Phibar cut (deg).
Exceptions
GException::invalid_argumentDRE cube has a non-positive Phibar bin size.
GException::invalid_valueNo BVC data available for pulsar selection

Compute DRE event cube for a COMPTEL observation.

Definition at line 265 of file GCOMDri.cpp.

References abs(), GCOMObservation::bvcs(), compute_tof_correction(), GPhases::contains(), GSkyMap::contains(), GPulsarEphemeris::dir(), GSkyMap::dir2pix(), GSkyDir::dist_deg(), GEbounds::emax(), GEbounds::emin(), GPulsar::ephemeris(), GObservation::events(), G_COMPUTE_DRE, GCOMOad::gcaz(), GCOMOad::gcel(), GEphemerides::geo2ssb(), GCOMOad::georad(), GCOMSelection::has_pulsar(), init_cube(), GCOMSelection::init_statistics(), init_statistics(), GCOMBvcs::is_empty(), GPhases::length(), m_dri, m_ebounds, m_gti, m_has_selection, m_num_skipped_superpackets, m_num_superpackets, m_num_used_superpackets, m_phasecor, m_phibin, m_phimin, m_selection, m_tstart, m_tstop, m_zetamin, nphibar(), GSkyMap::npix(), GCOMObservation::oads(), GPulsarEphemeris::phase(), GSkyMap::pix2inx(), GCOMOad::pos(), GCOMSelection::pulsar(), GCOMSelection::pulsar_phases(), GSkyDir::radec_deg(), GCOMTim::reduce(), GGti::reduce(), GCOMEventList::size(), GCOMOads::size(), GGti::size(), gammalib::str(), GCOMBvcs::tdelta(), GCOMObservation::tim(), GCOMOad::tstart(), GGti::tstart(), GCOMOad::tstop(), GGti::tstop(), GCOMSelection::use_event(), use_superpacket(), GPulsar::validity(), and gammalib::warning().

void GCOMDri::compute_drg ( const GCOMObservation obs,
const GCOMSelection select = GCOMSelection(),
const double &  zetamin = 5.0 
)
void GCOMDri::compute_drm ( const GCOMObservation obs,
const GModel model 
)

Compute DRM model.

Parameters
[in]obsCOMPTEL observation.
[in]modelModel.

Compute DRM model cube for a COMPTEL observation.

Definition at line 809 of file GCOMDri.cpp.

References GModel::eval(), m_dri, nchi(), nphibar(), and npsi().

void GCOMDri::compute_drx ( const GCOMObservation obs,
const GCOMSelection select = GCOMSelection() 
)

Compute DRX exposure map.

Parameters
[in]obsCOMPTEL observation.
[in]selectSelection set.
Exceptions
GException::invalid_valueNo BVC data available for pulsar selection

Compute DRX exposure map for a COMPTEL observation.

For a given superpacket, the exposure is computed using

\[ X_i(\theta_c) = 7 \pi r_1^2 \cos \theta_c \frac{1 - \exp \left( -\tau \ \cos \theta_c \right)} {1 - \exp \left( -\tau \right)} \]

where \(\tau=0.2\) is the typical thickness of a D1 module in radiation lengths, \(r_1=13.8\) cm is the radius of a D1 module, and \(\theta_c\) is the zenith angle in COMPTEL coordinates.

Definition at line 702 of file GCOMDri.cpp.

References GEbounds::clear(), cos(), d1_area, gammalib::deg2rad, exp(), GCOMSelection::has_pulsar(), init_cube(), init_statistics(), GSkyMap::inx2dir(), m_dri, m_ebounds, m_gti, m_num_skipped_superpackets, m_num_superpackets, m_num_used_superpackets, m_tstart, m_tstop, GSkyMap::npix(), GCOMObservation::oads(), GCOMSelection::pulsar(), GCOMTim::reduce(), GCOMOads::size(), size(), superpacket_duration, GCOMOad::theta(), GCOMObservation::tim(), use_superpacket(), and GPulsar::validity().

double GCOMDri::compute_geometry ( const int &  tjd,
const double &  theta,
const double &  phi,
const GCOMSelection select,
const GCOMStatus status 
) const
protected

Compute DRG geometry factor.

Parameters
[in]tjdTJD for module status
[in]thetaZenith angle in COMPTEL coordinates (deg).
[in]phiAzimuth angle in COMPTEL coordinates (deg).
[in]selectSelection set.
[in]statusD1 and D2 module status
Returns
Geometry factor.

Computes the DRG geometry factor as function of zenith and azimuth angles given in the COMPTEL coordinate system.

Definition at line 1506 of file GCOMDri.cpp.

References gammalib::acos(), gammalib::com_exd2r(), gammalib::com_exd2x(), gammalib::com_exd2y(), compute_overlap(), cos(), GCOMStatus::d1status(), GCOMStatus::d2status(), gammalib::deg2rad, GCOMSelection::fpmtflag(), gammalib::pi, r1, r1sq, sin(), sqrt(), tan(), GCOMSelection::use_d1(), and GCOMSelection::use_d2().

Referenced by compute_drg(), GCOMDris::compute_drws_energy(), GCOMDris::compute_drws_phibar(), GCOMDris::vetorate_generate(), and GCOMDris::vetorate_setup().

double GCOMDri::compute_overlap ( const double &  x1,
const double &  y1,
const double &  r1,
const double &  x2,
const double &  y2,
const double &  r2,
const double &  x3,
const double &  y3,
const double &  r3 
) const
protected

Compute overlap between three circles.

Parameters
[in]x1X position of D1 projection (cm).
[in]y1Y position of D1 projection (cm).
[in]r1Radius D1 module (cm).
[in]x2X position of D2 module (cm).
[in]y2Y position of D2 module (cm).
[in]r2Radius D2 module (cm).
[in]x3X position of dead PMT (cm).
[in]y3Y position of dead PMT (cm).
[in]r3Radius of dead PMT (cm).
Returns
Fractional overlap [0.0,...,1.0]

Compute fractional overlap between three circles, composed of projected D1 module, D2 module and failed PMT exclusion circle.

The method is a reimplementation of the COMPASS SKYDRS17.OVERLP function.

Definition at line 1747 of file GCOMDri.cpp.

References compute_surface(), gammalib::pi, r1, r1sq, size(), and sqrt().

Referenced by compute_geometry().

double GCOMDri::compute_surface ( const double &  x1,
const double &  y1,
const double &  r1,
const double &  x2,
const double &  y2,
const double &  r2 
) const
protected

Compute surface of overlap between two circles.

Parameters
[in]x1X position of D2 module (cm).
[in]y1Y position of D2 module (cm).
[in]r1Radius D2 module (cm).
[in]x2X position of dead PMT (cm).
[in]y2Y position of dead PMT (cm).
[in]r2Radius of dead PMT (cm).
Returns
Surface of overlap (cm^2)

Computes the surface of overlap in cm^2 between two circles, composed of D2 module and failed PMT exclusion circle.

The method is a reimplementation of the COMPASS SKYDRS17.COM2 function.

Definition at line 1673 of file GCOMDri.cpp.

References gammalib::acos(), gammalib::pi, r1, r1sq, and sqrt().

Referenced by compute_overlap().

void GCOMDri::compute_tof_correction ( void  )
protected

Compute ToF correction.

Compute the ToF correction according to COM-RP-ROL-DRG-057.

Definition at line 1863 of file GCOMDri.cpp.

References GEbounds::emax(), GEbounds::emin(), G_COMPUTE_TOF_CORRECTION, m_ebounds, m_selection, m_tofcor, GEnergy::MeV(), sqrt(), gammalib::str(), GCOMSelection::tof_max(), GCOMSelection::tof_min(), and gammalib::warning().

Referenced by compute_dre().

double GCOMDri::cone_content ( const GSkyDir dir,
const double &  armmin,
const double &  armmax 
) const

Compute content in cone.

Parameters
[in]dirSky direction of cone apex.
[in]armminMinimum Angular Resolution Measure (deg).
[in]armmaxMaximum Angular Resolution Measure (deg).
Returns
Content in cone.

Compute the sum of the DRI bins within an event cone with apex at a given sky direction. All bins with an Angular Resolution Measure comprised between armmin (inclusive) and armmax (exclusive) will be considered. The bin centres will be used for the computation of the Angular Resolution Measure. The Angular Resolution Measure is defined as phibar - phigeo.

Definition at line 845 of file GCOMDri.cpp.

References GSkyDir::dist_deg(), GSkyMap::extract(), GSkyMap::inx2dir(), m_dri, nphibar(), GSkyMap::npix(), phibin(), and phimin().

void GCOMDri::copy_members ( const GCOMDri dri)
protected
const GEbounds & GCOMDri::ebounds ( void  ) const
inline

Return energy boundaries of DRI cube.

Returns
Energy boundaries of DRI cube.

Definition at line 317 of file GCOMDri.hpp.

References m_ebounds.

Referenced by ebounds(), and GCOMEventCube::init_cube().

void GCOMDri::ebounds ( const GEbounds ebounds)
inline

Set energy boundaries of DRI cube.

Parameters
[in]eboundsEnergy boundaries of DRI cube.

Sets energy boundaries of DRI cube

Definition at line 331 of file GCOMDri.hpp.

References ebounds(), and m_ebounds.

void GCOMDri::free_members ( void  )
protected

Delete class members.

Definition at line 1171 of file GCOMDri.cpp.

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

const GGti & GCOMDri::gti ( void  ) const
inline

Return Good Time Intervals of DRI cube.

Returns
Good Time Intervals of DRI cube.

Definition at line 344 of file GCOMDri.hpp.

References m_gti.

Referenced by gti(), GCOMEventCube::init_cube(), and GCOMDris::vetorate_setup().

void GCOMDri::gti ( const GGti gti)
inline

Set Good Time Intervals of DRI cube.

Parameters
[in]gtiGood Time Intervals of DRI data.

Sets the Good Time Intervals of DRI cube.

Definition at line 358 of file GCOMDri.hpp.

References gti(), and m_gti.

void GCOMDri::init_cube ( void  )
protected

Initialise DRI cube.

Sets all DRI cube bins to zero.

Definition at line 1183 of file GCOMDri.cpp.

References size().

Referenced by compute_dre(), compute_drg(), and compute_drx().

void GCOMDri::init_statistics ( void  )
protected

Initialise computation statistics.

Definition at line 1198 of file GCOMDri.cpp.

References GTime::clear(), m_num_skipped_superpackets, m_num_superpackets, m_num_used_superpackets, m_tstart, and m_tstop.

Referenced by compute_dre(), compute_drg(), compute_drx(), and init_members().

void GCOMDri::load ( const GFilename filename)

Load COMPTEL Data Space from DRI FITS file.

Parameters
[in]filenameDRI FITS file name.

Definition at line 892 of file GCOMDri.cpp.

References GFits::close(), GFits::image(), m_drw_table, read(), and GFits::size().

Referenced by GCOMDri(), and GCOMModelDRM::read().

const std::string & GCOMDri::name ( void  ) const
inline

Return DRI cube name.

Returns
DRI cube name.

Definition at line 290 of file GCOMDri.hpp.

References m_name.

Referenced by name().

void GCOMDri::name ( const std::string &  name)
inline

Set DRI cube name.

Parameters
[in]nameDRI cube name.

Sets the name of the DRI cube.

Definition at line 304 of file GCOMDri.hpp.

References m_name, and name().

int GCOMDri::nchi ( void  ) const
inline
int GCOMDri::npsi ( void  ) const
inline
const int & GCOMDri::num_skipped_superpackets ( void  ) const
inline

Return number of superpackets skipped for DRI.

Returns
Number of superpackets skipped for DRI.

Returns the number of superpackets skipped for DRI.

Definition at line 517 of file GCOMDri.hpp.

References m_num_skipped_superpackets.

void GCOMDri::num_skipped_superpackets ( const int &  number)
inline

Set number of superpackets skipped for DRI.

Parameters
[in]numberNumber of superpackets skipped for DRI.

Set the number of superpackets skipped for DRI.

Definition at line 531 of file GCOMDri.hpp.

References m_num_skipped_superpackets, and gammalib::number().

const int & GCOMDri::num_superpackets ( void  ) const
inline

Return number of superpackets read for DRI.

Returns
Number of superpackets read for DRI.

Returns the number of superpackets read for DRI.

Definition at line 459 of file GCOMDri.hpp.

References m_num_superpackets.

void GCOMDri::num_superpackets ( const int &  number)
inline

Set number of superpackets read for DRI.

Parameters
[in]numberNumber of superpackets read for DRI.

Set the number of superpackets read for DRI.

Definition at line 473 of file GCOMDri.hpp.

References m_num_superpackets, and gammalib::number().

const int & GCOMDri::num_used_superpackets ( void  ) const
inline

Return number of superpackets used for DRI.

Returns
Number of superpackets used for DRI.

Returns the number of superpackets used for DRI.

Definition at line 488 of file GCOMDri.hpp.

References m_num_used_superpackets.

void GCOMDri::num_used_superpackets ( const int &  number)
inline

Set number of superpackets used for DRI.

Parameters
[in]numberNumber of superpackets used for DRI.

Set the number of superpackets used for DRI.

Definition at line 502 of file GCOMDri.hpp.

References m_num_used_superpackets, and gammalib::number().

GCOMDri & GCOMDri::operator= ( const GCOMDri dri)

Assignment operator.

Parameters
[in]driCOMPTEL Data Space.
Returns
COMPTEL Data Space.

Definition at line 197 of file GCOMDri.cpp.

References copy_members(), free_members(), and init_members().

double & GCOMDri::operator[] ( const int &  index)
inline

DRI bin access operators.

Parameters
[in]indexDRI bin index [0,...,size()-1].
Returns
Reference to DRI bin.

Definition at line 205 of file GCOMDri.hpp.

References m_dri, and GSkyMap::pixels().

const double & GCOMDri::operator[] ( const int &  index) const
inline

DRI bin access operators (const version)

Parameters
[in]indexDRI bin index [0,...,size()-1].
Returns
Reference to DRI bin.

Definition at line 218 of file GCOMDri.hpp.

References m_dri, and GSkyMap::pixels().

const double & GCOMDri::phase_correction ( void  ) const
inline

Return pulsar phase correction factor.

Returns
Pulsar phase correction factor.

Returns the pulsar phase correction factor that corrects for the phase selection for pulsar analysis.

Definition at line 429 of file GCOMDri.hpp.

References m_phasecor.

Referenced by GCOMResponse::irf(), GCOMResponse::irf_diffuse(), GCOMResponse::irf_elliptical(), GCOMResponse::irf_ptsrc(), GCOMResponse::irf_radial(), and GCOMEventCube::print().

void GCOMDri::phase_correction ( const double &  phasecor)
inline

Set pulsar phase correction factor.

Parameters
[in]phasecorPulsar phase correction factor.

Set the pulsar phase correction factor that corrects for the phase selection for pulsar analysis.

Definition at line 444 of file GCOMDri.hpp.

References m_phasecor.

const double & GCOMDri::phibin ( void  ) const
inline

Return Compton scatter angle bin of DRI cube.

Returns
Compton scatter angle bin of DRI cube (deg).

Definition at line 383 of file GCOMDri.hpp.

References m_phibin.

Referenced by GCOMDris::compute_drws_phibar(), cone_content(), GCOMDri(), and GCOMEventCube::set_scatter_angles().

const double & GCOMDri::phimin ( void  ) const
inline

Return minimum Compton scatter angle of DRI cube.

Returns
Minimum Compton scatter angle of DRI cube (deg).

Definition at line 371 of file GCOMDri.hpp.

References m_phimin.

Referenced by GCOMDris::compute_drws_phibar(), cone_content(), GCOMDri(), and GCOMEventCube::set_scatter_angles().

std::string GCOMDri::print ( const GChatter chatter = NORMAL) const
virtual
void GCOMDri::read ( const GFitsImage image)

Read COMPTEL Data Space from DRI FITS image.

Parameters
[in]imageDRI FITS image.

Definition at line 955 of file GCOMDri.cpp.

References clear(), gammalib::com_wcs_mer2car(), m_dri, GSkyMap::read(), and read_attributes().

Referenced by load(), GCOMObservation::load_drb(), GCOMObservation::load_drg(), GCOMObservation::load_drw(), GCOMObservation::load_drx(), and GCOMEventCube::read().

void GCOMDri::read_attributes ( const GFitsHDU hdu)
protected

Read DRI attributes from FITS HDU.

Parameters
[in]hduFITS HDU pointer.

Reads the time interval from the FITS header and sets the Phibar definiton and energy boundaries from the header keywords if they are provided.

Definition at line 1295 of file GCOMDri.cpp.

References GEbounds::clear(), gammalib::com_time(), GFitsHDU::has_card(), GFitsHDU::integer(), m_ebounds, m_gti, m_has_selection, m_num_skipped_superpackets, m_num_superpackets, m_num_used_superpackets, m_phasecor, m_phibin, m_phimin, m_selection, m_tofcor, m_zetamin, GCOMSelection::read(), and GFitsHDU::real().

Referenced by read().

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

Save COMPTEL Data Space into DRI FITS file.

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

Definition at line 926 of file GCOMDri.cpp.

References GFits::append(), GFits::close(), GFilename::extname(), gammalib::extname_dri, m_drw_table, GFitsTable::nrows(), GFits::saveto(), and write().

int GCOMDri::size ( void  ) const
inline

Return number of bins.

Returns
Number of bins.

Definition at line 230 of file GCOMDri.hpp.

References m_dri, GSkyMap::nmaps(), and GSkyMap::npix().

Referenced by compute_drg(), compute_drx(), compute_overlap(), init_cube(), GCOMEventCube::number(), and GCOMEventCube::size().

const double & GCOMDri::tof_correction ( void  ) const
inline

Return ToF correction factor.

Returns
ToF correction factor.

Returns the ToF correction factor that corrects for the event selection in a ToF window.

Definition at line 398 of file GCOMDri.hpp.

References m_tofcor.

Referenced by GCOMResponse::irf(), GCOMResponse::irf_diffuse(), GCOMResponse::irf_elliptical(), GCOMResponse::irf_ptsrc(), GCOMResponse::irf_radial(), and GCOMEventCube::print().

void GCOMDri::tof_correction ( const double &  tofcor)
inline

Set ToF correction factor.

Parameters
[in]tofcorToF correction factor.

Set the ToF correction factor that corrects for the event selection in a ToF window.

Definition at line 413 of file GCOMDri.hpp.

References m_tofcor.

bool GCOMDri::use_superpacket ( const GCOMOad oad,
const GCOMTim tim,
const GCOMSelection select 
)
protected

Check if superpacket should be used.

Parameters
[in]oadOrbit Aspect Data record (i.e. superpacket).
[in]timGood Time Intervals.
[in]selectSelection set.
Returns
True if superpacket should be used, false otherwise.

Checks if a superpacket should be used. A superpacket will be used if it is fully enclosed within the COMPTEL Good Time Intervals and the Good Time Intervals of the DRI dataset. In case that orbital phases are present in the selection set, the superpacket will be used when the start time is comprised in one of the orbital phases.

The method updates the superpacket statistics and selected time interval.

Definition at line 1228 of file GCOMDri.cpp.

References GPhases::contains(), GCOMTim::contains(), GGti::contains(), GPhases::is_empty(), m_gti, m_num_skipped_superpackets, m_num_superpackets, m_num_used_superpackets, m_tstart, m_tstop, GCOMSelection::orbital_phase(), GCOMSelection::orbital_phases(), GGti::size(), GCOMOad::tstart(), and GCOMOad::tstop().

Referenced by compute_dre(), compute_drg(), and compute_drx().

void GCOMDri::write ( GFits fits,
const std::string &  extname = "" 
) const

Write COMPTEL Data Space into FITS image.

Parameters
[in]fitsFITS file.
[in]extnameExtension name.

Definition at line 980 of file GCOMDri.cpp.

References m_dri, GSkyMap::write(), and write_attributes().

Referenced by save(), and GCOMEventCube::write().

Friends And Related Function Documentation

friend class GCOMDris
friend

Definition at line 64 of file GCOMDri.hpp.

Member Data Documentation

double GCOMDri::m_drw_e_fprompt
protected

DRW fprompt parameter error.

Definition at line 176 of file GCOMDri.hpp.

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

double GCOMDri::m_drw_fprompt
protected

DRW fitted fprompt parameter.

Definition at line 175 of file GCOMDri.hpp.

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

int GCOMDri::m_drw_iter
protected

DRW fitting iterations.

Definition at line 177 of file GCOMDri.hpp.

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

std::string GCOMDri::m_drw_method
protected

DRW method.

Definition at line 172 of file GCOMDri.hpp.

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

std::string GCOMDri::m_drw_status
protected

DRW fitting status.

Definition at line 174 of file GCOMDri.hpp.

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

GFitsBinTable GCOMDri::m_drw_table
protected

DRW binary table to append to the FITS file.

Definition at line 173 of file GCOMDri.hpp.

Referenced by copy_members(), init_members(), load(), save(), and GCOMDris::vetorate_finish().

GEbounds GCOMDri::m_ebounds
protected
GGti GCOMDri::m_gti
protected

Good Time Intervals of data cube.

Definition at line 158 of file GCOMDri.hpp.

Referenced by compute_dre(), compute_drg(), compute_drx(), copy_members(), gti(), init_members(), print(), read_attributes(), use_superpacket(), and write_attributes().

bool GCOMDri::m_has_selection
protected

Signal that selection was applied.

Definition at line 180 of file GCOMDri.hpp.

Referenced by compute_dre(), copy_members(), init_members(), print(), read_attributes(), and write_attributes().

std::string GCOMDri::m_name
protected

Data cube name.

Definition at line 155 of file GCOMDri.hpp.

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

int GCOMDri::m_num_skipped_superpackets
protected
int GCOMDri::m_num_superpackets
protected
int GCOMDri::m_num_used_superpackets
protected
double GCOMDri::m_phasecor
protected

Pulsar phase correction.

Definition at line 162 of file GCOMDri.hpp.

Referenced by compute_dre(), copy_members(), init_members(), phase_correction(), read_attributes(), and write_attributes().

double GCOMDri::m_phimin
protected
GCOMSelection GCOMDri::m_selection
protected
double GCOMDri::m_tofcor
protected
GTime GCOMDri::m_tstart
protected

Selection start time.

Definition at line 165 of file GCOMDri.hpp.

Referenced by compute_dre(), compute_drg(), GCOMDris::compute_drws(), compute_drx(), copy_members(), init_statistics(), and use_superpacket().

GTime GCOMDri::m_tstop
protected

Selection stop time.

Definition at line 166 of file GCOMDri.hpp.

Referenced by compute_dre(), compute_drg(), GCOMDris::compute_drws(), compute_drx(), copy_members(), init_statistics(), and use_superpacket().

double GCOMDri::m_zetamin
protected

Minimum zeta angle.

Definition at line 182 of file GCOMDri.hpp.

Referenced by compute_dre(), copy_members(), init_members(), read_attributes(), and write_attributes().


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