CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
RecAnalyzerMinbias Class Reference
Inheritance diagram for RecAnalyzerMinbias:
edm::EDAnalyzer edm::EDConsumerBase

Classes

struct  myInfo
 

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
virtual void beginJob ()
 
virtual void beginRun (edm::Run const &iRun, edm::EventSetup const &iSetup)
 
virtual void endJob ()
 
 RecAnalyzerMinbias (const edm::ParameterSet &)
 
 ~RecAnalyzerMinbias ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void analyzeHcal (const HBHERecHitCollection &, const HFRecHitCollection &, int, bool)
 

Private Attributes

int cells
 
std::map< DetId, double > corrFactor_
 
int depth
 
double eHighHB_
 
double eHighHE_
 
double eHighHF_
 
double eLowHB_
 
double eLowHE_
 
double eLowHF_
 
bool fillHist_
 
std::string fOutputFileName_
 
std::vector< unsigned int > hcalID_
 
std::map< HcalDetId, TH1D * > histHB_
 
std::map< HcalDetId, TH1D * > histHE_
 
std::map< HcalDetId, TH1D * > histHF_
 
std::vector< TH1D * > histo_
 
TFile * hOutputFile_
 
int ieta
 
bool ignoreL1_
 
bool init_
 
int iphi
 
float mom0_MB
 
float mom1_MB
 
float mom2_MB
 
float mom3_MB
 
float mom4_MB
 
std::map< std::pair< int,
HcalDetId >, myInfo
myMap_
 
int mysubd
 
TTree * myTree_
 
double rnnum_
 
double rnnumber
 
bool runNZS_
 
bool theRecalib_
 
edm::EDGetTokenT
< HBHERecHitCollection
tok_hbherecoMB_
 
edm::EDGetTokenT
< HFRecHitCollection
tok_hfrecoMB_
 
edm::EDGetTokenT
< L1GlobalTriggerObjectMapRecord
tok_hltL1GtMap_
 
int trigbit
 
std::vector< int > trigbit_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 34 of file RecAnalyzerMinbias.cc.

Constructor & Destructor Documentation

RecAnalyzerMinbias::RecAnalyzerMinbias ( const edm::ParameterSet iConfig)
explicit

Definition at line 78 of file RecAnalyzerMinbias.cc.

References funct::abs(), corrFactor_, depth, eHighHB_, eHighHE_, eHighHF_, eLowHB_, eLowHE_, eLowHF_, fillHist_, fOutputFileName_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), HcalBarrel, HcalEndcap, HcalForward, hcalID_, HcalOuter, ieta, ignoreL1_, EdgesToViz::infile, HLT_25ns14e33_v1_cff::InputTag, iphi, relval_steps::k, runNZS_, AlCaHLTBitMon_QueryRunRegistry::string, theRecalib_, tok_hbherecoMB_, tok_hfrecoMB_, tok_hltL1GtMap_, and trigbit_.

