Belle II Software  release-05-02-19
eclee5x5Algorithm.cc
1 
2 #include <ecl/calibration/eclee5x5Algorithm.h>
3 #include <ecl/dbobjects/ECLCrystalCalib.h>
4 
5 #include "TH2F.h"
6 #include "TFile.h"
7 #include "TF1.h"
8 #include "TROOT.h"
9 #include "TMatrixD.h"
10 #include "TMatrixDSym.h"
11 #include "TDecompLU.h"
12 
13 using namespace std;
14 using namespace Belle2;
15 using namespace ECL;
16 
18 eclee5x5Algorithm::eclee5x5Algorithm(): CalibrationAlgorithm("eclee5x5Collector")
19 {
21  "Perform energy calibration of ecl crystals by analyzing energy in 25-crystal sums from Bhabha events"
22  );
23 }
24 
26 {
29  TH1F* dummy;
30  dummy = (TH1F*)gROOT->FindObject("AverageExpECrys");
31  if (dummy) {delete dummy;}
32  dummy = (TH1F*)gROOT->FindObject("AverageElecCalib");
33  if (dummy) {delete dummy;}
34  dummy = (TH1F*)gROOT->FindObject("AverageInitCalib");
35  if (dummy) {delete dummy;}
36  dummy = (TH1F*)gROOT->FindObject("meanEnvsCrysID");
37  if (dummy) {delete dummy;}
38 
40  gROOT->SetBatch();
41 
44  B2INFO("eclee5x5Algorithm parameters:");
45  B2INFO("outputName = " << m_outputName);
46  B2INFO("minEntries = " << m_minEntries);
47  B2INFO("payloadName = " << m_payloadName);
48  B2INFO("storeConst = " << m_storeConst);
49  if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC" or m_payloadName == "None") {
50  B2INFO("fracLo = " << m_fracLo);
51  B2INFO("fracHiSym = " << m_fracHiSym);
52  B2INFO("fracHiASym = " << m_fracHiASym);
53  B2INFO("nsigLo = " << m_nsigLo);
54  B2INFO("nsigHiSym = " << m_nsigHiSym);
55  B2INFO("nsigHiASym = " << m_nsigHiASym);
56  }
57 
60  auto EnVsCrysID = getObjectPtr<TH2F>("EnVsCrysID");
61  auto RvsCrysID = getObjectPtr<TH1F>("RvsCrysID");
62  auto NRvsCrysID = getObjectPtr<TH1F>("NRvsCrysID");
63  auto Qmatrix = getObjectPtr<TH2F>("Qmatrix");
64  auto ElecCalibvsCrys = getObjectPtr<TH1F>("ElecCalibvsCrys");
65  auto ExpEvsCrys = getObjectPtr<TH1F>("ExpEvsCrys");
66  auto InitialCalibvsCrys = getObjectPtr<TH1F>("InitialCalibvsCrys");
67  auto CalibEntriesvsCrys = getObjectPtr<TH1F>("CalibEntriesvsCrys");
68  auto EntriesvsCrys = getObjectPtr<TH1F>("EntriesvsCrys");
69  auto dPhivsThetaID = getObjectPtr<TH2F>("dPhivsThetaID");
70 
74  TH1F* AverageExpECrys = new TH1F("AverageExpECrys", "Average expected E per crys from collector;Crystal ID;Energy (GeV)", 8736, 0,
75  8736);
76  TH1F* AverageElecCalib = new TH1F("AverageElecCalib", "Average electronics calib const vs crystal;Crystal ID;Calibration constant",
77  8736, 0, 8736);
78  TH1F* AverageInitCalib = new TH1F("AverageInitCalib", "Average initial calib const vs crystal;Crystal ID;Calibration constant",
79  8736, 0, 8736);
80  TH1F* meanEnvsCrysID = new TH1F("meanEnvsCrysID", "Mean normalized energy vs crystal;CrystalID;E/Eexp", 8736, 0, 8736);
81 
82  for (int cellID = 1; cellID <= 8736; cellID++) {
83  double TotEntries = CalibEntriesvsCrys->GetBinContent(cellID);
84  if (TotEntries > 0.) {
85  AverageElecCalib->SetBinContent(cellID, ElecCalibvsCrys->GetBinContent(cellID) / TotEntries);
86  AverageExpECrys->SetBinContent(cellID, ExpEvsCrys->GetBinContent(cellID) / TotEntries);
87  AverageInitCalib->SetBinContent(cellID, InitialCalibvsCrys->GetBinContent(cellID) / TotEntries);
88  }
89 
90  TH1D* En = EnVsCrysID->ProjectionY("En", cellID, cellID);
91  meanEnvsCrysID->SetBinContent(cellID, En->GetMean());
92  meanEnvsCrysID->SetBinError(cellID, En->GetStdDev());
93  }
94 
97  TString fName = m_outputName;
98  TFile* histfile = new TFile(fName, "recreate");
99  EnVsCrysID->Write();
100  RvsCrysID->Write();
101  NRvsCrysID->Write();
102  Qmatrix->Write();
103  AverageElecCalib->Write();
104  AverageExpECrys->Write();
105  AverageInitCalib->Write();
106  EntriesvsCrys->Write();
107  dPhivsThetaID->Write();
108  meanEnvsCrysID->Write();
109 
112  if (m_payloadName == "None") {
113  B2INFO("eclee5x5Algorithm has not been asked to find constants; copying input histograms and quitting");
114  histfile->Close();
115  return c_NotEnoughData;
116  }
117 
120  for (int cellID = 1; cellID <= 8736; cellID++) {
121 
123  if (AverageInitCalib->GetBinContent(cellID) > 0.) {
124  if (EntriesvsCrys->GetBinContent(cellID) < m_minEntries) {
125  histfile->Close();
126  B2INFO("eclee5x5Algorithm: insufficient data for cellID = " << cellID << " " << EntriesvsCrys->GetBinContent(cellID) << " entries");
127  return c_NotEnoughData;
128  }
129  }
130  }
131 
132 
136 
137 
140  bool foundConst = false;
141  TH1F* gVsCrysID = new TH1F("gVsCrysID", "Ratio of new to old calibration vs crystal ID;crystal ID;vector g", 8736, 0, 8736);
142  TH1F* CalibVsCrysID = new TH1F("CalibVsCrysID", "Calibration constant vs crystal ID;crystal ID;counts per GeV", 8736, 0, 8736);
143  TH1F* ExpEnergyperCrys = new TH1F("ExpEnergyperCrys", "Expected energy per crystal;Crystal ID;Energy in 25 crystal sum (GeV)", 8736,
144  0, 8736);
145  TString title = "dPhi cut per crystal " + m_payloadName + ";Crystal ID;dPhi requirement (deg)";
146  TH1F* dPhiperCrys = new TH1F("dPhiperCrys", title, 8736, 0, 8736);
147 
150  if (m_payloadName == "ECLExpee5x5E") {
151  for (int cellID = 1; cellID <= 8736; cellID++) {
152  float mean = meanEnvsCrysID->GetBinContent(cellID);
153  float stdDev = meanEnvsCrysID->GetBinError(cellID);
154  float inputE = AverageExpECrys->GetBinContent(cellID);
155  if (mean > 0. and stdDev > 0.) {
156  ExpEnergyperCrys->SetBinContent(cellID, mean * abs(inputE));
157  ExpEnergyperCrys->SetBinError(cellID, stdDev * abs(inputE));
158  } else {
159  ExpEnergyperCrys->SetBinContent(cellID, inputE);
160  ExpEnergyperCrys->SetBinError(cellID, 0.05 * inputE);
161  }
162  }
163 
164  //..Generate the payload, if requested
165  foundConst = true;
166  if (m_storeConst) {
167  std::vector<float> tempCalib;
168  std::vector<float> tempCalibStdDev;
169 
170  for (int cellID = 1; cellID <= 8736; cellID++) {
171  tempCalib.push_back(ExpEnergyperCrys->GetBinContent(cellID));
172  tempCalibStdDev.push_back(ExpEnergyperCrys->GetBinError(cellID));
173  }
174  ECLCrystalCalib* ExpectedE = new ECLCrystalCalib();
175  ExpectedE->setCalibVector(tempCalib, tempCalibStdDev);
176  saveCalibration(ExpectedE, "ECLExpee5x5E");
177  B2INFO("eclCosmicEAlgorithm: successfully stored expected energies ECLExpee5x5E");
178  }
179 
182  } else if (m_payloadName == "ECLCrystalEnergy5x5") {
183 
184  //..Create the Q matrix and the R vector
185  TMatrixDSym matrixQ(8736);
186  TVectorD vectorR(8736);
187  for (int ix = 1; ix <= 8736; ix++) {
188  vectorR[ix - 1] = RvsCrysID->GetBinContent(ix);
189  for (int iy = 1; iy <= 8736; iy++) {
190  matrixQ[ix - 1][iy - 1] = Qmatrix->GetBinContent(ix, iy);
191  }
192  }
193 
194  //..Crystals that are not being calibrated have no entries in NR. Adjust R and Q to get g=-1
195  int nNotCalibrated = 0;
196  for (int cellID = 1; cellID <= 8736; cellID++) {
197  if (NRvsCrysID->GetBinContent(cellID) < 0.5) {
198  for (int othercell = 1; othercell <= 8736; othercell++) {
199  matrixQ[cellID - 1][othercell - 1] = 0.;
200  matrixQ[othercell - 1][cellID - 1] = 0.;
201  }
202  matrixQ[cellID - 1][cellID - 1] = 1.;
203  vectorR[cellID - 1] = -1.;
204  nNotCalibrated++;
205  }
206  }
207  B2INFO("eclCosmicEAlgorithm: " << nNotCalibrated << " crystals will not be calibrated. ");
208 
209  //..Invert to solve Q g = R
210  TDecompLU lu(matrixQ);
211  bool solved;
212  TVectorD vectorg = lu.Solve(vectorR, solved);
213 
214  //..Fill histograms and check that there are no unexpected negative output values
215  if (solved) {
216  foundConst = true;
217  for (int cellID = 1; cellID <= 8736; cellID++) {
218  gVsCrysID->SetBinContent(cellID, vectorg[cellID - 1]);
219  gVsCrysID->SetBinError(cellID, 0.);
220  float newCalib = vectorg[cellID - 1] * abs(AverageInitCalib->GetBinContent(cellID));
221  CalibVsCrysID->SetBinContent(cellID, newCalib);
222  CalibVsCrysID->SetBinError(cellID, 0.);
223 
224  if (vectorg[cellID - 1] < 0. and NRvsCrysID->GetBinContent(cellID) > 0.) {foundConst = false;}
225  }
226  }
227 
228  //..Generate the payload if requested, and if the matrix inversion worked
229  if (m_storeConst and foundConst) {
230  std::vector<float> tempCalib;
231  std::vector<float> tempCalibStdDev;
232 
233  for (int cellID = 1; cellID <= 8736; cellID++) {
234  tempCalib.push_back(CalibVsCrysID->GetBinContent(cellID));
235  tempCalibStdDev.push_back(CalibVsCrysID->GetBinError(cellID));
236  }
237  ECLCrystalCalib* e5x5ECalib = new ECLCrystalCalib();
238  e5x5ECalib->setCalibVector(tempCalib, tempCalibStdDev);
239  saveCalibration(e5x5ECalib, "ECLCrystalEnergyee5x5");
240  B2INFO("eclCosmicEAlgorithm: successfully stored calibration ECLCrystalEnergyee5x5");
241  }
242 
245  } else if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC") {
246 
247  //..Find mean and sigma of Gaussian fit to ThetaID projections with sufficient statistics
248  float dPhiCenter[69] = {};
249  float dPhiHalfWidth[69] = {};
250  int firstID = 0;
251  int lastID = 0;
252  const int nbins = dPhivsThetaID->GetNbinsX();
253  for (int ib = 1; ib <= nbins; ib++) {
254  TH1D* proj = dPhivsThetaID->ProjectionY("proj", ib, ib);
255  if (proj->Integral() < 100) {continue;}
256  int thetaID = (int)dPhivsThetaID->GetXaxis()->GetBinCenter(ib);
257  if (firstID == 0) { firstID = thetaID; }
258  lastID = thetaID;
259 
260  //..Find the fit limits
261  double fracHi = m_fracHiSym;
262  double nsigHi = m_nsigHiSym;
263  if (thetaID <= m_lastLoThetaID) {
264  fracHi = m_fracHiASym;
265  nsigHi = m_nsigHiASym;
266  }
267 
268  //..Fit range includes all bins with entries>m_fracLo*peak on the low side of the
269  // peak, and entries>m_fracHi*peak on the high side.
270  // i.e. low edge of bin on low side, high edge on high side = low edge of bin+1
271  double peak = proj->GetMaximum();
272  int nPhiBins = proj->GetNbinsX();
273  int binLo = 1;
274  do {
275  binLo++;
276  } while (proj->GetBinContent(binLo) < m_fracLo * peak and binLo < nPhiBins);
277  double xfitLo = proj->GetBinLowEdge(binLo);
278 
279  //..Start search for the upper edge of the fit range at 175 deg to avoid gamma gamma
280  // and e gamma peaks
281  int binHi = proj->GetXaxis()->FindBin(175.01);
282  do {
283  binHi--;
284  } while (proj->GetBinContent(binHi)<fracHi* peak and binHi>1);
285  double xfitHi = proj->GetBinLowEdge(binHi + 1);
286 
287  //..Check that the fit region is sensible
288  int peakBin = proj->GetMaximumBin();
289  if (binLo >= peakBin or binHi <= peakBin) {
290  B2ERROR("Flawed dPhi fit range for thetaID = " << thetaID << " peakBin = " << peakBin << " binLo = " << binLo << "binHi = " <<
291  binHi);
292  }
293 
294  //..Now fit a Gaussian to the selected region
295  proj->Fit("gaus", "", "", xfitLo, xfitHi);
296 
297  //..Find mean, sigma, and selection range. Record center and half-width of range.
298  TF1* fitGaus = proj->GetFunction("gaus");
299  float mean = fitGaus->GetParameter(1);
300  float sigma = fitGaus->GetParameter(2);
301  float dPhiLo = mean - m_nsigLo * sigma;
302  float dPhiHi = mean + nsigHi * sigma;
303  dPhiCenter[thetaID] = 0.5 * (dPhiLo + dPhiHi);
304  dPhiHalfWidth[thetaID] = 0.5 * (dPhiHi - dPhiLo);
305  }
306 
307  //..Pad the thetaID's without fits with the first or last thetaID values
308  for (int thetaID = 0; thetaID < firstID; thetaID++) {
309  dPhiCenter[thetaID] = dPhiCenter[firstID];
310  dPhiHalfWidth[thetaID] = dPhiHalfWidth[firstID];
311  }
312  for (int thetaID = lastID + 1; thetaID < 69; thetaID++) {
313  dPhiCenter[thetaID] = dPhiCenter[lastID];
314  dPhiHalfWidth[thetaID] = dPhiHalfWidth[lastID];
315  }
316 
317  //..Now copy these to each crystal to generate the payload and fill the output histogram.
318  // We will use ECLNeighours to get the number of crystals in each theta ring
320  std::vector<float> tempCalib;
321  std::vector<float> tempCalibWidth;
322  tempCalib.resize(8736);
323  tempCalibWidth.resize(8736);
324  int crysID = 0;
325  for (int thetaID = 0; thetaID < 69; thetaID++) {
326  for (int ic = 0; ic < m_eclNeighbours5x5->getCrystalsPerRing(thetaID); ic++) {
327  tempCalib.at(crysID) = dPhiCenter[thetaID];
328  tempCalibWidth.at(crysID) = dPhiHalfWidth[thetaID];
329  dPhiperCrys->SetBinContent(crysID + 1, dPhiCenter[thetaID]);
330  dPhiperCrys->SetBinError(crysID + 1, dPhiHalfWidth[thetaID]);
331  crysID++;
332  }
333  }
334 
335  //..Store the payload, if requested
336  foundConst = true;
337  if (m_storeConst) {
338  ECLCrystalCalib* eedPhi = new ECLCrystalCalib();
339  eedPhi->setCalibVector(tempCalib, tempCalibWidth);
341  B2INFO("eclCosmicEAlgorithm: successfully stored calibration " << m_payloadName);
342  }
343 
346  } else {
347  B2ERROR("eclee5x5Algorithm: invalid payload name: m_payloadName = " << m_payloadName);
348  }
349 
352  histfile->cd();
353  if (m_payloadName == "ECLExpee5x5E") {
354  ExpEnergyperCrys->Write();
355  } else if (m_payloadName == "ECLCrystalEnergy5x5") {
356  gVsCrysID->Write();
357  CalibVsCrysID->Write();
358  } else if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC") {
359  dPhiperCrys->Write();
360  }
361  histfile->Close();
362 
363  dummy = (TH1F*)gROOT->FindObject("gVsCrysID"); delete dummy;
364  dummy = (TH1F*)gROOT->FindObject("CalibVsCrysID"); delete dummy;
365  dummy = (TH1F*)gROOT->FindObject("ExpEnergyperCrys"); delete dummy;
366  dummy = (TH1F*)gROOT->FindObject("dPhiperCrys"); delete dummy;
367 
370  if (!m_storeConst) {
371  if (foundConst) {
372  B2INFO("eclee5x5Algorithm successfully found constants but was not asked to store them");
373  } else {
374  B2INFO("eclee5x5Algorithm was not asked to store constants, and did not succeed in finding them");
375  }
376  return c_Failure;
377  } else if (!foundConst) {
378  if (m_payloadName == "ECLExpee5x5E") {
379  B2INFO("eclee5x5Algorithm: failed to store expected values");
380  } else if (m_payloadName == "ECLCrystalEnergy5x5") {
381  B2INFO("eclee5x5Algorithm: failed to store calibration constants");
382  } else if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC") {
383  B2INFO("eclee5x5Algorithm: failed to find dPhi* selection criteria");
384  }
385  return c_Failure;
386  }
387  return c_OK;
388 }
Belle2::ECL::eclee5x5Algorithm::m_fracHiASym
double m_fracHiASym
or fracHiASym*peak, at low values of thetaID
Definition: eclee5x5Algorithm.h:131
Belle2::ECL::eclee5x5Algorithm::m_fracLo
double m_fracLo
start dPhi fit where data is > fraclo*peak
Definition: eclee5x5Algorithm.h:129
Belle2::ECL::eclee5x5Algorithm::m_storeConst
bool m_storeConst
write payload to localdb if true
Definition: eclee5x5Algorithm.h:127
Belle2::ECL::eclee5x5Algorithm::m_fracHiSym
double m_fracHiSym
end dPhi fit where data is > fracHiSym*peak
Definition: eclee5x5Algorithm.h:130
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::ECL::eclee5x5Algorithm::m_outputName
std::string m_outputName
..Parameters to control job to find energy calibration using Bhabhas
Definition: eclee5x5Algorithm.h:123
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::ECLCrystalCalib
General DB object to store one calibration number per ECL crystal.
Definition: ECLCrystalCalib.h:34
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::ECL::eclee5x5Algorithm::m_lastLoThetaID
int m_lastLoThetaID
use asymmetric dPhi range for thetaID<= this value
Definition: eclee5x5Algorithm.h:135
Belle2::ECL::eclee5x5Algorithm::m_nsigHiASym
double m_nsigHiASym
or mean+nsigHiASym*sigma at low thetaID
Definition: eclee5x5Algorithm.h:134
Belle2::ECL::eclee5x5Algorithm::m_eclNeighbours5x5
ECL::ECLNeighbours * m_eclNeighbours5x5
Neighbours, used to get nCrys per ring.
Definition: eclee5x5Algorithm.h:128
Belle2::ECL::eclee5x5Algorithm::m_minEntries
int m_minEntries
all crystals to be calibrated must have this many entries
Definition: eclee5x5Algorithm.h:124
Belle2::ECL::ECLNeighbours
Class to get the neighbours for a given cell id.
Definition: ECLNeighbours.h:38
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CalibrationAlgorithm::c_Failure
@ c_Failure
Failed =3 in Python.
Definition: CalibrationAlgorithm.h:54
Belle2::ECL::eclee5x5Algorithm::m_nsigHiSym
double m_nsigHiSym
to mean + nsigHiSym*sigma
Definition: eclee5x5Algorithm.h:133
Belle2::ECL::eclee5x5Algorithm::calibrate
virtual EResult calibrate() override
..Run algorithm on events
Definition: eclee5x5Algorithm.cc:25
Belle2::ECL::eclee5x5Algorithm::m_nsigLo
double m_nsigLo
dPhi region is mean - nsigLo*sigma
Definition: eclee5x5Algorithm.h:132
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::ECL::eclee5x5Algorithm::m_payloadName
std::string m_payloadName
Name of the payload to be stored.
Definition: eclee5x5Algorithm.h:126
Belle2::ECLCrystalCalib::setCalibVector
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
Definition: ECLCrystalCalib.h:48
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::ECL::ECLNeighbours::getCrystalsPerRing
short int getCrystalsPerRing(const short int thetaid) const
return number of crystals in a given theta ring
Definition: ECLNeighbours.h:52