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

CTA observation class. More...

#include <GCTAObservation.hpp>

Inheritance diagram for GCTAObservation:
GObservation GBase

Public Member Functions

 GCTAObservation (void)
 Void constructor. More...
 
 GCTAObservation (const bool &dummy, const std::string &instrument)
 Instrument constructor. More...
 
 GCTAObservation (const GFilename &filename)
 Unbinned or binned analysis constructor. More...
 
 GCTAObservation (const GFilename &cntcube, const GFilename &expcube, const GFilename &psfcube, const GFilename &bkgcube)
 Stacked analysis constructor. More...
 
 GCTAObservation (const GFilename &cntcube, const GFilename &expcube, const GFilename &psfcube, const GFilename &edispcube, const GFilename &bkgcube)
 Stacked analysis with energy dispersion constructor. More...
 
 GCTAObservation (const GCTAObservation &obs)
 Copy constructor. More...
 
virtual ~GCTAObservation (void)
 Destructor. More...
 
GCTAObservationoperator= (const GCTAObservation &obs)
 Assignment operator. More...
 
virtual void clear (void)
 Clear CTA observation. More...
 
virtual GCTAObservationclone (void) const
 Clone CTA observation. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual void response (const GResponse &rsp)
 Set response function. More...
 
virtual const GCTAResponseresponse (void) const
 Return pointer to CTA response function. More...
 
virtual std::string instrument (void) const
 Return instrument name. 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 CTA observation information. More...
 
void instrument (const std::string &instrument)
 Set instrument name. More...
 
bool has_response (void) const
 Signal if CTA observation contains response information. More...
 
bool has_events (void) const
 Signal if CTA observation contains events. More...
 
void read (const GFits &fits)
 Read data from FITS file. More...
 
void write (GFits &fits, const std::string &evtname=gammalib::extname_cta_events, const std::string &gtiname=gammalib::extname_gti) const
 Write observation into FITS file. More...
 
void load (const GFilename &filename)
 Load unbinned or binned analysis data. More...
 
void load (const GFilename &cntcube, const GFilename &expcube, const GFilename &psfcube, const GFilename &bkgcube)
 Load stacked analysis data. More...
 
void load (const GFilename &cntcube, const GFilename &expcube, const GFilename &psfcube, const GFilename &edispcube, const GFilename &bkgcube)
 Load stacked analysis data with energy dispersion. More...
 
void save (const GFilename &filename, const bool &clobber=false) const
 Save CTA observation into FITS file. More...
 
void response (const std::string &rspname, const GCaldb &caldb)
 Set CTA response function. More...
 
void response (const GCTACubeExposure &expcube, const GCTACubePsf &psfcube, const GCTACubeBackground &bkgcube)
 Set CTA response function. More...
 
void response (const GCTACubeExposure &expcube, const GCTACubePsf &psfcube, const GCTACubeEdisp &edispcube, const GCTACubeBackground &bkgcube)
 Set CTA response function. More...
 
void pointing (const GCTAPointing &pointing)
 Set CTA pointing. More...
 
const GCTAPointingpointing (void) const
 Return CTA pointing. More...
 
void off_regions (const GSkyRegions &off_regions)
 Set sky off regions. More...
 
const GSkyRegionsoff_regions (void) const
 Return sky off regions. More...
 
GCTARoi roi (void) const
 Get Region of Interest. More...
 
GGti gti (void) const
 Get Good Time Intervals. More...
 
GEbounds ebounds (void) const
 Get energy boundaries. More...
 
void object (const std::string &object)
 Set CTA observation object name. More...
 
const std::string & object (void) const
 Return CTA observation object. More...
 
void ra_obj (const double &ra)
 Set CTA object Right Ascension. More...
 
const double & ra_obj (void) const
 Return CTA object Right Ascension. More...
 
void dec_obj (const double &dec)
 Set CTA object Declination. More...
 
const double & dec_obj (void) const
 Return CTA object Declination. 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. More...
 
void eventfile (const GFilename &filename)
 Set event file name. More...
 
const GFilenameeventfile (void) const
 Return event file name. More...
 
std::string eventtype (void) const
 Return event type string. More...
 
void dispose_events (void)
 Dispose events. More...
 
void lo_user_thres (const double &lo_user_thres)
 Set user low energy threshold. More...
 
const double & lo_user_thres (void) const
 Return user low energy threshold. More...
 
void hi_user_thres (const double &hi_user_thres)
 Set user high energy threshold. More...
 
const double & hi_user_thres (void) const
 Return user high energy threshold. More...
 
void n_tels (const int &tels)
 Set number of telescopes. More...
 