78  :
79  init_(false) {
80 
81  // get name of output file with histogramms
82  runNZS_ = iConfig.getParameter<bool>("RunNZS");
83  eLowHB_ = iConfig.getParameter<double>("ELowHB");
84  eHighHB_ = iConfig.getParameter<double>("EHighHB");
85  eLowHE_ = iConfig.getParameter<double>("ELowHE");
86  eHighHE_ = iConfig.getParameter<double>("EHighHE");
87  eLowHF_ = iConfig.getParameter<double>("ELowHF");
88  eHighHF_ = iConfig.getParameter<double>("EHighHF");
89  fOutputFileName_ = iConfig.getUntrackedParameter<std::string>("HistOutFile");
90  trigbit_ = iConfig.getUntrackedParameter<std::vector<int>>("TriggerBits");
91  ignoreL1_ = iConfig.getUntrackedParameter<bool>("IgnoreL1",false);
92  std::string cfile= iConfig.getUntrackedParameter<std::string>("CorrFile");
93  fillHist_ = iConfig.getUntrackedParameter<bool>("FillHisto",false);
94  std::vector<int> ieta = iConfig.getUntrackedParameter<std::vector<int>>("HcalIeta");
95  std::vector<int> iphi = iConfig.getUntrackedParameter<std::vector<int>>("HcalIphi");
96  std::vector<int> depth= iConfig.getUntrackedParameter<std::vector<int>>("HcalDepth");
97 
98  // get token names of modules, producing object collections
99  tok_hbherecoMB_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInputMB"));
100  tok_hfrecoMB_ = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInputMB"));
101  tok_hltL1GtMap_ = consumes<L1GlobalTriggerObjectMapRecord>(edm::InputTag("hltL1GtObjectMap"));
102 
103  // Read correction factors
104  std::ifstream infile(cfile);
105  if (!infile.is_open()) {
106  theRecalib_ = false;
107  edm::LogInfo("AnalyzerMB") << "Cannot open '" << cfile
108  << "' for the correction file";
109  } else {
110  unsigned int ndets(0), nrec(0);
111  while(1) {
112  unsigned int id;
113  double cfac;
114  infile >> id >> cfac;
115  if (!infile.good()) break;
116  HcalDetId detId(id);
117  nrec++;
118  std::map<DetId,double>::iterator itr = corrFactor_.find(detId);
119  if (itr == corrFactor_.end()) {
120  corrFactor_[detId] = cfac;
121  ndets++;
122  }
123  }
124  infile.close();
125  edm::LogInfo("AnalyzerMB") << "Reads " << nrec << " correction factors for "
126  << ndets << " detIds";
127  theRecalib_ = (ndets>0);
128  }
129 
130  edm::LogInfo("AnalyzerMB") << " Flags (ReCalib): " << theRecalib_
131  << " (IgnoreL1): " << ignoreL1_ << " (NZS) "
132  << runNZS_ << " and with " << ieta.size()
133  << " detId for full histogram";
134  edm::LogInfo("AnalyzerMB") << "Thresholds for HB " << eLowHB_ << ":"
135  << eHighHB_ << " for HE " << eLowHE_ << ":"
136  << eHighHE_ << " for HF " << eLowHF_ << ":"
137  << eHighHF_;
138  for (unsigned int k=0; k<ieta.size(); ++k) {
139  HcalSubdetector subd = ((std::abs(ieta[k]) > 29) ? HcalForward :
140  (std::abs(ieta[k]) > 16) ? HcalEndcap :
141  ((std::abs(ieta[k]) == 16) && (depth[k] == 3)) ? HcalEndcap :
142  (depth[k] == 4) ? HcalOuter : HcalBarrel);
143  unsigned int id = (HcalDetId(subd,ieta[k],iphi[k],depth[k])).rawId();
144  hcalID_.push_back(id);
145  edm::LogInfo("AnalyzerMB") << "DetId[" << k << "] " << HcalDetId(id);
146  }
147  edm::LogInfo("AnalyzerMB") << "Select on " << trigbit_.size()
148  << " L1 Trigger selection";
149  for (unsigned int k=0; k<trigbit_.size(); ++k)
150  edm::LogInfo("AnalyzerMB") << "Bit[" << k << "] " << trigbit_[k];
151 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > trigbit_
std::string fOutputFileName_
std::map< DetId, double > corrFactor_
edm::EDGetTokenT< HFRecHitCollection > tok_hfrecoMB_
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< unsigned int > hcalID_
edm::EDGetTokenT< HBHERecHitCollection > tok_hbherecoMB_
list infile
Definition: EdgesToViz.py:90
edm::EDGetTokenT< L1GlobalTriggerObjectMapRecord > tok_hltL1GtMap_
RecAnalyzerMinbias::~RecAnalyzerMinbias ( )

Definition at line 153 of file RecAnalyzerMinbias.cc.

153 {}

Member Function Documentation

void RecAnalyzerMinbias::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Implements edm::EDAnalyzer.

Definition at line 302 of file RecAnalyzerMinbias.cc.

References analyzeHcal(), spr::find(), edm::Event::getByToken(), ignoreL1_, edm::HandleBase::isValid(), convertSQLiteXML::ok, edm::Handle< T >::product(), mathSSE::return(), rnnum_, edm::Event::run(), runNZS_, benchmark_cfg::select, edm::SortedCollection< T, SORT >::size(), tok_hbherecoMB_, tok_hfrecoMB_, tok_hltL1GtMap_, and trigbit_.

302  {
303 
304  rnnum_ = (float)iEvent.run();
305 
307  iEvent.getByToken(tok_hbherecoMB_, hbheMB);
308  if (!hbheMB.isValid()) {
309  edm::LogWarning("AnalyzerMB") << "HcalCalibAlgos: Error! can't get hbhe product!";
310  return ;
311  }
312  const HBHERecHitCollection HithbheMB = *(hbheMB.product());
313  edm::LogInfo("AnalyzerMB") << "HBHE MB size of collection "<<HithbheMB.size();
314  if (HithbheMB.size() < 5100 && runNZS_) {
315  edm::LogWarning("AnalyzerMB") << "HBHE problem " << rnnum_ << " size "
316  << HithbheMB.size();
317  return;
318  }
319 
321  iEvent.getByToken(tok_hfrecoMB_, hfMB);
322  if (!hfMB.isValid()) {
323  edm::LogWarning("AnalyzerMB") << "HcalCalibAlgos: Error! can't get hbhe product!";
324  return ;
325  }
326  const HFRecHitCollection HithfMB = *(hfMB.product());
327  edm::LogInfo("AnalyzerMB") << "HF MB size of collection " << HithfMB.size();
328  if (HithfMB.size() < 1700 && runNZS_) {
329  edm::LogWarning("AnalyzerMB") << "HF problem " << rnnum_ << " size "
330  << HithfMB.size();
331  return;
332  }
333 
334  bool select(false);
335  if (trigbit_.size() > 0) {
337  iEvent.getByToken(tok_hltL1GtMap_, gtObjectMapRecord);
338  if (gtObjectMapRecord.isValid()) {
339  const std::vector<L1GlobalTriggerObjectMap>& objMapVec = gtObjectMapRecord->gtObjectMap();
340  for (std::vector<L1GlobalTriggerObjectMap>::const_iterator itMap = objMapVec.begin();
341  itMap != objMapVec.end(); ++itMap) {
342  bool resultGt = (*itMap).algoGtlResult();
343  if (resultGt) {
344  int algoBit = (*itMap).algoBitNumber();
345  if (std::find(trigbit_.begin(),trigbit_.end(),algoBit) !=
346  trigbit_.end()) {
347  select = true;
348  break;
349  }
350  }
351  }
352  }
353  }
354 
355  if (ignoreL1_ || ((trigbit_.size() > 0) && select)) {
356  analyzeHcal(HithbheMB, HithfMB, 1, true);
357  } else if ((!ignoreL1_) && (trigbit_.size() == 0)) {
359  iEvent.getByToken(tok_hltL1GtMap_, gtObjectMapRecord);
360  if (gtObjectMapRecord.isValid()) {
361  const std::vector<L1GlobalTriggerObjectMap>& objMapVec = gtObjectMapRecord->gtObjectMap();
362  bool ok(false);
363  for (std::vector<L1GlobalTriggerObjectMap>::const_iterator itMap = objMapVec.begin();
364  itMap != objMapVec.end(); ++itMap) {
365  bool resultGt = (*itMap).algoGtlResult();
366  if (resultGt) {
367  int algoBit = (*itMap).algoBitNumber();
368  analyzeHcal(HithbheMB, HithfMB, algoBit, (!ok));
369  ok = true;
370  }
371  }
372  if (!ok) {
373  edm::LogInfo("AnalyzerMB") << "No passed L1 Trigger found";
374  }
375  }
376  }
377 }
std::vector< int > trigbit_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
void analyzeHcal(const HBHERecHitCollection &, const HFRecHitCollection &, int, bool)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
return((rh^lh)&mask)
edm::EDGetTokenT< HFRecHitCollection > tok_hfrecoMB_
RunNumber_t run() const
Definition: Event.h:92
bool isValid() const
Definition: HandleBase.h:75
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< HBHERecHitCollection > tok_hbherecoMB_
size_type size() const
edm::EDGetTokenT< L1GlobalTriggerObjectMapRecord > tok_hltL1GtMap_
void RecAnalyzerMinbias::analyzeHcal ( const HBHERecHitCollection HithbheMB,
const HFRecHitCollection HithfMB,
int  algoBit,
bool  fill 
)
private

Definition at line 379 of file RecAnalyzerMinbias.cc.

References edm::SortedCollection< T, SORT >::begin(), corrFactor_, eHighHB_, eHighHE_, eHighHF_, eLowHB_, eLowHE_, eLowHF_, edm::SortedCollection< T, SORT >::end(), CaloRecHit::energy(), HcalEndcap, hcalID_, histHB_, histHE_, histHF_, histo_, i, info(), myMap_, DetId::rawId(), rnnum_, runNZS_, HcalDetId::subdet(), and theRecalib_.

Referenced by analyze().

381  {
382  // Signal part for HB HE
383  for (HBHERecHitCollection::const_iterator hbheItr=HithbheMB.begin();
384  hbheItr!=HithbheMB.end(); hbheItr++) {
385  // Recalibration of energy
386  DetId mydetid = hbheItr->id().rawId();
387  double icalconst(1.);
388  if (theRecalib_) {
389  std::map<DetId,double>::iterator itr = corrFactor_.find(mydetid);
390  if (itr != corrFactor_.end()) icalconst = itr->second;
391  }
392  HBHERecHit aHit(hbheItr->id(),hbheItr->energy()*icalconst,hbheItr->time());
393  double energyhit = aHit.energy();
394  DetId id = (*hbheItr).detid();
395  HcalDetId hid = HcalDetId(id);
396  double eLow = (hid.subdet() == HcalEndcap) ? eLowHE_ : eLowHB_;
397  double eHigh = (hid.subdet() == HcalEndcap) ? eHighHE_ : eHighHB_;
398  if (fill) {
399  for (unsigned int i = 0; i < hcalID_.size(); i++) {
400  if (hcalID_[i] == id.rawId()) {
401  histo_[i]->Fill(energyhit);
402  break;
403  }
404  }
405  std::map<HcalDetId,TH1D*>::iterator itr1 = histHB_.find(hid);
406  if (itr1 != histHB_.end()) itr1->second->Fill(energyhit);
407  std::map<HcalDetId,TH1D*>::iterator itr2 = histHE_.find(hid);
408  if (itr2 != histHE_.end()) itr2->second->Fill(energyhit);
409  }
410 
411  if (runNZS_ || (energyhit >= eLow && energyhit <= eHigh)) {
412  std::map<std::pair<int,HcalDetId>,myInfo>::iterator itr1 = myMap_.find(std::pair<int,HcalDetId>(algoBit,hid));
413  if (itr1 == myMap_.end()) {
414  myInfo info;
415  myMap_[std::pair<int,HcalDetId>(algoBit,hid)] = info;
416  itr1 = myMap_.find(std::pair<int,HcalDetId>(algoBit,hid));
417  }
418  itr1->second.theMB0++;
419  itr1->second.theMB1 += energyhit;
420  itr1->second.theMB2 += (energyhit*energyhit);
421  itr1->second.theMB3 += (energyhit*energyhit*energyhit);
422  itr1->second.theMB4 += (energyhit*energyhit*energyhit*energyhit);
423  itr1->second.runcheck = rnnum_;
424  }
425  } // HBHE_MB
426 
427  // Signal part for HF
428  for (HFRecHitCollection::const_iterator hfItr=HithfMB.begin();
429  hfItr!=HithfMB.end(); hfItr++) {
430  // Recalibration of energy
431  DetId mydetid = hfItr->id().rawId();
432  double icalconst(1.);
433  if (theRecalib_) {
434  std::map<DetId,double>::iterator itr = corrFactor_.find(mydetid);
435  if (itr != corrFactor_.end()) icalconst = itr->second;
436  }
437  HFRecHit aHit(hfItr->id(),hfItr->energy()*icalconst,hfItr->time());
438 
439  double energyhit = aHit.energy();
440  DetId id = (*hfItr).detid();
441  HcalDetId hid = HcalDetId(id);
442  if (fill) {
443  for (unsigned int i = 0; i < hcalID_.size(); i++) {
444  if (hcalID_[i] == id.rawId()) {
445  histo_[i]->Fill(energyhit);
446  break;
447  }
448  }
449  std::map<HcalDetId,TH1D*>::iterator itr1 = histHF_.find(hid);
450  if (itr1 != histHF_.end()) itr1->second->Fill(energyhit);
451  }
452 
453  //
454  // Remove PMT hits
455  //
456  if ((runNZS_ && fabs(energyhit) <= 40.) ||
457  (energyhit >= eLowHF_ && energyhit <= eHighHF_)) {
458  std::map<std::pair<int,HcalDetId>,myInfo>::iterator itr1 = myMap_.find(std::pair<int,HcalDetId>(algoBit,hid));
459  if (itr1 == myMap_.end()) {
460  myInfo info;
461  myMap_[std::pair<int,HcalDetId>(algoBit,hid)] = info;
462  itr1 = myMap_.find(std::pair<int,HcalDetId>(algoBit,hid));
463  }
464  itr1->second.theMB0++;
465  itr1->second.theMB1 += energyhit;
466  itr1->second.theMB2 += (energyhit*energyhit);
467  itr1->second.theMB3 += (energyhit*energyhit*energyhit);
468  itr1->second.theMB4 += (energyhit*energyhit*energyhit*energyhit);
469  itr1->second.runcheck = rnnum_;
470  }
471  }
472 }
std::map< HcalDetId, TH1D * > histHE_
int i
Definition: DBlmapReader.cc:9
static const TGPicture * info(bool iBackgroundIsBlack)
string fill
Definition: lumiContext.py:319
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:45
std::map< HcalDetId, TH1D * > histHB_
std::vector< HBHERecHit >::const_iterator const_iterator
std::map< DetId, double > corrFactor_
std::vector< TH1D * > histo_
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
float energy() const
Definition: CaloRecHit.h:17
const_iterator end() const
Definition: DetId.h:18
std::vector< unsigned int > hcalID_
std::map< std::pair< int, HcalDetId >, myInfo > myMap_
std::map< HcalDetId, TH1D * > histHF_
const_iterator begin() const
void RecAnalyzerMinbias::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 155 of file RecAnalyzerMinbias.cc.

References cells, depth, fOutputFileName_, AnalysisDataFormats_SUSYBSMObjects::hc, hcalID_, histo_, hOutputFile_, i, ieta, iphi, mom0_MB, mom1_MB, mom2_MB, mom4_MB, myMap_, mysubd, myTree_, mergeVDriftHistosByStation::name, rnnumber, AlCaHLTBitMon_QueryRunRegistry::string, indexGen::title, trigbit, SiStripMonitorClusterAlca_cfi::xmax, and SiStripMonitorClusterAlca_cfi::xmin.

155  {
156 
157  std::string hc[5] = {"Empty", "HB", "HE", "HO", "HF"};
158  char name[700], title[700];
159  for (unsigned int i=0; i<hcalID_.size(); i++) {
160  HcalDetId id = HcalDetId(hcalID_[i]);
161  int subdet = id.subdetId();
162  sprintf (name, "%s%d_%d_%d", hc[subdet].c_str(), id.ieta(), id.iphi(), id.depth());
163  sprintf (title, "Energy Distribution for %s ieta %d iphi %d depth %d", hc[subdet].c_str(), id.ieta(), id.iphi(), id.depth());
164  double xmin = (subdet == 4) ? -10 : -1;
165  double xmax = (subdet == 4) ? 90 : 9;
166  TH1D* hh = new TH1D(name, title, 50, xmin, xmax);
167  histo_.push_back(hh);
168  };
169 
170  hOutputFile_ = new TFile(fOutputFileName_.c_str(), "RECREATE");
171  myTree_ = new TTree("RecJet","RecJet Tree");
172  myTree_->Branch("cells", &cells, "cells/I");
173  myTree_->Branch("mysubd", &mysubd, "mysubd/I");
174  myTree_->Branch("depth", &depth, "depth/I");
175  myTree_->Branch("ieta", &ieta, "ieta/I");
176  myTree_->Branch("iphi", &iphi, "iphi/I");
177  myTree_->Branch("mom0_MB", &mom0_MB, "mom0_MB/F");
178  myTree_->Branch("mom1_MB", &mom1_MB, "mom1_MB/F");
179  myTree_->Branch("mom2_MB", &mom2_MB, "mom2_MB/F");
180  myTree_->Branch("mom3_MB", &mom2_MB, "mom3_MB/F");
181  myTree_->Branch("mom4_MB", &mom4_MB, "mom4_MB/F");
182  myTree_->Branch("trigbit", &trigbit, "trigbit/I");
183  myTree_->Branch("rnnumber", &rnnumber, "rnnumber/D");
184 
185  myMap_.clear();
186 }
int i
Definition: DBlmapReader.cc:9
std::string fOutputFileName_
std::vector< TH1D * > histo_
std::vector< unsigned int > hcalID_
susybsm::HSCParticleCollection hc
Definition: classes.h:25
std::map< std::pair< int, HcalDetId >, myInfo > myMap_
void RecAnalyzerMinbias::beginRun ( edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 189 of file RecAnalyzerMinbias.cc.

References depth, eHighHB_, eHighHE_, eHighHF_, eta, fillHist_, edm::EventSetup::get(), h, HcalBarrel, HcalEndcap, HcalForward, histHB_, histHE_, histHF_, init_, edm::ESHandleBase::isValid(), HcalTopology::maxDepthHB(), mergeVDriftHistosByStation::name, phi, edm::ESHandle< class >::product(), indexGen::title, and HcalTopology::valid().

189  {
190 
191  if (!init_) {
192  init_ = true;
193  if (fillHist_) {
195  iS.get<IdealGeometryRecord>().get(htopo);
196  if (htopo.isValid()) {
197  const HcalTopology* hcaltopology = htopo.product();
198 
199  char name[700], title[700];
200  // For HB
201  int maxDepthHB = hcaltopology->maxDepthHB();
202  int nbinHB = int(200*eHighHB_);
203  for (int eta = -50; eta < 50; eta++) {
204  for (int phi = 0; phi < 100; phi++) {
205  for (int depth = 1; depth < maxDepthHB; depth++) {
206  HcalDetId cell (HcalBarrel, eta, phi, depth);
207  if (hcaltopology->valid(cell)) {
208  sprintf (name, "HBeta%dphi%ddep%d", eta, phi, depth);
209  sprintf (title,"Energy (HB #eta %d #phi %d depth %d)", eta, phi, depth);
210  TH1D* h = new TH1D(name, title, nbinHB, 0, 2*eHighHB_);
211  histHB_[cell] = h;
212  }
213  }
214  }
215  }
216  // For HE
217  int maxDepthHE = hcaltopology->maxDepthHB();
218  int nbinHE = int(200*eHighHE_);
219  for (int eta = -50; eta < 50; eta++) {
220  for (int phi = 0; phi < 100; phi++) {
221  for (int depth = 1; depth < maxDepthHE; depth++) {
222  HcalDetId cell (HcalEndcap, eta, phi, depth);
223  if (hcaltopology->valid(cell)) {
224  sprintf (name, "HEeta%dphi%ddep%d", eta, phi, depth);
225  sprintf (title,"Energy (HE #eta %d #phi %d depth %d)", eta, phi, depth);
226  TH1D* h = new TH1D(name, title, nbinHE, 0, 2*eHighHE_);
227  histHE_[cell] = h;
228  }
229  }
230  }
231  }
232  // For HF
233  int maxDepthHF = 4;
234  int nbinHF = int(200*eHighHF_);
235  for (int eta = -50; eta < 50; eta++) {
236  for (int phi = 0; phi < 100; phi++) {
237  for (int depth = 1; depth < maxDepthHF; depth++) {
238  HcalDetId cell (HcalForward, eta, phi, depth);
239  if (hcaltopology->valid(cell)) {
240  sprintf (name, "HFeta%dphi%ddep%d", eta, phi, depth);
241  sprintf (title,"Energy (HF #eta %d #phi %d depth %d)", eta, phi, depth);
242  TH1D* h = new TH1D(name, title, nbinHF, 0, 2*eHighHF_);
243  histHF_[cell] = h;
244  }
245  }
246  }
247  }
248  }
249  }
250  }
251 }
std::map< HcalDetId, TH1D * > histHE_
std::map< HcalDetId, TH1D * > histHB_
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
T const * product() const
Definition: ESHandle.h:86
virtual bool valid(const DetId &id) const
int maxDepthHB() const
Definition: HcalTopology.h:128
bool isValid() const
Definition: ESHandle.h:47
std::map< HcalDetId, TH1D * > histHF_
void RecAnalyzerMinbias::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 255 of file RecAnalyzerMinbias.cc.

References cells, depth, histHB_, histHE_, histHF_, histo_, hOutputFile_, i, ieta, info(), iphi, python.multivaluedict::map(), mom0_MB, mom1_MB, mom2_MB, mom3_MB, mom4_MB, myMap_, mysubd, myTree_, rnnumber, RecAnalyzerMinbias::myInfo::runcheck, RecAnalyzerMinbias::myInfo::theMB0, RecAnalyzerMinbias::myInfo::theMB1, RecAnalyzerMinbias::myInfo::theMB2, RecAnalyzerMinbias::myInfo::theMB3, RecAnalyzerMinbias::myInfo::theMB4, and trigbit.

255  {
256 
257  cells = 0;
258  for (std::map<std::pair<int,HcalDetId>,myInfo>::const_iterator itr=myMap_.begin(); itr != myMap_.end(); ++itr) {
259  edm::LogInfo("AnalyzerMB") << "Fired trigger bit number "<<itr->first.first;
260  myInfo info = itr->second;
261  if (info.theMB0 > 0) {
262  mom0_MB = info.theMB0;
263  mom1_MB = info.theMB1;
264  mom2_MB = info.theMB2;
265  mom3_MB = info.theMB3;
266  mom4_MB = info.theMB4;
267  rnnumber = info.runcheck;
268  trigbit = itr->first.first;
269  mysubd = itr->first.second.subdet();
270  depth = itr->first.second.depth();
271  iphi = itr->first.second.iphi();
272  ieta = itr->first.second.ieta();
273  edm::LogInfo("AnalyzerMB") << " Result= " << trigbit << " " << mysubd
274  << " " << ieta << " " << iphi << " mom0 "
275  << mom0_MB << " mom1 " << mom1_MB << " mom2 "
276  << mom2_MB << " mom3 " << mom3_MB << " mom4 "
277  << mom4_MB;
278  myTree_->Fill();
279  cells++;
280  }
281  }
282  edm::LogInfo("AnalyzerMB") << "cells" << " " << cells;
283 
284  hOutputFile_->Write();
285  hOutputFile_->cd();
286  for (unsigned int i=0; i<histo_.size(); i++) histo_[i]->Write();
287  for (std::map<HcalDetId,TH1D*>::iterator itr=histHB_.begin();
288  itr != histHB_.end(); ++itr) itr->second->Write();
289  for (std::map<HcalDetId,TH1D*>::iterator itr=histHE_.begin();
290  itr != histHE_.end(); ++itr) itr->second->Write();
291  for (std::map<HcalDetId,TH1D*>::iterator itr=histHF_.begin();
292  itr != histHF_.end(); ++itr) itr->second->Write();
293  myTree_->Write();
294  hOutputFile_->Close() ;
295 }
std::map< HcalDetId, TH1D * > histHE_
int i
Definition: DBlmapReader.cc:9
static const TGPicture * info(bool iBackgroundIsBlack)
std::map< HcalDetId, TH1D * > histHB_
std::vector< TH1D * > histo_
std::map< std::pair< int, HcalDetId >, myInfo > myMap_
std::map< HcalDetId, TH1D * > histHF_

Member Data Documentation

int RecAnalyzerMinbias::cells
private

Definition at line 69 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

std::map<DetId,double> RecAnalyzerMinbias::corrFactor_
private

Definition at line 53 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), and RecAnalyzerMinbias().

int RecAnalyzerMinbias::depth
private
double RecAnalyzerMinbias::eHighHB_
private

Definition at line 51 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginRun(), and RecAnalyzerMinbias().

double RecAnalyzerMinbias::eHighHE_
private

Definition at line 51 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginRun(), and RecAnalyzerMinbias().

double RecAnalyzerMinbias::eHighHF_
private

Definition at line 52 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginRun(), and RecAnalyzerMinbias().

double RecAnalyzerMinbias::eLowHB_
private

Definition at line 51 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), and RecAnalyzerMinbias().

