CMS 3D CMS Logo

dEdxHitAnalyzer.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Loic Quertenmont
5  */
7 
12 
14 
18 
19 #include <string>
20 #include "TMath.h"
21 
23  : fullconf_( iConfig )
24  , conf_ (fullconf_.getParameter<edm::ParameterSet>("dEdxParameters") )
25  , doAllPlots_ ( conf_.getParameter<bool>("doAllPlots") )
26  , doDeDxPlots_ ( conf_.getParameter<bool>("doDeDxPlots") )
27  , genTriggerEventFlag_( new GenericTriggerEventFlag(conf_.getParameter<edm::ParameterSet>("genericTriggerEventPSet"),consumesCollector(), *this) )
28 {
29 
31  trackToken_ = consumes<reco::TrackCollection>(trackInputTag_);
32 
33  dEdxInputList_ = conf_.getParameter<std::vector<std::string> >("deDxHitProducers");
34  for (auto const& tag : dEdxInputList_) {
35  dEdxTokenList_.push_back(consumes<reco::DeDxHitInfoAss>(edm::InputTag(tag) ) );
36  }
37 
38  // parameters from the configuration
39  MEFolderName = conf_.getParameter<std::string>("FolderName");
40 
41  dEdxNHitBin = conf_.getParameter<int>( "dEdxNHitBin");
42  dEdxNHitMin = conf_.getParameter<double>("dEdxNHitMin");
43  dEdxNHitMax = conf_.getParameter<double>("dEdxNHitMax");
44 
45  dEdxStripBin = conf_.getParameter<int>( "dEdxStripBin");
46  dEdxStripMin = conf_.getParameter<double>("dEdxStripMin");
47  dEdxStripMax = conf_.getParameter<double>("dEdxStripMax");
48 
49  dEdxPixelBin = conf_.getParameter<int>( "dEdxPixelBin");
50  dEdxPixelMin = conf_.getParameter<double>("dEdxPixelMin");
51  dEdxPixelMax = conf_.getParameter<double>("dEdxPixelMax");
52 
53  dEdxHarm2Bin = conf_.getParameter<int>( "dEdxHarm2Bin");
54  dEdxHarm2Min = conf_.getParameter<double>("dEdxHarm2Min");
55  dEdxHarm2Max = conf_.getParameter<double>("dEdxHarm2Max");
56 }
57 
59 {
60 
62 }
63 
64 // -- BeginRun
65 //---------------------------------------------------------------------------------//
66 void dEdxHitAnalyzer::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
67 {
68  // Initialize the GenericTriggerEventFlag
69  if ( genTriggerEventFlag_->on() ) genTriggerEventFlag_->initRun( iRun, iSetup );
70 }
71 
72 
74  edm::Run const & iRun,
75  edm::EventSetup const & iSetup )
76 {
77 
79 
80  // book the Hit Property histograms
81  // ---------------------------------------------------------------------------------//
82 
83  if ( doDeDxPlots_ || doAllPlots_ ){
84  for(unsigned int i=0;i<dEdxInputList_.size();i++){
86  dEdxMEsVector.push_back(dEdxMEs() );
87 
88  histname = "Strip_dEdxPerCluster_";
90  dEdxMEsVector[i].ME_StripHitDeDx->setAxisTitle("dEdx of on-track strip cluster (ADC)");
91  dEdxMEsVector[i].ME_StripHitDeDx->setAxisTitle("Number of Strip clusters", 2);
92 
93  histname = "Pixel_dEdxPerCluster_";
95  dEdxMEsVector[i].ME_PixelHitDeDx->setAxisTitle("dEdx of on-track pixel cluster (ADC)");
96  dEdxMEsVector[i].ME_PixelHitDeDx->setAxisTitle("Number of Pixel clusters", 2);
97 
98  histname = "NumberOfdEdxHitsPerTrack_";
100  dEdxMEsVector[i].ME_NHitDeDx->setAxisTitle("Number of dEdxHits per Track");
101  dEdxMEsVector[i].ME_NHitDeDx->setAxisTitle("Number of Tracks", 2);
102 
103  histname = "Harm2_dEdxPerTrack_";
105  dEdxMEsVector[i].ME_Harm2DeDx->setAxisTitle("Harmonic2 dEdx estimator for each Track");
106  dEdxMEsVector[i].ME_Harm2DeDx->setAxisTitle("Number of Tracks", 2);
107  }
108  }
109 }
110 
112  if(!dedxHits)return -1;
113  std::vector<double> vect;
114  for(unsigned int h=0;h<dedxHits->size();h++){
115  DetId detid(dedxHits->detId(h));
116  double Norm = (detid.subdetId()<3)?3.61e-06:3.61e-06*265;
117  double ChargeOverPathlength = Norm*dedxHits->charge(h)/dedxHits->pathlength(h);
118  vect.push_back(ChargeOverPathlength); //save charge
119  }
120 
121  int size = vect.size();
122  if(size<=0)return -1;
123  double result=0;
124  double expo = -2;
125  for(int i = 0; i< size; i ++){
126  result+=pow(vect[i],expo);
127  }
128  return pow(result/size,1./expo);
129 }
130 
131 
132 // -- Analyse
133 // ---------------------------------------------------------------------------------//
135 {
136 
137  // Filter out events if Trigger Filtering is requested
138  if (genTriggerEventFlag_->on()&& ! genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
139 
140 
141  if ( doDeDxPlots_ || doAllPlots_ ){
142  edm::Handle<reco::TrackCollection> trackCollectionHandle;
143  iEvent.getByToken(trackToken_, trackCollectionHandle );
144  if(!trackCollectionHandle.isValid())return;
145 
146  for(unsigned int i=0;i<dEdxInputList_.size();i++){
147  edm::Handle<reco::DeDxHitInfoAss> dEdxObjectHandle;
148  iEvent.getByToken(dEdxTokenList_[i], dEdxObjectHandle );
149  if(!dEdxObjectHandle.isValid())continue;
150 
151  for(unsigned int t=0; t<trackCollectionHandle->size(); t++){
152  reco::TrackRef track = reco::TrackRef( trackCollectionHandle, t );
153 
154  if(track->quality(reco::TrackBase::highPurity) ) {
155  const reco::DeDxHitInfo* dedxHits = NULL;
156  if(!track.isNull()) {
157  reco::DeDxHitInfoRef dedxHitsRef = (*dEdxObjectHandle)[track];
158  if(!dedxHitsRef.isNull())dedxHits = &(*dedxHitsRef);
159  }
160  if(!dedxHits)continue;
161 
162  for(unsigned int h=0;h<dedxHits->size();h++){
163  DetId detid(dedxHits->detId(h));
164  if(detid.subdetId()>=3)dEdxMEsVector[i].ME_StripHitDeDx ->Fill(dedxHits->charge(h));
165  if(detid.subdetId()<3 )dEdxMEsVector[i].ME_PixelHitDeDx ->Fill(dedxHits->charge(h));
166  }
167  dEdxMEsVector[i].ME_NHitDeDx->Fill(dedxHits->size());
168  dEdxMEsVector[i].ME_Harm2DeDx->Fill(harmonic2(dedxHits));
169  }
170  }
171  }
172  }
173 }
174 
175 
176 
177 
178 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
179 void
181  //The following says we do not know what parameters are allowed so do no validation
182  // Please change this to state exactly what you do use, even if it is no parameters
184  desc.setUnknown();
185  descriptions.addDefault(desc);
186 }
187 
size
Write out results.
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
float charge(size_t i) const
Definition: DeDxHitInfo.h:42
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
#define NULL
Definition: scimark2.h:8
edm::EDGetTokenT< reco::TrackCollection > trackToken_
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
std::vector< std::string > dEdxInputList_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::string histname
int iEvent
Definition: GenABIO.cc:230
dEdxHitAnalyzer(const edm::ParameterSet &)
void addDefault(ParameterSetDescription const &psetDescription)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
edm::InputTag trackInputTag_
bool accept(const edm::Event &event, const edm::EventSetup &setup)
To be called from analyze/filter() methods.
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
GenericTriggerEventFlag * genTriggerEventFlag_
bool isValid() const
Definition: HandleBase.h:74
std::vector< edm::EDGetTokenT< reco::DeDxHitInfoAss > > dEdxTokenList_
bool isNull() const
Checks for null.
Definition: Ref.h:249
double harmonic2(const reco::DeDxHitInfo *dedxHits)
Definition: DetId.h:18
virtual void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
size_t size() const
Definition: DeDxHitInfo.h:41
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
HLT enums.
float pathlength(size_t i) const
Definition: DeDxHitInfo.h:43
void initRun(const edm::Run &run, const edm::EventSetup &setup)
To be called from beginRun() methods.
edm::ParameterSet conf_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::string MEFolderName
DetId detId(size_t i) const
Definition: DeDxHitInfo.h:44
Definition: Run.h:42
std::vector< dEdxMEs > dEdxMEsVector