CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentStats.cc
Go to the documentation of this file.
2 
10 
16 
19 
25 
26 using namespace std;
27 
29  src_(iConfig.getParameter<edm::InputTag>("src")),
30  overlapAM_(iConfig.getParameter<edm::InputTag>("OverlapAssoMap")),
31  keepTrackStats_(iConfig.getParameter<bool>("keepTrackStats")),
32  keepHitPopulation_(iConfig.getParameter<bool>("keepHitStats")),
33  statsTreeName_(iConfig.getParameter<string>("TrkStatsFileName")),
34  hitsTreeName_(iConfig.getParameter<string>("HitStatsFileName")),
35  prescale_(iConfig.getParameter<uint32_t>("TrkStatsPrescale")),
36  lastSetup_(nullptr)
37 {
38 
39  //sanity checks
40 
41  //init
42  outtree_=0;
43 
44 }//end constructor
45 
47  //
48 }
49 
50 void AlignmentStats::beginJob(){// const edm::EventSetup &iSetup
51 
52  //book track stats tree
53  treefile_=new TFile(statsTreeName_.c_str(),"RECREATE");
54  treefile_->cd();
55  outtree_=new TTree("AlignmentTrackStats","Statistics of Tracks used for Alignment");
56  // int nHitsinPXB[MAXTRKS_], nHitsinPXE[MAXTRKS_], nHitsinTEC[MAXTRKS_], nHitsinTIB[MAXTRKS_],nHitsinTOB[MAXTRKS_],nHitsinTID[MAXTRKS_];
57 
58 
59  outtree_->Branch("Ntracks" ,&ntracks ,"Ntracks/i");
60  outtree_->Branch("Run_" ,&run_ ,"Run_Nr/I");
61  outtree_->Branch("Event" ,&event_ ,"EventNr/I");
62  outtree_->Branch("Eta" ,&Eta ,"Eta[Ntracks]/F");
63  outtree_->Branch("Phi" ,&Phi ,"Phi[Ntracks]/F");
64  outtree_->Branch("P" ,&P ,"P[Ntracks]/F");
65  outtree_->Branch("Pt" ,&Pt ,"Pt[Ntracks]/F");
66  outtree_->Branch("Chi2n" ,&Chi2n ,"Chi2n[Ntracks]/F");
67  outtree_->Branch("Nhits" ,&Nhits ,"Nhits[Ntracks][7]/I");
68  /*
69  outtree_->Branch("NhitsPXB" , ,);
70  outtree_->Branch("NhitsPXE" , ,);
71  outtree_->Branch("NhitsTIB" , ,);
72  outtree_->Branch("NhitsTID" , ,);
73  outtree_->Branch("NhitsTOB" , ,);
74  outtree_->Branch("NhitsTOB" , ,);
75  */
76 
78  firstEvent_=true;
79 
80  //load the tracker geometry from the EventSetup
81  // iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeometry_);
82 
83 }//end beginJob
84 
85 
87 
88  lastSetup_ = &iSetup;
89 
90  //load list of detunits needed then in endJob
91 
92  if(firstEvent_){
93  edm::ESHandle<TrackerGeometry> tmpTkGeometry;
94  iSetup.get<TrackerDigiGeometryRecord>().get(tmpTkGeometry);
95  trackerGeometry_=&(*tmpTkGeometry);
96  firstEvent_=false;
97  }
98 
99 
100  //take trajectories and tracks to loop on
101  // edm::Handle<TrajTrackAssociationCollection> TrackAssoMap;
103  iEvent.getByLabel(src_, Tracks);
104 
105  //take overlap HitAssomap
107  iEvent.getByLabel(overlapAM_, hMap);
108  const AliClusterValueMap &OverlapMap=*hMap;
109 
110  // Initialise
111  run_=1;
112  event_=1;
113  ntracks=0;
114  run_=iEvent.id().run();
115  event_=iEvent.id().event();
116  ntracks=Tracks->size();
117  // if(ntracks>1) std::cout<<"~~~~~~~~~~~~\n For this event processing "<<ntracks<<" tracks"<<std::endl;
118 
119  unsigned int trk_cnt=0;
120 
121  for(int j=0;j<MAXTRKS_;j++){
122  Eta[j]=-9999.0;
123  Phi[j]=-8888.0;
124  P[j] =-7777.0;
125  Pt[j] =-6666.0;
126  Chi2n[j]=-2222.0;
127  for(int k=0;k<7;k++){
128  Nhits[j][k]=0;
129  }
130  }
131 
132  // int npxbhits=0;
133 
134  //loop on tracks
135  for(std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk; ++ittrk){
136  Eta[trk_cnt]=ittrk->eta();
137  Phi[trk_cnt]=ittrk->phi();
138  Chi2n[trk_cnt]=ittrk->normalizedChi2();
139  P[trk_cnt]=ittrk->p();
140  Pt[trk_cnt]=ittrk->pt();
141  Nhits[trk_cnt][0]=ittrk->numberOfValidHits();
142 
143  // if(ntracks>1)std::cout<<"Track #"<<trk_cnt+1<<" params: Eta="<< Eta[trk_cnt]<<" Phi="<< Phi[trk_cnt]<<" P="<<P[trk_cnt]<<" Nhits="<<Nhits[trk_cnt][0]<<std::endl;
144 
145  int nhit=0;
146  //loop on tracking rechits
147  //std::cout << " loop on hits of track #" << (itt - tracks->begin()) << std::endl;
148  for (trackingRecHit_iterator ith = ittrk->recHitsBegin(), edh = ittrk->recHitsEnd(); ith != edh; ++ith) {
149 
150  const TrackingRecHit *hit = ith->get(); // ith is an iterator on edm::Ref to rechit
151  if(! hit->isValid())continue;
152  DetId detid = hit->geographicalId();
153  int subDet = detid.subdetId();
154  uint32_t rawId = hit->geographicalId().rawId();
155 
156  // if(subDet==1)npxbhits++;
157 
158  //look if you find this detid in the map
159  DetHitMap::iterator mapiter;
160  mapiter= hitmap_.find(rawId);
161  if(mapiter!=hitmap_.end() ){//present, increase its value by one
162  // hitmap_[rawId]=hitmap_[rawId]+1;
163  ++(hitmap_[rawId]);
164  }
165  else{//not present, let's add this key to the map with value=1
166  hitmap_.insert(pair<uint32_t, uint32_t>(rawId, 1));
167  }
168 
169  AlignmentClusterFlag inval;
170 
171  bool hitInPixel= (subDet==PixelSubdetector::PixelBarrel)||(subDet==PixelSubdetector::PixelEndcap);
172  bool hitInStrip=(subDet==SiStripDetId::TIB) || (subDet==SiStripDetId::TID) ||(subDet==SiStripDetId::TOB) ||(subDet==SiStripDetId::TEC);
173 
174  if(!(hitInPixel||hitInStrip)){
175  //skip only this hit, don't stop everything throwing an exception
176  edm::LogError("AlignmentStats")<<"Hit not belonging neither to pixels nor to strips! Skipping it. SubDet= "<<subDet;
177  continue;
178  }
179 
180  //check also if the hit is an overlap. If yes fill a dedicated hitmap
181  if (hitInStrip){
182  const std::type_info &type = typeid(*hit);
183 
184  if (type == typeid(SiStripRecHit1D)) {
185  //Notice the difference respect to when one loops on Trajectories: the recHit is a TrackingRecHit and not a TransientTrackingRecHit
186  const SiStripRecHit1D* striphit=dynamic_cast<const SiStripRecHit1D*>(hit);
187  if(striphit!=0){
188  SiStripRecHit1D::ClusterRef stripclust(striphit->cluster());
189  inval = OverlapMap[stripclust];
190  }
191  else{
192  // edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit1D failed! TypeId of the RecHit: "<<className(*hit);
193  throw cms::Exception("NullPointerError")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit1D failed! TypeId of the RecHit: "<<className(*hit);
194  }
195  }//end if sistriprechit1D
196  if (type == typeid(SiStripRecHit2D)) {
197  const SiStripRecHit2D* striphit=dynamic_cast<const SiStripRecHit2D*>(hit);
198  if(striphit!=0){
199  SiStripRecHit2D::ClusterRef stripclust(striphit->cluster());
200  inval = OverlapMap[stripclust];
201  //cout<<"Taken the Strip Cluster with ProdId "<<stripclust.id() <<"; the Value in the map is "<<inval<<" (DetId is "<<hit->geographicalId().rawId()<<")"<<endl;
202  }
203  else{
204  throw cms::Exception("NullPointerError")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit2D failed! TypeId of the RecHit: "<<className(*hit);
205  // edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit2D failed! TypeId of the RecHit: "<<className(*hit);
206  }
207  }//end if sistriprechit2D
208 
209  }//end if hit in Strips
210  else{
211  const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
212  if(pixelhit!=0){
213  SiPixelClusterRefNew pixclust(pixelhit->cluster());
214  inval = OverlapMap[pixclust];
215  //cout<<"Taken the Pixel Cluster with ProdId "<<pixclust.id() <<"; the Value in the map is "<<inval<<" (DetId is "<<hit->geographicalId().rawId()<<")"<<endl;
216  }
217  else{
218  edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Pixel RecHit failed! TypeId of the RecHit: "<<className(*hit);
219  }
220  }//end else hit is in Pixel
221 
222 
223  bool isOverlapHit(inval.isOverlap());
224 
225  if( isOverlapHit ){
226  //cout<<"This hit is an overlap !"<<endl;
227  DetHitMap::iterator overlapiter;
228  overlapiter=overlapmap_.find(rawId);
229 
230  if(overlapiter!=overlapmap_.end() ){//the det already collected at least an overlap, increase its value by one
231  overlapmap_[rawId]=overlapmap_[rawId]+1;
232  }
233  else{//first overlap on det unit, let's add it to the map
234  overlapmap_.insert(pair<uint32_t, uint32_t>(rawId, 1));
235  }
236  }//end if the hit is an overlap
237 
238 
239 
240 
241  int subdethit= static_cast<int>(hit->geographicalId().subdetId());
242  // if(ntracks>1)std::cout<<"Hit in SubDet="<<subdethit<<std::endl;
243  Nhits[trk_cnt][subdethit]=Nhits[trk_cnt][subdethit]+1;
244  nhit++;
245  }//end loop on trackingrechits
246  trk_cnt++;
247 
248  }//end loop on tracks
249 
250  // cout<<"Total number of pixel hits is "<<npxbhits<<endl;
251 
252  tmpPresc_--;
253  if(tmpPresc_<1){
254  outtree_->Fill();
256  }
257  if(trk_cnt!=ntracks)edm::LogError("AlignmentStats")<<"\nERROR! trk_cnt="<<trk_cnt<<" ntracks="<<ntracks;
258  trk_cnt=0;
259 
260  return;
261 }
262 
264 
265  treefile_->cd();
266  edm::LogInfo("AlignmentStats")<<"Writing out the TrackStatistics in "<<gDirectory->GetPath()<<std::endl;
267  outtree_->Write();
268  delete outtree_;
269 
270  //create tree with hit maps (hitstree)
271  //book track stats tree
272  TFile *hitsfile=new TFile(hitsTreeName_.c_str(),"RECREATE");
273  hitsfile->cd();
274  TTree *hitstree=new TTree("AlignmentHitMap","Maps of Hits used for Alignment");
275 
276  unsigned int id=0,nhits=0,noverlaps=0;
277  float posX(-99999.0),posY(-77777.0),posZ(-88888.0);
278  float posEta(-6666.0),posPhi(-5555.0),posR(-4444.0);
279  int subdet=0;
280  unsigned int layer=0;
281  bool is2D=false,isStereo=false;
282  hitstree->Branch("DetId", &id , "DetId/i");
283  hitstree->Branch("Nhits", &nhits , "Nhits/i");
284  hitstree->Branch("Noverlaps",&noverlaps,"Noverlaps/i");
285  hitstree->Branch("SubDet", &subdet, "SubDet/I");
286  hitstree->Branch("Layer", &layer, "Layer/i");
287  hitstree->Branch("is2D" , &is2D, "is2D/B");
288  hitstree->Branch("isStereo", &isStereo, "isStereo/B");
289  hitstree->Branch("posX", &posX, "posX/F");
290  hitstree->Branch("posY", &posY, "posY/F");
291  hitstree->Branch("posZ", &posZ, "posZ/F");
292  hitstree->Branch("posR", &posR, "posR/F");
293  hitstree->Branch("posEta", &posEta, "posEta/F");
294  hitstree->Branch("posPhi", &posPhi, "posPhi/F");
295 
296  /*
297  TTree *overlapstree=new TTree("OverlapHitMap","Maps of Overlaps used for Alignment");
298  hitstree->Branch("DetId", &id , "DetId/i");
299  hitstree->Branch("NOverlaps", &nhits , "Nhits/i");
300  hitstree->Branch("SubDet", &subdet, "SubDet/I");
301  hitstree->Branch("Layer", &layer, "Layer/i");
302  hitstree->Branch("is2D" , &is2D, "is2D/B");
303  hitstree->Branch("isStereo",&isStereo,"isStereo/B");
304  hitstree->Branch("posX", &posX, "posX/F");
305  hitstree->Branch("posY", &posY, "posY/F");
306  hitstree->Branch("posZ", &posZ, "posZ/F");
307  hitstree->Branch("posR", &posR, "posR/F");
308  hitstree->Branch("posEta", &posEta, "posEta/F");
309  hitstree->Branch("posPhi", &posPhi, "posPhi/F");
310  */
311 
312  //Retrieve tracker topology from geometry
313  edm::ESHandle<TrackerTopology> tTopoHandle;
314  lastSetup_->get<IdealGeometryRecord>().get(tTopoHandle);
315  const TrackerTopology* const tTopo = tTopoHandle.product();
316 
317  AlignableTracker* theAliTracker=new AlignableTracker(&(*trackerGeometry_), tTopo);
318  const std::vector<Alignable*>& Detunitslist=theAliTracker->deepComponents();
319  int ndetunits=Detunitslist.size();
320  edm::LogInfo("AlignmentStats")<<"Number of DetUnits in the AlignableTracker: "<< ndetunits<<std::endl;
321 
322  for(int det_cnt=0;det_cnt< ndetunits;++det_cnt){
323 
324  //re-initialize for safety
325  id=0;
326  nhits=0;
327  noverlaps=0;
328  posX=-99999.0;
329  posY=-77777.0;
330  posZ=-88888.0;
331  posEta=-6666.0;
332  posPhi=-5555.0;
333  posR=-4444.0;
334  subdet=0;
335  layer=0;
336  is2D=false;
337  isStereo=false;
338 
339  //if detunit in vector is found also in the map, look for how many hits were collected
340  //and save in the tree this number
341  id=static_cast <uint32_t>( Detunitslist[det_cnt]->id() );
342  if( hitmap_.find(id) != hitmap_.end() ){
343  nhits=hitmap_[id];
344  }
345  //if not, save nhits=0
346  else{
347  nhits=0;
348  hitmap_.insert(pair<uint32_t, uint32_t>(id, 0));
349  }
350 
351  if( overlapmap_.find(id) != overlapmap_.end() ){
352  noverlaps=overlapmap_[id];
353  }
354  //if not, save nhits=0
355  else{
356  noverlaps=0;
357  overlapmap_.insert(pair<uint32_t, uint32_t>(id, 0));
358  }
359 
360  //take other geometrical infos from the det
361  posX= Detunitslist[det_cnt]->globalPosition().x();
362  posY= Detunitslist[det_cnt]->globalPosition().y();
363  posZ= Detunitslist[det_cnt]->globalPosition().z();
364 
365  align::GlobalVector vec(posX,posY,posZ);
366  posR = vec.perp();
367  posPhi = vec.phi();
368  posEta = vec.eta();
369  // posPhi = atan2(posY,posX);
370 
371  DetId detid(id);
372  subdet=detid.subdetId();
373 
374  //get layers, petals, etc...
375  if(subdet==PixelSubdetector::PixelBarrel){//PXB
376 
377  layer=tTopo->pxbLayer(id);
378  is2D=true;
379  isStereo=false;
380  }
381  else if(subdet==PixelSubdetector::PixelEndcap){
382 
383  layer=tTopo->pxfDisk(id);
384  is2D=true;
385  isStereo=false;
386  }
387  else if(subdet==SiStripDetId::TIB){
388 
389  layer=tTopo->tibLayer(id);
390  is2D=tTopo->tibIsDoubleSide(id);
391  isStereo=tTopo->tibIsStereo(id);
392  }
393  else if(subdet==SiStripDetId::TID){
394 
395  layer=tTopo->tidWheel(id);
396  is2D=tTopo->tidIsDoubleSide(id);
397  isStereo=tTopo->tidIsStereo(id);
398  }
399  else if(subdet==SiStripDetId::TOB){
400 
401  layer=tTopo->tobLayer(id);
402  is2D=tTopo->tobIsDoubleSide(id);
403  isStereo=tTopo->tobIsStereo(id);
404  }
405  else if(subdet==SiStripDetId::TEC){
406 
407  layer=tTopo->tecWheel(id);
408  is2D=tTopo->tecIsDoubleSide(id);
409  isStereo=tTopo->tecIsStereo(id);
410  }
411  else{
412  edm::LogError("AlignmentStats")<<"Detector not belonging neither to pixels nor to strips! Skipping it. SubDet= "<<subdet;
413  }
414 
415  //write in the hitstree
416  hitstree->Fill();
417  }//end loop over detunits
418 
419 
420  //save hitstree
421  hitstree->Write();
422  delete hitstree;
423  //delete Detunitslist;
424  hitmap_.clear();
425  overlapmap_.clear();
426  delete hitsfile;
427 }
428 // ========= MODULE DEF ==============
RunNumber_t run() const
Definition: EventID.h:42
float Eta[MAXTRKS_]
type
Definition: HCALResponse.h:21
EventNumber_t event() const
Definition: EventID.h:44
uint32_t prescale_
DetHitMap hitmap_
bool tecIsDoubleSide(const DetId &id) const
bool tobIsDoubleSide(const DetId &id) const
bool tibIsDoubleSide(const DetId &id) const
T perp() const
Definition: PV3DBase.h:72
unsigned int tibLayer(const DetId &id) const
const edm::EventSetup * lastSetup_
bool tobIsStereo(const DetId &id) const
edm::InputTag src_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
unsigned int pxfDisk(const DetId &id) const
DetHitMap overlapmap_
#define nullptr
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
float Chi2n[MAXTRKS_]
unsigned int tidWheel(const DetId &id) const
bool tidIsStereo(const DetId &id) const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
float P[MAXTRKS_]
bool tecIsStereo(const DetId &id) const
int iEvent
Definition: GenABIO.cc:230
uint32_t tmpPresc_
float Phi[MAXTRKS_]
int j
Definition: DBlmapReader.cc:9
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
int k[5][pyjets_maxn]
unsigned int pxbLayer(const DetId &id) const
Definition: DetId.h:18
const TrackerGeometry * trackerGeometry_
const Alignables & deepComponents() const
Definition: Alignable.h:71
std::string statsTreeName_
virtual void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup)
const T & get() const
Definition: EventSetup.h:55
bool tibIsStereo(const DetId &id) const
AlignmentStats(const edm::ParameterSet &iConfig)
bool tidIsDoubleSide(const DetId &id) const
T const * product() const
Definition: ESHandle.h:62
bool isValid() const
static const int MAXTRKS_
T eta() const
Definition: PV3DBase.h:76
unsigned int ntracks
edm::EventID id() const
Definition: EventBase.h:56
int Nhits[MAXTRKS_][7]
std::string hitsTreeName_
float Pt[MAXTRKS_]
DetId geographicalId() const
unsigned int tecWheel(const DetId &id) const
edm::InputTag overlapAM_
std::string className(const T &t)
Definition: ClassName.h:30
unsigned int tobLayer(const DetId &id) const
Pixel Reconstructed Hit.