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

Abstract observation base class. More...

#include <GObservation.hpp>

Inheritance diagram for GObservation:
GBase GCOMObservation GCTAObservation GCTAOnOffObservation GLATObservation GMWLObservation GSPIObservation

Classes

class  model_func
 
class  npred_func
 
class  npred_kern
 
class  npred_spec_kern
 

Public Member Functions

 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 void clear (void)=0
 Clear object. More...
 
virtual GObservationclone (void) const =0
 Clones object. More...
 
virtual std::string classname (void) const =0
 Return class name. More...
 
virtual void response (const GResponse &rsp)=0
 
virtual const GResponseresponse (void) const =0
 
virtual std::string instrument (void) const =0
 
virtual double ontime (void) const =0
 
virtual double livetime (void) const =0
 
virtual double deadc (const GTime &time=GTime()) const =0
 
virtual void read (const GXmlElement &xml)=0
 
virtual void write (GXmlElement &xml) const =0
 
virtual std::string print (const GChatter &chatter=NORMAL) const =0
 Print content of object. 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 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_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

Abstract observation base class.

This class provides an abstract interface for an observation. The observation collects information about the instrument, holds the measured events, and provides information about the analysis definition.

The response() method returns a pointer to the response function. The derived classes have to make sure that this method never returns NULL.

The method model() returns the probability for an event to be measured with a given instrument direction, a given energy and at a given time, given a source model and an instrument pointing direction. The method npred() returns the total number of expected events within the analysis region for a given source model and a given instrument pointing direction. The methods are defined as virtual and can be overloaded by derived classes that implement instrument specific observations in order to optimize the execution speed for data analysis.

Definition at line 70 of file GObservation.hpp.

Constructor & Destructor Documentation

GObservation::GObservation ( void  )

Void constructor.

Definition at line 80 of file GObservation.cpp.

References init_members().

GObservation::GObservation ( const GObservation obs)

Copy constructor.

Parameters
[in]obsObservation.

Instantiate the class by copying from an existing observation.

Definition at line 97 of file GObservation.cpp.

References copy_members(), and init_members().

GObservation::~GObservation ( void  )
virtual

Destructor.

Definition at line 113 of file GObservation.cpp.

References free_members().

Member Function Documentation

virtual std::string GObservation::classname ( void  ) const
pure virtual

Return class name.

Returns
String containing the class name.

Returns the class name for non-abstract classes in a human readable way.

Implements GBase.

Implemented in GCTAOnOffObservation, GCOMObservation, GCTAObservation, GSPIObservation, GMWLObservation, and GLATObservation.

virtual void GObservation::clear ( void  )
pure virtual

Clear object.

Sets the object to a clean initial state. After calling the method the object will be in the same state as it were if an empty instance of the object would have been created.

Implements GBase.

Implemented in GCTAOnOffObservation, GCOMObservation, GCTAObservation, GSPIObservation, GMWLObservation, and GLATObservation.

virtual GObservation* GObservation::clone ( void  ) const
pure virtual

Clones object.

Returns
Pointer to deep copy of object.

Creates a deep copy of the object and returns a pointer to the object.

Implements GBase.

Implemented in GCTAOnOffObservation, GCOMObservation, GCTAObservation, GSPIObservation, GMWLObservation, and GLATObservation.

Referenced by GObservations::append(), GObservations::insert(), and GObservations::set().

void GObservation::computed_gradient ( const GModel model,
const GModelPar par 
) const

Signals that an analytical gradient was computed for a model parameter.

Parameters
[in]modelModel.
[in]parModel parameter.

Signals that an analytical gradient was computed for a model parameter.

Definition at line 1235 of file GObservation.cpp.

References m_pars_with_gradients, GModel::name(), and GOptimizerPar::name().

Referenced by GModelData::eval(), GCTAResponseCube::irf_radial(), and GModelSky::signal_analytical_gradients().

void GObservation::copy_members ( const GObservation obs)
protected

Copy class members.

Parameters
[in]obsObservation.

Copy members from an observation.

Definition at line 1300 of file GObservation.cpp.

References GEvents::clone(), m_events, m_grad_step_size, m_id, m_name, m_pars_with_gradients, and m_statistic.

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

GEvents * GObservation::events ( void  )
virtual

Return events.

