CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Alignment/CommonAlignmentMonitor/plugins/AlignmentStats.cc

Go to the documentation of this file.
00001 #include "Alignment/CommonAlignmentMonitor/plugins/AlignmentStats.h"
00002 
00003 #include "DataFormats/Common/interface/View.h"
00004 #include "DataFormats/DetId/interface/DetId.h"
00005 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00006 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00007 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00008 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00009 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00010 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00011 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00012 #include "Alignment/CommonAlignment/interface/Alignable.h"
00013 #include "Alignment/CommonAlignment/interface/Utilities.h"
00014 
00015 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00016 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00017 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00018 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
00019 #include "RecoTracker/TransientTrackingRecHit/interface/TSiPixelRecHit.h"
00020 
00021 #include "DataFormats/Alignment/interface/AlignmentClusterFlag.h"
00022 #include "DataFormats/Alignment/interface/AliClusterValueMap.h"
00023 
00024 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00025 #include "DataFormats/TrackReco/interface/Track.h"
00026 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00027 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00028 #include "Utilities/General/interface/ClassName.h"
00029 
00030 //using namespace edm;
00031 
00032 AlignmentStats::AlignmentStats(const edm::ParameterSet &iConfig) :
00033   src_(iConfig.getParameter<edm::InputTag>("src")),
00034   overlapAM_(iConfig.getParameter<edm::InputTag>("OverlapAssoMap")),
00035   keepTrackStats_(iConfig.getParameter<bool>("keepTrackStats")),
00036   keepHitPopulation_(iConfig.getParameter<bool>("keepHitStats")),
00037   statsTreeName_(iConfig.getParameter<string>("TrkStatsFileName")),
00038   hitsTreeName_(iConfig.getParameter<string>("HitStatsFileName")),
00039   prescale_(iConfig.getParameter<uint32_t>("TrkStatsPrescale"))
00040 {
00041 
00042   //sanity checks
00043 
00044   //init
00045   outtree_=0;
00046 
00047 }//end constructor
00048 
00049 AlignmentStats::~AlignmentStats(){
00050                                //
00051 }
00052 
00053 void AlignmentStats::beginJob(){// const edm::EventSetup &iSetup
00054 
00055   //book track stats tree
00056   treefile_=new TFile(statsTreeName_.c_str(),"RECREATE");
00057   treefile_->cd();
00058   outtree_=new TTree("AlignmentTrackStats","Statistics of Tracks used for Alignment");
00059   // int nHitsinPXB[MAXTRKS_], nHitsinPXE[MAXTRKS_], nHitsinTEC[MAXTRKS_], nHitsinTIB[MAXTRKS_],nHitsinTOB[MAXTRKS_],nHitsinTID[MAXTRKS_];
00060   
00061 
00062   outtree_->Branch("Ntracks"    ,&ntracks ,"Ntracks/i");
00063   outtree_->Branch("Run_"        ,&run_     ,"Run_Nr/I");
00064   outtree_->Branch("Event"      ,&event_   ,"EventNr/I");
00065   outtree_->Branch("Eta"        ,&Eta     ,"Eta[Ntracks]/F");
00066   outtree_->Branch("Phi"        ,&Phi     ,"Phi[Ntracks]/F");
00067   outtree_->Branch("P"          ,&P       ,"P[Ntracks]/F");
00068   outtree_->Branch("Pt"         ,&Pt      ,"Pt[Ntracks]/F");
00069   outtree_->Branch("Chi2n"      ,&Chi2n   ,"Chi2n[Ntracks]/F");
00070   outtree_->Branch("Nhits"      ,&Nhits   ,"Nhits[Ntracks][7]/I");
00071   /*
00072     outtree_->Branch("NhitsPXB"       , ,);
00073     outtree_->Branch("NhitsPXE"       , ,);
00074     outtree_->Branch("NhitsTIB"       , ,);
00075     outtree_->Branch("NhitsTID"       , ,);
00076     outtree_->Branch("NhitsTOB"       , ,);
00077     outtree_->Branch("NhitsTOB"       , ,);
00078   */
00079 
00080   tmpPresc_=prescale_;
00081   firstEvent_=true;
00082 
00083   //load the tracker geometry from the EventSetup
00084   //  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeometry_);
00085  
00086 }//end beginJob
00087 
00088 
00089 void AlignmentStats::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup){
00090 
00091   //load list of detunits needed then in endJob
00092  
00093   if(firstEvent_){
00094     edm::ESHandle<TrackerGeometry> tmpTkGeometry;
00095     iSetup.get<TrackerDigiGeometryRecord>().get(tmpTkGeometry);
00096     trackerGeometry_=&(*tmpTkGeometry);
00097     firstEvent_=false;
00098   }
00099 
00100 
00101   //take trajectories and tracks to loop on
00102   // edm::Handle<TrajTrackAssociationCollection> TrackAssoMap;
00103   edm::Handle<reco::TrackCollection> Tracks;
00104   iEvent.getByLabel(src_, Tracks);
00105 
00106   //take overlap HitAssomap
00107   edm::Handle<AliClusterValueMap> hMap;
00108   iEvent.getByLabel(overlapAM_, hMap);
00109   const AliClusterValueMap &OverlapMap=*hMap;
00110 
00111   // Initialise
00112   run_=1;
00113   event_=1;
00114   ntracks=0;
00115   run_=iEvent.id().run();
00116   event_=iEvent.id().event();
00117   ntracks=Tracks->size();
00118   //  if(ntracks>1) std::cout<<"~~~~~~~~~~~~\n For this event processing "<<ntracks<<" tracks"<<std::endl;
00119 
00120   unsigned int trk_cnt=0;
00121   
00122   for(int j=0;j<MAXTRKS_;j++){
00123     Eta[j]=-9999.0;
00124     Phi[j]=-8888.0;
00125     P[j]  =-7777.0;
00126     Pt[j] =-6666.0;
00127     Chi2n[j]=-2222.0;
00128     for(int k=0;k<7;k++){
00129       Nhits[j][k]=0;
00130     }
00131   }
00132   
00133   // int npxbhits=0;
00134 
00135   //loop on tracks
00136   for(std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk; ++ittrk){
00137     Eta[trk_cnt]=ittrk->eta();
00138     Phi[trk_cnt]=ittrk->phi();
00139     Chi2n[trk_cnt]=ittrk->normalizedChi2();
00140     P[trk_cnt]=ittrk->p();
00141     Pt[trk_cnt]=ittrk->pt();
00142     Nhits[trk_cnt][0]=ittrk->numberOfValidHits();
00143  
00144     //   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;
00145    
00146     int nhit=0;
00147     //loop on tracking rechits
00148     //std::cout << "   loop on hits of track #" << (itt - tracks->begin()) << std::endl;
00149     for (trackingRecHit_iterator ith = ittrk->recHitsBegin(), edh = ittrk->recHitsEnd(); ith != edh; ++ith) {
00150 
00151       const TrackingRecHit *hit = ith->get(); // ith is an iterator on edm::Ref to rechit
00152       if(! hit->isValid())continue;
00153       DetId detid = hit->geographicalId();
00154       int subDet = detid.subdetId();
00155       uint32_t rawId = hit->geographicalId().rawId();
00156 
00157       //  if(subDet==1)npxbhits++;
00158 
00159       //look if you find this detid in the map
00160       DetHitMap::iterator mapiter;
00161       mapiter= hitmap_.find(rawId);
00162       if(mapiter!=hitmap_.end() ){//present, increase its value by one
00163         //      hitmap_[rawId]=hitmap_[rawId]+1;
00164         ++(hitmap_[rawId]);
00165       }
00166       else{//not present, let's add this key to the map with value=1
00167         hitmap_.insert(pair<uint32_t, uint32_t>(rawId, 1));
00168       }
00169       
00170       AlignmentClusterFlag inval;
00171 
00172       bool hitInPixel= (subDet==PixelSubdetector::PixelBarrel)||(subDet==PixelSubdetector::PixelEndcap);
00173       bool hitInStrip=(subDet==SiStripDetId::TIB) || (subDet==SiStripDetId::TID) ||(subDet==SiStripDetId::TOB) ||(subDet==SiStripDetId::TEC);
00174 
00175       if(!(hitInPixel||hitInStrip)){
00176         //skip only this hit, don't stop everything throwing an exception
00177         edm::LogError("AlignmentStats")<<"Hit not belonging neither to pixels nor to strips! Skipping it. SubDet= "<<subDet;
00178         continue;
00179       }
00180 
00181       //check also if the hit is an overlap. If yes fill a dedicated hitmap
00182       if (hitInStrip){
00183         const std::type_info &type = typeid(*hit);
00184 
00185         if (type == typeid(SiStripRecHit1D)) {
00186           //Notice the difference respect to when one loops on Trajectories: the recHit is a TrackingRecHit and not a TransientTrackingRecHit
00187           const SiStripRecHit1D* striphit=dynamic_cast<const  SiStripRecHit1D*>(hit);
00188           if(striphit!=0){
00189             SiStripRecHit1D::ClusterRef stripclust(striphit->cluster());
00190             inval = OverlapMap[stripclust];
00191           }
00192           else{
00193             //  edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit1D failed!   TypeId of the RecHit: "<<className(*hit);
00194             throw cms::Exception("NullPointerError")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit1D failed!   TypeId of the RecHit: "<<className(*hit);
00195           }
00196         }//end if sistriprechit1D
00197         if (type == typeid(SiStripRecHit2D)) {
00198           const SiStripRecHit2D* striphit=dynamic_cast<const  SiStripRecHit2D*>(hit);
00199           if(striphit!=0){
00200             SiStripRecHit2D::ClusterRef stripclust(striphit->cluster());
00201             inval = OverlapMap[stripclust];
00202             //cout<<"Taken the Strip Cluster with ProdId "<<stripclust.id() <<"; the Value in the map is "<<inval<<"  (DetId is "<<hit->geographicalId().rawId()<<")"<<endl;
00203           }
00204           else{
00205             throw cms::Exception("NullPointerError")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit2D failed!   TypeId of the RecHit: "<<className(*hit);
00206             //    edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit2D failed!   TypeId of the RecHit: "<<className(*hit);
00207           }
00208         }//end if sistriprechit2D
00209 
00210       }//end if hit in Strips
00211       else{
00212         const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
00213         if(pixelhit!=0){
00214           SiPixelClusterRefNew pixclust(pixelhit->cluster());
00215           inval = OverlapMap[pixclust];
00216           //cout<<"Taken the Pixel Cluster with ProdId "<<pixclust.id() <<"; the Value in the map is "<<inval<<"  (DetId is "<<hit->geographicalId().rawId()<<")"<<endl;
00217         }
00218         else{
00219           edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Pixel RecHit failed!   TypeId of the RecHit: "<<className(*hit);
00220         }
00221       }//end else hit is in Pixel
00222 
00223 
00224       bool isOverlapHit(inval.isOverlap());
00225       
00226       if( isOverlapHit ){
00227         //cout<<"This hit is an overlap !"<<endl;
00228         DetHitMap::iterator overlapiter;
00229         overlapiter=overlapmap_.find(rawId);
00230         
00231         if(overlapiter!=overlapmap_.end() ){//the det already collected at least an overlap, increase its value by one
00232             overlapmap_[rawId]=overlapmap_[rawId]+1;
00233         }
00234         else{//first overlap on det unit, let's add it to the map
00235           overlapmap_.insert(pair<uint32_t, uint32_t>(rawId, 1));
00236         }
00237       }//end if the hit is an overlap
00238       
00239 
00240 
00241 
00242       int subdethit= static_cast<int>(hit->geographicalId().subdetId());
00243       // if(ntracks>1)std::cout<<"Hit in SubDet="<<subdethit<<std::endl;
00244       Nhits[trk_cnt][subdethit]=Nhits[trk_cnt][subdethit]+1;
00245       nhit++;
00246     }//end loop on trackingrechits
00247     trk_cnt++;
00248 
00249   }//end loop on tracks
00250 
00251   //  cout<<"Total number of pixel hits is "<<npxbhits<<endl;
00252 
00253   tmpPresc_--;
00254   if(tmpPresc_<1){
00255     outtree_->Fill();
00256     tmpPresc_=prescale_;
00257   }
00258   if(trk_cnt!=ntracks)edm::LogError("AlignmentStats")<<"\nERROR! trk_cnt="<<trk_cnt<<"   ntracks="<<ntracks;
00259   trk_cnt=0;
00260 
00261   return;
00262 }
00263 
00264 void AlignmentStats::endJob(){
00265 
00266   treefile_->cd();
00267   edm::LogInfo("AlignmentStats")<<"Writing out the TrackStatistics in "<<gDirectory->GetPath()<<std::endl;
00268   outtree_->Write();
00269   delete outtree_;
00270 
00271   //create tree with hit maps (hitstree)
00272   //book track stats tree
00273   TFile *hitsfile=new TFile(hitsTreeName_.c_str(),"RECREATE");
00274   hitsfile->cd();
00275   TTree *hitstree=new TTree("AlignmentHitMap","Maps of Hits used for Alignment");
00276 
00277   unsigned int id=0,nhits=0,noverlaps=0;
00278   float posX(-99999.0),posY(-77777.0),posZ(-88888.0);
00279   float posEta(-6666.0),posPhi(-5555.0),posR(-4444.0);
00280   int subdet=0;
00281   unsigned int layer=0; 
00282   bool is2D=false,isStereo=false;
00283   hitstree->Branch("DetId",    &id ,      "DetId/i");
00284   hitstree->Branch("Nhits",    &nhits ,   "Nhits/i");
00285   hitstree->Branch("Noverlaps",&noverlaps,"Noverlaps/i");
00286   hitstree->Branch("SubDet",   &subdet,   "SubDet/I");
00287   hitstree->Branch("Layer",    &layer,    "Layer/i");
00288   hitstree->Branch("is2D" ,    &is2D,     "is2D/B");
00289   hitstree->Branch("isStereo", &isStereo, "isStereo/B");
00290   hitstree->Branch("posX",     &posX,     "posX/F");
00291   hitstree->Branch("posY",     &posY,     "posY/F");
00292   hitstree->Branch("posZ",     &posZ,     "posZ/F");
00293   hitstree->Branch("posR",     &posR,     "posR/F");
00294   hitstree->Branch("posEta",   &posEta,   "posEta/F");
00295   hitstree->Branch("posPhi",   &posPhi,   "posPhi/F");
00296 
00297   /*
00298   TTree *overlapstree=new TTree("OverlapHitMap","Maps of Overlaps used for Alignment");
00299   hitstree->Branch("DetId",   &id ,     "DetId/i");
00300   hitstree->Branch("NOverlaps",   &nhits ,  "Nhits/i");
00301   hitstree->Branch("SubDet",  &subdet,  "SubDet/I");
00302   hitstree->Branch("Layer",   &layer,   "Layer/i");
00303   hitstree->Branch("is2D" ,   &is2D,    "is2D/B");
00304   hitstree->Branch("isStereo",&isStereo,"isStereo/B");
00305   hitstree->Branch("posX",    &posX,    "posX/F");
00306   hitstree->Branch("posY",    &posY,    "posY/F");
00307   hitstree->Branch("posZ",    &posZ,    "posZ/F");
00308   hitstree->Branch("posR",    &posR,    "posR/F");
00309   hitstree->Branch("posEta",  &posEta,  "posEta/F");
00310   hitstree->Branch("posPhi",  &posPhi,  "posPhi/F");
00311   */
00312 
00313 
00314   AlignableTracker* theAliTracker=new AlignableTracker(&(*trackerGeometry_));
00315   const std::vector<Alignable*>& Detunitslist=theAliTracker->deepComponents();
00316   int ndetunits=Detunitslist.size();
00317   edm::LogInfo("AlignmentStats")<<"Number of DetUnits in the AlignableTracker: "<< ndetunits<<std::endl;
00318 
00319   for(int det_cnt=0;det_cnt< ndetunits;++det_cnt){
00320     
00321     //re-initialize for safety
00322     id=0;
00323     nhits=0;
00324     noverlaps=0;
00325     posX=-99999.0;
00326     posY=-77777.0;
00327     posZ=-88888.0;
00328     posEta=-6666.0;
00329     posPhi=-5555.0;
00330     posR=-4444.0;
00331     subdet=0;
00332     layer=0; 
00333     is2D=false;
00334     isStereo=false;
00335 
00336     //if detunit in vector is found also in the map, look for how many hits were collected 
00337     //and save in the tree this number
00338     id=static_cast <uint32_t>( Detunitslist[det_cnt]->id() );
00339     if( hitmap_.find(id) != hitmap_.end() ){
00340       nhits=hitmap_[id];
00341     }
00342     //if not, save nhits=0
00343     else{
00344       nhits=0;
00345       hitmap_.insert(pair<uint32_t, uint32_t>(id, 0));
00346     }
00347 
00348     if( overlapmap_.find(id) != overlapmap_.end() ){
00349       noverlaps=overlapmap_[id];
00350     }
00351     //if not, save nhits=0
00352     else{
00353       noverlaps=0;
00354       overlapmap_.insert(pair<uint32_t, uint32_t>(id, 0));
00355     }
00356 
00357   //take other geometrical infos from the det
00358    posX= Detunitslist[det_cnt]->globalPosition().x();
00359    posY= Detunitslist[det_cnt]->globalPosition().y();
00360    posZ= Detunitslist[det_cnt]->globalPosition().z();
00361 
00362    align::GlobalVector vec(posX,posY,posZ);
00363    posR = vec.perp();
00364    posPhi = vec.phi();
00365    posEta = vec.eta();
00366    //   posPhi = atan2(posY,posX);      
00367 
00368    DetId detid(id);
00369    subdet=detid.subdetId();
00370 
00371    //get layers, petals, etc...
00372    if(subdet==PixelSubdetector::PixelBarrel){//PXB
00373      PXBDetId pxbdet=(id);
00374      layer=pxbdet.layer();
00375      is2D=true;
00376      isStereo=false;
00377    }
00378    else if(subdet==PixelSubdetector::PixelEndcap){
00379      PXFDetId pxfdet(id);
00380      layer=pxfdet.disk();
00381      is2D=true;
00382      isStereo=false;
00383    }
00384    else if(subdet==SiStripDetId::TIB){
00385      TIBDetId tibdet(id);
00386      layer=tibdet.layerNumber();
00387      is2D=tibdet.isDoubleSide();
00388      isStereo=tibdet.isStereo();
00389    }
00390    else if(subdet==SiStripDetId::TID){
00391      TIDDetId tiddet(id);
00392      layer=tiddet.diskNumber();
00393      is2D=tiddet.isDoubleSide();
00394      isStereo=tiddet.isStereo();
00395    }
00396    else if(subdet==SiStripDetId::TOB){
00397      TOBDetId tobdet(id);
00398      layer=tobdet.layerNumber();
00399      is2D=tobdet.isDoubleSide();
00400      isStereo=tobdet.isStereo();
00401    }
00402    else if(subdet==SiStripDetId::TEC){
00403      TECDetId tecdet(id);
00404      layer=tecdet.wheelNumber();
00405      is2D=tecdet.isDoubleSide();
00406      isStereo=tecdet.isStereo();
00407    }
00408    else{
00409      edm::LogError("AlignmentStats")<<"Detector not belonging neither to pixels nor to strips! Skipping it. SubDet= "<<subdet;
00410    }
00411 
00412   //write in the hitstree
00413     hitstree->Fill();
00414   }//end loop over detunits
00415 
00416 
00417   //save hitstree
00418   hitstree->Write();
00419   delete hitstree;
00420   //delete Detunitslist;
00421   hitmap_.clear();
00422   overlapmap_.clear();
00423   delete hitsfile;
00424 }
00425 // ========= MODULE DEF ==============
00426 #include "FWCore/PluginManager/interface/ModuleDef.h"
00427 #include "FWCore/Framework/interface/MakerMacros.h"
00428 DEFINE_FWK_MODULE(AlignmentStats);