CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
cms::MinBias Class Reference

#include <MinBias.h>

Inheritance diagram for cms::MinBias:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void endJob ()
 
 MinBias (const edm::ParameterSet &)
 
 ~MinBias ()
 
- 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 Attributes

bool allowMissingInputs_
 
int depth
 
float eta
 
std::string fOutputFileName
 
const CaloGeometrygeo
 
std::string hbheLabel_
 
edm::EDGetTokenT
< HBHERecHitCollection
hbheToken_
 
std::string hfLabel_
 
edm::EDGetTokenT
< HFRecHitCollection
hfToken_
 
std::string hoLabel_
 
edm::EDGetTokenT
< HORecHitCollection
hoToken_
 
TFile * hOutputFile
 
int ieta
 
int ievent
 
int iphi
 
float mom1
 
float mom2
 
float mom3
 
float mom4
 
int mydet
 
int mysubd
 
TTree * myTree
 
float occup
 
float phi
 
std::map< DetId, double > theFillDetMap0
 
std::map< DetId, double > theFillDetMap1
 
std::map< DetId, double > theFillDetMap2
 
std::map< DetId, double > theFillDetMap3
 
std::map< DetId, double > theFillDetMap4
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- 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 42 of file MinBias.h.

Constructor & Destructor Documentation

cms::MinBias::MinBias ( const edm::ParameterSet iConfig)
explicit

Definition at line 9 of file MinBias.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), HLT_25ns14e33_v1_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

