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