49#define G_LIKELIHOOD "GObservation::likelihood(GModels&, GVector*,"\
50 " GMatrixSparse*, double*)"
51#define G_MODEL1 "GObservation::model(GModels&, GEvent&, GVector*)"
52#define G_MODEL2 "GObservation::model(GModels&, GMatrixSparse*)"
53#define G_EVENTS "GObservation::events()"
54#define G_NPRED "GObservation::npred(GModel&)"
55#define G_NPRED_SPEC "GObservation::npred_spec(GModel&, GTime&)"
64#define G_LN_ENERGY_INT
176 std::string msg =
"No events allocated for observation.";
197 std::string msg =
"No events allocated for observation.";
270 std::string msg =
"Invalid statistic \""+
statistic+
"\". Unbinned "
271 "optimization requires \"POISSON\" or \"CSTAT\" "
299 std::string msg =
"Invalid statistic \""+
statistic+
"\". Binned "
300 "optimization requires \"POISSON\", \"CSTAT\", "
301 "\"GAUSSIAN\" or \"CHI2\" statistic.";
343 bool use_grad = (gradients != NULL);
351 grad_size = gradients->
size();
360 for (
int i = 0; i < models.
size(); ++i) {
363 const GModel* mptr = models[i];
371 model += mptr->
eval(event, *
this, use_grad);
377 #if defined(G_RANGE_CHECK)
378 if (igrad+mptr->
size() > grad_size) {
379 std::string msg =
"Vector has not enough elements to "
380 "store the model parameter gradients. "+
382 "requested while vector only contains "+
397 for (
int index = 0; index < npars; ++index) {
408 (*gradients)[igrad+ipar] =
412 (*gradients)[igrad+ipar] =
417 (*gradients)[igrad+ipar] = 0.0;
428 for (
int ipar = 0; ipar < mptr->
size(); ++ipar) {
436 (*gradients)[igrad+ipar] =
440 (*gradients)[igrad+ipar] =
445 (*gradients)[igrad+ipar] = 0.0;
457 igrad += mptr->
size();
490 bool use_grad = (gradients != NULL);
505 int ncolumns = (nevents > 0) ? models.
npars() : 0;
506 if (gradients->
columns() != ncolumns) {
508 "columns in gradient matrix differs from expected "
510 "specify compatible arguments.";
515 int nrows = (ncolumns > 0) ? nevents : 0;
516 if (gradients->
rows() != nrows) {
518 "rows in gradient matrix differs from expected "
520 "specify compatible arguments.";
527 for (
int i = 0; i < models.
size(); ++i) {
530 const GModel* mptr = models[i];
540 values += mptr->
eval(*
this, NULL);
545 else if (nevents > 0) {
554 values += mptr->
eval(*
this, &grads);
568 for (
int index = 0; index < npars; ++index) {
582 grad = grads.
column(ipar);
589 gradients->
column(igrad+ipar, grad);
601 for (
int ipar = 0; ipar < mptr->
size(); ++ipar) {
612 grad = grads.
column(ipar);
613 #if defined(G_DEBUG_VECTOR_MODEL)
616 for (
int k = 0; k < nevents; ++k) {
617 if (grad[k] != 0.0 || num_grad[k] != 0.0) {
618 std::cout <<
"ipar=" << ipar;
619 std::cout <<
" par=" << par.
name();
620 std::cout <<
" k=" << k;
621 std::cout <<
" ana=" << grad[k];
622 std::cout <<
" num=" << num_grad[k];
623 std::cout << std::endl;
634 gradients->
column(igrad+ipar, grad);
647 igrad += mptr->
size();
702 #if defined(G_RANGE_CHECK)
703 if (gradients != NULL) {
704 if (models.
npars() != gradients->
size()) {
705 std::string msg =
"Model parameter gradients vector size "+
708 " model parameters. Please specify a gradients "
709 "vector with the appropriate size.";
720 if (gradients != NULL) {
725 for (
int i = 0; i < models.
size(); ++i) {
728 const GModel* mptr = models[i];
739 if (gradients != NULL) {
740 for (
int k = 0; k < mptr->
size(); ++k) {
742 (*gradients)[igrad+k] =
npred_grad(*mptr, par);
749 igrad += mptr->
size();
809 if (
model.is_constant()) {
828 for (
int i = 0; i <
events()->
gti().size(); ++i) {
831 double tstart =
events()->
gti().tstart(i).secs();
832 double tstop =
events()->
gti().tstop(i).secs();
835 if (tstop <= tstart) {
836 std::string msg =
"Start time "+
838 " is past the stop time "+
840 " for Good Time Interval "+
842 "that all Good Time Intervals have start "
843 "times that are earlier than the stop "
881 const GEvent& event)
const
1011 gradients = (fs1 - fs2) / h;
1028 gradients = (fs1 - fs2) / h;
1043 gradients = 0.5 * (fs1 - fs2) / h;
1180 std::string model_name =
id() +
":" +
name;
1208 std::string
id =
model.name() +
":" + par.
name();
1238 bool in_stack =
false;
1241 std::string
id =
model.name() +
":" + par.
name();
1370 double* npred)
const
1377 int npars = gradients->
size();
1380 int* inx =
new int[npars];
1381 double* values =
new double[npars];
1386 double npred_value = this->
npred(models, &wrk_grad);
1389 value += npred_value;
1390 *
npred += npred_value;
1391 *gradients += wrk_grad;
1394 GVector model_vector = this->
model(models, &wrk_matrix);
1400 for (
int i = 0; i < nevents; ++i) {
1408 double model = model_vector[i];
1421 for (
int ipar = 0; ipar < npars; ++ipar) {
1431 value -= std::log(
model);
1439 double fb = 1.0 /
model;
1440 double fa = fb /
model;
1441 for (
int jdev = 0; jdev < ndev; ++jdev) {
1444 register int jpar = inx[jdev];
1445 double g = wrk_grad[jpar];
1446 double fa_i = fa * g;
1449 (*gradients)[jpar] -= fb * g;
1452 register int* ipar = inx;
1454 for (
register int idev = 0; idev < ndev; ++idev, ++ipar) {
1455 values[idev] = fa_i * wrk_grad[*ipar];
1466 if (values != NULL)
delete [] values;
1467 if (inx != NULL)
delete [] inx;
1503 double* npred)
const
1509 #if defined(G_OPT_DEBUG)
1512 int n_small_model = 0;
1513 int n_zero_data = 0;
1514 int n_exclude_data = 0;
1515 double sum_data = 0.0;
1516 double sum_model = 0.0;
1517 double init_value = value;
1522 int npars = gradients->
size();
1525 int* inx =
new int[npars];
1526 double* values =
new double[npars];
1530 GVector model_vector = this->
model(models, &wrk_matrix);
1536 for (
int i = 0; i < nevents; ++i) {
1544 #if defined(G_OPT_DEBUG)
1553 double data = bin->
counts();
1557 #if defined(G_OPT_DEBUG)
1564 double model = model_vector[i] * bin->
size();
1568 #if defined(G_OPT_DEBUG)
1575 #if defined(G_OPT_DEBUG)
1590 for (
int ipar = 0; ipar < npars; ++ipar) {
1615 double fb = data /
model;
1616 double fc = (1.0 - fb);
1617 double fa = fb /
model;
1620 for (
int jdev = 0; jdev < ndev; ++jdev) {
1623 register int jpar = inx[jdev];
1624 double g = wrk_grad[jpar];
1625 double fa_i = fa * g;
1628 (*gradients)[jpar] += fc * g;
1631 register int* ipar = inx;
1632 for (
register int idev = 0; idev < ndev; ++idev, ++ipar) {
1633 values[idev] = fa_i * wrk_grad[*ipar];
1647 #if defined(G_OPT_DEBUG)
1661 register int* ipar = inx;
1662 for (
register int idev = 0; idev < ndev; ++idev, ++ipar) {
1663 (*gradients)[*ipar] += wrk_grad[*ipar];
1671 if (values != NULL)
delete [] values;
1672 if (inx != NULL)
delete [] inx;
1675 #if defined(G_OPT_DEBUG)
1676 std::cout <<
"Number of bins: " << n_bins << std::endl;
1677 std::cout <<
"Number of bins used for computation: " << n_used << std::endl;
1678 std::cout <<
"Number of bins excluded from data: " << n_exclude_data << std::endl;
1679 std::cout <<
"Number of bins excluded due to small model: " << n_small_model << std::endl;
1680 std::cout <<
"Number of bins with zero data: " << n_zero_data << std::endl;
1681 std::cout <<
"Sum of data: " << sum_data << std::endl;
1682 std::cout <<
"Sum of model: " << sum_model << std::endl;
1683 std::cout <<
"Initial statistic: " << init_value << std::endl;
1684 std::cout <<
"Statistic: " << value-init_value << std::endl;
1722 double* npred)
const
1729 int npars = gradients->
size();
1732 int* inx =
new int[npars];
1733 double* values =
new double[npars];
1737 GVector model_vector = this->
model(models, &wrk_matrix);
1743 for (
int i = 0; i < nevents; ++i) {
1755 double data = bin->
counts();
1763 double sigma = bin->
error();
1771 double model = model_vector[i] * bin->
size();
1787 for (
int ipar = 0; ipar < npars; ++ipar) {
1796 double weight = 1.0 / (sigma * sigma);
1799 double fa = data -
model;
1800 value += 0.5 * (fa * fa * weight);
1808 for (
int jdev = 0; jdev < ndev; ++jdev) {
1811 register int jpar = inx[jdev];
1812 double fa_i = wrk_grad[jpar] * weight;
1815 (*gradients)[jpar] -= fa * fa_i;
1818 register int* ipar = inx;
1819 for (
register int idev = 0; idev < ndev; ++idev, ++ipar) {
1820 values[idev] = fa_i * wrk_grad[*ipar];
1831 if (values != NULL)
delete [] values;
1832 if (inx != NULL)
delete [] inx;
1894 m_par->factor_value(x);
1897 double npred = m_parent->npred(*m_model);
1916 return (m_parent->npred_spec(*m_model, time));
1957 const GTime& obsTime)
const
1960 static const int iter = 8;
1963 double result = 0.0;
1969 for (
int i = 0; i < ebounds.
size(); ++i) {
1972 double emin = ebounds.
emin(i).
MeV();
1973 double emax = ebounds.
emax(i).
MeV();
1986 #if defined(G_LN_ENERGY_INT)
1987 emin = std::log(emin);
1988 emax = std::log(emax);
1990 result += integral.
romberg(emin, emax);
1997 #if defined(G_NAN_CHECK)
1999 std::cout <<
"*** ERROR: GObservation::npred_spec:";
2000 std::cout <<
" NaN/Inf encountered";
2001 std::cout <<
" (result=" << result;
2002 std::cout <<
")" << std::endl;
2024 #if defined(G_LN_ENERGY_INT)
2025 double expx = std::exp(x);
2032 double value = m_model->npred(eng, *m_time, *m_parent);
2035 #if defined(G_NAN_CHECK)
2036 double value_out = value;
2040 #if defined(G_LN_ENERGY_INT)
2045 #if defined(G_NAN_CHECK)
2047 std::cout <<
"*** ERROR: GObservation::npred_spec_kern::eval";
2048 std::cout <<
"(x=" << x <<
"): ";
2049 std::cout <<
" NaN/Inf encountered";
2050 std::cout <<
" (value=" << value;
2051 std::cout <<
" (value_out=" << value_out;
2052 #if defined(G_LN_ENERGY_INT)
2053 std::cout <<
" exp(x)=" << expx;
2055 std::cout <<
")" << std::endl;
GDerivative class interface definition.
Abstract event bin base class definition.
Abstract event bin container class interface definition.
Abstract event atom container class interface definition.
Exception handler interface definition.
Integration class interface definition.
Sparse matrix class definition.
Abstract data model base class interface definition.
Model parameter class interface definition.
Sky model class interface definition.
Abstract model base class interface definition.
Model container class definition.
const double minerr
Minimum statistical error.
const double minmod
Minimum model value.
Abstract observation base class interface definition.
Abstract response base class definition.
Vector class interface definition.
Numerical derivatives class.
double left_difference(const double &x, const double &h)
Returns gradient computed from left-sided function difference.
double right_difference(const double &x, const double &h)
Returns gradient computed from right-sided function difference.
double difference(const double &x, const double &h)
Returns gradient computed from symmetric function difference.
Energy boundaries container class.
const GEnergy & emax(void) const
Return maximum energy of all intervals.
int size(void) const
Return number of energy boundaries.
const GEnergy & emin(void) const
Return minimum energy of all intervals.
Class that handles energies in a unit independent way.
double MeV(void) const
Return energy in MeV.
Abstract interface for the event bin class.
virtual double error(void) const =0
Return error in number of counts.
virtual double counts(void) const =0
virtual double size(void) const =0
Abstract event bin container class.
Abstract event atom container class.
Abstract interface for the event classes.
Abstract event container class.
void gti(const GGti >i)
Set Good Time Intervals.
void ebounds(const GEbounds &ebounds)
Set energy boundaries.
virtual int size(void) const =0
virtual int number(void) const =0
virtual GEvents * clone(void) const =0
Clones object.
GIntegral class interface definition.
double romberg(std::vector< double > bounds, const int &order=5)
Perform Romberg integration.
void fixed_iter(const int &iter)
Set fixed number of iterations.
const int & rows(void) const
Return number of matrix rows.
const int & columns(void) const
Return number of matrix columns.
Sparse matrix class interface definition.
virtual GVector column(const int &column) const
Extract column as vector from matrix.
GMatrixSparse transpose(void) const
Return transposed matrix.
void stack_destroy(void)
Destroy matrix stack.
virtual void add_to_column(const int &column, const GVector &vector)
Add vector column into matrix.
void stack_init(const int &size=0, const int &entries=0)
Initialises matrix filling stack.
const bool & has_eval_indices(void) const
Signals that parameter indices updated by eval() method were set.
bool is_valid(const std::string &instruments, const std::string &ids) const
Verifies if model is valid for a given instrument and identifier.
int size(void) const
Return number of parameters in model.
virtual double eval(const GEvent &event, const GObservation &obs, const bool &gradients=false) const =0
const std::vector< int > & eval_indices(void) const
Return vector of parameter indices updated by eval() method.
int size(void) const
Return number of models in container.
int npars(void) const
Return number of model parameters in container.
GModelPar * m_par
Model parameter.
const GModel * m_model
Model.
const GObservation * m_parent
Observation.
double eval(const double &x)
Model function evaluation for gradient computation.
const GEvent * m_event
Event.
double eval(const double &x)
Npred function evaluation for gradient computation.
double eval(const double &x)
Integration kernel for npred() method.
double eval(const double &x)
Integration kernel for npred_spec() method.
Abstract observation base class.
std::string m_statistic
Optimizer statistic.
virtual double ontime(void) const =0
const std::string & statistic(void) const
Return optimizer statistic.
virtual double npred(const GModels &models, GVector *gradients=NULL) const
Return total number (and optionally gradients) of predicted counts for all models.
virtual double model(const GModels &models, const GEvent &event, GVector *gradients=NULL) const
Return model value and (optionally) gradients.
void copy_members(const GObservation &obs)
Copy class members.
virtual int nobserved(void) const
Return total number of observed events.
const std::string & id(void) const
Return observation identifier.
virtual double likelihood(const GModels &models, GVector *gradients, GMatrixSparse *curvature, double *npred) const
Compute likelihood function.
const std::string & name(void) const
Return observation name.
virtual bool use_event_for_likelihood(const int &index) const
Check whether bin should be used for likelihood analysis.
virtual ~GObservation(void)
Destructor.
GObservation(void)
Void constructor.
double m_grad_step_size
Gradient step size.
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 ar...
virtual GEvents * events(void)
Return events.
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 arr...
void computed_gradient(const GModel &model, const GModelPar &par) const
Signals that an analytical gradient was computed for a model parameter.
virtual const GResponse * response(void) const =0
bool has_gradient(const GModel &model, const GModelPar &par) const
Check whether a model parameter has an analytical gradient.
void init_members(void)
Initialise class members.
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 arra...
virtual double npred_spec(const GModel &model, const GTime &obsTime) const
Integrates spatially integrated Npred kernel spectrally.
virtual GObservation & operator=(const GObservation &obs)
Assignment operator.
std::vector< std::string > m_pars_with_gradients
GEvents * m_events
Pointer to event container.
virtual std::string instrument(void) const =0
void free_members(void)
Delete class members.
virtual void remove_response_cache(const std::string &name)
Response cache removal hook.
std::string m_name
Observation name.
virtual double npred_grad(const GModel &model, const GModelPar &par) const
Returns parameter gradient of Npred.
virtual double model_grad(const GModel &model, const GModelPar &par, const GEvent &event) const
Returns parameter gradient of model for a given event.
std::string m_id
Observation identifier.
const double & factor_max(void) const
Return parameter maximum factor boundary.
const double & factor_value(void) const
Return parameter factor value.
bool is_free(void) const
Signal if parameter is free.
bool has_min(void) const
Signal if parameter has minimum boundary.
const double & factor_gradient(void) const
Return parameter factor gradient.
bool has_max(void) const
Signal if parameter has maximum boundary.
const double & factor_min(void) const
Return parameter minimum factor boundary.
const std::string & name(void) const
Return parameter name.
Abstract instrument response base class.
const double & secs(void) const
Return time in seconds in native reference (TT)
const int & size(void) const
Return size of vector.
bool is_infinite(const double &x)
Signal if argument is infinite.
bool is_notanumber(const double &x)
Signal if argument is not a number.
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
std::string toupper(const std::string &s)
Convert string to upper case.