CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/DQM/SiPixelMonitorTrack/src/SiPixelHitEfficiencySource.cc

Go to the documentation of this file.
00001 // Package:    SiPixelMonitorTrack
00002 // Class:      SiPixelHitEfficiencySource
00003 // 
00004 // class SiPixelHitEfficiencyModule SiPixelHitEfficiencyModule.cc 
00005 //       DQM/SiPixelMonitorTrack/src/SiPixelHitEfficiencyModule.cc
00006 //
00007 // Description: SiPixel hit efficiency data quality monitoring modules
00008 // Implementation: prototype -> improved -> never final - end of the 1st step 
00009 //
00010 // Original Authors: Romain Rougny & Luca Mucibello
00011 //         Created: Mar Nov 10 13:29:00 CET 2009
00012 
00013 
00014 #include <iostream>
00015 #include <map>
00016 #include <string>
00017 #include <vector>
00018 #include <utility>
00019 
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022 
00023 #include "DataFormats/Common/interface/Handle.h"
00024 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00025 
00026 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00027 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00028 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00029 
00030 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
00031 
00032 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00033 
00034 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
00035 
00036 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00037 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00038 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00039 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00040 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00041 
00042 #include "DQMServices/Core/interface/DQMStore.h"
00043 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
00044 #include "DQM/SiPixelMonitorTrack/interface/SiPixelHitEfficiencySource.h"
00045 
00046 
00047 using namespace std;
00048 using namespace edm;
00049 
00050 
00051 SiPixelHitEfficiencySource::SiPixelHitEfficiencySource(const edm::ParameterSet& pSet) :
00052   pSet_(pSet),
00053   modOn( pSet.getUntrackedParameter<bool>("modOn",true) ),
00054   ladOn( pSet.getUntrackedParameter<bool>("ladOn",false) ), 
00055   layOn( pSet.getUntrackedParameter<bool>("layOn",false) ), 
00056   phiOn( pSet.getUntrackedParameter<bool>("phiOn",false) ), 
00057   ringOn( pSet.getUntrackedParameter<bool>("ringOn",false) ), 
00058   bladeOn( pSet.getUntrackedParameter<bool>("bladeOn",false) ), 
00059   diskOn( pSet.getUntrackedParameter<bool>("diskOn",false) )
00060   //updateEfficiencies( pSet.getUntrackedParameter<bool>("updateEfficiencies",false) )
00061  { 
00062    pSet_ = pSet; 
00063    debug_ = pSet_.getUntrackedParameter<bool>("debug", false); 
00064    tracksrc_ = pSet_.getParameter<edm::InputTag>("trajectoryInput");
00065    applyEdgeCut_ = pSet_.getUntrackedParameter<bool>("applyEdgeCut");
00066    nSigma_EdgeCut_ = pSet_.getUntrackedParameter<double>("nSigma_EdgeCut");
00067    dbe_ = edm::Service<DQMStore>().operator->();
00068 
00069   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource constructor" << endl;
00070   LogInfo ("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/" 
00071             << layOn << "/" << phiOn << std::endl;
00072   LogInfo ("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/" 
00073             << ringOn << std::endl;
00074 }
00075 
00076 
00077 SiPixelHitEfficiencySource::~SiPixelHitEfficiencySource() {
00078   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource destructor" << endl;
00079 
00080   std::map<uint32_t,SiPixelHitEfficiencyModule*>::iterator struct_iter;
00081   for (struct_iter = theSiPixelStructure.begin() ; struct_iter != theSiPixelStructure.end() ; struct_iter++){
00082     delete struct_iter->second;
00083     struct_iter->second = 0;
00084   }
00085 }
00086 
00087 void SiPixelHitEfficiencySource::beginJob() {
00088   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource beginJob()" << endl;
00089   firstRun = true;
00090 }
00091 
00092 void SiPixelHitEfficiencySource::beginRun(const edm::Run& r, edm::EventSetup const& iSetup) {
00093   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource beginRun()" << endl;
00094   
00095   if(firstRun){
00096   // retrieve TrackerGeometry for pixel dets
00097   edm::ESHandle<TrackerGeometry> TG;
00098   iSetup.get<TrackerDigiGeometryRecord>().get(TG);
00099   if (debug_) LogVerbatim("PixelDQM") << "TrackerGeometry "<< &(*TG) <<" size is "<< TG->dets().size() << endl;
00100  
00101   // build theSiPixelStructure with the pixel barrel and endcap dets from TrackerGeometry
00102   for (TrackerGeometry::DetContainer::const_iterator pxb = TG->detsPXB().begin();  
00103        pxb!=TG->detsPXB().end(); pxb++) {
00104     if (dynamic_cast<PixelGeomDetUnit*>((*pxb))!=0) {
00105       SiPixelHitEfficiencyModule* module = new SiPixelHitEfficiencyModule((*pxb)->geographicalId().rawId());
00106       theSiPixelStructure.insert(pair<uint32_t, SiPixelHitEfficiencyModule*>((*pxb)->geographicalId().rawId(), module));
00107     }
00108   }
00109   for (TrackerGeometry::DetContainer::const_iterator pxf = TG->detsPXF().begin(); 
00110        pxf!=TG->detsPXF().end(); pxf++) {
00111     if (dynamic_cast<PixelGeomDetUnit*>((*pxf))!=0) {
00112       SiPixelHitEfficiencyModule* module = new SiPixelHitEfficiencyModule((*pxf)->geographicalId().rawId());
00113       theSiPixelStructure.insert(pair<uint32_t, SiPixelHitEfficiencyModule*>((*pxf)->geographicalId().rawId(), module));
00114     }
00115   }
00116   LogInfo("PixelDQM") << "SiPixelStructure size is " << theSiPixelStructure.size() << endl;
00117 
00118   // book residual histograms in theSiPixelFolder - one (x,y) pair of histograms per det
00119   SiPixelFolderOrganizer theSiPixelFolder;
00120   for (std::map<uint32_t, SiPixelHitEfficiencyModule*>::iterator pxd = theSiPixelStructure.begin(); 
00121        pxd!=theSiPixelStructure.end(); pxd++) {
00122 
00123     if(modOn){
00124       if (theSiPixelFolder.setModuleFolder((*pxd).first)) (*pxd).second->book(pSet_);
00125       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Folder Creation Failed! "; 
00126     }
00127     if(ladOn){
00128       if (theSiPixelFolder.setModuleFolder((*pxd).first,1)) (*pxd).second->book(pSet_,1);
00129       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource ladder Folder Creation Failed! "; 
00130     }
00131     if(layOn){
00132       if (theSiPixelFolder.setModuleFolder((*pxd).first,2)) (*pxd).second->book(pSet_,2);
00133       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource layer Folder Creation Failed! "; 
00134     }
00135     if(phiOn){
00136       if (theSiPixelFolder.setModuleFolder((*pxd).first,3)) (*pxd).second->book(pSet_,3);
00137       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource phi Folder Creation Failed! "; 
00138     }
00139     if(bladeOn){
00140       if (theSiPixelFolder.setModuleFolder((*pxd).first,4)) (*pxd).second->book(pSet_,4);
00141       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Blade Folder Creation Failed! "; 
00142     }
00143     if(diskOn){
00144       if (theSiPixelFolder.setModuleFolder((*pxd).first,5)) (*pxd).second->book(pSet_,5);
00145       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Disk Folder Creation Failed! "; 
00146     }
00147     if(ringOn){
00148       if (theSiPixelFolder.setModuleFolder((*pxd).first,6)) (*pxd).second->book(pSet_,6);
00149       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Ring Folder Creation Failed! "; 
00150     }
00151   }
00152   
00153   nvalid=0;
00154   nmissing=0;
00155   
00156   firstRun = false;
00157   }
00158 }
00159 
00160 
00161 void SiPixelHitEfficiencySource::endJob(void) {
00162   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource endJob()";
00163 
00164   // save the residual histograms to an output root file
00165   bool saveFile = pSet_.getUntrackedParameter<bool>("saveFile", true);
00166   if (saveFile) { 
00167     std::string outputFile = pSet_.getParameter<std::string>("outputFile");
00168     LogInfo("PixelDQM") << " - saving histograms to "<< outputFile.data();
00169     dbe_->save(outputFile);
00170   } 
00171   LogInfo("PixelDQM") << endl; // dbe_->showDirStructure();
00172   
00173   //std::cout<< "********** SUMMARY **********"<<std::endl;
00174   //std::cout<< "number of valid hits: "<<nvalid<<std::endl;
00175   //std::cout<< "number of missing hits: "<<nmissing<<std::endl;
00176 }
00177 
00178 
00179 void SiPixelHitEfficiencySource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00180 
00181   //Get the geometry
00182   ESHandle<TrackerGeometry> TG;
00183   iSetup.get<TrackerDigiGeometryRecord>().get(TG);
00184   const TrackerGeometry* theTrackerGeometry = TG.product();
00185   
00186   //get the map
00187   edm::Handle<TrajTrackAssociationCollection> match;
00188   iEvent.getByLabel(tracksrc_,match);  
00189   const TrajTrackAssociationCollection ttac = *(match.product());
00190 
00191   if(debug_){
00192     //std::cout << "Trajectories\t : " << trajColl.size() << std::endl;
00193     //std::cout << "recoTracks  \t : " << trackColl.size() << std::endl;
00194     std::cout << "+++ NEW EVENT +++"<< std::endl;
00195     std::cout << "Map entries \t : " << ttac.size() << std::endl;
00196   }
00197 
00198   std::set<SiPixelCluster> clusterSet;
00199   TrajectoryStateCombiner tsoscomb;
00200 
00201   //Loop over map entries
00202   for(TrajTrackAssociationCollection::const_iterator it =  ttac.begin();it !=  ttac.end(); ++it){
00203     const edm::Ref<std::vector<Trajectory> > traj_iterator = it->key;  
00204     // Trajectory Map, extract Trajectory for this track
00205     reco::TrackRef trackref = it->val;
00206     //tracks++;
00207     bool isBpixtrack = false, isFpixtrack = false;
00208     std::vector<TrajectoryMeasurement> tmeasColl =traj_iterator->measurements();
00209     std::vector<TrajectoryMeasurement>::const_iterator tmeasIt;
00210     //loop on measurements to find out whether there are bpix and/or fpix hits
00211     for(tmeasIt = tmeasColl.begin();tmeasIt!=tmeasColl.end();tmeasIt++){
00212       //if(! tmeasIt->updatedState().isValid()) continue; NOT NECESSARY (I HOPE)
00213       TransientTrackingRecHit::ConstRecHitPointer testhit = tmeasIt->recHit();
00214       if(testhit->geographicalId().det() != DetId::Tracker) continue; 
00215       uint testSubDetID = (testhit->geographicalId().subdetId()); 
00216       if(testSubDetID==PixelSubdetector::PixelBarrel) isBpixtrack = true;
00217       if(testSubDetID==PixelSubdetector::PixelEndcap) isFpixtrack = true;
00218     }
00219     if(isBpixtrack || isFpixtrack){
00220     
00221       if(debug_){
00222         std::cout << "isBpixtrack : " << isBpixtrack << std::endl;
00223         std::cout << "isFpixtrack : " << isFpixtrack << std::endl;
00224       }
00225       for(std::vector<TrajectoryMeasurement>::const_iterator tmeasIt = tmeasColl.begin(); tmeasIt!=tmeasColl.end(); tmeasIt++){   
00226         //if(! tmeasIt->updatedState().isValid()) continue; 
00227         
00228         TrajectoryStateOnSurface tsos = tsoscomb( tmeasIt->forwardPredictedState(), tmeasIt->backwardPredictedState() );
00229         TransientTrackingRecHit::ConstRecHitPointer hit = tmeasIt->recHit();
00230         if(hit->geographicalId().det() != DetId::Tracker )
00231           continue; 
00232         else {
00233           
00234 //        //residual
00235           const DetId & hit_detId = hit->geographicalId();
00236           //uint IntRawDetID = (hit_detId.rawId());
00237           uint IntSubDetID = (hit_detId.subdetId());
00238           
00239           if(IntSubDetID == 0 ){
00240             if(debug_) std::cout << "NO IntSubDetID\n";
00241             continue;
00242           }
00243           if(IntSubDetID!=PixelSubdetector::PixelBarrel && IntSubDetID!=PixelSubdetector::PixelEndcap)
00244             continue;
00245           
00246           
00247               //check wether hit is valid or missing using track algo flag
00248           bool isHitValid   =hit->hit()->getType()==TrackingRecHit::valid;
00249           bool isHitMissing =hit->hit()->getType()==TrackingRecHit::missing;
00250           //std::cout<<"------ New Hit"<<std::endl;
00251           //std::cout<<(hit->hit()->getType()==TrackingRecHit::missing)<<std::endl;
00252           
00253           // get the enclosed persistent hit
00254           //const TrackingRecHit *persistentHit = hit->hit();
00255           // check if it's not null, and if it's a valid pixel hit
00256           //if ((persistentHit != 0) && (typeid(*persistentHit) == typeid(SiPixelRecHit))) {
00257           
00258             if(debug_) std::cout << "the hit is persistent\n";
00259             
00260             // tell the C++ compiler that the hit is a pixel hit
00261             //const SiPixelRecHit* pixhit = dynamic_cast<const SiPixelRecHit*>( hit->hit() );
00262          
00263             //define tracker and pixel geometry and topology
00264             const TrackerGeometry& theTracker(*theTrackerGeometry);
00265             const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*> (theTracker.idToDet(hit_detId) );
00266             //test if PixelGeomDetUnit exists
00267             if(theGeomDet == 0) {
00268               if(debug_) std::cout << "NO THEGEOMDET\n";
00269               continue;
00270             }
00271                       
00272               //const RectangularPixelTopology * topol = dynamic_cast<const RectangularPixelTopology*>(&(theGeomDet->specificTopology()));
00273 
00274               std::map<uint32_t, SiPixelHitEfficiencyModule*>::iterator pxd = theSiPixelStructure.find((*hit).geographicalId().rawId());
00275 
00276               // calculate alpha and beta from cluster position
00277               LocalTrajectoryParameters ltp = tsos.localParameters();
00278               //LocalVector localDir = ltp.momentum()/ltp.momentum().mag();
00279               
00280               //float clust_alpha = atan2(localDir.z(), localDir.x());
00281               //float clust_beta = atan2(localDir.z(), localDir.y());
00282               
00283               
00284               
00285               // THE CUTS
00286               int nrows = theGeomDet->specificTopology().nrows();
00287               int ncols = theGeomDet->specificTopology().ncolumns();
00288               //
00289               //std::pair<float,float> pitchTest = theGeomDet->specificTopology().pitch();
00290               //RectangularPixelTopology rectTopolTest = RectangularPixelTopology(nrows, ncols, pitch.first, pitch.second);
00291               //std::pair<float,float> pixelTest = rectTopol.pixel(tsos.localPosition());
00292               //
00293               std::pair<float,float> pixel = theGeomDet->specificTopology().pixel(tsos.localPosition());
00294 
00295               //*************** Edge cut ********************
00296               if(applyEdgeCut_){
00297                 bool passedEdgeCut = false; 
00298                 double nsigma = nSigma_EdgeCut_;
00299               
00300                 //transforming localX,Y coordinates into a number of rows/columns
00301                 LocalPoint  actual = LocalPoint( tsos.localPosition().x(),tsos.localPosition().y() );
00302                 std::pair<float,float> pixelActual = theGeomDet->specificTopology().pixel(actual);
00303               
00304                 //Adding/substracting nsigma times the error position to the acual traj prediction position
00305                 //Depending on the signes of localX,Y => 4 "quadrants"
00306                 LocalPoint  exagerated;
00307                 LocalError tsosErr = tsos.localError().positionError();  
00308                 if ( tsos.localPosition().x()>0 && tsos.localPosition().y()>0 )
00309                   exagerated = LocalPoint(tsos.localPosition().x()+nsigma*std::sqrt(tsosErr.xx()),tsos.localPosition().y()+nsigma*std::sqrt(tsosErr.yy()) );
00310                 if ( tsos.localPosition().x()<0 && tsos.localPosition().y()>0)
00311                   exagerated = LocalPoint(tsos.localPosition().x()-nsigma*std::sqrt(tsosErr.xx()),tsos.localPosition().y()+nsigma*std::sqrt(tsosErr.yy())  );
00312                 if ( tsos.localPosition().x()<0 && tsos.localPosition().y()<0)
00313                   exagerated = LocalPoint(tsos.localPosition().x()-nsigma*std::sqrt(tsosErr.xx()),tsos.localPosition().y()-nsigma*std::sqrt(tsosErr.yy()) );
00314                 if ( tsos.localPosition().x()>0 && tsos.localPosition().y()<0)
00315                   exagerated = LocalPoint(tsos.localPosition().x()+nsigma*std::sqrt(tsosErr.xx()),tsos.localPosition().y()-nsigma*std::sqrt(tsosErr.yy()) );
00316               
00317                 //transforming localX,Y coordinates into a number of rows/columns
00318                 std::pair<float,float> pixelExagerated = theGeomDet->specificTopology().pixel(exagerated);
00319          
00320                 //if the modified traj prediction falls out of the module, we cut
00321                 if( pixelExagerated.first>nrows  || pixelExagerated.first<0)
00322                   passedEdgeCut = false;
00323                 if( pixelExagerated.second>ncols || pixelExagerated.second<0)
00324                   passedEdgeCut = false;
00325               
00326                 if(!passedEdgeCut)
00327                   continue;
00328               }
00329               //*************** Telescope cut ********************
00330               //not needed for collisions !!
00331               
00332               
00333               
00334               
00335               
00336               
00337               
00338               
00339               if(debug_){
00340                 std::cout << "Ready to add hit in histogram:\n";
00341                 //std::cout << "detid: "<<hit_detId<<std::endl;
00342                 std::cout << "isHitValid: "<<isHitValid<<std::endl;
00343                 std::cout << "isHitMissing: "<<isHitMissing<<std::endl;
00344                 //std::cout << "passedEdgeCut: "<<passedEdgeCut<<std::endl;             
00345               }
00346               
00347               
00348               if(pxd!=theSiPixelStructure.end() && isHitValid)
00349                 ++nvalid;
00350               if(pxd!=theSiPixelStructure.end() && isHitMissing)
00351                 ++nmissing;
00352                 
00353               if (pxd!=theSiPixelStructure.end() && (isHitValid || isHitMissing))
00354                 (*pxd).second->fill(ltp, isHitValid, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);      
00355 
00356           //}//end if (persistent hit exists and is pixel hit)
00357           
00358         }//end of else 
00359         
00360         
00361       }//end for (all traj measurements of pixeltrack)
00362     }//end if (is pixeltrack)
00363     else
00364       if(debug_) std::cout << "no pixeltrack:\n";
00365     
00366   }//end loop on map entries
00367 }
00368 
00369 
00370 DEFINE_FWK_MODULE(SiPixelHitEfficiencySource); // define this as a plug-in