GammaLib  2.0.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GObservations Class Reference

Observation container class. More...

#include <GObservations.hpp>

Inheritance diagram for GObservations:
GContainer GBase

Classes

class  likelihood
 

Public Member Functions

 GObservations (void)
 Void constructor. More...
 
 GObservations (const GObservations &obs)
 Copy constructor. More...
 
 GObservations (const GFilename &filename)
 Load constructor. More...
 
virtual ~GObservations (void)
 Destructor. More...
 
GObservationsoperator= (const GObservations &obs)
 Assignment operator. More...
 
GObservationoperator[] (const int &index)
 Return pointer to observation. More...
 
const GObservationoperator[] (const int &index) const
 Return pointer to observation (const version) More...
 
void clear (void)
 Clear observations. More...
 
GObservationsclone (void) const
 Clone observations. More...
 
std::string classname (void) const
 Return class name. More...
 
GObservationat (const int &index)
 Return pointer to observation. More...
 
const GObservationat (const int &index) const
 Return pointer to observation (const version) More...
 
int size (void) const
 Return number of observations in container. More...
 
bool is_empty (void) const
 Signals if there are no observations in container. More...
 
GObservationset (const int &index, const GObservation &obs)
 Set observation in container. More...
 
GObservationappend (const GObservation &obs)
 Append observation to container. More...
 
GObservationinsert (const int &index, const GObservation &obs)
 Insert observation into container. More...
 
void remove (const int &index)
 Remove observation from container. More...
 
void reserve (const int &num)
 Reserves space for observations in container. More...
 
void extend (const GObservations &obs)
 Append observations from observation container. More...
 
bool contains (const std::string &instrument, const std::string &id) const
 Signals if observation exists. More...
 
void load (const GFilename &filename)
 Load observations from XML file. More...
 
void save (const GFilename &filename) const
 Save observations into XML file. More...
 
void read (const GXml &xml)
 Read observations from XML document. More...
 
void write (GXml &xml) const
 Write observations into XML document. More...
 
void models (const GModels &models)
 Set model container. More...
 
void models (const GFilename &filename)
 Load models from XML file. More...
 
const GModelsmodels (void) const
 Return model container. More...
 
void optimize (GOptimizer &opt)
 Optimize model parameters using optimizer. More...
 
void errors (GOptimizer &opt)
 Computes parameter errors using optimizer. More...
 
void errors_hessian (void)
 Computes parameter errors using hessian matrix and optimizer. More...
 
void eval (void)
 Evaluate function. More...
 
double logL (void) const
 Return log-likelihood of models. More...
 
int nobserved (void) const
 Return total number of observed events. More...
 
double npred (void) const
 Return total number of predicted events in models. More...
 
double npred (const std::string &name, const std::string &instrument="") const
 Return number of predicted events for model with a given name. More...
 
double npred (const GModel &model, const std::string &instrument="") const
 Return number of predicted events for a given model. More...
 
void remove_response_cache (void)
 Remove response cache for all models. More...
 
void remove_response_cache (const std::string &name)
 Remove response cache for model. More...
 
std::string print (const GChatter &chatter=NORMAL) const
 Print observation list information. More...
 
const GObservations::likelihoodfunction (void) const
 Return likelihood function. More...
 
- Public Member Functions inherited from GContainer
virtual ~GContainer (void)
 Destructor. 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 GObservations &obs)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
int get_index (const std::string &instrument, const std::string &id) const
 Return observation index by instrument and identifier. More...
 

Protected Attributes

std::vector< GObservation * > m_obs
 List of observations. More...
 
GModels m_models
 List of models. More...
 
GObservations::likelihood m_fct
 Optimizer function. More...
 

Detailed Description

Observation container class.

This is the main user interface class that provides high-level access to gamma-ray observations and manages their scientific analysis. For a given instruments, the identifiers of all observations in the container have to be unique, i.e. every identifier can occur only once. This allows for accessing the observations by instrument and by identifier.

GObservations holds a list of pointers to GObservation objects which implement gamma-ray observations. The class derives from the abstract container interface class GContainer, and provides methods to manage the list of observations. The append() appends an observation to the container, the insert() method inserts an observation at a specific position in the list (though the order of the observations in the list is not relevant for the scientific analysis). The remove() method removes a specific observation from the list and the extend() method extends the list of observations by appending another list of observations to it. The size() method provides the number of observations in the list.