const int & n_tels (void) const
 Return number of telescopes. 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 npred (const GModel &model) const
 Return total number of predicted counts for one model. 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 GCTAObservation &obs)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void read_attributes (const GFitsHDU &hdu)
 Read observation attributes. More...
 
void write_attributes (GFitsHDU &hdu) const
 Write observation attributes. 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 bool use_event_for_likelihood (const int &index) const
 Check whether bin should be used for likelihood analysis. 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...
 
GFilename m_eventfile
 Event filename. More...
 
GCTAResponsem_response
 Pointer to instrument response functions. More...
 
GCTAPointing m_pointing
 Pointing direction. More...
 
GSkyRegions m_off_regions
 Off regions. More...
 
double m_ontime
 Ontime (seconds) More...
 
double m_livetime
 Livetime (seconds) More...
 
double m_deadc
 Deadtime correction (livetime/ontime) More...
 
double m_ra_obj
 Right Ascension of object (degrees) More...
 
double m_dec_obj
 Declination of object (degrees) More...
 
double m_lo_user_thres
 User defined lower energy threshold. More...
 
double m_hi_user_thres
 User defined upper energy boundary. More...
 
int m_n_tels
 Number of telescopes. More...
 
std::string m_object
 Object name. More...
 
std::string m_bgdfile
 Background filename. 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
 

Friends

class GCTAModelCubeBackground
 

Detailed Description

CTA observation class.

This class implements a CTA observation.

Definition at line 61 of file GCTAObservation.hpp.

Constructor & Destructor Documentation

GCTAObservation::GCTAObservation ( void  )

Void constructor.

Constructs an empty CTA observation.

Definition at line 92 of file GCTAObservation.cpp.

References init_members().

Referenced by clone().

GCTAObservation::GCTAObservation ( const bool &  dummy,
const std::string &  instrument 
)

Instrument constructor.

Parameters
[in]dummyDummy flag.
[in]instrumentInstrument name.

Constructs an empty CTA observation for a given instrument.

This method is specifically used allocation of global class instances for observation registry. By specifying explicit instrument names it is possible to use the "CTA" module are for other Imaging Air Cherenkov Telescopes. So far, the following instrument codes are supported: "CTA", "HESS", "VERITAS", "MAGIC", "ASTRI".

Definition at line 116 of file GCTAObservation.cpp.

References init_members(), instrument(), and m_instrument.

GCTAObservation::GCTAObservation ( const GFilename filename)
explicit

Unbinned or binned analysis constructor.

Parameters
[in]filenameEvent list or counts cube FITS file name.

Constructs a CTA observation from either an event list of a counts cube.

Definition at line 138 of file GCTAObservation.cpp.

References init_members(), and load().

GCTAObservation::GCTAObservation ( const GFilename cntcube,
const GFilename expcube,
const GFilename psfcube,
const GFilename bkgcube 
)

Stacked analysis constructor.

Parameters
[in]cntcubeCounts cube file name.
[in]expcubeExposure cube file name.
[in]psfcubePoint spread function cube file name.
[in]bkgcubeBackground cube file name.

Constructs a CTA observation from a counts cube, an exposure cube, a point spread function cube, and a background cube.

Definition at line 162 of file GCTAObservation.cpp.

References init_members(), and load().

GCTAObservation::GCTAObservation ( const GFilename cntcube,
const GFilename expcube,
const GFilename psfcube,
const GFilename edispcube,
const GFilename bkgcube 
)

Stacked analysis with energy dispersion constructor.

Parameters
[in]cntcubeCounts cube file name.
[in]expcubeExposure cube file name.
[in]psfcubePoint spread function cube file name.
[in]edispcubeEnergy dispersion cube file name.
[in]bkgcubeBackground cube file name.

Constructs a CTA observation from a counts cube, an exposure cube, a point spread function cube, an energy dispersion cube, and a background cube.

Definition at line 190 of file GCTAObservation.cpp.

References init_members(), and load().

GCTAObservation::GCTAObservation ( const GCTAObservation obs)

Copy constructor.

Parameters
[in]obsCTA observation.

Constructs a CTA observation by copying an existing CTA observation.

Definition at line 214 of file GCTAObservation.cpp.

References copy_members(), and init_members().

GCTAObservation::~GCTAObservation ( void  )
virtual

Destructor.

Destructs CTA observation.

Definition at line 232 of file GCTAObservation.cpp.

References free_members().

Member Function Documentation

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

Return class name.

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

Implements GObservation.

Definition at line 192 of file GCTAObservation.hpp.

void GCTAObservation::clear ( void  )
virtual

Clear CTA observation.

Clear CTA observation.

Implements GObservation.

Definition at line 290 of file GCTAObservation.cpp.

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

Referenced by read().

