CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes
hcalCalib Class Reference

#include <hcalCalib.h>

Inheritance diagram for hcalCalib:

Public Member Functions

void Begin (TTree *tree) override
 
void GetCoefFromMtrxInvOfAve ()
 
Int_t GetEntry (Long64_t entry, Int_t getall=0) override
 
TList * GetOutputList () const override
 
 hcalCalib (TTree *=0)
 
void Init (TTree *tree) override
 
void makeTextFile ()
 
Bool_t Notify () override
 
Bool_t Process (Long64_t entry) override
 
Bool_t ReadPhiSymCor ()
 
void SetApplyPhiSymCorFlag (Bool_t b)
 
void SetCalibAbsIEtaMax (Int_t i)
 
void SetCalibAbsIEtaMin (Int_t i)
 
void SetCalibMethod (const TString &s)
 
void SetCalibType (const TString &s)
 
void SetCaloGeometry (const CaloGeometry *g, const HcalTopology *topo)
 
void SetCombinePhiFlag (Bool_t b)
 
void SetConeMaxDist (Float_t d)
 
void SetHbClusterSize (Int_t i)
 
void SetHeClusterSize (Int_t i)
 
void SetHistoFileName (const TString &filename)
 
void SetInputList (TList *input) override
 
void SetMaxEOverP (Float_t e)
 
void SetMaxEtThirdJet (Float_t et)
 
void SetMaxProbeJetEmFrac (Float_t f)
 
void SetMaxTagJetAbsEta (Float_t e)
 
void SetMaxTagJetEmFrac (Float_t f)
 
void SetMaxTargetE (Float_t e)
 
void SetMaxTrkEmE (Float_t e)
 
void SetMinCellE (Float_t e)
 
void SetMinDPhiDiJets (Float_t dphi)
 
void SetMinEOverP (Float_t e)
 
void SetMinProbeJetAbsEta (Float_t e)
 
void SetMinTagJetEt (Float_t e)
 
void SetMinTargetE (Float_t e)
 
void SetObject (TObject *obj) override
 
void SetOption (const char *option) override
 
void SetOutputCorCoefFileName (const TString &filename)
 
void SetPhiSymCorFileName (const TString &filename)
 
void SetSumDepthsFlag (Bool_t b)
 
void SetSumSmallDepthsFlag (Bool_t b)
 
void SetUseConeClustering (Bool_t b)
 
void Terminate () override
 
Int_t Version () const override
 
 ~hcalCalib () override
 

Public Attributes

Bool_t APPLY_PHI_SYM_COR_FLAG
 
TBranch * b_cells
 
TBranch * b_emEnergy
 
TBranch * b_etVetoJet
 
TBranch * b_eventNumber
 
TBranch * b_iEtaHit
 
TBranch * b_iPhiHit
 
TBranch * b_probeJetEmFrac
 
TBranch * b_probeJetP4
 
TBranch * b_runNumber
 
TBranch * b_tagJetEmFrac
 
TBranch * b_tagJetP4
 
TBranch * b_targetE
 
TBranch * b_xTrkEcal
 
TBranch * b_xTrkHcal
 
TBranch * b_yTrkEcal
 
TBranch * b_yTrkHcal
 
TBranch * b_zTrkEcal
 
TBranch * b_zTrkHcal
 
Int_t CALIB_ABS_IETA_MAX
 
Int_t CALIB_ABS_IETA_MIN
 
TString CALIB_METHOD
 
TString CALIB_TYPE
 
std::vector< std::vector< Float_t > > cellEnergies
 
std::vector< std::vector< UInt_t > > cellIds
 
TClonesArray * cells
 
Bool_t COMBINE_PHI
 
Float_t emEnergy
 
Float_t etVetoJet
 
UInt_t eventNumber
 pointer to the analyzed TTree or TChain More...
 
TTree * fChain
 
Int_t HB_CLUSTER_SIZE
 
Int_t HE_CLUSTER_SIZE
 
TString HISTO_FILENAME
 
std::map< Int_t, Float_t > iEtaCoefMap
 
Int_t iEtaHit
 
UInt_t iPhiHit
 
Float_t MAX_CONE_DIST
 
Float_t MAX_EOVERP
 
Float_t MAX_ET_THIRD_JET
 
Float_t MAX_PROBEJET_EMFRAC
 
Float_t MAX_TAGJET_ABSETA
 
Float_t MAX_TAGJET_EMFRAC
 
Float_t MAX_TARGET_E
 
Float_t MAX_TRK_EME
 
Float_t MIN_CELL_E
 
Float_t MIN_DPHI_DIJETS
 
Float_t MIN_EOVERP
 
Float_t MIN_PROBEJET_ABSETA
 
Float_t MIN_TAGJET_ET
 
