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

Interface class for COMPTEL observations. More...

#include <GCOMObservation.hpp>

Inheritance diagram for GCOMObservation:
GObservation GBase

Public Member Functions

 GCOMObservation (void)
 Void constructor. More...
 
 GCOMObservation (const GXmlElement &xml)
 XML constructor. More...
 
 GCOMObservation (const GCOMDri &dre, const GCOMDri &drb, const GCOMDri &drg, const GCOMDri &drx)
 Binned observation DRI constructor. More...
 
 GCOMObservation (const GCOMDri &dre, const GCOMDri &drb, const GCOMDri &drw, const GCOMDri &drg, const GCOMDri &drx)
 Binned observation DRI constructor. More...
 
 GCOMObservation (const GFilename &drename, const GFilename &drbname, const GFilename &drwname, const GFilename &drgname, const GFilename &drxname)
 Binned observation filename constructor. More...
 
 GCOMObservation (const GFilename &evpname, const GFilename &timname, const std::vector< GFilename > &oadnames, const std::vector< GFilename > &hkdnames=std::vector< GFilename >(), const GFilename &bvcname="")
 Unbinned observation constructor. More...
 
 GCOMObservation (const GCOMObservation &obs)
 Copy constructor. More...
 
virtual ~GCOMObservation (void)
 Destructor. More...
 
virtual GCOMObservationoperator= (const GCOMObservation &obs)
 Assignment operator. More...
 
virtual void clear (void)
 Clear COMPTEL observation. More...
 
virtual GCOMObservationclone (void) const
 Clone COMPTEL observation. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual void response (const GResponse &rsp)
 Set response function. More...
 
virtual const GCOMResponseresponse (void) const
 Return response function. More...
 
virtual std::string instrument (void) const
 Return instrument. More...
 
virtual double ontime (void) const
 Return ontime. More...
 
virtual double livetime (void) const
 Return livetime. More...
 
virtual double deadc (const GTime &time=GTime()) const
 Return deadtime correction factor. More...
 
virtual void read (const GXmlElement &xml)
 Read observation from XML element. More...
 
virtual void write (GXmlElement &xml) const
 Write observation into XML element. More...
 
virtual std::string print (const GChatter &chatter=NORMAL) const
 Print observation information. More...
 
virtual double npred (const GModel &model) const
 Return total number of predicted counts for one model. More...
 
bool is_unbinned (void) const
 Check whether observation is unbinned. More...
 
bool is_binned (void) const
 Check whether observation is binned. More...
 
void load (const GFilename &drename, const GFilename &drbname, const GFilename &drwname, const GFilename &drgname, const GFilename &drxname)
 Load data for a binned observation. More...
 
void load (const GFilename &evpname, const GFilename &timname, const std::vector< GFilename > &oadnames, const std::vector< GFilename > &hkdnames=std::vector< GFilename >(), const GFilename &bvcname="")
 Load data for an unbinned observation. More...
 
void response (const GCaldb &caldb, const std::string &rspname)
 Set response function. More...
 
void response (const GCOMResponse &response)
 Set response function. More...
 
void obs_id (const double &id)
 Set observation ID. More...
 
void ontime (const double &ontime)
 Set ontime. More...
 
void livetime (const double &livetime)
 Set livetime. More...
 
void deadc (const double &deadc)
 Set deadtime correction factor. More...
 
void ewidth (const double &ewidth)
 Set energy width. More...
 
const double & obs_id (void) const
 Return observation ID. More...
 
const double & ewidth (void) const
 Return energy width. More...
 
const GCOMDridrb (void) const
 Return background model. More...
 
const GCOMDridrw (void) const
 Return weighting cube. More...
 
const GCOMDridrg (void) const
 Return geometry factors. More...
 
const GCOMDridrx (void) const
 Return exposure. More...
 
GCOMDri drm (const GModels &models) const
 Compute DRM cube. More...
 
const GCOMTimtim (void) const
 Return COMPTEL Good Time Intervals. More...
 
void tim (const GCOMTim &tim)
 Set COMPTEL Good Time Intervals. More...
 
const GCOMOadsoads (void) const
 Return Orbit Aspect Data. More...
 
void oads (const GCOMOads &oads)
 Set Orbit Aspect Data. More...
 
const GCOMHkdshkds (void) const
 Return Housekeeping Data collection. More...
 
void hkds (const GCOMHkds &hkds)
 Set Housekeeping Data collection. More...
 
const GCOMBvcsbvcs (void) const
 Return Solar System Barycentre Data. More...
 
void bvcs (const GCOMBvcs &bvcs)
 Set Solar System Barycentre Data. More...
 
const GFilenamedrename (void) const
 Return DRE filename. More...
 
const GFilenamedrbname (void) const
 Return DRB filename. More...
 
const GFilenamedrwname (void) const
 Return DRW filename. More...
 
const GFilenamedrgname (void) const
 Return DRG filename. More...
 
const GFilenamedrxname (void) const
 Return DRX filename. More...
 
const GFilenamerspname (void) const
 Return response cache filename. More...
 
const int & phi_first (void) const
 Return index of first Phibar layer to be used for likelihood fitting. More...
 
const int & phi_last (void) const
 Return index of last Phibar layer to be used for likelihood fitting. More...
 
void drename (const GFilename &drename)
 Set DRE filename. More...
 
void drbname (const GFilename &drbname)
 Set DRB filename. More...
 
void drwname (const GFilename &drwname)
 Set DRW filename. More...
 
void drgname (const GFilename &drgname)
 Set DRG filename. More...
 
void drxname (const GFilename &drxname)
 Set DRX filename. More...
 
void rspname (const GFilename &rspname)
 Set response cache filename filename. More...
 
void phi_first (const int &phi_first)
 Set index of first Phibar layer to be used for likelihood fitting. More...
 
void phi_last (const int &phi_last)
 Set index of last Phibar layer to be used for likelihood fitting. More...
 
void compute_drb (const std::string &method, const GCOMDri &drm, const int &nrunav=3, const int &navgr=3, const int &nincl=15, const int &nexcl=0)
 Compute DRB cube. More...
 
- Public Member Functions inherited from GObservation
 GObservation (void)
 Void constructor. More...
 
 GObservation (const GObservation &obs)
 Copy constructor. More...
 
virtual ~GObservation (void)
 Destructor. More...
 
virtual GObservationoperator= (const GObservation &obs)
 Assignment operator. More...
 
virtual GEventsevents (void)
 Return events. More...
 
virtual const GEventsevents (void) const
 Return events (const version) More...
 
virtual void events (const GEvents &events)
 Set event container. More...
 
virtual double likelihood (const GModels &models, GVector *gradients, GMatrixSparse *curvature, double *npred) const
 Compute likelihood function. More...
 
virtual double model (const GModels &models, const GEvent &event, GVector *gradients=NULL) const
 Return model value and (optionally) gradients. More...
 
virtual GVector model (const GModels &models, GMatrixSparse *gradients=NULL) const
 Return vector of model values and (optionally) gradients. More...
 