Exceptions
GException::no_eventsNo events allocated for observation.

Returns pointer to events.

Definition at line 172 of file GObservation.cpp.

References G_EVENTS, and m_events.

Referenced by GCOMObservation::compute_drb_bgdlixa(), GCOMObservation::compute_drb_bgdlixe(), GCOMObservation::compute_drb_bgdlixf(), GCOMObservation::compute_drb_phinor(), GCOMDri::compute_dre(), GCOMDris::compute_drws(), GResponse::convolve(), gammalib::cta_event_cube(), gammalib::cta_event_list(), GModelData::eval(), GModels::eval(), GResponse::eval_probs(), GCTAModelSpatialLookup::fill_buffer(), GSPIModelDataSpace::GSPIModelDataSpace(), GCOMResponse::irf(), GSPIResponse::irf(), GResponse::irf_composite(), GCOMResponse::irf_diffuse(), GResponse::irf_diffuse(), GCOMResponse::irf_elliptical(), GResponse::irf_elliptical(), GCOMResponse::irf_ptsrc(), GResponse::irf_ptsrc(), GCOMResponse::irf_radial(), GCTAResponseCube::irf_radial(), GResponse::irf_radial(), GCTAResponseCube::irf_spatial(), GResponse::irf_spatial(), likelihood(), likelihood_gaussian_binned(), likelihood_poisson_binned(), likelihood_poisson_unbinned(), GCTAObservation::load(), GLATObservation::load_binned(), GLATObservation::load_unbinned(), GCTAModelRadialAcceptance::mc(), model(), model_grad(), npred(), GCOMObservation::npred(), npred_spec(), GSPIObservation::read(), GCTAObservation::read(), GLATMeanPsf::set(), GSPIResponse::set(), GCTAOnOffObservation::set(), GLATMeanPsf::set_map_corrections(), GSPIModelDataSpace::setup_model(), GCTAObservation::write(), and GCTAObservation::write_attributes().

const GEvents * GObservation::events ( void  ) const
virtual

Return events (const version)

Exceptions
GException::no_eventsNo events allocated for observation.

Returns const pointer to events.

Definition at line 193 of file GObservation.cpp.

References G_EVENTS, and m_events.

void GObservation::events ( const GEvents events)
virtual

Set event container.

Parameters
[in]eventsEvent container.

Set the event container for this observation by cloning the container specified in the argument.

Definition at line 214 of file GObservation.cpp.

References GEvents::clone(), and m_events.

void GObservation::free_members ( void  )
protected
const double & GObservation::grad_step_size ( void  ) const
inlinevirtual

Return gradient step size.

Returns
Gradient step size.

Definition at line 333 of file GObservation.hpp.

References m_grad_step_size.

bool GObservation::has_events ( void  ) const
inline

Signal if observation has events.

Returns
True if observation contains events.

Definition at line 240 of file GObservation.hpp.

References m_events.

bool GObservation::has_gradient ( const GModel model,
const GModelPar par 
) const

Check whether a model parameter has an analytical gradient.

Parameters
[in]modelModel.
[in]parModel parameter.
Returns
True if model parameter is free and has an analytical gradient.

Checks whether a model parameter is free and has an analytical gradient.

Definition at line 1199 of file GObservation.cpp.

References GOptimizerPar::is_free(), m_pars_with_gradients, GModel::name(), and GOptimizerPar::name().

Referenced by model().

void GObservation::init_members ( void  )
protected

Initialise class members.

The step size for gradient computation is fixed by default to 0.0002 degrees, which corresponds to about 1 arcsec for parameters that are given in degrees. The reasoning behind this value is that parameters that use numerical gradients are typically angles, such as for example the position, and we want to achieve arcsec precision with this method.

Definition at line 1276 of file GObservation.cpp.

References m_events, m_grad_step_size, m_id, m_name, m_pars_with_gradients, and m_statistic.

Referenced by GLATObservation::clear(), GMWLObservation::clear(), GSPIObservation::clear(), GCTAObservation::clear(), GCOMObservation::clear(), GObservation(), and operator=().

double GObservation::likelihood ( const GModels models,
GVector gradients,
GMatrixSparse curvature,
double *  npred 
) const
virtual

Compute likelihood function.

