GammaLib 2.0.0
Loading...
Searching...
No Matches
GCTASupport.cpp
Go to the documentation of this file.
1/***************************************************************************
2 * GCTASupport.cpp - CTA support functions *
3 * ----------------------------------------------------------------------- *
4 * copyright (C) 2011-2020 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 GCTASupport.cpp
23 * @brief Implementation of support function used by CTA classes
24 * @author Juergen Knoedlseder
25 */
26
27/* __ Includes ___________________________________________________________ */
28#ifdef HAVE_CONFIG_H
29#include <config.h>
30#endif
31#include <typeinfo>
32#include <cmath>
33#include <iostream>
34#include <string>
35#include "GCTASupport.hpp"
36#include "GTools.hpp"
37#include "GMath.hpp"
38#include "GFitsHDU.hpp"
39#include "GEbounds.hpp"
40#include "GPhases.hpp"
41#include "GResponse.hpp"
42#include "GModelData.hpp"
43#include "GModelSpectral.hpp"
44#include "GModelTemporal.hpp"
45#include "GCTAObservation.hpp"
46#include "GCTAPointing.hpp"
47#include "GCTAResponse.hpp"
48#include "GCTAResponseIrf.hpp"
49#include "GCTAResponseCube.hpp"
50#include "GCTAAeff.hpp"
51#include "GCTABackground.hpp"
52#include "GCTAEventList.hpp"
53#include "GCTAEventCube.hpp"
54#include "GCTARoi.hpp"
55#include "GCTAInstDir.hpp"
61
62/* __ Method name definitions ____________________________________________ */
63#define G_READ_DS_ROI "gammalib::read_ds_roi(GFitsHDU&)"
64#define G_READ_DS_EBOUNDS "gammalib::read_ds_ebounds(GFitsHDU&)"
65#define G_READ_DS_PHASE "gammalib::read_ds_phase(GFitsHDU&)"
66#define G_READ_DS_GTI "gammalib::read_ds_gti_extname(GFitsHDU&)"
67#define G_CTA_MODEL_NROI "gammalib::cta_model_nroi(GObservation*, "\
68 "GModelData*, GEnergy&, GTime&)"
69
70/* __ Coding definitions _________________________________________________ */
71
72/* __ Debug definitions __________________________________________________ */
73
74
75/***********************************************************************//**
76 * @brief Extract ROI from data sub-space keywords
77 *
78 * @param[in] hdu FITS HDU
79 *
80 * @exception GException::invalid_value
81 * Invalid ROI data sub-space encountered
82 *
83 * Reads the ROI data sub-space keywords by searching for a DSTYPx keyword
84 * named "POS(RA,DEC)". The data sub-space information is expected to be in
85 * the format "CIRCLE(267.0208,-24.78,4.5)", where the 3 arguments are Right
86 * Ascension, Declination and radius in units of degrees. No detailed syntax
87 * checking is performed.
88 *
89 * If no ROI information has been found, an GCTARoi object with initial
90 * values will be returned.
91 ***************************************************************************/
93{
94 // Initialise ROI
95 GCTARoi roi;
96
97 // Get number of data sub-space keywords (default to 0 if keyword is
98 // not found)
99 int ndskeys = (hdu.has_card("NDSKEYS")) ? hdu.integer("NDSKEYS") : 0;
100
101 // Loop over all data selection keys
102 for (int i = 1; i <= ndskeys; ++i) {
103
104 // Set data sub-space key strings
105 std::string type_key = "DSTYP"+gammalib::str(i);
106 //std::string unit_key = "DSUNI"+gammalib::str(i);
107 std::string value_key = "DSVAL"+gammalib::str(i);
108
109 // Continue only if type_key is found and if this key is POS(RA,DEC)
110 if (hdu.has_card(type_key) && hdu.string(type_key) == "POS(RA,DEC)") {
111
112 // ...
113 //std::string unit = gammalib::toupper(hdu.string(unit_key));
114 std::string value = hdu.string(value_key);
115 std::string value_proc = gammalib::strip_chars(value, "CIRCLE(");
116 value_proc = gammalib::strip_chars(value_proc, ")");
117 std::vector<std::string> args = gammalib::split(value_proc, ",");
118 if (args.size() == 3) {
119
120 // Extract sky direction and radius
121 double ra = gammalib::todouble(args[0]);
122 double dec = gammalib::todouble(args[1]);
123 double rad = gammalib::todouble(args[2]);
124 GSkyDir skydir;
125 skydir.radec_deg(ra, dec);
126
127 // If the header contains pointing information then set a
128 // full instrument direction, including DETX and DETY
129 // coordinates.
130 if (hdu.has_card("RA_PNT") && hdu.has_card("DEC_PNT")) {
131 double ra_pnt = hdu.real("RA_PNT");
132 double dec_pnt = hdu.real("DEC_PNT");
133 GSkyDir dir_pnt;
134 dir_pnt.radec_deg(ra_pnt, dec_pnt);
135 GCTAPointing pnt(dir_pnt);
136 GCTAInstDir dir = pnt.instdir(skydir);
137 roi.centre(dir);
138 }
139
140 // ... otherwise just set the sky direction
141 else {
142 GCTAInstDir dir(skydir);
143 roi.centre(dir);
144 }
145
146 // Set RoI radius
147 roi.radius(rad);
148 }
149 else {
150 std::string msg = "Invalid acceptance cone value \""+value+
151 "\" encountered in data sub-space "
152 "key \""+value_key+"\".";
154 }
155
156 } // endif: POS(RA,DEC) type found
157
158 } // endfor: looped over data sub-space keys
159
160 // Return roi
161 return roi;
162}
163
164
165/***********************************************************************//**
166 * @brief Read energy boundary data sub-space keywords
167 *
168 * @param[in] hdu FITS HDU
169 *
170 * @exception GException::invalid_value
171 * Invalid energy data sub-space encountered
172 *
173 * Reads the energy boundary data sub-space keywords by searching for a
174 * DSTYPx keyword named "ENERGY". The data sub-space information is expected
175 * to be in the format "200:50000", where the 2 arguments are the minimum
176 * and maximum energy. The energy unit is given by the keyword DSUNIx, which
177 * supports keV, MeV, GeV and TeV (case independent). No detailed syntax
178 * checking is performed.
179 ***************************************************************************/
181{
182 // Initialise energy boundaries
183 GEbounds ebounds;
184
185 // Get number of data sub-space keywords (default to 0 if keyword is
186 // not found)
187 int ndskeys = (hdu.has_card("NDSKEYS")) ? hdu.integer("NDSKEYS") : 0;
188
189 // Loop over all data sub-space keys
190 for (int i = 1; i <= ndskeys; ++i) {
191
192 // Set data sub-space key strings
193 std::string type_key = "DSTYP"+gammalib::str(i);
194 std::string unit_key = "DSUNI"+gammalib::str(i);
195 std::string value_key = "DSVAL"+gammalib::str(i);
196
197 // Continue only if type_key is found and if this key is ENERGY
198 if (hdu.has_card(type_key) && hdu.string(type_key) == "ENERGY") {
199
200 // Extract energy boundaries
201 std::string value = hdu.string(value_key);
202 std::string unit = hdu.string(unit_key);
203 std::vector<std::string> values = gammalib::split(value, ":");
204 if (values.size() == 2) {
205 double emin = gammalib::todouble(values[0]);
206 double emax = gammalib::todouble(values[1]);
207 GEnergy e_min(emin, unit);
208 GEnergy e_max(emax, unit);
209 ebounds.append(e_min, e_max);
210 }
211 else {
212 std::string msg = "Invalid energy value \""+value+
213 "\" encountered in data selection key \""+
214 value_key+"\"";
216 }
217
218 } // endif: ENERGY type_key found
219
220 } // endfor: looped over data selection keys
221
222 // Return
223 return ebounds;
224}
225
226
227/***********************************************************************//**
228 * @brief Read phase boundary data sub-space keywords
229 *
230 * @param[in] hdu FITS HDU
231 * @return Phase intervals
232 *
233 * @exception GException::invalid_value
234 * Invalid phase data sub-space encountered
235 *
236 * Reads the phase boundary data sub-space keywords by searching for a
237 * DSTYPx keyword named "PHASE". The data sub-space information is expected
238 * to be in the format "0.1:0.3,0.5:0.7", where each subset of numbers
239 * represents the minimum and maximum phase. No detailed syntax
240 * checking is performed.
241 ***************************************************************************/
243{
244 // Initialise phase intervals
245 GPhases phases;
246
247 // Get number of data sub-space keywords (default to 0 if keyword is
248 // not found)
249 int ndskeys = (hdu.has_card("NDSKEYS")) ? hdu.integer("NDSKEYS") : 0;
250
251 // Loop over all data sub-space keys
252 for (int i = 1; i <= ndskeys; ++i) {
253
254 // Set data sub-space key strings
255 std::string type_key = "DSTYP"+gammalib::str(i);
256 //std::string unit_key = "DSUNI"+gammalib::str(i);
257 std::string value_key = "DSVAL"+gammalib::str(i);
258
259 // Continue only if type_key is found and if this key is PHASE
260 if (hdu.has_card(type_key) && hdu.string(type_key) == "PHASE") {
261
262 // Extract phase boundaries
263 std::string value = hdu.string(value_key);
264 std::vector<std::string> intervals = gammalib::split(value, ",");
265 for (int j = 0; j < intervals.size(); ++j) {
266 std::vector<std::string> values = gammalib::split(intervals[j], ":");
267 if (values.size() == 2) {
268 double phmin = gammalib::todouble(values[0]);
269 double phmax = gammalib::todouble(values[1]);
270 phases.append(phmin, phmax);
271 }
272 else {
273 std::string msg = "Invalid phase value \""+value+
274 "\" encountered in data selection key \""+
275 value_key+"\"";
277 }
278
279 } //endfor: looped over phase intervals
280
281 } // endif: PHASE type_key found
282
283 } // endfor: looped over data selection keys
284
285 // Return
286 return phases;
287}
288
289
290/***********************************************************************//**
291 * @brief Return Good Time Intervals extension name from data sub-space
292 * keywords
293 *
294 * @param[in] hdu FITS HDU
295 * @return Good Time Interval extension name
296 *
297 * @exception GException::invalid_value
298 * Invalid Good Time Intervals data sub-space encountered
299 *
300 * Returns the name of the FITS extension that contains the Good Time
301 * Intervals by screening the data sub-space keywords that are present in
302 * the FITS header. The method searches for a DSTYPx keyword named "TIME"
303 * and a corresponding DSVALx keyword named "TABLE", and the extension name
304 * is extracted from the corresponding DSREFx keyword. Note that by
305 * convention the extension name is preceeded by a colon, which is stripped
306 * by this method.
307 ***************************************************************************/
309{
310 // Initialise extension name
311 std::string extname;
312
313 // Get number of data sub-space keys (default to 0 if "NDSKEYS" keyword
314 // is not found)
315 int ndskeys = (hdu.has_card("NDSKEYS")) ? hdu.integer("NDSKEYS") : 0;
316
317 // Loop over all data sub-space keys
318 for (int i = 1; i <= ndskeys; ++i) {
319
320 // Set data sub-space key strings
321 std::string type_key = "DSTYP"+gammalib::str(i);
322 //std::string unit_key = "DSUNI"+gammalib::str(i);
323 std::string value_key = "DSVAL"+gammalib::str(i);
324 std::string ref_key = "DSREF"+gammalib::str(i);
325
326 // Skip if DSTYPi keyword does not exist, or if it is not "TIME"
327 if (!hdu.has_card(type_key) || hdu.string(type_key) != "TIME") {
328 continue;
329 }
330
331 // If DSVALi keyword does not exist or if it's value is not
332 // "TABLE" then throw an exception
333 if (!hdu.has_card(value_key)) {
334 std::string msg = "Keyword \""+value_key+"\" missing for data "
335 "sub-space selection of type "+type_key+
336 "=\"TIME\". Please correct FITS header.";
338 }
339 if (hdu.string(value_key) != "TABLE") {
340 std::string msg = "Cannot interpret keyword \""+value_key+"\" "
341 "value \""+hdu.string(value_key)+"\". Only "
342 "the value \"TABLE\" is supported.";
344 }
345
346 // If DSREFi keyword does not exist then throw an exception
347 if (!hdu.has_card(ref_key)) {
348 std::string msg = "Keyword \""+ref_key+"\" missing for data "
349 "sub-space selection of type "+type_key+
350 "=\"TIME\". Please correct FITS header.";
352 }
353
354 // Get extension name (strip any leading and trailing colons)
355 extname = gammalib::strip_whitespace(gammalib::strip_chars(hdu.string(ref_key), ":"));
356
357 } // endfor: looped over data sub-space keywords
358
359 // Return
360 return extname;
361}
362
363
364/***********************************************************************//**
365 * @brief Return extension name for GADF response table of given HDU class 4
366 *
367 * @param[in] fits FITS file.
368 * @param[in] hduclas4 HDU class 4.
369 * @return Extension name.
370 *
371 * Returns the extension name for GADF response table of given HDU class 4.
372 * If the response table is not found, an empty extension name is returned.
373 ***************************************************************************/
374std::string gammalib::gadf_hduclas4(const GFits& fits,
375 const std::string& hduclas4)
376{
377 // Initialise extension name
378 std::string extname;
379
380 // Loop over all FITS HDUs
381 for (int i = 0; i < fits.size(); ++i) {
382
383 // Get FITS HDU
384 const GFitsHDU* hdu = fits[i];
385
386 // Skip HDU if it has no "HDUCLASS", "HDUCLAS1" and "HDUCLAS4" keywords
387 if (!(hdu->has_card("HDUCLASS") &&
388 hdu->has_card("HDUCLAS1") &&
389 hdu->has_card("HDUCLAS4"))) {
390 continue;
391 }
392
393 // Skip HDU if "HDUCLASS" is not "GADF", "HDUCLAS1" is not response
394 // and "HDUCLAS4" is not hduclas4
395 if (!((gammalib::strip_whitespace(hdu->string("HDUCLASS")) == "GADF") &&
396 (gammalib::strip_whitespace(hdu->string("HDUCLAS1")) == "RESPONSE") &&
397 (gammalib::strip_whitespace(hdu->string("HDUCLAS4")) == hduclas4))) {
398 continue;
399 }
400
401 // Extract extension name and break
402 extname = hdu->extname();
403
404 } // endfor: loop over headers
405
406 // Return
407 return extname;
408}
409
410
411/***********************************************************************//**
412 * @brief Determine number of radial Romberg iterations
413 *
414 * @param[in] rho_max Maximum radial offset (radians).
415 * @param[in] resolution Requested angular resolution (radians).
416 * @param[in] iter_min Minimum number of iterations.
417 * @param[in] iter_max Maximum number of iterations.
418 * @return Number of radial Romberg iterations.
419 *
420 * Determines the number of radial Romberg iterations using the formula
421 *
422 * \f[
423 * iter = \log_2 \left( \frac{\rho_{\rm max}}{resolution} \right) + 1
424 * \f]
425 *
426 * where
427 * \f$\rho_{\rm max}\f$ is the maximum radial offset and
428 * \f$resolution\f$ is the required angular resolution.
429 *
430 * The result will be constrained to the interval [@p iter_min,@p iter_max].
431 ***************************************************************************/
432int gammalib::iter_rho(const double& rho_max,
433 const double& resolution,
434 const int& iter_min,
435 const int& iter_max)
436{
437 // Initialise number of iterations
438 int iter = iter_min;
439
440 // Continue only if radial offset and resolution are positive
441 if ((rho_max > 0.0) && (resolution > 0.0)) {
442 double arg = rho_max / resolution;
443 iter = int(std::log(arg) * gammalib::inv_ln2) + 1;
444 if (iter < iter_min) {
445 iter = iter_min;
446 }
447 else if (iter > iter_max) {
448 iter = iter_max;
449 }
450 }
451
452 // Return number of iterations
453 return iter;
454}
455
456
457/***********************************************************************//**
458 * @brief Determine number of azimuthal Romberg iterations
459 *
460 * @param[in] rho Radial offset (radians).
461 * @param[in] resolution Requested angular resolution (radians).
462 * @param[in] iter_min Minimum number of iterations.
463 * @param[in] iter_max Maximum number of iterations.
464 * @return Number of azimuthal Romberg iterations.
465 *
466 * Determines the number of azimuthal Romberg iterations using the formula
467 *
468 * \f[
469 * iter = \log_2 \left( \frac{2\pi \rho}{resolution} \right) + 1
470 * \f]
471 *
472 * where
473 * \f$\rho\f$ is the radial offset and
474 * \f$resolution\f$ is the required angular resolution.
475 *
476 * The result will be constrained to the interval [@p iter_min,@p iter_max].
477 ***************************************************************************/
478int gammalib::iter_phi(const double& rho,
479 const double& resolution,
480 const int& iter_min,
481 const int& iter_max)
482{
483 // Initialise number of iterations
484 int iter = iter_min;
485
486 // Continue only if radial offset and resolution are positive
487 if ((rho > 0.0) && (resolution > 0.0)) {
488 double arg = gammalib::twopi * rho / resolution;
489 iter = int(std::log(arg) * gammalib::inv_ln2) + 1;
490 if (iter < iter_min) {
491 iter = iter_min;
492 }
493 else if (iter > iter_max) {
494 iter = iter_max;
495 }
496 }
497
498 // Return number of iterations
499 return iter;
500}
501
502
503/***********************************************************************//**
504 * @brief Retrieve CTA observation from generic observation
505 *
506 * @param[in] origin Method asking for pointer retrieval.
507 * @param[in] obs Generic observation.
508 * @return Reference to CTA observation.
509 *
510 * @exception GException::invalid_argument
511 * Observation @p obs is not a CTA observations.
512 *
513 * Dynamically casts generic observation into a CTA observation. If the
514 * generic observation is not a CTA observation, an exception is thrown.
515 ***************************************************************************/
516const GCTAObservation& gammalib::cta_obs(const std::string& origin,
517 const GObservation& obs)
518{
519 // Get pointer on CTA observation
520 const GCTAObservation* cta = dynamic_cast<const GCTAObservation*>(&obs);
521
522 // If pointer is not valid then throw an exception
523 if (cta == NULL) {
524 std::string cls = std::string(typeid(&obs).name());
525 std::string msg = "Invalid observation type \""+cls+"\" provided with "
526 "name \""+obs.name()+"\" (ID="+obs.id()+"). "
527 "Please specify a \"GCTAObservation\" instance as "
528 "argument.";
529 throw GException::invalid_argument(origin, msg);
530 }
531
532 // Return reference
533 return (*cta);
534}
535
536
537/***********************************************************************//**
538 * @brief Retrieve CTA pointing from generic observation
539 *
540 * @param[in] origin Method asking for pointer retrieval.
541 * @param[in] obs Generic observation.
542 * @return Reference to CTA pointing.
543 *
544 * Extract CTA pointing from a CTA observation.
545 ***************************************************************************/
546const GCTAPointing& gammalib::cta_pnt(const std::string& origin,
547 const GObservation& obs)
548{
549 // Retrieve CTA observation
550 const GCTAObservation& cta = gammalib::cta_obs(origin, obs);
551
552 // Return CTA pointing
553 return (cta.pointing());
554}
555
556
557/***********************************************************************//**
558 * @brief Retrieve CTA response from generic observation
559 *
560 * @param[in] origin Method asking for pointer retrieval.
561 * @param[in] rsp Generic response.
562 * @return Pointer to CTA response.
563 *
564 * @exception GException::invalid_argument
565 * Response @p rsp is not a CTA response.
566 *
567 * Extract CTA response from a CTA observation.
568 ***************************************************************************/
569const GCTAResponse* gammalib::cta_rsp(const std::string& origin,
570 const GResponse& rsp)
571{
572 // Retrieve CTA response
573 const GCTAResponse* ptr = dynamic_cast<const GCTAResponse*>(&rsp);
574
575 // If pointer is not valid then throw an exception
576 if (ptr == NULL) {
577 std::string cls = std::string(typeid(&rsp).name());
578 std::string msg = "Invalid response type \""+cls+"\" provided. Please "
579 "specify a \"GCTAResponse\" as argument.";
580 throw GException::invalid_argument(origin, msg);
581 }
582
583 // Return CTA response
584 return (ptr);
585}
586
587
588/***********************************************************************//**
589 * @brief Retrieve CTA IRF response from generic observation
590 *
591 * @param[in] origin Method asking for pointer retrieval.
592 * @param[in] obs Generic observation.
593 * @return Reference to CTA IRF response.
594 *
595 * @exception GException::invalid_argument
596 * Observation @p obs does not contain a CTA IRF response.
597 *
598 * Extract CTA IRF response from a CTA observation.
599 ***************************************************************************/
600const GCTAResponseIrf& gammalib::cta_rsp_irf(const std::string& origin,
601 const GObservation& obs)
602{
603 // Retrieve CTA observation
604 const GCTAObservation& cta = gammalib::cta_obs(origin, obs);
605
606 // Get pointer on CTA IRF response
607 const GCTAResponseIrf* rsp = dynamic_cast<const GCTAResponseIrf*>(cta.response());
608
609 // If pointer is not valid then throw an exception
610 if (rsp == NULL) {
611 std::string cls = std::string(typeid(cta.response()).name());
612 std::string msg = "Invalid response type \""+cls+"\" provided in "
613 "CTA observation \""+obs.name()+"\" (ID="+obs.id()+"). "
614 "Please specify a CTA observation containing an IRF "
615 "response as argument.";
616 throw GException::invalid_argument(origin, msg);
617 }
618
619 // Return CTA IRF response
620 return (*rsp);
621}
622
623
624/***********************************************************************//**
625 * @brief Retrieve CTA cube response from generic observation
626 *
627 * @param[in] origin Method asking for pointer retrieval.
628 * @param[in] obs Generic observation.
629 * @return Reference to CTA cube response.
630 *
631 * @exception GException::invalid_argument
632 * Observation @p obs does not contain a CTA cube response.
633 *
634 * Extract CTA cube response from a CTA observation.
635 ***************************************************************************/
636const GCTAResponseCube& gammalib::cta_rsp_cube(const std::string& origin,
637 const GObservation& obs)
638{
639 // Retrieve CTA observation
640 const GCTAObservation& cta = gammalib::cta_obs(origin, obs);
641
642 // Get pointer on CTA response cube
643 const GCTAResponseCube* rsp = dynamic_cast<const GCTAResponseCube*>(cta.response());
644
645 // If pointer is not valid then throw an exception
646 if (rsp == NULL) {
647 std::string cls = std::string(typeid(cta.response()).name());
648 std::string msg = "Invalid response type \""+cls+"\" provided. Please "
649 "specify a CTA observation containing a response cube "
650 "as argument.";
651 throw GException::invalid_argument(origin, msg);
652 }
653
654 // Return CTA response cube
655 return (*rsp);
656}
657
658
659/***********************************************************************//**
660 * @brief Retrieve CTA effective area response from generic observation
661 *
662 * @param[in] origin Method asking for pointer retrieval.
663 * @param[in] obs Generic observation.
664 * @return Reference to CTA effective area response.
665 *
666 * @exception GException::invalid_argument
667 * Observation @p obs does not contain a CTA effective area
668 * response.
669 *
670 * Extract CTA effective area response from a CTA observation.
671 ***************************************************************************/
672const GCTAAeff& gammalib::cta_rsp_aeff(const std::string& origin,
673 const GObservation& obs)
674{
675 // Retrieve CTA IRF response
676 const GCTAResponseIrf& rsp = gammalib::cta_rsp_irf(origin, obs);
677
678 // Get pointer on CTA effective area
679 const GCTAAeff* aeff = rsp.aeff();
680 if (aeff == NULL) {
681 std::string msg = "Specified observation does not contain a valid "
682 "effective area response. Please specify a CTA "
683 "observation containing an effective area response.";
684 throw GException::invalid_argument(origin, msg);
685 }
686
687 // Return CTA effective area response
688 return (*aeff);
689}
690
691
692/***********************************************************************//**
693 * @brief Retrieve CTA background response from generic observation
694 *
695 * @param[in] origin Method asking for pointer retrieval.
696 * @param[in] obs Generic observation.
697 * @return Reference to CTA background response.
698 *
699 * @exception GException::invalid_argument
700 * Observation @p obs does not contain a CTA background response.
701 *
702 * Extract CTA background response from a CTA observation.
703 ***************************************************************************/
704const GCTABackground& gammalib::cta_rsp_bkg(const std::string& origin,
705 const GObservation& obs)
706{
707 // Retrieve CTA IRF response
708 const GCTAResponseIrf& rsp = gammalib::cta_rsp_irf(origin, obs);
709
710 // Get pointer on CTA background response
711 const GCTABackground* bkg = rsp.background();
712 if (bkg == NULL) {
713 std::string msg = "Specified observation does not contain a valid "
714 "background response. Please specify a CTA "
715 "observation containing a background response.";
716 throw GException::invalid_argument(origin, msg);
717 }
718
719 // Return CTA background response
720 return (*bkg);
721}
722
723
724/***********************************************************************//**
725 * @brief Retrieve CTA event list from generic observation
726 *
727 * @param[in] origin Method asking for pointer retrieval.
728 * @param[in] obs Generic observation.
729 * @return Reference to CTA event list.
730 *
731 * @exception GException::invalid_argument
732 * Observation @p obs does not contain a CTA event list.
733 *
734 * Extract CTA event list from a CTA observation.
735 ***************************************************************************/
736const GCTAEventList& gammalib::cta_event_list(const std::string& origin,
737 const GObservation& obs)
738{
739 // Retrieve CTA observation
740 const GCTAObservation& cta = gammalib::cta_obs(origin, obs);
741
742 // Get pointer on CTA events list
743 const GCTAEventList* events = dynamic_cast<const GCTAEventList*>(cta.events());
744
745 // If pointer is not valid then throw an exception
746 if (events == NULL) {
747 std::string msg = "Specified observation does not contain a CTA event "
748 "list. Please specify a CTA observation containing "
749 "a CTA event list as argument.";
750 throw GException::invalid_argument(origin, msg);
751 }
752
753 // Return CTA event list
754 return (*events);
755}
756
757
758/***********************************************************************//**
759 * @brief Retrieve CTA event cube from generic observation
760 *
761 * @param[in] origin Method asking for pointer retrieval.
762 * @param[in] obs Generic observation.
763 * @return Reference to CTA event cube.
764 *
765 * @exception GException::invalid_argument
766 * Observation @p obs does not contain a CTA event cube.
767 *
768 * Extract CTA event cube from a CTA observation.
769 ***************************************************************************/
770const GCTAEventCube& gammalib::cta_event_cube(const std::string& origin,
771 const GObservation& obs)
772{
773 // Retrieve CTA observation
774 const GCTAObservation& cta = gammalib::cta_obs(origin, obs);
775
776 // Get pointer on CTA event cube
777 const GCTAEventCube* events = dynamic_cast<const GCTAEventCube*>(cta.events());
778
779 // If pointer is not valid then throw an exception
780 if (events == NULL) {
781 std::string msg = "Specified observation does not contain a CTA event "
782 "cube. Please specify a CTA observation containing "
783 "a CTA event cube as argument.";
784 throw GException::invalid_argument(origin, msg);
785 }
786
787 // Return CTA event cube
788 return (*events);
789}
790
791
792/***********************************************************************//**
793 * @brief Retrieve CTA instrument direction from generic event
794 *
795 * @param[in] origin Method asking for pointer retrieval.
796 * @param[in] event Generic event.
797 * @return Reference to CTA RoI.
798 *
799 * @exception GException::invalid_argument
800 * @p event does not contain a CTA instrument direction.
801 *
802 * Extract CTA Instrument Direction from an event.
803 ***************************************************************************/
804const GCTAInstDir& gammalib::cta_dir(const std::string& origin,
805 const GEvent& event)
806{
807 // Get pointer on CTA instrument direction
808 const GCTAInstDir* dir = dynamic_cast<const GCTAInstDir*>(&(event.dir()));
809
810 // If pointer is not valid then throw an exception
811 if (dir == NULL) {
812 std::string cls = std::string(typeid(&event).name());
813 std::string msg = "Invalid event type \""+cls+"\" provided. Please "
814 "specify a \"GCTAEventAtom\" or \"GCTAEventBin\" "
815 "instance as argument.";
816 throw GException::invalid_argument(origin, msg);
817 }
818
819 // Return reference
820 return (*dir);
821}
822
823
824/***********************************************************************//**
825 * @brief Retrieve spectral component from CTA background model
826 *
827 * @param[in] model Data model.
828 * @return Pointer to spectral component.
829 *
830 * Retrieves the spectral component from a CTA background model. If the
831 * background model has no spectral component, a NULL pointer will be
832 * returned.
833 ***************************************************************************/
835{
836 // Initialise pointer on spectral component
837 GModelSpectral* spectral = NULL;
838
839 // Cast pointer to possible CTA background models
840 const GCTAModelBackground* cta = dynamic_cast<const GCTAModelBackground*>(&model);
841 const GCTAModelAeffBackground* aeff = dynamic_cast<const GCTAModelAeffBackground*>(&model);
842 const GCTAModelIrfBackground* irf = dynamic_cast<const GCTAModelIrfBackground*>(&model);
843 const GCTAModelCubeBackground* cube = dynamic_cast<const GCTAModelCubeBackground*>(&model);
844 const GCTAModelRadialAcceptance* rad = dynamic_cast<const GCTAModelRadialAcceptance*>(&model);
845
846 // Extract spectral component
847 if (cta != NULL) {
848 spectral = cta->spectral();
849 }
850 else if (aeff != NULL) {
851 spectral = aeff->spectral();
852 }
853 else if (irf != NULL) {
854 spectral = irf->spectral();
855 }
856 else if (cube != NULL) {
857 spectral = cube->spectral();
858 }
859 else if (rad != NULL) {
860 spectral = rad->spectral();
861 }
862
863 // Return pointer to spectral component
864 return spectral;
865}
866
867
868/***********************************************************************//**
869 * @brief Retrieve temporal component from CTA background model
870 *
871 * @param[in] model Data model.
872 * @return Pointer to temporal component.
873 *
874 * Retrieves the temporal component from a CTA background model. If the
875 * background model has no temporal component, a NULL pointer will be
876 * returned.
877 ***************************************************************************/
879{
880 // Initialise pointer on temporal component
881 GModelTemporal* temporal = NULL;
882
883 // Cast pointer to possible CTA background models
884 const GCTAModelBackground* cta = dynamic_cast<const GCTAModelBackground*>(&model);
885 const GCTAModelAeffBackground* aeff = dynamic_cast<const GCTAModelAeffBackground*>(&model);
886 const GCTAModelIrfBackground* irf = dynamic_cast<const GCTAModelIrfBackground*>(&model);
887 const GCTAModelCubeBackground* cube = dynamic_cast<const GCTAModelCubeBackground*>(&model);
888 const GCTAModelRadialAcceptance* rad = dynamic_cast<const GCTAModelRadialAcceptance*>(&model);
889
890 // Extract temporal component
891 if (cta != NULL) {
892 temporal = cta->temporal();
893 }
894 else if (aeff != NULL) {
895 temporal = aeff->temporal();
896 }
897 else if (irf != NULL) {
898 temporal = irf->temporal();
899 }
900 else if (cube != NULL) {
901 temporal = cube->temporal();
902 }
903 else if (rad != NULL) {
904 temporal = rad->temporal();
905 }
906
907 // Return pointer to temporal component
908 return temporal;
909}
CTA effective area base class definition.
CTA background model base class definition.
CTA event bin container class interface definition.
CTA event list class interface definition.
CTA instrument direction class interface definition.
CTA Aeff background model class definition.
Background model class interface definition.
CTA cube background model class interface definition.
CTA IRF background model class definition.
Radial acceptance model class interface definition.
CTA observation class interface definition.
CTA pointing class interface definition.
CTA cube-style response function class definition.
CTA instrument response function class definition.
CTA response abstract base class definition.
CTA region of interest class interface definition.
#define G_READ_DS_GTI
#define G_READ_DS_EBOUNDS
#define G_READ_DS_ROI
#define G_READ_DS_PHASE
Definition of support function used by CTA classes.
Energy boundaries class interface definition.
Abstract FITS extension base class definition.
Mathematical function definitions.
Abstract data model base class interface definition.
Abstract spectral model base class interface definition.
Abstract temporal model base class interface definition.
Phase intervals class interface definition.
Abstract response base class definition.
Gammalib tools definition.
Abstract base class for the CTA effective area.
Definition GCTAAeff.hpp:54
Abstract base class for the CTA background model.
CTA event bin container class.
CTA event list class.
CTA instrument direction class.
GModelTemporal * temporal(void) const
Return temporal model component.
GModelSpectral * spectral(void) const
Return spectral model component.
Background model class.
GModelSpectral * spectral(void) const
Return spectral model component.
GModelTemporal * temporal(void) const
Return temporal model component.
CTA cube background model class.
GModelSpectral * spectral(void) const
Return spectral model component.
GModelTemporal * temporal(void) const
Return temporal model component.
CTA IRF background model class.
GModelTemporal * temporal(void) const
Return temporal model component.
GModelSpectral * spectral(void) const
Return spectral model component.
Radial acceptance model class.
GModelSpectral * spectral(void) const
Return spectral model component.
GModelTemporal * temporal(void) const
Return temporal model component.
CTA observation class.
void pointing(const GCTAPointing &pointing)
Set CTA pointing.
virtual void response(const GResponse &rsp)
Set response function.
CTA pointing class.
GCTAInstDir instdir(const GSkyDir &skydir) const
Get instrument direction from sky direction.
CTA cube-style response function class.
CTA instrument response function class.
const GCTAAeff * aeff(void) const
Return pointer to effective area response.
const GCTABackground * background(void) const
Return pointer to background model.
CTA instrument response function class.
Interface for the CTA region of interest class.
Definition GCTARoi.hpp:49
const double & radius(void) const
Returns radius of region of interest in degrees.
Definition GCTARoi.hpp:125
const GCTAInstDir & centre(void) const
Returns region of interest centre.
Definition GCTARoi.hpp:111
Energy boundaries container class.
Definition GEbounds.hpp:60
void append(const GEnergy &emin, const GEnergy &emax)
Append energy interval.
Definition GEbounds.cpp:301
Class that handles energies in a unit independent way.
Definition GEnergy.hpp:48
Abstract interface for the event classes.
Definition GEvent.hpp:71
Abstract FITS extension base class.
Definition GFitsHDU.hpp:51
bool has_card(const int &cardno) const
Check existence of header card.
Definition GFitsHDU.hpp:233
const std::string & extname(void) const
Return extension name.
Definition GFitsHDU.hpp:162
double real(const std::string &keyname) const
Return card value as double precision.
Definition GFitsHDU.hpp:423
std::string string(const std::string &keyname) const
Return card value as string.
Definition GFitsHDU.hpp:410
int integer(const std::string &keyname) const
Return card value as integer.
Definition GFitsHDU.hpp:436
FITS file class.
Definition GFits.hpp:63
int size(void) const
Return number of HDUs in FITS file.
Definition GFits.hpp:237
Abstract data model class.
Abstract spectral model base class.
Abstract temporal model base class.
Abstract observation base class.
virtual GEvents * events(void)
Return events.
void name(const std::string &name)
Set observation name.
void id(const std::string &id)
Set observation identifier.
Phase Intervals class.
Definition GPhases.hpp:42
void append(const double &pmin, const double &pmax)
Append phase interval.
Definition GPhases.cpp:226
Abstract instrument response base class.
Definition GResponse.hpp:77
Sky direction class.
Definition GSkyDir.hpp:62
void radec_deg(const double &ra, const double &dec)
Set equatorial sky direction (degrees)
Definition GSkyDir.cpp:224
const GCTAObservation & cta_obs(const std::string &origin, const GObservation &obs)
Retrieve CTA observation from generic observation.
std::string read_ds_gti_extname(const GFitsHDU &hdu)
Return Good Time Intervals extension name from data sub-space keywords.
const GCTAPointing & cta_pnt(const std::string &origin, const GObservation &obs)
Retrieve CTA pointing from generic observation.
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
Definition GTools.cpp:489
double todouble(const std::string &arg)
Convert string into double precision value.
Definition GTools.cpp:926
const GCTAInstDir & cta_dir(const std::string &origin, const GEvent &event)
Retrieve CTA instrument direction from generic event.
const double inv_ln2
Definition GMath.hpp:47
const GCTABackground & cta_rsp_bkg(const std::string &origin, const GObservation &obs)
Retrieve CTA background response from generic observation.
GCTARoi read_ds_roi(const GFitsHDU &hdu)
Extract ROI from data sub-space keywords.
std::string strip_whitespace(const std::string &arg)
Strip leading and trailing whitespace from string.
Definition GTools.cpp:80
const GCTAEventList & cta_event_list(const std::string &origin, const GObservation &obs)
Retrieve CTA event list from generic observation.
int iter_phi(const double &rho, const double &resolution, const int &iter_min, const int &iter_max)
Determine number of azimuthal Romberg iterations.
double resolution(const GModelSpatial *model)
Determine resolution of spatial model.
const GCTAResponse * cta_rsp(const std::string &origin, const GResponse &rsp)
Retrieve CTA response from generic observation.
const GModelSpectral * cta_model_spectral(const GModelData &model)
Retrieve spectral component from CTA background model.
std::string gadf_hduclas4(const GFits &fits, const std::string &hduclas4)
Return extension name for GADF response table of given HDU class 4.
const double twopi
Definition GMath.hpp:36
const GCTAResponseIrf & cta_rsp_irf(const std::string &origin, const GObservation &obs)
Retrieve CTA IRF response from generic observation.
GEbounds read_ds_ebounds(const GFitsHDU &hdu)
Read energy boundary data sub-space keywords.
const GCTAAeff & cta_rsp_aeff(const std::string &origin, const GObservation &obs)
Retrieve CTA effective area response from generic observation.
int iter_rho(const double &rho_max, const double &resolution, const int &iter_min, const int &iter_max)
Determine number of radial Romberg iterations.
std::vector< std::string > split(const std::string &s, const std::string &sep)
Split string.
Definition GTools.cpp:983
GPhases read_ds_phase(const GFitsHDU &hdu)
Read phase boundary data sub-space keywords.
const GCTAEventCube & cta_event_cube(const std::string &origin, const GObservation &obs)
Retrieve CTA event cube from generic observation.
const GCTAResponseCube & cta_rsp_cube(const std::string &origin, const GObservation &obs)
Retrieve CTA cube response from generic observation.
std::string strip_chars(const std::string &arg, const std::string &chars)
Strip leading and trailing character from string.
Definition GTools.cpp:94
const GModelTemporal * cta_model_temporal(const GModelData &model)
Retrieve temporal component from CTA background model.