The operator[] provides access to an observation in the list by index, returning a reference to the observation. The operator does not perform index range checking. If range checking is required, use the at() method.

The observation information may be stored into a XML file. The save() method saves all information into a XML file while the load() method loads the information into the list. Alternatively, a XML file constructor can be used to construct an instance of GObservations by loading observation information from a XML file. In addition, the read() and write() method handle the reading and writing of observation information from and into a XML object of type GXml.

The class also holds a list of models that is implemented by the GModels class. The models() method allows setting and retrieving the list of models, as well as loading of models from an XML file.

Based on the list of models, the optimize() method will optimize all model parameters that are marked as free in the list of models. The npred() method returns the total number of events that are predicted by all models after optimization.

GObservations also provides an optimizer class that is derived from the abstract GOptimizerFunction base class. The GObservations::optimizer class is the object that is used for model parameter optimization.

Definition at line 92 of file GObservations.hpp.

Constructor & Destructor Documentation

GObservations::GObservations ( void  )

Void constructor.

Definition at line 67 of file GObservations.cpp.

References init_members().

Referenced by clone().

GObservations::GObservations ( const GObservations obs)

Copy constructor.

Parameters
obsObservation container.

Definition at line 82 of file GObservations.cpp.

References copy_members(), and init_members().

GObservations::GObservations ( const GFilename filename)
explicit

Load constructor.

Parameters
[in]filenameXML filename.

Construct the observation container by loading all relevant information from a XML file. Please refer to the read() method for more information about the structure of the XML file.

Definition at line 104 of file GObservations.cpp.

References init_members(), and load().

GObservations::~GObservations ( void  )
virtual

Destructor.

Definition at line 120 of file GObservations.cpp.

References free_members().

Member Function Documentation

GObservation * GObservations::append ( const GObservation obs)

Append observation to container.

Parameters
[in]obsObservation.
Returns
Pointer to deep copy of observation.
Exceptions
GException::invalid_valueObservation with same instrument and identifier already exists in container.

Appends a deep copy of an observation to the container.

Definition at line 309 of file GObservations.cpp.

References GObservation::clone(), G_APPEND, get_index(), GObservation::id(), GObservation::instrument(), m_obs, and gammalib::str().

Referenced by read().

GObservation * GObservations::at ( const int &  index)

Return pointer to observation.