10 {
11  // get names of modules, producing object collections
12  hbheLabel_= iConfig.getParameter<std::string>("hbheInput");
13  hoLabel_= iConfig.getParameter<std::string>("hoInput");
14  hfLabel_= iConfig.getParameter<std::string>("hfInput");
15  hbheToken_= mayConsume<HBHERecHitCollection>(edm::InputTag(hbheLabel_));
16  hoToken_=mayConsume<HORecHitCollection>(edm::InputTag(hoLabel_));
17  hfToken_=mayConsume<HFRecHitCollection>(edm::InputTag(hfLabel_));
18  allowMissingInputs_=iConfig.getUntrackedParameter<bool>("AllowMissingInputs",false);
19  // get name of output file with histogramms
20  fOutputFileName = iConfig.getUntrackedParameter<std::string>("HistOutFile");
21 
22 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< HFRecHitCollection > hfToken_
Definition: MinBias.h:58
std::string hoLabel_
Definition: MinBias.h:55
edm::EDGetTokenT< HORecHitCollection > hoToken_
Definition: MinBias.h:57
bool allowMissingInputs_
Definition: MinBias.h:62
edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
Definition: MinBias.h:56
std::string hbheLabel_
Definition: MinBias.h:55
std::string hfLabel_
Definition: MinBias.h:55
std::string fOutputFileName
Definition: MinBias.h:61
cms::MinBias::~MinBias ( )

Definition at line 25 of file MinBias.cc.

26 {
27 
28  // do anything here that needs to be done at desctruction time
29  // (e.g. close files, deallocate resources etc.)
30 
31 }

Member Function Documentation

void cms::MinBias::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 107 of file MinBias.cc.

References gather_cfg::cout, cond::rpcobgas::detid, edm::EventSetup::get(), edm::Event::getByToken(), DetId::Hcal, edm::HandleBase::isValid(), funct::pow(), and edm::ESHandle< class >::product().

108 {
109 
110  using namespace edm;
111  if(ievent == 0 ){
113  iSetup.get<CaloGeometryRecord>().get(pG);
114  geo = pG.product();
115  std::vector<DetId> did = geo->getValidDetIds();
116 
117  for(std::vector<DetId>::const_iterator id=did.begin(); id != did.end(); id++)
118  {
119  if( (*id).det() == DetId::Hcal ) {
120  theFillDetMap0[*id] = 0.;
121  theFillDetMap1[*id] = 0.;
122  theFillDetMap2[*id] = 0.;
123  theFillDetMap3[*id] = 0.;
124  theFillDetMap4[*id] = 0.;
125  }
126  }
127  }
128 
129 
130 
131  if (!hbheLabel_.empty()) {
133  iEvent.getByToken(hbheToken_,hbhe);
134  if (!hbhe.isValid()) {
135  // can't find it!
136  if (!allowMissingInputs_) {
137  *hbhe; // will throw the proper exception
138  }
139  } else {
140  for(HBHERecHitCollection::const_iterator hbheItr = (*hbhe).begin();
141  hbheItr != (*hbhe).end(); ++hbheItr)
142  {
143  DetId id = (hbheItr)->detid();
144  if( (*hbheItr).energy() > 0. ) std::cout<<" Energy = "<<(*hbheItr).energy()<<std::endl;
145  theFillDetMap0[id] = theFillDetMap0[id]+ 1.;
146  theFillDetMap1[id] = theFillDetMap1[id]+(*hbheItr).energy();
147  theFillDetMap2[id] = theFillDetMap2[id]+pow((*hbheItr).energy(),2);
148  theFillDetMap3[id] = theFillDetMap3[id]+pow((*hbheItr).energy(),3);
149  theFillDetMap4[id] = theFillDetMap4[id]+pow((*hbheItr).energy(),4);
150  }
151  }
152  }
153 
154  if (!hoLabel_.empty()) {
156  iEvent.getByToken(hoToken_,ho);
157  if (!ho.isValid()) {
158  // can't find it!
159  if (!allowMissingInputs_) {
160  *ho; // will throw the proper exception
161  }
162  } else {
163  for(HORecHitCollection::const_iterator hoItr = (*ho).begin();
164  hoItr != (*ho).end(); ++hoItr)
165  {
166  DetId id = (hoItr)->detid();
167  theFillDetMap0[id] = theFillDetMap0[id]+ 1.;
168  theFillDetMap1[id] = theFillDetMap1[id]+(*hoItr).energy();
169  theFillDetMap2[id] = theFillDetMap2[id]+pow((*hoItr).energy(),2);
170  theFillDetMap3[id] = theFillDetMap3[id]+pow((*hoItr).energy(),3);
171  theFillDetMap4[id] = theFillDetMap4[id]+pow((*hoItr).energy(),4);
172  }
173  }
174  }
175 
176  if (!hfLabel_.empty()) {
178  iEvent.getByToken(hfToken_,hf);
179  if (!hf.isValid()) {
180  // can't find it!
181  if (!allowMissingInputs_) {
182  *hf; // will throw the proper exception
183  }
184  } else {
185  for(HFRecHitCollection::const_iterator hfItr = (*hf).begin();
186  hfItr != (*hf).end(); ++hfItr)
187  {
188  DetId id = (hfItr)->detid();
189  theFillDetMap0[id] = theFillDetMap0[id]+ 1.;
190  theFillDetMap1[id] = theFillDetMap1[id]+(*hfItr).energy();
191  theFillDetMap2[id] = theFillDetMap2[id]+pow((*hfItr).energy(),2);
192  theFillDetMap3[id] = theFillDetMap3[id]+pow((*hfItr).energy(),3);
193  theFillDetMap4[id] = theFillDetMap4[id]+pow((*hfItr).energy(),4);
194  }
195  }
196  }
197 
198 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
int ievent
Definition: MinBias.h:68
edm::EDGetTokenT< HFRecHitCollection > hfToken_
Definition: MinBias.h:58
std::vector< HBHERecHit >::const_iterator const_iterator
std::string hoLabel_
Definition: MinBias.h:55
std::map< DetId, double > theFillDetMap4
Definition: MinBias.h:77
edm::EDGetTokenT< HORecHitCollection > hoToken_
Definition: MinBias.h:57
bool allowMissingInputs_
Definition: MinBias.h:62
std::map< DetId, double > theFillDetMap1
Definition: MinBias.h:74
const CaloGeometry * geo
Definition: MinBias.h:71
std::map< DetId, double > theFillDetMap0
Definition: MinBias.h:73
std::map< DetId, double > theFillDetMap2
Definition: MinBias.h:75
edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
Definition: MinBias.h:56
bool isValid() const
Definition: HandleBase.h:75
std::map< DetId, double > theFillDetMap3
Definition: MinBias.h:76
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:90
std::string hbheLabel_
Definition: MinBias.h:55
tuple cout
Definition: gather_cfg.py:121
std::string hfLabel_
Definition: MinBias.h:55
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void cms::MinBias::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 33 of file MinBias.cc.

References HLT_25ns14e33_v1_cff::depth, eta(), and phi.

34 {
35  hOutputFile = new TFile( fOutputFileName.c_str(), "RECREATE" ) ;
36  myTree = new TTree("RecJet","RecJet Tree");
37  myTree->Branch("mydet", &mydet, "mydet/I");
38  myTree->Branch("mysubd", &mysubd, "mysubd/I");
39  myTree->Branch("depth", &depth, "depth/I");
40  myTree->Branch("ieta", &ieta, "ieta/I");
41  myTree->Branch("iphi", &iphi, "iphi/I");
42  myTree->Branch("eta", &eta, "eta/F");
43  myTree->Branch("phi", &phi, "phi/F");
44  myTree->Branch("mom1", &mom1, "mom1/F");
45  myTree->Branch("mom2", &mom2, "mom2/F");
46  myTree->Branch("mom3", &mom3, "mom3/F");
47  myTree->Branch("mom4", &mom4, "mom4/F");
48 
49  ievent = 0;
50 
51 }
float mom2
Definition: MinBias.h:70
float mom1
Definition: MinBias.h:70
float mom4
Definition: MinBias.h:70
int ievent
Definition: MinBias.h:68
int mydet
Definition: MinBias.h:68
int ieta
Definition: MinBias.h:68
int depth
Definition: MinBias.h:68
float eta
Definition: MinBias.h:69
TTree * myTree
Definition: MinBias.h:66
float phi
Definition: MinBias.h:69
float mom3
Definition: MinBias.h:70
int iphi
Definition: MinBias.h:68
TFile * hOutputFile
Definition: MinBias.h:64
int mysubd
Definition: MinBias.h:68
std::string fOutputFileName
Definition: MinBias.h:61
void cms::MinBias::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 53 of file MinBias.cc.

References gather_cfg::cout, HcalDetId::depth(), HLT_25ns14e33_v1_cff::depth, PV3DBase< T, PVType, FrameType >::eta(), eta(), DetId::Hcal, i, HcalDetId::ieta(), HcalDetId::iphi(), phi, PV3DBase< T, PVType, FrameType >::phi(), funct::pow(), and HcalDetId::subdet().

54 {
55 
56  const std::vector<DetId>& did = geo->getSubdetectorGeometry( DetId::Hcal, 1 )->getValidDetIds() ;
57  int i=0;
58  for(std::vector<DetId>::const_iterator id=did.begin(); id != did.end(); id++)
59  {
60 // if( (*id).det() == DetId::Hcal ) {
61  GlobalPoint pos = geo->getPosition(*id);
62  mydet = ((*id).rawId()>>28)&0xF;
63  mysubd = ((*id).rawId()>>25)&0x7;
64  depth = HcalDetId(*id).depth();
65  ieta = HcalDetId(*id).ieta();
66  iphi = HcalDetId(*id).iphi();
67  phi = pos.phi();
68  eta = pos.eta();
69  if ( theFillDetMap0[*id] > 0. )
70  {
74  2.*pow(mom2,3);
75  mom4 = (theFillDetMap4[*id]-4.*mom1*theFillDetMap3[*id]+6.*pow(mom1,2)*theFillDetMap2[*id])/theFillDetMap0[*id]-3.*pow(mom1,4);
76 
77  } else
78  {
79  mom1 = 0.; mom2 = 0.; mom3 = 0.; mom4 = 0.;
80  }
81  std::cout<<" Detector "<<(*id).rawId()<<" mydet "<<mydet<<" "<<mysubd<<" "<<depth<<" "<<
82  HcalDetId(*id).subdet()<<" "<<ieta<<" "<<iphi<<" "<<pos.eta()<<" "<<pos.phi()<<std::endl;
83  std::cout<<" Energy "<<mom1<<" "<<mom2<<std::endl;
84  myTree->Fill();
85  i++;
86 // }
87  }
88  std::cout<<" The number of CaloDet records "<<did.size()<<std::endl;
89  std::cout<<" The number of Hcal records "<<i<<std::endl;
90 
91 
92  std::cout << "===== Start writing user histograms =====" << std::endl;
93  hOutputFile->SetCompressionLevel(2);
94  hOutputFile->cd();
95  myTree->Write();
96  hOutputFile->Close() ;
97  std::cout << "===== End writing user histograms =======" << std::endl;
98 }
float mom2
Definition: MinBias.h:70
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
int i
Definition: DBlmapReader.cc:9
float mom1
Definition: MinBias.h:70
float mom4
Definition: MinBias.h:70
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:30
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
std::map< DetId, double > theFillDetMap4
Definition: MinBias.h:77
virtual const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const
Get a list of valid detector ids (for the given subdetector)
int mydet
Definition: MinBias.h:68
int ieta
Definition: MinBias.h:68
int depth() const
get the tower depth
Definition: HcalDetId.h:40
std::map< DetId, double > theFillDetMap1
Definition: MinBias.h:74
const CaloGeometry * geo
Definition: MinBias.h:71
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
std::map< DetId, double > theFillDetMap0
Definition: MinBias.h:73
std::map< DetId, double > theFillDetMap2
Definition: MinBias.h:75
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
int depth
Definition: MinBias.h:68
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
std::map< DetId, double > theFillDetMap3
Definition: MinBias.h:76
float eta
Definition: MinBias.h:69
TTree * myTree
Definition: MinBias.h:66
float phi
Definition: MinBias.h:69
float mom3
Definition: MinBias.h:70
T eta() const
Definition: PV3DBase.h:76
int iphi
Definition: MinBias.h:68
TFile * hOutputFile
Definition: MinBias.h:64
tuple cout
Definition: gather_cfg.py:121
int mysubd
Definition: MinBias.h:68
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40

Member Data Documentation

bool cms::MinBias::allowMissingInputs_
private

Definition at line 62 of file MinBias.h.

int cms::MinBias::depth
private
float cms::MinBias::eta
private
std::string cms::MinBias::fOutputFileName
private

Definition at line 61 of file MinBias.h.

const CaloGeometry* cms::MinBias::geo
private

Definition at line 71 of file MinBias.h.

std::string cms::MinBias::hbheLabel_
private

Definition at line 55 of file MinBias.h.

edm::EDGetTokenT<HBHERecHitCollection> cms::MinBias::hbheToken_
private

Definition at line 56 of file MinBias.h.

std::string cms::MinBias::hfLabel_
private

Definition at line 55 of file MinBias.h.

edm::EDGetTokenT<HFRecHitCollection> cms::MinBias::hfToken_
private

Definition at line 58 of file MinBias.h.

std::string cms::MinBias::hoLabel_
private

Definition at line 55 of file MinBias.h.

edm::EDGetTokenT<HORecHitCollection> cms::MinBias::hoToken_
private

Definition at line 57 of file MinBias.h.

TFile* cms::MinBias::hOutputFile
private

Definition at line 64 of file MinBias.h.

int cms::MinBias::ieta
private

Definition at line 68 of file MinBias.h.

int cms::MinBias::ievent
private

Definition at line 68 of file MinBias.h.

int cms::MinBias::iphi
private

Definition at line 68 of file MinBias.h.

float cms::MinBias::mom1
private

Definition at line 70 of file MinBias.h.

float cms::MinBias::mom2
private

Definition at line 70 of file MinBias.h.

float cms::MinBias::mom3
private

Definition at line 70 of file MinBias.h.

float cms::MinBias::mom4
private

Definition at line 70 of file MinBias.h.

int cms::MinBias::mydet
private

Definition at line 68 of file MinBias.h.

int cms::MinBias::mysubd
private

Definition at line 68 of file MinBias.h.

TTree* cms::MinBias::myTree
private

Definition at line 66 of file MinBias.h.

float cms::MinBias::occup
private

Definition at line 70 of file MinBias.h.

float cms::MinBias::phi
private

Definition at line 69 of file MinBias.h.

Referenced by Particle.Particle::__str__().

std::map<DetId,double> cms::MinBias::theFillDetMap0
private

Definition at line 73 of file MinBias.h.

std::map<DetId,double> cms::MinBias::theFillDetMap1
private

Definition at line 74 of file MinBias.h.

std::map<DetId,double> cms::MinBias::theFillDetMap2
private

Definition at line 75 of file MinBias.h.

std::map<DetId,double> cms::MinBias::theFillDetMap3
private

Definition at line 76 of file MinBias.h.

std::map<DetId,double> cms::MinBias::theFillDetMap4
private

Definition at line 77 of file MinBias.h.