GCTAObservation * GCTAObservation::clone ( void  ) const
virtual

Clone CTA observation.

Returns
Pointer to deep copy of CTA observation.

Returns a pointer to a deep copy of a CTA observation.

Implements GObservation.

Definition at line 312 of file GCTAObservation.cpp.

References GCTAObservation().

void GCTAObservation::copy_members ( const GCTAObservation obs)
protected

Copy class members.

Parameters
[in]obsCTA observation.

Definition at line 1555 of file GCTAObservation.cpp.

References GCTAResponse::clone(), m_bgdfile, m_deadc, m_dec_obj, m_eventfile, m_hi_user_thres, m_instrument, m_livetime, m_lo_user_thres, m_n_tels, m_object, m_off_regions, m_ontime, m_pointing, m_ra_obj, and m_response.

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

double GCTAObservation::deadc ( const GTime time = GTime()) const
inlinevirtual

Return deadtime correction factor.

Parameters
[in]timeTime (default: GTime()).
Returns
Deadtime correction factor.

Returns the deadtime correction factor. Optionally, this method takes a time argument that takes provision for returning the deadtime correction factor as function of time.

The deadtime correction factor is defined as the livetime divided by the ontime.

Implements GObservation.

Definition at line 261 of file GCTAObservation.hpp.

References m_deadc.

Referenced by deadc(), response(), GCTAOnOffObservation::set(), and write_attributes().

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

Set deadtime correction.

Parameters
[in]deadcDeadtime correction.

Definition at line 448 of file GCTAObservation.hpp.

References deadc(), and m_deadc.

void GCTAObservation::dec_obj ( const double &  dec)
inline

Set CTA object Declination.

Parameters
[in]decObject Declination.

Definition at line 397 of file GCTAObservation.hpp.

References m_dec_obj.

const double & GCTAObservation::dec_obj ( void  ) const
inline

Return CTA object Declination.

Returns
Object Declination.

Definition at line 410 of file GCTAObservation.hpp.

References m_dec_obj.

Referenced by write_attributes().

void GCTAObservation::dispose_events ( void  )

Dispose events.

Disposes the events to save memory. The method only applies to event lists. If does nothing in case that the observation does not contain an event list.

Definition at line 1416 of file GCTAObservation.cpp.

References GCTAEventList::dispose(), and GObservation::m_events.

GEbounds GCTAObservation::ebounds ( void  ) const

Get energy boundaries.

Returns
Energy boundaries.
Exceptions
GException::invalid_valueObservation does not contain events.

Extract the energy boundaries from the events. An exception is thrown if the observation does not contain events.

Definition at line 550 of file GCTAObservation.cpp.

References GEvents::ebounds(), G_EBOUNDS, and GObservation::m_events.

Referenced by GCTAOnOffObservation::apply_ebounds(), GCTACubeBackground::fill(), GCTACubeExposure::fill_cube(), GCTACubePsf::fill_cube(), GCTACubeEdisp::fill_cube(), read(), and GCTAOnOffObservation::set().

void GCTAObservation::eventfile ( const GFilename filename)
inline

Set event file name.

Parameters
[in]filenameEvent file name.

Definition at line 461 of file GCTAObservation.hpp.

References m_eventfile.

const GFilename & GCTAObservation::eventfile ( void  ) const
inline

Return event file name.

Returns
Event file name.

Definition at line 474 of file GCTAObservation.hpp.

References m_eventfile.

Referenced by print().

std::string GCTAObservation::eventtype ( void  ) const

Return event type string.

Returns
Event type string.

Returns "EventList" if the observation contains an event list, "CountsCube" if it contains a counts cube, "Events" if it container an unknown type of events (which should never occur), and an empty string if no events have been allocated.

Definition at line 1375 of file GCTAObservation.cpp.

References GObservation::m_events.

Referenced by GCTACubeBackground::fill(), GCTACubeExposure::fill(), GCTAModelSpatialLookup::fill(), GCTACubeEdisp::fill(), GCTACubePsf::fill(), GCTACubeExposure::fill_cube(), GCTACubePsf::fill_cube(), GCTACubeEdisp::fill_cube(), print(), and write().

void GCTAObservation::free_members ( void  )
protected

Delete class members.

Definition at line 1584 of file GCTAObservation.cpp.

References m_response.

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

GGti GCTAObservation::gti ( void  ) const

Get Good Time Intervals.

Returns
Good Time Intervals.
Exceptions
GException::invalid_valueObservation does not contain events.

Extracts the Good Time Intervals from the events. An exception is thrown if the observation does not contain events.

Definition at line 525 of file GCTAObservation.cpp.

References G_GTI, GEvents::gti(), and GObservation::m_events.

