CMS 3D CMS Logo

TkAlCaRecoMonitor.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  */
5 
7 
15 
20 
22 
23 #include <string>
24 //#include <sstream>
25 #include "TLorentzVector.h"
26 
28  conf_ = iConfig;
29  trackProducer_ = consumes<reco::TrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer"));
31  consumes<reco::TrackCollection>(conf_.getParameter<edm::InputTag>("ReferenceTrackProducer"));
32  jetCollection_ = mayConsume<reco::CaloJetCollection>(conf_.getParameter<edm::InputTag>("CaloJetCollection"));
33 }
34 
36 
38  std::string histname; // for naming the histograms according to algorithm used
39 
40  std::string AlgoName = conf_.getParameter<std::string>("AlgoName");
41  std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
42 
43  daughterMass_ = conf_.getParameter<double>("daughterMass");
44 
45  maxJetPt_ = conf_.getParameter<double>("maxJetPt");
46 
47  iBooker.setCurrentFolder(MEFolderName + "/TkAlignmentSpecific");
48  fillInvariantMass_ = conf_.getParameter<bool>("fillInvariantMass");
49  runsOnReco_ = conf_.getParameter<bool>("runsOnReco");
50  useSignedR_ = conf_.getParameter<bool>("useSignedR");
51  fillRawIdMap_ = conf_.getParameter<bool>("fillRawIdMap");
52 
53  //
54  unsigned int MassBin = conf_.getParameter<unsigned int>("MassBin");
55  double MassMin = conf_.getParameter<double>("MassMin");
56  double MassMax = conf_.getParameter<double>("MassMax");
57 
58  if (fillInvariantMass_) {
59  histname = "InvariantMass_";
60  invariantMass_ = iBooker.book1D(histname + AlgoName, histname + AlgoName, MassBin, MassMin, MassMax);
61  invariantMass_->setAxisTitle("invariant Mass / GeV");
62  } else {
63  invariantMass_ = nullptr;
64  }
65 
66  unsigned int TrackPtPositiveBin = conf_.getParameter<unsigned int>("TrackPtBin");
67  double TrackPtPositiveMin = conf_.getParameter<double>("TrackPtMin");
68  double TrackPtPositiveMax = conf_.getParameter<double>("TrackPtMax");
69 
70  histname = "TrackPtPositive_";
71  TrackPtPositive_ = iBooker.book1D(
72  histname + AlgoName, histname + AlgoName, TrackPtPositiveBin, TrackPtPositiveMin, TrackPtPositiveMax);
73  TrackPtPositive_->setAxisTitle("p_{T} of tracks charge > 0");
74 
75  unsigned int TrackPtNegativeBin = conf_.getParameter<unsigned int>("TrackPtBin");
76  double TrackPtNegativeMin = conf_.getParameter<double>("TrackPtMin");
77  double TrackPtNegativeMax = conf_.getParameter<double>("TrackPtMax");
78 
79  histname = "TrackPtNegative_";
80  TrackPtNegative_ = iBooker.book1D(
81  histname + AlgoName, histname + AlgoName, TrackPtNegativeBin, TrackPtNegativeMin, TrackPtNegativeMax);
82  TrackPtNegative_->setAxisTitle("p_{T} of tracks charge < 0");
83 
84  histname = "TrackQuality_";
85  TrackQuality_ = iBooker.book1D(
86  histname + AlgoName, histname + AlgoName, reco::TrackBase::qualitySize, -0.5, reco::TrackBase::qualitySize - 0.5);
87  TrackQuality_->setAxisTitle("quality");
88  for (int i = 0; i < reco::TrackBase::qualitySize; ++i) {
89  TrackQuality_->getTH1()->GetXaxis()->SetBinLabel(
91  }
92 
93  unsigned int SumChargeBin = conf_.getParameter<unsigned int>("SumChargeBin");
94  double SumChargeMin = conf_.getParameter<double>("SumChargeMin");
95  double SumChargeMax = conf_.getParameter<double>("SumChargeMax");
96 
97  histname = "SumCharge_";
98  sumCharge_ = iBooker.book1D(histname + AlgoName, histname + AlgoName, SumChargeBin, SumChargeMin, SumChargeMax);
99  sumCharge_->setAxisTitle("#SigmaCharge");
100 
101  unsigned int TrackCurvatureBin = conf_.getParameter<unsigned int>("TrackCurvatureBin");
102  double TrackCurvatureMin = conf_.getParameter<double>("TrackCurvatureMin");
103  double TrackCurvatureMax = conf_.getParameter<double>("TrackCurvatureMax");
104 
105  histname = "TrackCurvature_";
107  iBooker.book1D(histname + AlgoName, histname + AlgoName, TrackCurvatureBin, TrackCurvatureMin, TrackCurvatureMax);
108  TrackCurvature_->setAxisTitle("#kappa track");
109 
110  if (runsOnReco_) {
111  unsigned int JetPtBin = conf_.getParameter<unsigned int>("JetPtBin");
112  double JetPtMin = conf_.getParameter<double>("JetPtMin");
113  double JetPtMax = conf_.getParameter<double>("JetPtMax");
114 
115  histname = "JetPt_";
116  jetPt_ = iBooker.book1D(histname + AlgoName, histname + AlgoName, JetPtBin, JetPtMin, JetPtMax);
117  jetPt_->setAxisTitle("jet p_{T} / GeV");
118 
119  unsigned int MinJetDeltaRBin = conf_.getParameter<unsigned int>("MinJetDeltaRBin");
120  double MinJetDeltaRMin = conf_.getParameter<double>("MinJetDeltaRMin");
121  double MinJetDeltaRMax = conf_.getParameter<double>("MinJetDeltaRMax");
122 
123  histname = "MinJetDeltaR_";
124  minJetDeltaR_ =
125  iBooker.book1D(histname + AlgoName, histname + AlgoName, MinJetDeltaRBin, MinJetDeltaRMin, MinJetDeltaRMax);
126  minJetDeltaR_->setAxisTitle("minimal Jet #DeltaR / rad");
127  } else {
128  jetPt_ = nullptr;
129  minJetDeltaR_ = nullptr;
130  }
131 
132  unsigned int MinTrackDeltaRBin = conf_.getParameter<unsigned int>("MinTrackDeltaRBin");
133  double MinTrackDeltaRMin = conf_.getParameter<double>("MinTrackDeltaRMin");
134  double MinTrackDeltaRMax = conf_.getParameter<double>("MinTrackDeltaRMax");
135 
136  histname = "MinTrackDeltaR_";
138  iBooker.book1D(histname + AlgoName, histname + AlgoName, MinTrackDeltaRBin, MinTrackDeltaRMin, MinTrackDeltaRMax);
139  minTrackDeltaR_->setAxisTitle("minimal Track #DeltaR / rad");
140 
141  unsigned int TrackEfficiencyBin = conf_.getParameter<unsigned int>("TrackEfficiencyBin");
142  double TrackEfficiencyMin = conf_.getParameter<double>("TrackEfficiencyMin");
143  double TrackEfficiencyMax = conf_.getParameter<double>("TrackEfficiencyMax");
144 
145  histname = "AlCaRecoTrackEfficiency_";
147  histname + AlgoName, histname + AlgoName, TrackEfficiencyBin, TrackEfficiencyMin, TrackEfficiencyMax);
148  Labels l_tp, l_rtp;
149  labelsForToken(referenceTrackProducer_, l_rtp);
150  labelsForToken(trackProducer_, l_tp);
151  AlCaRecoTrackEfficiency_->setAxisTitle("n(" + std::string(l_tp.module) + ") / n(" + std::string(l_rtp.module) + ")");
152 
153  int zBin = conf_.getParameter<unsigned int>("HitMapsZBin"); // 300
154  double zMax = conf_.getParameter<double>("HitMapZMax"); // 300.0; //cm
155 
156  int rBin = conf_.getParameter<unsigned int>("HitMapsRBin"); // 120;
157  double rMax = conf_.getParameter<double>("HitMapRMax"); // 120.0; //cm
158 
159  histname = "Hits_ZvsR_";
160  double rMin = 0.0;
161  if (useSignedR_)
162  rMin = -rMax;
163 
164  Hits_ZvsR_ = iBooker.book2D(histname + AlgoName, histname + AlgoName, zBin, -zMax, zMax, rBin, rMin, rMax);
165 
166  histname = "Hits_XvsY_";
167  Hits_XvsY_ = iBooker.book2D(histname + AlgoName, histname + AlgoName, rBin, -rMax, rMax, rBin, -rMax, rMax);
168 
169  if (fillRawIdMap_) {
170  histname = "Hits_perDetId_";
171 
172  // leads to differences in axsis between samples??
173  // int nModules = binByRawId_.size();
174  // Hits_perDetId_ = iBooker.book1D(histname+AlgoName, histname+AlgoName,
175  // nModules, static_cast<double>(nModules) -0.5,
176  // static_cast<double>(nModules) -0.5);
177  Hits_perDetId_ = iBooker.book1D(histname + AlgoName, histname + AlgoName, 16601, -0.5, 16600.5);
178  Hits_perDetId_->setAxisTitle("rawId Bins");
179 
181  // std::stringstream binLabel;
182  // for( std::map<int,int>::iterator it = binByRawId_.begin(); it !=
183  // binByRawId_.end(); ++it ){
184  // binLabel.str() = "";
185  // binLabel << (*it).first;
186  // Hits_perDetId_->getTH1()->GetXaxis()->SetBinLabel( (*it).second +1,
187  // binLabel.str().c_str());
188  // }
189  }
190 }
191 //
192 // -- Analyse
193 //
196  iEvent.getByToken(trackProducer_, trackCollection);
197  if (!trackCollection.isValid()) {
198  edm::LogError("Alignment") << "invalid trackcollection encountered!";
199  return;
200  }
201 
202  edm::Handle<reco::TrackCollection> referenceTrackCollection;
203  iEvent.getByToken(referenceTrackProducer_, referenceTrackCollection);
204  if (!trackCollection.isValid()) {
205  edm::LogError("Alignment") << "invalid reference track-collection encountered!";
206  return;
207  }
208 
210  iSetup.get<TrackerDigiGeometryRecord>().get(geometry);
211  if (!geometry.isValid()) {
212  edm::LogError("Alignment") << "invalid geometry found in event setup!";
213  }
214 
216  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
217  if (!magneticField.isValid()) {
218  edm::LogError("Alignment") << "invalid magnetic field configuration encountered!";
219  return;
220  }
221 
223  if (runsOnReco_) {
224  iEvent.getByToken(jetCollection_, jets);
225  if (!jets.isValid()) {
226  edm::LogError("Alignment") << "no jet collection found in event!";
227  }
228  }
229  // fill only once - not yet in beginJob since no access to geometry
230  if (fillRawIdMap_ && binByRawId_.empty())
231  this->fillRawIdMap(*geometry);
232 
233  AlCaRecoTrackEfficiency_->Fill(static_cast<double>((*trackCollection).size()) / (*referenceTrackCollection).size());
234 
235  double sumOfCharges = 0;
236  for (reco::TrackCollection::const_iterator track = (*trackCollection).begin(); track < (*trackCollection).end();
237  ++track) {
238  double dR = 0;
239  if (runsOnReco_) {
240  double minJetDeltaR = 10; // some number > 2pi
241  for (reco::CaloJetCollection::const_iterator itJet = jets->begin(); itJet != jets->end(); ++itJet) {
242  jetPt_->Fill((*itJet).pt());
243  dR = deltaR((*track), (*itJet));
244  if ((*itJet).pt() > maxJetPt_ && dR < minJetDeltaR)
245  minJetDeltaR = dR;
246 
247  // edm::LogInfo("Alignment") <<"> isolated: "<< isolated << " jetPt "<<
248  // (*itJet).pt() <<" deltaR: "<< deltaR(*(*it),(*itJet)) ;
249  }
250  minJetDeltaR_->Fill(minJetDeltaR);
251  }
252 
253  double minTrackDeltaR = 10; // some number > 2pi
254  for (reco::TrackCollection::const_iterator track2 = (*trackCollection).begin(); track2 < (*trackCollection).end();
255  ++track2) {
256  dR = deltaR((*track), (*track2));
257  if (dR < minTrackDeltaR && dR > 1e-6)
258  minTrackDeltaR = dR;
259  }
260 
261  for (int i = 0; i < reco::TrackBase::qualitySize; ++i) {
262  if ((*track).quality(reco::TrackBase::TrackQuality(i))) {
263  TrackQuality_->Fill(i);
264  }
265  }
266 
267  GlobalPoint gPoint((*track).vx(), (*track).vy(), (*track).vz());
268  double B = magneticField->inTesla(gPoint).z();
269  double curv = -(*track).charge() * 0.002998 * B / (*track).pt();
270  TrackCurvature_->Fill(curv);
271 
272  if ((*track).charge() > 0)
273  TrackPtPositive_->Fill((*track).pt());
274  if ((*track).charge() < 0)
275  TrackPtNegative_->Fill((*track).pt());
276 
277  minTrackDeltaR_->Fill(minTrackDeltaR);
278  fillHitmaps(*track, *geometry);
279  sumOfCharges += (*track).charge();
280  }
281 
282  sumCharge_->Fill(sumOfCharges);
283 
284  if (fillInvariantMass_) {
285  if ((*trackCollection).size() == 2) {
286  TLorentzVector track0(
287  (*trackCollection).at(0).px(),
288  (*trackCollection).at(0).py(),
289  (*trackCollection).at(0).pz(),
290  sqrt(((*trackCollection).at(0).p() * (*trackCollection).at(0).p()) + daughterMass_ * daughterMass_));
291  TLorentzVector track1(
292  (*trackCollection).at(1).px(),
293  (*trackCollection).at(1).py(),
294  (*trackCollection).at(1).pz(),
295  sqrt(((*trackCollection).at(1).p() * (*trackCollection).at(1).p()) + daughterMass_ * daughterMass_));
296  TLorentzVector mother = track0 + track1;
297 
298  invariantMass_->Fill(mother.M());
299  } else {
300  edm::LogInfo("Alignment") << "wrong number of tracks trackcollection encountered: " << (*trackCollection).size();
301  }
302  }
303 }
304 
306  for (trackingRecHit_iterator iHit = track.recHitsBegin(); iHit != track.recHitsEnd(); ++iHit) {
307  if ((*iHit)->isValid()) {
308  const TrackingRecHit *hit = (*iHit);
309  const DetId geoId(hit->geographicalId());
310  const GeomDet *gd = geometry.idToDet(geoId);
311  // since 2_1_X local hit positions are transient. taking center of the hit
312  // module for now. The alternative would be the coarse estimation or a
313  // refit.
314  // const GlobalPoint globP( gd->toGlobal( hit->localPosition() ) );
315  const GlobalPoint globP(gd->toGlobal(Local3DPoint(0., 0., 0.)));
316  double r = sqrt(globP.x() * globP.x() + globP.y() * globP.y());
317  if (useSignedR_)
318  r *= globP.y() / fabs(globP.y());
319 
320  Hits_ZvsR_->Fill(globP.z(), r);
321  Hits_XvsY_->Fill(globP.x(), globP.y());
322 
323  if (fillRawIdMap_)
324  Hits_perDetId_->Fill(binByRawId_[geoId.rawId()]);
325  }
326  }
327 }
328 
330  std::vector<int> sortedRawIds;
331  for (std::vector<DetId>::const_iterator iDetId = geometry.detUnitIds().begin(); iDetId != geometry.detUnitIds().end();
332  ++iDetId) {
333  sortedRawIds.push_back((*iDetId).rawId());
334  }
335  std::sort(sortedRawIds.begin(), sortedRawIds.end());
336 
337  int i = 0;
338  for (std::vector<int>::iterator iRawId = sortedRawIds.begin(); iRawId != sortedRawIds.end(); ++iRawId) {
339  binByRawId_[*iRawId] = i;
340  ++i;
341  }
342 }
343 
const DetIdContainer & detUnitIds() const override
Returm a vector of all GeomDetUnit DetIds.
std::map< int, int > binByRawId_
T getParameter(std::string const &) const
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:571
void analyze(const edm::Event &, const edm::EventSetup &) override
void fillRawIdMap(const TrackerGeometry &geometry)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
TrackQuality
track quality
Definition: TrackBase.h:151
MonitorElement * TrackQuality_
TH1 * getTH1() const
~TkAlCaRecoMonitor() override
MonitorElement * minJetDeltaR_
MonitorElement * Hits_XvsY_
MonitorElement * Hits_ZvsR_
void Fill(long long x)
MonitorElement * AlCaRecoTrackEfficiency_
MonitorElement * sumCharge_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
edm::EDGetTokenT< reco::TrackCollection > trackProducer_
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
T z() const
Definition: PV3DBase.h:64
MonitorElement * invariantMass_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * TrackPtNegative_
void fillHitmaps(const reco::Track &track, const TrackerGeometry &geometry)
edm::EDGetTokenT< reco::CaloJetCollection > jetCollection_
static const std::string B
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:109
bool isValid() const
Definition: HandleBase.h:74
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Definition: DetId.h:18
MonitorElement * Hits_perDetId_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
MonitorElement * TrackCurvature_
MonitorElement * jetPt_
edm::ParameterSet conf_
ESHandle< TrackerGeometry > geometry
T get() const
Definition: EventSetup.h:71
MonitorElement * minTrackDeltaR_
const TrackerGeomDet * idToDet(DetId) const override
DetId geographicalId() const
bool isValid() const
Definition: ESHandle.h:44
edm::EDGetTokenT< reco::TrackCollection > referenceTrackProducer_
MonitorElement * TrackPtPositive_
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Definition: Run.h:45
TkAlCaRecoMonitor(const edm::ParameterSet &)
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:114