Parameters
[in]indexObservation index [0,...,size()[.
Returns
Pointer to observation.
Exceptions
GException::out_of_rangeOperation index is out of range.

Returns a pointer to the observation with specified index.

Definition at line 208 of file GObservations.cpp.

References G_AT, m_obs, and size().

const GObservation * GObservations::at ( const int &  index) const

Return pointer to observation (const version)

Parameters
[in]indexObservation index [0,...,size()[.
Returns
Pointer to observation.
Exceptions
GException::out_of_rangeOperation index is out of range.

Returns a const pointer to the observation with specified index.

Definition at line 232 of file GObservations.cpp.

References G_AT, m_obs, and size().

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

Return class name.

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

Implements GBase.

Definition at line 210 of file GObservations.hpp.

void GObservations::clear ( void  )
virtual

Clear observations.

Implements GBase.

Definition at line 172 of file GObservations.cpp.

References free_members(), and init_members().

Referenced by load().

GObservations * GObservations::clone ( void  ) const
virtual

Clone observations.

Returns
Pointer to deep copy of observation container.

Implements GBase.

Definition at line 190 of file GObservations.cpp.

References GObservations().

Referenced by extend().

bool GObservations::contains ( const std::string &  instrument,
const std::string &  id 
) const

Signals if observation exists.

Parameters
[in]instrumentInstrument.
[in]idObservation identifier.
Returns
True if observation with specified instrument and identifier id exists.

Searches all observations for a match with the specified instrument and identifier. If the specified attributes have been found, true is returned.

Definition at line 493 of file GObservations.cpp.

References get_index().

void GObservations::copy_members ( const GObservations obs)
protected

Copy class members.

Parameters
[in]obsObservation container.

Copies all class members. Deep copies are created for the observations, which is what is needed to have a real physical copy of the members.

Definition at line 1086 of file GObservations.cpp.

References m_fct, m_models, m_obs, and GObservations::likelihood::set().

Referenced by GObservations(), operator=(), and GObservations::likelihood::operator=().

void GObservations::errors ( GOptimizer opt)

Computes parameter errors using optimizer.

Parameters
[in]optOptimizer.

Calculates the errors of the free parameters of the models by using the optimizer that has been provided by the opt argument.

Definition at line 749 of file GObservations.cpp.

References GOptimizer::errors(), m_fct, m_models, and GModels::pars().

void GObservations::errors_hessian ( void  )

Computes parameter errors using hessian matrix and optimizer.

Calculates the errors of the free model parameters as the square roots of the diagonal elements of the inverse Hessian matrix. The Hessian matrix is computed using finite differences.

Definition at line 769 of file GObservations.cpp.

References GMatrixSparse::cholesky_decompose(), GMatrixSparse::cholesky_solver(), GObservations::likelihood::hessian(), m_fct, m_models, GModels::pars(), GOptimizerPars::size(), and sqrt().

void GObservations::eval ( void  )

Evaluate function.

Evaluates the likelihood funtion at the actual set of parameters.

Definition at line 853 of file GObservations.cpp.

References GObservations::likelihood::eval(), m_fct, m_models, and GModels::pars().

Referenced by GObservations::likelihood::hessian().

void GObservations::extend ( const GObservations obs)

Append observations from observation container.

Parameters
[in]obsObservations.
Exceptions
GException::invalid_valueObservation with same instrument and identifier already exists in container.

Appends deep copies of observations to the container.

Definition at line 438 of file GObservations.cpp.

References clone(), G_EXTEND, get_index(), is_empty(), m_obs, reserve(), size(), and gammalib::str().

void GObservations::free_members ( void  )
protected

Delete class members.

Deallocates all observations. Since container classes that hold pointers need to handle the proper deallocation of memory, we loop here over all pointers and make sure that we deallocate all observations in the container.

Definition at line 1115 of file GObservations.cpp.

References m_obs.

Referenced by clear(), operator=(), GObservations::likelihood::operator=(), ~GObservations(), and GObservations::likelihood::~likelihood().

const GObservations::likelihood & GObservations::function ( void  ) const
inline

Return likelihood function.

Returns
Reference to likelihood function.

Returns a reference to the likelihood function.

Definition at line 357 of file GObservations.hpp.

References m_fct.

int GObservations::get_index ( const std::string &  instrument,
const std::string &  id 
) const
protected

Return observation index by instrument and identifier.

Parameters
[in]instrumentInstrument.
[in]idObservation identifier.
Returns
Observation index (-1 of instrument and observation identifier has not been found)

Returns observation index based on the specified instrument and observation identifier id. If no observation with the specified attributes has been found, the method returns -1.

Definition at line 1140 of file GObservations.cpp.

References m_obs, and size().

Referenced by append(), contains(), extend(), insert(), and set().

void GObservations::init_members ( void  )
protected

Initialise class members.

< Makes sure that optimizer points to this instance

Definition at line 1066 of file GObservations.cpp.

References GModels::clear(), m_fct, m_models, m_obs, and GObservations::likelihood::set().

Referenced by clear(), GObservations(), operator=(), and GObservations::likelihood::operator=().

GObservation * GObservations::insert ( const int &  index,
const GObservation obs 
)

Insert observation into container.

Parameters
[in]indexObservation index [0,...,size()[.
[in]obsObservation.
Returns
Pointer to deep copy of observation.
Exceptions
GException::out_of_rangeObservation index is out of range.
GException::invalid_valueObservation with same instrument and identifier already exists in container.

Inserts a deep copy of an observation into the container before the observation with the specified index.

Definition at line 352 of file GObservations.cpp.

References GObservation::clone(), G_INSERT, get_index(), GObservation::id(), GObservation::instrument(), is_empty(), m_obs, size(), and gammalib::str().

bool GObservations::is_empty ( void  ) const
inlinevirtual

Signals if there are no observations in container.

Returns
True if container is empty, false otherwise.

Signals if the observation container does not contain any observation.

Implements GContainer.

Definition at line 267 of file GObservations.hpp.

References m_obs.

Referenced by extend(), and insert().

void GObservations::load ( const GFilename filename)

Load observations from XML file.

Parameters
[in]filenameXML filename.

Loads observation from a XML file into the container. Please refer to the read() method for more information about the structure of the XML file.

Definition at line 512 of file GObservations.cpp.

References clear(), read(), and GFilename::url().

Referenced by GObservations().

double GObservations::logL ( void  ) const
inline

Return log-likelihood of models.

Returns
Log-likelihood of models.

Returns the log-likelihood of models. This value if computed following a call of the optimize() or eval() methods. If optimize() or eval() have never been called, the method returns 0.

Definition at line 327 of file GObservations.hpp.

References m_fct, and GObservations::likelihood::value().

void GObservations::models ( const GModels models)
inline

Set model container.

Parameters
[in]modelsModel container.

Sets the model container for the observations.

Definition at line 296 of file GObservations.hpp.

References m_models, and models().

Referenced by GCTACubeBackground::fill().

void GObservations::models ( const GFilename filename)

Load models from XML file.

Parameters
[in]filenameXML filename.

Loads the models from a XML file. Please refer to the GModels::read() method for more information about the expected structure of the XML file.

Definition at line 710 of file GObservations.cpp.

References GModels::load(), and m_models.

const GModels & GObservations::models ( void  ) const
inline

Return model container.

Returns
Model container.

Returns the model container of the observations.

Definition at line 311 of file GObservations.hpp.

References m_models.

Referenced by models().

int GObservations::nobserved ( void  ) const

Return total number of observed events.

Returns
Total number of observed events.

Returns the total number of observed events that is container in the observation container.

Definition at line 874 of file GObservations.cpp.

References size().

Referenced by print().

double GObservations::npred ( void  ) const
inline

Return total number of predicted events in models.

Returns
Total number of predicted events in models.

Returns the total number of events that is predicted by the models. This number if computed following a call of the optimize() or eval() methods. If optimize() or eval() have never been called, the method returns 0.

Definition at line 343 of file GObservations.hpp.

References m_fct, and GObservations::likelihood::npred().

Referenced by npred(), and print().

double GObservations::npred ( const std::string &  name,
const std::string &  instrument = "" 
) const

Return number of predicted events for model with a given name.

Parameters
[in]nameModel name.
[in]instrumentInstrument name.
Returns
Total number of predicted events for model.
Exceptions
GException::invalid_argumentModel with name not found.

Returns the total number of predicted events for the model with a given name.

Definition at line 902 of file GObservations.cpp.

References GModels::contains(), G_NPRED, m_models, and npred().

double GObservations::npred ( const GModel model,
const std::string &  instrument = "" 
) const

Return number of predicted events for a given model.

Parameters
[in]modelModel.
[in]instrumentInstrument name.
Returns
Total number of predicted events for model.
Exceptions
GException::invalid_argumentModel with name not found.

Returns the total number of predicted events for the model with a given name.

Definition at line 941 of file GObservations.cpp.

References GModel::is_valid(), npred(), and size().

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

Assignment operator.

Parameters
[in]obsObservation container.
Returns
Observation container.

Definition at line 142 of file GObservations.cpp.

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

GObservation * GObservations::operator[] ( const int &  index)
inline

Return pointer to observation.

Parameters
[in]indexObservation index [0,...,size()-1].
Returns
Observation.

Returns a pointer to the observation with specified index.

Definition at line 225 of file GObservations.hpp.

References m_obs.

const GObservation * GObservations::operator[] ( const int &  index) const
inline

Return pointer to observation (const version)

Parameters
[in]indexObservation index [0,...,size()-1].

Returns a const pointer to the observation with specified index.

Definition at line 239 of file GObservations.hpp.

References m_obs.

void GObservations::optimize ( GOptimizer opt)

Optimize model parameters using optimizer.

Parameters
[in]optOptimizer.

Optimizes the free parameters of the models by using the optimizer that has been provided by the opt argument.

Definition at line 728 of file GObservations.cpp.

References m_fct, m_models, GOptimizer::optimize(), and GModels::pars().

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

Print observation list information.

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

Implements GBase.

Definition at line 1013 of file GObservations.cpp.

References EXPLICIT, m_models, nobserved(), npred(), gammalib::parformat(), GModels::print(), gammalib::reduce(), SILENT, size(), GModels::size(), gammalib::str(), and VERBOSE.

void GObservations::read ( const GXml xml)

Read observations from XML document.

Parameters
[in]xmlXML document.
Exceptions
GException::invalid_valueInvalid instrument encountered in XML file.

Reads observations from the first observation list that is found in the XML document. The decoding of the instrument specific observation definition is done within the observation's GObservation::read() method. The following file structure is expected:

<observation_list title="observation library">
  <observation name="..." id="..." instrument="...">
    ...
  </observation>
  <observation name="..." id="..." instrument="...">
    ...
  </observation>
  ...
</observation_list>

The name and id attributes allow for a unique identification of an observation within the observation container. The instrument attributes specifies the instrument to which the observation applies. This attribute will be used to allocate the appropriate instrument specific GObservation class variant by using the GObservationRegistry class.

The structure within the observation tag is defined by the instrument specific GObservation class.

Todo:
Observation names and IDs are not verified so far for uniqueness. This would be required to achieve an unambiguous update of parameters in an already existing XML file when using the write method.

Definition at line 589 of file GObservations.cpp.

References GObservationRegistry::alloc(), append(), GXmlElement::attribute(), GRegistry::content(), GXmlNode::element(), GXml::element(), GXmlNode::elements(), G_READ, GObservation::id(), GObservation::name(), and GObservation::read().

Referenced by load().

void GObservations::remove ( const int &  index)
virtual

Remove observation from container.

Parameters
[in]indexObservation index [0,...,size()[.
Exceptions
GException::out_of_rangeObservation index is out of range.

Removes observation of specified index from the container.

Implements GContainer.

Definition at line 406 of file GObservations.cpp.

References G_REMOVE, m_obs, and size().

void GObservations::remove_response_cache ( void  )

Remove response cache for all models.

Remove response cache for all models in the container.

Definition at line 975 of file GObservations.cpp.

References m_models, and GModels::size().

void GObservations::remove_response_cache ( const std::string &  name)

Remove response cache for model.

Parameters
[in]nameModel name.

Remove response cache for model name in the container.

Definition at line 994 of file GObservations.cpp.

References size().

void GObservations::reserve ( const int &  num)
inlinevirtual

Reserves space for observations in container.

Parameters
[in]numNumber of observations.

Reserves space for num observations in the container.

Implements GContainer.

Definition at line 281 of file GObservations.hpp.

References m_obs.

Referenced by extend().

void GObservations::save ( const GFilename filename) const

Save observations into XML file.

Parameters
[in]filenameXML filename.

Saves observations into a XML file. Please refer to the read() method for more information about the structure of the XML file.

Definition at line 536 of file GObservations.cpp.

References GXml::save(), and write().

GObservation * GObservations::set ( const int &  index,
const GObservation obs 
)

Set observation in container.

Parameters
[in]indexObservation index [0,...,size()[.
[in]obsObservation.
Returns
Pointer to deep copy of observation.
Exceptions
GException::out_of_rangeObservation index is out of range.
GException::invalid_valueObservation with same instrument and identifier already exists in container.

Set a deep copy and observation obs at the specified index in the container.

Definition at line 261 of file GObservations.cpp.

References GObservation::clone(), G_SET, get_index(), GObservation::id(), GObservation::instrument(), m_obs, size(), and gammalib::str().

int GObservations::size ( void  ) const
inlinevirtual

Return number of observations in container.

Returns
Number of observations in container.

Returns the number of observations in the observation container.

Implements GContainer.

Definition at line 253 of file GObservations.hpp.

References m_obs.

Referenced by at(), extend(), GCTACubeBackground::fill(), GCTACubeExposure::fill(), GCTAModelSpatialLookup::fill(), GCTACubePsf::fill(), GCTACubeEdisp::fill(), GCTAOnOffObservation::GCTAOnOffObservation(), get_index(), insert(), nobserved(), npred(), print(), remove(), remove_response_cache(), GObservations::likelihood::save_csv(), GObservations::likelihood::save_fits(), set(), and write().

void GObservations::write ( GXml xml) const

Write observations into XML document.

Parameters
[in]xmlXML document.

Write observations into the first observation library that is found in the XML document. In case that no observation library exists, one is added to the document. Please refer to the read() method for more information about the structure of the XML document.

Definition at line 654 of file GObservations.cpp.

References GXmlNode::append(), GXml::append(), GXmlElement::attribute(), GXml::element(), GXml::elements(), m_obs, and size().

Referenced by save().

Member Data Documentation

GObservations::likelihood GObservations::m_fct
protected

Optimizer function.

Definition at line 200 of file GObservations.hpp.

Referenced by copy_members(), errors(), errors_hessian(), eval(), function(), init_members(), logL(), npred(), and optimize().

GModels GObservations::m_models
protected
std::vector<GObservation*> GObservations::m_obs
protected

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