CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
hcalCalib.cc
Go to the documentation of this file.
1 #define hcalCalib_cxx
2 
3 // TSelector-based code for getting the HCAL resp. correction
4 // from physics events. Works for DiJet and IsoTrack calibration.
5 //
6 // Anton Anastassov (Northwestern)
7 // Email: aa@fnal.gov
8 //
9 //
10 
13 
14 #include <TH2.h>
15 #include <TStyle.h>
16 #include "TFile.h"
17 
18 #include <iostream>
19 #include <fstream>
20 
21 #include <sstream>
22 #include <string>
23 
24 #include <map>
25 #include <numeric>
26 #include <algorithm>
27 #include <set>
28 
30 
31 #include "TMatrixF.h"
32 #include "TMatrixD.h"
33 #include "TDecompSVD.h"
34 #include "TDecompQRH.h"
35 
38 
39 
40 using namespace std;
41 
42 UInt_t nEvents;
43 
44 TFile* histoFile;
45 
46 // sanity check histograms
47 TH1F* h1_trkP;
48 TH1F* h1_allTrkP;
49 
51 
52 TH1F* h1_rawSumE;
53 TH1F* h1_rawResp;
54 TH1F* h1_corResp;
60 
63 
64 // histograms based on iEta, iPhi of refPosition forthe cluster (at the moment: hottest tower)
65 // expect range |iEta|<=24 (to do: add flexibility for arbitrary range)
66 TH1F* h1_corRespIEta[48];
67 
68 
69 void hcalCalib::Begin(TTree * /*tree*/) {
70  TString option = GetOption();
71 
72  nEvents = 0;
73 
74  if (APPLY_PHI_SYM_COR_FLAG && !ReadPhiSymCor()) {
75  cout << "\nERROR: Failed to read the phi symmetry corrections." << endl;
76  cout << "Check if the filename is correct. If the corrections are not needed, set the corresponding flag to \"false\"\n" << endl;
77 
78  cout << "\nThe program will be terminated\n" << endl;
79 
80  exit(1);
81 
82  }
83 
84 
85  // cellEnergies.reserve(1000000);
86  // cellIds.reserve(1000000);
87  // targetEnergies.reserve(1000000);
88 
89  histoFile = new TFile(HISTO_FILENAME.Data(), "RECREATE");
90 
91 
92  h1_trkP = new TH1F("h1_trkP", "Track momenta; p_{trk} (GeV); Number of tracks", 100, 0, 200);
93  h1_allTrkP = new TH1F("h1_allTrkP", "Track momenta - all tracks; p_{trk} (GeV); Number of tracks", 100, 0, 200);
94 
95  h1_selTrkP_iEta10 = new TH1F("h1_selTrkP_iEta10", "Track momenta - tracks with |iEta|<10; p_{trk} (GeV); Number of tracks", 100, 0, 200);
96 
97 
98 
99  if (CALIB_TYPE=="ISO_TRACK")
100  h1_rawSumE = new TH1F("h1_rawSumE", "Cluster Energy; E_{cl} (GeV); Number of tracks", 100, 0, 200);
101  else
102  h1_rawSumE = new TH1F("h1_rawSumE", "Cluster Energy; E_{cl} (GeV); Number of tracks", 1000, 0, 2000);
103 
104 
105 
106  h1_rawResp = new TH1F("h1_rawResp", "Uncorrected response: |iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
107  h1_corResp = new TH1F("h1_corResp", "Corrected response: |iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
108 
109  h1_rawRespBarrel = new TH1F("h1_rawRespBarrel", "Uncorrected response: |iEta|<15; E_{had}/p; Number of tracks", 300, 0, 3);
110  h1_corRespBarrel = new TH1F("h1_corRespBarrel", "Corrected response: |iEta|<15; E_{had}/p; Number of tracks", 300, 0, 3);
111 
112  h1_rawRespEndcap = new TH1F("h1_rawRespEndcap", "Uncorrected response: 17<|iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
113  h1_corRespEndcap = new TH1F("h1_corRespEndcap", "Corrected response: 17<|iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
114 
115  h1_numEventsTwrIEta = new TH1F("h1_numEventsTwrIEta", "h1_numEventsTwrIEta", 80, -40, 40);
116 
117  h2_dHitRefBarrel = new TH2F("h2_dHitRefBarrel","{#Delta}i{#phi} vs {#Delta}i{#eta} of hit and most energetic tower(|i{#eta}|<16);{#Delta}i{#eta}; {#Delta}i{#phi}", 10, -5, 5, 10, -5, 5);
118  h2_dHitRefEndcap = new TH2F("h2_dHitRefEndcap","{#Delta}i{#phi} vs {#Delta}i{#eta} of hit and most energetic tower (16<|i{#eta}|<25) ;{#Delta}i{#eta}; {#Delta}i{#phi}", 10, -5, 5, 10, -5, 5);
119 
120 
121 
122  TString histoName = "isoTrack_";
123 
124  for (Int_t i=0; i<48; ++i) {
125  Long_t iEta;
126  if (i<24) iEta = i-24;
127  else iEta = i-23;
128  TString hn = histoName + iEta;
129  h1_corRespIEta[i] = new TH1F(hn, hn, 300, 0, 3.0);
130  }
131 
132 
133 } // end of Begin()
134 
135 
136 
137 //void hcalCalib::SlaveBegin(TTree * /*tree*/) {
138 // TString option = GetOption();
139 //}
140 
141 
142 Bool_t hcalCalib::Process(Long64_t entry) {
143 
144  // fChain->GetTree()->GetEntry(entry);
145  GetEntry(entry);
146 
147 
148  set<UInt_t> uniqueIds; // for testing: check if there are duplicate cells (AA)
149 
150 
151  Bool_t acceptEvent = kTRUE;
152 
153  ++nEvents;
154 
155  if (!(nEvents%100000)) cout << "event: " << nEvents << endl;
156 
157 
158  h1_allTrkP->Fill(targetE);
159 
160  if (targetE < MIN_TARGET_E || targetE > MAX_TARGET_E) return kFALSE;;
161 
162 
163  // make local copy as the cells may be modified due to phi/depth sum, phi corrections etc
164  vector<TCell> selectCells;
165 
166 
167  if (cells->GetSize()==0) return kFALSE;
168 
169  if (CALIB_TYPE=="DI_JET" && probeJetEmFrac > 0.999) return kTRUE;
170 
171 
172  for (Int_t i=0; i<cells->GetSize(); ++i) {
173  TCell* thisCell = (TCell*) cells->At(i);
174 
175  if (HcalDetId(thisCell->id()).subdet()== HcalOuter) continue; // reject HO, make a switch!
176 
177  if (HcalDetId(thisCell->id()).subdet() != HcalBarrel &&
178  HcalDetId(thisCell->id()).subdet() != HcalEndcap &&
179  HcalDetId(thisCell->id()).subdet() != HcalForward) {
180 
181  cout << "Unknown or wrong hcal subdetector: " << HcalDetId(thisCell->id()).subdet() << endl;
182 
183  }
184 
185 
186  // Apply phi symmetry corrections if the flag is set
187  if (APPLY_PHI_SYM_COR_FLAG) thisCell->SetE(phiSymCor[thisCell->id()] * thisCell->e());
188 
189  if (thisCell->e() > MIN_CELL_E) selectCells.push_back(*thisCell);
190  }
191 
192 
193  if (selectCells.size()==0) {
194  cout << "NO CELLS ABOVE THRESHOLD FOUND FOR TARGET!!!" << endl;
195  }
196 
197 
198  if (SUM_DEPTHS) sumDepths(selectCells);
199  else if (SUM_SMALL_DEPTHS) sumSmallDepths(selectCells); // depth 1,2 in twrs 15,16
200 
201 
202 
203  // most energetic tower (IsoTracks) or centroid of probe jet (DiJets)
204  pair<Int_t, UInt_t> refPos;
205 
206  Int_t dEtaHitRef = 999;
207  Int_t dPhiHitRef = 999;
208 
209  if (CALIB_TYPE=="ISO_TRACK") {
210 
211  Int_t iEtaMaxE; // filled by reference in getIEtaIPhiForHighestE
212  UInt_t iPhiMaxE; //
213 
214  getIEtaIPhiForHighestE(selectCells, iEtaMaxE, iPhiMaxE);
215 
216  dEtaHitRef = iEtaMaxE - iEtaHit;
217  dPhiHitRef = iPhiMaxE - iPhiHit;
218 
219  if (dPhiHitRef < -36) dPhiHitRef += 72;
220  if (dPhiHitRef > 36) dPhiHitRef -= 72;
221 
222  if (iEtaHit*iEtaMaxE < 0) {
223  if (dEtaHitRef<0) dEtaHitRef += 1;
224  if (dEtaHitRef>0) dEtaHitRef -= 1;
225  }
226 
227 
228  if (abs(iEtaHit)<16) h2_dHitRefBarrel->Fill(dEtaHitRef, dPhiHitRef);
229  if (abs(iEtaHit)>16 && abs(iEtaHit)<25) h2_dHitRefEndcap->Fill(dEtaHitRef, dPhiHitRef);
230 
231  // --------------------------------------------------
232  // Choice of cluster definition
233  //
234  // fixed size NxN clusters as specified in to config file
235  if (!USE_CONE_CLUSTERING) {
236  if (abs(iEtaMaxE)<16 && HB_CLUSTER_SIZE==3) filterCells3x3(selectCells, iEtaMaxE, iPhiMaxE);
237  if (abs(iEtaMaxE)>15 && HE_CLUSTER_SIZE==3) filterCells3x3(selectCells, iEtaMaxE, iPhiMaxE);
238 
239  if (abs(iEtaMaxE)<16 && HB_CLUSTER_SIZE==5) filterCells5x5(selectCells, iEtaMaxE, iPhiMaxE);
240  if (abs(iEtaMaxE)>15 && HE_CLUSTER_SIZE==5) filterCells5x5(selectCells, iEtaMaxE, iPhiMaxE);
241  }
242  else {
243  // calculate distance at hcal surface
244  const GlobalPoint hitPositionHcal(xTrkHcal, yTrkHcal, zTrkHcal);
245  filterCellsInCone(selectCells, hitPositionHcal, MAX_CONE_DIST, theCaloGeometry);
246  }
247 
248 
249 
250  refPos.first = iEtaMaxE;
251  refPos.second = iPhiMaxE;
252 
253  }
254  else if (CALIB_TYPE=="DI_JET") { // Apply selection cuts on DiJet events here
255 
256  if (etVetoJet>MAX_ET_THIRD_JET) acceptEvent = kFALSE;
257 
258  Float_t jetsDPhi = probeJetP4->DeltaPhi(*tagJetP4);
259  if (fabs(jetsDPhi * 180.0/M_PI) < MIN_DPHI_DIJETS) acceptEvent = kFALSE;
260 
261  if (probeJetEmFrac>MAX_PROBEJET_EMFRAC) acceptEvent = kFALSE;
262  if (fabs(probeJetP4->Eta())<MIN_PROBEJET_ABSETA) acceptEvent = kFALSE;
263  if (fabs(tagJetP4->Eta())>MAX_TAGJET_ABSETA) acceptEvent = kFALSE;
264  if (fabs(tagJetP4->Et())< MIN_TAGJET_ET) acceptEvent = kFALSE;
265 
266  if (acceptEvent) {
267 
268  Int_t iEtaMaxE; // filled by reference in getIEtaIPhiForHighestE
269  UInt_t iPhiMaxE; //
270 
271  getIEtaIPhiForHighestE(selectCells, iEtaMaxE, iPhiMaxE);
272 
273  // The ref position for the jet is not used in the minimization at this time.
274  // It will be needed if we attempt to do matrix inversion: then the question is
275  // which value is better suited: the centroid of the jet or the hottest tower...
276 
277  // refPos.first = iEtaHit;
278  // refPos.second = iPhiHit;
279 
280  refPos.first = iEtaMaxE;
281  refPos.second = iPhiMaxE;
282 
283  if (abs(iEtaMaxE)>40) acceptEvent = kFALSE; // for testing : set as parameter (AA)
284 
285  }
286 
287  }
288 
289 
290  if (COMBINE_PHI) combinePhi(selectCells);
291 
292  // fill the containers for the minimization prcedures
293  vector<Float_t> energies;
294  vector<UInt_t> ids;
295 
296  for (vector<TCell>::iterator i_it = selectCells.begin(); i_it!=selectCells.end(); ++i_it) {
297 
298  // for testing : fill only unique id's
299 
300  if (uniqueIds.insert(i_it->id()).second) {
301  energies.push_back(i_it->e());
302  ids.push_back(i_it->id());
303  }
304 
305  }
306 
307  if (CALIB_TYPE=="ISO_TRACK") {
308  if (accumulate(energies.begin(), energies.end(), 0.0) / targetE < MIN_EOVERP) acceptEvent=kFALSE;
309  if (accumulate(energies.begin(), energies.end(), 0.0) / targetE > MAX_EOVERP) acceptEvent=kFALSE;
310 
311  if (emEnergy > MAX_TRK_EME) acceptEvent=kFALSE;
312 
313  if (abs(dEtaHitRef)>1 || abs(dPhiHitRef)>1) acceptEvent=kFALSE;
314 
315  // Have to check if for |iEta|>20 (and neighboring region) the dPhiHitRef
316  // should be relaxed to 2. The neighboring towers have dPhi=2...
317 
318 
319  }
320 
321 
322  h1_rawSumE->Fill(accumulate(energies.begin(), energies.end(), 0.0));
323 
324 
325  // here we fill the information for the minimization procedure
326  if (acceptEvent) {
327  cellEnergies.push_back(energies);
328  cellIds.push_back(ids);
329  targetEnergies.push_back(targetE);
330  refIEtaIPhi.push_back(refPos);
331 
332  if (abs(refPos.first)<=10)
333  h1_selTrkP_iEta10->Fill(targetE);
334 
335 
336  }
337 
338 
339 
340 
341  // Clean up
342  energies.clear();
343  ids.clear();
344  selectCells.clear();
345 
346  return kTRUE;
347 }
348 
349 //void hcalCalib::SlaveTerminate() {}
350 
352 
353  cout << "\n\nFinished reading the events.\n";
354  cout << "Number of input objects: " << cellIds.size() << endl;
355  cout << "Performing minimization: depending on selected method can take some time...\n\n";
356 
357  for (vector<pair<Int_t, UInt_t> >::iterator it_rp=refIEtaIPhi.begin(); it_rp!=refIEtaIPhi.end(); ++it_rp ) {
358  Float_t weight = (abs(it_rp->first)<21)? 1.0/72.0 : 1.0/36.0;
359  h1_numEventsTwrIEta->Fill(it_rp->first, weight);
360  }
361 
362  if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
363  GetCoefFromMtrxInvOfAve();
364  }
365  else if (CALIB_METHOD=="L3" || CALIB_METHOD=="L3_AND_MTRX_INV") {
366  int eventWeight = 2; // 2 is the default (try at some point 0,1,2,3)
367  MinL3AlgoUniv<UInt_t>* thisL3Algo = new MinL3AlgoUniv<UInt_t>(eventWeight);
368  int numIterations = 10; // default 10
369 
370  solution = thisL3Algo->iterate(cellEnergies, cellIds, targetEnergies, numIterations);
371 
372  // in order to handle the case where sumDepths="false", but the flag to sum depths 1,2 in HB towers 15, 16
373  // is set (sumSmallDepths) we create entries in "solution" to create equal correction here
374  // for each encountered coef in depth one.
375 
376  if (!SUM_DEPTHS && SUM_SMALL_DEPTHS) {
377  vector<UInt_t> idForSummedCells;
378 
379  for (map<UInt_t, Float_t>::iterator m_it = solution.begin(); m_it != solution.end(); ++m_it) {
380  if (HcalDetId(m_it->first).ietaAbs()!=15 && HcalDetId(m_it->first).ietaAbs()!=16) continue;
381  if (HcalDetId(m_it->first).subdet()!=HcalBarrel) continue;
382  if (HcalDetId(m_it->first).depth()==1)
383  idForSummedCells.push_back(HcalDetId(m_it->first));
384  }
385 
386  for (vector<UInt_t>::iterator v_it=idForSummedCells.begin(); v_it!=idForSummedCells.end(); ++v_it) {
387  UInt_t addCoefId = HcalDetId(HcalBarrel, HcalDetId(*v_it).ieta(), HcalDetId(*v_it).iphi(), 2);
388  solution[addCoefId] = solution[*v_it];
389  }
390 
391  } // end of special treatment for "sumSmallDepths" mode
392 
393 
394  if (CALIB_METHOD=="L3_AND_MTRX_INV") {
395  GetCoefFromMtrxInvOfAve();
396 
397  // loop over the solution from L3 and multiply by the additional factor from
398  // the matrix inversion. Set coef outside of the valid calibration region =1.
399  for (map<UInt_t, Float_t>::iterator it_s=solution.begin(); it_s!=solution.end(); ++it_s) {
400  Int_t iEtaSol = HcalDetId(it_s->first).ieta();
401  if (abs(iEtaSol)<CALIB_ABS_IETA_MIN || abs(iEtaSol)> CALIB_ABS_IETA_MAX) it_s->second = 1.0;
402  else it_s->second *= iEtaCoefMap[iEtaSol];
403  }
404 
405  } // if CALIB_METHOD=="L3_AND_MTRX_INV"
406 
407  } // end of L3 or L3 + mtrx inversion minimization
408 
409  // done getting the constants -> write the formatted file
410  makeTextFile();
411 
412  // fill some histograms
413  Float_t rawResp = 0;
414  Float_t corResp = 0;
415  Int_t maxIEta = 0;
416  Int_t minIEta = 999;
417 
418 
419  for (unsigned int i=0; i<cellEnergies.size(); ++i) {
420  Int_t iEta;
421  for (unsigned int j=0; j<(cellEnergies[i]).size(); ++j) {
422  iEta = HcalDetId(cellIds[i][j]).ieta();
423  rawResp += (cellEnergies[i])[j];
424 
425  if (CALIB_METHOD=="L3_AND_MTRX_INV") {
426  corResp += (solution[cellIds[i][j]] * cellEnergies[i][j]);
427  }
428  else if (CALIB_METHOD=="L3") {
429  corResp += (solution[cellIds[i][j]] * (cellEnergies[i][j]));
430  }
431  else if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
432  corResp += (iEtaCoefMap[iEta]* cellEnergies[i][j]);
433  }
434 
435  if (maxIEta < abs(iEta) ) maxIEta = abs(iEta);
436  if (minIEta > abs(iEta) ) minIEta = abs(iEta);
437  }
438 
439  rawResp /= targetEnergies[i];
440  corResp /= targetEnergies[i];
441 
442 
443  // fill histograms based on iEta on ref point of the cluster (for now: hottest tower)
444  // expect range |iEta|<=24 (to do: add flexibility for arbitrary range)
445 
446  if (CALIB_TYPE=="ISO_TRACK") {
447  Int_t ind = refIEtaIPhi[i].first;
448  (ind<0) ? (ind +=24) : (ind +=23);
449  if (ind>=0 && ind<48) {
450  h1_corRespIEta[ind]->Fill(corResp);
451  }
452 
453  // fill histograms for cases where all towers are in the barrel or endcap
454  if (maxIEta<25) {
455  h1_rawResp->Fill(rawResp);
456  h1_corResp->Fill(corResp);
457  }
458  if (maxIEta<15) {
459  h1_rawRespBarrel->Fill(rawResp);
460  h1_corRespBarrel->Fill(corResp);
461  }
462  else if (maxIEta<25 && minIEta>16) {
463  h1_rawRespEndcap->Fill(rawResp);
464  h1_corRespEndcap->Fill(corResp);
465  }
466  } // histograms for isotrack calibration
467 
468 
469  else {
470 
471  // put jet plots here
472 
473  }
474 
475 
476 
477  rawResp = 0;
478  corResp = 0;
479  maxIEta = 0;
480  minIEta = 999;
481  }
482 
483 
484  // save the histograms
485  h1_trkP->Write();
486  h1_allTrkP->Write();
487 
488 h1_selTrkP_iEta10->Write();
489 
490  h1_rawSumE->Write();
491  h1_rawResp->Write();
492  h1_corResp->Write();
493  h1_rawRespBarrel->Write();
494  h1_corRespBarrel->Write();
495  h1_rawRespEndcap->Write();
496  h1_corRespEndcap->Write();
497  h1_numEventsTwrIEta->Write();
498  h2_dHitRefBarrel->Write();
499  h2_dHitRefEndcap->Write();
500  for (Int_t i=0; i<48; ++i) {
501  h1_corRespIEta[i]->Write();
502  }
503 
504 
505 
506 
507 
508  histoFile->Write();
509  histoFile->Close();
510 
511  cout << "\n Finished calibration.\n " << endl;
512 
513 
514 } // end of Terminate()
515 
516 
517 
518 
519 //**************************************************************************************************************
520 
522 
523  // these maps are keyed by iEta
524  map<Int_t, Float_t> aveTargetE;
525  map<Int_t, Int_t> nEntries; // count hits
526 
527  // iEtaRef iEtaCell, energy
528  map<Int_t, map<Int_t, Float_t> > aveHitE; // add energies in the loop, normalize after that
529 
530  for (unsigned int i=0; i<cellEnergies.size(); ++i) {
531  Int_t iEtaRef = refIEtaIPhi[i].first;
532  aveTargetE[iEtaRef] += targetEnergies[i];
533  nEntries[iEtaRef]++;
534 
535  // for hybrid method: matrix inv of averages preceeded by L3
536  if (CALIB_METHOD=="L3_AND_MTRX_INV") {
537  for (unsigned int j=0; j<(cellEnergies[i]).size(); ++j) {
538  aveHitE[iEtaRef][HcalDetId(cellIds[i][j]).ieta()] += (solution[cellIds[i][j]] * cellEnergies[i][j]);
539  }
540  }
541  else if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
542  for (unsigned int j=0; j<(cellEnergies[i]).size(); ++j) {
543  aveHitE[iEtaRef][HcalDetId(cellIds[i][j]).ieta()] += cellEnergies[i][j];
544  }
545  }
546 
547  } // end of loop of entries
548 
549 
550  // scale by number of entries to get the averages
551  Float_t norm = 1.0;
552  for (map<Int_t, Float_t>::iterator m_it = aveTargetE.begin(); m_it!=aveTargetE.end(); ++m_it){
553  Int_t iEta = m_it->first;
554  norm = (nEntries[iEta] > 0)? 1.0/(nEntries[iEta]) : 1.0;
555  aveTargetE[iEta] *= norm;
556 
557  map<Int_t, Float_t>::iterator n_it = (aveHitE[iEta]).begin();
558 
559  Float_t sumRawE = 0;
560  for (; n_it!=(aveHitE[iEta]).end(); ++n_it) {
561  (n_it->second) *= norm;
562  sumRawE += (n_it->second);
563  }
564 
565  } // end of scaling by number of entries
566 
567  Int_t ONE_SIDE_IETA_RANGE = CALIB_ABS_IETA_MAX - CALIB_ABS_IETA_MIN +1;
568 
569  // conversion from iEta to index for the linear system
570  // contains elements only in the valid range for *matrix inversion*
571  vector<Int_t> iEtaList;
572 
573  for (Int_t i=-CALIB_ABS_IETA_MAX; i<=CALIB_ABS_IETA_MAX; ++i) {
574  if (abs(i)<CALIB_ABS_IETA_MIN) continue;
575  iEtaList.push_back(i);
576  }
577 
578  TMatrixD A(2*ONE_SIDE_IETA_RANGE, 2*ONE_SIDE_IETA_RANGE);
579  TMatrixD b(2*ONE_SIDE_IETA_RANGE, 1);
580  TMatrixD x(2*ONE_SIDE_IETA_RANGE, 1);
581 
582  for (Int_t i=0; i<2*ONE_SIDE_IETA_RANGE; ++i) {
583  for (Int_t j=0; j<2*ONE_SIDE_IETA_RANGE; ++j) {
584  A(i, j) = 0.0;
585  }
586  }
587 
588  for (UInt_t i=0; i<iEtaList.size(); ++i) {
589  b(i, 0) = aveTargetE[iEtaList[i]];
590 
591  map<Int_t, Float_t>::iterator n_it = aveHitE[iEtaList[i]].begin();
592  for (; n_it!=aveHitE[iEtaList[i]].end(); ++n_it){
593 
594  if (fabs(n_it->first)>CALIB_ABS_IETA_MAX || fabs(n_it->first)<CALIB_ABS_IETA_MIN) continue;
595  Int_t j = Int_t(find(iEtaList.begin(),iEtaList.end(), n_it->first) - iEtaList.begin());
596  A(i, j) = n_it->second;
597 
598  }
599 
600  }
601 
602  TMatrixD coef = b;
603  TDecompQRH qrh(A);
604  Bool_t hasSolution = qrh.MultiSolve(coef);
605 
606  if (hasSolution && (CALIB_METHOD=="L3_AND_MTRX_INV" || CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE")) {
607  for (UInt_t i=0; i<iEtaList.size(); ++i) {
608  iEtaCoefMap[iEtaList[i]] = coef(i, 0);
609  }
610  }
611 }
612 
613 
615 
616  std::ifstream phiSymFile(PHI_SYM_COR_FILENAME.Data());
617 
618  if (!phiSymFile) {
619  cout << "\nERROR: Can not find file with phi symmetry constants \"" << PHI_SYM_COR_FILENAME.Data() << "\"" << endl;
620  return kFALSE;
621  }
622 
623  // assumes the format used in CSA08, first line is a comment
624 
625  Int_t iEta;
626  UInt_t iPhi;
627  UInt_t depth;
628  TString sdName;
629  UInt_t detId;
630 
631  Float_t value;
633 
635 
636  while (getline(phiSymFile, line)) {
637 
638  if(!line.size() || line[0]=='#') continue;
639 
640  std::istringstream linestream(line);
641  linestream >> iEta >> iPhi >> depth >> sdName >> value >> hex >> detId;
642 
643  if (sdName=="HB") sd = HcalBarrel;
644  else if (sdName=="HE") sd = HcalEndcap;
645  else if (sdName=="HO") sd = HcalOuter;
646  else if (sdName=="HF") sd = HcalForward;
647  else {
648  cout << "\nInvalid detector name in phi symmetry constants file: " << sdName.Data() << endl;
649  cout << "Check file and rerun!\n" << endl;
650  return kFALSE;
651  }
652 
653  // check if the data is consistent
654 
655  if (HcalDetId(sd, iEta, iPhi, depth) != HcalDetId(detId)) {
656  cout << "\nInconsistent info in phi symmetry file: subdet, iEta, iPhi, depth do not match rawId!\n" << endl;
657  return kFALSE;
658  }
659  HcalDetId hId(detId);
660  if (!topo_->valid(hId)) {
661  cout << "\nInvalid DetId from: iEta=" << iEta << " iPhi=" << iPhi << " depth=" << depth
662  << " subdet=" << sdName.Data() << " detId=" << detId << endl << endl;
663  return kFALSE;
664  }
665 
666  phiSymCor[HcalDetId(sd, iEta, iPhi, depth)] = value;
667  }
668 
669  return kTRUE;
670 }
671 
672 
674 
675  //{ HcalEmpty=0, HcalBarrel=1, HcalEndcap=2, HcalOuter=3, HcalForward=4, HcalTriggerTower=5, HcalOther=7 };
676 
677  TString sdName[8] = {"EMPTY", "HB", "HE", "HO", "HF", "TRITWR", "UNDEF", "OTHER"};
678 
679  FILE *constFile = fopen(OUTPUT_COR_COEF_FILENAME.Data(), "w+");
680 
681  // header of the constants file
682  fprintf(constFile, "%1s%16s%16s%16s%16s%9s%11s\n",
683  "#", "eta", "phi", "depth", "det", "value", "DetId");
684 
685  // Order loops to produce sequence of constants as for phi symmetry
686 
687  for(Int_t sd=1; sd<=4; sd++) {
688 
689  for(Int_t e=1; e<=50; e++) {
690  Int_t eta;
691 
692  for(Int_t side=0; side<2; side++) {
693  eta = (side==0)? -e : e; //ta *= (-1);
694 
695  for(Int_t phi=1; phi<=72; phi++) {
696 
697  for(Int_t d=1; d<5; d++) {
698 
699  if (!topo_->valid(HcalDetId(HcalSubdetector(sd), eta, phi, d))) continue;
700  HcalDetId id(HcalSubdetector(sd), eta, phi, d);
701  Float_t corrFactor = 1.0;
702 
703 
704  if (abs(eta)>=CALIB_ABS_IETA_MIN && abs(eta)<=CALIB_ABS_IETA_MAX && HcalSubdetector(sd)!=HcalOuter) {
705  // if (abs(eta)>=CALIB_ABS_IETA_MIN && abs(eta)<=22 && HcalSubdetector(sd)!=HcalOuter) {
706 
707 
708  // need some care when depths were summed for iEta=16 =>
709  // the coeficients are saved in depth 1 of HB: affects
710  Int_t subdetInd = sd;
711 
712  if (abs(eta)==16 && HcalSubdetector(sd)==HcalEndcap && SUM_DEPTHS) {
713  subdetInd = 1;
714  }
715 
716  if (CALIB_METHOD=="L3" || CALIB_METHOD=="L3_AND_MTRX_INV") {
717 
718 
719  if (SUM_DEPTHS && COMBINE_PHI) corrFactor = solution[HcalDetId(HcalSubdetector(subdetInd), eta, 1, 1)];
720  else if (SUM_DEPTHS) corrFactor = solution[HcalDetId(HcalSubdetector(subdetInd), eta, phi, 1)];
721  else if (COMBINE_PHI) corrFactor = solution[HcalDetId(HcalSubdetector(sd), eta, 1, d)];
722  else corrFactor = solution[HcalDetId(HcalSubdetector(sd), eta, phi, d)];
723 
724 
725  // Remark: a new case was added (sumSmallDepths) where the first two depths in towers 15,16
726  // are summed and stored in depth 1.
727  // For now we create the correction coef for depth 2 (set equal to depth 1)
728  // after the call to the L3 minimizer so that this case is also handled without modifying the
729  // logic above. Probably it is better to move it here?
730 
731 
732  }// L3
733 
734  else if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
735  corrFactor = iEtaCoefMap[eta];
736  }
737 
738  else if (CALIB_METHOD=="ISO_TRK_PHI_SYM") {
739  corrFactor = solution[HcalDetId(HcalSubdetector(sd), eta, phi, d)];
740  }
741 
742 
743  } // if within the calibration range
744 
745  fprintf(constFile, "%17i%16i%16i%16s%9.5f%11X\n",
746  eta, phi, d, sdName[sd].Data(), corrFactor, id.rawId());
747 
748  }
749  }
750  }
751  }
752  }
753 
754  return;
755 }
756 
void filterCellsInCone(std::vector< TCell > &selectCells, const GlobalPoint hitPositionHcal, Float_t maxConeDist, const CaloGeometry *theCaloGeometry)
TH1F * h1_corRespBarrel
Definition: hcalCalib.cc:56
TFile * histoFile
Definition: hcalCalib.cc:44
TH1F * h1_corRespIEta[48]
Definition: hcalCalib.cc:66
int i
Definition: DBlmapReader.cc:9
virtual Bool_t Process(Long64_t entry)
Definition: hcalCalib.cc:142
TH1F * h1_rawRespBarrel
Definition: hcalCalib.cc:55
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:45
tuple sumSmallDepths
Definition: diJetCalib.py:58
void GetCoefFromMtrxInvOfAve()
Definition: hcalCalib.cc:521
TH2F * h2_dHitRefEndcap
Definition: hcalCalib.cc:62
TH2F * h2_dHitRefBarrel
Definition: hcalCalib.cc:61
TH1F * h1_allTrkP
Definition: hcalCalib.cc:48
Float_t e()
Definition: TCell.h:26
Definition: TCell.h:15
TH1F * h1_numEventsTwrIEta
Definition: hcalCalib.cc:59
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
TH1F * h1_rawSumE
Definition: hcalCalib.cc:52
virtual void Begin(TTree *tree)
Definition: hcalCalib.cc:69
U second(std::pair< T, U > const &p)
tuple d
Definition: ztail.py:151
T x() const
Cartesian x coordinate.
TH1F * h1_selTrkP_iEta10
Definition: hcalCalib.cc:50
int depth() const
get the tower depth
Definition: HcalDetId.h:55
UInt_t id()
Definition: TCell.h:27
int ieta() const
get the cell ieta
Definition: HcalDetId.h:51
tuple combinePhi
Definition: diJetCalib.py:54
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
void SetE(Float_t e)
Definition: TCell.h:29
void filterCells3x3(std::vector< TCell > &selectCells, Int_t iEta, UInt_t iPhi)
#define M_PI
tuple sumDepths
Definition: diJetCalib.py:51
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:49
TH1F * h1_corRespEndcap
Definition: hcalCalib.cc:58
TH1F * h1_rawResp
Definition: hcalCalib.cc:53
int iphi() const
get the cell iphi
Definition: HcalDetId.h:53
TH1F * h1_rawRespEndcap
Definition: hcalCalib.cc:57
TH1F * h1_trkP
Definition: hcalCalib.cc:47
virtual void Terminate()
Definition: hcalCalib.cc:351
TH1F * h1_corResp
Definition: hcalCalib.cc:54
double sd
double b
Definition: hdecay.h:120
Geom::Phi< T > phi() const
void makeTextFile()
Definition: hcalCalib.cc:673
Bool_t ReadPhiSymCor()
Definition: hcalCalib.cc:614
#define begin
Definition: vmac.h:30
void filterCells5x5(std::vector< TCell > &selectCells, Int_t iEta, UInt_t iPhi)
tuple cout
Definition: gather_cfg.py:121
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)
Definition: MinL3AlgoUniv.h:81
int weight
Definition: histoStyle.py:50
UInt_t nEvents
Definition: hcalCalib.cc:42
void getIEtaIPhiForHighestE(std::vector< TCell > &selectCells, Int_t &iEta, UInt_t &iPhi)
tuple size
Write out results.