virtual int nobserved (void) const
 Return total number of observed events. More...
 
virtual double npred (const GModels &models, GVector *gradients=NULL) const
 Return total number (and optionally gradients) of predicted counts for all models. More...
 
virtual double model_grad (const GModel &model, const GModelPar &par, const GEvent &event) const
 Returns parameter gradient of model for a given event. More...
 
virtual GVector model_grad (const GModel &model, const GModelPar &par) const
 Returns parameter gradients of model for all events. More...
 
virtual double npred_grad (const GModel &model, const GModelPar &par) const
 Returns parameter gradient of Npred. More...
 
virtual void remove_response_cache (const std::string &name)
 Response cache removal hook. More...
 
virtual const double & grad_step_size (void) const
 Return gradient step size. More...
 
bool has_events (void) const
 Signal if observation has events. More...
 
bool has_gradient (const GModel &model, const GModelPar &par) const
 Check whether a model parameter has an analytical gradient. More...
 
void name (const std::string &name)
 Set observation name. More...
 
void id (const std::string &id)
 Set observation identifier. More...
 
void statistic (const std::string &statistic)
 Set optimizer statistic. More...
 
const std::string & name (void) const
 Return observation name. More...
 
const std::string & id (void) const
 Return observation identifier. More...
 
const std::string & statistic (void) const
 Return optimizer statistic. More...
 
void computed_gradient (const GModel &model, const GModelPar &par) const
 Signals that an analytical gradient was computed for a model parameter. 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 GCOMObservation &obs)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void load_dre (const GFilename &drename)
 Load event cube data from DRE file. More...
 
void load_drb (const GFilename &drbname)
 Load background model from DRB file. More...
 
void load_drw (const GFilename &drwname)
 Load weighting cube from DRW file. More...
 
void load_drg (const GFilename &drgname)
 Load geometry factors from DRG file. More...
 
void load_drx (const GFilename &drxname)
 Load exposure from DRX file. More...
 
bool check_dri (const GCOMDri &map) const
 Check if DRI is compatible with event cube. More...
 
void read_attributes (const GFitsHDU *hdu)
 Read observation attributes. More...
 
void write_attributes (GFitsHDU *hdu) const
 Write observation attributes. More...
 
void compute_drb_phinor (const GCOMDri &drm)
 Compute DRB cube using PHINOR method. More...
 
void compute_drb_bgdlixa (const GCOMDri &drm, const int &nrunav=3, const int &navgr=3, const int &nincl=15, const int &nexcl=0)
 Compute DRB cube using BGDLIXA method. More...
 
void compute_drb_bgdlixe (const GCOMDri &drm, const int &nrunav=3, const int &navgr=3, const int &nincl=15, const int &nexcl=0)
 Compute DRB cube using BGDLIXE method. More...
 
void compute_drb_bgdlixf (const GCOMDri &drm, const int &nrunav=3, const int &navgr=3, const int &nincl=15, const int &nexcl=0)
 Compute DRB cube using BGDLIXF method. More...
 
GSkyMap get_weighted_drg_map (void) const
 Return weighted DRG map. More...
 
void get_bgdlixa_phibar_indices (const int &iphibar, const int &nincl, const int &nexcl, int *isel1, int *iex1, int *iex2, int *isel2) const
 Compute Phibar index range for BGDLIXA background method. More...
 
virtual bool use_event_for_likelihood (const int &index) const
 Check whether bin should be used for likelihood analysis. More...
 
- Protected Member Functions inherited from GObservation
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GObservation &obs)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
virtual double likelihood_poisson_unbinned (const GModels &models, GVector *gradients, GMatrixSparse *curvature, double *npred) const
 Evaluate log-likelihood function for Poisson statistic and unbinned analysis (version with working arrays) More...
 
virtual double likelihood_poisson_binned (const GModels &models, GVector *gradients, GMatrixSparse *curvature, double *npred) const
 Evaluate log-likelihood function for Poisson statistic and binned analysis (version with working arrays) More...
 
virtual double likelihood_gaussian_binned (const GModels &models, GVector *gradients, GMatrixSparse *curvature, double *npred) const
 Evaluate log-likelihood function for Gaussian statistic and binned analysis (version with working arrays) More...
 
virtual double npred_spec (const GModel &model, const GTime &obsTime) const
 Integrates spatially integrated Npred kernel spectrally. More...
 

Protected Attributes

std::string m_instrument
 Instrument name. More...
 
GCOMResponse m_response
 Response functions. More...
 
double m_obs_id
 Observation ID. More...
 
double m_ontime
 Ontime (sec) More...
 
double m_livetime
 Livetime (sec) More...
 
double m_deadc
 Deadtime correction. More...
 
GFilename m_drename
 DRE filename. More...
 
GFilename m_drbname
 DRB filename. More...
 
GFilename m_drwname
 DRW filename. More...
 
GFilename m_drgname
 DRG filename. More...
 
GFilename m_drxname
 DRX filename. More...
 
GFilename m_rspname
 Response cache filename. More...
 
GCOMDri m_drb
 Background model. More...
 
GCOMDri m_drw
 Weighting cube. More...
 
GCOMDri m_drg
 Geometry factors. More...
 
GCOMDri m_drx
 Exposure map. More...
 
double m_ewidth
 Energy width (MeV) More...
 
int m_phi_first
 First Phibar layer to use for likelihood. More...
 
int m_phi_last
 Last Phibar layer to use for likelihood. More...
 
GFilename m_evpname
 EVP filename. More...
 
GFilename m_timname
 TIM filename. More...
 
std::vector< GFilenamem_oadnames
 OAD filenames. More...
 
std::vector< GFilenamem_hkdnames
 HKD filenames. More...
 
GFilename m_bvcname
 BVC filename. More...
 
GCOMTim m_tim
 COMPTEL Good Time Intervals. More...
 
GCOMOads m_oads
 Orbit Aspect Data. More...
 
GCOMHkds m_hkds
 Housekeeping Data. More...
 
GCOMBvcs m_bvcs
 Solar System Barycentre Data. More...
 
- Protected Attributes inherited from GObservation
std::string m_name
 Observation name. More...
 
std::string m_id
 Observation identifier. More...
 
std::string m_statistic
 Optimizer statistic. More...
 
GEventsm_events
 Pointer to event container. More...
 
double m_grad_step_size
 Gradient step size. More...
 
std::vector< std::string > m_pars_with_gradients
 

Detailed Description

Interface class for COMPTEL observations.

This class implements a COMPTEL observation. Each COMPTEL observation is defined for a given energy range, and is composed of a DRE, DRB, DRG and DRX file. The DRE file contains the event data, the DRB file contains a background model, the DRG file contains geometry factors, and the DRX file contains the exposure.

Definition at line 66 of file GCOMObservation.hpp.

Constructor & Destructor Documentation

GCOMObservation::GCOMObservation ( void  )

Void constructor.

Creates an empty COMPTEL observation.