Referenced by GCTACubeBackground::fill(), GCTACubeExposure::fill_cube(), and read().

bool GCTAObservation::has_events ( void  ) const
inline

Signal if CTA observation contains events.

Returns
True if CTA observation contains events.

Definition at line 285 of file GCTAObservation.hpp.

References GFilename::length(), m_eventfile, and GObservation::m_events.

bool GCTAObservation::has_response ( void  ) const
inline

Signal if CTA observation contains response information.

Returns
True if CTA observation contains response information.

Definition at line 273 of file GCTAObservation.hpp.

References GCTAResponse::is_valid(), and m_response.

void GCTAObservation::hi_user_thres ( const double &  hi_user_thres)
inline

Set user high energy threshold.

Parameters
[in]hi_user_thresUser high energy threshold.

Definition at line 511 of file GCTAObservation.hpp.

References hi_user_thres(), and m_hi_user_thres.

const double & GCTAObservation::hi_user_thres ( void  ) const
inline

Return user high energy threshold.

Returns
User high energy threshold.

Definition at line 524 of file GCTAObservation.hpp.

References m_hi_user_thres.

Referenced by hi_user_thres().

void GCTAObservation::init_members ( void  )
protected
std::string GCTAObservation::instrument ( void  ) const
inlinevirtual
void GCTAObservation::instrument ( const std::string &  instrument)
inline

Set instrument name.

Parameters
[in]instrumentInstrument name.

Definition at line 204 of file GCTAObservation.hpp.

References instrument(), and m_instrument.

double GCTAObservation::livetime ( void  ) const
inlinevirtual
void GCTAObservation::livetime ( const double &  livetime)
inline

Set livetime.

Parameters
[in]livetimeLivetime.

Definition at line 435 of file GCTAObservation.hpp.

References livetime(), and m_livetime.

void GCTAObservation::lo_user_thres ( const double &  lo_user_thres)
inline

Set user low energy threshold.

Parameters
[in]lo_user_thresUser low energy threshold.

Definition at line 486 of file GCTAObservation.hpp.

References lo_user_thres(), and m_lo_user_thres.

const double & GCTAObservation::lo_user_thres ( void  ) const
inline

Return user low energy threshold.

Returns
User low energy threshold.

Definition at line 499 of file GCTAObservation.hpp.

References m_lo_user_thres.

Referenced by lo_user_thres().

void GCTAObservation::load ( const GFilename filename)

Load unbinned or binned analysis data.

Parameters
[in]filenameEvent list or counts cube FITS file name.

Loads either an event list or a counts cube from a FITS file.

Definition at line 1187 of file GCTAObservation.cpp.

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

Referenced by GCTAObservation(), load(), and read().

void GCTAObservation::load ( const GFilename cntcube,
const GFilename expcube,
const GFilename psfcube,
const GFilename bkgcube 
)

Load stacked analysis data.

Parameters
[in]cntcubeCounts cube file name.
[in]expcubeExposure cube file name.
[in]psfcubePoint spread function cube file name.
[in]bkgcubeBackground cube file name.

Loads a counts cube, an exposure cube, a point spread function cube, and a background cube for stacked analysis.

Definition at line 1220 of file GCTAObservation.cpp.

References GObservation::events(), G_LOAD1, load(), response(), and GFilename::url().

void GCTAObservation::load ( const GFilename cntcube,
const GFilename expcube,
const GFilename psfcube,
const GFilename edispcube,
const GFilename bkgcube 
)

Load stacked analysis data with energy dispersion.

Parameters
[in]cntcubeCounts cube file name.
[in]expcubeExposure cube file name.
[in]psfcubePoint spread function cube file name.
[in]edispcubeEnergy dispersion cube file name.
[in]bkgcubeBackground cube file name.

Loads a counts cube, an exposure cube, a point spread function cube, a energy dispersion cube and a background cube for stacked analysis.

Definition at line 1265 of file GCTAObservation.cpp.

References GObservation::events(), G_LOAD2, load(), response(), and GFilename::url().

void GCTAObservation::n_tels ( const int &  tels)
inline

Set number of telescopes.

Parameters
[in]telsNumber of telescopes.

Definition at line 536 of file GCTAObservation.hpp.

References m_n_tels.

const int & GCTAObservation::n_tels ( void  ) const
inline

Return number of telescopes.

Returns
Number of telescopes.

Definition at line 549 of file GCTAObservation.hpp.

References m_n_tels.

Referenced by write_attributes().

void GCTAObservation::object ( const std::string &  object)
inline

Set CTA observation object name.

Parameters
[in]objectObject name.

Definition at line 347 of file GCTAObservation.hpp.

References m_object, and object().