double RecAnalyzerMinbias::eLowHE_
private

Definition at line 51 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), and RecAnalyzerMinbias().

double RecAnalyzerMinbias::eLowHF_
private

Definition at line 52 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), and RecAnalyzerMinbias().

bool RecAnalyzerMinbias::fillHist_
private

Definition at line 50 of file RecAnalyzerMinbias.cc.

Referenced by beginRun(), and RecAnalyzerMinbias().

std::string RecAnalyzerMinbias::fOutputFileName_
private

Definition at line 49 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and RecAnalyzerMinbias().

std::vector<unsigned int> RecAnalyzerMinbias::hcalID_
private

Definition at line 54 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginJob(), and RecAnalyzerMinbias().

std::map<HcalDetId,TH1D*> RecAnalyzerMinbias::histHB_
private

Definition at line 58 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginRun(), and endJob().

std::map<HcalDetId,TH1D*> RecAnalyzerMinbias::histHE_
private

Definition at line 58 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginRun(), and endJob().

std::map<HcalDetId,TH1D*> RecAnalyzerMinbias::histHF_
private

Definition at line 58 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginRun(), and endJob().

std::vector<TH1D*> RecAnalyzerMinbias::histo_
private

Definition at line 57 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginJob(), and endJob().

TFile* RecAnalyzerMinbias::hOutputFile_
private

