CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/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/TrackerCommon/interface/TrackerTopology.h"
00027 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00028 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00029 #include "DataFormats/DetId/interface/DetId.h"
00030 #include "DataFormats/SiPixelDetId/interface/PixelBarrelName.h"
00031 #include "DataFormats/SiPixelDetId/interface/PixelEndcapName.h"
00032 
00033 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
00034 #include "DataFormats/VertexReco/interface/Vertex.h"
00035 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00036 
00037 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00038 
00039 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
00040 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
00041 //#include "RecoLocalTracker/SiPixelRecHits/plugins/PixelCPEGenericESProducer.h"
00042 
00043 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00044 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00045 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00046 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00047 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00048 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00049 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00050 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00051 #include "TrackingTools/TrackAssociator/interface/TrackDetectorAssociator.h"
00052 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00053 
00054 #include "DQMServices/Core/interface/DQMStore.h"
00055 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
00056 #include "DQM/SiPixelMonitorTrack/interface/SiPixelHitEfficiencySource.h"
00057 
00058 
00059 using namespace std;
00060 using namespace edm;
00061 
00062 
00063 SiPixelHitEfficiencySource::SiPixelHitEfficiencySource(const edm::ParameterSet& pSet) :
00064   pSet_(pSet),
00065   modOn( pSet.getUntrackedParameter<bool>("modOn",true) ),
00066   ladOn( pSet.getUntrackedParameter<bool>("ladOn",false) ), 
00067   layOn( pSet.getUntrackedParameter<bool>("layOn",false) ), 
00068   phiOn( pSet.getUntrackedParameter<bool>("phiOn",false) ), 
00069   ringOn( pSet.getUntrackedParameter<bool>("ringOn",false) ), 
00070   bladeOn( pSet.getUntrackedParameter<bool>("bladeOn",false) ), 
00071   diskOn( pSet.getUntrackedParameter<bool>("diskOn",false) )
00072   //updateEfficiencies( pSet.getUntrackedParameter<bool>("updateEfficiencies",false) )
00073  { 
00074    pSet_ = pSet; 
00075    debug_ = pSet_.getUntrackedParameter<bool>("debug", false); 
00076    tracksrc_ = pSet_.getParameter<edm::InputTag>("trajectoryInput");
00077    applyEdgeCut_ = pSet_.getUntrackedParameter<bool>("applyEdgeCut");
00078    nSigma_EdgeCut_ = pSet_.getUntrackedParameter<double>("nSigma_EdgeCut");
00079    dbe_ = edm::Service<DQMStore>().operator->();
00080 
00081   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource constructor" << endl;
00082   LogInfo ("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/" 
00083             << layOn << "/" << phiOn << std::endl;
00084   LogInfo ("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/" 
00085             << ringOn << std::endl;
00086 }
00087 
00088 
00089 SiPixelHitEfficiencySource::~SiPixelHitEfficiencySource() {
00090   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource destructor" << endl;
00091 
00092   std::map<uint32_t,SiPixelHitEfficiencyModule*>::iterator struct_iter;
00093   for (struct_iter = theSiPixelStructure.begin() ; struct_iter != theSiPixelStructure.end() ; struct_iter++){
00094     delete struct_iter->second;
00095     struct_iter->second = 0;
00096   }
00097 }
00098 
00099 void SiPixelHitEfficiencySource::beginJob() {
00100   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource beginJob()" << endl;
00101   firstRun = true;
00102 }
00103 
00104 void SiPixelHitEfficiencySource::beginRun(const edm::Run& r, edm::EventSetup const& iSetup) {
00105   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource beginRun()" << endl;
00106   
00107   if(firstRun){
00108   // retrieve TrackerGeometry for pixel dets
00109   edm::ESHandle<TrackerGeometry> TG;
00110   iSetup.get<TrackerDigiGeometryRecord>().get(TG);
00111   if (debug_) LogVerbatim("PixelDQM") << "TrackerGeometry "<< &(*TG) <<" size is "<< TG->dets().size() << endl;
00112  
00113   // build theSiPixelStructure with the pixel barrel and endcap dets from TrackerGeometry
00114   for (TrackerGeometry::DetContainer::const_iterator pxb = TG->detsPXB().begin();  
00115        pxb!=TG->detsPXB().end(); pxb++) {
00116     if (dynamic_cast<PixelGeomDetUnit*>((*pxb))!=0) {
00117       SiPixelHitEfficiencyModule* module = new SiPixelHitEfficiencyModule((*pxb)->geographicalId().rawId());
00118       theSiPixelStructure.insert(pair<uint32_t, SiPixelHitEfficiencyModule*>((*pxb)->geographicalId().rawId(), module));
00119     }
00120   }
00121   for (TrackerGeometry::DetContainer::const_iterator pxf = TG->detsPXF().begin(); 
00122        pxf!=TG->detsPXF().end(); pxf++) {
00123     if (dynamic_cast<PixelGeomDetUnit*>((*pxf))!=0) {
00124       SiPixelHitEfficiencyModule* module = new SiPixelHitEfficiencyModule((*pxf)->geographicalId().rawId());
00125       theSiPixelStructure.insert(pair<uint32_t, SiPixelHitEfficiencyModule*>((*pxf)->geographicalId().rawId(), module));
00126     }
00127   }
00128   LogInfo("PixelDQM") << "SiPixelStructure size is " << theSiPixelStructure.size() << endl;
00129 
00130   // book residual histograms in theSiPixelFolder - one (x,y) pair of histograms per det
00131   SiPixelFolderOrganizer theSiPixelFolder;
00132   for (std::map<uint32_t, SiPixelHitEfficiencyModule*>::iterator pxd = theSiPixelStructure.begin(); 
00133        pxd!=theSiPixelStructure.end(); pxd++) {
00134 
00135     if(modOn){
00136       if (theSiPixelFolder.setModuleFolder((*pxd).first)) (*pxd).second->book(pSet_);
00137       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Folder Creation Failed! "; 
00138     }
00139     if(ladOn){
00140       if (theSiPixelFolder.setModuleFolder((*pxd).first,1)) (*pxd).second->book(pSet_,1);
00141       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource ladder Folder Creation Failed! "; 
00142     }
00143     if(layOn){
00144       if (theSiPixelFolder.setModuleFolder((*pxd).first,2)) (*pxd).second->book(pSet_,2);
00145       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource layer Folder Creation Failed! "; 
00146     }
00147     if(phiOn){
00148       if (theSiPixelFolder.setModuleFolder((*pxd).first,3)) (*pxd).second->book(pSet_,3);
00149       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource phi Folder Creation Failed! "; 
00150     }
00151     if(bladeOn){
00152       if (theSiPixelFolder.setModuleFolder((*pxd).first,4)) (*pxd).second->book(pSet_,4);
00153       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Blade Folder Creation Failed! "; 
00154     }
00155     if(diskOn){
00156       if (theSiPixelFolder.setModuleFolder((*pxd).first,5)) (*pxd).second->book(pSet_,5);
00157       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Disk Folder Creation Failed! "; 
00158     }
00159     if(ringOn){
00160       if (theSiPixelFolder.setModuleFolder((*pxd).first,6)) (*pxd).second->book(pSet_,6);
00161       else throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Ring Folder Creation Failed! "; 
00162     }
00163   }
00164   
00165   nvalid=0;
00166   nmissing=0;
00167   
00168   firstRun = false;
00169   }
00170 }
00171 
00172 
00173 void SiPixelHitEfficiencySource::endJob(void) {
00174   LogInfo("PixelDQM") << "SiPixelHitEfficiencySource endJob()";
00175 
00176   // save the residual histograms to an output root file
00177   bool saveFile = pSet_.getUntrackedParameter<bool>("saveFile", true);
00178   if (saveFile) { 
00179     std::string outputFile = pSet_.getParameter<std::string>("outputFile");
00180     LogInfo("PixelDQM") << " - saving histograms to "<< outputFile.data();
00181     dbe_->save(outputFile);
00182   } 
00183   LogInfo("PixelDQM") << endl; // dbe_->showDirStructure();
00184   
00185   //std::cout<< "********** SUMMARY **********"<<std::endl;
00186   //std::cout<< "number of valid hits: "<<nvalid<<std::endl;
00187   //std::cout<< "number of missing hits: "<<nmissing<<std::endl;
00188 }
00189 
00190 
00191 void SiPixelHitEfficiencySource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00192 
00193   edm::Handle<reco::VertexCollection> vertexCollectionHandle;
00194   iEvent.getByLabel("offlinePrimaryVertices", vertexCollectionHandle);
00195   if(!vertexCollectionHandle.isValid()) return;
00196   nvtx_=0;
00197   vtxntrk_=-9999;
00198   vtxD0_=-9999.; vtxX_=-9999.; vtxY_=-9999.; vtxZ_=-9999.; vtxndof_=-9999.; vtxchi2_=-9999.;
00199   const reco::VertexCollection & vertices = *vertexCollectionHandle.product();
00200   reco::VertexCollection::const_iterator bestVtx=vertices.end();
00201   for(reco::VertexCollection::const_iterator it=vertices.begin(); it!=vertices.end(); ++it){
00202     if(!it->isValid()) continue;
00203     if(vtxntrk_==-9999 ||
00204        vtxntrk_<int(it->tracksSize()) ||
00205        (vtxntrk_==int(it->tracksSize()) && fabs(vtxZ_)>fabs(it->z()))){
00206       vtxntrk_=it->tracksSize();
00207       vtxD0_=it->position().rho();
00208       vtxX_=it->x();
00209       vtxY_=it->y();
00210       vtxZ_=it->z();
00211       vtxndof_=it->ndof();
00212       vtxchi2_=it->chi2();
00213       bestVtx=it;
00214     }
00215     if(fabs(it->z())<=20. && fabs(it->position().rho())<=2. && it->ndof()>4) nvtx_++;
00216   }
00217   if(nvtx_<1) return;
00218     
00219   //Get the geometry
00220   ESHandle<TrackerGeometry> TG;
00221   iSetup.get<TrackerDigiGeometryRecord>().get(TG);
00222   const TrackerGeometry* theTrackerGeometry = TG.product();
00223   
00224   //get the map
00225   edm::Handle<TrajTrackAssociationCollection> match;
00226   iEvent.getByLabel(tracksrc_,match);  
00227   const TrajTrackAssociationCollection ttac = *(match.product());
00228 
00229   if(debug_){
00230     //std::cout << "Trajectories\t : " << trajColl.size() << std::endl;
00231     //std::cout << "recoTracks  \t : " << trackColl.size() << std::endl;
00232     std::cout << "+++ NEW EVENT +++"<< std::endl;
00233     std::cout << "Map entries \t : " << ttac.size() << std::endl;
00234   }
00235 
00236   std::set<SiPixelCluster> clusterSet;
00237   TrajectoryStateCombiner tsoscomb;
00238 
00239   //Loop over map entries
00240   for(TrajTrackAssociationCollection::const_iterator it =  ttac.begin();it !=  ttac.end(); ++it){
00241     const edm::Ref<std::vector<Trajectory> > traj_iterator = it->key;  
00242     // Trajectory Map, extract Trajectory for this track
00243     reco::TrackRef trackref = it->val;
00244     //tracks++;
00245     bool isBpixtrack = false, isFpixtrack = false;
00246     int nStripHits=0; int L1hits=0; int L2hits=0; int L3hits=0; int D1hits=0; int D2hits=0;
00247     std::vector<TrajectoryMeasurement> tmeasColl =traj_iterator->measurements();
00248     std::vector<TrajectoryMeasurement>::const_iterator tmeasIt;
00249     //loop on measurements to find out what kind of hits there are
00250     for(tmeasIt = tmeasColl.begin();tmeasIt!=tmeasColl.end();tmeasIt++){
00251       //if(! tmeasIt->updatedState().isValid()) continue; NOT NECESSARY (I HOPE)
00252       TransientTrackingRecHit::ConstRecHitPointer testhit = tmeasIt->recHit();
00253       if(testhit->geographicalId().det() != DetId::Tracker) continue; 
00254       uint testSubDetID = (testhit->geographicalId().subdetId()); 
00255       const DetId & hit_detId = testhit->geographicalId();
00256       if(testSubDetID==PixelSubdetector::PixelBarrel){
00257         isBpixtrack = true;
00258         int layer = PixelBarrelName(hit_detId).layerName();
00259         if(layer==1) L1hits++;
00260         if(layer==2) L2hits++;
00261         if(layer==3) L3hits++;
00262       }
00263       if(testSubDetID==PixelSubdetector::PixelEndcap){
00264         isFpixtrack = true;
00265         int disk = PixelEndcapName(hit_detId).diskName();
00266         if(disk==1) D1hits++;
00267         if(disk==2) D2hits++;
00268       }
00269       if(testSubDetID==StripSubdetector::TIB) nStripHits++;
00270       if(testSubDetID==StripSubdetector::TOB) nStripHits++;
00271       if(testSubDetID==StripSubdetector::TID) nStripHits++;
00272       if(testSubDetID==StripSubdetector::TEC) nStripHits++;
00273     }
00274     if(isBpixtrack || isFpixtrack){
00275       if(trackref->pt()<0.6 ||
00276          nStripHits<11 ||
00277          fabs(trackref->dxy(bestVtx->position()))>0.05 ||
00278          fabs(trackref->dz(bestVtx->position()))>0.5) continue;
00279     
00280       if(debug_){
00281         std::cout << "isBpixtrack : " << isBpixtrack << std::endl;
00282         std::cout << "isFpixtrack : " << isFpixtrack << std::endl;
00283       }
00284       //std::cout<<"This tracks has so many hits: "<<tmeasColl.size()<<std::endl;
00285       for(std::vector<TrajectoryMeasurement>::const_iterator tmeasIt = tmeasColl.begin(); tmeasIt!=tmeasColl.end(); tmeasIt++){   
00286         //if(! tmeasIt->updatedState().isValid()) continue; 
00287         
00288         TrajectoryStateOnSurface tsos = tsoscomb( tmeasIt->forwardPredictedState(), tmeasIt->backwardPredictedState() );
00289         TransientTrackingRecHit::ConstRecHitPointer hit = tmeasIt->recHit();
00290         if(hit->geographicalId().det() != DetId::Tracker )
00291           continue; 
00292         else {
00293           
00294 //        //residual
00295           const DetId & hit_detId = hit->geographicalId();
00296           //uint IntRawDetID = (hit_detId.rawId());
00297           uint IntSubDetID = (hit_detId.subdetId());
00298           
00299           if(IntSubDetID == 0 ){
00300             if(debug_) std::cout << "NO IntSubDetID\n";
00301             continue;
00302           }
00303           if(IntSubDetID!=PixelSubdetector::PixelBarrel && IntSubDetID!=PixelSubdetector::PixelEndcap)
00304             continue;
00305           
00306           int disk=0; int layer=0; int panel=0; int module=0; bool isHalfModule=false;
00307           if(IntSubDetID==PixelSubdetector::PixelBarrel){ // it's a BPIX hit
00308             layer = PixelBarrelName(hit_detId).layerName();
00309             isHalfModule = PixelBarrelName(hit_detId).isHalfModule();
00310           }else if(IntSubDetID==PixelSubdetector::PixelEndcap){ // it's an FPIX hit
00311             disk = PixelEndcapName(hit_detId).diskName();
00312             panel = PixelEndcapName(hit_detId).pannelName();
00313             module = PixelEndcapName(hit_detId).plaquetteName();
00314           }
00315           if(layer==1){
00316             if(fabs(trackref->dxy(bestVtx->position()))>0.01 ||
00317                fabs(trackref->dz(bestVtx->position()))>0.1) continue;
00318             if(!(L2hits>0&&L3hits>0) && !(L2hits>0&&D1hits>0) && !(D1hits>0&&D2hits>0)) continue;
00319           }else if(layer==2){
00320             if(fabs(trackref->dxy(bestVtx->position()))>0.02 ||
00321                fabs(trackref->dz(bestVtx->position()))>0.1) continue;
00322             if(!(L1hits>0&&L3hits>0) && !(L1hits>0&&D1hits>0)) continue;
00323           }else if(layer==3){
00324             if(fabs(trackref->dxy(bestVtx->position()))>0.02 ||
00325                fabs(trackref->dz(bestVtx->position()))>0.1) continue;
00326             if(!(L1hits>0&&L2hits>0)) continue;
00327           }else if(disk==1){
00328             if(fabs(trackref->dxy(bestVtx->position()))>0.05 ||
00329                fabs(trackref->dz(bestVtx->position()))>0.5) continue;
00330             if(!(L1hits>0&&D2hits>0) && !(L2hits>0&&D2hits>0)) continue;
00331           }else if(disk==2){
00332             if(fabs(trackref->dxy(bestVtx->position()))>0.05 ||
00333                fabs(trackref->dz(bestVtx->position()))>0.5) continue;
00334             if(!(L1hits>0&&D1hits>0)) continue;
00335           }
00336           
00337               //check wether hit is valid or missing using track algo flag
00338           bool isHitValid   =hit->hit()->getType()==TrackingRecHit::valid;
00339           bool isHitMissing =hit->hit()->getType()==TrackingRecHit::missing;
00340           //std::cout<<"------ New Hit"<<std::endl;
00341           //std::cout<<(hit->hit()->getType()==TrackingRecHit::missing)<<std::endl;
00342           
00343           // get the enclosed persistent hit
00344           //const TrackingRecHit *persistentHit = hit->hit();
00345           // check if it's not null, and if it's a valid pixel hit
00346           //if ((persistentHit != 0) && (typeid(*persistentHit) == typeid(SiPixelRecHit))) {
00347           
00348             if(debug_) std::cout << "the hit is persistent\n";
00349             
00350             // tell the C++ compiler that the hit is a pixel hit
00351             //const SiPixelRecHit* pixhit = dynamic_cast<const SiPixelRecHit*>( hit->hit() );
00352          
00353             //define tracker and pixel geometry and topology
00354             const TrackerGeometry& theTracker(*theTrackerGeometry);
00355             const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*> (theTracker.idToDet(hit_detId) );
00356             //test if PixelGeomDetUnit exists
00357             if(theGeomDet == 0) {
00358               if(debug_) std::cout << "NO THEGEOMDET\n";
00359               continue;
00360             }
00361                       
00362               //const RectangularPixelTopology * topol = dynamic_cast<const RectangularPixelTopology*>(&(theGeomDet->specificTopology()));
00363 
00364               std::map<uint32_t, SiPixelHitEfficiencyModule*>::iterator pxd = theSiPixelStructure.find((*hit).geographicalId().rawId());
00365 
00366               // calculate alpha and beta from cluster position
00367               LocalTrajectoryParameters ltp = tsos.localParameters();
00368               //LocalVector localDir = ltp.momentum()/ltp.momentum().mag();
00369               
00370               //float clust_alpha = atan2(localDir.z(), localDir.x());
00371               //float clust_beta = atan2(localDir.z(), localDir.y());
00372               
00373               
00374               
00375               // THE CUTS
00376               //int nrows = theGeomDet->specificTopology().nrows();
00377               //int ncols = theGeomDet->specificTopology().ncolumns();
00378               //
00379               //std::pair<float,float> pitchTest = theGeomDet->specificTopology().pitch();
00380               //RectangularPixelTopology rectTopolTest = RectangularPixelTopology(nrows, ncols, pitch.first, pitch.second);
00381               //std::pair<float,float> pixelTest = rectTopol.pixel(tsos.localPosition());
00382               //
00383               
00384               
00385               //*************** Edge cut ********************
00386                //double glx=tsos.globalPosition().x();
00387                //double gly=tsos.globalPosition().y();
00388                //double glz=tsos.globalPosition().z();
00389                double lx=tsos.localPosition().x();
00390                double ly=tsos.localPosition().y();
00391                //double lz=tsos.localPosition().z();
00392                //double lx_err=tsos.localError().positionError().xx();
00393                //double ly_err=tsos.localError().positionError().yy();
00394                //int telescope=0; int telescope_valid=0; 
00395                if(fabs(lx)>0.55 || fabs(ly)>3.0) continue;
00396                //LocalTrajectoryParameters predTrajParam=tsos.localParameters();
00397                //LocalVector dir=predTrajParam.momentum()/predTrajParam.momentum().mag();
00398                //double alpha=atan2(dir.z(), dir.x());
00399                //double beta=atan2(dir.z(), dir.y());
00400                bool passedFiducial=true;
00401                // Module fiducials:
00402                if(IntSubDetID==PixelSubdetector::PixelBarrel && fabs(ly)>=3.1) passedFiducial=false;
00403                if(IntSubDetID==PixelSubdetector::PixelEndcap &&
00404                   !((panel==1 &&
00405                     ((module==1 && fabs(ly)<0.7) ||
00406                      ((module==2 && fabs(ly)<1.1) &&
00407                       !(disk==-1 && ly>0.8 && lx>0.2) &&
00408                       !(disk==1 && ly<-0.7 && lx>0.2) &&
00409                       !(disk==2 && ly<-0.8)) ||
00410                      ((module==3 && fabs(ly)<1.5) &&
00411                       !(disk==-2 && lx>0.1 && ly>1.0) &&
00412                       !(disk==2 && lx>0.1 && ly<-1.0)) ||
00413                      ((module==4 && fabs(ly)<1.9) &&
00414                       !(disk==-2 && ly>1.5) &&
00415                       !(disk==2 && ly<-1.5)))) ||
00416                     (panel==2 &&
00417                     ((module==1 && fabs(ly)<0.7) ||
00418                      (module==2 && fabs(ly)<1.2 &&
00419                       !(disk>0 && ly>1.1) &&
00420                       !(disk<0 && ly<-1.1)) ||
00421                      (module==3 && fabs(ly)<1.6 &&
00422                       !(disk>0 && ly>1.5) &&
00423                       !(disk<0 && ly<-1.5)))))) passedFiducial=false;
00424                if(IntSubDetID==PixelSubdetector::PixelEndcap &&
00425                   ((panel==1 && (module==1 || (module>=3 && abs(disk)==1))) ||
00426                    (panel==2 && ((module==1 && abs(disk)==2) ||
00427                                  (module==3 && abs(disk)==1))))) passedFiducial=false;
00428                // ROC fiducials:
00429                double ly_mod = fabs(ly);
00430                if(IntSubDetID==PixelSubdetector::PixelEndcap && (panel+module)%2==1) ly_mod=ly_mod+0.405;
00431                float d_rocedge = fabs(fmod(ly_mod,0.81)-0.405);
00432                if(d_rocedge<=0.0625) passedFiducial=false;
00433                if(!( (IntSubDetID==PixelSubdetector::PixelBarrel &&
00434                       ((!isHalfModule && fabs(lx)<0.6) ||
00435                        (isHalfModule && lx>-0.3 && lx<0.2))) ||
00436                      (IntSubDetID==PixelSubdetector::PixelEndcap &&
00437                       ((panel==1 &&
00438                        ((module==1 && fabs(lx)<0.2) ||
00439                         (module==2 &&
00440                          ((fabs(lx)<0.55 && abs(disk)==1) ||
00441                           (lx>-0.5 && lx<0.2 && disk==-2) ||
00442                           (lx>-0.5 && lx<0.0 && disk==2))) ||
00443                         (module==3 && lx>-0.6 && lx<0.5) ||
00444                         (module==4 && lx>-0.3 && lx<0.15))) ||
00445                        (panel==2 &&
00446                         ((module==1 && fabs(lx)<0.6) ||
00447                          (module==2 &&
00448                           ((fabs(lx)<0.55 && abs(disk)==1) ||
00449                            (lx>-0.6 && lx<0.5 && abs(disk)==2))) ||
00450                          (module==3 && fabs(lx)<0.5))))))) passedFiducial=false;
00451                if(((IntSubDetID==PixelSubdetector::PixelBarrel && !isHalfModule) || 
00452                    (IntSubDetID==PixelSubdetector::PixelEndcap && !(panel==1 && (module==1 || module==4)))) &&
00453                    fabs(lx)<0.06) passedFiducial=false;
00454                
00455                
00456               //*************** find closest clusters ********************
00457               float dx_cl[2]; float dy_cl[2]; dx_cl[0]=dx_cl[1]=dy_cl[0]=dy_cl[1]=-9999.;
00458               ESHandle<PixelClusterParameterEstimator> cpEstimator;
00459               iSetup.get<TkPixelCPERecord>().get("PixelCPEGeneric", cpEstimator);
00460               if(cpEstimator.isValid()){
00461                 const PixelClusterParameterEstimator &cpe(*cpEstimator);
00462                 edm::ESHandle<TrackerGeometry> tracker;
00463                 iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
00464                 if(tracker.isValid()){
00465                   const TrackerGeometry *tkgeom=&(*tracker);
00466                   edm::Handle<edmNew::DetSetVector<SiPixelCluster> > clusterCollectionHandle;
00467                   iEvent.getByLabel("siPixelClusters", clusterCollectionHandle);
00468                   if(clusterCollectionHandle.isValid()){
00469                     const edmNew::DetSetVector<SiPixelCluster>& clusterCollection=*clusterCollectionHandle;
00470                     edmNew::DetSetVector<SiPixelCluster>::const_iterator itClusterSet=clusterCollection.begin();
00471                     float minD[2]; minD[0]=minD[1]=10000.;
00472                     for( ; itClusterSet!=clusterCollection.end(); itClusterSet++){
00473                       DetId detId(itClusterSet->id());
00474                       if(detId.rawId()!=hit->geographicalId().rawId()) continue;
00475                       //unsigned int sdId=detId.subdetId();
00476                       const PixelGeomDetUnit *pixdet=(const PixelGeomDetUnit*) tkgeom->idToDetUnit(detId);
00477                       edmNew::DetSet<SiPixelCluster>::const_iterator itCluster=itClusterSet->begin();
00478                       for( ; itCluster!=itClusterSet->end(); ++itCluster){
00479                         LocalPoint lp(itCluster->x(), itCluster->y(), 0.);
00480                         PixelClusterParameterEstimator::LocalValues params=cpe.localParameters(*itCluster,*pixdet);
00481                         lp=params.first;
00482                         float D = sqrt((lp.x()-lx)*(lp.x()-lx)+(lp.y()-ly)*(lp.y()-ly));
00483                         if(D<minD[0]){
00484                           minD[1]=minD[0];
00485                           dx_cl[1]=dx_cl[0];
00486                           dy_cl[1]=dy_cl[0];
00487                           minD[0]=D;
00488                           dx_cl[0]=lp.x();
00489                           dy_cl[0]=lp.y();
00490                         }else if(D<minD[1]){
00491                           minD[1]=D;
00492                           dx_cl[1]=lp.x();
00493                           dy_cl[1]=lp.y();
00494                         }
00495                       }  // loop on clusterSets 
00496                     } // loop on clusterCollection
00497                     for(size_t i=0; i<2; i++){
00498                       if(minD[i]<9999.){
00499                         dx_cl[i]=fabs(dx_cl[i]-lx);
00500                         dy_cl[i]=fabs(dy_cl[i]-ly);
00501                       }
00502                     }
00503                   } // valid clusterCollectionHandle
00504                 } // valid tracker
00505               } // valid cpEstimator
00506               // distance of hit from closest cluster!
00507               float d_cl[2]; d_cl[0]=d_cl[1]=-9999.;
00508               if(dx_cl[0]!=-9999. && dy_cl[0]!=-9999.) d_cl[0]=sqrt(dx_cl[0]*dx_cl[0]+dy_cl[0]*dy_cl[0]);
00509               if(dx_cl[1]!=-9999. && dy_cl[1]!=-9999.) d_cl[1]=sqrt(dx_cl[1]*dx_cl[1]+dy_cl[1]*dy_cl[1]);
00510               if(isHitMissing && (d_cl[0]<0.05 || d_cl[1]<0.05)){ isHitMissing=0; isHitValid=1; }
00511               
00512               
00513               if(debug_){
00514                 std::cout << "Ready to add hit in histogram:\n";
00515                 //std::cout << "detid: "<<hit_detId<<std::endl;
00516                 std::cout << "isHitValid: "<<isHitValid<<std::endl;
00517                 std::cout << "isHitMissing: "<<isHitMissing<<std::endl;
00518                 //std::cout << "passedEdgeCut: "<<passedFiducial<<std::endl;            
00519               }
00520               
00521               
00522               if(pxd!=theSiPixelStructure.end() && isHitValid && passedFiducial)
00523                 ++nvalid;
00524               if(pxd!=theSiPixelStructure.end() && isHitMissing && passedFiducial)
00525                 ++nmissing;
00526                 
00527               if (pxd!=theSiPixelStructure.end() && passedFiducial && (isHitValid || isHitMissing))
00528                 (*pxd).second->fill(ltp, isHitValid, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);      
00529 
00530           //}//end if (persistent hit exists and is pixel hit)
00531           
00532         }//end of else 
00533         
00534         
00535       }//end for (all traj measurements of pixeltrack)
00536     }//end if (is pixeltrack)
00537     else
00538       if(debug_) std::cout << "no pixeltrack:\n";
00539     
00540   }//end loop on map entries
00541 }
00542 
00543 
00544 DEFINE_FWK_MODULE(SiPixelHitEfficiencySource); // define this as a plug-in