Definition at line 89 of file GCOMObservation.cpp.

References init_members().

Referenced by clone().

GCOMObservation::GCOMObservation ( const GXmlElement xml)
explicit

XML constructor.

Parameters
[in]xmlXML element.

Constructs a COMPTEL observation from the information that is found in an XML element.

Definition at line 107 of file GCOMObservation.cpp.

References init_members(), and read().

GCOMObservation::GCOMObservation ( const GCOMDri dre,
const GCOMDri drb,
const GCOMDri drg,
const GCOMDri drx 
)

Binned observation DRI constructor.

Parameters
[in]dreEvent cube.
[in]drbBackground cube.
[in]drgGeometry cube.
[in]drxExposure map.

Creates COMPTEL observation from DRI instances.

The method fixes the deadtime correction factor deadc to 0.965.

Definition at line 132 of file GCOMObservation.cpp.

References drb(), drg(), drx(), GEvents::emax(), GEvents::emin(), GEvents::gti(), init_members(), m_deadc, m_drb, m_drbname, m_drename, m_drg, m_drgname, m_drwname, m_drx, m_drxname, GObservation::m_events, m_ewidth, m_livetime, GObservation::m_name, m_obs_id, m_ontime, and GEnergy::MeV().

GCOMObservation::GCOMObservation ( const GCOMDri dre,
const GCOMDri drb,
const GCOMDri drw,
const GCOMDri drg,
const GCOMDri drx 
)

Binned observation DRI constructor.

Parameters
[in]dreEvent cube.
[in]drbBackground cube.
[in]drwWeighting cube.
[in]drgGeometry cube.
[in]drxExposure map.

Creates COMPTEL observation from DRI instances.

The method fixes the deadtime correction factor deadc to 0.965.

Definition at line 177 of file GCOMObservation.cpp.

References drb(), drg(), drw(), drx(), GEvents::emax(), GEvents::emin(), GEvents::gti(), init_members(), m_deadc, m_drb, m_drbname, m_drename, m_drg, m_drgname, m_drw, m_drwname, m_drx, m_drxname, GObservation::m_events, m_ewidth, m_livetime, GObservation::m_name, m_obs_id, m_ontime, and GEnergy::MeV().

GCOMObservation::GCOMObservation ( const GFilename drename,
const GFilename drbname,
const GFilename drwname,
const GFilename drgname,
const GFilename drxname 
)

Binned observation filename constructor.

Parameters
[in]drenameEvent cube name.
[in]drbnameBackground cube name.
[in]drwnameWeighting cube name.
[in]drgnameGeometry cube name.
[in]drxnameExposure map name.

Creates COMPTEL observation by loading the following FITS files:

DRE - Events cube
DRB - Background model cube
DRW - Weighting cube
DRG - Geometry factors cube
DRX - Exposure map

Each of the four files is mandatory.

Definition at line 230 of file GCOMObservation.cpp.

References init_members(), and load().

GCOMObservation::GCOMObservation ( const GFilename evpname,
const GFilename timname,
const std::vector< GFilename > &  oadnames,
const std::vector< GFilename > &  hkdnames = std::vector<GFilename> (),
const GFilename bvcname = "" 
)

Unbinned observation constructor.

Parameters
[in]evpnameEvent list FITS file name.
[in]timnameGood Time Intervals FITS file name.
[in]oadnamesList of Orbit Aspect Data FITS file names.
[in]hkdnamesList of Housekeeping Data FITS file names.
[in]bvcnameSolar System Barycentre Data FITS file name.

Creates a COMPTEL unbinned observation by loading the event list, Good Time Interval, Orbit Aspect Data and optionally the Solar System Barycentre Data from FITS files.

Except of the Housekeeping Data and the Solar System Barycentre Data all files are mandatory. The Housekeeping Data and the Solar System Barycentre Data will only be loaded if the file names are not empty.

See load() method for more information.

Definition at line 266 of file GCOMObservation.cpp.

References init_members(), and load().

GCOMObservation::GCOMObservation ( const GCOMObservation obs)

Copy constructor.

Parameters
[in]obsCOMPTEL observation.

Creates COMPTEL observation by copying an existing COMPTEL observation.

Definition at line 290 of file GCOMObservation.cpp.

References copy_members(), and init_members().

GCOMObservation::~GCOMObservation ( void  )
virtual

Destructor.

Definition at line 306 of file GCOMObservation.cpp.

References free_members().

Member Function Documentation

const GCOMBvcs & GCOMObservation::bvcs ( void  ) const
inline

Return Solar System Barycentre Data.

Returns
Solar System Barycentre Data

Definition at line 568 of file GCOMObservation.hpp.

References m_bvcs.

Referenced by bvcs(), and GCOMDri::compute_dre().

void GCOMObservation::bvcs ( const GCOMBvcs bvcs)
inline

Set Solar System Barycentre Data.

Parameters
[in]bvcsSolar System Barycentre Data.

Definition at line 581 of file GCOMObservation.hpp.

References bvcs(), and m_bvcs.

bool GCOMObservation::check_dri ( const GCOMDri dri) const
protected

Check if DRI is compatible with event cube.

Parameters
[in]driDRI.
Returns
True if DRI is compatible, false otherwise.

Compares the dimension and the WCS definition of a DRI to that of the event cube. If both are identical, true is returned, false otherwise.

Definition at line 1491 of file GCOMObservation.cpp.

References GObservation::m_events, GCOMDri::map(), GSkyMap::nmaps(), GSkyMap::nx(), GSkyMap::ny(), and GSkyMap::projection().

Referenced by compute_drb_bgdlixa(), compute_drb_bgdlixe(), compute_drb_bgdlixf(), compute_drb_phinor(), load_drb(), load_drg(), and load_drw().

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

Return class name.

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

Implements GObservation.

Definition at line 255 of file GCOMObservation.hpp.

void GCOMObservation::clear ( void  )
virtual

Clear COMPTEL observation.

Implements GObservation.

Definition at line 363 of file GCOMObservation.cpp.

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

Referenced by load(), and read().

GCOMObservation * GCOMObservation::clone ( void  ) const
virtual

Clone COMPTEL observation.

Returns
Pointer to deep copy of COMPTEL observation.

Implements GObservation.

Definition at line 383 of file GCOMObservation.cpp.

References GCOMObservation().

void GCOMObservation::compute_drb ( const std::string &  method,
const GCOMDri drm,
const int &  nrunav = 3,
const int &  navgr = 3,
const int &  nincl = 15,
const int &  nexcl = 0 
)

Compute DRB cube.

Parameters
[in]methodBackground method (PHINOR, BGDLIXA, BGDLIXE or BGDLIXF).
[in]drmDRM cube.
[in]nrunavBGDLIXn: number of bins used for running average.
[in]navgrBGDLIXn: number of bins used for averaging.
[in]ninclBGDLIXn: number of Phibar layers to include.
[in]nexclBGDLIXn: number of Phibar layers to exclude.