const std::string & GCTAObservation::object ( void  ) const
inline

Return CTA observation object.

Returns
Object name.

Definition at line 360 of file GCTAObservation.hpp.

References m_object.

Referenced by object().

void GCTAObservation::off_regions ( const GSkyRegions off_regions)
inline

Set sky off regions.

Parameters
[in]off_regionsSky off regions.

Definition at line 334 of file GCTAObservation.hpp.

References m_off_regions, and off_regions().

const GSkyRegions & GCTAObservation::off_regions ( void  ) const
inline

Return sky off regions.

Returns
Sky off regions

Definition at line 322 of file GCTAObservation.hpp.

References m_off_regions.

Referenced by off_regions(), and print().

double GCTAObservation::ontime ( void  ) const
inlinevirtual

Return ontime.

Returns
Ontime in seconds.

Implements GObservation.

Definition at line 229 of file GCTAObservation.hpp.

References m_ontime.

Referenced by GCTACubeBackground::fill(), ontime(), print(), read(), response(), GCTAOnOffObservation::set(), and write_attributes().

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

Set ontime.

Parameters
[in]ontimeOntime.

Definition at line 422 of file GCTAObservation.hpp.

References m_ontime, and ontime().

GCTAObservation & GCTAObservation::operator= ( const GCTAObservation obs)

Assignment operator.

Parameters
[in]obsCTA observation.

Assign CTA observation to this object.

Definition at line 255 of file GCTAObservation.cpp.

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

const GCTAPointing & GCTAObservation::pointing ( void  ) const
inline

Return CTA pointing.

Returns
CTA pointing

Definition at line 297 of file GCTAObservation.hpp.

References m_pointing.

Referenced by pointing(), and print().

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

Print CTA observation information.

Parameters
[in]chatterChattiness.
Returns
String containing observation information.

Implements GObservation.

Definition at line 1437 of file GCTAObservation.cpp.

References eventfile(), eventtype(), instrument(), livetime(), m_deadc, GObservation::m_events, m_hi_user_thres, m_lo_user_thres, m_response, GObservation::name(), off_regions(), ontime(), gammalib::parformat(), pointing(), GCTAResponse::print(), GEvents::print(), gammalib::reduce(), SILENT, GObservation::statistic(), and gammalib::str().

void GCTAObservation::ra_obj ( const double &  ra)
inline

Set CTA object Right Ascension.

Parameters
[in]raObject Right Ascension.

Definition at line 372 of file GCTAObservation.hpp.

References m_ra_obj.

const double & GCTAObservation::ra_obj ( void  ) const
inline

Return CTA object Right Ascension.

Returns
Object Right Ascension.

Definition at line 385 of file GCTAObservation.hpp.

References m_ra_obj.

Referenced by write_attributes().

void GCTAObservation::read ( const GXmlElement xml)
virtual

Read observation from XML element.

Parameters
[in]xmlXML element.
Exceptions
GException::invalid_valueInvalid number of parameters found in XML element.
GException::xml_invalid_parnamesInvalid parameter names found in XML element.
GException::invalid_valueInvalid statistic attribute encountered

Reads information for a CTA observation from an XML element. This method handles two variants: a first where an event list of counts cube is given and a second where the observation definition information is provided by parameter tags.

The XML format for an event list is

<observation name="..." id="..." instrument="...">
  <parameter name="EventList" file="..."/>
  ...
</observation>

and for a counts cube is

<observation name="..." id="..." instrument="...">
  <parameter name="CountsCube" file="..."/>
  ...
</observation>

The second variant without event information has the following format

<observation name="..." id="..." instrument="...">
  <parameter name="Pointing" ra="..." dec="..."/>
  <parameter name="Deadtime" deadc="..."/>
  ...
</observation>

In addition, calibration information can be specified using the format

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

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

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

In case that a CountsCube is handled, the stacked response can also be provided in the format

 <observation name="..." id="..." instrument="...">
   ...
   <parameter name="ExposureCube" file="..."/>
   <parameter name="PsfCube"      file="..."/>
   <parameter name="EdispCube"    file="..."/>
   <parameter name="BkgCube"      file="..."/>
 </observation>

Note that the EdispCube is an optional parameter.

Optional user energy thresholds can be specified by adding the emin and emax attributes to the observation tag. Also the statistic used for maximum likelihood fitting can be specified:

<observation name="..." id="..." instrument="..." emin="..." emax="..." statistic="...">

The method does not load the events into memory but stores the file name of the event file. The events are only loaded when required. This reduces the memory needs for an CTA observation object and allows for loading of event information upon need.

Implements GObservation.

Definition at line 645 of file GCTAObservation.cpp.

