25 from datetime
import datetime
28 from cscripts
import calutils
36 Generates IRFs in CALDB from a ROOT offaxis performance file
45 self._init_cscript(self.__class__.__name__, ctools.__version__, argv)
54 Get parameters from parfile
57 self[
'infile'].filename()
58 self[
'outdir'].string()
61 self[
'analysis'].string()
63 self[
'azimuth'].real()
64 self[
'version'].string()
65 self[
'psftype'].string()
66 self[
'split'].boolean()
67 self[
'norm1d'].boolean()
68 self[
'rebin'].boolean()
69 self[
'eascale'].real()
70 self[
'bgdscale'].real()
71 self[
'bgdoversample'].integer()
72 self[
'bgdinfill'].boolean()
75 self._log_parameters(gammalib.TERSE)
82 Initialise dictionary containing the IRF metadata
87 IRF metadata dictionary
93 cal_name =
'NAME('+self[
'id'].string()+
')'
94 cal_version =
'VERSION('+self[
'version'].string()+
')'
95 cal_cut =
'CLASS(BEST)'
96 cal_analysis =
'ANALYSIS('+self[
'analysis'].string()+
')'
97 cal_zenith =
'ZENITH(%.3f)deg' % self[
'zenith'].real()
98 cal_azimuth =
'AZIMUTH(%.3f)deg' % self[
'azimuth'].real()
99 cal_bounds = [cal_name, cal_version, cal_cut, cal_analysis, \
100 cal_zenith, cal_azimuth]
103 irf[
'CAL_TEL'] =
'CTA'
104 irf[
'CAL_INST'] = self[
'inst'].string().upper()
105 irf[
'CAL_OBSID'] = self[
'id'].string()
106 irf[
'CAL_DET'] =
'NONE'
107 irf[
'CAL_FLT'] =
'NONE'
108 irf[
'CAL_CLASS'] =
'BCF'
109 irf[
'CAL_TYPE'] =
'DATA'
111 irf[
'CAL_DATE'] =
'14/01/30'
112 irf[
'VAL_DATE'] =
'2014-01-30'
113 irf[
'VAL_TIME'] =
'00:00:00'
114 irf[
'REF_TIME'] = 51544.0
115 irf[
'EA_NAME'] =
'EFF_AREA'
116 irf[
'EA_HDUCLAS3'] =
'FULL-ENCLOSURE'
117 irf[
'EA_HDUCLAS4'] =
'2D'
118 irf[
'EA_DOC'] =
'CAL/GEN/92-019'
119 irf[
'EA_BOUNDS'] = copy.deepcopy(cal_bounds)
120 irf[
'EA_DESC'] =
'CTA effective area'
121 irf[
'PSF_NAME'] =
'RPSF'
122 irf[
'PSF_HDUCLAS3'] =
'FULL-ENCLOSURE'
123 irf[
'PSF_HDUCLAS4'] =
'3GAUSS'
124 irf[
'PSF_DOC'] =
'CAL/GEN/92-020'
125 irf[
'PSF_BOUNDS'] = copy.deepcopy(cal_bounds)
126 irf[
'PSF_DESC'] =
'CTA point spread function'
127 irf[
'EDISP_NAME'] =
'EDISP'
128 irf[
'EDISP_HDUCLAS3'] =
'FULL-ENCLOSURE'
129 irf[
'EDISP_HDUCLAS4'] =
'2D'
130 irf[
'EDISP_DOC'] =
'???'
131 irf[
'EDISP_BOUNDS'] = copy.deepcopy(cal_bounds)
132 irf[
'EDISP_DESC'] =
'CTA energy dispersion'
133 irf[
'BKG_NAME'] =
'BKG'
134 irf[
'BKG_HDUCLAS3'] =
'FULL-ENCLOSURE'
135 irf[
'BKG_HDUCLAS4'] =
'3D'
136 irf[
'BKG_DOC'] =
'???'
137 irf[
'BKG_BOUNDS'] = copy.deepcopy(cal_bounds)
138 irf[
'BKG_DESC'] =
'CTA background'
141 if self[
'psftype'].string() ==
'King':
142 irf[
'PSF_HDUCLAS4'] =
'KING'
144 irf[
'PSF_HDUCLAS4'] =
'3GAUSS'
151 Generate CALDB directory structure for one observation identifier
153 The structure is given by
155 data/<tel>/<inst>/bcf/<obsid>
157 where <tel> is "cta" and <inst> is the instrument specified in the
158 CALDB constructor (the instrument may be used for different array
166 IRF metadata dictionary
171 Directory structure information
178 base_dir +=
'/'+irf[
'CAL_TEL'].lower()
179 base_dir +=
'/'+irf[
'CAL_INST'].lower()
182 obs_dir = base_dir+
'/bcf/'+irf[
'CAL_OBSID']
183 ds[
'EA_DIR'] = obs_dir
184 ds[
'PSF_DIR'] = obs_dir
185 ds[
'EDISP_DIR'] = obs_dir
186 ds[
'BKG_DIR'] = obs_dir
190 outdir = self[
'outdir'].string()
192 ds[
'BASE_PATH'] = outdir+
'/'+base_dir
193 ds[
'EA_PATH'] = outdir+
'/'+ds[
'EA_DIR']
194 ds[
'PSF_PATH'] = outdir+
'/'+ds[
'PSF_DIR']
195 ds[
'EDISP_PATH'] = outdir+
'/'+ds[
'EDISP_DIR']
196 ds[
'BKG_PATH'] = outdir+
'/'+ds[
'BKG_DIR']
198 ds[
'BASE_PATH'] = base_dir
199 ds[
'EA_PATH'] = ds[
'EA_DIR']
200 ds[
'PSF_PATH'] = ds[
'PSF_DIR']
201 ds[
'EDISP_PATH'] = ds[
'EDISP_DIR']
202 ds[
'BKG_PATH'] = ds[
'BKG_DIR']
206 if self[
'split'].boolean():
207 ds[
'EA_FILE'] =
'ea_'+version+
'.fits'
208 ds[
'PSF_FILE'] =
'psf_'+version+
'.fits'
209 ds[
'EDISP_FILE'] =
'edisp_'+version+
'.fits'
210 ds[
'BKG_FILE'] =
'bgd_'+version+
'.fits'
212 ds[
'EA_FILE'] =
'irf_'+version+
'.fits'
213 ds[
'PSF_FILE'] =
'irf_'+version+
'.fits'
214 ds[
'EDISP_FILE'] =
'irf_'+version+
'.fits'
215 ds[
'BKG_FILE'] =
'irf_'+version+
'.fits'
218 if not os.path.isdir(ds[
'EA_PATH']):
219 os.makedirs(ds[
'EA_PATH'])
220 if not os.path.isdir(ds[
'PSF_PATH']):
221 os.makedirs(ds[
'PSF_PATH'])
222 if not os.path.isdir(ds[
'EDISP_PATH']):
223 os.makedirs(ds[
'EDISP_PATH'])
224 if not os.path.isdir(ds[
'BKG_PATH']):
225 os.makedirs(ds[
'BKG_PATH'])
232 Open existing or create new calibration
234 The actual version will put all calibrations in the same file, although
235 each part of the response function will have its own logical name. We
236 can thus easily modify the script to put each calibration information
242 IRF metadata dictionary
244 Directory structure dictionary
247 ds[
'CIF'] = gammalib.GFits(ds[
'BASE_PATH']+
'/caldb.indx',
True)
250 if not ds[
'CIF'].contains(
'CIF'):
251 ds[
'CIF'].append(calutils.create_cif_table())
252 ds[
'HDU_CIF'] = ds[
'CIF'].table(
'CIF')
255 ea_filename = ds[
'EA_PATH']+
'/'+ds[
'EA_FILE']
256 psf_filename = ds[
'PSF_PATH']+
'/'+ds[
'PSF_FILE']
257 edisp_filename = ds[
'EDISP_PATH']+
'/'+ds[
'EDISP_FILE']
258 bgd_filename = ds[
'BKG_PATH']+
'/'+ds[
'BKG_FILE']
261 if self[
'split'].boolean():
262 ds[
'EA_FITS'] = gammalib.GFits(ea_filename,
True)
263 ds[
'PSF_FITS'] = gammalib.GFits(psf_filename,
True)
264 ds[
'EDISP_FITS'] = gammalib.GFits(edisp_filename,
True)
265 ds[
'BKG_FITS'] = gammalib.GFits(bgd_filename,
True)
267 ds[
'EA_FITS'] = gammalib.GFits(ea_filename,
True)
268 ds[
'PSF_FITS'] = ds[
'EA_FITS']
269 ds[
'EDISP_FITS'] = ds[
'EA_FITS']
270 ds[
'BKG_FITS'] = ds[
'EA_FITS']
273 ds[
'HDU_EA'] = self.
_open_hdu(ds[
'EA_FITS'],
"EFFECTIVE AREA",
274 irf[
'EA_NAME'], irf[
'EA_HDUCLAS3'],
275 irf[
'EA_HDUCLAS4'], irf[
'EA_DOC'],
277 ds[
'HDU_PSF'] = self.
_open_hdu(ds[
'PSF_FITS'],
"POINT SPREAD FUNCTION",
278 irf[
'PSF_NAME'], irf[
'PSF_HDUCLAS3'],
279 irf[
'PSF_HDUCLAS4'], irf[
'PSF_DOC'],
281 ds[
'HDU_EDISP'] = self.
_open_hdu(ds[
'EDISP_FITS'],
"ENERGY DISPERSION",
282 irf[
'EDISP_NAME'], irf[
'EDISP_HDUCLAS3'],
283 irf[
'EDISP_HDUCLAS4'], irf[
'EDISP_DOC'],
285 ds[
'HDU_BKG'] = self.
_open_hdu(ds[
'BKG_FITS'],
"BACKGROUND",
286 irf[
'BKG_NAME'], irf[
'BKG_HDUCLAS3'],
287 irf[
'BKG_HDUCLAS4'], irf[
'BKG_DOC'],
295 Close calibration FITS files
300 IRF metadata dictionary
302 Directory structure dictionary
316 ds[
'EA_FITS'].save(
True)
317 ds[
'EA_FITS'].close()
318 if self[
'split'].boolean():
319 ds[
'PSF_FITS'].save(
True)
320 ds[
'EDISP_FITS'].save(
True)
321 ds[
'BKG_FITS'].save(
True)
322 ds[
'PSF_FITS'].close()
323 ds[
'EDISP_FITS'].close()
324 ds[
'BKG_FITS'].close()
329 def _open_hdu(self, fits, extname, name, hduclas3, hduclas4, doc, irf):
333 Opens a FITS binary table with given "extname". If HDU does not exist
334 in the FITS file it will be created and appended to the FITS file.
338 fits : `~gammalib.GFits`
351 IRF metadata dictionary
355 table : `~gammalib.GFitsBinTable`
359 if not fits.contains(extname):
362 table = gammalib.GFitsBinTable()
365 table.extname(extname)
369 [
'RESPONSE', name, hduclas3, hduclas4],
376 return fits.table(extname)
380 Set standard OGIP keywords for extension
384 hdu : `~gammalig.GFitsHDU`
387 Documentation reference string
388 hduclas : list of str
389 List of HDUCLAS fields
394 utc = datetime.utcnow().isoformat()[:19]
397 hdu.card(
'ORIGIN',
'IRAP',
'Name of organization making this file')
398 hdu.card(
'DATE', utc,
'File creation date (YYYY-MM-DDThh:mm:ss UTC)')
399 hdu.card(
'TELESCOP', irf[
'CAL_TEL'],
'Name of telescope')
400 hdu.card(
'INSTRUME', irf[
'CAL_INST'],
'Name of instrument')
401 hdu.card(
'DETNAM', irf[
'CAL_DET'],
'Name of detector')
402 hdu.card(
'HDUCLASS',
'OGIP',
'HDU class')
403 hdu.card(
'HDUDOC', hdudoc,
'HDU documentation')
404 for i, item
in enumerate(hduclas):
405 key =
'HDUCLAS%d' % (i+1)
406 hdu.card(key, item,
'HDU class')
407 hdu.card(
'HDUVERS',
'0.2',
'HDU version')
417 Add information to CIF extension
422 IRF metadata dictionary
424 Directory structure dictionary
427 names = [
'EA',
'PSF',
'EDISP',
'BKG']
430 row = ds[
'HDU_CIF'].nrows()
433 ds[
'HDU_CIF'].append_rows(len(names))
439 key_dir =
'%s_DIR' % name
440 key_file =
'%s_FILE' % name
441 key_name =
'%s_NAME' % name
442 key_desc =
'%s_DESC' % name
443 key_bounds =
'%s_BOUNDS' % name
446 ds[
'HDU_CIF'][
'TELESCOP'][row] = irf[
'CAL_TEL']
447 ds[
'HDU_CIF'][
'INSTRUME'][row] = irf[
'CAL_INST']
448 ds[
'HDU_CIF'][
'DETNAM'][row] = irf[
'CAL_DET']
449 ds[
'HDU_CIF'][
'FILTER'][row] = irf[
'CAL_FLT']
450 ds[
'HDU_CIF'][
'CAL_DEV'][row] =
'ONLINE'
451 ds[
'HDU_CIF'][
'CAL_CLAS'][row] = irf[
'CAL_CLASS']
452 ds[
'HDU_CIF'][
'CAL_DTYP'][row] = irf[
'CAL_TYPE']
453 ds[
'HDU_CIF'][
'CAL_VSD'][row] = irf[
'VAL_DATE']
454 ds[
'HDU_CIF'][
'CAL_VST'][row] = irf[
'VAL_TIME']
455 ds[
'HDU_CIF'][
'REF_TIME'][row] = irf[
'REF_TIME']
456 ds[
'HDU_CIF'][
'CAL_QUAL'][row] = irf[
'CAL_QUAL']
457 ds[
'HDU_CIF'][
'CAL_DATE'][row] = irf[
'CAL_DATE']
460 ds[
'HDU_CIF'][
'CAL_DIR'][row] = ds[key_dir]
461 ds[
'HDU_CIF'][
'CAL_FILE'][row] = ds[key_file]
462 ds[
'HDU_CIF'][
'CAL_CNAM'][row] = irf[key_name]
463 ds[
'HDU_CIF'][
'CAL_DESC'][row] = irf[key_desc]
464 ds[
'HDU_CIF'][
'CAL_XNO'][row] = 1
466 if i >= len(irf[key_bounds]):
467 ds[
'HDU_CIF'][
'CAL_CBD'][row,i] =
'NONE'
469 ds[
'HDU_CIF'][
'CAL_CBD'][row,i] = irf[key_bounds][i]
479 Set standard CIF keywords for extension
483 hdu : `~gammalib.GFitsHDU`
488 Calibration boundaries
490 Calibration description
492 IRF metadata dictionary
495 hdu.card(
'CSYSNAME',
'XMA_POL',
'')
496 hdu.card(
'CCLS0001', irf[
'CAL_CLASS'],
'Calibration class')
497 hdu.card(
'CDTP0001', irf[
'CAL_TYPE'],
'Calibration type')
498 hdu.card(
'CCNM0001', name,
'Calibration name')
502 keyname =
'CBD%d0001' % (i+1)
507 hdu.card(keyname, value,
'Calibration boundary')
510 hdu.card(
'CVSD0001', irf[
'VAL_DATE'],
511 'Calibration validity start date (UTC)')
512 hdu.card(
'CVST0001', irf[
'VAL_TIME'],
513 'Calibration validity start time (UTC)')
514 hdu.card(
'CDEC0001', desc,
'Calibration description')
523 Translate ROOT to CALDB information
528 IRF metadata dictionary
530 Directory structure dictionary
533 from ROOT
import TFile
536 tfile = TFile(self[
'infile'].filename().url())
555 Translate ROOT to CALDB effective area extension
557 The following ROOT histograms are used:
558 - EffectiveAreaEtrueNoTheta2cut_offaxis -> EFFAREA
559 - EffectiveAreaEtrueNoTheta2cut (for normalization)
563 tfile : `~ROOT.TFile`
566 IRF metadata dictionary
568 Directory structure dictionary
571 self._log_header1(gammalib.TERSE,
'Generate effective area extension')
574 etrue = tfile.Get(
'EffectiveAreaEtrueNoTheta2cut_offaxis')
579 if self[
'norm1d'].boolean():
580 etrue_1D = tfile.Get(
'EffectiveAreaEtrueNoTheta2cut')
585 if self[
'rebin'].boolean():
587 neng = etrue.GetXaxis().GetNbins()
588 noffset = etrue.GetYaxis().GetNbins()
589 for ioff
in range(noffset):
590 for ieng
in range(neng):
591 value = etrue.GetBinContent(ieng+1,ioff+1) / 10.0
592 etrue.SetBinContent(ieng+1,ioff+1,value)
596 eascale = self[
'eascale'].real()
600 neng = etrue.GetXaxis().GetNbins()
601 noffset = etrue.GetYaxis().GetNbins()
602 for ioff
in range(noffset):
603 for ieng
in range(neng):
604 value = etrue.GetBinContent(ieng+1,ioff+1) * eascale
605 etrue.SetBinContent(ieng+1,ioff+1,value)
609 irf[
'EA_BOUNDS'], irf[
'EA_DESC'], irf)
612 if self[
'emin'].is_valid():
613 emin = self[
'emin'].real()
614 ds[
'HDU_EA'].card(
'LO_THRES', emin,
'[TeV] Low energy threshold')
615 if self[
'emax'].is_valid():
616 emin = self[
'emax'].real()
617 ds[
'HDU_EA'].card(
'HI_THRES', emin,
'[TeV] High energy threshold')
620 self.
_make_2D(etrue, ds[
'HDU_EA'],
'EFFAREA',
'm**2')
627 Translate ROOT to CALDB point spread function extension
631 tfile : `~ROOT.TFile`
634 IRF metadata dictionary
636 Directory structure dictionary
639 if self[
'psftype'].string() ==
'King':
651 Generates 68% and 80% containment radii histograms from
652 AngularPSF2DEtrue_offaxis histogram
656 tfile : `~ROOT.TFile`
661 r68, r80 : Tuple of `~ROOT.THF2`
662 68% and 80% containment radius ROOT histograms
665 psf = tfile.Get(
'AngularPSF2DEtrue_offaxis')
668 energies = psf.GetXaxis()
669 dists = psf.GetYaxis()
670 offsets = psf.GetZaxis()
673 neng = energies.GetNbins()
674 ndist = dists.GetNbins()
675 noffset = offsets.GetNbins()
678 from ROOT
import TH2F
681 r68 = TH2F(
'AngRes68Etrue_offaxis',
682 'Angular resolution (68% containment)',
683 neng, energies.GetBinLowEdge(1), energies.GetBinUpEdge(neng),
684 noffset, offsets.GetBinLowEdge(1), offsets.GetBinUpEdge(noffset))
685 r80 = TH2F(
'AngRes80Etrue_offaxis',
686 'Angular resolution (80% containment)',
687 neng, energies.GetBinLowEdge(1), energies.GetBinUpEdge(neng),
688 noffset, offsets.GetBinLowEdge(1), offsets.GetBinUpEdge(noffset))
691 for ieng
in range(neng):
692 for ioff
in range(noffset):
696 for idist
in range(ndist):
697 pdf_integral += psf.GetBinContent(ieng+1,idist+1,ioff+1)
702 for idist
in range(ndist):
703 pdf_68 += psf.GetBinContent(ieng+1,idist+1,ioff+1)
704 if pdf_68 > 0.68*pdf_integral:
705 r68_angle = dists.GetBinLowEdge(idist+1)
711 for idist
in range(ndist):
712 pdf_80 += psf.GetBinContent(ieng+1,idist+1,ioff+1)
713 if pdf_80 > 0.80*pdf_integral:
714 r80_angle = dists.GetBinLowEdge(idist+1)
718 r68.SetBinContent(ieng+1,ioff+1,r68_angle)
719 r80.SetBinContent(ieng+1,ioff+1,r80_angle)
726 Translate ROOT to CALDB point spread function extension
728 The following ROOT histograms are used:
729 - 1/(2*pi*SIGMA_1) -> SCALE
730 - AngRes_offaxis -> SIGMA_1 (scaling: 1/0.8)
738 tfile : `~ROOT.TFile`
741 IRF metadata dictionary
743 Directory structure dictionary
746 self._log_header1(gammalib.TERSE,
747 'Generate Gaussian point spread function extension')
753 neng = r68.GetXaxis().GetNbins()
754 noffset = r68.GetYaxis().GetNbins()
757 r68_to_sigma = 0.6624305
758 for ioff
in range(noffset):
759 for ieng
in range(neng):
760 sigma = r68.GetBinContent(ieng+1,ioff+1) * r68_to_sigma
761 r68.SetBinContent(ieng+1,ioff+1,sigma)
765 for ioff
in range(noffset):
766 for ieng
in range(neng):
767 integral = 2.0 * math.pi * r68.GetBinContent(ieng+1,ioff+1)
769 value = 1.0 / integral
772 scale.SetBinContent(ieng+1,ioff+1,value)
776 for ioff
in range(noffset):
777 for ieng
in range(neng):
778 zero.SetBinContent(ieng+1,ioff+1,0.0)
781 irf[
'PSF_BOUNDS'].append(
'PSF(GAUSS)')
785 irf[
'PSF_BOUNDS'], irf[
'PSF_DESC'], irf)
788 self.
_make_2D(scale, ds[
'HDU_PSF'],
'SCALE',
'')
791 self.
_make_2D(r68, ds[
'HDU_PSF'],
'SIGMA_1',
'deg')
794 self.
_make_2D(zero, ds[
'HDU_PSF'],
'AMPL_2',
'')
797 self.
_make_2D(zero, ds[
'HDU_PSF'],
'SIGMA_2',
'deg')
800 self.
_make_2D(zero, ds[
'HDU_PSF'],
'AMPL_3',
'')
803 self.
_make_2D(zero, ds[
'HDU_PSF'],
'SIGMA_3',
'deg')
810 Translate ROOT to CALDB point spread function extension
812 The following ROOT histograms are used:
818 tfile : `~ROOT.TFile`
821 IRF metadata dictionary
823 Directory structure dictionary
826 self._log_header1(gammalib.TERSE,
827 'Generate King point spread function extension')
833 neng = r68.GetXaxis().GetNbins()
834 noffset = r68.GetYaxis().GetNbins()
837 gamma2D = r68.Clone()
838 sigma2D = r68.Clone()
841 for ioff
in range(noffset):
848 for ieng
in range(neng):
851 r_68 = r68.GetBinContent(ieng+1,ioff+1)
852 r_80 = r80.GetBinContent(ieng+1,ioff+1)
859 if r_68 > 0
and r_80 > 0:
864 c = r_68*r_68/(r_80*r_80)
870 f1 = (math.pow(a,x1) - 1.0)/(math.pow(b,x1) - 1.0) - c
871 f2 = (math.pow(a,x2) - 1.0)/(math.pow(b,x2) - 1.0) - c
873 x = x1 - f1 * (x1-x2)/(f1-f2)
874 f = (math.pow(a,x) - 1.0)/(math.pow(b,x) - 1.0) - c
890 denominator = 2.0 * gamma * (math.pow(a, x) - 1.0)
891 if denominator > 0.0:
892 sigma = r_68 * math.sqrt(1.0/denominator)
894 denominator = 2.0 * gamma * (math.pow(b, x) - 1.0)
895 if denominator > 0.0:
896 sigma = r_80 * math.sqrt(1.0/denominator)
904 if gamma == 0.0
and sigma == 0.0:
907 status =
' (use last)'
912 self._log_value(gammalib.EXPLICIT,
913 'ieng=%d ioff=%d%s' % (ieng,ioff,status),
914 'r68=%f r80=%f gamma=%f sigma=%f' %
915 (r_68, r_80, gamma, sigma))
922 gamma2D.SetBinContent(ieng+1,ioff+1,gamma)
923 sigma2D.SetBinContent(ieng+1,ioff+1,sigma)
926 irf[
'PSF_BOUNDS'].append(
'PSF(KING)')
930 irf[
'PSF_BOUNDS'], irf[
'PSF_DESC'], irf)
933 self.
_make_2D(gamma2D, ds[
'HDU_PSF'],
'GAMMA',
'')
936 self.
_make_2D(sigma2D, ds[
'HDU_PSF'],
'SIGMA',
'deg')
943 Generates migration matrix histogram from AngularPSF2DEtrue_offaxis
948 tfile : `~ROOT.TFile`
957 matrix = tfile.Get(
'MigMatrixNoTheta2cut_offaxis')
960 ereco = matrix.GetXaxis()
961 etrue = matrix.GetYaxis()
962 offsets = matrix.GetZaxis()
965 nereco = ereco.GetNbins()
966 netrue = etrue.GetNbins()
967 noffset = offsets.GetNbins()
973 from ROOT
import TH3F
976 migra = TH3F(
'EestOverEtrue_offaxis',
978 netrue, etrue.GetBinLowEdge(1), etrue.GetBinUpEdge(netrue),
980 noffset, offsets.GetBinLowEdge(1), offsets.GetBinUpEdge(noffset))
983 for ietrue
in range(netrue):
986 e_true = math.pow(10.0, migra.GetXaxis().GetBinCenter(ietrue+1))
989 for imigra
in range(nmigra):
992 migration = migra.GetYaxis().GetBinCenter(imigra+1)
995 e_reco = e_true * migration
1001 ireco = ereco.FindFixBin(math.log10(e_reco))
1004 if ireco > 0
and ireco <= nereco:
1007 for ioff
in range(noffset):
1010 value = matrix.GetBinContent(ireco,ietrue+1,ioff+1)
1013 migra.SetBinContent(ietrue+1,imigra+1,ioff+1,value)
1020 Translate ROOT to CALDB energy dispersion extension
1022 The following ROOT histograms are used:
1023 - EestOverEtrue_offaxis -> MATRIX
1027 tfile : `~ROOT.TFile`
1030 IRF metadata dictionary
1032 Directory structure dictionary
1035 self._log_header1(gammalib.TERSE,
'Generate energy dispersion extension')
1042 irf[
'EDISP_BOUNDS'], irf[
'EDISP_DESC'], irf)
1052 Translate ROOT to CALDB background extension.
1054 The following ROOT histograms are used:
1055 - BGRatePerSqDeg_offaxis -> BKG
1059 tfile : `~ROOT.TFile`
1062 IRF metadata dictionary
1064 Directory structure dictionary
1067 self._log_header1(gammalib.TERSE,
'Generate 3D background extension')
1070 array = tfile.Get(
'BGRatePerSqDeg_offaxis')
1076 if self[
'bgdinfill'].boolean():
1082 if self[
'norm1d'].boolean():
1083 array_1D = tfile.Get(
'BGRatePerSqDeg')
1088 irf[
'BKG_BOUNDS'], irf[
'BKG_DESC'], irf)
1091 self.
_make_3D(array, ds[
'HDU_BKG'],
'BKG',
'1/(MeV s sr)')
1098 Append column of ROOT axis values to HDU
1102 hdu : `~gammalib.GFitsHDU`
1108 axis : `~ROOT.TAxis`
1110 log : bool, optional
1114 if not hdu.contains(name):
1117 self._log_header3(gammalib.TERSE,
'Append axis column "%s"' % name)
1120 nbins = axis.GetNbins()
1123 hdu.append(gammalib.GFitsTableFloatCol(name, 1, nbins))
1124 hdu[name].unit(unit)
1127 if name[-2:] ==
'LO':
1133 for i
in range(nbins):
1135 value = axis.GetBinLowEdge(i+1)
1137 value = axis.GetBinUpEdge(i+1)
1139 value = pow(10.0, value)
1140 hdu[name][0,i] = value
1143 self._log_value(gammalib.NORMAL,
'Number of axis bins', nbins)
1144 self._log_value(gammalib.NORMAL,
'Unit', unit)
1151 Append column of values to HDU
1155 hdu : `~gammalib.GFitsHDU`
1161 values : list of float
1165 if not hdu.contains(name):
1168 self._log_header3(gammalib.TERSE,
'Append value column "%s"' % name)
1174 hdu.append(gammalib.GFitsTableFloatCol(name, 1, nbins))
1175 hdu[name].unit(unit)
1178 for i
in range(nbins):
1179 hdu[name][0,i] = values[i]
1182 self._log_value(gammalib.NORMAL,
'Number of values', nbins)
1183 self._log_value(gammalib.NORMAL,
'Unit', unit)
1190 Renormalise 2D histogram (energy,offset) on 1D histogram (energy)
1192 This method makes sure that a 2D histogram has the same onaxis values
1193 as the corresponding 1D histogram.
1197 hist2D : `~ROOT.TH2F`
1199 hist1D : `~ROOT.TH1F`
1203 neng = hist2D.GetXaxis().GetNbins()
1204 noffset = hist2D.GetYaxis().GetNbins()
1208 if neng != hist1D.GetXaxis().GetNbins():
1211 self._log_value(gammalib.TERSE,
'On-axis normalisation',
1212 'Impossible since energy binning of 1D histogram '
1213 '"%s" does not match that of 2D histogram "%s".' %
1214 (hist1D.GetName(), hist2D.GetName()))
1217 self._log_value(gammalib.TERSE,
'On-axis normalisation',
1218 'On-axis values of 1D histogram "%s" imposed on '
1219 '2D histogram "%s".' %
1220 (hist1D.GetName(), hist2D.GetName()))
1223 hist2D_copy = hist2D.Clone()
1227 for ieng
in range(neng):
1228 onaxis = hist1D.GetBinContent(ieng+1)
1230 for ioff
in range(noffset):
1231 value = hist2D.GetBinContent(ieng+1,ioff+1) / onaxis
1232 hist2D.SetBinContent(ieng+1,ioff+1,value)
1235 for ieng
in range(neng):
1236 for ioff
in range(noffset):
1238 value = (2.0 * hist2D.GetBinContent(ieng+1,ioff+1) +
1239 1.0 * hist2D.GetBinContent(ieng+2,ioff+1)) / 3.0
1240 elif ieng == neng-1:
1241 value = (2.0 * hist2D.GetBinContent(ieng+1,ioff+1) +
1242 1.0 * hist2D.GetBinContent(ieng+0,ioff+1)) / 3.0
1244 value = (2.0 * hist2D.GetBinContent(ieng+1,ioff+1) +
1245 1.0 * hist2D.GetBinContent(ieng+0,ioff+1) +
1246 1.0 * hist2D.GetBinContent(ieng+2,ioff+1)) / 4.0
1247 hist2D_copy.SetBinContent(ieng+1,ioff+1,value)
1250 for ieng
in range(neng):
1251 onaxis = hist2D_copy.GetBinContent(ieng+1,1)
1253 for ioff
in range(noffset):
1254 value = hist2D_copy.GetBinContent(ieng+1,ioff+1) / onaxis
1255 hist2D_copy.SetBinContent(ieng+1,ioff+1,value)
1258 for ieng
in range(neng):
1259 onaxis = hist1D.GetBinContent(ieng+1)
1260 for ioff
in range(noffset):
1261 value = hist2D_copy.GetBinContent(ieng+1,ioff+1) * onaxis
1262 hist2D.SetBinContent(ieng+1,ioff+1,value)
1269 Replace background rate values by power law values
1273 hist2D : `~ROOT.TH2F`
1279 energies = hist2D.GetXaxis()
1280 neng = energies.GetNbins()
1281 offsets = hist2D.GetYaxis()
1282 noffset = offsets.GetNbins()
1285 if ethres < math.pow(10.0, energies.GetBinUpEdge(neng)):
1289 for ieng
in range(neng):
1290 energy = 0.5 * (math.pow(10.0, energies.GetBinLowEdge(ieng+1)) +
1291 math.pow(10.0, energies.GetBinUpEdge(ieng+1)))
1295 for ioff
in range(noffset):
1298 offset = 0.5 * (offsets.GetBinLowEdge(ioff+1) +
1299 offsets.GetBinUpEdge(ioff+1))
1306 for ieng
in range(neng):
1308 rate = hist2D.GetBinContent(ieng+1,ioff+1)
1321 if coeff[
'm'] != 0.0
and coeff[
't'] != 0.0:
1322 for ieng
in range(neng):
1326 hist2D.SetBinContent(ieng+1,ioff+1,value)
1329 self._log_value(gammalib.NORMAL,
1330 'Power law replacement',
1331 '%e cts/(s deg**2) (E=%8.4f TeV (%d) '
1332 'Off=%.2f deg (%d)' %
1333 (value,engs[ieng],ieng,offset,ioff))
1340 Determine the energy range over which the power law is valid
1344 energies : list of float
1346 rates : list of float
1347 Background rate values
1352 for ieng, energy
in enumerate(energies):
1353 if rates[ieng] > max_rate:
1354 max_rate = rates[ieng]
1359 for ieng
in range(len(energies)-1,-1,-1):
1360 if rates[ieng] > 0.0:
1361 max_energy = energies[ieng]
1365 return (min_energy, max_energy)
1369 Determines the power law coefficients
1373 energies : list of float
1375 rates : list of float
1376 Background rate values
1378 Minimum energy (TeV)
1380 Maximum energy (TeV)
1387 for ieng, energy
in enumerate(energies):
1388 if abs(energy-emin) < delta_emin:
1389 delta_emin = abs(energy-emin)
1391 if abs(energy-emax) < delta_emax:
1392 delta_emax = abs(energy-emax)
1396 while rates[i_emax] == 0.0
and i_emax > i_emin:
1398 if rates[i_emax] > 0.0:
1399 x1 = math.log10(energies[i_emin])
1400 x2 = math.log10(energies[i_emax])
1401 y1 = math.log10(rates[i_emin])
1402 y2 = math.log10(rates[i_emax])
1410 return {
'm': m,
't': t}
1414 Determine the power law value
1419 Dictionary of coefficients
1424 value = math.pow(10.0, coeff[
'm'] * math.log10(energy) + coeff[
't'])
1431 Fill all empty bins in 2D histogram by preceeding values in energy
1435 hist2D : `~ROOT.TH2F`
1439 energies = hist2D.GetXaxis()
1440 neng = energies.GetNbins()
1441 offsets = hist2D.GetYaxis()
1442 noffset = offsets.GetNbins()
1446 for ieng
in range(neng):
1447 energy = 0.5 * (math.pow(10.0, energies.GetBinLowEdge(ieng+1)) +
1448 math.pow(10.0, energies.GetBinUpEdge(ieng+1)))
1452 for ioff
in range(noffset):
1455 offset = 0.5 * (offsets.GetBinLowEdge(ioff+1) +
1456 offsets.GetBinUpEdge(ioff+1))
1464 for ieng
in range(neng):
1465 if hist2D.GetBinContent(ieng+1,ioff+1) > 0.0:
1467 value = hist2D.GetBinContent(ieng+1,ioff+1)
1469 for ieng
in range(neng-1, -1, -1):
1470 if hist2D.GetBinContent(ieng+1,ioff+1) > 0.0:
1475 if i_start > -1
and i_stop > -1:
1476 for ieng
in range(i_start, i_stop+1):
1477 if hist2D.GetBinContent(ieng+1,ioff+1) == 0.0:
1478 hist2D.SetBinContent(ieng+1,ioff+1,value)
1481 self._log_value(gammalib.NORMAL,
1482 'Background infill',
1483 '%e cts/(s deg**2) (E=%8.4f TeV (%d) '
1484 'Off=%.2f deg (%d)' %
1485 (value,engs[ieng],ieng,offset,ioff))
1488 value = hist2D.GetBinContent(ieng+1,ioff+1)
1495 Make 2D data column as function of energy and offset angle
1497 If the HDU has already the energy and offset angle columns, this method
1498 will simply add another data column. If name==None, the method will not
1499 append any data column.
1503 array : `~ROOT.TH2F`
1505 hdu : `~gammalib.GFitsHDU`
1511 scale : float, optional
1512 Scaling factor for histogram values
1515 self._log_header3(gammalib.TERSE,
'Make 2D data column "%s"' % name)
1518 energies = array.GetXaxis()
1519 offsets = array.GetYaxis()
1520 neng = energies.GetNbins()
1521 noffset = offsets.GetNbins()
1524 self._log_value(gammalib.NORMAL,
'Number of energies', neng)
1525 self._log_value(gammalib.NORMAL,
'Number of offsets', noffset)
1534 if name !=
None and not hdu.contains(name):
1535 hdu.append(gammalib.GFitsTableFloatCol(name, 1, neng*noffset))
1536 hdu[name].unit(unit)
1537 hdu[name].dim([neng, noffset])
1538 for ioff
in range(noffset):
1539 for ieng
in range(neng):
1540 index = ieng + ioff * neng
1541 value = array.GetBinContent(ieng+1,ioff+1)
1542 hdu[name][0,index] = value * scale
1545 bd_eng =
'ENERG(%.4f-%.2f)TeV' % \
1546 (pow(10.0, energies.GetBinLowEdge(1)), \
1547 pow(10.0, energies.GetBinUpEdge(neng)))
1548 bd_off =
'THETA(%.2f-%.2f)deg' % \
1549 (offsets.GetBinLowEdge(1), \
1550 offsets.GetBinUpEdge(noffset))
1551 bd_phi =
'PHI(0-360)deg'
1552 bounds = [bd_eng, bd_off, bd_phi]
1559 Make 3D data column as function of DETX, DETY and energy
1561 If the HDU has already the energy and offset angle columns, this method
1562 will simply add another data column. If name==None, the method will not
1563 append any data column.
1567 array : `~ROOT.TH2F`
1569 hdu : `~gammalib.GFitsHDU`
1577 self._log_header3(gammalib.TERSE,
'Make 3D data column "%s"' % name)
1580 scale = self[
'bgdscale'].real()
1581 oversample = self[
'bgdoversample'].integer()
1584 deg2sr = 0.01745329*0.01745329
1587 energies = array.GetXaxis()
1588 neng = energies.GetNbins()
1589 offsets = array.GetYaxis()
1590 noffset = offsets.GetNbins()
1591 theta_max = offsets.GetBinUpEdge(noffset)
1593 for ieng
in range(neng):
1594 ewidth.append(pow(10.0, energies.GetBinUpEdge(ieng+1) +6.0) -
1595 pow(10.0, energies.GetBinLowEdge(ieng+1)+6.0))
1598 self._log_value(gammalib.NORMAL,
'Scale', scale)
1599 self._log_value(gammalib.NORMAL,
'Oversample', oversample)
1600 self._log_value(gammalib.NORMAL,
'Number of energies', neng)
1601 self._log_value(gammalib.NORMAL,
'Number of offsets', noffset)
1602 self._log_value(gammalib.NORMAL,
'Maximum offset', theta_max)
1605 ndet = array.GetYaxis().GetNbins()
1606 ndets = 2*ndet*oversample
1607 det_max = array.GetYaxis().GetBinUpEdge(ndet)
1608 det_bin = det_max/float(ndet*oversample)
1612 for i
in range(ndets):
1613 det_lo = -det_max + i*det_bin
1614 det_hi = det_lo + det_bin
1615 det_val = 0.5*(det_lo + det_hi)
1616 dets_lo.append(det_lo)
1617 dets_hi.append(det_hi)
1618 dets2.append(det_val*det_val)
1621 self._log_value(gammalib.NORMAL,
'Number of DETXY bins', ndets)
1622 self._log_value(gammalib.NORMAL,
'Maximum DETXY', det_max)
1623 self._log_value(gammalib.NORMAL,
'DETXY binsize', det_bin)
1636 if name !=
None and not hdu.contains(name):
1637 hdu.append(gammalib.GFitsTableFloatCol(name, 1, ndets*ndets*neng))
1638 hdu[name].unit(unit)
1639 hdu[name].dim([ndets,ndets,neng])
1640 for ix
in range(ndets):
1641 for iy
in range(ndets):
1642 for ieng
in range(neng):
1643 index = ix + (iy + ieng * ndets) * ndets
1644 theta = math.sqrt(dets2[ix] + dets2[iy])
1645 if theta < theta_max:
1646 binsq = array.Interpolate(energies.GetBinCenter(ieng+1),
1648 value = binsq / deg2sr / (ewidth[ieng])
1649 hdu[name][0,index] = value * scale
1652 bd_detx =
'DETX(%.2f-%.2f)deg' % (-det_max, det_max)
1653 bd_dety =
'DETY(%.2f-%.2f)deg' % (-det_max, det_max)
1654 bd_eng =
'ENERG(%.4f-%.2f)TeV' % \
1655 (pow(10.0, energies.GetBinLowEdge(1)),
1656 pow(10.0, energies.GetBinUpEdge(neng)))
1657 bounds = [bd_detx, bd_dety, bd_eng]
1664 Make 3D data column as function of ENERG, MIGRA and THETA
1666 If the HDU has already the energy and offset angle columns, this method
1667 will simply add another data column. If name==None, the method will not
1668 append any data column.
1672 array : `~ROOT.TH3F`
1674 hdu : `~gammalib.GFitsHDU`
1680 scale : float, optional
1681 Scaling factor for histogram values
1684 self._log_header3(gammalib.TERSE,
1685 'Make 3D migration matrix data column "%s"' % name)
1688 etrue = array.GetXaxis()
1689 netrue = etrue.GetNbins()
1690 migra = array.GetYaxis()
1691 nmigra = migra.GetNbins()
1692 offsets = array.GetZaxis()
1693 noffset = offsets.GetNbins()
1695 for ieng
in range(netrue):
1696 ewidth.append(pow(10.0, etrue.GetBinUpEdge(ieng+1) +6.0) - \
1697 pow(10.0, etrue.GetBinLowEdge(ieng+1)+6.0))
1700 self._log_value(gammalib.NORMAL,
'Number of energies', netrue)
1701 self._log_value(gammalib.NORMAL,
'Number of migrations', nmigra)
1702 self._log_value(gammalib.NORMAL,
'Number of offsets', noffset)
1713 if name !=
None and not hdu.contains(name):
1714 hdu.append(gammalib.GFitsTableFloatCol(name, 1, netrue*nmigra*noffset))
1715 hdu[name].unit(unit)
1716 hdu[name].dim([netrue, nmigra, noffset])
1717 for ioff
in range(noffset):
1718 for imigra
in range(nmigra):
1719 for ieng
in range(netrue):
1720 index = ieng + (imigra + ioff * nmigra) * netrue
1721 value = array.GetBinContent(ieng+1,imigra+1,ioff+1)
1722 hdu[name][0,index] = value * scale
1725 bd_eng =
'ENERG(%.4f-%.2f)TeV' % \
1726 (pow(10.0, etrue.GetBinLowEdge(1)),
1727 pow(10.0, etrue.GetBinUpEdge(netrue)))
1728 bd_migra =
'MIGRA(%.3f-%.3f)' % \
1729 (migra.GetBinLowEdge(1),
1730 migra.GetBinUpEdge(nmigra))
1731 bd_off =
'THETA(%.2f-%.2f)deg' % \
1732 (offsets.GetBinLowEdge(1),
1733 offsets.GetBinUpEdge(noffset))
1734 bd_phi =
'PHI(0-360)deg'
1735 bounds = [bd_eng, bd_migra, bd_off, bd_phi]
1748 from ROOT
import TFile, TH2F, TH3F
1755 gammalib.warning(
'csroot2caldb',
'ROOT Python module not present, '
1756 'script will create empty IRF.')
1784 if __name__ ==
'__main__':
def _append_column_values
def _root2psf_get_r68_and_r80
def _root2edisp_get_migra