GammaLib
2.1.0.dev
|
CTA On/Off observation class. More...
#include <GCTAOnOffObservation.hpp>
Public Member Functions | |
GCTAOnOffObservation (void) | |
Void constructor. More... | |
GCTAOnOffObservation (const bool &dummy, const std::string &instrument) | |
Instrument constructor. More... | |
GCTAOnOffObservation (const GObservations &obs) | |
CTA On/Off observation stacking constructor. More... | |
GCTAOnOffObservation (const GPha &pha_on, const GPha &pha_off, const GArf &arf, const GRmf &rmf) | |
CTA observation constructor. More... | |
GCTAOnOffObservation (const GCTAObservation &obs, const GModels &models, const std::string &srcname, const GEbounds &etrue, const GEbounds &ereco, const GSkyRegions &on, const GSkyRegions &off, const bool &use_model_bkg=true) | |
CTA observation constructor (same On and Off observation) More... | |
GCTAOnOffObservation (const GCTAObservation &obs_on, const GCTAObservation &obs_off, const GModels &models, const std::string &srcname, const GEbounds &etrue, const GEbounds &ereco, const GSkyRegions &on, const GSkyRegions &off, const bool &use_model_bkg=true) | |
CTA observation constructor (different On and Off observations) More... | |
GCTAOnOffObservation (const GCTAOnOffObservation &obs) | |
Copy constructor. More... | |
virtual | ~GCTAOnOffObservation (void) |
Destructor. More... | |
GCTAOnOffObservation & | operator= (const GCTAOnOffObservation &obs) |
Assignment operator. More... | |
virtual void | clear (void) |
Clear instance. More... | |
virtual GCTAOnOffObservation * | clone (void) const |
Clone instance. More... | |
virtual std::string | classname (void) const |
Return class name. More... | |
virtual void | response (const GResponse &rsp) |
Set response function. More... | |
virtual const GCTAResponse * | response (void) const |
Return pointer to CTA response function. More... | |
virtual std::string | instrument (void) const |
Return instrument. 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 On/Off observation from an XML element. More... | |
virtual void | write (GXmlElement &xml) const |
write observation to an xml element More... | |
virtual std::string | print (const GChatter &chatter=NORMAL) const |
Print On/Off observation information. More... | |
virtual double | likelihood (const GModels &models, GVector *gradient, GMatrixSparse *curvature, double *npred) const |
Evaluate log-likelihood function for On/Off analysis. More... | |
virtual int | nobserved (void) const |
Return number of observed events. More... | |
void | instrument (const std::string &instrument) |
Set instrument. More... | |
bool | has_response (void) const |
Signal if observation contains response information. More... | |
const GPha & | on_spec (void) const |
Return On spectrum. More... | |
const GPha & | off_spec (void) const |
Return Off spectrum. More... | |
const GArf & | arf (void) const |
Return Auxiliary Response File. More... | |
const GRmf & | rmf (void) const |
Return Redistribution Matrix File. More... | |
GPha | model_gamma (const GModels &models) const |
GPha | model_background (const GModels &models) const |
Public Member Functions inherited from GObservation | |
GObservation (void) | |
Void constructor. More... | |
GObservation (const GObservation &obs) | |
Copy constructor. More... | |
virtual | ~GObservation (void) |
Destructor. More... | |
virtual GObservation & | operator= (const GObservation &obs) |
Assignment operator. More... | |
virtual GEvents * | events (void) |
Return events. More... | |
virtual const GEvents * | events (void) const |
Return events (const version) More... | |
virtual void | events (const GEvents &events) |
Set event container. 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 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 GCTAOnOffObservation &obs) |
Copy class members. More... | |
void | free_members (void) |
Delete class members. More... | |
void | set_exposure (void) |
Set ontime, livetime and deadtime correction factor. More... | |
void | check_consistency (const std::string &method) const |
Check consistency of data members. More... | |
GArf | arf_stacked (const GArf &arf, const GEnergy &emin, const GEnergy &emax) const |
Return ARF for stacking. More... | |
GRmf | rmf_stacked (const GRmf &rmf, const GEnergy &emin, const GEnergy &emax) const |
Return RMF for stacking. More... | |
void | set (const GCTAObservation &obs_on, const GCTAObservation &obs_off, const GModels &models, const std::string &srcname, const GSkyRegions &on, const GSkyRegions &off, const bool &use_model_bkg) |
Set On/Off observation from a CTA observation. More... | |
void | compute_arf (const GCTAObservation &obs, const GModelSpatial &spatial, const GSkyRegions &on) |
Compute ARF of On/Off observation. More... | |
void | compute_arf_cut (const GCTAObservation &obs, const GModelSpatial &spatial, const GSkyRegions &on) |
Compute ARF of On/Off observation for a IRF with radius cut. More... | |
void | compute_bgd (const GCTAObservation &obs, const GSkyRegions &off, const GModels &models, const bool &use_model_bkg) |
Compute background rate in Off regions. More... | |
void | compute_alpha (const GCTAObservation &obs_on, const GCTAObservation &obs_off, const GSkyRegions &on, const GSkyRegions &off, const GModels &models, const bool &use_model_bkg) |
Compute vector of alpha parameters. More... | |
void | compute_rmf (const GCTAObservation &obs, const GSkyRegions &on) |
Compute RMF of On/Off observation. More... | |
void | apply_ebounds (const GCTAObservation &obs) |
Apply CTA observation energy boundaries to On/Off observation. More... | |
double | arf_rad_max (const GCTAObservation &obs, const GSkyRegions &on) const |
Check if effective area IRF has a radius cut. More... | |
double | N_gamma (const GModels &models, const int &ibin, GVector *grad) const |
double | N_bgd (const GModels &models, const int &ibin, GVector *grad) const |
virtual double | likelihood_cstat (const GModels &models, GVector *gradient, GMatrixSparse *curvature, double *npred) const |
Evaluate log-likelihood function for On/Off analysis in the case of Poisson signal with modeled Poisson background. More... | |
virtual double | likelihood_wstat (const GModels &models, GVector *gradient, GMatrixSparse *curvature, double *npred) const |
Evaluate log-likelihood function for On/Off analysis in the case of Poisson signal with measured Poisson background. More... | |
virtual double | wstat_value (const double &non, const double &noff, const double &alpha, const double &ngam, double &nonpred, double &nbgd, double &dlogLdsky, double &d2logLdsky2) const |
Evaluate log-likelihood value in energy bin for On/Off analysis by profiling over the number of background counts. 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... | |
GCTAResponse * | m_response |
Pointer to IRFs. More... | |
double | m_ontime |
Ontime of On observation (seconds) More... | |
double | m_livetime |
Livetime of On observation (seconds) More... | |
double | m_deadc |
Deadtime correction (livetime/ontime) More... | |
GPha | m_on_spec |
On counts spectrum. More... | |
GPha | m_off_spec |
Off counts spectrum. More... | |
GArf | m_arf |
Auxiliary Response Function vector. More... | |
GRmf | m_rmf |
Redistribution matrix. 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... | |
GEvents * | m_events |
Pointer to event container. More... | |
double | m_grad_step_size |
Gradient step size. More... | |
std::vector< std::string > | m_pars_with_gradients |
CTA On/Off observation class.
This class defines a CTA On/Off observation. An On/Off observation is defined by two spectra, one for an On region including the source of interest, and one for an Off region including only background. The response of an On/Off observation is given by the Auxiliary Response File (ARF) and the Redistribution Matrix File (RMF).
The class uses GPha objects to store the On and Off spectra, and GArf and GRmf objects to store the response information. On and Off regions are defined via GSkyRegionMap objects which are essentially finely-pixellized skymaps that allow handling arbitrarily complex shapes
Definition at line 67 of file GCTAOnOffObservation.hpp.
GCTAOnOffObservation::GCTAOnOffObservation | ( | void | ) |
Void constructor.
Constructs empty On/Off observation.
Definition at line 127 of file GCTAOnOffObservation.cpp.
References init_members().
Referenced by clone().
GCTAOnOffObservation::GCTAOnOffObservation | ( | const bool & | dummy, |
const std::string & | instrument | ||
) |
Instrument constructor.
[in] | dummy | Dummy flag. |
[in] | instrument | Instrument name. |
Constructs an empty CTA On/Off 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: "CTAOnOff", "HESSOnOff", "VERITASOnOff", "MAGICOnOff", "ASTRIOnOff".
Definition at line 151 of file GCTAOnOffObservation.cpp.
References init_members(), instrument(), and m_instrument.
GCTAOnOffObservation::GCTAOnOffObservation | ( | const GObservations & | obs | ) |
CTA On/Off observation stacking constructor.
[in] | obs | Observation container. |
GException::invalid_value Incompatible On/Off observation definition.
Constructs On/Off observation by stacking all On/Off observation in the observation container into a single On/Off observation.
The constructor uses the following formulae:
\[ N^{\rm on}(E_{\rm reco}) = \sum_i N^{\rm on}_i(E_{\rm reco}) \]
\[ N^{\rm off}(E_{\rm reco}) = \sum_i N^{\rm off}_i(E_{\rm reco}) \]
\[ \alpha(E_{\rm reco}) = \frac{\sum_i \alpha_i(E_{\rm reco}) B_i(E_{\rm reco}) \tau_off_i} {\sum_i B_i(E_{\rm reco}) \tau_off_i} \]
\[ B(E_{\rm reco}) = \frac{\sum_i B_i(E_{\rm reco}) \tau_off_i}{\sum_i \tau_off_i} \]
\[ ARF(E_{\rm true}) = \frac{\sum_i ARF_i(E_{\rm true}) \tau_on_i} {\sum_i \tau_i} \]
\[ RMF(E_{\rm true}, E_{\rm reco}) = \frac{\sum_i RMF_i(E_{\rm true}, E_{\rm reco}) ARF_i(E_{\rm true}) \tau_on_i} {\sum_i ARF_i(E_{\rm true}) \tau_on_i} \]
where \(N^{\rm on}_i(E_{\rm reco})\) is the On spectrum of observation \(i\), \(N^{\rm off}_i(E_{\rm reco})\) is the Off spectrum of observation \(i\), \(\alpha_i(E_{\rm reco})\) is the background scaling of observation \(i\), \(B_i(E_{\rm reco})\) is the background rate for observation \(i\), \(ARF_i(E_{\rm true})\) is the Auxiliary Response File of observation \(i\), \(RMF_i(E_{\rm true}, E_{\rm reco})\) is the Redistribution Matrix File of observation \(i\), and \(\tau_on_i\) is the livetime of On observation \(i\). \(\tau_off_i\) is the livetime of Off observation \(i\).
The method extracts the instrument name from the first On/Off observation in the observation container and only stacks subsequent On/Off observations with the same instrument name.
Definition at line 391 of file GCTAOnOffObservation.cpp.
References arf(), arf_stacked(), GPha::backscal(), check_consistency(), GArf::ebounds(), GPha::ebounds(), GPha::emax_obs(), GRmf::emeasured(), GPha::emin_obs(), GRmf::etrue(), GPha::exposure(), G_CONSTRUCTOR2, init_members(), instrument(), livetime(), m_arf, m_instrument, m_livetime, m_off_spec, m_on_spec, m_ontime, m_rmf, GRmf::nmeasured(), norm(), GRmf::ntrue(), off_spec(), on_spec(), ontime(), rmf(), rmf_stacked(), GPha::size(), and GObservations::size().
GCTAOnOffObservation::GCTAOnOffObservation | ( | const GPha & | pha_on, |
const GPha & | pha_off, | ||
const GArf & | arf, | ||
const GRmf & | rmf | ||
) |
CTA observation constructor.
[in] | pha_on | On spectrum. |
[in] | pha_off | Off spectrum. |
[in] | arf | Auxiliary Response File. |
[in] | rmf | Redistribution Matrix File. |
Constructs On/Off observation from On and Off spectra, an Auxiliary Response File and a Redistribution Matrix File.
Definition at line 196 of file GCTAOnOffObservation.cpp.
References arf(), check_consistency(), G_CONSTRUCTOR1, init_members(), m_arf, m_off_spec, m_on_spec, m_rmf, rmf(), and set_exposure().
GCTAOnOffObservation::GCTAOnOffObservation | ( | const GCTAObservation & | obs, |
const GModels & | models, | ||
const std::string & | srcname, | ||
const GEbounds & | etrue, | ||
const GEbounds & | ereco, | ||
const GSkyRegions & | on, | ||
const GSkyRegions & | off, | ||
const bool & | use_model_bkg = true |
||
) |
CTA observation constructor (same On and Off observation)
[in] | obs | CTA observation. |
[in] | models | Models including source and background model. |
[in] | srcname | Name of source in models. |
[in] | etrue | True energy boundaries. |
[in] | ereco | Reconstructed energy boundaries. |
[in] | on | On regions. |
[in] | off | Off regions. |
[in] | use_model_bkg | Use model background. |
Constructs On/Off observation from one osbervation by filling the On and Off spectra and computing the Auxiliary Response File (ARF) and Redistribution Matrix File (RMF).
The method requires the specification of the true and reconstructed energy boundaries, as well as the definition of On and Off regions.
The use_model_bkg
flag controls whether a background model should be used for the computations or not. This impacts the computations of the BACKSCAL
column in the On spectrum and the BACKRESP
column in the Off spectrum. See the compute_alpha() and compute_bgd() methods for more information on the applied formulae.
Definition at line 246 of file GCTAOnOffObservation.cpp.
References init_members(), m_arf, m_off_spec, m_on_spec, m_rmf, and set().
GCTAOnOffObservation::GCTAOnOffObservation | ( | const GCTAObservation & | obs_on, |
const GCTAObservation & | obs_off, | ||
const GModels & | models, | ||
const std::string & | srcname, | ||
const GEbounds & | etrue, | ||
const GEbounds & | ereco, | ||
const GSkyRegions & | on, | ||
const GSkyRegions & | off, | ||
const bool & | use_model_bkg = true |
||
) |
CTA observation constructor (different On and Off observations)
[in] | obs_on | On CTA observation. |
[in] | obs_off | Off CTA observation. |
[in] | models | Models including source and background models. |
[in] | srcname | Name of source in models. |
[in] | etrue | True energy boundaries. |
[in] | ereco | Reconstructed energy boundaries. |
[in] | on | On regions. |
[in] | off | Off regions. |
[in] | use_model_bkg | Use model background. |
Constructs On/Off observation from two osbervations by filling the On and Off spectra and computing the Auxiliary Response File (ARF) and Redistribution Matrix File (RMF).
The method requires the specification of the true and reconstructed energy boundaries, as well as the definition of On and Off regions.
The use_model_bkg
flag controls whether a background model should be used for the computations or not. This impacts the computations of the BACKSCAL
column in the On spectrum and the BACKRESP
column in the Off spectrum. See the compute_alpha() and compute_bgd() methods for more information on the applied formulae.
Definition at line 301 of file GCTAOnOffObservation.cpp.
References init_members(), m_arf, m_off_spec, m_on_spec, m_rmf, and set().
GCTAOnOffObservation::GCTAOnOffObservation | ( | const GCTAOnOffObservation & | obs | ) |
Copy constructor.
[in] | obs | On/Off observation. |
Definition at line 171 of file GCTAOnOffObservation.cpp.
References copy_members(), and init_members().
|
virtual |
|
protected |
Apply CTA observation energy boundaries to On/Off observation.
[in] | obs | CTA observation. |
Applies CTA energy boundaries to all histograms used for the On/Off analysis in.
For the PHA On and Off spectra, all bins are set to zero that do not fully overlap with the CTA observation energy boundaries. Specifically, partially overlapping bins are excluded. Some margin is applied that effectively reduces the width of the PHA energy bin, which should cope with any rounding errors.
For the background response, stored in the Off spectrum, all bins are set to zero that do not fully overlap with the CTA observation energy boundaries.
For the RMF, all reconstructued energy bins are set to zero that do not fully overlap with the CTA observation energy boundaries. True energy bins remain unchanged to properly account for energy migration.
Definition at line 2243 of file GCTAOnOffObservation.cpp.
References GEbounds::contains(), GCTAObservation::ebounds(), GEbounds::emax(), GRmf::emeasured(), GEbounds::emin(), GRmf::etrue(), m_off_spec, m_on_spec, m_rmf, and GEbounds::size().
Referenced by set().
|
inline |
Return Auxiliary Response File.
Definition at line 325 of file GCTAOnOffObservation.hpp.
References m_arf.
Referenced by arf_stacked(), GCTAOnOffObservation(), N_gamma(), and read().
|
protected |
Check if effective area IRF has a radius cut.
[in] | obs | CTA observation. |
[in] | on | On regions. |
GException::invalid_argument | IRF has a radius cut that is incompatible with the On regions |
Checks if the effective area IRF has a radius cut. If a radius cut is found the cut value is checked against the radii of the On regions. If the On regions are not circular, or if they have a radius that differs from the IRF cut value, an exception is thrown.
Definition at line 2293 of file GCTAOnOffObservation.cpp.
References abs(), GCTAResponseIrf::aeff(), G_ARF_RADIUS_CUT, GCTAAeff2D::rad_max(), GSkyRegionCircle::radius(), GCTAObservation::response(), GSkyRegions::size(), and gammalib::str().
Referenced by set().
|
protected |
Return ARF for stacking.
[in] | arf | Auxiliary Response File. |
[in] | emin | Minimum observation energy. |
[in] | emax | Maximum observation energy. |
Returns an Auxiliary Response File for stacking that has all elements with true energies outside the interval [emin
, emax
] clipped to zero. This prevents spill over from one observation into another.
Definition at line 1233 of file GCTAOnOffObservation.cpp.
References arf(), GArf::ebounds(), GEbounds::emax(), GEbounds::emin(), and GEbounds::size().
Referenced by GCTAOnOffObservation().
|
protected |
Check consistency of data members.
[in] | method | Calling method. |
GException::invalid_value | Inconsistency found. |
Checks the consistency of data members and throw exceptions in case that inconsistencies are found.
Definition at line 1194 of file GCTAOnOffObservation.cpp.
References GArf::ebounds(), GPha::ebounds(), GRmf::emeasured(), GRmf::etrue(), m_arf, m_off_spec, m_on_spec, and m_rmf.
Referenced by GCTAOnOffObservation(), and read().
|
inlinevirtual |
Return class name.
Implements GObservation.
Definition at line 208 of file GCTAOnOffObservation.hpp.
|
virtual |
Clear instance.
Clears the On/Off observation. All class members will be set to the initial state. Any information that was present in the object before will be lost.
Implements GObservation.
Definition at line 689 of file GCTAOnOffObservation.cpp.
References free_members(), and init_members().
Referenced by read().
|
virtual |
Clone instance.
Returns a pointer to a deep copy of an On/Off observation.
Implements GObservation.
Definition at line 709 of file GCTAOnOffObservation.cpp.
References GCTAOnOffObservation().
|
protected |
Compute vector of alpha parameters.
[in] | obs_on | On CTA observation. |
[in] | obs_off | Off CTA observation. |
[in] | on | On regions. |
[in] | off | Off regions. |
[in] | models | Models including background model. |
[in] | use_model_bkg | Use model background. |
GException::invalid_argument | Observation does not contain relevant response or background information |
Compute the \(\alpha\) parameters for all reconstructed energy bins. The \(\alpha\) parameter gives the ratio between the On and Off region background acceptance multiplied by the On and Off region solid angles.
If use_model_bkg
is true
, the IRF background template will be used for the computation, and \(\alpha(E_{\rm reco})\) is given by
\[ \alpha(E_{\rm reco}) = \frac{\sum_{\rm on} \sum_p {\tt BKG}_{{\rm on},p}(E_{\rm reco}) \times \Omega_{{\rm on},p}} {\sum_{\rm off} \sum_p {\tt BKG}_{{\rm off},p}(E_{\rm reco}) \times \Omega_{{\rm off},p}} \]
where the nominator sums over the On regions, indicated by the index \({\rm on}\), and the denominator sums over the Off regions, indicated by the index \({\rm off}\). Each On or Off region is defined by a region sky map of type GSkyRegionMap, and the pixels of these maps are index by \(p\). \({\tt BKG}_{{\rm on/off},p}(E_{\rm reco})\) is the background rate as provided by the IRF background template in units of \({\rm events} \, {\rm MeV}^{-1} \, {\rm s}^{-1} \, {\rm sr}^{-1}\) for a reconstructed energy \(E_{\rm reco}\) and a pixel index \(p\) in the On or Off region \({\rm on/off}\). \(\Omega_{{\rm on/off},p}\) is the solid angle in units of \({\rm sr}\) of the pixel index \(p\) in the On or Off region \({\rm on/off}\).
If use_model_bkg
is false
, the background acceptance is assumed constant and hence cancels out, which makes the \(\alpha\) parameter independent of reconstructed energy \(E_{\rm reco}\). The \(\alpha\) parameter is then given by
\[ \alpha = \frac{\sum_{\rm on} \sum_p \Omega_{{\rm on},p}} {\sum_{\rm off} \sum_p \Omega_{{\rm off},p}} \]
Definition at line 1922 of file GCTAOnOffObservation.cpp.
References GPha::backscal(), GEbounds::elogmean(), GRmf::emeasured(), GModels::eval(), GCTAPointing::instdir(), GCTAObservation::livetime(), m_on_spec, m_rmf, GSkyRegionMap::map(), GSkyRegionMap::nonzero_indices(), GCTAObservation::pointing(), GSkyRegions::size(), and GEbounds::size().
Referenced by set().
|
protected |
Compute ARF of On/Off observation.
[in] | obs | CTA observation. |
[in] | spatial | Spatial source model. |
[in] | on | On regions. |
GException::invalid_argument | No CTA response found in CTA observation. |
Computes the ARF for an On/Off observation by integration over the IRF for the specified spatial
source model over the on
regions.
All On regions contained in on
are expected to be sky region maps, i.e. of type GSkyRegionMap.
Definition at line 1524 of file GCTAOnOffObservation.cpp.
References GCTAResponseIrf::apply_edisp(), gammalib::cta_rsp_irf(), GCTAEventBin::dir(), GArf::ebounds(), GEbounds::elogmean(), G_COMPUTE_ARF, GResponse::irf_spatial(), m_arf, m_deadc, GSkyRegionMap::map(), GSkyRegionMap::nonzero_indices(), GSkyRegions::size(), and GEbounds::size().
Referenced by set().
|
protected |
Compute ARF of On/Off observation for a IRF with radius cut.
[in] | obs | CTA observation. |
[in] | spatial | Spatial source model. |
[in] | on | On regions. |
GException::invalid_argument | No CTA response found in CTA observation. |
Computes the ARF for an On/Off observation by computing the average effective area over the on
regions for the specified spatial
source model.
All On regions contained in on
are expected to be sky region maps, i.e. of type GSkyRegionMap.
Definition at line 1624 of file GCTAOnOffObservation.cpp.
References GCTAResponseIrf::aeff(), GCTAPointing::azimuth(), gammalib::cta_rsp_irf(), GCTAPointing::dir(), GSkyDir::dist(), GArf::ebounds(), GEbounds::elogmean(), G_COMPUTE_ARF, GEnergy::log10TeV(), m_arf, GSkyRegionMap::map(), GSkyRegionMap::nonzero_indices(), GCTAObservation::pointing(), GSkyDir::posang(), GSkyRegions::size(), GEbounds::size(), sum(), and GCTAPointing::zenith().
Referenced by set().
|
protected |
Compute background rate in Off regions.
[in] | obs | CTA observation. |
[in] | off | Off regions. |
[in] | models | Models including background model. |
[in] | use_model_bkg | Use model background. |
Computes the background rate in units of \({\rm events} \, {\rm MeV}^{-1} \, {\rm s}^{-1}\) in the Off regions and stores the result as additional column with name BACKRESP
in the Off spectrum.
All Off regions contained in off
are expected to be sky region maps, i.e. of type GSkyRegionMap.
If use_model_bkg
is true
, the IRF background template will be used for the computation, and BACKRESP
will be computed using
\[ {\tt BACKRESP}(E_{\rm reco}) = \sum_{\rm off} \sum_p {\tt BKG}_{{\rm off},p}(E_{\rm reco}) \times \Omega_{{\rm off},p} \]
where \({\rm off}\) is the index of the Off region and \(p\) is the pixel in the Off region (note that each Off region is transformed into a region map and \(p\) indicates the pixels of this map). \({\tt BKG}_{{\rm off},p}(E_{\rm reco})\) is the background rate as provided by the IRF background template in units of \({\rm events} \, {\rm MeV}^{-1} \, {\rm s}^{-1} \, {\rm sr}^{-1}\) for a reconstructed energy \(E_{\rm reco}\) and a pixel index \(p\) in the Off region \({\rm off}\). \(\Omega_{{\rm off},p}\) is the solid angle in units of \({\rm sr}\) of the pixel index \(p\) in the Off region \({\rm off}\).
If use_model_bkg
is false
, BACKRESP
will be computed using
\[ {\tt BACKRESP}(E_{\rm reco}) = \sum_{\rm off} \sum_p \Omega_{{\rm off},p} \]
which actually assumes that the background rate is constant over the field of view.
Definition at line 1755 of file GCTAOnOffObservation.cpp.
References GPha::append(), GPha::ebounds(), GEbounds::elogmean(), GModels::eval(), GCTAPointing::instdir(), m_off_spec, GSkyRegionMap::map(), GSkyRegionMap::nonzero_indices(), GCTAObservation::pointing(), GSkyRegions::size(), and GEbounds::size().
Referenced by set().
|
protected |
Compute RMF of On/Off observation.
[in] | obs | CTA observation. |
[in] | on | On regions. |
GException::invalid_argument | Observation does not contain IRF response |
Compute the energy redistribution matrix for an On/Off observation. The method requires that the RMF energy axes have been defined before.
Definition at line 2123 of file GCTAOnOffObservation.cpp.
References GCTAResponseIrf::aeff(), GCTAPointing::azimuth(), gammalib::cta_rsp_irf(), GCTAPointing::dir(), GSkyDir::dist(), GCTAResponseIrf::edisp(), GEbounds::elogmean(), GEbounds::emax(), GRmf::emeasured(), GEbounds::emin(), GRmf::etrue(), G_COMPUTE_RMF, GEnergy::log10TeV(), m_rmf, GSkyRegionMap::map(), GSkyRegionMap::nonzero_indices(), GCTAObservation::pointing(), GSkyDir::posang(), GCTAEdisp::prob_erecobin(), GSkyRegions::size(), GEbounds::size(), and GCTAPointing::zenith().
Referenced by set().
|
protected |
Copy class members.
[in] | obs | CTA On/Off observation. |
Definition at line 1137 of file GCTAOnOffObservation.cpp.
References GCTAResponse::clone(), m_arf, m_deadc, m_instrument, m_livetime, m_off_spec, m_on_spec, m_ontime, m_response, and m_rmf.
Referenced by GCTAOnOffObservation(), and operator=().
Return deadtime correction factor.
[in] | time | Time. |
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 264 of file GCTAOnOffObservation.hpp.
References m_deadc.
|
protected |
Delete class members.
Definition at line 1161 of file GCTAOnOffObservation.cpp.
Referenced by clear(), operator=(), and ~GCTAOnOffObservation().
|
inline |
Signal if observation contains response information.
Definition at line 276 of file GCTAOnOffObservation.hpp.
References m_response.
|
protected |
Initialise class members.
Definition at line 1114 of file GCTAOnOffObservation.cpp.
References GRmf::clear(), GArf::clear(), GPha::clear(), m_arf, m_deadc, m_instrument, m_livetime, m_off_spec, m_on_spec, m_ontime, m_response, and m_rmf.
Referenced by clear(), GCTAOnOffObservation(), and operator=().
|
inlinevirtual |
Return instrument.
Implements GObservation.
Definition at line 220 of file GCTAOnOffObservation.hpp.
References m_instrument.
Referenced by GCTAOnOffObservation(), instrument(), N_gamma(), print(), and set().
|
inline |
Set instrument.
[in] | instrument | Instrument. |
Definition at line 288 of file GCTAOnOffObservation.hpp.
References instrument(), and m_instrument.
|
virtual |
Evaluate log-likelihood function for On/Off analysis.
[in] | models | Models. |
[in,out] | gradient | Pointer to gradients. |
[in,out] | curvature | Pointer to curvature matrix. |
[in,out] | npred | Pointer to Npred value. |
GException::invalid_value | Invalid statistic encountered. |
Reimplemented from GObservation.
Definition at line 914 of file GCTAOnOffObservation.cpp.
References G_LIKELIHOOD, likelihood_cstat(), likelihood_wstat(), GObservation::statistic(), and gammalib::toupper().
|
protectedvirtual |
Evaluate log-likelihood function for On/Off analysis in the case of Poisson signal with modeled Poisson background.
[in] | models | Models. |
[in,out] | gradient | Pointer to gradients. |
[in,out] | curvature | Pointer to curvature matrix. |
[in,out] | npred | Pointer to Npred value. |
GException::invalid_value | There are no model parameters. |
Computes the log-likelihood value for the On/Off observation. The method loops over the energy bins to update the function value, its derivatives and the curvature matrix. The number of On counts \(N_{\rm on}\) and Off counts \(N_{\rm off}\) are taken from the On and Off spectra, the expected number of gamma-ray events \(N_{\gamma}\) and background events \(N_{\rm bgd}\) are computed from the spectral models of the relevant components in the model container (spatial and temporal components are ignored so far). See the N_gamma() and N_bgd() methods for details about the model computations.
The log-likelihood is given by
\[ \ln L = \sum_i \ln L_i \]
where the sum is taken over all energy bins \(i\) and
\[ \ln L_i = - N_{\rm on} \ln N_{\rm pred} + N_{\rm pred} - N_{\rm off} \ln N_{\rm bgd} + N_{\rm bgd} \]
with
\[ N_{\rm pred} = N_{\gamma} + \alpha N_{\rm bgd} \]
being the total number of predicted events for an energy bin in the On region, \(N_{\rm on}\) is the total number of observed events for an energy bin in the On region, \(N_{\rm off}\) is the total number of observed events for an energy bin in the Off region, and \(N_{\rm bgd}\) is the predicted number of background events for an energy bin in the Off region.
The log-likelihood gradient with respect to sky model parameters \(p_{\rm sky}\) is given by
\[ \left( \frac{\partial \ln L_i}{\partial p_{\rm sky}} \right) = \left( 1 - \frac{N_{\rm on}}{N_{\rm pred}} \right) \left( \frac{\partial N_{\gamma}}{\partial p_{\rm sky}} \right) \]
and with respect to background model parameters \(p_{\rm bgd}\) is given by
\[ \left( \frac{\partial \ln L_i}{\partial p_{\rm bgd}} \right) = \left( 1 + \alpha - \frac{N_{\rm off}}{N_{\rm bgd}} - \frac{\alpha N_{\rm on}}{N_{\rm pred}} \right) \left( \frac{\partial N_{\rm bgd}}{\partial p_{\rm bgd}} \right) \]
The curvature matrix elements are given by
\[ \left( \frac{\partial^2 \ln L_i}{\partial^2 p_{\rm sky}} \right) = \left( \frac{N_{\rm on}}{N_{\rm pred}^2} \right) \left( \frac{\partial N_{\gamma}}{\partial p_{\rm sky}} \right)^2 \]
\[ \left( \frac{\partial^2 \ln L_i}{\partial p_{\rm sky} \partial p_{\rm bgd}} \right) = \left( \frac{\alpha N_{\rm on}}{N_{\rm pred}^2} \right) \left( \frac{\partial N_{\gamma}}{\partial p_{\rm sky}} \right) \left( \frac{\partial N_{\rm bgd}}{\partial p_{\rm bgd}} \right) \]
\[ \left( \frac{\partial^2 \ln L_i}{\partial p_{\rm bgd} \partial p_{\rm sky}} \right) = \left( \frac{\alpha N_{\rm on}}{N_{\rm pred}^2} \right) \left( \frac{\partial N_{\rm bgd}}{\partial p_{\rm bgd}} \right) \left( \frac{\partial N_{\gamma}}{\partial p_{\rm sky}} \right) \]
\[ \left( \frac{\partial^2 \ln L_i}{\partial^2 p_{\rm bgd}} \right) = \left( \frac{N_{\rm off}}{N_{\rm bgd}^2} + \frac{\alpha^2 N_{\rm on}}{N_{\rm pred}^2} \right) \left( \frac{\partial N_{\rm bgd}}{\partial p_{\rm bgd}} \right)^2 \]
Definition at line 2744 of file GCTAOnOffObservation.cpp.
References GMatrixSparse::add_to_column(), GPha::backscal(), gammalib::is_infinite(), log(), m_off_spec, m_on_spec, minmod, N_bgd(), N_gamma(), GModels::npars(), and GPha::size().
Referenced by likelihood().
|
protectedvirtual |
Evaluate log-likelihood function for On/Off analysis in the case of Poisson signal with measured Poisson background.
[in] | models | Models. |
[in,out] | gradient | Pointer to gradients. |
[in,out] | curvature | Pointer to curvature matrix. |
[in,out] | npred | Pointer to Npred value. |
Computes the log-likelihood value for the On/Off observation. The method loops over the energy bins to update the function value, its derivatives and the curvature matrix.
The number of On counts \(n_{\rm on}\) and Off counts \(n_{\rm off}\) are taken from the On and Off spectra, the expected number of gamma-ray events \(\mu_s\) is computed from the spectral models of the relevant components in the model container (spatial and temporal components are ignored so far). See the N_gamma() method for details about the model computations.
The estimated number of background counts \(\mu_b\) is derived based on the measurement in the Off region by analytical minimization of the Poisson likelihood, i.e., it is treated as a nuisance parameter. See the wstat_value function for the derivation.
Definition at line 2983 of file GCTAOnOffObservation.cpp.
References GMatrixSparse::add_to_column(), GPha::backscal(), gammalib::is_infinite(), m_off_spec, m_on_spec, N_gamma(), GModels::npars(), GPha::size(), and wstat_value().
Referenced by likelihood().
|
inlinevirtual |
Return livetime.
Implements GObservation.
Definition at line 244 of file GCTAOnOffObservation.hpp.
References m_livetime.
Referenced by GCTAOnOffObservation(), and print().
Definition at line 992 of file GCTAOnOffObservation.cpp.
References GPha::backscal(), GPha::ebounds(), GEbounds::emean(), GPha::fill(), G_MODEL_BACKGROUND, m_off_spec, m_on_spec, N_bgd(), N_gamma(), GModels::npars(), GPha::size(), GObservation::statistic(), gammalib::toupper(), and wstat_value().
Definition at line 957 of file GCTAOnOffObservation.cpp.
References GPha::ebounds(), GEbounds::emean(), GPha::fill(), m_on_spec, N_gamma(), GModels::npars(), and GPha::size().
|
protected |
Definition at line 2550 of file GCTAOnOffObservation.cpp.
References gammalib::cta_model_spectral(), GPha::ebounds(), GEbounds::elogmean(), GModelSpectral::eval(), GEbounds::ewidth(), GPha::exposure(), GOptimizerPar::factor_gradient(), GObservation::id(), GOptimizerPar::is_free(), GModel::is_valid(), m_off_spec, m_on_spec, GEnergy::MeV(), norm(), GModels::npars(), GPha::size(), GModel::size(), and GModels::size().
Referenced by likelihood_cstat(), and model_background().
|
protected |
Definition at line 2381 of file GCTAOnOffObservation.cpp.
References arf(), GArf::ebounds(), GEbounds::elogmean(), GEbounds::emax(), GEbounds::emin(), GModelSpectral::eval(), GEbounds::ewidth(), GPha::exposure(), GOptimizerPar::factor_gradient(), GModelSpectral::flux(), GModelSpectral::has_par(), GModel::has_scales(), instrument(), GOptimizerPar::is_free(), GModel::is_valid(), m_arf, m_on_spec, m_rmf, GEnergy::MeV(), GOptimizerPar::name(), GModels::npars(), rmf(), GOptimizerPar::scale(), GModel::scale(), GArf::size(), GPha::size(), GModel::size(), GModels::size(), GModelSky::spectral(), and GOptimizerPar::value().
Referenced by likelihood_cstat(), likelihood_wstat(), model_background(), and model_gamma().
|
inlinevirtual |
Return number of observed events.
Reimplemented from GObservation.
Definition at line 349 of file GCTAOnOffObservation.hpp.
References GPha::counts(), and m_on_spec.
|
inline |
Return Off spectrum.
Definition at line 313 of file GCTAOnOffObservation.hpp.
References m_off_spec.
Referenced by GCTAOnOffObservation().
|
inline |
Return On spectrum.
Definition at line 301 of file GCTAOnOffObservation.hpp.
References m_on_spec.
Referenced by GCTAOnOffObservation().
|
inlinevirtual |
Return ontime.
Implements GObservation.
Definition at line 232 of file GCTAOnOffObservation.hpp.
References m_ontime.
Referenced by GCTAOnOffObservation(), and print().
GCTAOnOffObservation & GCTAOnOffObservation::operator= | ( | const GCTAOnOffObservation & | obs | ) |
Assignment operator.
[in] | obs | On/Off observation. |
Assigns one On/Off observation to another On/Off observation object.
Definition at line 652 of file GCTAOnOffObservation.cpp.
References copy_members(), free_members(), init_members(), and GObservation::operator=().
Print On/Off observation information.
[in] | chatter | Chattiness. |
Implements GObservation.
Definition at line 1070 of file GCTAOnOffObservation.cpp.
References instrument(), livetime(), m_arf, m_deadc, GObservation::m_id, GObservation::m_name, m_off_spec, m_on_spec, m_rmf, ontime(), gammalib::parformat(), GRmf::print(), GArf::print(), GPha::print(), gammalib::reduce(), SILENT, GObservation::statistic(), and gammalib::str().
|
virtual |
Read On/Off observation from an XML element.
[in] | xml | XML element. |
GException::invalid_value | Invalid statistic attribute encountered |
Reads information for an On/Off observation from an XML element. The expected format of the XML element is
<observation name="..." id="..." instrument="..."> <parameter name="Pha_on" file="..."/> <parameter name="Pha_off" file="..."/> <parameter name="Arf" file="..."/> <parameter name="Rmf" file="..."/> </observation>
Optionally, the statistic used for maximum likelihood fitting can be specified:
<observation name="..." id="..." instrument="..." statistic="...">
Implements GObservation.
Definition at line 794 of file GCTAOnOffObservation.cpp.
References arf(), GXmlElement::attribute(), check_consistency(), clear(), G_READ, GRmf::load(), GArf::load(), GPha::load(), m_arf, m_instrument, m_off_spec, m_on_spec, m_rmf, rmf(), set_exposure(), GObservation::statistic(), gammalib::toupper(), gammalib::xml_file_expand(), and gammalib::xml_get_attr().
|
virtual |
Set response function.
[in] | rsp | Response function. |
GException::invalid_argument | Invalid response class specified. |
Sets the response function for the On/Off observation.
Implements GObservation.
Definition at line 725 of file GCTAOnOffObservation.cpp.
References GCTAResponse::clone(), gammalib::cta_rsp(), G_RESPONSE_SET, and m_response.
|
virtual |
Return pointer to CTA response function.
GException::invalid_value | No 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 756 of file GCTAOnOffObservation.cpp.
References G_RESPONSE_GET, and m_response.
|
inline |
Return Redistribution Matrix File.
Definition at line 337 of file GCTAOnOffObservation.hpp.
References m_rmf.
Referenced by GCTAOnOffObservation(), N_gamma(), read(), and rmf_stacked().
|
protected |
Return RMF for stacking.
[in] | rmf | Redistribution Matrix File. |
[in] | emin | Minimum observation energy. |
[in] | emax | Maximum observation energy. |
Returns a Redistribution Matrix File for stacking that has all matrix elements with true energies outside the interval [emin
, emax
] clipped to zero. This prevents spill over from one observation into another.
To correct for the missing events at the edges, the Redistribution Matrix File is renormalized so that for each reconstructed energy the sum over the true energies is the same as before the clipping.
Definition at line 1272 of file GCTAOnOffObservation.cpp.
References GEbounds::emax(), GRmf::emeasured(), GEbounds::emin(), GRmf::etrue(), norm(), rmf(), GEbounds::size(), and sum().
Referenced by GCTAOnOffObservation().
|
protected |
Set On/Off observation from a CTA observation.
[in] | obs_on | CTA On observation. |
[in] | obs_off | CTA Off observation. |
[in] | models | Models including source and background model. |
[in] | srcname | Name of soucre in models. |
[in] | on | On regions. |
[in] | off | Off regions. |
[in] | use_model_bkg | Use model background. |
GException::invalid_value | No CTA event list found in CTA observation. |
Sets an On/Off observation from a CTA observation by filling the events that fall in the On and Off regions into the PHA spectra and by computing the corresponding ARF and RMF response functions.
The use_model_bkg
flags controls whether the background model should be used for computations or not. This impacts the computations of the BACKSCAL
column in the On spectrum and the BACKRESP
column in the Off spectrum. See the compute_alpha() and compute_bgd() methods for more information on the applied formulae.
Definition at line 1343 of file GCTAOnOffObservation.cpp.
References GSkyRegions::append(), GFitsHeader::append(), GModels::append(), apply_ebounds(), arf_rad_max(), GCTAResponseIrf::caldb(), GModels::classname(), compute_alpha(), compute_arf(), compute_arf_cut(), compute_bgd(), compute_rmf(), GSkyRegions::contains(), GCTAObservation::deadc(), GCTAEventAtom::dir(), GCTAInstDir::dir(), GCTAObservation::ebounds(), GEbounds::emax(), GPha::emax_obs(), GEbounds::emin(), GPha::emin_obs(), GCTAEventAtom::energy(), GObservation::events(), GPha::exposure(), GPha::fill(), G_SET, GRmf::header(), GArf::header(), GPha::header(), GObservation::id(), GCTAObservation::instrument(), instrument(), GCTAObservation::livetime(), m_arf, m_deadc, m_instrument, m_livetime, m_off_spec, m_on_spec, m_ontime, m_rmf, GObservation::model(), GObservation::name(), GCTAObservation::ontime(), GCTAObservation::response(), GCTAResponseIrf::rspname(), GSkyRegions::size(), GCTAEventList::size(), GModels::size(), GModelSky::spatial(), and gammalib::toupper().
Referenced by GCTAOnOffObservation().
|
protected |
Set ontime, livetime and deadtime correction factor.
Definition at line 1171 of file GCTAOnOffObservation.cpp.
References GPha::exposure(), m_deadc, m_livetime, m_on_spec, and m_ontime.
Referenced by GCTAOnOffObservation(), and read().
|
virtual |
write observation to an xml element
[in] | xml | XML element. |
Writes information for an On/Off observation into an XML element. The expected format of the XML element is
<observation name="..." id="..." instrument="..." statistic="..."> <parameter name="Pha_on" file="..."/> <parameter name="Pha_off" file="..."/> <parameter name="Arf" file="..."/> <parameter name="Rmf" file="..."/> </observation>
The actual files described in the XML elements are not written.
Implements GObservation.
Definition at line 871 of file GCTAOnOffObservation.cpp.
References GXmlElement::attribute(), GRmf::filename(), GArf::filename(), GPha::filename(), G_WRITE, m_arf, m_off_spec, m_on_spec, m_rmf, GObservation::statistic(), gammalib::xml_file_reduce(), and gammalib::xml_need_par().
|
protectedvirtual |
Evaluate log-likelihood value in energy bin for On/Off analysis by profiling over the number of background counts.
[in] | non | number of On counts |
[in] | noff | number of Off counts |
[in] | alpha | background scaling rate |
[in] | ngam | number of predicted gamma-ray counts |
[in,out] | nonpred | number of predicted On counts |
[in,out] | nbgd | number of predicted background counts |
[in,out] | dlogLdsky | first derivative of log-like w.r.t. sky pars |
[in,out] | d2logLdsky2 | second derivative of log-like w.r.t. sky pars |
Computes the log-likelihood value for the On/Off observation in an energy bin by treating the number of background counts as nuisance parameter. The method performs an analytical minimisation of the Poisson likelihood and updates the relevant values. In the general case, the log-likelihood function is computed using
\[ L = \mu_s + (1+\alpha) \mu_b - n_{\rm on} - n_{\rm off} - n_{\rm on} \left( \ln(\mu_s + \alpha \mu_b) - \ln(n_{\rm on}) \right) - n_{\rm off} \left( \ln(\mu_b) - \ln(n_{\rm off}) \right) \]
where
\[ \mu_b = \frac{C+D}{2\alpha(1+\alpha)} \]
are the estimated number of background counts with
\[ C = \alpha (n_{\rm on} + n_{\rm off}) - (1 + \alpha) \mu_s \]
and
\[ D = \sqrt{C^2 + 4 (1 + \alpha) \, \alpha \, n_{\rm off} \, \mu_s} \]
\
The first derivative of the log-likelihood function is given by
\[ \frac{\delta L}{\delta \mu_s} = 1 + (1 + \alpha) \frac{\delta \mu_b}{\delta \mu_s} - \frac{n_{\rm on}}{\mu_s + \alpha \mu_b} \left( 1 + \alpha \frac{\delta \mu_b}{\delta \mu_s} \right) - \frac{n_{\rm off}}{\mu_b} \frac{\delta \mu_b}{\delta \mu_s} \]
with
\[ \frac{\delta \mu_b}{\delta \mu_s} = \frac{n_{\rm off} - C}{D} - \frac{1}{2 \alpha} \]
The second derivative of the log-likelihood function is given by
\[ \frac{\delta^2 L}{\delta \mu_s^2} = \frac{n_{\rm on}}{(\mu_s + \alpha \mu_b)^2} \left( 1 + \alpha \frac{\delta \mu_b}{\delta \mu_s} \right) + \frac{\delta^2 \mu_b}{\delta \mu_s^2} \left( (1 + \alpha) - \frac{\alpha n_{\rm on}}{\mu_s + \alpha \mu_b} - \frac{n_{\rm off}}{\mu_b} \right) + \frac{\delta \mu_b}{\delta \mu_s} \left( \frac{\alpha n_{\rm on}}{(\mu_s + \alpha \mu_b)^2} \left( 1 + \alpha \frac{\delta \mu_b}{\delta \mu_s} \right) + \frac{n_{\rm off}}{\mu_b^2} \frac{\delta \mu_b}{\delta \mu_s} \right) \]
with
\[ \frac{\delta^2 \mu_b}{\delta \mu_s^2} = \frac{-1}{2 \alpha} \left( \frac{1}{D} \frac{\delta C}{\delta \mu_s} + \frac{2 \alpha \, n_{\rm off} - C}{D^2} \frac{\delta D}{\delta \mu_s} \right) \]
where
\[ \frac{\delta C}{\delta \mu_s} = -(1 + \alpha) \]
and
\[ \frac{\delta D}{\delta \mu_s} = \frac{4 (1 + \alpha) \, \alpha \, n_{\rm off} - 2 \, C \, (1 + \alpha)} {2D} \]
In the special case \(n_{\rm on}=n_{\rm off}=0\) the formal background estimate becomes negative, hence we set \(\mu_b=0\) and the log-likelihood function becomes
\[ L = \mu_s \]
the first derivative
\[ \frac{\delta L}{\delta \mu_s} = 1 \]
and the second derivative
\[ \frac{\delta^2 L}{\delta \mu_s^2} = 0 \]
In the special case \(n_{\rm on}=0\) and \(n_{\rm off}>0\) the log-likelihood function becomes
\[ L = \mu_s + n_{\rm off} \ln(1 + \alpha) \]
the background estimate
\[ \mu_b = \frac{n_{\rm off}}{1+\alpha} \]
the first derivative
\[ \frac{\delta L}{\delta \mu_s} = 1 \]
and the second derivative
\[ \frac{\delta^2 L}{\delta \mu_s^2} = 0 \]
In the special case \(n_{\rm on}>0\) and \(n_{\rm off}=0\) the background estimate becomes
\[ \mu_b = \frac{n_{\rm on}}{1+\alpha} - \frac{\mu_s}{\alpha} \]
which is positive for
\[ \mu_s < n_{\rm on} \frac{\alpha}{1+\alpha} \]
For positive \(\mu_b\) the log-likelihood function is given by
\[ L = -\frac{\mu_s}{\alpha} - n_{\rm on} \ln \left(\frac{\alpha}{1 + \alpha} \right) \]
the first derivative
\[ \frac{\delta L}{\delta \mu_s} = -\frac{1}{\alpha} \]
and the second derivative
\[ \frac{\delta^2 L}{\delta \mu_s^2} = 0 \]
For negative \(\mu_b\) we set \(\mu_b=0\) and the log-likelihood function becomes
\[ L = \mu_s - n_{\rm on} - n_{\rm on} \left( \ln(\mu_s) - \ln(n_{\rm on}) \right) \]
the first derivative
\[ \frac{\delta L}{\delta \mu_s} = 1 - \frac{n_{\rm on}}{\mu_s} \]
and the second derivative
\[ \frac{\delta^2 L}{\delta \mu_s^2} = \frac{1}{\mu_s^2} \]
Note that the fit results may be biased and the statistical errors overestimated if for some bins \(n_{\rm on}=0\) and/or \(n_{\rm off}=0\) (i.e. if the special cases are encountered).
Definition at line 3343 of file GCTAOnOffObservation.cpp.
References gammalib::is_infinite(), gammalib::is_notanumber(), log(), and sqrt().
Referenced by likelihood_wstat(), and model_background().
|
protected |
Auxiliary Response Function vector.
Definition at line 197 of file GCTAOnOffObservation.hpp.
Referenced by arf(), check_consistency(), compute_arf(), compute_arf_cut(), copy_members(), GCTAOnOffObservation(), init_members(), N_gamma(), print(), read(), set(), and write().
|
protected |
Deadtime correction (livetime/ontime)
Definition at line 194 of file GCTAOnOffObservation.hpp.
Referenced by compute_arf(), copy_members(), deadc(), init_members(), print(), set(), and set_exposure().
|
protected |
Instrument name.
Definition at line 190 of file GCTAOnOffObservation.hpp.
Referenced by copy_members(), GCTAOnOffObservation(), init_members(), instrument(), read(), and set().
|
protected |
Livetime of On observation (seconds)
Definition at line 193 of file GCTAOnOffObservation.hpp.
Referenced by copy_members(), GCTAOnOffObservation(), init_members(), livetime(), set(), and set_exposure().
|
protected |
Off counts spectrum.
Definition at line 196 of file GCTAOnOffObservation.hpp.
Referenced by apply_ebounds(), check_consistency(), compute_bgd(), copy_members(), GCTAOnOffObservation(), init_members(), likelihood_cstat(), likelihood_wstat(), model_background(), N_bgd(), off_spec(), print(), read(), set(), and write().
|
protected |
On counts spectrum.
Definition at line 195 of file GCTAOnOffObservation.hpp.
Referenced by apply_ebounds(), check_consistency(), compute_alpha(), copy_members(), GCTAOnOffObservation(), init_members(), likelihood_cstat(), likelihood_wstat(), model_background(), model_gamma(), N_bgd(), N_gamma(), nobserved(), on_spec(), print(), read(), set(), set_exposure(), and write().
|
protected |
Ontime of On observation (seconds)
Definition at line 192 of file GCTAOnOffObservation.hpp.
Referenced by copy_members(), GCTAOnOffObservation(), init_members(), ontime(), set(), and set_exposure().
|
protected |
Pointer to IRFs.
Definition at line 191 of file GCTAOnOffObservation.hpp.
Referenced by copy_members(), has_response(), init_members(), and response().
|
protected |
Redistribution matrix.
Definition at line 198 of file GCTAOnOffObservation.hpp.
Referenced by apply_ebounds(), check_consistency(), compute_alpha(), compute_rmf(), copy_members(), GCTAOnOffObservation(), init_members(), N_gamma(), print(), read(), rmf(), set(), and write().