References GXmlElement::attribute(), clear(), GEvents::ebounds(), ebounds(), GXmlNode::element(), GXmlNode::elements(), GObservation::events(), G_READ, GEvents::gti(), gti(), livetime(), load(), m_bgdfile, m_deadc, m_hi_user_thres, m_instrument, m_lo_user_thres, m_off_regions, m_pointing, m_response, ontime(), GGti::ontime(), GCTARoi::read(), GCTAPointing::read(), GCTAResponse::read(), GGti::read(), GEbounds::read(), GCTAEventList::roi(), roi(), GObservation::statistic(), gammalib::str(), gammalib::todouble(), gammalib::toupper(), gammalib::xml_check_par(), gammalib::xml_file_expand(), gammalib::xml_get_par(), and gammalib::xml_has_par().

Referenced by load().

void GCTAObservation::read ( const GFits fits)

Read data from FITS file.

Parameters
[in]fitsFITS file.

Reads event data from a FITS file and sets the observation attributes.

The method automatically switches between an event list and a counts cube, depending of the information provided in the FITS file. If an extension name is specified, the method checks whether the extension exists and loads the extension as event list. Otherwise, it checks whether the file contains an EVENTS extension and loads the extension as event list. If none of the above are satistified, the method loads a counts cube.

Definition at line 1058 of file GCTAObservation.cpp.

References GFits::at(), GFits::contains(), GObservation::events(), GFilename::extname(), gammalib::extname_cta_events, GFits::filename(), GObservation::m_events, GCTAEventList::read(), GCTAEventCube::read(), and read_attributes().

void GCTAObservation::read_attributes ( const GFitsHDU hdu)
protected

Read observation attributes.

Parameters
[in]hduFITS HDU.

Reads CTA observation attributes from HDU. Mandatory attributes are

RA_PNT   - Right Ascension of pointing
DEC_PNT  - Declination of pointing
ONTIME   - Exposure time
LIVETIME - Livetime

and optional attributes are

OBJECT   - Name of observed object
DEADC    - Deadtime correction
RA_OBJ   - Right Ascension of observed object,
DEC_OBJ  - Declination of observed object,
OBS_ID   - Observation identifier
ALT_PNT  - Altitude of pointing above horizon
AZ_PNT   - Azimuth of pointing
TELESCOP - Telescope name
N_TELS   - Number of telescopes

Based on RA_PNT and DEC_PNT, the CTA pointing direction is set. Note that DEADC is computed using DEADC=LIVETIME/ONTIME

Todo:
The actual reader is a minimal reader to accomodate as many different datasets as possible. Once the CTA data format is fixed the reader should have more mandatory attributes.

Definition at line 1628 of file GCTAObservation.cpp.

References GCTAPointing::azimuth(), GCTAPointing::dir(), GFitsHDU::has_card(), GObservation::id(), GFitsHDU::integer(), m_deadc, m_dec_obj, m_instrument, m_livetime, m_n_tels, m_object, m_ontime, m_pointing, m_ra_obj, GSkyDir::radec_deg(), GFitsHDU::real(), GFitsHDU::string(), and GCTAPointing::zenith().

Referenced by read().

void GCTAObservation::response ( const GResponse rsp)
virtual

Set response function.

Parameters
[in]rspResponse function.

Sets the response function for the observation.

Implements GObservation.

Definition at line 325 of file GCTAObservation.cpp.

References GCTAResponse::clone(), gammalib::cta_rsp(), G_RESPONSE_SET, and m_response.

Referenced by GCTAOnOffObservation::arf_rad_max(), gammalib::cta_rsp_cube(), gammalib::cta_rsp_irf(), GCTACubeExposure::fill_cube(), GCTACubePsf::fill_cube(), GCTACubeEdisp::fill_cube(), and GCTAOnOffObservation::set().

const GCTAResponse * GCTAObservation::response ( void  ) const
virtual

Return pointer to CTA response function.

Returns
Pointer to CTA response function.
Exceptions
GException::invalid_valueNo valid response found in CTA observation.

Returns a pointer to the CTA response function. An exception is thrown if the pointer is not valid, hence the user does not need to verify the validity of the pointer.

Implements GObservation.

Definition at line 356 of file GCTAObservation.cpp.

References G_RESPONSE_GET, and m_response.

Referenced by load().

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

Set CTA response function.

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

Sets the CTA response function by specifying a response name and a calibration database. This method also loads the response function so that it is available for data analysis.

Definition at line 380 of file GCTAObservation.cpp.

References GCTAResponseIrf::caldb(), GCTAResponseIrf::load(), and m_response.

void GCTAObservation::response ( const GCTACubeExposure expcube,
const GCTACubePsf psfcube,
const GCTACubeBackground bkgcube 
)