Definition at line 55 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

int RecAnalyzerMinbias::ieta
private

Definition at line 69 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), endJob(), and RecAnalyzerMinbias().

bool RecAnalyzerMinbias::ignoreL1_
private

Definition at line 50 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), and RecAnalyzerMinbias().

bool RecAnalyzerMinbias::init_
private

Definition at line 50 of file RecAnalyzerMinbias.cc.

Referenced by beginRun().

int RecAnalyzerMinbias::iphi
private

Definition at line 69 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), endJob(), and RecAnalyzerMinbias().

float RecAnalyzerMinbias::mom0_MB
private

Definition at line 70 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

float RecAnalyzerMinbias::mom1_MB
private

Definition at line 70 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

float RecAnalyzerMinbias::mom2_MB
private

Definition at line 70 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

float RecAnalyzerMinbias::mom3_MB
private

Definition at line 70 of file RecAnalyzerMinbias.cc.

Referenced by endJob().

float RecAnalyzerMinbias::mom4_MB
private

Definition at line 70 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

std::map<std::pair<int,HcalDetId>,myInfo> RecAnalyzerMinbias::myMap_
private

Definition at line 71 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), beginJob(), and endJob().

int RecAnalyzerMinbias::mysubd
private

Definition at line 69 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

TTree* RecAnalyzerMinbias::myTree_
private