Float_t MIN_TARGET_E
 
TString OUTPUT_COR_COEF_FILENAME
 
TString PHI_SYM_COR_FILENAME
 
std::map< UInt_t, Float_t > phiSymCor
 
Float_t probeJetEmFrac
 
TLorentzVector * probeJetP4
 
std::vector< std::pair< Int_t, UInt_t > > refIEtaIPhi
 
UInt_t runNumber
 
std::map< UInt_t, Float_t > solution
 
Bool_t SUM_DEPTHS
 
Bool_t SUM_SMALL_DEPTHS
 
Float_t tagJetEmFrac
 
TLorentzVector * tagJetP4
 
Float_t targetE
 
std::vector< Float_t > targetEnergies
 
const CaloGeometrytheCaloGeometry
 
const HcalTopologytopo_
 
Bool_t USE_CONE_CLUSTERING
 
Float_t xTrkEcal
 
Float_t xTrkHcal
 
Float_t yTrkEcal
 
Float_t yTrkHcal
 
Float_t zTrkEcal
 
Float_t zTrkHcal
 

Detailed Description

Definition at line 39 of file hcalCalib.h.

Constructor & Destructor Documentation

hcalCalib::hcalCalib ( TTree *  = 0)
inline

Definition at line 89 of file hcalCalib.h.

89 { }
hcalCalib::~hcalCalib ( )
inlineoverride

Definition at line 90 of file hcalCalib.h.

90 { }

Member Function Documentation

void hcalCalib::Begin ( TTree *  tree)
override

Definition at line 69 of file hcalCalib.cc.

References gather_cfg::cout, cmsRelvalreport::exit, h1_allTrkP, h1_corResp, h1_corRespBarrel, h1_corRespEndcap, h1_corRespIEta, h1_numEventsTwrIEta, h1_rawResp, h1_rawRespBarrel, h1_rawRespEndcap, h1_rawSumE, h1_selTrkP_iEta10, h1_trkP, h2_dHitRefBarrel, h2_dHitRefEndcap, histoFile, mps_fire::i, nEvents, and TSGForRoadSearch_cfi::option.

Referenced by Version().