Parameters
[in]modelsModels.
[in,out]gradientsPointer to gradients.
[in,out]curvaturePointer to curvature matrix.
[in,out]npredPointer to Npred value.
Returns
Likelihood.

Computes the likelihood for a specified set of models. The method also returns the gradients, the curvature matrix, and the number of events that are predicted by all models.

Reimplemented in GCTAOnOffObservation.

Definition at line 243 of file GObservation.cpp.

References events(), G_LIKELIHOOD, likelihood_gaussian_binned(), likelihood_poisson_binned(), likelihood_poisson_unbinned(), statistic(), and gammalib::toupper().

double GObservation::likelihood_gaussian_binned ( const GModels models,
GVector gradients,
GMatrixSparse curvature,
double *  npred 
) const
protectedvirtual

Evaluate log-likelihood function for Gaussian statistic and binned analysis (version with working arrays)

Parameters
[in]modelsModels.
[in,out]gradientsGradient.
[in,out]curvatureCurvature matrix.
[in,out]npredNumber of predicted events.
Returns
Likelihood value.

This method evaluates the -(log-likelihood) function for parameter optimisation using binned analysis and Poisson statistic. The -(log-likelihood) function is given by

\[ L = 1/2 \sum_i (n_i - e_i)^2 \sigma_i^{-2} \]

where the sum is taken over all data space bins, \(n_i\) is the observed number of counts, \(e_i\) is the model and \(\sigma_i\) is the statistical uncertainty. This method also computes the parameter gradients \(\delta L/dp\) and the curvature matrix \(\delta^2 L/dp_1 dp_2\) and also updates the total number of predicted events m_npred.

Definition at line 1719 of file GObservation.cpp.

References GMatrixSparse::add_to_column(), GMatrixSparse::column(), GEventBin::counts(), GEventBin::error(), events(), gammalib::is_infinite(), minerr, minmod, model(), GEventBin::size(), GEvents::size(), GVector::size(), GMatrixSparse::transpose(), and use_event_for_likelihood().

Referenced by likelihood().

double GObservation::likelihood_poisson_binned ( const GModels models,
GVector gradients,
GMatrixSparse curvature,
double *  npred 
) const
protectedvirtual

Evaluate log-likelihood function for Poisson statistic and binned analysis (version with working arrays)

Parameters
[in]modelsModels.
[in,out]gradientsGradient.
[in,out]curvatureCurvature matrix.
[in,out]npredNumber of predicted events.
Returns
Likelihood value.

This method evaluates the -(log-likelihood) function for parameter optimisation using binned analysis and Poisson statistic. The -(log-likelihood) function is given by

\[ L=-\sum_i n_i \log e_i - e_i \]

where the sum is taken over all data space bins, \(n_i\) is the observed number of counts and \(e_i\) is the model. This method also computes the parameter gradients \(\delta L/dp\) and the curvature matrix \(\delta^2 L/dp_1 dp_2\) and also updates the total number of predicted events m_npred.

Definition at line 1500 of file GObservation.cpp.

References GMatrixSparse::add_to_column(), GMatrixSparse::column(), GEventBin::counts(), events(), gammalib::is_infinite(), log(), minmod, model(), GEventBin::size(), GEvents::size(), GVector::size(), GMatrixSparse::transpose(), and use_event_for_likelihood().

Referenced by likelihood().

double GObservation::likelihood_poisson_unbinned ( const GModels models,
GVector gradients,
GMatrixSparse curvature,
double *  npred 
) const
protectedvirtual

Evaluate log-likelihood function for Poisson statistic and unbinned analysis (version with working arrays)

Parameters
[in]modelsModels.
[in,out]gradientsGradient.
[in,out]curvatureCurvature matrix.
[in,out]npredNumber of predicted events.
Returns
Likelihood value.

This method evaluates the -(log-likelihood) function for parameter optimisation using unbinned analysis and Poisson statistic. The -(log-likelihood) function is given by

\[ L = N_{\rm pred} - \sum_i \log e_i \]

where \(N_{\rm pred}\) is the number of events predicted by the model, and the sum is taken over all events. This method also computes the parameter gradients \(\delta L/dp\) and the curvature matrix \(\delta^2 L/dp_1 dp_2\).

Definition at line 1367 of file GObservation.cpp.

