GammaLib 2.1.0.dev
Loading...
Searching...
No Matches
GObservations.cpp
Go to the documentation of this file.
1/***************************************************************************
2 * GObservations.cpp - Observation container class *
3 * ----------------------------------------------------------------------- *
4 * copyright (C) 2009-2024 by Juergen Knoedlseder *
5 * ----------------------------------------------------------------------- *
6 * *
7 * This program is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * This program is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have received a copy of the GNU General Public License *
18 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
19 * *
20 ***************************************************************************/
21/**
22 * @file GObservations.cpp
23 * @brief Observations container class implementation
24 * @author Juergen Knoedlseder
25 */
26
27/* __ Includes ___________________________________________________________ */
28#ifdef HAVE_CONFIG_H
29#include <config.h>
30#endif
31#include "GTools.hpp"
32#include "GException.hpp"
33#include "GFilename.hpp"
34#include "GXml.hpp"
35#include "GModel.hpp"
36#include "GObservations.hpp"
38#include "GMatrixSparse.hpp"
39#include "GOptimizer.hpp"
40
41/* __ Method name definitions ____________________________________________ */
42#define G_AT "GObservations::at(int&)"
43#define G_SET "GObservations::set(int&, GObservation&)"
44#define G_APPEND "GObservations::append(GObservation&)"
45#define G_INSERT "GObservations::insert(int&, GObservation&)"
46#define G_REMOVE "GObservations::remove(int&)"
47#define G_EXTEND "GObservations::extend(GObservations&)"
48#define G_READ "GObservations::read(GXml&)"
49#define G_NPRED "GObservations::npred(std::string&, std::string&)"
50
51/* __ Macros _____________________________________________________________ */
52
53/* __ Coding definitions _________________________________________________ */
54
55/* __ Debug definitions __________________________________________________ */
56
57
58/*==========================================================================
59 = =
60 = Constructors/destructors =
61 = =
62 ==========================================================================*/
63
64/***********************************************************************//**
65 * @brief Void constructor
66 ***************************************************************************/
68{
69 // Initialise members
71
72 // Return
73 return;
74}
75
76
77/***********************************************************************//**
78 * @brief Copy constructor
79 *
80 * @param obs Observation container.
81 ***************************************************************************/
83{
84 // Initialise members
86
87 // Copy members
88 copy_members(obs);
89
90 // Return
91 return;
92}
93
94
95/***********************************************************************//**
96 * @brief Load constructor
97 *
98 * @param[in] filename XML filename.
99 *
100 * Construct the observation container by loading all relevant information
101 * from a XML file. Please refer to the read() method for more information
102 * about the structure of the XML file.
103 ***************************************************************************/
105{
106 // Initialise members
107 init_members();
108
109 // Load XML file
110 load(filename);
111
112 // Return
113 return;
114}
115
116
117/***********************************************************************//**
118 * @brief Destructor
119 ***************************************************************************/
121{
122 // Free members
123 free_members();
124
125 // Return
126 return;
127}
128
129
130/*==========================================================================
131 = =
132 = Operators =
133 = =
134 ==========================================================================*/
135
136/***********************************************************************//**
137 * @brief Assignment operator
138 *
139 * @param[in] obs Observation container.
140 * @return Observation container.
141 ***************************************************************************/
143{
144 // Execute only if object is not identical
145 if (this != &obs) {
146
147 // Free members
148 free_members();
149
150 // Initialise members
151 init_members();
152
153 // Copy members
154 copy_members(obs);
155
156 } // endif: object was not identical
157
158 // Return this object
159 return *this;
160}
161
162
163/*==========================================================================
164 = =
165 = Public methods =
166 = =
167 ==========================================================================*/
168
169/***********************************************************************//**
170 * @brief Clear observations
171 ***************************************************************************/
173{
174 // Free members
175 free_members();
176
177 // Initialise members
178 init_members();
179
180 // Return
181 return;
182}
183
184
185/***********************************************************************//**
186 * @brief Clone observations
187 *
188 * @return Pointer to deep copy of observation container.
189 ***************************************************************************/
191{
192 // Clone observations
193 return new GObservations(*this);
194}
195
196
197/***********************************************************************//**
198 * @brief Return pointer to observation
199 *
200 * @param[in] index Observation index [0,...,size()[.
201 * @return Pointer to observation.
202 *
203 * @exception GException::out_of_range
204 * Operation index is out of range.
205 *
206 * Returns a pointer to the observation with specified @p index.
207 ***************************************************************************/
209{
210 // Raise exception if index is out of range
211 if (index < 0 || index >= size()) {
212 throw GException::out_of_range(G_AT, "Observation index",
213 index, size());
214 }
215
216 // Return pointer
217 return m_obs[index];
218}
219
220
221/***********************************************************************//**
222 * @brief Return pointer to observation (const version)
223 *
224 * @param[in] index Observation index [0,...,size()[.
225 * @return Pointer to observation.
226 *
227 * @exception GException::out_of_range
228 * Operation index is out of range.
229 *
230 * Returns a const pointer to the observation with specified @p index.
231 ***************************************************************************/
232const GObservation* GObservations::at(const int& index) const
233{
234 // Raise exception if index is out of range
235 if (index < 0 || index >= size()) {
236 throw GException::out_of_range(G_AT, "Observation index",
237 index, size());
238 }
239
240 // Return pointer
241 return m_obs[index];
242}
243
244
245/***********************************************************************//**
246 * @brief Set observation in container
247 *
248 * @param[in] index Observation index [0,...,size()[.
249 * @param[in] obs Observation.
250 * @return Pointer to deep copy of observation.
251 *
252 * @exception GException::out_of_range
253 * Observation index is out of range.
254 * @exception GException::invalid_value
255 * Observation with same instrument and identifier already
256 * exists in container.
257 *
258 * Set a deep copy and observation @p obs at the specified @p index in the
259 * container.
260 ***************************************************************************/
261GObservation* GObservations::set(const int& index, const GObservation& obs)
262{
263 // Compile option: raise exception if index is out of range
264 #if defined(G_RANGE_CHECK)
265 if (index < 0 || index >= size()) {
266 throw GException::out_of_range(G_SET, "Observation index",
267 index, size());
268 }
269 #endif
270
271 // Raise an exception if an observation with specified instrument and
272 // identifier already exists
273 int inx = get_index(obs.instrument(), obs.id());
274 if (inx != -1 && inx != index) {
275 std::string msg =
276 "Attempt to set \""+obs.instrument()+"\" observation with"
277 " identifier \""+obs.id()+"\" in observation container at"
278 " index "+gammalib::str(index)+", but an observation with the"
279 " same attributes exists already at index "+gammalib::str(inx)+
280 " in the container.\n"
281 "Every observation for a given instrument in the observation"
282 " container needs a unique identifier.";
284 }
285
286 // Delete existing observation
287 if (m_obs[index] != NULL) delete m_obs[index];
288
289 // Store pointer to a deep copy of the observation
290 m_obs[index] = obs.clone();
291
292 // Return pointer
293 return m_obs[index];
294}
295
296
297/***********************************************************************//**
298 * @brief Append observation to container
299 *
300 * @param[in] obs Observation.
301 * @return Pointer to deep copy of observation.
302 *
303 * @exception GException::invalid_value
304 * Observation with same instrument and identifier already
305 * exists in container.
306 *
307 * Appends a deep copy of an observation to the container.
308 ***************************************************************************/
310{
311 // Raise an exception if an observation with specified instrument and
312 // identifier already exists
313 int inx = get_index(obs.instrument(), obs.id());
314 if (inx != -1) {
315 std::string msg =
316 "Attempt to append \""+obs.instrument()+"\" observation with"
317 " identifier \""+obs.id()+"\" to observation container, but an"
318 " observation with the same attributes exists already at"
319 " index "+gammalib::str(inx)+" in the container.\n"
320 "Every observation for a given instrument in the observation"
321 " container needs a unique identifier.";
323 }
324
325 // Clone observation
326 GObservation* ptr = obs.clone();
327
328 // Append to list
329 m_obs.push_back(ptr);
330
331 // Return pointer
332 return ptr;
333}
334
335
336/***********************************************************************//**
337 * @brief Insert observation into container
338 *
339 * @param[in] index Observation index [0,...,size()[.
340 * @param[in] obs Observation.
341 * @return Pointer to deep copy of observation.
342 *
343 * @exception GException::out_of_range
344 * Observation index is out of range.
345 * @exception GException::invalid_value
346 * Observation with same instrument and identifier already
347 * exists in container.
348 *
349 * Inserts a deep copy of an observation into the container before the
350 * observation with the specified @p index.
351 ***************************************************************************/
352GObservation* GObservations::insert(const int& index, const GObservation& obs)
353{
354 // Compile option: raise exception if index is out of range
355 #if defined(G_RANGE_CHECK)
356 if (is_empty()) {
357 if (index > 0) {
358 throw GException::out_of_range(G_INSERT, "Observation index",
359 index, size());
360 }
361 }
362 else {
363 if (index < 0 || index >= size()) {
364 throw GException::out_of_range(G_INSERT, "Observation index",
365 index, size());
366 }
367 }
368 #endif
369
370 // Raise an exception if an observation with specified instrument and
371 // identifier already exists
372 int inx = get_index(obs.instrument(), obs.id());
373 if (inx != -1 && inx != index) {
374 std::string msg =
375 "Attempt to insert \""+obs.instrument()+"\" observation with"
376 " identifier \""+obs.id()+"\" in observation container before"
377 " index "+gammalib::str(index)+", but an observation with the"
378 " same attributes exists already at index "+gammalib::str(inx)+
379 " in the container.\n"
380 "Every observation for a given instrument in the observation"
381 " container needs a unique identifier.";
383 }
384
385 // Clone observation
386 GObservation* ptr = obs.clone();
387
388 // Clone observation and insert into list
389 m_obs.insert(m_obs.begin()+index, ptr);
390
391 // Return pointer
392 return ptr;
393}
394
395
396/***********************************************************************//**
397 * @brief Remove observation from container
398 *
399 * @param[in] index Observation index [0,...,size()[.
400 *
401 * @exception GException::out_of_range
402 * Observation index is out of range.
403 *
404 * Removes observation of specified @p index from the container.
405 ***************************************************************************/
406void GObservations::remove(const int& index)
407{
408 // Compile option: If index is outside boundary then raise exception
409 #if defined(G_RANGE_CHECK)
410 if (index < 0 || index >= size()) {
411 throw GException::out_of_range(G_REMOVE, "Observation index",
412 index, size());
413 }
414 #endif
415
416 // Delete observation
417 delete m_obs[index];
418
419 // Erase observation component from container
420 m_obs.erase(m_obs.begin() + index);
421
422 // Return
423 return;
424}
425
426
427/***********************************************************************//**
428 * @brief Append observations from observation container
429 *
430 * @param[in] obs Observations.
431 *
432 * @exception GException::invalid_value
433 * Observation with same instrument and identifier already
434 * exists in container.
435 *
436 * Appends deep copies of observations to the container.
437 ***************************************************************************/
439{
440 // Do nothing if observation container is empty
441 if (!obs.is_empty()) {
442
443 // Get size. Note that we extract the size first to avoid an
444 // endless loop that arises when a container is appended to
445 // itself.
446 int num = obs.size();
447
448 // Reserve enough space
449 reserve(size() + num);
450
451 // Loop over all observations, clone them and append them to the
452 // list
453 for (int i = 0; i < num; ++i) {
454
455 // Raise an exception if an observation with specified
456 // instrument and identifier already exists
457 int inx = get_index(obs[i]->instrument(), obs[i]->id());
458 if (inx != -1) {
459 std::string msg =
460 "Attempt to append \""+obs[i]->instrument()+"\""
461 " observation with identifier \""+obs[i]->id()+"\""
462 " to observation container, but an observation with the"
463 " same attributes exists already at index "+
464 gammalib::str(inx)+" in the container.\n"
465 "Every observation for a given instrument in the"
466 " observation container needs a unique identifier.";
468 }
469
470 // Append observation to container
471 m_obs.push_back(obs[i]->clone());
472 }
473
474 } // endif: observation container was not empty
475
476 // Return
477 return;
478}
479
480
481/***********************************************************************//**
482 * @brief Signals if observation exists
483 *
484 * @param[in] instrument Instrument.
485 * @param[in] id Observation identifier.
486 * @return True if observation with specified @p instrument and identifier
487 * @p id exists.
488 *
489 * Searches all observations for a match with the specified @p instrument
490 * and @p identifier. If the specified attributes have been found, true is
491 * returned.
492 ***************************************************************************/
493bool GObservations::contains(const std::string& instrument,
494 const std::string& id) const
495{
496 // Get observation index
497 int index = get_index(instrument, id);
498
499 // Return
500 return (index != -1);
501}
502
503
504/***********************************************************************//**
505 * @brief Load observations from XML file
506 *
507 * @param[in] filename XML filename.
508 *
509 * Loads observation from a XML file into the container. Please refer to the
510 * read() method for more information about the structure of the XML file.
511 ***************************************************************************/
512void GObservations::load(const GFilename& filename)
513{
514 // Clear any existing observations
515 clear();
516
517 // Load XML document
518 GXml xml(filename.url());
519
520 // Read observations from XML document
521 read(xml);
522
523 // Return
524 return;
525}
526
527
528/***********************************************************************//**
529 * @brief Save observations into XML file
530 *
531 * @param[in] filename XML filename.
532 *
533 * Saves observations into a XML file. Please refer to the read() method for
534 * more information about the structure of the XML file.
535 ***************************************************************************/
536void GObservations::save(const GFilename& filename) const
537{
538 // Declare empty XML document
539 GXml xml;
540
541 // Write observations into XML file
542 write(xml);
543
544 // Save XML document
545 xml.save(filename);
546
547 // Return
548 return;
549}
550
551
552/***********************************************************************//**
553 * @brief Read observations from XML document
554 *
555 * @param[in] xml XML document.
556 *
557 * @exception GException::invalid_value
558 * Invalid instrument encountered in XML file.
559 *
560 * Reads observations from the first observation list that is found in the
561 * XML document. The decoding of the instrument specific observation
562 * definition is done within the observation's GObservation::read() method.
563 * The following file structure is expected:
564 *
565 * <observation_list title="observation library">
566 * <observation name="..." id="..." instrument="...">
567 * ...
568 * </observation>
569 * <observation name="..." id="..." instrument="...">
570 * ...
571 * </observation>
572 * ...
573 * </observation_list>
574 *
575 * The @p name and @p id attributes allow for a unique identification of an
576 * observation within the observation container. The @p instrument
577 * attributes specifies the instrument to which the observation applies.
578 * This attribute will be used to allocate the appropriate instrument
579 * specific GObservation class variant by using the GObservationRegistry
580 * class.
581 *
582 * The structure within the @p observation tag is defined by the instrument
583 * specific GObservation class.
584 *
585 * @todo Observation names and IDs are not verified so far for uniqueness.
586 * This would be required to achieve an unambiguous update of parameters
587 * in an already existing XML file when using the write method.
588 ***************************************************************************/
589void GObservations::read(const GXml& xml)
590{
591 // Get pointer on observation library
592 const GXmlElement* lib = xml.element("observation_list", 0);
593
594 // Loop over all observations
595 int n = lib->elements("observation");
596 for (int i = 0; i < n; ++i) {
597
598 // Get pointer on observation
599 const GXmlElement* obs = lib->element("observation", i);
600
601 // Get attributes
602 std::string name = obs->attribute("name");
603 std::string id = obs->attribute("id");
604 std::string instrument = obs->attribute("instrument");
605
606 // Allocate instrument specific observation
607 GObservationRegistry registry;
608 GObservation* ptr = registry.alloc(instrument);
609
610 // If observation is valid then read its definition from XML file
611 if (ptr != NULL) {
612
613 // Read definition
614 ptr->read(*obs);
615
616 // Set attributes
617 ptr->name(name);
618 ptr->id(id);
619
620 } // endif: observation was valid
621
622 // ... otherwise throw an exception
623 else {
624 std::string msg = "Instrument \""+instrument+"\" unknown. The "
625 "following instruments are available: "+
626 registry.content()+". Please specify one of "
627 "the available instruments.";
629 }
630
631 // Append observation to container
632 append(*ptr);
633
634 // Free observation (the append method made a deep copy)
635 delete ptr;
636
637 } // endfor: looped over all observations
638
639 // Return
640 return;
641}
642
643
644/***********************************************************************//**
645 * @brief Write observations into XML document
646 *
647 * @param[in] xml XML document.
648 *
649 * Write observations into the first observation library that is found in the
650 * XML document. In case that no observation library exists, one is added to
651 * the document. Please refer to the read() method for more information
652 * about the structure of the XML document.
653 ***************************************************************************/
655{
656 // If there is no observation library then append one
657 if (xml.elements("observation_list") == 0) {
658 xml.append(GXmlElement("observation_list title=\"observation list\""));
659 }
660
661 // Get pointer on observation library
662 GXmlElement* lib = xml.element("observation_list", 0);
663
664 // Write all observations into library
665 for (int i = 0; i < size(); ++i) {
666
667 // Initialise pointer on observation
668 GXmlElement* obs = NULL;
669
670 // Search corresponding observation
671 int n = xml.elements("observation");
672 for (int k = 0; k < n; ++k) {
673 GXmlElement* element = xml.element("observation", k);
674 if (element->attribute("name") == m_obs[i]->name() &&
675 element->attribute("id") == m_obs[i]->id() &&
676 element->attribute("instrument") == m_obs[i]->instrument()) {
677 obs = element;
678 break;
679 }
680 }
681
682 // If no observation with corresponding name, ID and instrument was
683 // found then append one now
684 if (obs == NULL) {
685 obs = lib->append("observation");
686 obs->attribute("name", m_obs[i]->name());
687 obs->attribute("id", m_obs[i]->id());
688 obs->attribute("instrument", m_obs[i]->instrument());
689 }
690
691 // Write now observation
692 m_obs[i]->write(*obs);
693
694 } // endfor: looped over all observaitons
695
696 // Return
697 return;
698}
699
700
701/***********************************************************************//**
702 * @brief Load models from XML file
703 *
704 * @param[in] filename XML filename.
705 *
706 * Loads the models from a XML file. Please refer to the GModels::read()
707 * method for more information about the expected structure of the XML
708 * file.
709 ***************************************************************************/
710void GObservations::models(const GFilename& filename)
711{
712 // Load models
713 m_models.load(filename);
714
715 // Return
716 return;
717}
718
719
720/***********************************************************************//**
721 * @brief Optimize model parameters using optimizer
722 *
723 * @param[in] opt Optimizer.
724 *
725 * Optimizes the free parameters of the models by using the optimizer
726 * that has been provided by the @p opt argument.
727 ***************************************************************************/
729{
730 // Call model setup hook for all observations
731 int nobs = size();
732 for (int i = 0; i < nobs; ++i) {
733 m_models.setup(*(m_obs[i]));
734 }
735
736 // Extract optimizer parameter container from model container
738
739 // Optimize model parameters
740 opt.optimize(m_fct, pars);
741
742 // Return
743 return;
744}
745
746
747/***********************************************************************//**
748 * @brief Computes parameter errors using optimizer
749 *
750 * @param[in] opt Optimizer.
751 *
752 * Calculates the errors of the free parameters of the models by using
753 * the optimizer that has been provided by the @p opt argument.
754 ***************************************************************************/
756{
757 // Extract optimizer parameter container from model container
759
760 // Compute model parameter errors
761 opt.errors(m_fct, pars);
762
763 // Return
764 return;
765}
766
767
768/***********************************************************************//**
769 * @brief Computes parameter errors using hessian matrix and optimizer
770 *
771 * Calculates the errors of the free model parameters as the square roots
772 * of the diagonal elements of the inverse Hessian matrix. The Hessian
773 * matrix is computed using finite differences.
774 ***************************************************************************/
776{
777 // Extract optimizer parameter container from model container
779
780 // Get number of parameters
781 int npars = pars.size();
782
783 // Compute Hessian matrix
784 GMatrixSparse hessian = m_fct.hessian(pars);
785
786 // Signal no diagonal element loading
787 bool diag_loaded = false;
788
789 // Loop over error computation (maximum 2 passes)
790 for (int i = 0; i < 2; ++i) {
791
792 // Solve: hessian * X = unit
793 try {
794 GMatrixSparse decomposition = hessian.cholesky_decompose(true);
795 GVector unit(npars);
796 for (int ipar = 0; ipar < npars; ++ipar) {
797 unit[ipar] = 1.0;
798 GVector x = decomposition.cholesky_solver(unit, true);
799 if (x[ipar] >= 0.0) {
800 pars[ipar]->factor_error(sqrt(x[ipar]));
801 }
802 else {
803 pars[ipar]->factor_error(0.0);
804 }
805 unit[ipar] = 0.0;
806 }
807 }
808 catch (GException::runtime_error &e) {
809
810 // Load diagonal if this has not yet been tried
811 if (!diag_loaded) {
812
813 // Flag errors as inaccurate
814 std::cout << "Non-Positive definite hessian matrix encountered."
815 << std::endl;
816 std::cout << "Load diagonal elements with 1e-10."
817 << " Fit errors may be inaccurate."
818 << std::endl;
819
820 // Try now with diagonal loaded matrix
821 for (int ipar = 0; ipar < npars; ++ipar) {
822 hessian(ipar,ipar) += 1.0e-10;
823 }
824
825 // Signal loading
826 diag_loaded = true;
827
828 // Try again
829 continue;
830
831 } // endif: diagonal has not yet been loaded
832
833 // ... otherwise signal an error
834 else {
835 std::cout << "Non-Positive definite hessian matrix encountered,"
836 << " even after diagonal loading." << std::endl;
837 break;
838 }
839 }
840 catch (std::exception &e) {
841 throw;
842 }
843
844 // If no error occured then break now
845 break;
846
847 } // endfor: looped over error computation
848
849 // Return
850 return;
851}
852
853
854/***********************************************************************//**
855 * @brief Evaluate function
856 *
857 * Evaluates the likelihood funtion at the actual set of parameters.
858 ***************************************************************************/
860{
861 // Extract optimizer parameter container from model container
863
864 // Compute log-likelihood
865 m_fct.eval(pars);
866
867 // Return
868 return;
869}
870
871
872/***********************************************************************//**
873 * @brief Return total number of observed events
874 *
875 * @return Total number of observed events.
876 *
877 * Returns the total number of observed events that is container in the
878 * observation container.
879 ***************************************************************************/
880double GObservations::nobserved(void) const
881{
882 // Initialise number of observed events
883 double nobserved = 0;
884
885 // Compute number of observed events
886 for (int i = 0; i < size(); ++i) {
887 nobserved += (*this)[i]->nobserved();
888 }
889
890 // Return number of observed events
891 return nobserved;
892}
893
894
895/***********************************************************************//**
896 * @brief Return number of predicted events for model with a given @p name
897 *
898 * @param[in] name Model name.
899 * @param[in] instrument Instrument name.
900 * @return Total number of predicted events for model.
901 *
902 * @exception GException::invalid_argument
903 * Model with @p name not found.
904 *
905 * Returns the total number of predicted events for the model with a given
906 * @p name.
907 ***************************************************************************/
908double GObservations::npred(const std::string& name,
909 const std::string& instrument) const
910{
911 // Initialise number of predicted events
912 double npred = 0.0;
913
914 // Throw an exception if model does not exist
915 if (!m_models.contains(name)) {
916 std::string msg = "Model \""+name+"\" not found in observation "
917 "container. Please specify a valid model name.";
919 }
920
921 // Get pointer to model
922 const GModel* model = m_models[name];
923
924 // Continue only if pointer is valid
925 if (model != NULL) {
926 npred = this->npred(*model, instrument);
927 }
928
929 // Return number of predicted events
930 return npred;
931}
932
933
934/***********************************************************************//**
935 * @brief Return number of predicted events for a given model
936 *
937 * @param[in] model Model.
938 * @param[in] instrument Instrument name.
939 * @return Total number of predicted events for model.
940 *
941 * @exception GException::invalid_argument
942 * Model with @p name not found.
943 *
944 * Returns the total number of predicted events for the model with a given
945 * @p name.
946 ***************************************************************************/
947double GObservations::npred(const GModel& model,
948 const std::string& instrument) const
949{
950 // Initialise number of predicted events
951 double npred = 0.0;
952
953 // Compute number of predicted events
954 for (int i = 0; i < size(); ++i) {
955
956 // Skip observation if it does not correspond to the specified
957 // instrument
958 if (!instrument.empty()) {
959 if ((*this)[i]->instrument() != instrument) {
960 continue;
961 }
962 }
963
964 // Add only relevant events
965 if (model.is_valid(instrument, "")) {
966 npred += (*this)[i]->npred(model);
967 }
968
969 } // endfor: looped over observations
970
971 // Return number of predicted events
972 return npred;
973}
974
975
976/***********************************************************************//**
977 * @brief Remove response cache for all models
978 *
979 * Remove response cache for all models in the container.
980 ***************************************************************************/
982{
983 // Loop over all models and remove response cache for each model
984 for (int i = 0; i < m_models.size(); ++i) {
986 }
987
988 // Return
989 return;
990}
991
992
993/***********************************************************************//**
994 * @brief Remove response cache for model
995 *
996 * @param[in] name Model name.
997 *
998 * Remove response cache for model @p name in the container.
999 ***************************************************************************/
1000void GObservations::remove_response_cache(const std::string& name)
1001{
1002 // Loop over all observations and remove response cache for the
1003 // specified source
1004 for (int i = 0; i < size(); ++i) {
1005 (*this)[i]->remove_response_cache(name);
1006 }
1007
1008 // Return
1009 return;
1010}
1011
1012
1013/***********************************************************************//**
1014 * @brief Print observation list information
1015 *
1016 * @param[in] chatter Chattiness.
1017 * @return String containing observation list information
1018 ***************************************************************************/
1019std::string GObservations::print(const GChatter& chatter) const
1020{
1021 // Initialise result string
1022 std::string result;
1023
1024 // Continue only if chatter is not silent
1025 if (chatter != SILENT) {
1026
1027 // Append header
1028 result.append("=== GObservations ===");
1029
1030 // Append information
1031 result.append("\n"+gammalib::parformat("Number of observations"));
1032 result.append(gammalib::str(size()));
1033 result.append("\n"+gammalib::parformat("Number of models"));
1034 result.append(gammalib::str(m_models.size()));
1035 result.append("\n"+gammalib::parformat("Number of observed events"));
1036 result.append(gammalib::str(nobserved()));
1037 result.append("\n"+gammalib::parformat("Number of predicted events"));
1038 result.append(gammalib::str(npred()));
1039
1040 // Get reduced chatter level for detailed information
1041 GChatter reduced_chatter = gammalib::reduce(chatter);
1042
1043 // Append observations for chatter >= EXPLICIT
1044 if (chatter >= EXPLICIT) {
1045 for (int i = 0; i < size(); ++i) {
1046 result.append("\n");
1047 result.append((*this)[i]->print(reduced_chatter));
1048 }
1049 }
1050
1051 // Append models for chatter >= VERBOSE
1052 if (chatter >= VERBOSE) {
1053 result.append("\n"+m_models.print(reduced_chatter));
1054 }
1055
1056 } // endif: chatter was not silent
1057
1058 // Return result
1059 return result;
1060}
1061
1062
1063/*==========================================================================
1064 = =
1065 = Private methods =
1066 = =
1067 ==========================================================================*/
1068
1069/***********************************************************************//**
1070 * @brief Initialise class members
1071 ***************************************************************************/
1073{
1074 // Initialise members
1075 m_obs.clear();
1076 m_models.clear();
1077 m_fct.set(this); //!< Makes sure that optimizer points to this instance
1078
1079 // Return
1080 return;
1081}
1082
1083
1084/***********************************************************************//**
1085 * @brief Copy class members
1086 *
1087 * @param[in] obs Observation container.
1088 *
1089 * Copies all class members. Deep copies are created for the observations,
1090 * which is what is needed to have a real physical copy of the members.
1091 ***************************************************************************/
1093{
1094 // Copy attributes
1095 m_models = obs.m_models;
1096 m_fct = obs.m_fct;
1097
1098 // Copy observations
1099 m_obs.clear();
1100 for (int i = 0; i < obs.m_obs.size(); ++i) {
1101 m_obs.push_back((obs.m_obs[i]->clone()));
1102 }
1103
1104 // Set likelihood function pointer to this instance. This makes sure
1105 // that the optimizer points to this instance
1106 m_fct.set(this);
1107
1108 // Return
1109 return;
1110}
1111
1112
1113/***********************************************************************//**
1114 * @brief Delete class members
1115 *
1116 * Deallocates all observations. Since container classes that hold pointers
1117 * need to handle the proper deallocation of memory, we loop here over all
1118 * pointers and make sure that we deallocate all observations in the
1119 * container.
1120 ***************************************************************************/
1122{
1123 // Free observations
1124 for (int i = 0; i < m_obs.size(); ++i) {
1125 delete m_obs[i];
1126 m_obs[i] = NULL;
1127 }
1128
1129 // Return
1130 return;
1131}
1132
1133
1134/***********************************************************************//**
1135 * @brief Return observation index by instrument and identifier
1136 *
1137 * @param[in] instrument Instrument.
1138 * @param[in] id Observation identifier.
1139 * @return Observation index (-1 of instrument and observation identifier
1140 * has not been found)
1141 *
1142 * Returns observation index based on the specified @p instrument and
1143 * observation identifier @p id. If no observation with the specified
1144 * attributes has been found, the method returns -1.
1145 ***************************************************************************/
1146int GObservations::get_index(const std::string& instrument,
1147 const std::string& id) const
1148{
1149 // Initialise index
1150 int index = -1;
1151
1152 // Search observation with specified instrument and id
1153 for (int i = 0; i < size(); ++i) {
1154 if ((m_obs[i]->instrument() == instrument) &&
1155 (m_obs[i]->id() == id)) {
1156 index = i;
1157 break;
1158 }
1159 }
1160
1161 // Return index
1162 return index;
1163}
#define G_AT
#define G_READ
#define G_APPEND
#define G_EXTEND
#define G_SET
Definition GEnergies.cpp:46
Exception handler interface definition.
Filename class interface definition.
#define G_INSERT
#define G_REMOVE
Sparse matrix class definition.
#define G_NPRED
Abstract model base class interface definition.
Observation registry class definition.
Observations container class interface definition.
Abstract optimizer abstract base class interface definition.
Gammalib tools definition.
GChatter
Definition GTypemaps.hpp:33
@ EXPLICIT
Definition GTypemaps.hpp:37
@ SILENT
Definition GTypemaps.hpp:34
@ VERBOSE
Definition GTypemaps.hpp:38
GVector sqrt(const GVector &vector)
Computes square root of vector elements.
Definition GVector.cpp:1426
XML class interface definition.
Filename class.
Definition GFilename.hpp:62
std::string url(void) const
Return Uniform Resource Locator (URL)
Sparse matrix class interface definition.
GMatrixSparse cholesky_decompose(const bool &compress=true) const
Return Cholesky decomposition.
GVector cholesky_solver(const GVector &vector, const bool &compress=true) const
Cholesky solver.
Abstract model class.
Definition GModel.hpp:100
bool is_valid(const std::string &instruments, const std::string &ids) const
Verifies if model is valid for a given instrument and identifier.
Definition GModel.cpp:585
bool contains(const std::string &name) const
Signals if model name exists.
Definition GModels.cpp:670
void setup(const GObservation &obs) const
Setup all models for a given observation.
Definition GModels.cpp:972
void clear(void)
Clear object.
Definition GModels.cpp:233
void load(const GFilename &filename)
Load models from XML file.
Definition GModels.cpp:688
GOptimizerPars pars(void) const
Return optimizer parameter container.
Definition GModels.cpp:878
int size(void) const
Return number of models in container.
Definition GModels.hpp:260
std::string print(const GChatter &chatter=NORMAL) const
Print models.
Definition GModels.cpp:992
Interface definition for the observation registry class.
GObservation * alloc(const std::string &name) const
Allocate observation of given name.
Abstract observation base class.
virtual void read(const GXmlElement &xml)=0
virtual GObservation * clone(void) const =0
Clones object.
virtual std::string instrument(void) const =0
void name(const std::string &name)
Set observation name.
void id(const std::string &id)
Set observation identifier.
void set(GObservations *obs)
Set observation container.
virtual void eval(const GOptimizerPars &pars)
Evaluate log-likelihood function.
GMatrixSparse hessian(const GOptimizerPars &pars)
Compute Hessian matrix.
Observation container class.
bool contains(const std::string &instrument, const std::string &id) const
Signals if observation exists.
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.
double nobserved(void) const
Return total number of observed events.
GModels m_models
List of models.
int size(void) const
Return number of observations in container.
GObservations * clone(void) const
Clone observations.
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.
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 parameter container class.
int size(void) const
Return number of parameters in container.
Abstract optimizer abstract base class.
virtual void optimize(GOptimizerFunction &fct, GOptimizerPars &pars)=0
virtual void errors(GOptimizerFunction &fct, GOptimizerPars &pars)=0
std::string content(void) const
Return list of names in registry.
Definition GRegistry.cpp:54
Vector class.
Definition GVector.hpp:46
XML element node class.
const GXmlAttribute * attribute(const int &index) const
Return attribute.
virtual GXmlNode * append(const GXmlNode &node)
Append XML child node.
Definition GXmlNode.cpp:287
virtual GXmlElement * element(const int &index)
Return pointer to GXMLElement child.
Definition GXmlNode.cpp:640
virtual int elements(void) const
Return number of GXMLElement children of node.
Definition GXmlNode.cpp:586
XML class.
Definition GXml.hpp:172
void save(const GFilename &filename) const
Save XML document into file.
Definition GXml.cpp:581
GXmlElement * element(const int &index)
Return pointer to child element.
Definition GXml.cpp:419
GXmlNode * append(const GXmlNode &node)
Append child node to XML document root.
Definition GXml.cpp:279
int elements(void) const
Return number of child elements in XML document root.
Definition GXml.cpp:384
std::string parformat(const std::string &s, const int &indent=0)
Convert string in parameter format.
Definition GTools.cpp:1162
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
Definition GTools.cpp:508
GChatter reduce(const GChatter &chatter)
Reduce chattiness by one level.
Definition GTypemaps.hpp:65