69  {
70  TString option = GetOption();
71 
72  nEvents = 0;
73 
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()
TH1F * h1_corRespBarrel
Definition: hcalCalib.cc:56
TFile * histoFile
Definition: hcalCalib.cc:44
TH1F * h1_corRespIEta[48]
Definition: hcalCalib.cc:66
TH1F * h1_rawRespBarrel
Definition: hcalCalib.cc:55
TH2F * h2_dHitRefEndcap
Definition: hcalCalib.cc:62
TH2F * h2_dHitRefBarrel
Definition: hcalCalib.cc:61
TString CALIB_TYPE
Definition: hcalCalib.h:137
TH1F * h1_allTrkP
Definition: hcalCalib.cc:48
TH1F * h1_numEventsTwrIEta
Definition: hcalCalib.cc:59
TH1F * h1_rawSumE
Definition: hcalCalib.cc:52
TH1F * h1_selTrkP_iEta10
Definition: hcalCalib.cc:50
Bool_t APPLY_PHI_SYM_COR_FLAG
Definition: hcalCalib.h:141
TH1F * h1_corRespEndcap
Definition: hcalCalib.cc:58
TH1F * h1_rawResp
Definition: hcalCalib.cc:53
TH1F * h1_rawRespEndcap
Definition: hcalCalib.cc:57
TH1F * h1_trkP
Definition: hcalCalib.cc:47
TH1F * h1_corResp
Definition: hcalCalib.cc:54
Bool_t ReadPhiSymCor()
Definition: hcalCalib.cc:614
TString HISTO_FILENAME
Definition: hcalCalib.h:144
UInt_t nEvents
Definition: hcalCalib.cc:42
void hcalCalib::GetCoefFromMtrxInvOfAve ( )

Definition at line 521 of file hcalCalib.cc.

References patCaloMETCorrections_cff::A, funct::abs(), b, begin, end, spr::find(), mps_fire::i, HcalDetId::ieta(), and findQualityFiles::size.

Referenced by SetCaloGeometry().

521  {
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 }
size
Write out results.
std::map< Int_t, Float_t > iEtaCoefMap
Definition: hcalCalib.h:201
std::vector< std::pair< Int_t, UInt_t > > refIEtaIPhi
Definition: hcalCalib.h:195
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
Int_t CALIB_ABS_IETA_MIN
Definition: hcalCalib.h:128
std::vector< Float_t > targetEnergies
Definition: hcalCalib.h:196
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define end
Definition: vmac.h:39
std::map< UInt_t, Float_t > solution
Definition: hcalCalib.h:200
TString CALIB_METHOD
Definition: hcalCalib.h:138
double b
Definition: hdecay.h:120
Int_t CALIB_ABS_IETA_MAX
Definition: hcalCalib.h:127
#define begin
Definition: vmac.h:32
std::vector< std::vector< Float_t > > cellEnergies
Definition: hcalCalib.h:193
std::vector< std::vector< UInt_t > > cellIds
Definition: hcalCalib.h:194
Int_t hcalCalib::GetEntry ( Long64_t  entry,
Int_t  getall = 0 
)
inlineoverride

Definition at line 97 of file hcalCalib.h.

97 { return fChain ? fChain->GetTree()->GetEntry(entry, getall) : 0; }
TTree * fChain
Definition: hcalCalib.h:41
TList* hcalCalib::GetOutputList ( ) const
inlineoverride

Definition at line 101 of file hcalCalib.h.

References Terminate().

101 { return fOutput; }
void hcalCalib::Init ( TTree *  tree)
override

Referenced by Version().

void hcalCalib::makeTextFile ( )

Definition at line 673 of file hcalCalib.cc.

References funct::abs(), edmIntegrityCheck::d, MillePedeFileConverter_cfg::e, PVValHelper::eta, HcalEndcap, HcalOuter, triggerObjects_cff::id, and sd.

Referenced by SetCaloGeometry().

673  {
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 
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 }
std::map< Int_t, Float_t > iEtaCoefMap
Definition: hcalCalib.h:201
bool valid(const DetId &id) const override
Bool_t SUM_DEPTHS
Definition: hcalCalib.h:117
Int_t CALIB_ABS_IETA_MIN
Definition: hcalCalib.h:128
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const HcalTopology * topo_
Definition: hcalCalib.h:148
std::map< UInt_t, Float_t > solution
Definition: hcalCalib.h:200
Bool_t COMBINE_PHI
Definition: hcalCalib.h:119
TString CALIB_METHOD
Definition: hcalCalib.h:138
TString OUTPUT_COR_COEF_FILENAME
Definition: hcalCalib.h:143
double sd
Int_t CALIB_ABS_IETA_MAX
Definition: hcalCalib.h:127
Bool_t hcalCalib::Notify ( )
override

Referenced by Version().

Bool_t hcalCalib::Process ( Long64_t  entry)
override

Definition at line 142 of file hcalCalib.cc.

References funct::abs(), combinePhi(), gather_cfg::cout, TCell::e(), filterCells3x3(), filterCells5x5(), filterCellsInCone(), getIEtaIPhiForHighestE(), h1_allTrkP, h1_rawSumE, h1_selTrkP_iEta10, h2_dHitRefBarrel, h2_dHitRefEndcap, HcalBarrel, HcalEndcap, HcalForward, HcalOuter, mps_fire::i, TCell::id(), electrons_cff::ids, M_PI, nEvents, edm::second(), TCell::SetE(), sumDepths(), and sumSmallDepths().

Referenced by Version().

142  {
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.empty()) {
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 }
void filterCellsInCone(std::vector< TCell > &selectCells, const GlobalPoint hitPositionHcal, Float_t maxConeDist, const CaloGeometry *theCaloGeometry)
Float_t MAX_TARGET_E
Definition: hcalCalib.h:107
TH2F * h2_dHitRefEndcap
Definition: hcalCalib.cc:62
TH2F * h2_dHitRefBarrel
Definition: hcalCalib.cc:61
Float_t xTrkHcal
Definition: hcalCalib.h:52
TString CALIB_TYPE
Definition: hcalCalib.h:137
TLorentzVector * tagJetP4
Definition: hcalCalib.h:59
TH1F * h1_allTrkP
Definition: hcalCalib.cc:48
Float_t zTrkHcal
Definition: hcalCalib.h:54
std::map< UInt_t, Float_t > phiSymCor
Definition: hcalCalib.h:198
Float_t e()
Definition: TCell.h:26
Float_t MAX_EOVERP
Definition: hcalCalib.h:111
Float_t targetE
Definition: hcalCalib.h:49
Bool_t SUM_DEPTHS
Definition: hcalCalib.h:117
Definition: TCell.h:15
Float_t MIN_CELL_E
Definition: hcalCalib.h:109
TLorentzVector * probeJetP4
Definition: hcalCalib.h:60
std::vector< std::pair< Int_t, UInt_t > > refIEtaIPhi
Definition: hcalCalib.h:195
TH1F * h1_rawSumE
Definition: hcalCalib.cc:52
std::vector< Float_t > targetEnergies
Definition: hcalCalib.h:196
U second(std::pair< T, U > const &p)
void combinePhi(std::vector< TCell > &selectCells)
TH1F * h1_selTrkP_iEta10
Definition: hcalCalib.cc:50
Int_t HE_CLUSTER_SIZE
Definition: hcalCalib.h:122
UInt_t id()
Definition: TCell.h:27
Int_t iEtaHit
Definition: hcalCalib.h:45
Float_t MAX_TRK_EME
Definition: hcalCalib.h:112
Float_t MIN_DPHI_DIJETS
Definition: hcalCalib.h:115
Int_t HB_CLUSTER_SIZE
Definition: hcalCalib.h:121
Float_t MIN_TAGJET_ET
Definition: hcalCalib.h:133
Bool_t APPLY_PHI_SYM_COR_FLAG
Definition: hcalCalib.h:141
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Float_t emEnergy
Definition: hcalCalib.h:48
Float_t MIN_EOVERP
Definition: hcalCalib.h:110
void SetE(Float_t e)
Definition: TCell.h:29
void filterCells3x3(std::vector< TCell > &selectCells, Int_t iEta, UInt_t iPhi)
Float_t MIN_PROBEJET_ABSETA
Definition: hcalCalib.h:135
#define M_PI
Bool_t COMBINE_PHI
Definition: hcalCalib.h:119
Float_t MAX_ET_THIRD_JET
Definition: hcalCalib.h:114
UInt_t iPhiHit
Definition: hcalCalib.h:46
Bool_t USE_CONE_CLUSTERING
Definition: hcalCalib.h:124
void sumSmallDepths(std::vector< TCell > &selectCells)
Float_t etVetoJet
Definition: hcalCalib.h:50
TClonesArray * cells
Definition: hcalCalib.h:47
void filterCells5x5(std::vector< TCell > &selectCells, Int_t iEta, UInt_t iPhi)
Bool_t SUM_SMALL_DEPTHS
Definition: hcalCalib.h:118
Int_t GetEntry(Long64_t entry, Int_t getall=0) override
Definition: hcalCalib.h:97
Float_t MAX_TAGJET_ABSETA
Definition: hcalCalib.h:132
Float_t MAX_PROBEJET_EMFRAC
Definition: hcalCalib.h:130
const CaloGeometry * theCaloGeometry
Definition: hcalCalib.h:147
void sumDepths(std::vector< TCell > &selectCells)
std::vector< std::vector< Float_t > > cellEnergies
Definition: hcalCalib.h:193
UInt_t nEvents
Definition: hcalCalib.cc:42
std::vector< std::vector< UInt_t > > cellIds
Definition: hcalCalib.h:194
void getIEtaIPhiForHighestE(std::vector< TCell > &selectCells, Int_t &iEta, UInt_t &iPhi)
Float_t MAX_CONE_DIST
Definition: hcalCalib.h:125
Float_t yTrkHcal
Definition: hcalCalib.h:53
Float_t probeJetEmFrac
Definition: hcalCalib.h:63
Bool_t hcalCalib::ReadPhiSymCor ( )

Definition at line 614 of file hcalCalib.cc.

References gather_cfg::cout, particleFlowClusterECALTimeSelected_cfi::depth, HcalBarrel, HcalEndcap, HcalForward, HcalOuter, geometryCSVtoXML::line, sd, AlCaHLTBitMon_QueryRunRegistry::string, and relativeConstraints::value.

Referenced by SetCaloGeometry().

614  {
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.empty() || 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 }
TString PHI_SYM_COR_FILENAME
Definition: hcalCalib.h:140
bool valid(const DetId &id) const override
std::map< UInt_t, Float_t > phiSymCor
Definition: hcalCalib.h:198
HcalSubdetector
Definition: HcalAssistant.h:31
const HcalTopology * topo_
Definition: hcalCalib.h:148
double sd
void hcalCalib::SetApplyPhiSymCorFlag ( Bool_t  b)
inline

Definition at line 172 of file hcalCalib.h.

References b.

Referenced by HcalCalibrator::endJob().

Bool_t APPLY_PHI_SYM_COR_FLAG
Definition: hcalCalib.h:141
double b
Definition: hdecay.h:120
void hcalCalib::SetCalibAbsIEtaMax ( Int_t  i)
inline

Definition at line 168 of file hcalCalib.h.

References mps_fire::i.

Referenced by HcalCalibrator::endJob().

168 { CALIB_ABS_IETA_MAX = i; }
Int_t CALIB_ABS_IETA_MAX
Definition: hcalCalib.h:127
void hcalCalib::SetCalibAbsIEtaMin ( Int_t  i)
inline

Definition at line 169 of file hcalCalib.h.

References mps_fire::i.

Referenced by HcalCalibrator::endJob().

169 { CALIB_ABS_IETA_MIN = i; }
Int_t CALIB_ABS_IETA_MIN
Definition: hcalCalib.h:128
void hcalCalib::SetCalibMethod ( const TString &  s)
inline

Definition at line 161 of file hcalCalib.h.

References alignCSCRings::s.

Referenced by HcalCalibrator::endJob().

161 { CALIB_METHOD = s; }
TString CALIB_METHOD
Definition: hcalCalib.h:138
void hcalCalib::SetCalibType ( const TString &  s)
inline

Definition at line 160 of file hcalCalib.h.

References alignCSCRings::s.

Referenced by HcalCalibrator::endJob().

160 { CALIB_TYPE = s; }
TString CALIB_TYPE
Definition: hcalCalib.h:137
void hcalCalib::SetCaloGeometry ( const CaloGeometry g,
const HcalTopology topo 
)
inline

Definition at line 183 of file hcalCalib.h.

References g, GetCoefFromMtrxInvOfAve(), makeTextFile(), and ReadPhiSymCor().

Referenced by HcalCalibrator::endJob().

183 { theCaloGeometry = g; topo_=topo; }
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
const HcalTopology * topo_
Definition: hcalCalib.h:148
const CaloGeometry * theCaloGeometry
Definition: hcalCalib.h:147
void hcalCalib::SetCombinePhiFlag ( Bool_t  b)
inline

Definition at line 155 of file hcalCalib.h.

References b.

Referenced by HcalCalibrator::endJob().

155 { COMBINE_PHI = b; }
Bool_t COMBINE_PHI
Definition: hcalCalib.h:119
double b
Definition: hdecay.h:120
void hcalCalib::SetConeMaxDist ( Float_t  d)
inline

Definition at line 166 of file hcalCalib.h.

References edmIntegrityCheck::d.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetHbClusterSize ( Int_t  i)
inline

Definition at line 162 of file hcalCalib.h.

References mps_fire::i.

Referenced by HcalCalibrator::endJob().

162 { HB_CLUSTER_SIZE = i; }
Int_t HB_CLUSTER_SIZE
Definition: hcalCalib.h:121
void hcalCalib::SetHeClusterSize ( Int_t  i)
inline

Definition at line 163 of file hcalCalib.h.

References mps_fire::i.

Referenced by HcalCalibrator::endJob().

163 { HE_CLUSTER_SIZE = i; }
Int_t HE_CLUSTER_SIZE
Definition: hcalCalib.h:122
void hcalCalib::SetHistoFileName ( const TString &  filename)
inline

Definition at line 180 of file hcalCalib.h.

References corrVsCorr::filename.

Referenced by HcalCalibrator::endJob().

TString HISTO_FILENAME
Definition: hcalCalib.h:144
void hcalCalib::SetInputList ( TList *  input)
inlineoverride

Definition at line 100 of file hcalCalib.h.

References input.

100 { fInput = input; }
static std::string const input
Definition: EdmProvDump.cc:44
void hcalCalib::SetMaxEOverP ( Float_t  e)
inline

Definition at line 158 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetMaxEtThirdJet ( Float_t  et)
inline

Definition at line 170 of file hcalCalib.h.

References stringResolutionProvider_cfi::et.

Referenced by HcalCalibrator::endJob().

170 { MAX_ET_THIRD_JET = et; }
Float_t MAX_ET_THIRD_JET
Definition: hcalCalib.h:114
et
define resolution functions of each parameter
void hcalCalib::SetMaxProbeJetEmFrac ( Float_t  f)
inline

Definition at line 174 of file hcalCalib.h.

References f.

Referenced by HcalCalibrator::endJob().

174 { MAX_PROBEJET_EMFRAC = f; }
double f[11][100]
Float_t MAX_PROBEJET_EMFRAC
Definition: hcalCalib.h:130
void hcalCalib::SetMaxTagJetAbsEta ( Float_t  e)
inline
void hcalCalib::SetMaxTagJetEmFrac ( Float_t  f)
inline

Definition at line 175 of file hcalCalib.h.

References f.

Referenced by HcalCalibrator::endJob().

175 { MAX_TAGJET_EMFRAC = f; }
double f[11][100]
Float_t MAX_TAGJET_EMFRAC
Definition: hcalCalib.h:131
void hcalCalib::SetMaxTargetE ( Float_t  e)
inline

Definition at line 152 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetMaxTrkEmE ( Float_t  e)
inline

Definition at line 159 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetMinCellE ( Float_t  e)
inline

Definition at line 156 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetMinDPhiDiJets ( Float_t  dphi)
inline

Definition at line 171 of file hcalCalib.h.

Referenced by HcalCalibrator::endJob().

171 { MIN_DPHI_DIJETS = dphi; }
Float_t MIN_DPHI_DIJETS
Definition: hcalCalib.h:115
void hcalCalib::SetMinEOverP ( Float_t  e)
inline

Definition at line 157 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetMinProbeJetAbsEta ( Float_t  e)
inline
void hcalCalib::SetMinTagJetEt ( Float_t  e)
inline

Definition at line 177 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetMinTargetE ( Float_t  e)
inline

Definition at line 151 of file hcalCalib.h.

References MillePedeFileConverter_cfg::e.

Referenced by HcalCalibrator::endJob().

void hcalCalib::SetObject ( TObject *  obj)
inlineoverride

Definition at line 99 of file hcalCalib.h.

References GetRecoTauVFromDQM_MC_cff::obj.

void hcalCalib::SetOption ( const char *  option)
inlineoverride

Definition at line 98 of file hcalCalib.h.

References TSGForRoadSearch_cfi::option.

98 { fOption = option; }
void hcalCalib::SetOutputCorCoefFileName ( const TString &  filename)
inline

Definition at line 179 of file hcalCalib.h.

References corrVsCorr::filename.

Referenced by HcalCalibrator::endJob().

TString OUTPUT_COR_COEF_FILENAME
Definition: hcalCalib.h:143
void hcalCalib::SetPhiSymCorFileName ( const TString &  filename)
inline

Definition at line 173 of file hcalCalib.h.

References corrVsCorr::filename.

Referenced by HcalCalibrator::endJob().

TString PHI_SYM_COR_FILENAME
Definition: hcalCalib.h:140
void hcalCalib::SetSumDepthsFlag ( Bool_t  b)
inline

Definition at line 153 of file hcalCalib.h.

References b.

Referenced by HcalCalibrator::endJob().

153 { SUM_DEPTHS = b; }
Bool_t SUM_DEPTHS
Definition: hcalCalib.h:117
double b
Definition: hdecay.h:120
void hcalCalib::SetSumSmallDepthsFlag ( Bool_t  b)
inline

Definition at line 154 of file hcalCalib.h.

References b.

Referenced by HcalCalibrator::endJob().

154 { SUM_SMALL_DEPTHS = b; }
double b
Definition: hdecay.h:120
Bool_t SUM_SMALL_DEPTHS
Definition: hcalCalib.h:118
void hcalCalib::SetUseConeClustering ( Bool_t  b)
inline

Definition at line 165 of file hcalCalib.h.

References b.

Referenced by HcalCalibrator::endJob().

165 { USE_CONE_CLUSTERING = b; }
Bool_t USE_CONE_CLUSTERING
Definition: hcalCalib.h:124
double b
Definition: hdecay.h:120
void hcalCalib::Terminate ( )
override

Definition at line 351 of file hcalCalib.cc.

References funct::abs(), gather_cfg::cout, HcalDetId::depth(), h1_allTrkP, h1_corResp, h1_corRespBarrel, h1_corRespEndcap, h1_corRespIEta, h1_numEventsTwrIEta, h1_rawResp, h1_rawRespBarrel, h1_rawRespEndcap, h1_rawSumE, h1_selTrkP_iEta10, h1_trkP, h2_dHitRefBarrel, h2_dHitRefEndcap, HcalBarrel, histoFile, mps_fire::i, HcalDetId::ieta(), HcalDetId::ietaAbs(), HcalDetId::iphi(), MinL3AlgoUniv< IDdet >::iterate(), findQualityFiles::size, and HcalDetId::subdet().

Referenced by GetOutputList().

351  {
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") {
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") {
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()
size
Write out results.
TH1F * h1_corRespBarrel
Definition: hcalCalib.cc:56
TFile * histoFile
Definition: hcalCalib.cc:44
TH1F * h1_corRespIEta[48]
Definition: hcalCalib.cc:66
std::map< Int_t, Float_t > iEtaCoefMap
Definition: hcalCalib.h:201
TH1F * h1_rawRespBarrel
Definition: hcalCalib.cc:55
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
void GetCoefFromMtrxInvOfAve()
Definition: hcalCalib.cc:521
TH2F * h2_dHitRefEndcap
Definition: hcalCalib.cc:62
TH2F * h2_dHitRefBarrel
Definition: hcalCalib.cc:61
TString CALIB_TYPE
Definition: hcalCalib.h:137
TH1F * h1_allTrkP
Definition: hcalCalib.cc:48
Definition: weight.py:1
Bool_t SUM_DEPTHS
Definition: hcalCalib.h:117
TH1F * h1_numEventsTwrIEta
Definition: hcalCalib.cc:59
std::vector< std::pair< Int_t, UInt_t > > refIEtaIPhi
Definition: hcalCalib.h:195
TH1F * h1_rawSumE
Definition: hcalCalib.cc:52
Int_t CALIB_ABS_IETA_MIN
Definition: hcalCalib.h:128
std::vector< Float_t > targetEnergies
Definition: hcalCalib.h:196
TH1F * h1_selTrkP_iEta10
Definition: hcalCalib.cc:50
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< UInt_t, Float_t > solution
Definition: hcalCalib.h:200
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.cc:98
TH1F * h1_corRespEndcap
Definition: hcalCalib.cc:58
TH1F * h1_rawResp
Definition: hcalCalib.cc:53
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
TString CALIB_METHOD
Definition: hcalCalib.h:138
TH1F * h1_rawRespEndcap
Definition: hcalCalib.cc:57
TH1F * h1_trkP
Definition: hcalCalib.cc:47
TH1F * h1_corResp
Definition: hcalCalib.cc:54
Int_t CALIB_ABS_IETA_MAX
Definition: hcalCalib.h:127
void makeTextFile()
Definition: hcalCalib.cc:673
Bool_t SUM_SMALL_DEPTHS
Definition: hcalCalib.h:118
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
std::vector< std::vector< Float_t > > cellEnergies
Definition: hcalCalib.h:193
std::vector< std::vector< UInt_t > > cellIds
Definition: hcalCalib.h:194
Int_t hcalCalib::Version ( ) const
inlineoverride

Definition at line 91 of file hcalCalib.h.

References Begin(), mps_splice::entry, Init(), Notify(), and Process().

91 { return 2; }

Member Data Documentation

Bool_t hcalCalib::APPLY_PHI_SYM_COR_FLAG

Definition at line 141 of file hcalCalib.h.

TBranch* hcalCalib::b_cells

Definition at line 70 of file hcalCalib.h.

TBranch* hcalCalib::b_emEnergy

Definition at line 71 of file hcalCalib.h.

TBranch* hcalCalib::b_etVetoJet

Definition at line 73 of file hcalCalib.h.

TBranch* hcalCalib::b_eventNumber

Definition at line 66 of file hcalCalib.h.

TBranch* hcalCalib::b_iEtaHit

Definition at line 68 of file hcalCalib.h.

TBranch* hcalCalib::b_iPhiHit

Definition at line 69 of file hcalCalib.h.

TBranch* hcalCalib::b_probeJetEmFrac

Definition at line 83 of file hcalCalib.h.

TBranch* hcalCalib::b_probeJetP4

Definition at line 86 of file hcalCalib.h.

TBranch* hcalCalib::b_runNumber

Definition at line 67 of file hcalCalib.h.

TBranch* hcalCalib::b_tagJetEmFrac

Definition at line 82 of file hcalCalib.h.

TBranch* hcalCalib::b_tagJetP4

Definition at line 85 of file hcalCalib.h.

TBranch* hcalCalib::b_targetE

Definition at line 72 of file hcalCalib.h.

TBranch* hcalCalib::b_xTrkEcal

Definition at line 78 of file hcalCalib.h.

TBranch* hcalCalib::b_xTrkHcal

Definition at line 75 of file hcalCalib.h.

TBranch* hcalCalib::b_yTrkEcal

Definition at line 79 of file hcalCalib.h.

TBranch* hcalCalib::b_yTrkHcal

Definition at line 76 of file hcalCalib.h.

TBranch* hcalCalib::b_zTrkEcal

Definition at line 80 of file hcalCalib.h.

TBranch* hcalCalib::b_zTrkHcal

Definition at line 77 of file hcalCalib.h.

Int_t hcalCalib::CALIB_ABS_IETA_MAX

Definition at line 127 of file hcalCalib.h.

Int_t hcalCalib::CALIB_ABS_IETA_MIN

Definition at line 128 of file hcalCalib.h.

TString hcalCalib::CALIB_METHOD

Definition at line 138 of file hcalCalib.h.

TString hcalCalib::CALIB_TYPE

Definition at line 137 of file hcalCalib.h.

std::vector< std::vector<Float_t> > hcalCalib::cellEnergies

Definition at line 193 of file hcalCalib.h.

std::vector< std::vector<UInt_t> > hcalCalib::cellIds

Definition at line 194 of file hcalCalib.h.

TClonesArray* hcalCalib::cells

Definition at line 47 of file hcalCalib.h.

Bool_t hcalCalib::COMBINE_PHI

Definition at line 119 of file hcalCalib.h.

Float_t hcalCalib::emEnergy

Definition at line 48 of file hcalCalib.h.

Float_t hcalCalib::etVetoJet

Definition at line 50 of file hcalCalib.h.

UInt_t hcalCalib::eventNumber

pointer to the analyzed TTree or TChain

Definition at line 43 of file hcalCalib.h.

Referenced by Vispa.Plugins.EdmBrowser.EdmDataAccessor.EdmDataAccessor::setFilterBranches().

TTree* hcalCalib::fChain

Definition at line 41 of file hcalCalib.h.

Int_t hcalCalib::HB_CLUSTER_SIZE

Definition at line 121 of file hcalCalib.h.

Int_t hcalCalib::HE_CLUSTER_SIZE

Definition at line 122 of file hcalCalib.h.

TString hcalCalib::HISTO_FILENAME

Definition at line 144 of file hcalCalib.h.

std::map<Int_t, Float_t> hcalCalib::iEtaCoefMap

Definition at line 201 of file hcalCalib.h.

Int_t hcalCalib::iEtaHit

Definition at line 45 of file hcalCalib.h.

UInt_t hcalCalib::iPhiHit

Definition at line 46 of file hcalCalib.h.

Float_t hcalCalib::MAX_CONE_DIST

Definition at line 125 of file hcalCalib.h.

Float_t hcalCalib::MAX_EOVERP

Definition at line 111 of file hcalCalib.h.

Float_t hcalCalib::MAX_ET_THIRD_JET

Definition at line 114 of file hcalCalib.h.

Float_t hcalCalib::MAX_PROBEJET_EMFRAC

Definition at line 130 of file hcalCalib.h.

Float_t hcalCalib::MAX_TAGJET_ABSETA

Definition at line 132 of file hcalCalib.h.

Float_t hcalCalib::MAX_TAGJET_EMFRAC

Definition at line 131 of file hcalCalib.h.

Float_t hcalCalib::MAX_TARGET_E

Definition at line 107 of file hcalCalib.h.

Float_t hcalCalib::MAX_TRK_EME

Definition at line 112 of file hcalCalib.h.

Float_t hcalCalib::MIN_CELL_E

Definition at line 109 of file hcalCalib.h.

Float_t hcalCalib::MIN_DPHI_DIJETS

Definition at line 115 of file hcalCalib.h.

Float_t hcalCalib::MIN_EOVERP

Definition at line 110 of file hcalCalib.h.

Float_t hcalCalib::MIN_PROBEJET_ABSETA

Definition at line 135 of file hcalCalib.h.

Float_t hcalCalib::MIN_TAGJET_ET

Definition at line 133 of file hcalCalib.h.

Float_t hcalCalib::MIN_TARGET_E

Definition at line 106 of file hcalCalib.h.

TString hcalCalib::OUTPUT_COR_COEF_FILENAME

Definition at line 143 of file hcalCalib.h.

TString hcalCalib::PHI_SYM_COR_FILENAME

Definition at line 140 of file hcalCalib.h.

std::map<UInt_t, Float_t> hcalCalib::phiSymCor

Definition at line 198 of file hcalCalib.h.

Float_t hcalCalib::probeJetEmFrac

Definition at line 63 of file hcalCalib.h.

TLorentzVector* hcalCalib::probeJetP4

Definition at line 60 of file hcalCalib.h.

std::vector< std::pair<Int_t, UInt_t> > hcalCalib::refIEtaIPhi

Definition at line 195 of file hcalCalib.h.

UInt_t hcalCalib::runNumber

Definition at line 44 of file hcalCalib.h.

Referenced by TH2PolyOfflineMaps.TH2PolyOfflineMaps::PrintTrackerMaps().

std::map<UInt_t, Float_t> hcalCalib::solution

Definition at line 200 of file hcalCalib.h.

Bool_t hcalCalib::SUM_DEPTHS

Definition at line 117 of file hcalCalib.h.

Bool_t hcalCalib::SUM_SMALL_DEPTHS

Definition at line 118 of file hcalCalib.h.

Float_t hcalCalib::tagJetEmFrac

Definition at line 62 of file hcalCalib.h.

TLorentzVector* hcalCalib::tagJetP4

Definition at line 59 of file hcalCalib.h.

Float_t hcalCalib::targetE

Definition at line 49 of file hcalCalib.h.

std::vector< Float_t> hcalCalib::targetEnergies

Definition at line 196 of file hcalCalib.h.

const CaloGeometry* hcalCalib::theCaloGeometry

Definition at line 147 of file hcalCalib.h.

const HcalTopology* hcalCalib::topo_

Definition at line 148 of file hcalCalib.h.

Bool_t hcalCalib::USE_CONE_CLUSTERING

Definition at line 124 of file hcalCalib.h.

Float_t hcalCalib::xTrkEcal

Definition at line 55 of file hcalCalib.h.

Float_t hcalCalib::xTrkHcal

Definition at line 52 of file hcalCalib.h.

Float_t hcalCalib::yTrkEcal

Definition at line 56 of file hcalCalib.h.

Float_t hcalCalib::yTrkHcal

Definition at line 53 of file hcalCalib.h.

Float_t hcalCalib::zTrkEcal

Definition at line 57 of file hcalCalib.h.

Float_t hcalCalib::zTrkHcal

Definition at line 54 of file hcalCalib.h.