References GMatrixSparse::add_to_column(), GMatrixSparse::column(), events(), gammalib::is_infinite(), log(), minmod, model(), npred(), GEvents::size(), GVector::size(), GMatrixSparse::transpose(), and use_event_for_likelihood().

Referenced by likelihood().

virtual double GObservation::livetime ( void  ) const
pure virtual
double GObservation::model ( const GModels models,
const GEvent event,
GVector gradients = NULL 
) const
virtual

Return model value and (optionally) gradients.

Parameters
[in]modelsModel container.
[in]eventObserved event.
[out]gradientsPointer to gradient vector (optional).
Returns
Model value.
Exceptions
GException::invalid_valueDimension of gradient vector mismatches number of parameters.

Implements generic model and gradient evaluation for an observation. The model gives the probability for an event to occur with a given instrument direction, at a given energy and at a given time. The gradient is the parameter derivative of this probability. If NULL is passed for the gradient vector, then gradients will not be computed.

The method will only operate on models for which the list of instruments and observation identifiers matches those of the observation. Models that do not match will be skipped.

Definition at line 333 of file GObservation.cpp.

References GModel::eval(), GModel::eval_indices(), GOptimizerPar::factor_gradient(), G_MODEL1, GModel::has_eval_indices(), has_gradient(), instrument(), GOptimizerPar::is_free(), GModel::is_valid(), m_pars_with_gradients, model_grad(), GModels::npars(), GVector::size(), GModel::size(), GModels::size(), and gammalib::str().

Referenced by GCOMObservation::drm(), likelihood_gaussian_binned(), likelihood_poisson_binned(), likelihood_poisson_unbinned(), model_grad(), GCOMObservation::npred(), npred_grad(), and GCTAOnOffObservation::set().

GVector GObservation::model ( const GModels models,
GMatrixSparse gradients = NULL 
) const
virtual

Return vector of model values and (optionally) gradients.

Parameters
[in]modelsModel container.
[out]gradientsPointer to sparse gradient matrix.
Returns
Vector of model values.
Exceptions
GException::invalid_argumentGradient matrix mismatches number of events or model parameters.

Returns the model values for each event in the observation.

If gradients is not NULL, the matrix contains on output the model factor gradients for all events. Each row of the gradients matrix corresponds to one event, the columns correspond to the parameters of the models container.

Definition at line 485 of file GObservation.cpp.

References GMatrixSparse::column(), GMatrixBase::columns(), GModel::eval(), GModel::eval_indices(), events(), G_MODEL2, GModel::has_eval_indices(), has_gradient(), instrument(), GOptimizerPar::is_free(), GModel::is_valid(), m_pars_with_gradients, model_grad(), GOptimizerPar::name(), GModels::npars(), GMatrixBase::rows(), GEvents::size(), GModel::size(), GModels::size(), GMatrixSparse::stack_init(), and gammalib::str().

double GObservation::model_grad ( const GModel model,
const GModelPar par,
const GEvent event 
) const
virtual

Returns parameter gradient of model for a given event.

Parameters
[in]modelModel.
[in]parModel parameter.
[in]eventEvent.

This method uses a robust but simple difference method to estimate parameter gradients that have not been provided by the model. We use here a simple method as this method is likely used for spatial model parameters, and the spatial model may eventually be noisy due to numerical integration limits.

Definition at line 879 of file GObservation.cpp.

References GOptimizerPar::factor_max(), GOptimizerPar::factor_min(), GOptimizerPar::factor_value(), GOptimizerPar::has_max(), GOptimizerPar::has_min(), GOptimizerPar::is_free(), m_grad_step_size, and model().

Referenced by model().

GVector GObservation::model_grad ( const GModel model,
const GModelPar par 
) const
virtual

Returns parameter gradients of model for all events.

Parameters
[in]modelModel.
[in]parModel parameter.

This method uses a robust but simple difference method to estimate parameter gradients that have not been provided by the model. We use here a simple method as this method is likely used for spatial model parameters, and the spatial model may eventually be noisy due to numerical integration limits.

Definition at line 961 of file GObservation.cpp.

References GModel::eval(), events(), GOptimizerPar::factor_max(), GOptimizerPar::factor_min(), GOptimizerPar::factor_value(), GOptimizerPar::has_max(), GOptimizerPar::has_min(), GOptimizerPar::is_free(), and m_grad_step_size.

