CMS 3D CMS Logo

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

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void beginRun (edm::Run const &, edm::EventSetup const &) override
 
void endJob () override
 
 MinBias (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Attributes

bool allowMissingInputs_
 
int depth
 
float eta
 
const CaloGeometrygeo_
 
std::string hbheLabel_
 
edm::EDGetTokenT< HBHERecHitCollectionhbheToken_
 
std::string hfLabel_
 
edm::EDGetTokenT< HFRecHitCollectionhfToken_
 
std::string hoLabel_
 
edm::EDGetTokenT< HORecHitCollectionhoToken_
 
int ieta
 
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
 
- 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 &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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 35 of file MinBias.h.

Constructor & Destructor Documentation

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

Definition at line 10 of file MinBias.cc.

References allowMissingInputs_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hbheLabel_, hbheToken_, hfLabel_, hfToken_, hoLabel_, hoToken_, and AlCaHLTBitMon_QueryRunRegistry::string.

10  : geo_(nullptr) {
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 
20 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< HFRecHitCollection > hfToken_
Definition: MinBias.h:51
std::string hoLabel_
Definition: MinBias.h:48
edm::EDGetTokenT< HORecHitCollection > hoToken_
Definition: MinBias.h:50
const CaloGeometry * geo_
Definition: MinBias.h:61
bool allowMissingInputs_
Definition: MinBias.h:53
edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
Definition: MinBias.h:49
std::string hbheLabel_
Definition: MinBias.h:48
std::string hfLabel_
Definition: MinBias.h:48

Member Function Documentation

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

Definition at line 107 of file MinBias.cc.

References allowMissingInputs_, edm::Event::getByToken(), photonIsolationHIProducer_cfi::hbhe, hbheLabel_, hbheToken_, photonIsolationHIProducer_cfi::hf, hfLabel_, hfToken_, photonIsolationHIProducer_cfi::ho, hoLabel_, hoToken_, triggerObjects_cff::id, edm::HandleBase::isValid(), funct::pow(), theFillDetMap0_, theFillDetMap1_, theFillDetMap2_, theFillDetMap3_, and theFillDetMap4_.

107  {
108 
109  if (!hbheLabel_.empty()) {
111  iEvent.getByToken(hbheToken_,hbhe);
112  if (!hbhe.isValid()) {
113  // can't find it!
114  if (!allowMissingInputs_) {
115  *hbhe; // will throw the proper exception
116  }
117  } else {
118  for (auto const& hbheItr : (HBHERecHitCollection)(*hbhe)) {
119  DetId id = (hbheItr).detid();
120  if (hbheItr.energy() > 0. )
121  edm::LogWarning("MinBias")<<" Energy = "<<hbheItr.energy();
122  theFillDetMap0_[id] += 1.;
123  theFillDetMap1_[id] += hbheItr.energy();
124  theFillDetMap2_[id] += pow(hbheItr.energy(),2);
125  theFillDetMap3_[id] += pow(hbheItr.energy(),3);
126  theFillDetMap4_[id] += pow(hbheItr.energy(),4);
127  }
128  }
129  }
130 
131  if (!hoLabel_.empty()) {
133  iEvent.getByToken(hoToken_,ho);
134  if (!ho.isValid()) {
135  // can't find it!
136  if (!allowMissingInputs_) {
137  *ho; // will throw the proper exception
138  }
139  } else {
140  for (auto const& hoItr : (HORecHitCollection)(*ho)) {
141  DetId id = hoItr.detid();
142  theFillDetMap0_[id] += 1.;
143  theFillDetMap1_[id] += hoItr.energy();
144  theFillDetMap2_[id] += pow(hoItr.energy(),2);
145  theFillDetMap3_[id] += pow(hoItr.energy(),3);
146  theFillDetMap4_[id] += pow(hoItr.energy(),4);
147  }
148  }
149  }
150 
151  if (!hfLabel_.empty()) {
153  iEvent.getByToken(hfToken_,hf);
154  if (!hf.isValid()) {
155  // can't find it!
156  if (!allowMissingInputs_) {
157  *hf; // will throw the proper exception
158  }
159  } else {
160  for (auto const hfItr : (HFRecHitCollection)(*hf)) {
161  DetId id = hfItr.detid();
162  theFillDetMap0_[id] += 1.;
163  theFillDetMap1_[id] += hfItr.energy();
164  theFillDetMap2_[id] += pow(hfItr.energy(),2);
165  theFillDetMap3_[id] += pow(hfItr.energy(),3);
166  theFillDetMap4_[id] += pow(hfItr.energy(),4);
167  }
168  }
169  }
170 
171 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< HFRecHitCollection > hfToken_
Definition: MinBias.h:51
std::string hoLabel_
Definition: MinBias.h:48
std::map< DetId, double > theFillDetMap4_
Definition: MinBias.h:67
edm::EDGetTokenT< HORecHitCollection > hoToken_
Definition: MinBias.h:50
bool allowMissingInputs_
Definition: MinBias.h:53
std::map< DetId, double > theFillDetMap3_
Definition: MinBias.h:66
edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
Definition: MinBias.h:49
bool isValid() const
Definition: HandleBase.h:74
Definition: DetId.h:18
std::map< DetId, double > theFillDetMap2_
Definition: MinBias.h:65
std::string hbheLabel_
Definition: MinBias.h:48
std::map< DetId, double > theFillDetMap0_
Definition: MinBias.h:63
std::map< DetId, double > theFillDetMap1_
Definition: MinBias.h:64
std::string hfLabel_
Definition: MinBias.h:48
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void cms::MinBias::beginJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 22 of file MinBias.cc.

References depth, eta, ieta, iphi, TFileService::make(), mom1, mom2, mom3, mom4, mydet, mysubd, myTree_, and phi.

22  {
24  myTree_ = fs->make<TTree>("RecJet","RecJet Tree");
25  myTree_->Branch("mydet", &mydet, "mydet/I");
26  myTree_->Branch("mysubd", &mysubd, "mysubd/I");
27  myTree_->Branch("depth", &depth, "depth/I");
28  myTree_->Branch("ieta", &ieta, "ieta/I");
29  myTree_->Branch("iphi", &iphi, "iphi/I");
30  myTree_->Branch("eta", &eta, "eta/F");
31  myTree_->Branch("phi", &phi, "phi/F");
32  myTree_->Branch("mom1", &mom1, "mom1/F");
33  myTree_->Branch("mom2", &mom2, "mom2/F");
34  myTree_->Branch("mom3", &mom3, "mom3/F");
35  myTree_->Branch("mom4", &mom4, "mom4/F");
36 
37 }
float mom2
Definition: MinBias.h:60
float mom1
Definition: MinBias.h:60
float mom4
Definition: MinBias.h:60
TTree * myTree_
Definition: MinBias.h:56
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
int mydet
Definition: MinBias.h:58
int ieta
Definition: MinBias.h:58
int depth
Definition: MinBias.h:58
float eta
Definition: MinBias.h:59
float phi
Definition: MinBias.h:59
float mom3
Definition: MinBias.h:60
int iphi
Definition: MinBias.h:58
int mysubd
Definition: MinBias.h:58
void cms::MinBias::beginRun ( edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 39 of file MinBias.cc.

References geo_, edm::EventSetup::get(), CaloGeometry::getValidDetIds(), DetId::Hcal, triggerObjects_cff::id, edm::ESHandle< T >::product(), theFillDetMap0_, theFillDetMap1_, theFillDetMap2_, theFillDetMap3_, and theFillDetMap4_.

39  {
41  iSetup.get<CaloGeometryRecord>().get(pG);
42  geo_ = pG.product();
43  std::vector<DetId> did = geo_->getValidDetIds();
44 
45  for (auto const& id : did) {
46  if( (id).det() == DetId::Hcal ) {
47  theFillDetMap0_[id] = 0.;
48  theFillDetMap1_[id] = 0.;
49  theFillDetMap2_[id] = 0.;
50  theFillDetMap3_[id] = 0.;
51  theFillDetMap4_[id] = 0.;
52  }
53  }
54 }
std::map< DetId, double > theFillDetMap4_
Definition: MinBias.h:67
const CaloGeometry * geo_
Definition: MinBias.h:61
std::map< DetId, double > theFillDetMap3_
Definition: MinBias.h:66
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
std::map< DetId, double > theFillDetMap2_
Definition: MinBias.h:65
std::map< DetId, double > theFillDetMap0_
Definition: MinBias.h:63
std::map< DetId, double > theFillDetMap1_
Definition: MinBias.h:64
T const * product() const
Definition: ESHandle.h:86
void cms::MinBias::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 56 of file MinBias.cc.

References depth, HcalDetId::depth(), eta, PV3DBase< T, PVType, FrameType >::eta(), geo_, HcalGeometry::getPosition(), CaloGeometry::getSubdetectorGeometry(), HcalGeometry::getValidDetIds(), DetId::Hcal, mps_fire::i, triggerObjects_cff::id, ieta, HcalDetId::ieta(), createfilelist::int, iphi, HcalDetId::iphi(), mom1, mom2, mom3, mom4, mydet, mysubd, myTree_, phi, PV3DBase< T, PVType, FrameType >::phi(), funct::pow(), theFillDetMap0_, theFillDetMap1_, theFillDetMap2_, theFillDetMap3_, and theFillDetMap4_.

56  {
57  const HcalGeometry* hgeo = static_cast<const HcalGeometry*>(geo_->getSubdetectorGeometry(DetId::Hcal,1));
58  const std::vector<DetId>& did = hgeo->getValidDetIds() ;
59  int i=0;
60  for (const auto& id : did) {
61  // if( id.det() == DetId::Hcal ) {
62  GlobalPoint pos = hgeo->getPosition(id);
63  mydet = (int)(id.det());
64  mysubd = (id.subdetId());
65  depth = HcalDetId(id).depth();
66  ieta = HcalDetId(id).ieta();
67  iphi = HcalDetId(id).iphi();
68  phi = pos.phi();
69  eta = pos.eta();
70  if ( theFillDetMap0_[id] > 0. ) {
74  2.*pow(mom2,3);
76 
77  } else {
78  mom1 = 0.; mom2 = 0.; mom3 = 0.; mom4 = 0.;
79  }
80  edm::LogWarning("MinBias")<<" Detector "<< id.rawId()<<" mydet "<<mydet
81  <<" "<<mysubd<<" "<<depth<<" "<<ieta<<" "
82  <<iphi<<" "<<pos.eta()<<" "<<pos.phi();
83  edm::LogWarning("MinBias")<<" Energy "<<mom1<<" "<<mom2<<std::endl;
84  myTree_->Fill();
85  i++;
86  // }
87  }
88  edm::LogWarning("MinBias")<<" The number of CaloDet records "<<did.size();
89  edm::LogWarning("MinBias")<<" 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  */
99 }
float mom2
Definition: MinBias.h:60
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
float mom1
Definition: MinBias.h:60
float mom4
Definition: MinBias.h:60
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
TTree * myTree_
Definition: MinBias.h:56
std::map< DetId, double > theFillDetMap4_
Definition: MinBias.h:67
int mydet
Definition: MinBias.h:58
int ieta
Definition: MinBias.h:58
const CaloGeometry * geo_
Definition: MinBias.h:61
int depth() const
get the tower depth
Definition: HcalDetId.h:166
std::map< DetId, double > theFillDetMap3_
Definition: MinBias.h:66
int ieta() const
get the cell ieta
Definition: HcalDetId.h:159
GlobalPoint getPosition(const DetId &id) const
int depth
Definition: MinBias.h:58
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
float eta
Definition: MinBias.h:59
float phi
Definition: MinBias.h:59
float mom3
Definition: MinBias.h:60
T eta() const
Definition: PV3DBase.h:76
int iphi
Definition: MinBias.h:58
std::map< DetId, double > theFillDetMap2_
Definition: MinBias.h:65
std::map< DetId, double > theFillDetMap0_
Definition: MinBias.h:63
int mysubd
Definition: MinBias.h:58
std::map< DetId, double > theFillDetMap1_
Definition: MinBias.h:64
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 53 of file MinBias.h.

Referenced by analyze(), and MinBias().

int cms::MinBias::depth
private
float cms::MinBias::eta
private
const CaloGeometry* cms::MinBias::geo_
private

Definition at line 61 of file MinBias.h.

Referenced by beginRun(), and endJob().

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

Definition at line 48 of file MinBias.h.

Referenced by analyze(), and MinBias().

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

Definition at line 49 of file MinBias.h.

Referenced by analyze(), and MinBias().

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

Definition at line 48 of file MinBias.h.

Referenced by analyze(), and MinBias().

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

Definition at line 51 of file MinBias.h.

Referenced by analyze(), and MinBias().

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

Definition at line 48 of file MinBias.h.

Referenced by analyze(), and MinBias().

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

Definition at line 50 of file MinBias.h.

Referenced by analyze(), and MinBias().

int cms::MinBias::ieta
private

Definition at line 58 of file MinBias.h.

Referenced by beginJob(), and endJob().

int cms::MinBias::iphi
private

Definition at line 58 of file MinBias.h.

Referenced by beginJob(), and endJob().

float cms::MinBias::mom1
private

Definition at line 60 of file MinBias.h.

Referenced by beginJob(), and endJob().

float cms::MinBias::mom2
private

Definition at line 60 of file MinBias.h.

Referenced by beginJob(), and endJob().

float cms::MinBias::mom3
private

Definition at line 60 of file MinBias.h.

Referenced by beginJob(), and endJob().

float cms::MinBias::mom4
private

Definition at line 60 of file MinBias.h.

Referenced by beginJob(), and endJob().

int cms::MinBias::mydet
private

Definition at line 58 of file MinBias.h.

Referenced by beginJob(), and endJob().

int cms::MinBias::mysubd
private

Definition at line 58 of file MinBias.h.

Referenced by beginJob(), and endJob().

TTree* cms::MinBias::myTree_
private

Definition at line 56 of file MinBias.h.

Referenced by beginJob(), and endJob().

float cms::MinBias::occup
private

Definition at line 60 of file MinBias.h.

float cms::MinBias::phi
private
std::map<DetId,double> cms::MinBias::theFillDetMap0_
private

Definition at line 63 of file MinBias.h.

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

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

Definition at line 64 of file MinBias.h.

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

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

Definition at line 65 of file MinBias.h.

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

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

Definition at line 66 of file MinBias.h.

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

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

Definition at line 67 of file MinBias.h.

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