29 TString
option = GetOption();
34 edm::LogError(
"HcalCalib") <<
"\nERROR: Failed to read the phi symmetry corrections.\n" 35 <<
"Check if the filename is correct. If the corrections are not needed, set the " 36 "corresponding flag to \"false\"\n" 37 <<
"\nThe program will be terminated\n";
48 h1_trkP =
new TH1F(
"h1_trkP",
"Track momenta; p_{trk} (GeV); Number of tracks", 100, 0, 200);
49 h1_allTrkP =
new TH1F(
"h1_allTrkP",
"Track momenta - all tracks; p_{trk} (GeV); Number of tracks", 100, 0, 200);
52 "h1_selTrkP_iEta10",
"Track momenta - tracks with |iEta|<10; p_{trk} (GeV); Number of tracks", 100, 0, 200);
55 h1_rawSumE =
new TH1F(
"h1_rawSumE",
"Cluster Energy; E_{cl} (GeV); Number of tracks", 100, 0, 200);
57 h1_rawSumE =
new TH1F(
"h1_rawSumE",
"Cluster Energy; E_{cl} (GeV); Number of tracks", 1000, 0, 2000);
59 h1_rawResp =
new TH1F(
"h1_rawResp",
"Uncorrected response: |iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
60 h1_corResp =
new TH1F(
"h1_corResp",
"Corrected response: |iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
63 new TH1F(
"h1_rawRespBarrel",
"Uncorrected response: |iEta|<15; E_{had}/p; Number of tracks", 300, 0, 3);
65 new TH1F(
"h1_corRespBarrel",
"Corrected response: |iEta|<15; E_{had}/p; Number of tracks", 300, 0, 3);
68 new TH1F(
"h1_rawRespEndcap",
"Uncorrected response: 17<|iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
70 new TH1F(
"h1_corRespEndcap",
"Corrected response: 17<|iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
75 "{#Delta}i{#phi} vs {#Delta}i{#eta} of hit and most energetic " 76 "tower(|i{#eta}|<16);{#Delta}i{#eta}; {#Delta}i{#phi}",
84 "{#Delta}i{#phi} vs {#Delta}i{#eta} of hit and most energetic tower (16<|i{#eta}|<25) " 85 ";{#Delta}i{#eta}; {#Delta}i{#phi}",
95 for (Int_t
i = 0;
i < 48; ++
i) {
115 std::set<UInt_t> uniqueIds;
117 Bool_t acceptEvent = kTRUE;
131 std::vector<TCell> selectCells;
133 if (
cells->GetSize() == 0)
139 for (Int_t
i = 0;
i <
cells->GetSize(); ++
i) {
155 selectCells.push_back(*thisCell);
158 if (selectCells.empty()) {
159 edm::LogWarning(
"HcalCalib") <<
"NO CELLS ABOVE THRESHOLD FOUND FOR TARGET!!!";
168 std::pair<Int_t, UInt_t> refPos;
170 Int_t dEtaHitRef = 999;
171 Int_t dPhiHitRef = 999;
179 dEtaHitRef = iEtaMaxE -
iEtaHit;
180 dPhiHitRef = iPhiMaxE -
iPhiHit;
182 if (dPhiHitRef < -36)
219 refPos.first = iEtaMaxE;
220 refPos.second = iPhiMaxE;
225 acceptEvent = kFALSE;
229 acceptEvent = kFALSE;
232 acceptEvent = kFALSE;
234 acceptEvent = kFALSE;
236 acceptEvent = kFALSE;
238 acceptEvent = kFALSE;
253 refPos.first = iEtaMaxE;
254 refPos.second = iPhiMaxE;
256 if (
abs(iEtaMaxE) > 40)
257 acceptEvent = kFALSE;
265 std::vector<Float_t> energies;
266 std::vector<UInt_t> ids;
268 for (std::vector<TCell>::iterator i_it = selectCells.begin(); i_it != selectCells.end(); ++i_it) {
271 if (uniqueIds.insert(i_it->id()).
second) {
272 energies.push_back(i_it->e());
273 ids.push_back(i_it->id());
279 acceptEvent = kFALSE;
281 acceptEvent = kFALSE;
284 acceptEvent = kFALSE;
286 if (
abs(dEtaHitRef) > 1 ||
abs(dPhiHitRef) > 1)
287 acceptEvent = kFALSE;
293 h1_rawSumE->Fill(accumulate(energies.begin(), energies.end(), 0.0));
302 if (
abs(refPos.first) <= 10)
318 <<
"Number of input objects: " <<
cellIds.size()
319 <<
"\nPerforming minimization: depending on selected method can take some time...\n\n";
323 Float_t
weight = (
abs(it_rp->first) < 21) ? 1.0 / 72.0 : 1.0 / 36.0;
332 int numIterations = 10;
341 std::vector<UInt_t> idForSummedCells;
343 for (std::map<UInt_t, Float_t>::iterator m_it =
solution.begin(); m_it !=
solution.end(); ++m_it) {
349 idForSummedCells.push_back(
HcalDetId(m_it->first));
352 for (std::vector<UInt_t>::iterator v_it = idForSummedCells.begin(); v_it != idForSummedCells.end(); ++v_it) {
364 for (std::map<UInt_t, Float_t>::iterator it_s =
solution.begin(); it_s !=
solution.end(); ++it_s) {
413 (ind < 0) ? (ind += 24) : (ind += 23);
414 if (ind >= 0 && ind < 48) {
426 }
else if (maxIEta < 25 && minIEta > 16) {
458 for (Int_t
i = 0;
i < 48; ++
i) {
473 std::map<Int_t, Float_t> aveTargetE;
474 std::map<Int_t, Int_t> nEntries;
477 std::map<Int_t, std::map<Int_t, Float_t> > aveHitE;
499 for (std::map<Int_t, Float_t>::iterator m_it = aveTargetE.begin(); m_it != aveTargetE.end(); ++m_it) {
500 Int_t
iEta = m_it->first;
501 norm = (nEntries[
iEta] > 0) ? 1.0 / (nEntries[
iEta]) : 1.0;
502 aveTargetE[
iEta] *= norm;
504 std::map<Int_t, Float_t>::iterator n_it = (aveHitE[
iEta]).begin();
507 for (; n_it != (aveHitE[
iEta]).
end(); ++n_it) {
508 (n_it->second) *= norm;
509 sumRawE += (n_it->second);
518 std::vector<Int_t> iEtaList;
523 iEtaList.push_back(
i);
526 TMatrixD
A(2 * ONE_SIDE_IETA_RANGE, 2 * ONE_SIDE_IETA_RANGE);
527 TMatrixD
b(2 * ONE_SIDE_IETA_RANGE, 1);
528 TMatrixD
x(2 * ONE_SIDE_IETA_RANGE, 1);
530 for (Int_t
i = 0;
i < 2 * ONE_SIDE_IETA_RANGE; ++
i) {
531 for (Int_t
j = 0;
j < 2 * ONE_SIDE_IETA_RANGE; ++
j) {
536 for (UInt_t
i = 0;
i < iEtaList.size(); ++
i) {
537 b(
i, 0) = aveTargetE[iEtaList[
i]];
539 std::map<Int_t, Float_t>::iterator n_it = aveHitE[iEtaList[
i]].begin();
540 for (; n_it != aveHitE[iEtaList[
i]].end(); ++n_it) {
543 Int_t
j = Int_t(
find(iEtaList.begin(), iEtaList.end(), n_it->first) - iEtaList.begin());
544 A(
i,
j) = n_it->second;
550 Bool_t hasSolution = qrh.MultiSolve(coef);
553 for (UInt_t
i = 0;
i < iEtaList.size(); ++
i) {
563 edm::LogWarning(
"HcalCalib") <<
"\nERROR: Can not find file with phi symmetry constants \"" 581 while (getline(phiSymFile,
line)) {
585 std::istringstream linestream(
line);
597 edm::LogWarning(
"HcalCalib") <<
"\nInvalid detector name in phi symmetry constants file: " <<
sdName.Data()
598 <<
"\nCheck file and rerun!\n";
606 <<
"\nInconsistent info in phi symmetry file: subdet, iEta, iPhi, depth do not match rawId!\n";
612 <<
" subdet=" <<
sdName.Data() <<
" detId=" << detId <<
"\n";
625 TString
sdName[8] = {
"EMPTY",
"HB",
"HE",
"HO",
"HF",
"TRITWR",
"UNDEF",
"OTHER"};
630 fprintf(constFile,
"%1s%16s%16s%16s%16s%9s%11s\n",
"#",
"eta",
"phi",
"depth",
"det",
"value",
"DetId");
634 for (Int_t sd = 1; sd <= 4; sd++) {
635 for (Int_t
e = 1;
e <= 50;
e++) {
638 for (Int_t side = 0; side < 2; side++) {
639 eta = (side == 0) ? -
e :
e;
642 for (Int_t
d = 1;
d < 5;
d++) {
646 Float_t corrFactor = 1.0;
653 Int_t subdetInd = sd;
687 fprintf(constFile,
"%17i%16i%16i%16s%9.5f%11X\n",
eta,
phi,
d,
sdName[sd].Data(), corrFactor,
id.rawId());
void filterCellsInCone(std::vector< TCell > &selectCells, const GlobalPoint hitPositionHcal, Float_t maxConeDist, const CaloGeometry *theCaloGeometry)
std::vector< std::pair< Int_t, UInt_t > > refIEtaIPhi
Log< level::Info, true > LogVerbatim
std::vector< Float_t > targetEnergies
std::map< Int_t, Float_t > iEtaCoefMap
UInt_t eventNumber
pointer to the analyzed TTree or TChain
std::vector< std::vector< UInt_t > > cellIds
TString PHI_SYM_COR_FILENAME
void GetCoefFromMtrxInvOfAve()
void Begin(TTree *tree) override
constexpr int ietaAbs() const
get the absolute value of the cell ieta
TBranch * b_probeJetEmFrac
bool valid(const DetId &id) const override
TLorentzVector * tagJetP4
std::map< UInt_t, Float_t > phiSymCor
TLorentzVector * probeJetP4
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
U second(std::pair< T, U > const &p)
void combinePhi(std::vector< TCell > &selectCells)
constexpr HcalSubdetector subdet() const
get the subdetector
constexpr int ieta() const
get the cell ieta
void Terminate() override
Bool_t APPLY_PHI_SYM_COR_FLAG
Abs< T >::type abs(const T &t)
TH1F * h1_corRespIEta[48]
TH1F * h1_numEventsTwrIEta
const HcalTopology * topo_
std::map< UInt_t, Float_t > solution
std::vector< std::vector< Float_t > > cellEnergies
void filterCells3x3(std::vector< TCell > &selectCells, Int_t iEta, UInt_t iPhi)
Float_t MIN_PROBEJET_ABSETA
Bool_t USE_CONE_CLUSTERING
TString OUTPUT_COR_COEF_FILENAME
void Init(TTree *tree) override
void sumSmallDepths(std::vector< TCell > &selectCells)
Bool_t Process(Long64_t entry) override
void filterCells5x5(std::vector< TCell > &selectCells, Int_t iEta, UInt_t iPhi)
Int_t GetEntry(Long64_t entry, Int_t getall=0) override
IDmap iterate(const std::vector< std::vector< float > > &eventMatrix, const std::vector< std::vector< IDdet > > &idMatrix, const std::vector< float > &energyVector, const int &nIter, const bool &normalizeFlag=false)
Float_t MAX_TAGJET_ABSETA
Float_t MAX_PROBEJET_EMFRAC
const CaloGeometry * theCaloGeometry
void sumDepths(std::vector< TCell > &selectCells)
Log< level::Warning, false > LogWarning
void getIEtaIPhiForHighestE(std::vector< TCell > &selectCells, Int_t &iEta, UInt_t &iPhi)
constexpr int iphi() const
get the cell iphi
constexpr int depth() const
get the tower depth