Computes a COMPTEL DRB cube using either the PHINOR or BGDLIX method. There are three variants of the BGFLIX method (BGDLIXA, BGDLIXE and BGDLIXF). See the protected methods compute_drb_phinor(), compute_drb_bgdlixa(), compute_drb_bgdlixe(), and compute_drb_bgdlixf() for more information.

Definition at line 1027 of file GCOMObservation.cpp.

References compute_drb_bgdlixa(), compute_drb_bgdlixe(), compute_drb_bgdlixf(), compute_drb_phinor(), and G_COMPUTE_DRB.

void GCOMObservation::compute_drb_bgdlixa ( const GCOMDri drm,
const int &  nrunav = 3,
const int &  navgr = 3,
const int &  nincl = 15,
const int &  nexcl = 0 
)
protected

Compute DRB cube using BGDLIXA method.

Parameters
[in]drmDRM cube.
[in]nrunavNumber of bins used for running average.
[in]navgrNumber of bins used for averaging.
[in]ninclNumber of Phibar layers to include.
[in]nexclNumber of Phibar layers to exclude.
Exceptions
GException::invalid_valueObservation does not contain an event cube
GException::invalid_argumentDRM cube is incompatible with DRE

Computes a DRB cube using the BGDLIXA method that is documented in Rob van Dijk's PhD thesis. The revelant equations from the thesis that are implemented here are 3.12, 3.12 and 3.14.

Definition at line 1665 of file GCOMObservation.cpp.

References check_dri(), GCOMEventCube::dre(), GObservation::events(), G_COMPUTE_DRB_BGDLIXA, get_bgdlixa_phibar_indices(), get_weighted_drg_map(), GObservation::id(), m_drb, m_drg, GCOMDri::map(), GObservation::name(), GCOMDri::nchi(), norm(), GCOMDri::nphibar(), and GCOMDri::npsi().

Referenced by compute_drb().

void GCOMObservation::compute_drb_bgdlixe ( const GCOMDri drm,
const int &  nrunav = 3,
const int &  navgr = 3,
const int &  nincl = 15,
const int &  nexcl = 0 
)
protected

Compute DRB cube using BGDLIXE method.

Parameters
[in]drmDRM cube.
[in]nrunavNumber of bins used for running average.
[in]navgrNumber of bins used for averaging.
[in]ninclNumber of Phibar layers to include.
[in]nexclNumber of Phibar layers to exclude.
Exceptions
GException::invalid_valueObservation does not contain an event cube
GException::invalid_argumentDRM cube is incompatible with DRE

Computes a DRB cube using the BGDLIXE method. This method differs from the BGDLIXA method in the last step.

Definition at line 2016 of file GCOMObservation.cpp.

References check_dri(), GCOMEventCube::dre(), GObservation::events(), G_COMPUTE_DRB_BGDLIXE, get_bgdlixa_phibar_indices(), get_weighted_drg_map(), GObservation::id(), m_drb, m_drg, GCOMDri::map(), GObservation::name(), GCOMDri::nchi(), norm(), GCOMDri::nphibar(), and GCOMDri::npsi().

Referenced by compute_drb().

void GCOMObservation::compute_drb_bgdlixf ( const GCOMDri drm,
const int &  nrunav = 3,
const int &  navgr = 3,
const int &  nincl = 15,
const int &  nexcl = 0 
)
protected

Compute DRB cube using BGDLIXF method.

Parameters
[in]drmDRM cube.
[in]nrunavNumber of bins used for running average.
[in]navgrNumber of bins used for averaging.
[in]ninclNumber of Phibar layers to include.
[in]nexclNumber of Phibar layers to exclude.
Exceptions
GException::invalid_valueObservation does not contain an event cube
GException::invalid_argumentDRM cube is incompatible with DRE

Computes a DRB cube using the BGDLIXF method. This method differs from the BGDLIXE method in the DRW instead of the DRG is used for background model computation.

Definition at line 2227 of file GCOMObservation.cpp.

References check_dri(), GCOMEventCube::dre(), GObservation::events(), G_COMPUTE_DRB_BGDLIXE, G_COMPUTE_DRB_BGDLIXF, get_bgdlixa_phibar_indices(), GObservation::id(), m_drb, m_drw, GCOMDri::map(), GObservation::name(), GCOMDri::nchi(), norm(), GCOMDri::nphibar(), and GCOMDri::npsi().

Referenced by compute_drb().

void GCOMObservation::compute_drb_phinor ( const GCOMDri drm)
protected

Compute DRB cube using PHINOR method.

Parameters
[in]drmDRM cube.
Exceptions
GException::invalid_valueObservation does not contain an event cube
GException::invalid_argumentDRM cube is incompatible with DRE

Definition at line 1577 of file GCOMObservation.cpp.

References check_dri(), GCOMEventCube::dre(), GObservation::events(), G_COMPUTE_DRB_PHINOR, get_weighted_drg_map(), GObservation::id(), m_drb, m_drg, GCOMDri::map(), GObservation::name(), norm(), GCOMDri::nphibar(), and GSkyMap::npix().

Referenced by compute_drb().

void GCOMObservation::copy_members ( const GCOMObservation obs)
protected
double GCOMObservation::deadc ( const GTime time = GTime()) const
inlinevirtual

Return deadtime correction factor.

Parameters
[in]timeTime.
Returns
Deadtime correction factor.

Implements GObservation.

Definition at line 334 of file GCOMObservation.hpp.

References m_deadc.

Referenced by deadc(), GCOMResponse::irf_diffuse(), GCOMResponse::irf_elliptical(), GCOMResponse::irf_ptsrc(), and GCOMResponse::irf_radial().

void GCOMObservation::deadc ( const double &  deadc)
inline

Set deadtime correction factor.

Parameters
[in]deadcDeadtime correction factor.

Definition at line 386 of file GCOMObservation.hpp.

References deadc(), and m_deadc.

const GCOMDri & GCOMObservation::drb ( void  ) const
inline

Return background model.

Returns
Background model.

Definition at line 438 of file GCOMObservation.hpp.

References m_drb.

Referenced by GCOMModelDRBPhibarBins::eval(), GCOMModelDRBPhibarNodes::eval(), and GCOMObservation().

const GFilename & GCOMObservation::drbname ( void  ) const
inline

Return DRB filename.

Returns
DRB filename.

Definition at line 631 of file GCOMObservation.hpp.

References m_drbname.

Referenced by drbname(), load_drb(), and read().

void GCOMObservation::drbname ( const GFilename drbname)
inline

Set DRB filename.

Parameters
[in]drbnameDRB filename.

Definition at line 735 of file GCOMObservation.hpp.

References drbname(), and m_drbname.

const GFilename & GCOMObservation::drename ( void  ) const
inline

Return DRE filename.

Returns
DRE filename.

Definition at line 618 of file GCOMObservation.hpp.

References m_drename.

Referenced by drename(), load_dre(), and read().

void GCOMObservation::drename ( const GFilename drename)
inline

Set DRE filename.

