ATLAS Offline Software
Loading...
Searching...
No Matches
ZDCNLCalibration.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
6// This class has been automatically generated on
7// Fri Aug 26 05:51:09 2016 by ROOT version 6.06/02
8// from TTree zdcTree/ZDC Tree
9// found on file: ../user.steinber.data15_hi.00287931.calibration_zdcCalib.merge.AOD.c932_m1533.42testNL3_output.root
11
12#ifndef ZDCNLCalibration_h
13#define ZDCNLCalibration_h
14
15#include <TROOT.h>
16#include <TTree.h>
17#include <TFile.h>
18
19
20#include <array>
21#include <map>
22#include <algorithm>
23#include <iostream>
24#include <vector>
25#include <string>
26
27#include <TVectorD.h>
28#include <TMatrixD.h>
29#include <TDecompLU.h>
30
31class TH1;
32class TH1D;
33
34struct CalibData
35{
36 std::array<std::vector<float>, 4> weights;
37 unsigned int LBStart;
38 unsigned int LBEnd;
39
41 {
42 std::vector<float> unity(1, 1);
43 weights = {unity, unity, unity, unity};
44 LBStart = 0;
45 LBEnd = 10000;
46 };
47
48CalibData(int inLBStart, int inLBEnd, const std::array<std::vector<float>, 4>& inWeights) :
49 LBStart(inLBStart), LBEnd(inLBEnd), weights(inWeights)
50 {};
51};
52
53// Header file for the classes stored in the TTree if any.
54
56{
57 TFile* m_TFile;
58 TTree* m_tree;
59
63
64
65
66 const float m_SNEnergy;
67 const std::vector<float> m_HEFraction;
68
69 // Declaration of leaf types
70 UInt_t runNumber;
72 UInt_t lumiBlock;
73 UInt_t bcid;
74 UInt_t passBits;
75
76 Float_t zdc_ZdcAmp[2];
77
79 Float_t zdc_ZdcModuleAmp[2][4];
80
81 Bool_t L1_ZDC_A;
82 Bool_t L1_ZDC_C;
83 Bool_t L1_ZDC_AND;
84 Bool_t L1_ZDC_A_C;
85
86 // List of branches
87 TBranch *b_runNumber;
88 TBranch *b_eventNumber;
89 TBranch *b_lumiBlock;
90 TBranch *b_bcid;
91
92 TBranch *b_passBits;
93
94 TBranch *b_zdc_ZdcAmp;
95
98
99 TBranch *b_L1_ZDC_A;
100 TBranch *b_L1_ZDC_C;
101 TBranch *b_L1_ZDC_AND;
102 TBranch *b_L1_ZDC_A_C;
103
104 // Map that keeps lumi block and event associations to optimize
105 // processing of lumi block ranges
106 //
107 typedef std::multimap<unsigned int, std::pair<unsigned int, unsigned int> > LBEvtMap;
109
110public:
111
112 std::array<std::map<std::string, CalibData>, 2> m_calibrations;
113
116 std::array<TH1D*, 4> m_testCalibHEFracHist;
117 std::array<TH1D*, 4> m_testCalibEnergyHist;
118 TTree* m_testTree{};
119
120
121 ZDCNLCalibration(const std::string & file, int maxNLPower = 3, bool useGRL = true, int debugLevel = 0);
122 virtual ~ZDCNLCalibration() {}
123
124 std::array<float, 4> FindSNPeaks(size_t LBLow, size_t LBHigh, size_t side);
125
126 std::pair<float, float> FindSNRange(size_t LBLow, size_t LBHigh, size_t side);
127
128 std::pair<std::pair<float, float>,std::pair<float, float> > FindSNTwoNRanges(size_t LBLow, size_t LBHigh, size_t side);
129
130 void SetDefaultCalibration(size_t side, const CalibData& calib) {
131 AddCalibration(side, "default", calib);
132 }
133
134 void Calibrate(size_t side, const std::string & calibInput, const std::string & calibOutput,
135 size_t LBLow, size_t LBHigh, std::array<int, 4> maxPowerModule,
136 const std::vector<std::pair<double, double> >& nNeutERange,
137 bool excludeHE, float heSumThresh, float HEDeweight);
138
139 void TestCalibration(int side, const std::string& calibName);
140
142 TH1* GetTestFracHist(size_t module) {return (m_haveTest ? m_testCalibHEFracHist.at(module) : 0);}
143
144 TTree* GetTestTree() {return m_testTree;}
145
146private:
147
148 void FillLumiBlockEvtMap();
149
150 void FillMinimizationData(TMatrixD& minimMatrix, TVectorD& minimVector,
151 std::array<int, 4> maxPowerModule, float HEDeweight,
152 const std::vector<std::vector<double> >& sums1DVec, const std::vector<double>& sumsHE,
153 const std::vector<std::vector<double> >& sums2DVec, const std::vector<double>& sumsHE2D);
154
155 void AddCalibration(size_t side, const std::string & tag, const CalibData& calib);
156 CalibData GetCalibration(size_t side, const std::string & tag);
157
158 // Add the four amplitudes from a given event to the set of sums used in the NL weight fitting
159 // we need to make all combinations of amplitude products and for all combinations of powers
160 //
161 void AddToSums(std::vector<double>& sums1D, std::vector<double>& sums2D, float* amps)
162 {
163 size_t index1D = 0;
164 size_t index2D = 0;
165
166 for (size_t module1 : {0, 1, 2, 3}) {
167 double amp1 = amps[module1];
168 double amp1Pow = amp1;
169
170 for (size_t power1 = 0; power1 < m_maxNLPower; power1++) {
171 for (size_t module2 : {0, 1, 2, 3}) {
172 double amp2 = amps[module2];
173 double amp2Pow = amp2;
174
175 for (size_t power2 = 0; power2 < m_maxNLPower; power2++) {
176 sums2D.at(index2D++) += amp1Pow*amp2Pow;
177 amp2Pow *= amp2;
178 }
179 }
180
181 sums1D.at(index1D++) += amp1Pow;
182 amp1Pow *= amp1;
183 }
184 }
185 }
186
187 static double CalculateEnergy(const float* moduleAmps, const CalibData& calib)
188 {
189 double energy = 0;
190
191 for (size_t module : {0, 1, 2, 3}) {
192 float amp = moduleAmps[module];
193 float ampPow = amp;
194
195 for (size_t power = 0; power < calib.weights[module].size(); power++) {
196 energy += calib.weights[module][power]*ampPow;
197 ampPow *= amp;
198 }
199 }
200
201 return energy;
202 }
203
204};
205
206#endif
207
208#ifdef ZDCNLCalibration_cxx
209ZDCNLCalibration::ZDCNLCalibration(const std::string & file, int maxNLPower, bool useGRL, int debugLevel) :
210 m_TFile(0), m_tree(0), m_maxNLPower(maxNLPower),
211 m_useGRL(useGRL), m_debugLevel(debugLevel),
212 m_SNEnergy(2510), m_HEFraction({0.31, 0.27, 0.21, 0.21}),
213 m_haveTest(false)
214{
215 std::cout << "Initializing ZDCNLCalibration with debug level " << m_debugLevel << std::endl;
216
217 TFile* temp_p = new TFile(file.c_str());
218 if (!temp_p->IsOpen()) {
219 std::cout << "Error opening input file " << file << std::endl;
220 return;
221 }
222
223 m_TFile = temp_p;
224 m_tree = static_cast<TTree*>(m_TFile->GetObjectChecked("zdcTree","TTree"));
225 if (!m_tree) {
226 std::cout << "Error reading tree from input file " << file << std::endl;
227 return;
228 }
229
230 m_tree->SetBranchAddress("runNumber", &runNumber, &b_runNumber);
231 m_tree->SetBranchAddress("eventNumber", &eventNumber, &b_eventNumber);
232 m_tree->SetBranchAddress("lumiBlock", &lumiBlock, &b_lumiBlock);
233 m_tree->SetBranchAddress("bcid", &bcid, &b_bcid);
234 m_tree->SetBranchAddress("passBits", &passBits, &b_passBits);
235
236 if (m_tree->FindBranch("zdc_ZdcModuleMask")) {
237 m_tree->SetBranchAddress("zdc_ZdcAmp", zdc_ZdcAmp, &b_zdc_ZdcAmp);
238 m_tree->SetBranchAddress("zdc_ZdcModuleMask", &zdc_ZdcModuleMask, &b_zdc_ZdcModuleMask);
239 m_tree->SetBranchAddress("zdc_ZdcModuleAmp", zdc_ZdcModuleAmp, &b_zdc_ZdcModuleAmp);
240 }
241 else if (m_tree->FindBranch("zdc_ModuleMask")) {
242 m_tree->SetBranchAddress("zdc_OptSumAmp", zdc_ZdcAmp, &b_zdc_ZdcAmp);
243 m_tree->SetBranchAddress("zdc_ModuleMask", (int*) &zdc_ZdcModuleMask, &b_zdc_ZdcModuleMask);
244 m_tree->SetBranchAddress("zdc_OptAmp", zdc_ZdcModuleAmp, &b_zdc_ZdcModuleAmp);
245 }
246 else throw std::runtime_error("ZDCNLCalibration::ZDCNLCalibration valid branch not found");
247
248 m_tree->SetBranchAddress("L1_ZDC_A", &L1_ZDC_A, &b_L1_ZDC_A);
249 m_tree->SetBranchAddress("L1_ZDC_C", &L1_ZDC_C, &b_L1_ZDC_C);
250 m_tree->SetBranchAddress("L1_ZDC_AND", &L1_ZDC_AND, &b_L1_ZDC_AND);
251 m_tree->SetBranchAddress("L1_ZDC_A_C", &L1_ZDC_A_C, &b_L1_ZDC_A_C);
252
253 if (m_debugLevel > 0) {
254 std::cout << "Filling LumiBlock-event map" << std::endl;
255 }
256
257 FillLumiBlockEvtMap();
258
259 // Provide "null" default calibration
260 //
261 AddCalibration(0, "default", CalibData());
262 AddCalibration(1, "default", CalibData());
263}
264
265#endif // #ifdef ZDCNLCalibration_cxx
void TestCalibration(int side, const std::string &calibName)
TBranch * b_zdc_ZdcModuleAmp
virtual ~ZDCNLCalibration()
const std::vector< float > m_HEFraction
std::array< TH1D *, 4 > m_testCalibHEFracHist
void AddCalibration(size_t side, const std::string &tag, const CalibData &calib)
CalibData GetCalibration(size_t side, const std::string &tag)
std::array< std::map< std::string, CalibData >, 2 > m_calibrations
ZDCNLCalibration(const std::string &file, int maxNLPower=3, bool useGRL=true, int debugLevel=0)
void FillMinimizationData(TMatrixD &minimMatrix, TVectorD &minimVector, std::array< int, 4 > maxPowerModule, float HEDeweight, const std::vector< std::vector< double > > &sums1DVec, const std::vector< double > &sumsHE, const std::vector< std::vector< double > > &sums2DVec, const std::vector< double > &sumsHE2D)
TBranch * b_zdc_ZdcModuleMask
static double CalculateEnergy(const float *moduleAmps, const CalibData &calib)
TH1 * GetTestFracHist(size_t module)
std::array< TH1D *, 4 > m_testCalibEnergyHist
std::multimap< unsigned int, std::pair< unsigned int, unsigned int > > LBEvtMap
std::array< float, 4 > FindSNPeaks(size_t LBLow, size_t LBHigh, size_t side)
Float_t zdc_ZdcModuleAmp[2][4]
void SetDefaultCalibration(size_t side, const CalibData &calib)
const float m_SNEnergy
void AddToSums(std::vector< double > &sums1D, std::vector< double > &sums2D, float *amps)
std::pair< float, float > FindSNRange(size_t LBLow, size_t LBHigh, size_t side)
void Calibrate(size_t side, const std::string &calibInput, const std::string &calibOutput, size_t LBLow, size_t LBHigh, std::array< int, 4 > maxPowerModule, const std::vector< std::pair< double, double > > &nNeutERange, bool excludeHE, float heSumThresh, float HEDeweight)
std::pair< std::pair< float, float >, std::pair< float, float > > FindSNTwoNRanges(size_t LBLow, size_t LBHigh, size_t side)
unsigned int LBEnd
unsigned int LBStart
CalibData(int inLBStart, int inLBEnd, const std::array< std::vector< float >, 4 > &inWeights)
std::array< std::vector< float >, 4 > weights
TFile * file