Set CTA response function.

Parameters
[in]expcubeExposure cube.
[in]psfcubePsf cube.
[in]bkgcubeBackground cube.

Sets the CTA response function fur cube analysis by specifying the exposure cube, the Psf cube and the background cube. The method also copies over the ontime, the livetime and the deadtime correction factor from the exposure cube.

Definition at line 415 of file GCTAObservation.cpp.

References GCTACubeExposure::deadc(), deadc(), GCTACubeExposure::livetime(), livetime(), m_response, GCTACubeExposure::ontime(), and ontime().

void GCTAObservation::response ( const GCTACubeExposure expcube,
const GCTACubePsf psfcube,
const GCTACubeEdisp edispcube,
const GCTACubeBackground bkgcube 
)

Set CTA response function.

Parameters
[in]expcubeExposure cube.
[in]psfcubePsf cube.
[in]edispcubeEdisp cube.
[in]bkgcubeBackground cube.

Sets the CTA response function fur cube analysis by specifying the exposure cube, the Psf cube, the exposure cube and the background cube. The method also copies over the ontime, the livetime and the deadtime correction factor from the exposure cube.

Definition at line 452 of file GCTAObservation.cpp.

References GCTACubeExposure::deadc(), deadc(), GCTACubeExposure::livetime(), livetime(), m_response, GCTACubeExposure::ontime(), and ontime().

GCTARoi GCTAObservation::roi ( void  ) const

Get Region of Interest.

Returns
Region of Interest.
Exceptions
GException::invalid_valueObservation does not contain events. Observation does not contain an event list.

Extracts the Region of Interest from the event list. An exception is thrown if the observation does not contain an event list.

Definition at line 492 of file GCTAObservation.cpp.

References G_ROI, GObservation::m_events, and GCTAEventList::roi().

Referenced by GCTACubeBackground::fill(), GCTACubeExposure::fill_cube(), GCTACubePsf::fill_cube(), GCTACubeEdisp::fill_cube(), GCTAModelSpatial::mc(), GCTAModelSpatialGradient::mc_max_value(), and read().

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

Save CTA observation into FITS file.

Parameters
[in]filenameFITS filename.
[in]clobberOverwrite existing FITS file? (default: false)

Saves the CTA observation into a FITS file.

If the CTA observation contains an event list, the method will write an events and a Good Time Intervals extension to the FITS file. The names for both extension can be optionally specified in the filename using the format

 filename[EVENTS;GTI]

where the string in the squared bracket defines the extension names. The part before the semi-colon is the events extension name and the part after the semi-colon is the Good Time Intervals extension name.

If the CTA observation contains an event cube, the method will write the cube into the primary image, followed by binary tables containing the energy boundaries and the Good Time Intervals. The extension names of these binary tables are EBOUNDS and GTI, and cannot be modified.

Definition at line 1326 of file GCTAObservation.cpp.

References GFilename::extname(), gammalib::extname_cta_events, gammalib::extname_gti, GFilename::has_extname(), gammalib::split(), gammalib::strip_whitespace(), GFilename::url(), and write().

void GCTAObservation::write ( GXmlElement xml) const
virtual

Write observation into XML element.

Parameters
[in]xmlXML element.
Exceptions
GException::invalid_valueNo valid events found in observation.

Writes information for a CTA observation into an XML element.

Dependent on the content of the CTA observation, different XML formats will be written. If the CTA observation contains an event list that has been loaded from disk, the method will write the file name of the event list using the format

<observation name="..." id="..." instrument="..." statistic="...">
  <parameter name="EventList" file="..."/>
  ...
</observation>

If the CTA observation contains a counts cube that has been loaded from disk, the method will write the file name of the counts cube using the format

<observation name="..." id="..." instrument="..." statistic="...">
  <parameter name="CountsCube" file="..."/>
  ...
</observation>

In all other cases the method will only write the metadata information for the CTA observation using the format

<observation name="..." id="..." instrument="..." statistic="...">
  <parameter name="Pointing" ra="..." dec="..."/>
  <parameter name="Deadtime" deadc="..."/>
  ...
</observation>

In case that response information is available and that the response information has been either loaded or saved to disk, the method will write the reponse information into the XML file. If the response for an unbinned or binned observation has been loaded from the calibration database, the format

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

will be written. If the response has been loaded from individual files, the format

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

will be written. If the observation contains information about the response to a stacked analysis, the format

<observation name="..." id="..." instrument="..." statistic="...">
  ...
  <parameter name="ExposureCube" file="..."/>
  <parameter name="PsfCube"      file="..."/>
  <parameter name="EdispCube"    file="..."/>
  <parameter name="BkgCube"      file="..."/>