Definition at line 56 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

double RecAnalyzerMinbias::rnnum_
private

Definition at line 60 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), and analyzeHcal().

double RecAnalyzerMinbias::rnnumber
private

Definition at line 68 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

bool RecAnalyzerMinbias::runNZS_
private

Definition at line 50 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), analyzeHcal(), and RecAnalyzerMinbias().

bool RecAnalyzerMinbias::theRecalib_
private

Definition at line 50 of file RecAnalyzerMinbias.cc.

Referenced by analyzeHcal(), and RecAnalyzerMinbias().

edm::EDGetTokenT<HBHERecHitCollection> RecAnalyzerMinbias::tok_hbherecoMB_
private

Definition at line 72 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), and RecAnalyzerMinbias().

edm::EDGetTokenT<HFRecHitCollection> RecAnalyzerMinbias::tok_hfrecoMB_
private

Definition at line 73 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), and RecAnalyzerMinbias().

edm::EDGetTokenT<L1GlobalTriggerObjectMapRecord> RecAnalyzerMinbias::tok_hltL1GtMap_
private

Definition at line 74 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), and RecAnalyzerMinbias().

int RecAnalyzerMinbias::trigbit
private

Definition at line 69 of file RecAnalyzerMinbias.cc.

Referenced by beginJob(), and endJob().

std::vector<int> RecAnalyzerMinbias::trigbit_
private

Definition at line 59 of file RecAnalyzerMinbias.cc.

Referenced by analyze(), and RecAnalyzerMinbias().