27#ifndef GOBSERVATIONS_HPP
28#define GOBSERVATIONS_HPP
112 int size(
void)
const;
117 void remove(
const int& index);
120 bool contains(
const std::string& instrument,
121 const std::string&
id)
const;
133 double logL(
void)
const;
135 double npred(
void)
const;
136 double npred(
const std::string& name,
137 const std::string& instrument =
"")
const;
139 const std::string& instrument =
"")
const;
158 virtual double value(
void)
const;
164 double npred(
void)
const;
194 int get_index(
const std::string& instrument,
195 const std::string&
id)
const;
212 return (
"GObservations");
255 return (
int)
m_obs.size();
269 return (
m_obs.empty());
Definition of interface for container classes.
Model container class definition.
Abstract observation base class interface definition.
Optimizer function abstract base class.
Interface class for container classes.
Sparse matrix class interface definition.
Abstract observation base class.
void set(GObservations *obs)
Set observation container.
void save(const GFilename &filename) const
Save likelihood fit results into a CSV or FITS file.
GObservations * m_this
Pointer to GObservations object.
GVector * m_gradient
Pointer to gradient vector.
~likelihood(void)
Destructor.
double m_value
Function value.
virtual GVector * gradient(void)
Return pointer to gradient vector.
virtual double value(void) const
Return log-likelihood value of optimizer function.
virtual GMatrixSparse * curvature(void)
Return pointer to curvature matrix.
void save_fits(const GFilename &filename) const
Save likelihood fit results into FITS file.
GMatrixSparse hessian(const GOptimizerPars &pars)
Compute Hessian matrix.
void init_members(void)
Initialise class members.
void free_members(void)
Delete class members.
GMatrixSparse covariance(void) const
Compute covariance matrix.
void save_csv(const GFilename &filename) const
Save likelihood fit results into CSV file.
likelihood & operator=(const likelihood &fct)
Assignment operator.
double npred(void) const
Return total number of predicted events.
double m_npred
Total number of predicted events.
std::vector< std::string > covariance_names(void) const
Return covariance matrix row and column names.
likelihood(void)
Void constructor.
GMatrixSparse * m_curvature
Pointer to curvature matrix.
void copy_members(const likelihood &fct)
Copy class members.
Observation container class.
GObservation * operator[](const int &index)
Return pointer to observation.
int nobserved(void) const
Return total number of observed events.
bool contains(const std::string &instrument, const std::string &id) const
Signals if observation exists.
double logL(void) const
Return log-likelihood of models.
void init_members(void)
Initialise class members.
GObservations & operator=(const GObservations &obs)
Assignment operator.
std::vector< GObservation * > m_obs
List of observations.
void reserve(const int &num)
Reserves space for observations in container.
void copy_members(const GObservations &obs)
Copy class members.
void extend(const GObservations &obs)
Append observations from observation container.
void clear(void)
Clear observations.
GModels m_models
List of models.
int size(void) const
Return number of observations in container.
GObservations * clone(void) const
Clone observations.
std::string classname(void) const
Return class name.
GObservations::likelihood m_fct
Optimizer function.
void read(const GXml &xml)
Read observations from XML document.
bool is_empty(void) const
Signals if there are no observations in container.
void remove(const int &index)
Remove observation from container.
GObservation * set(const int &index, const GObservation &obs)
Set observation in container.
std::string print(const GChatter &chatter=NORMAL) const
Print observation list information.
const GObservations::likelihood & function(void) const
Return likelihood function.
virtual ~GObservations(void)
Destructor.
GObservations(void)
Void constructor.
void remove_response_cache(void)
Remove response cache for all models.
const GModels & models(void) const
Return model container.
void eval(void)
Evaluate function.
int get_index(const std::string &instrument, const std::string &id) const
Return observation index by instrument and identifier.
double npred(void) const
Return total number of predicted events in models.
void save(const GFilename &filename) const
Save observations into XML file.
void optimize(GOptimizer &opt)
Optimize model parameters using optimizer.
void errors(GOptimizer &opt)
Computes parameter errors using optimizer.
GObservation * insert(const int &index, const GObservation &obs)
Insert observation into container.
void write(GXml &xml) const
Write observations into XML document.
void errors_hessian(void)
Computes parameter errors using hessian matrix and optimizer.
void load(const GFilename &filename)
Load observations from XML file.
GObservation * append(const GObservation &obs)
Append observation to container.
void free_members(void)
Delete class members.
GObservation * at(const int &index)
Return pointer to observation.
Optimizer function abstract base class.
Optimizer parameter container class.
Abstract optimizer abstract base class.