</observation>

will be written. Note that the EdispCube parameter will only be written in case that an energy dispersion cube had been defined (energy dispersion is optional).

If user energy thresholds are defined (i.e. threshold values are >0) the additional emin and emax attributes will be written to the observation tag:

<observation name="..." id="..." instrument="..." statistic="..." emin="..." emax="...">

Implements GObservation.

Definition at line 966 of file GCTAObservation.cpp.

References GXmlElement::attribute(), GEvents::ebounds(), GObservation::events(), eventtype(), GSkyRegions::filename(), G_WRITE, GEvents::gti(), GFilename::is_empty(), m_deadc, m_eventfile, m_hi_user_thres, m_lo_user_thres, m_off_regions, m_pointing, m_response, GCTAEventList::roi(), GObservation::statistic(), gammalib::str(), GFilename::url(), GCTAPointing::write(), GCTAResponse::write(), gammalib::xml_file_reduce(), and gammalib::xml_need_par().

Referenced by save().

void GCTAObservation::write ( GFits fits,
const std::string &  evtname = gammalib::extname_cta_events,
const std::string &  gtiname = gammalib::extname_gti 
) const

Write observation into FITS file.

Parameters
[in]fitsFITS file.
[in]evtnameEvents FITS extension name.
[in]gtinameGood Time Intervals FITS extension name.

Writes the observation into a FITS file.

If the observation contains an event list, the event list and Good Time Intervals will be written to the FITS file. The FITS extension name for the event list and the Good Time Intervals can be optionally specified thru the evtname and gtiname arguments.

If the observation contains an event cube, the event cube will be written into the primary extension of the FITS file.

This method does nothing if no events are in the CTA observation.

Definition at line 1132 of file GCTAObservation.cpp.

References GFits::at(), GFits::contains(), GObservation::events(), GFits::remove(), GCTAEventList::write(), GCTAEventCube::write(), and write_attributes().

Friends And Related Function Documentation

friend class GCTAModelCubeBackground
friend

Definition at line 64 of file GCTAObservation.hpp.

Member Data Documentation

std::string GCTAObservation::m_bgdfile
protected

Background filename.

Definition at line 182 of file GCTAObservation.hpp.

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

double GCTAObservation::m_deadc
protected

Deadtime correction (livetime/ontime)

Definition at line 172 of file GCTAObservation.hpp.

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

double GCTAObservation::m_dec_obj
protected

Declination of object (degrees)

Definition at line 174 of file GCTAObservation.hpp.

Referenced by copy_members(), dec_obj(), init_members(), and read_attributes().

GFilename GCTAObservation::m_eventfile
protected

Event filename.

Definition at line 166 of file GCTAObservation.hpp.

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

double GCTAObservation::m_hi_user_thres
protected

User defined upper energy boundary.

Definition at line 176 of file GCTAObservation.hpp.

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

std::string GCTAObservation::m_instrument
protected

Instrument name.

Definition at line 165 of file GCTAObservation.hpp.

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

double GCTAObservation::m_livetime
protected

Livetime (seconds)

Definition at line 171 of file GCTAObservation.hpp.

Referenced by copy_members(), init_members(), livetime(), and read_attributes().

double GCTAObservation::m_lo_user_thres
protected

User defined lower energy threshold.

Definition at line 175 of file GCTAObservation.hpp.

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

int GCTAObservation::m_n_tels
protected

Number of telescopes.

Definition at line 177 of file GCTAObservation.hpp.

Referenced by copy_members(), init_members(), n_tels(), and read_attributes().

std::string GCTAObservation::m_object
protected

Object name.

Definition at line 181 of file GCTAObservation.hpp.

Referenced by copy_members(), init_members(), object(), and read_attributes().

GSkyRegions GCTAObservation::m_off_regions
protected

Off regions.

Definition at line 169 of file GCTAObservation.hpp.

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

double GCTAObservation::m_ontime
protected

Ontime (seconds)

Definition at line 170 of file GCTAObservation.hpp.

Referenced by copy_members(), init_members(), ontime(), and read_attributes().

GCTAPointing GCTAObservation::m_pointing
protected

Pointing direction.

Definition at line 168 of file GCTAObservation.hpp.

Referenced by copy_members(), init_members(), pointing(), read(), read_attributes(), write(), and write_attributes().

double GCTAObservation::m_ra_obj
protected

Right Ascension of object (degrees)

Definition at line 173 of file GCTAObservation.hpp.

Referenced by copy_members(), init_members(), ra_obj(), and read_attributes().

GCTAResponse* GCTAObservation::m_response
protected

Pointer to instrument response functions.

Definition at line 167 of file GCTAObservation.hpp.

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


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