Belle II Software development
eclMergingCrystalEAlgorithm.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <ecl/calibration/eclMergingCrystalEAlgorithm.h>
11
12/* ECL headers. */
13#include <ecl/dbobjects/ECLCrystalCalib.h>
14
15/* ROOT headers. */
16#include <TDirectory.h>
17#include <TFile.h>
18#include <TH1F.h>
19#include <TString.h>
20
21using namespace std;
22using namespace Belle2;
23using namespace ECL;
24using namespace Calibration;
25
28 m_existingGammaGamma("ECLCrystalEnergyGammaGamma"), m_existingMuMu("ECLCrystalEnergyMuMu"),
29 m_existing5x5("ECLCrystalEnergyee5x5"), m_existing(m_payloadName)
30{
32 "Perform energy calibration of ecl crystals by combining previous values from the DB for different calibrations."
33 );
34}
35
37{
38 //------------------------------------------------------------------------
39 // Get the input run list (should be only 1) for us to use to update the DBObjectPtrs
40 auto runs = getRunList();
41 // Take the first run
42 ExpRun chosenRun = runs.front();
43 B2INFO("merging using the ExpRun (" << chosenRun.second << "," << chosenRun.first << ")");
44 // After here your DBObjPtrs are correct
45 updateDBObjPtrs(1, chosenRun.second, chosenRun.first);
46
47 //------------------------------------------------------------------------
48 // Test the DBObjects we want to exist and fail if not all of them do.
49 bool allObjectsFound = true;
50 // Test that the DBObjects are valid
52 allObjectsFound = false;
53 B2ERROR("No valid DBObject found for 'ECLCrystalEnergyGammaGamma'");
54 }
55 if (!m_existingMuMu) {
56 allObjectsFound = false;
57 B2ERROR("No valid DBObject found for 'ECLCrystalEnergyMuMu'");
58 }
59 if (!m_existing5x5) {
60 allObjectsFound = false;
61 B2ERROR("No valid DBObject found for 'ECLCrystalEnergyee5x5'");
62 }
63 if (!m_existing) {
64 allObjectsFound = false;
65 B2ERROR("No valid DBObject found for 'ECLCrystalEnergy'");
66 }
67
68 if (allObjectsFound) {
69 B2INFO("Valid objects found for both 'ECLCrystalEnergy' and 'ECLCrystalEnergyGammaGamma'");
70 } else {
71 B2INFO("Exiting with failure");
72 return c_Failure;
73 }
74
75 //------------------------------------------------------------------------
77 vector<float> gammaGammaCalib = m_existingGammaGamma->getCalibVector();
78 vector<float> gammaGammaCalibUnc = m_existingGammaGamma->getCalibUncVector();
79
80 vector<float> existingCalib = m_existing->getCalibVector();
81 vector<float> existingCalibUnc = m_existing->getCalibUncVector();
82
83 //------------------------------------------------------------------------
85 //..For now, use Gamma Gamma if available; otherwise, use existing value
86 vector<float> newCalib(m_numCrystals);
87 vector<float> newCalibUnc(m_numCrystals);
88
89 for (int ic = 0; ic < m_numCrystals; ic++) {
90 if (gammaGammaCalib[ic] > 0.) {
91 newCalib[ic] = gammaGammaCalib[ic];
92 newCalibUnc[ic] = gammaGammaCalibUnc[ic];
93 } else {
94 newCalib[ic] = existingCalib[ic];
95 newCalibUnc[ic] = existingCalibUnc[ic];
96 }
97 }
98
99 //------------------------------------------------------------------------
101 for (int ic = 0; ic < 9000; ic += 1000) {
102 B2INFO(ic + 1 << " " << existingCalib[ic] << " " << existingCalibUnc[ic] << " "
103 << gammaGammaCalib[ic] << " " << gammaGammaCalibUnc[ic] << " "
104 << newCalib[ic] << " " << newCalibUnc[ic]);
105 }
106
107 //------------------------------------------------------------------------
108 //..Histograms of existing calibration, new calibration, and ratio new/old
109
110 // Just in case, we remember the current TDirectory so we can return to it
111 TDirectory* executeDir = gDirectory;
112
113 TString fname = m_payloadName;
114 fname += ".root";
115 TFile hfile(fname, "recreate");
116
117 TString htitle = m_payloadName;
118 htitle += " existing values;cellID";
119 TH1F* existingPayload = new TH1F("existingPayload", htitle, m_numCrystals, 1, 8737);
120
121 htitle = m_payloadName;
122 htitle += " new values;cellID";
123 TH1F* newPayload = new TH1F("newPayload", htitle, m_numCrystals, 1, 8737);
124
125 htitle = m_payloadName;
126 htitle += " ratio new/old;cellID";
127 TH1F* payloadRatioVsCellID = new TH1F("payloadRatioVsCellID", htitle, m_numCrystals, 1, 8737);
128
129 htitle = m_payloadName;
130 htitle += " ratio new/old";
131 TH1F* payloadRatio = new TH1F("payloadRatio", htitle, 200, 0.95, 1.05);
132
133 for (int cellID = 1; cellID <= m_numCrystals; cellID++) {
134 existingPayload->SetBinContent(cellID, existingCalib[cellID - 1]);
135 existingPayload->SetBinError(cellID, existingCalibUnc[cellID - 1]);
136
137 newPayload->SetBinContent(cellID, newCalib[cellID - 1]);
138 newPayload->SetBinError(cellID, newCalibUnc[cellID - 1]);
139
140 float ratio = 1.;
141 float ratioUnc = 0.;
142 if (abs(existingCalib[cellID - 1]) > 1.0e-12) {
143 ratio = newCalib[cellID - 1] / existingCalib[cellID - 1];
144 float rUnc0 = existingCalibUnc[cellID - 1] / existingCalib[cellID - 1];
145 float rUnc1 = 0.;
146 if (abs(newCalib[cellID - 1]) > 1.0e-12) {rUnc1 = newCalibUnc[cellID - 1] / newCalib[cellID - 1];}
147 ratioUnc = ratio * sqrt(rUnc0 * rUnc0 + rUnc1 * rUnc1);
148 }
149
150 payloadRatioVsCellID->SetBinContent(cellID, ratio);
151 payloadRatioVsCellID->SetBinError(cellID, ratioUnc);
152
153 payloadRatio->Fill(ratio);
154 }
155 hfile.cd();
156 hfile.Write();
157 hfile.Close();
158 B2INFO("Debugging histograms written to " << fname);
159 // Go back to original TDirectory
160 executeDir->cd();
161
162 ECLCrystalCalib* newCrystalEnergy = new ECLCrystalCalib();
163 newCrystalEnergy->setCalibVector(newCalib, newCalibUnc);
164 saveCalibration(newCrystalEnergy, m_payloadName);
165 return c_OK;
166}
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_Failure
Failed =3 in Python.
General DB object to store one calibration number per ECL crystal.
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
DBObjPtr< ECLCrystalCalib > m_existing5x5
ECLCrystalEnergyee5x5 payloads that we want to read from the DB.
DBObjPtr< ECLCrystalCalib > m_existing
ECLCrystalEnergy payloads that we want to read from the DB.
const std::string m_payloadName
The output payload name.
DBObjPtr< ECLCrystalCalib > m_existingGammaGamma
ECLCrystalEnergyGammaGamma payloads that we want to read from the DB.
const int m_numCrystals
Number of Crystals expected.
virtual EResult calibrate() override
..Run algorithm
DBObjPtr< ECLCrystalCalib > m_existingMuMu
ECLCrystalEnergyMuMu payloads that we want to read from the DB.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.
Struct containing exp number and run number.
Definition: Splitter.h:51