Belle II Software  release-05-01-25
eclEdgeAlgorithm.cc
1 #include <ecl/calibration/eclEdgeAlgorithm.h>
2 #include <ecl/dbobjects/ECLCrystalCalib.h>
3 #include "TH1F.h"
4 #include "TString.h"
5 #include "TFile.h"
6 #include "TDirectory.h"
7 
8 using namespace std;
9 using namespace Belle2;
10 using namespace ECL;
11 using namespace Calibration;
12 
13 
15 eclEdgeAlgorithm::eclEdgeAlgorithm(): CalibrationAlgorithm("eclEdgeCollector")
16 {
18  "Generate payloads ECLCrystalThetaEdge and ECLCrystalPhiEdge found by eclEdgeCollector"
19  );
20 }
21 
22 
25 {
26 
27  //-----------------------------------------------------------------------------------
28  //..Read in histograms created by the collector
29  auto eclCrystalX = getObjectPtr<TH1F>("eclCrystalX");
30  auto eclCrystalY = getObjectPtr<TH1F>("eclCrystalY");
31  auto eclCrystalZ = getObjectPtr<TH1F>("eclCrystalZ");
32  auto eclCrystalR = getObjectPtr<TH1F>("eclCrystalR");
33  auto eclCrystalTheta = getObjectPtr<TH1F>("eclCrystalTheta");
34  auto eclCrystalPhi = getObjectPtr<TH1F>("eclCrystalPhi");
35  auto eclCrystalDirTheta = getObjectPtr<TH1F>("eclCrystalDirTheta");
36  auto eclCrystalDirPhi = getObjectPtr<TH1F>("eclCrystalDirPhi");
37  auto eclCrystalEdgeTheta = getObjectPtr<TH1F>("eclCrystalEdgeTheta");
38  auto eclCrystalEdgePhi = getObjectPtr<TH1F>("eclCrystalEdgePhi");
39 
40  //..And write them to disk
41  TFile* histfile = new TFile("eclEdgeAlgorithm.root", "recreate");
42  eclCrystalX->Write();
43  eclCrystalY->Write();
44  eclCrystalZ->Write();
45  eclCrystalR->Write();
46  eclCrystalTheta->Write();
47  eclCrystalPhi->Write();
48  eclCrystalDirTheta->Write();
49  eclCrystalDirPhi->Write();
50  eclCrystalEdgeTheta->Write();
51  eclCrystalEdgePhi->Write();
52  histfile->Close();
53 
54 
55  //-----------------------------------------------------------------------------------
56  //..Store the payloads
57 
58  //..Edges of the crystals in theta
59  std::vector<float> tempTheta;
60  std::vector<float> tempUnc(8736, 0.);
61  for (int cellID = 1; cellID <= 8736; cellID++) {
62  tempTheta.push_back(eclCrystalEdgeTheta->GetBinContent(cellID));
63  }
64  ECLCrystalCalib* crystalTheta = new ECLCrystalCalib();
65  crystalTheta->setCalibVector(tempTheta, tempUnc);
66  saveCalibration(crystalTheta, "ECLCrystalThetaEdge");
67  B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalThetaEdge");
68 
69  //..Edges of the crystals in phi
70  std::vector<float> tempPhi;
71  for (int cellID = 1; cellID <= 8736; cellID++) {
72  tempPhi.push_back(eclCrystalEdgePhi->GetBinContent(cellID));
73  }
74  ECLCrystalCalib* crystalPhi = new ECLCrystalCalib();
75  crystalPhi->setCalibVector(tempPhi, tempUnc);
76  saveCalibration(crystalPhi, "ECLCrystalPhiEdge");
77  B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalPhiEdge");
78 
79  return c_OK;
80 }
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::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::eclEdgeAlgorithm::calibrate
virtual EResult calibrate() override
..Run algorithm
Definition: eclEdgeAlgorithm.cc:24
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
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