Parameters
[in]drenameDRE filename.

Definition at line 722 of file GCOMObservation.hpp.

References drename(), and m_drename.

const GCOMDri & GCOMObservation::drg ( void  ) const
inline
const GFilename & GCOMObservation::drgname ( void  ) const
inline

Return DRG filename.

Returns
DRG filename.

Definition at line 657 of file GCOMObservation.hpp.

References m_drgname.

Referenced by drgname(), load_drg(), and read().

void GCOMObservation::drgname ( const GFilename drgname)
inline

Set DRG filename.

Parameters
[in]drgnameDRG filename.

Definition at line 761 of file GCOMObservation.hpp.

References drgname(), and m_drgname.

GCOMDri GCOMObservation::drm ( const GModels models) const

Compute DRM cube.

Parameters
[in]modelsModel container.
Returns
DRM cube (units of counts)
Exceptions
GException::invalid_valueObservation does not contain an event cube.

Computes a COMPTEL DRM cube from the information provided in a model container. The values of the DRM cube are in units of counts.

Definition at line 968 of file GCOMObservation.cpp.

References GCOMEventBin::counts(), GCOMEventCube::dre(), GModels::eval(), G_DRM, GObservation::m_events, GObservation::model(), GObservation::name(), GCOMEventBin::size(), and GCOMEventCube::size().

const GCOMDri & GCOMObservation::drw ( void  ) const
inline

Return weighting cube.

Returns
Weighting cube.

Definition at line 451 of file GCOMObservation.hpp.

References m_drw.

Referenced by GCOMObservation().

const GFilename & GCOMObservation::drwname ( void  ) const
inline

Return DRW filename.

Returns
DRW filename.

Definition at line 644 of file GCOMObservation.hpp.

References m_drwname.

Referenced by drwname(), load_drw(), and read().

void GCOMObservation::drwname ( const GFilename drwname)
inline

Set DRW filename.

Parameters
[in]drwnameDRW filename.

Definition at line 748 of file GCOMObservation.hpp.

References drwname(), and m_drwname.

const GCOMDri & GCOMObservation::drx ( void  ) const
inline
const GFilename & GCOMObservation::drxname ( void  ) const
inline

Return DRX filename.

Returns
DRX filename.

Definition at line 670 of file GCOMObservation.hpp.

References m_drxname.

Referenced by drxname(), load_drx(), and read().

void GCOMObservation::drxname ( const GFilename drxname)
inline

Set DRX filename.

Parameters
[in]drxnameDRX filename.

Definition at line 774 of file GCOMObservation.hpp.

References drxname(), and m_drxname.

void GCOMObservation::ewidth ( const double &  ewidth)
inline

Set energy width.

Parameters
[in]ewidthEnergy width (MeV).

Definition at line 399 of file GCOMObservation.hpp.

References ewidth(), and m_ewidth.

const double & GCOMObservation::ewidth ( void  ) const
inline

Return energy width.

Returns
Energy width (MeV).

Definition at line 425 of file GCOMObservation.hpp.

References m_ewidth.

Referenced by ewidth().

void GCOMObservation::free_members ( void  )
protected

Delete class members.

Definition at line 1253 of file GCOMObservation.cpp.

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

void GCOMObservation::get_bgdlixa_phibar_indices ( const int &  iphibar,
const int &  nincl,
const int &  nexcl,
int *  isel1,
int *  iex1,
int *  iex2,
int *  isel2 
) const
protected

Compute Phibar index range for BGDLIXA background method.

Parameters
[in]iphibarPhibar layer index.
[in]ninclNumber of Phibar layers to include.
[in]nexclNumber of Phibar layers to exclude.
[out]isel1Start index for first sum.
[out]iex1Stop index for first sum.
[out]iex2Start index for second sum.
[out]isel2Stop index for second sum.

This method is a helper method for the BGDLIXA background method. It computes the Phibar index range for the third step of the background model computation. The third step corresponds to Equation (3.14) in Rob van Dijk's PhD thesis.

The Phibar sum will be taken over the index ranges [isel1, iex1] and [ixe2, isel2].

Definition at line 2469 of file GCOMObservation.cpp.

References m_drg, and GCOMDri::nphibar().

Referenced by compute_drb_bgdlixa(), compute_drb_bgdlixe(), and compute_drb_bgdlixf().

GSkyMap GCOMObservation::get_weighted_drg_map ( void  ) const
protected

Return weighted DRG map.

Returns
Weighted DRG map.

Returns a DRG as sky map where each pixel of the map is multiplied by the solidangle of the pixel.

Definition at line 2428 of file GCOMObservation.cpp.

References drg(), m_drg, GCOMDri::map(), GSkyMap::nmaps(), GSkyMap::npix(), and GSkyMap::solidangle().

Referenced by compute_drb_bgdlixa(), compute_drb_bgdlixe(), and compute_drb_phinor().

const GCOMHkds & GCOMObservation::hkds ( void  ) const
inline

Return Housekeeping Data collection.

Returns
Housekeeping Data collection

Definition at line 542 of file GCOMObservation.hpp.

References m_hkds.

Referenced by hkds(), and GCOMDris::vetorate_setup().

void GCOMObservation::hkds ( const GCOMHkds hkds)
inline

Set Housekeeping Data collection.

Parameters
[in]hkdsHousekeeping Data collection.

Definition at line 555 of file GCOMObservation.hpp.

References hkds(), and m_hkds.

std::string GCOMObservation::instrument ( void  ) const
inlinevirtual

Return instrument.

Returns
Instrument name.

Implements GObservation.

Definition at line 293 of file GCOMObservation.hpp.

References m_instrument.

Referenced by print().

bool GCOMObservation::is_binned ( void  ) const
inline

Check whether observation is binned.

Returns
True if observation is unbinned.

Definition at line 606 of file GCOMObservation.hpp.

References GObservation::m_events.

Referenced by write().

bool GCOMObservation::is_unbinned ( void  ) const
inline

Check whether observation is unbinned.

Returns
True if observation is unbinned.

Definition at line 594 of file GCOMObservation.hpp.

References GObservation::m_events.

Referenced by write().

double GCOMObservation::livetime ( void  ) const
inlinevirtual

Return livetime.

Returns
Livetime (seconds).

Implements GObservation.

Definition at line 319 of file GCOMObservation.hpp.

References m_livetime.

Referenced by livetime(), and print().

void GCOMObservation::livetime ( const double &  livetime)
inline

Set livetime.

Parameters
[in]livetimeLivetime.

Definition at line 373 of file GCOMObservation.hpp.

References livetime(), and m_livetime.

void GCOMObservation::load ( const GFilename drename,
const GFilename drbname,
const GFilename drwname,
const GFilename drgname,
const GFilename drxname 
)

Load data for a binned observation.

Parameters
[in]drenameEvent cube name.
[in]drbnameBackground cube name.
[in]drwnameWeighting cube name.
[in]drgnameGeometry cube name.
[in]drxnameExposure map name.