int GObservation::nobserved ( void  ) const
virtual

Return total number of observed events.

Returns
Total number of observed events.

Returns the total number of observed events. If the observation does not contain any events the method returns zero.

Reimplemented in GCTAOnOffObservation.

Definition at line 666 of file GObservation.cpp.

References m_events, and GEvents::number().

double GObservation::npred ( const GModels models,
GVector gradients = NULL 
) const
virtual

Return total number (and optionally gradients) of predicted counts for all models.

Parameters
[in]modelsModels.
[out]gradientsModel parameter gradients (optional).
Exceptions
GException::invalid_argumentDimension of gradient vector mismatches number of parameters.

Returns the total number of predicted counts within the analysis region. If NULL is passed for the gradient vector then gradients will not be computed.

The method will only operate on models for which the list of instruments and observation identifiers matches those of the observation. Models that do not match will be skipped.

Definition at line 699 of file GObservation.cpp.

References G_NPRED, instrument(), GModel::is_valid(), GModels::npars(), npred_grad(), GVector::size(), GModel::size(), GModels::size(), and gammalib::str().

Referenced by GObservation::npred_func::eval(), likelihood_poisson_unbinned(), and npred().

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

Return total number of predicted counts for one model.

Parameters
[in]modelGamma-ray source model.
Exceptions
GException::gti_invalidGood Time Interval is invalid.

Computes

