CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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"));
30  referenceTrackProducer_ = consumes<reco::TrackCollection>(conf_.getParameter<edm::InputTag>("ReferenceTrackProducer"));
31  jetCollection_ = mayConsume<reco::CaloJetCollection>(conf_.getParameter<edm::InputTag>("CaloJetCollection"));
32 }
33 
35 
37 
38  std::string histname; //for naming the histograms according to algorithm used
39 
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 
59  histname = "InvariantMass_";
60  invariantMass_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, MassBin, MassMin, MassMax);
61  invariantMass_->setAxisTitle("invariant Mass / GeV");
62  }else{
63  invariantMass_ = 0;
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(histname+AlgoName, histname+AlgoName, TrackPtPositiveBin, TrackPtPositiveMin, TrackPtPositiveMax);
72  TrackPtPositive_->setAxisTitle("p_{T} of tracks charge > 0");
73 
74  unsigned int TrackPtNegativeBin = conf_.getParameter<unsigned int>("TrackPtBin");
75  double TrackPtNegativeMin = conf_.getParameter<double>("TrackPtMin");
76  double TrackPtNegativeMax = conf_.getParameter<double>("TrackPtMax");
77 
78  histname = "TrackPtNegative_";
79  TrackPtNegative_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, TrackPtNegativeBin, TrackPtNegativeMin, TrackPtNegativeMax);
80  TrackPtNegative_->setAxisTitle("p_{T} of tracks charge < 0");
81 
82  histname = "TrackQuality_";
83  TrackQuality_ = iBooker.book1D(histname+AlgoName, histname+AlgoName,
85  TrackQuality_->setAxisTitle("quality");
86  for ( int i = 0; i<reco::TrackBase::qualitySize; ++i){
87  TrackQuality_->getTH1()->GetXaxis()->SetBinLabel(i+1,
89  }
90 
91  unsigned int SumChargeBin = conf_.getParameter<unsigned int>("SumChargeBin");
92  double SumChargeMin = conf_.getParameter<double>("SumChargeMin");
93  double SumChargeMax = conf_.getParameter<double>("SumChargeMax");
94 
95  histname = "SumCharge_";
96  sumCharge_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, SumChargeBin, SumChargeMin, SumChargeMax);
97  sumCharge_->setAxisTitle("#SigmaCharge");
98 
99  unsigned int TrackCurvatureBin = conf_.getParameter<unsigned int>("TrackCurvatureBin");
100  double TrackCurvatureMin = conf_.getParameter<double>("TrackCurvatureMin");
101  double TrackCurvatureMax = conf_.getParameter<double>("TrackCurvatureMax");
102 
103  histname = "TrackCurvature_";
104  TrackCurvature_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, TrackCurvatureBin, TrackCurvatureMin, TrackCurvatureMax);
105  TrackCurvature_->setAxisTitle("#kappa track");
106 
107  if( runsOnReco_ ){
108 
109  unsigned int JetPtBin = conf_.getParameter<unsigned int>("JetPtBin");
110  double JetPtMin = conf_.getParameter<double>("JetPtMin");
111  double JetPtMax = conf_.getParameter<double>("JetPtMax");
112 
113  histname = "JetPt_";
114  jetPt_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, JetPtBin, JetPtMin, JetPtMax);
115  jetPt_->setAxisTitle("jet p_{T} / GeV");
116 
117  unsigned int MinJetDeltaRBin = conf_.getParameter<unsigned int>("MinJetDeltaRBin");
118  double MinJetDeltaRMin = conf_.getParameter<double>("MinJetDeltaRMin");
119  double MinJetDeltaRMax = conf_.getParameter<double>("MinJetDeltaRMax");
120 
121  histname = "MinJetDeltaR_";
122  minJetDeltaR_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, MinJetDeltaRBin, MinJetDeltaRMin, MinJetDeltaRMax);
123  minJetDeltaR_->setAxisTitle("minimal Jet #DeltaR / rad");
124  }else{
125  jetPt_ = NULL;
127  }
128 
129  unsigned int MinTrackDeltaRBin = conf_.getParameter<unsigned int>("MinTrackDeltaRBin");
130  double MinTrackDeltaRMin = conf_.getParameter<double>("MinTrackDeltaRMin");
131  double MinTrackDeltaRMax = conf_.getParameter<double>("MinTrackDeltaRMax");
132 
133  histname = "MinTrackDeltaR_";
134  minTrackDeltaR_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, MinTrackDeltaRBin, MinTrackDeltaRMin, MinTrackDeltaRMax);
135  minTrackDeltaR_->setAxisTitle("minimal Track #DeltaR / rad");
136 
137  unsigned int TrackEfficiencyBin = conf_.getParameter<unsigned int>("TrackEfficiencyBin");
138  double TrackEfficiencyMin = conf_.getParameter<double>("TrackEfficiencyMin");
139  double TrackEfficiencyMax = conf_.getParameter<double>("TrackEfficiencyMax");
140 
141  histname = "AlCaRecoTrackEfficiency_";
142  AlCaRecoTrackEfficiency_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, TrackEfficiencyBin, TrackEfficiencyMin, TrackEfficiencyMax);
143  Labels l_tp, l_rtp;
147 
148  int zBin = conf_.getParameter<unsigned int>("HitMapsZBin"); //300
149  double zMax = conf_.getParameter<double>("HitMapZMax"); //300.0; //cm
150 
151  int rBin = conf_.getParameter<unsigned int>("HitMapsRBin");//120;
152  double rMax = conf_.getParameter<double>("HitMapRMax"); //120.0; //cm
153 
154  histname = "Hits_ZvsR_";
155  double rMin = 0.0;
156  if( useSignedR_ )
157  rMin = -rMax;
158 
159  Hits_ZvsR_ = iBooker.book2D(histname+AlgoName, histname+AlgoName, zBin, -zMax, zMax, rBin, rMin, rMax);
160 
161  histname = "Hits_XvsY_";
162  Hits_XvsY_ = iBooker.book2D(histname+AlgoName, histname+AlgoName, rBin, -rMax, rMax, rBin, -rMax, rMax);
163 
164  if( fillRawIdMap_){
165  histname = "Hits_perDetId_";
166 
167  //leads to differences in axsis between samples??
168  //int nModules = binByRawId_.size();
169  //Hits_perDetId_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, nModules, static_cast<double>(nModules) -0.5, static_cast<double>(nModules) -0.5);
170  Hits_perDetId_ = iBooker.book1D(histname+AlgoName, histname+AlgoName, 16601,-0.5, 16600.5 );
171  Hits_perDetId_->setAxisTitle("rawId Bins");
172 
174  // std::stringstream binLabel;
175  // for( std::map<int,int>::iterator it = binByRawId_.begin(); it != binByRawId_.end(); ++it ){
176  // binLabel.str() = "";
177  // binLabel << (*it).first;
178  // Hits_perDetId_->getTH1()->GetXaxis()->SetBinLabel( (*it).second +1, binLabel.str().c_str());
179  // }
180  }
181 }
182 //
183 // -- Analyse
184 //
186 
188  iEvent.getByToken(trackProducer_, trackCollection);
189  if (!trackCollection.isValid()){
190  edm::LogError("Alignment")<<"invalid trackcollection encountered!";
191  return;
192  }
193 
194  edm::Handle<reco::TrackCollection> referenceTrackCollection;
195  iEvent.getByToken(referenceTrackProducer_, referenceTrackCollection);
196  if (!trackCollection.isValid()){
197  edm::LogError("Alignment")<<"invalid reference track-collection encountered!";
198  return;
199  }
200 
202  iSetup.get<TrackerDigiGeometryRecord>().get(geometry);
203  if(! geometry.isValid()){
204  edm::LogError("Alignment")<<"invalid geometry found in event setup!";
205  }
206 
208  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
209  if (!magneticField.isValid()){
210  edm::LogError("Alignment")<<"invalid magnetic field configuration encountered!";
211  return;
212  }
213 
215  if(runsOnReco_){
216  iEvent.getByToken(jetCollection_, jets);
217  if(! jets.isValid()){
218  edm::LogError("Alignment")<<"no jet collection found in event!";
219  }
220  }
221  // fill only once - not yet in beginJob since no access to geometry
222  if (fillRawIdMap_ && binByRawId_.empty()) this->fillRawIdMap(*geometry);
223 
224  AlCaRecoTrackEfficiency_->Fill( static_cast<double>((*trackCollection).size()) / (*referenceTrackCollection).size() );
225 
226  double sumOfCharges = 0;
227  for( reco::TrackCollection::const_iterator track = (*trackCollection).begin(); track < (*trackCollection).end(); ++track ){
228  double dR = 0;
229  if(runsOnReco_){
230  double minJetDeltaR = 10; // some number > 2pi
231  for(reco::CaloJetCollection::const_iterator itJet = jets->begin(); itJet != jets->end() ; ++itJet) {
232  jetPt_->Fill( (*itJet).pt() );
233  dR = deltaR((*track),(*itJet));
234  if((*itJet).pt() > maxJetPt_ && dR < minJetDeltaR)
235  minJetDeltaR = dR;
236 
237  //edm::LogInfo("Alignment") <<"> isolated: "<< isolated << " jetPt "<< (*itJet).pt() <<" deltaR: "<< deltaR(*(*it),(*itJet)) ;
238  }
239  minJetDeltaR_->Fill( minJetDeltaR );
240  }
241 
242  double minTrackDeltaR = 10; // some number > 2pi
243  for( reco::TrackCollection::const_iterator track2 = (*trackCollection).begin(); track2 < (*trackCollection).end(); ++track2 ){
244  dR = deltaR((*track),(*track2));
245  if(dR < minTrackDeltaR && dR > 1e-6)
246  minTrackDeltaR = dR;
247  }
248 
249  for ( int i = 0; i<reco::TrackBase::qualitySize; ++i){
250  if( (*track).quality( reco::TrackBase::TrackQuality(i) ) ){
251  TrackQuality_->Fill( i );
252  }
253  }
254 
255  GlobalPoint gPoint((*track).vx(), (*track).vy(), (*track).vz());
256  double B = magneticField->inTesla(gPoint).z();
257  double curv = -(*track).charge()*0.002998*B/(*track).pt();
258  TrackCurvature_->Fill( curv );
259 
260  if( (*track).charge() > 0 )
261  TrackPtPositive_->Fill( (*track).pt() );
262  if( (*track).charge() < 0 )
263  TrackPtNegative_->Fill( (*track).pt() );
264 
265  minTrackDeltaR_->Fill( minTrackDeltaR );
266  fillHitmaps( *track, *geometry );
267  sumOfCharges += (*track).charge();
268  }
269 
270  sumCharge_->Fill( sumOfCharges );
271 
272  if(fillInvariantMass_){
273  if((*trackCollection).size() == 2){
274  TLorentzVector track0((*trackCollection).at(0).px(),(*trackCollection).at(0).py(),(*trackCollection).at(0).pz(),
275  sqrt(((*trackCollection).at(0).p()*(*trackCollection).at(0).p())+daughterMass_*daughterMass_));
276  TLorentzVector track1((*trackCollection).at(1).px(),(*trackCollection).at(1).py(),(*trackCollection).at(1).pz(),
277  sqrt(((*trackCollection).at(1).p()*(*trackCollection).at(1).p())+daughterMass_*daughterMass_));
278  TLorentzVector mother = track0+track1;
279 
280  invariantMass_->Fill( mother.M() );
281  }else{
282  edm::LogInfo("Alignment")<<"wrong number of tracks trackcollection encountered: "<<(*trackCollection).size();
283  }
284  }
285 }
286 
288 {
289  for (trackingRecHit_iterator iHit = track.recHitsBegin(); iHit != track.recHitsEnd(); ++iHit) {
290  if( (*iHit)->isValid() ){
291  const TrackingRecHit *hit = (*iHit);
292  const DetId geoId(hit->geographicalId());
293  const GeomDet* gd = geometry.idToDet(geoId);
294  // since 2_1_X local hit positions are transient. taking center of the hit module for now.
295  // The alternative would be the coarse estimation or a refit.
296  //const GlobalPoint globP( gd->toGlobal( hit->localPosition() ) );
297  const GlobalPoint globP( gd->toGlobal( Local3DPoint(0.,0.,0.) ) );
298  double r = sqrt( globP.x()*globP.x() + globP.y()*globP.y() );
299  if( useSignedR_ )
300  r*= globP.y() / fabs( globP.y() );
301 
302  Hits_ZvsR_->Fill( globP.z(), r );
303  Hits_XvsY_->Fill( globP.x(), globP.y() );
304 
305  if( fillRawIdMap_)
306  Hits_perDetId_->Fill( binByRawId_[ geoId.rawId() ]);
307  }
308  }
309 }
310 
312 {
313  std::vector<int> sortedRawIds;
314  for (std::vector<DetId>::const_iterator iDetId = geometry.detUnitIds().begin();
315  iDetId != geometry.detUnitIds().end(); ++iDetId) {
316  sortedRawIds.push_back((*iDetId).rawId());
317  }
318  std::sort(sortedRawIds.begin(), sortedRawIds.end());
319 
320  int i = 0;
321  for (std::vector<int>::iterator iRawId = sortedRawIds.begin();
322  iRawId != sortedRawIds.end(); ++iRawId){
323  binByRawId_[*iRawId] = i;
324  ++i;
325  }
326 }
327 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:495
virtual 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:457
TrackQuality
track quality
Definition: TrackBase.h:149
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * TrackQuality_
#define NULL
Definition: scimark2.h:8
MonitorElement * minJetDeltaR_
std::map< int, int > binByRawId_
MonitorElement * Hits_XvsY_
MonitorElement * Hits_ZvsR_
void Fill(long long x)
MonitorElement * AlCaRecoTrackEfficiency_
MonitorElement * sumCharge_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::TrackCollection > trackProducer_
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
MonitorElement * invariantMass_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * TrackPtNegative_
TH1 * getTH1(void) const
void fillHitmaps(const reco::Track &track, const TrackerGeometry &geometry)
edm::EDGetTokenT< reco::CaloJetCollection > jetCollection_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
bool isValid() const
Definition: HandleBase.h:75
char const * module
Definition: ProductLabels.h:5
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
Definition: DetId.h:18
MonitorElement * Hits_perDetId_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
const T & get() const
Definition: EventSetup.h:55
MonitorElement * TrackCurvature_
MonitorElement * jetPt_
edm::ParameterSet conf_
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
ESHandle< TrackerGeometry > geometry
MonitorElement * minTrackDeltaR_
virtual const DetIdContainer & detUnitIds() const
Returm a vector of all GeomDetUnit DetIds.
DetId geographicalId() const
bool isValid() const
Definition: ESHandle.h:47
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)
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
Definition: Run.h:41
TkAlCaRecoMonitor(const edm::ParameterSet &)
virtual const TrackerGeomDet * idToDet(DetId) const
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109