Load event cube from DRE file, background model from DRB file, weigthing cube from DRW file, geometry factors from DRG file and the exposure map from the DRX file. All files are mandatory.

Definition at line 814 of file GCOMObservation.cpp.

References load_drb(), load_dre(), load_drg(), load_drw(), and load_drx().

Referenced by GCOMObservation(), and read().

void GCOMObservation::load ( const GFilename evpname,
const GFilename timname,
const std::vector< GFilename > &  oadnames,
const std::vector< GFilename > &  hkdnames = std::vector<GFilename> (),
const GFilename bvcname = "" 
)

Load data for an unbinned observation.

Parameters
[in]evpnameEvent list FITS file name.
[in]timnameGood Time Intervals FITS file name.
[in]oadnamesList of Orbit Aspect Data FITS file names.
[in]hkdnamesList of Housekeeping Data FITS file names.
[in]bvcnameSolar System Barycentre Data FITS file name.

Loads the event list, Good Time Interval, Orbit Aspect Data and optionally Housekeeping data and Solar System Barycentre Data for an unbinned observation.

Except of the Housekeeping Data and the Solar System Barycentre Data all files are mandatory. The Housekeeping Data and the Solar System Barycentre Data will only be loaded if the file names are not empty.

The method fixes the deadtime correction factor deadc to 0.965.

Definition at line 859 of file GCOMObservation.cpp.

References clear(), GFits::close(), GCOMOads::extend(), GCOMHkds::extend(), GCOMTim::gti(), GCOMTim::load(), GCOMBvcs::load(), m_bvcname, m_bvcs, m_deadc, GObservation::m_events, m_evpname, m_hkdnames, m_hkds, m_livetime, m_oadnames, m_oads, m_ontime, m_tim, m_timname, oads(), GGti::ontime(), read_attributes(), GCOMOads::size(), and GCOMHkds::size().

void GCOMObservation::load_drb ( const GFilename drbname)
protected

Load background model from DRB file.

Parameters
[in]drbnameDRB filename.
Exceptions
GException::invalid_valueDRB data space incompatible with DRE data space.

Load the background model from the primary image of the specified FITS file. Since a DRB file is optional the method does nothing if the DRB filename is empty.

Definition at line 1326 of file GCOMObservation.cpp.

References check_dri(), GFits::close(), gammalib::com_wcs_mer2car(), drbname(), G_LOAD_DRB, GFits::image(), GFilename::is_empty(), m_drb, m_drbname, m_drename, GCOMDri::map(), and GCOMDri::read().

Referenced by load().

void GCOMObservation::load_dre ( const GFilename drename)
protected

Load event cube data from DRE file.

Parameters
[in]drenameDRE filename.

Loads the event cube from a DRE file.

The ontime is extracted from the Good Time Intervals. The deadtime correction factor deadc is fixed to 0.965. The livetime is computed by multiplying the deadtime correction by the ontime, i.e. LIVETIME = ONTIME * DEADC.

Definition at line 1272 of file GCOMObservation.cpp.

References GFits::close(), drename(), GEvents::emax(), GEvents::emin(), GEvents::gti(), m_deadc, m_drename, GObservation::m_events, m_ewidth, m_livetime, m_ontime, GEnergy::MeV(), GEvents::read(), and read_attributes().

Referenced by load().

void GCOMObservation::load_drg ( const GFilename drgname)
protected

Load geometry factors from DRG file.

Parameters
[in]drgnameDRG filename.
Exceptions
GException::invalid_valueDRG data space incompatible with DRE data space.

Load the geometry factors from the primary image of the specified FITS file.

Definition at line 1418 of file GCOMObservation.cpp.

References check_dri(), GFits::close(), gammalib::com_wcs_mer2car(), drgname(), G_LOAD_DRG, GFits::image(), m_drename, m_drg, m_drgname, GCOMDri::map(), and GCOMDri::read().

Referenced by load().

void GCOMObservation::load_drw ( const GFilename drwname)
protected

Load weighting cube from DRW file.

Parameters
[in]drwnameDRW filename.
Exceptions
GException::invalid_valueDRW data space incompatible with DRE data space.

Load the weighting cube from the primary image of the specified FITS file.

Definition at line 1373 of file GCOMObservation.cpp.

References check_dri(), GFits::close(), drwname(), G_LOAD_DRW, GFits::image(), GFilename::is_empty(), m_drename, m_drw, m_drwname, and GCOMDri::read().

Referenced by load().

void GCOMObservation::load_drx ( const GFilename drxname)
protected

Load exposure from DRX file.

Parameters
[in]drxnameDRX filename.

Load the exposure map from the primary image of the specified FITS file.

Definition at line 1457 of file GCOMObservation.cpp.

References GFits::close(), gammalib::com_wcs_mer2car(), drxname(), GFits::image(), m_drx, m_drxname, GCOMDri::map(), and GCOMDri::read().

Referenced by load().

double GCOMObservation::npred ( const GModel model) const
virtual

Return total number of predicted counts for one model.

Parameters
[in]modelGamma-ray source model.
Returns
Total number of predicted counts in model.

Returns the total number of predicted counts in a given model component.

Reimplemented from GObservation.

Definition at line 452 of file GCOMObservation.cpp.

References GModels::append(), GObservation::events(), GObservation::model(), GEventBin::size(), GVector::size(), and use_event_for_likelihood().

const GCOMOads & GCOMObservation::oads ( void  ) const
inline
void GCOMObservation::oads ( const GCOMOads oads)
inline

Set Orbit Aspect Data.

Parameters
[in]oadsOrbit Aspect Data.

Definition at line 529 of file GCOMObservation.hpp.

References m_oads, and oads().

void GCOMObservation::obs_id ( const double &  id)
inline

Set observation ID.

Parameters
[in]idObservation ID.

Definition at line 347 of file GCOMObservation.hpp.

References GObservation::id(), and m_obs_id.

const double & GCOMObservation::obs_id ( void  ) const
inline

Return observation ID.

Returns
Observation ID.

Definition at line 412 of file GCOMObservation.hpp.

References m_obs_id.

double GCOMObservation::ontime ( void  ) const
inlinevirtual

Return ontime.

Returns
Ontime (seconds).

Implements GObservation.

Definition at line 306 of file GCOMObservation.hpp.

References m_ontime.

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

void GCOMObservation::ontime ( const double &  ontime)
inline

Set ontime.

Parameters
[in]ontimeOntime.

Definition at line 360 of file GCOMObservation.hpp.

References m_ontime, and ontime().

GCOMObservation & GCOMObservation::operator= ( const GCOMObservation obs)
virtual

Assignment operator.

Parameters
[in]obsCOMPTEL observation.
Returns
COMPTEL observation.

Assign COMPTEL observation to this object.

Definition at line 330 of file GCOMObservation.cpp.

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

const int & GCOMObservation::phi_first ( void  ) const
inline

Return index of first Phibar layer to be used for likelihood fitting.

Returns
Index of first Phibar layer.