\[ N_{\rm pred} = \int_{\rm GTI} \int_{E_{\rm bounds}} \int_{\rm ROI} P(p',E',t') \, dp' \, dE' \, dt' \]

of the event probability

\[ P(p',E',t') = \int \int \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \, dE \, dt \]

where \(S(p,E,t)\) is the source model, \(R(p',E',t'|p,E,t)\) is the instrument response function, \(p'\) is the measured photon direction, \(E'\) is the measured photon energy, \(t'\) is the measured photon arrival time, \(p\) is the true photon arrival direction, \(E\) is the true photon energy, and \(t\) is the true photon arrival time.

This method performs the time integration over the Good Time Intervals \({\rm GTI}\). Note that MET is used for the time integration interval. This, however, is no specialisation since npred_grad_kern_spec::eval() converts the argument back in a GTime object by assuming that the argument is in MET, hence the correct time system will be used at the end by the method.

Reimplemented in GCOMObservation.

Definition at line 799 of file GObservation.cpp.

References events(), G_NPRED, GEvents::gti(), instrument(), GModel::is_constant(), GModel::is_valid(), npred(), npred_spec(), ontime(), GIntegral::romberg(), and gammalib::str().

double GObservation::npred_grad ( const GModel model,
const GModelPar par 
) const
virtual

Returns parameter gradient of Npred.

Parameters
[in]modelGamma-ray source model.
[in]parModel parameter.

Computes

\[ \frac{{\rm d} N_{\rm pred}}{{\rm d} a_i} \]

where

\[ N_{\rm pred} = \int_{\rm GTI} \int_{E_{\rm bounds}} \int_{\rm ROI} P(p',E',t') \, dp' \, dE' \, dt' \]

is the integral of the event probability

\[ P(p',E',t') = \int \int \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \, dE \, dt \]

and \(a_i\) is the model parameter \(i\). Furthermore \(S(p,E,t)\) is the source model, \(R(p',E',t'|p,E,t)\) is the instrument response function, \(p'\) is the measured photon direction, \(E'\) is the measured photon energy, \(t'\) is the measured photon arrival time, \(p\) is the true photon arrival direction, \(E\) is the true photon energy, and \(t\) is the true photon arrival time.

This method uses a robust but simple difference method to estimate parameter gradients that have not been provided by the model. We use here a simple method as this method is likely used for spatial model parameters, and the spatial model may eventually be noisy due to numerical integration limits.

Definition at line 1102 of file GObservation.cpp.

References GOptimizerPar::factor_max(), GOptimizerPar::factor_min(), GOptimizerPar::factor_value(), GOptimizerPar::has_max(), GOptimizerPar::has_min(), GOptimizerPar::is_free(), m_grad_step_size, and model().

Referenced by npred().

double GObservation::npred_spec ( const GModel model,
const GTime obsTime 
) const
protectedvirtual

Integrates spatially integrated Npred kernel spectrally.

Parameters
[in]modelGamma-ray source model.
[in]obsTimeMeasured photon arrival time.
Exceptions
GException::erange_invalidEnergy range is invalid.

Computes

\[ N_{\rm pred}(t') = \int_{E_{\rm bounds}} \int_{\rm ROI} P(p',E',t') \, dp' \, dE' \]

of the event probability

\[ P(p',E',t') = \int \int \int S(p,E,t) \times R(p',E',t'|p,E,t) \, dp \, dE \, dt \]

where \(S(p,E,t)\) is the source model, \(R(p',E',t'|p,E,t)\) is the instrument response function, \(p'\) is the measured photon direction, \(E'\) is the measured photon energy, \(t'\) is the measured photon arrival time, \(p\) is the true photon arrival direction, \(E\) is the true photon energy, and \(t\) is the true photon arrival time.

This method performs the energy intergration over the energy boundaries \(E_{\rm bounds}\).

Definition at line 1956 of file GObservation.cpp.

References GEvents::ebounds(), GEbounds::emax(), GEbounds::emin(), events(), GIntegral::fixed_iter(), gammalib::is_infinite(), gammalib::is_notanumber(), log(), GEnergy::MeV(), GIntegral::romberg(), and GEbounds::size().

Referenced by npred().

virtual double GObservation::ontime ( void  ) const
pure virtual
GObservation & GObservation::operator= ( const GObservation obs)
virtual

Assignment operator.

Parameters
[in]obsObservation.
Returns
Observation.

Assign observation.

Definition at line 137 of file GObservation.cpp.

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

Referenced by GLATObservation::operator=(), GMWLObservation::operator=(), GSPIObservation::operator=(), GCTAObservation::operator=(), GCOMObservation::operator=(), and GCTAOnOffObservation::operator=().

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

Print content of object.

Parameters
[in]chatterChattiness (defaults to NORMAL).
Returns
String containing the content of the object.

Formats the content in a standard way and puts this content in a C++ string that is returned.

Implements GBase.

Implemented in GCTAOnOffObservation, GCOMObservation, GCTAObservation, GSPIObservation, GMWLObservation, and GLATObservation.

virtual void GObservation::read ( const GXmlElement xml)
pure virtual
void GObservation::remove_response_cache ( const std::string &  name)
virtual

Response cache removal hook.

Parameters
[in]nameModel name.

Remove response cache for model name from response cache.

Definition at line 1177 of file GObservation.cpp.

References id(), name(), and response().

virtual void GObservation::response ( const GResponse rsp)
pure virtual
virtual const GResponse* GObservation::response ( void  ) const
pure virtual
void GObservation::statistic ( const std::string &  statistic)
inline

Set optimizer statistic.

Parameters
[in]statisticOptimizer statistic.

Set optimizer statistic for the observation.

Definition at line 284 of file GObservation.hpp.

References m_statistic, and statistic().

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

Check whether bin should be used for likelihood analysis.

Parameters
[in]indexEvent index.
Returns
True.

This is a dummy virtual method that allows implementation of a hook for event selection in the likelihood computation. The dummy method always returns true.

Reimplemented in GCOMObservation.

Definition at line 1849 of file GObservation.cpp.

Referenced by likelihood_gaussian_binned(), likelihood_poisson_binned(), and likelihood_poisson_unbinned().

virtual void GObservation::write ( GXmlElement xml) const
pure virtual

Member Data Documentation

double GObservation::m_grad_step_size
protected

Gradient step size.

Definition at line 227 of file GObservation.hpp.

Referenced by copy_members(), grad_step_size(), init_members(), model_grad(), and npred_grad().

std::string GObservation::m_id
protected

Observation identifier.

Definition at line 224 of file GObservation.hpp.

Referenced by copy_members(), id(), init_members(), GCTAOnOffObservation::print(), and GCOMObservation::read().

std::string GObservation::m_name
protected
std::vector<std::string> GObservation::m_pars_with_gradients
mutableprotected

Definition at line 230 of file GObservation.hpp.

Referenced by computed_gradient(), copy_members(), has_gradient(), init_members(), and model().

std::string GObservation::m_statistic
protected

Optimizer statistic.

Definition at line 225 of file GObservation.hpp.

Referenced by copy_members(), GMWLObservation::init_members(), init_members(), and statistic().


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