Definition at line 696 of file GCOMObservation.hpp.

References m_phi_first.

Referenced by phi_first(), and print().

void GCOMObservation::phi_first ( const int &  phi_first)
inline

Set index of first Phibar layer to be used for likelihood fitting.

Parameters
[in]phi_firstIndex of first Phibar layer.

Definition at line 800 of file GCOMObservation.hpp.

References m_phi_first, and phi_first().

const int & GCOMObservation::phi_last ( void  ) const
inline

Return index of last Phibar layer to be used for likelihood fitting.

Returns
Index of last Phibar layer.

Definition at line 709 of file GCOMObservation.hpp.

References m_phi_last.

Referenced by phi_last(), and print().

void GCOMObservation::phi_last ( const int &  phi_last)
inline

Set index of last Phibar layer to be used for likelihood fitting.

Parameters
[in]phi_lastIndex of last Phibar layer.

Definition at line 813 of file GCOMObservation.hpp.

References m_phi_last, and phi_last().

void GCOMObservation::read ( const GXmlElement xml)
virtual

Read observation from XML element.

Parameters
[in]xmlXML element.

Reads information for a COMPTEL observation from an XML element. The method supports both an unbinned and a binned observation.

For an unbinned observation the XML format is

<observation name="Crab" id="000001" instrument="COM">
  <parameter name="EVP" file="m16992_evp.fits"/>
  <parameter name="TIM" file="m10695_tim.fits"/>
  <parameter name="OAD" file="m20039_oad.fits"/>
  <parameter name="OAD" file="m20041_oad.fits"/>
  <parameter name="HKD" file="m20035_hkd.fits"/>
  <parameter name="HKD" file="m20037_hkd.fits"/>
  <parameter name="BVC" file="s10150_bvc.fits"/>
  ...
</observation>

where the observation can contain an arbitrary number of OAD file parameters. The file attribute provide either absolute or relative file names. If a file name includes no access path it is assumed that the file resides in the same location as the XML file. The HKD and BVC files are optional and do not need to be specified.

For a binned observation the XML format is

<observation name="Crab" id="000001" instrument="COM">
  <parameter name="DRE" file="m50438_dre.fits"/>
  <parameter name="DRB" file="m34997_drg.fits"/>
  <parameter name="DRW" file="m34997_drw.fits"/>
  <parameter name="DRG" file="m34997_drg.fits"/>
  <parameter name="DRX" file="m32171_drx.fits"/>
  <parameter name="IAQ" value="UNH(1.0-3.0)MeV"/>
</observation>

Note that the DRW parameter is optional.

Implements GObservation.

Definition at line 532 of file GCOMObservation.cpp.

References GXmlElement::attribute(), clear(), drbname(), drename(), drgname(), drwname(), drxname(), GXmlNode::element(), GXmlNode::elements(), GFilename::exists(), G_READ, GXmlElement::has_attribute(), GFilename::is_fits(), load(), GCOMResponse::load_cache(), GObservation::m_id, m_instrument, GObservation::m_name, m_phi_first, m_phi_last, m_response, m_rspname, response(), gammalib::toint(), gammalib::xml_file_expand(), gammalib::xml_get_attr(), and gammalib::xml_has_par().

Referenced by GCOMObservation().

void GCOMObservation::read_attributes ( const GFitsHDU hdu)
protected

Read observation attributes.

Parameters
[in]hduFITS HDU pointer

Reads optional attributes are

OBS_ID   - Observation identifier
OBJECT   - Object

Nothing is done if the HDU pointer is NULL.

Definition at line 1530 of file GCOMObservation.cpp.

References GFitsHDU::has_card(), GObservation::m_name, m_obs_id, GFitsHDU::real(), and GFitsHDU::string().

Referenced by load(), and load_dre().

void GCOMObservation::response ( const GResponse rsp)
virtual

Set response function.

Parameters
[in]rspResponse function.
Exceptions
GException::invalid_argumentSpecified response is not a COMPTEL response.

Sets the response function for the observation.

Implements GObservation.

Definition at line 399 of file GCOMObservation.cpp.

References G_RESPONSE, m_response, and GObservation::name().

const GCOMResponse * GCOMObservation::response ( void  ) const
inlinevirtual

Return response function.

Returns
Response function.

Implements GObservation.

Definition at line 267 of file GCOMObservation.hpp.

References m_response.

Referenced by print(), read(), and response().

void GCOMObservation::response ( const GCaldb caldb,
const std::string &  rspname 
)

Set response function.

Parameters
[in]caldbCalibration database.
[in]rspnameName of COMPTEL response.

Sets the response function by loading the response information from the calibration database.

Definition at line 428 of file GCOMObservation.cpp.

References GCOMResponse::caldb(), GCOMResponse::clear(), GCOMResponse::load(), and m_response.

void GCOMObservation::response ( const GCOMResponse response)
inline

Set response function.

Parameters
[in]responseResponse function.

Definition at line 280 of file GCOMObservation.hpp.

References m_response, and response().

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

Return response cache filename.

Returns
Response cache filename.

Definition at line 683 of file GCOMObservation.hpp.

References m_rspname.

Referenced by print(), and rspname().

void GCOMObservation::rspname ( const GFilename rspname)
inline

Set response cache filename filename.

Parameters
[in]cachenameResponse cache filename.

Definition at line 787 of file GCOMObservation.hpp.

References m_rspname, and rspname().

const GCOMTim & GCOMObservation::tim ( void  ) const
inline

Return COMPTEL Good Time Intervals.

Returns
COMPTEL Good Time Intervals.

Definition at line 490 of file GCOMObservation.hpp.

References m_tim.

Referenced by GCOMDri::compute_dre(), GCOMDri::compute_drg(), GCOMDris::compute_drws_energy(), GCOMDris::compute_drws_phibar(), GCOMDri::compute_drx(), tim(), GCOMDris::vetorate_generate(), and GCOMDris::vetorate_setup().

void GCOMObservation::tim ( const GCOMTim tim)
inline

Set COMPTEL Good Time Intervals.

Parameters
[in]timCOMPTEL Good Time Intervals.

Definition at line 503 of file GCOMObservation.hpp.

References m_tim, and tim().

bool GCOMObservation::use_event_for_likelihood ( const int &  index) const
protectedvirtual

Check whether bin should be used for likelihood analysis.

Parameters
[in]indexEvent index.
Returns
True if event with specified index should be used.

Implements the Phibar event selection.

Reimplemented from GObservation.

Definition at line 2535 of file GCOMObservation.cpp.

References m_drg, m_phi_first, m_phi_last, GCOMDri::map(), and GSkyMap::npix().

Referenced by npred().

void GCOMObservation::write ( GXmlElement xml) const
virtual

Write observation into XML element.

Parameters
[in]xmlXML element.

Writes information for a COMPTEL observation into an XML element. The method supports both an unbinned and a binned observation.

For an unbinned observation the XML format is

<observation name="Crab" id="000001" instrument="COM">
  <parameter name="EVP" file="m16992_evp.fits"/>
  <parameter name="TIM" file="m10695_tim.fits"/>
  <parameter name="OAD" file="m20039_oad.fits"/>
  <parameter name="OAD" file="m20041_oad.fits"/>
  <parameter name="HKD" file="m20035_hkd.fits"/>
  <parameter name="HKD" file="m20037_hkd.fits"/>
  <parameter name="BVC" file="s10150_bvc.fits"/>
  ...
</observation>

where the observation can contain an arbitrary number of OAD file parameters. The file attribute provide either absolute or relative file names. If a file name includes no access path it is assumed that the file resides in the same location as the XML file. The HKD and BVC files are optional and are only written if HKD and BVC information is contained in the observation.

For a binned observation the XML format is

<observation name="Crab" id="000001" instrument="COM">
  <parameter name="DRE" file="m50438_dre.fits"/>
  <parameter name="DRB" file="m34997_drg.fits"/>
  <parameter name="DRW" file="m34997_drw.fits"/>
  <parameter name="DRG" file="m34997_drg.fits"/>
  <parameter name="DRX" file="m32171_drx.fits"/>
  <parameter name="IAQ" value="UNH(1.0-3.0)MeV"/>
</observation>

Note that the DRW parameter is optional.

Implements GObservation.

Definition at line 694 of file GCOMObservation.cpp.

References GXmlNode::append(), GXmlElement::attribute(), GXmlNode::element(), GXmlNode::elements(), GFilename::exists(), G_WRITE, is_binned(), GCOMBvcs::is_empty(), GFilename::is_empty(), is_unbinned(), m_bvcname, m_bvcs, m_drbname, m_drename, m_drgname, m_drwname, m_drxname, m_evpname, m_hkdnames, m_oadnames, m_phi_first, m_phi_last, m_response, m_rspname, m_timname, GXmlNode::remove(), GCOMResponse::rspname(), GCOMResponse::save_cache(), gammalib::str(), gammalib::xml_file_reduce(), and gammalib::xml_need_par().

void GCOMObservation::write_attributes ( GFitsHDU hdu) const
protected

Write observation attributes.

Parameters
[in]hduFITS HDU pointer

Nothing is done if the HDU pointer is NULL.

Todo:
Implement method.

Definition at line 1555 of file GCOMObservation.cpp.

Member Data Documentation

GFilename GCOMObservation::m_bvcname
protected

BVC filename.

Definition at line 241 of file GCOMObservation.hpp.

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

GCOMBvcs GCOMObservation::m_bvcs
protected

Solar System Barycentre Data.

Definition at line 245 of file GCOMObservation.hpp.

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

double GCOMObservation::m_deadc
protected

Deadtime correction.

Definition at line 219 of file GCOMObservation.hpp.

Referenced by copy_members(), deadc(), GCOMObservation(), init_members(), load(), load_dre(), and print().

GCOMDri GCOMObservation::m_drb
protected
GFilename GCOMObservation::m_drbname
protected

DRB filename.

Definition at line 223 of file GCOMObservation.hpp.

Referenced by copy_members(), drbname(), GCOMObservation(), init_members(), load_drb(), and write().

GFilename GCOMObservation::m_drename
protected
GFilename GCOMObservation::m_drgname
protected

DRG filename.

Definition at line 225 of file GCOMObservation.hpp.

Referenced by copy_members(), drgname(), GCOMObservation(), init_members(), load_drg(), and write().

GCOMDri GCOMObservation::m_drw
protected

Weighting cube.

Definition at line 229 of file GCOMObservation.hpp.

Referenced by compute_drb_bgdlixf(), copy_members(), drw(), GCOMObservation(), init_members(), and load_drw().

GFilename GCOMObservation::m_drwname
protected

DRW filename.

Definition at line 224 of file GCOMObservation.hpp.

Referenced by copy_members(), drwname(), GCOMObservation(), init_members(), load_drw(), and write().

GCOMDri GCOMObservation::m_drx
protected

Exposure map.

Definition at line 231 of file GCOMObservation.hpp.

Referenced by copy_members(), drx(), GCOMObservation(), init_members(), and load_drx().

GFilename GCOMObservation::m_drxname
protected

DRX filename.

Definition at line 226 of file GCOMObservation.hpp.

Referenced by copy_members(), drxname(), GCOMObservation(), init_members(), load_drx(), and write().

GFilename GCOMObservation::m_evpname
protected

EVP filename.

Definition at line 237 of file GCOMObservation.hpp.

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

double GCOMObservation::m_ewidth
protected

Energy width (MeV)

Definition at line 232 of file GCOMObservation.hpp.

Referenced by copy_members(), ewidth(), GCOMObservation(), init_members(), and load_dre().

std::vector<GFilename> GCOMObservation::m_hkdnames
protected

HKD filenames.

Definition at line 240 of file GCOMObservation.hpp.

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

GCOMHkds GCOMObservation::m_hkds
protected

Housekeeping Data.

Definition at line 244 of file GCOMObservation.hpp.

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

std::string GCOMObservation::m_instrument
protected

Instrument name.

Definition at line 214 of file GCOMObservation.hpp.

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

double GCOMObservation::m_livetime
protected

Livetime (sec)

Definition at line 218 of file GCOMObservation.hpp.

Referenced by copy_members(), GCOMObservation(), init_members(), livetime(), load(), and load_dre().

std::vector<GFilename> GCOMObservation::m_oadnames
protected

OAD filenames.

Definition at line 239 of file GCOMObservation.hpp.

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

GCOMOads GCOMObservation::m_oads
protected

Orbit Aspect Data.

Definition at line 243 of file GCOMObservation.hpp.

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

double GCOMObservation::m_obs_id
protected

Observation ID.

Definition at line 216 of file GCOMObservation.hpp.

Referenced by copy_members(), GCOMObservation(), init_members(), obs_id(), and read_attributes().

double GCOMObservation::m_ontime
protected

Ontime (sec)

Definition at line 217 of file GCOMObservation.hpp.

Referenced by copy_members(), GCOMObservation(), init_members(), load(), load_dre(), and ontime().

int GCOMObservation::m_phi_first
protected

First Phibar layer to use for likelihood.

Definition at line 233 of file GCOMObservation.hpp.

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

int GCOMObservation::m_phi_last
protected

Last Phibar layer to use for likelihood.

Definition at line 234 of file GCOMObservation.hpp.

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

GCOMResponse GCOMObservation::m_response
protected

Response functions.

Definition at line 215 of file GCOMObservation.hpp.

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

GFilename GCOMObservation::m_rspname
protected

Response cache filename.

Definition at line 227 of file GCOMObservation.hpp.

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

GCOMTim GCOMObservation::m_tim
protected

COMPTEL Good Time Intervals.

Definition at line 242 of file GCOMObservation.hpp.

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

GFilename GCOMObservation::m_timname
protected

TIM filename.

Definition at line 238 of file GCOMObservation.hpp.

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


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