00001
00002
00003
00004
00005
00006
00008
00009
00010 #include <memory>
00011 #include <string>
00012 #include <iostream>
00013
00014 #include "FWCore/Framework/interface/Frameworkfwd.h"
00015 #include "FWCore/Framework/interface/EDAnalyzer.h"
00016 #include "FWCore/Framework/interface/Event.h"
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "CalibTracker/SiStripHitEfficiency/interface/HitEff.h"
00021 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00022 #include "DataFormats/Common/interface/Handle.h"
00023 #include "FWCore/Framework/interface/ESHandle.h"
00024 #include "FWCore/Framework/interface/EventSetup.h"
00025 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00026 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00027 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00028 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00029 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00030 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00031 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00032 #include "DataFormats/TrackReco/interface/Track.h"
00033 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00034 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00035 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00036 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00037 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00038 #include "CalibTracker/SiStripHitEfficiency/interface/TrajectoryAtInvalidHit.h"
00039 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00040 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00041 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00042 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00043 #include "RecoLocalTracker/ClusterParameterEstimator/interface/StripClusterParameterEstimator.h"
00044 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00045 #include "DataFormats/TrackReco/interface/DeDxData.h"
00046 #include "DataFormats/DetId/interface/DetIdCollection.h"
00047 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00048 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00049
00050 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00051 #include "AnalysisDataFormats/SiStripClusterInfo/interface/SiStripClusterInfo.h"
00052 #include "CalibTracker/Records/interface/SiStripDetCablingRcd.h"
00053 #include "CalibFormats/SiStripObjects/interface/SiStripDetCabling.h"
00054 #include "CalibTracker/Records/interface/SiStripQualityRcd.h"
00055 #include "CalibFormats/SiStripObjects/interface/SiStripQuality.h"
00056 #include "DataFormats/SiStripDetId/interface/SiStripSubStructure.h"
00057 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00058 #include "DataFormats/Common/interface/DetSetVector.h"
00059 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00060 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00061
00062 #include "DataFormats/MuonReco/interface/Muon.h"
00063 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00064
00065
00066 #include "TMath.h"
00067 #include "TH1F.h"
00068
00069
00070
00071
00072
00073 using namespace std;
00074 HitEff::HitEff(const edm::ParameterSet& conf) :
00075 conf_(conf)
00076 {
00077 layers =conf_.getParameter<int>("Layer");
00078 DEBUG = conf_.getParameter<bool>("Debug");
00079 }
00080
00081
00082 HitEff::~HitEff() { }
00083
00084 void HitEff::beginJob(){
00085
00086 edm::Service<TFileService> fs;
00087
00088 traj = fs->make<TTree>("traj","tree of trajectory positions");
00089 traj->Branch("TrajGlbX",&TrajGlbX,"TrajGlbX/F");
00090 traj->Branch("TrajGlbY",&TrajGlbY,"TrajGlbY/F");
00091 traj->Branch("TrajGlbZ",&TrajGlbZ,"TrajGlbZ/F");
00092 traj->Branch("TrajLocX",&TrajLocX,"TrajLocX/F");
00093 traj->Branch("TrajLocY",&TrajLocY,"TrajLocY/F");
00094 traj->Branch("TrajLocErrX",&TrajLocErrX,"TrajLocErrX/F");
00095 traj->Branch("TrajLocErrY",&TrajLocErrY,"TrajLocErrY/F");
00096 traj->Branch("TrajLocAngleX",&TrajLocAngleX,"TrajLocAngleX/F");
00097 traj->Branch("TrajLocAngleY",&TrajLocAngleY,"TrajLocAngleY/F");
00098 traj->Branch("ClusterLocX",&ClusterLocX,"ClusterLocX/F");
00099 traj->Branch("ClusterLocY",&ClusterLocY,"ClusterLocY/F");
00100 traj->Branch("ClusterLocErrX",&ClusterLocErrX,"ClusterLocErrX/F");
00101 traj->Branch("ClusterLocErrY",&ClusterLocErrY,"ClusterLocErrY/F");
00102 traj->Branch("ClusterStoN",&ClusterStoN,"ClusterStoN/F");
00103 traj->Branch("ResX",&ResX,"ResX/F");
00104 traj->Branch("ResXSig",&ResXSig,"ResXSig/F");
00105 traj->Branch("ModIsBad",&ModIsBad,"ModIsBad/i");
00106 traj->Branch("SiStripQualBad",&SiStripQualBad,"SiStripQualBad/i");
00107 traj->Branch("withinAcceptance",&withinAcceptance,"withinAcceptance/O");
00108 traj->Branch("nHits",&nHits,"nHits/I");
00109 traj->Branch("nLostHits",&nLostHits,"nLostHits/I");
00110 traj->Branch("chi2",&chi2,"chi2/F");
00111 traj->Branch("p",&p,"p/F");
00112 traj->Branch("pT",&pT,"pT/F");
00113 traj->Branch("trajHitValid", &trajHitValid, "trajHitValid/i");
00114 traj->Branch("Id",&Id,"Id/i");
00115 traj->Branch("run",&run,"run/i");
00116 traj->Branch("event",&event,"event/i");
00117 traj->Branch("layer",&whatlayer,"layer/i");
00118 traj->Branch("timeDT",&timeDT,"timeDT/F");
00119 traj->Branch("timeDTErr",&timeDTErr,"timeDTErr/F");
00120 traj->Branch("timeDTDOF",&timeDTDOF,"timeDTDOF/I");
00121 traj->Branch("timeECAL",&timeECAL,"timeECAL/F");
00122 traj->Branch("dedx",&dedx,"dedx/F");
00123 traj->Branch("dedxNOM",&dedxNOM,"dedxNOM/I");
00124 traj->Branch("tquality",&tquality,"tquality/I");
00125 traj->Branch("istep",&istep,"istep/I");
00126 traj->Branch("bunchx",&bunchx,"bunchx/I");
00127
00128 events = 0;
00129 EventTrackCKF = 0;
00130
00131 }
00132
00133
00134 void HitEff::analyze(const edm::Event& e, const edm::EventSetup& es){
00135
00136
00137
00138 if (DEBUG) cout << "beginning analyze from HitEff" << endl;
00139
00140 using namespace edm;
00141 using namespace reco;
00142
00143
00144 int run_nr = e.id().run();
00145 int ev_nr = e.id().event();
00146 int bunch_nr = e.bunchCrossing();
00147
00148
00149 edm::Handle<reco::TrackCollection> trackCollectionCKF;
00150 edm::InputTag TkTagCKF = conf_.getParameter<edm::InputTag>("combinatorialTracks");
00151 e.getByLabel(TkTagCKF,trackCollectionCKF);
00152
00153 edm::Handle<std::vector<Trajectory> > TrajectoryCollectionCKF;
00154 edm::InputTag TkTrajCKF = conf_.getParameter<edm::InputTag>("trajectories");
00155 e.getByLabel(TkTrajCKF,TrajectoryCollectionCKF);
00156
00157
00158
00159 edm::Handle< edmNew::DetSetVector<SiStripCluster> > theClusters;
00160 e.getByLabel("siStripClusters", theClusters);
00161
00162
00163 edm::ESHandle<TrackerGeometry> tracker;
00164 es.get<TrackerDigiGeometryRecord>().get(tracker);
00165 const TrackerGeometry * tkgeom=&(* tracker);
00166
00167
00168
00169 edm::ESHandle<StripClusterParameterEstimator> parameterestimator;
00170 es.get<TkStripCPERecord>().get("StripCPEfromTrackAngle", parameterestimator);
00171 const StripClusterParameterEstimator &stripcpe(*parameterestimator);
00172
00173
00174 edm::ESHandle<SiStripQuality> SiStripQuality_;
00175 try {
00176 es.get<SiStripQualityRcd>().get("forCluster",SiStripQuality_);
00177 }
00178 catch (...) {
00179 es.get<SiStripQualityRcd>().get(SiStripQuality_);
00180 }
00181
00182 edm::ESHandle<MagneticField> magFieldHandle;
00183 es.get<IdealMagneticFieldRecord>().get(magFieldHandle);
00184 const MagneticField* magField_ = magFieldHandle.product();
00185
00186
00187 edm::Handle< DetIdCollection > fedErrorIds;
00188 e.getByLabel("siStripDigis", fedErrorIds );
00189
00190 ESHandle<MeasurementTracker> measurementTrackerHandle;
00191 es.get<CkfComponentsRecord>().get(measurementTrackerHandle);
00192
00193 edm::ESHandle<Chi2MeasurementEstimatorBase> est;
00194 es.get<TrackingComponentsRecord>().get("Chi2",est);
00195
00196 edm::ESHandle<Propagator> prop;
00197 es.get<TrackingComponentsRecord>().get("PropagatorWithMaterial",prop);
00198 const Propagator* thePropagator = prop.product();
00199
00200 events++;
00201
00202
00203 const edmNew::DetSetVector<SiStripCluster>& input = *theClusters;
00204
00205
00206
00207
00208 for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
00209
00210 unsigned int ClusterId = DSViter->id();
00211 DetId ClusterDetId(ClusterId);
00212 const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
00213
00214 edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
00215 edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end();
00216 for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
00217
00218 StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
00219
00220 const Surface* surface;
00221 surface = &(tracker->idToDet(ClusterDetId)->surface());
00222 LocalPoint lp = parameters.first;
00223 GlobalPoint gp = surface->toGlobal(lp);
00224
00225
00226 }
00227 }
00228
00229
00230 const reco::TrackCollection *tracksCKF=trackCollectionCKF.product();
00231 if (DEBUG) cout << "number ckf tracks found = " << tracksCKF->size() << endl;
00232
00233 if (tracksCKF->size() > 0 && tracksCKF->size()<100) {
00234 if (DEBUG) cout << "starting checking good event with < 100 tracks" << endl;
00235
00236 EventTrackCKF++;
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255 Handle<MuonCollection> muH;
00256 if(e.getByLabel("muonsWitht0Correction",muH)){
00257 const MuonCollection & muonsT0 = *muH.product();
00258 if(muonsT0.size()!=0) {
00259 MuonTime mt0 = muonsT0[0].time();
00260 timeDT = mt0.timeAtIpInOut;
00261 timeDTErr = mt0.timeAtIpInOutErr;
00262 timeDTDOF = mt0.nDof;
00263
00264 bool hasCaloEnergyInfo = muonsT0[0].isEnergyValid();
00265 if (hasCaloEnergyInfo) timeECAL = muonsT0[0].calEnergy().ecal_time;
00266 }
00267 } else {
00268 timeDT = -999.0; timeDTErr = -999.0; timeDTDOF = -999; timeECAL = -999.0;
00269 }
00270
00271
00272
00273 for (vector<Trajectory>::const_iterator itraj = TrajectoryCollectionCKF.product()->begin();
00274 itraj != TrajectoryCollectionCKF.product()->end();
00275 itraj++) {
00276
00277
00278 nHits = itraj->foundHits();
00279 nLostHits = itraj->lostHits();
00280 chi2 = (itraj->chiSquared()/itraj->ndof());
00281 pT = sqrt( ( itraj->lastMeasurement().updatedState().globalMomentum().x() *
00282 itraj->lastMeasurement().updatedState().globalMomentum().x()) +
00283 ( itraj->lastMeasurement().updatedState().globalMomentum().y() *
00284 itraj->lastMeasurement().updatedState().globalMomentum().y()) );
00285 p = itraj->lastMeasurement().updatedState().globalMomentum().mag();
00286
00287
00288
00289
00290 std::vector<TrajectoryMeasurement> TMeas=itraj->measurements();
00291 vector<TrajectoryMeasurement>::iterator itm;
00292 double xloc = 0.;
00293 double yloc = 0.;
00294 double xErr = 0.;
00295 double yErr = 0.;
00296 double angleX = -999.;
00297 double angleY = -999.;
00298 double xglob,yglob,zglob;
00299
00300 for (itm=TMeas.begin();itm!=TMeas.end();itm++){
00301 ConstReferenceCountingPointer<TransientTrackingRecHit> theInHit;
00302 theInHit = (*itm).recHit();
00303
00304 if(DEBUG) cout << "theInHit is valid = " << theInHit->isValid() << endl;
00305
00306 unsigned int iidd = theInHit->geographicalId().rawId();
00307
00308 StripSubdetector strip=StripSubdetector(iidd);
00309 unsigned int TKlayers = checkLayer(iidd);
00310 if (DEBUG) cout << "TKlayer from trajectory: " << TKlayers << " from module = " << iidd << " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
00311
00312
00313 if ( TKlayers == 10 || TKlayers == 22 ) {
00314 if (DEBUG) cout << "skipping original TM for TOB 6 or TEC 9" << endl;
00315 continue;
00316 }
00317
00318
00319 std::vector<TrajectoryAtInvalidHit> TMs;
00320
00321
00322 AnalyticalPropagator propagator(magField_,anyDirection);
00323
00324
00325
00326
00327 if (isDoubleSided(iidd) && ((iidd & 0x3)==0) ) {
00328
00329
00330 TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 1));
00331 TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 2));
00332 } else if ( isDoubleSided(iidd) && (!check2DPartner(iidd, TMeas)) ) {
00333
00334
00335 TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 1));
00336 TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 2));
00337 if (DEBUG) cout << " found a hit with a missing partner" << endl;
00338 } else {
00339
00340 TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator));
00341 }
00342
00344
00345
00346
00347
00348
00349 bool isValid = theInHit->isValid();
00350 bool isLast = (itm==(TMeas.end()-1));
00351 bool isLastTOB5 = true;
00352 if (!isLast) {
00353 if ( checkLayer((++itm)->recHit()->geographicalId().rawId()) == 9 ) isLastTOB5 = false;
00354 else isLastTOB5 = true;
00355 --itm;
00356 }
00357
00358 if ( TKlayers==9 && isValid && isLastTOB5 ) {
00359
00360
00361 std::vector< BarrelDetLayer*> barrelTOBLayers = measurementTrackerHandle->geometricSearchTracker()->tobLayers() ;
00362 const DetLayer* tob6 = barrelTOBLayers[barrelTOBLayers.size()-1];
00363 const MeasurementEstimator* estimator = est.product();
00364 const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(&*measurementTrackerHandle);
00365 const TrajectoryStateOnSurface tsosTOB5 = itm->updatedState();
00366 vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements(*tob6, tsosTOB5, *thePropagator, *estimator);
00367
00368 if ( !tmp.empty()) {
00369 if (DEBUG) cout << "size of TM from propagation = " << tmp.size() << endl;
00370
00371
00372
00373
00374 TrajectoryMeasurement tob6TM(tmp.back());
00375 ConstReferenceCountingPointer<TransientTrackingRecHit> tob6Hit;
00376 tob6Hit = tob6TM.recHit();
00377
00378 if (tob6Hit->geographicalId().rawId()!=0) {
00379 if (DEBUG) cout << "tob6 hit actually being added to TM vector" << endl;
00380 TMs.push_back(TrajectoryAtInvalidHit(tob6TM,tkgeom, propagator));
00381 }
00382 }
00383 }
00384
00385 bool isLastTEC8 = true;
00386 if (!isLast) {
00387 if ( checkLayer((++itm)->recHit()->geographicalId().rawId()) == 21 ) isLastTEC8 = false;
00388 else isLastTEC8 = true;
00389 --itm;
00390 }
00391
00392 if ( TKlayers==21 && isValid && isLastTEC8 ) {
00393
00394 std::vector< ForwardDetLayer*> posTecLayers = measurementTrackerHandle->geometricSearchTracker()->posTecLayers() ;
00395 const DetLayer* tec9pos = posTecLayers[posTecLayers.size()-1];
00396 std::vector< ForwardDetLayer*> negTecLayers = measurementTrackerHandle->geometricSearchTracker()->negTecLayers() ;
00397 const DetLayer* tec9neg = negTecLayers[negTecLayers.size()-1];
00398
00399 const MeasurementEstimator* estimator = est.product();
00400 const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(&*measurementTrackerHandle);
00401 const TrajectoryStateOnSurface tsosTEC9 = itm->updatedState();
00402
00403
00404 if (!iidd == StripSubdetector::TEC) cout << "there is a problem with TEC 9 extrapolation" << endl;
00405 TECDetId tecdetid(iidd);
00406
00407 vector<TrajectoryMeasurement> tmp;
00408 if ( tecdetid.side() == 1 ) {
00409 tmp = theLayerMeasurements->measurements(*tec9neg, tsosTEC9, *thePropagator, *estimator);
00410
00411 }
00412 if ( tecdetid.side() == 2 ) {
00413 tmp = theLayerMeasurements->measurements(*tec9pos, tsosTEC9, *thePropagator, *estimator);
00414
00415 }
00416
00417 if ( !tmp.empty()) {
00418
00419
00420
00421 TrajectoryMeasurement tec9TM(tmp.back());
00422 ConstReferenceCountingPointer<TransientTrackingRecHit> tec9Hit;
00423 tec9Hit = tec9TM.recHit();
00424
00425 unsigned int tec9id = tec9Hit->geographicalId().rawId();
00426 if (DEBUG) cout << "tec9id = " << tec9id << " is Double sided = " << isDoubleSided(tec9id) << " and 0x3 = " << (tec9id & 0x3) << endl;
00427
00428 if (tec9Hit->geographicalId().rawId()!=0) {
00429 if (DEBUG) cout << "tec9 hit actually being added to TM vector" << endl;
00430
00431
00432 if (isDoubleSided(tec9id)) {
00433 TMs.push_back(TrajectoryAtInvalidHit(tec9TM,tkgeom, propagator, 1));
00434 TMs.push_back(TrajectoryAtInvalidHit(tec9TM,tkgeom, propagator, 2));
00435 } else
00436 TMs.push_back(TrajectoryAtInvalidHit(tec9TM,tkgeom, propagator));
00437 }
00438 }
00439 }
00440
00442
00443
00444
00445 for(std::vector<TrajectoryAtInvalidHit>::const_iterator TM=TMs.begin();TM!=TMs.end();++TM) {
00446
00447
00448 iidd = TM->monodet_id();
00449 if (DEBUG) cout << "setting iidd = " << iidd << " before checking efficiency and ";
00450
00451 xloc = TM->localX();
00452 yloc = TM->localY();
00453
00454 xErr = TM->localErrorX();
00455 yErr = TM->localErrorY();
00456
00457 angleX = atan( TM->localDxDz() );
00458 angleY = atan( TM->localDyDz() );
00459
00460 xglob = TM->globalX();
00461 yglob = TM->globalY();
00462 zglob = TM->globalZ();
00463 withinAcceptance = TM->withinAcceptance();
00464
00465 trajHitValid = TM->validHit();
00466
00467
00468 TKlayers = checkLayer(iidd);
00469
00470 if ( (layers == TKlayers) || (layers == 0) ) {
00471 whatlayer = TKlayers;
00472 if (DEBUG) cout << "Looking at layer under study" << endl;
00473 TrajGlbX = 0.0; TrajGlbY = 0.0; TrajGlbZ = 0.0; ModIsBad = 2; Id = 0; SiStripQualBad = 0;
00474 run = 0; event = 0; TrajLocX = 0.0; TrajLocY = 0.0; TrajLocErrX = 0.0; TrajLocErrY = 0.0;
00475 TrajLocAngleX = -999.0; TrajLocAngleY = -999.0; ResX = 0.0; ResXSig = 0.0;
00476 ClusterLocX = 0.0; ClusterLocY = 0.0; ClusterLocErrX = 0.0; ClusterLocErrY = 0.0; ClusterStoN = 0.0;
00477 bunchx = 0;
00478
00479
00480
00481 if (input.size() > 0 ) {
00482 if (DEBUG) cout << "Checking clusters with size = " << input.size() << endl;
00483 int nClusters = 0;
00484 std::vector< std::vector<float> > VCluster_info;
00485 for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
00486
00487
00488 unsigned int ClusterId = DSViter->id();
00489 if (ClusterId == iidd) {
00490 if (DEBUG) cout << "found (ClusterId == iidd) with ClusterId = " << ClusterId << " and iidd = " << iidd << endl;
00491 DetId ClusterDetId(ClusterId);
00492 const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
00493
00494 for(edmNew::DetSet<SiStripCluster>::const_iterator iter=DSViter->begin();iter!=DSViter->end();++iter) {
00495
00496 StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
00497 float res = (parameters.first.x() - xloc);
00498 float sigma = checkConsistency(parameters , xloc, xErr);
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509 std::vector< float > cluster_info;
00510 cluster_info.push_back(res);
00511 cluster_info.push_back(sigma);
00512 cluster_info.push_back(parameters.first.x());
00513 cluster_info.push_back(sqrt(parameters.second.xx()));
00514 cluster_info.push_back(parameters.first.y());
00515 cluster_info.push_back(sqrt(parameters.second.yy()));
00516
00517
00518
00519
00520 VCluster_info.push_back(cluster_info);
00521 nClusters++;
00522 if (DEBUG) cout << "Have ID match. residual = " << VCluster_info.back()[0] << " res sigma = " << VCluster_info.back()[1] << endl;
00523 if (DEBUG) cout << "trajectory measurement compatability estimate = " << (*itm).estimate() << endl;
00524 if (DEBUG) cout << "hit position = " << parameters.first.x() << " hit error = " << sqrt(parameters.second.xx()) << " trajectory position = " << xloc << " traj error = " << xErr << endl;
00525 }
00526 }
00527 }
00528 float FinalResSig = 1000.0;
00529 float FinalCluster[7]= {1000.0, 1000.0, 0.0, 0.0, 0.0, 0.0, 0.0};
00530 if (nClusters > 0) {
00531 if (DEBUG) cout << "found clusters > 0" << endl;
00532 if (nClusters > 1) {
00533
00534 vector< vector<float> >::iterator ires;
00535 for (ires=VCluster_info.begin(); ires!=VCluster_info.end(); ires++){
00536 if ( abs((*ires)[1]) < abs(FinalResSig)) {
00537 FinalResSig = (*ires)[1];
00538 for (unsigned int i = 0; i<ires->size(); i++) {
00539 if (DEBUG) cout << "filling final cluster. i = " << i << " before fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
00540 FinalCluster[i] = (*ires)[i];
00541 if (DEBUG) cout << "filling final cluster. i = " << i << " after fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
00542 }
00543 }
00544 if (DEBUG) cout << "iresidual = " << (*ires)[0] << " isigma = " << (*ires)[1] << " and FinalRes = " << FinalCluster[0] << endl;
00545 }
00546 }
00547 else {
00548 FinalResSig = VCluster_info.at(0)[1];
00549 for (unsigned int i = 0; i<VCluster_info.at(0).size(); i++) {
00550 FinalCluster[i] = VCluster_info.at(0)[i];
00551 }
00552 }
00553 nClusters=0;
00554 VCluster_info.clear();
00555 }
00556
00557 if (DEBUG) cout << "Final residual in X = " << FinalCluster[0] << "+-" << (FinalCluster[0]/FinalResSig) << endl;
00558 if (DEBUG) cout << "Checking location of trajectory: abs(yloc) = " << abs(yloc) << " abs(xloc) = " << abs(xloc) << endl;
00559 if (DEBUG) cout << "Checking location of cluster hit: yloc = " << FinalCluster[4] << "+-" << FinalCluster[5] << " xloc = " << FinalCluster[2] << "+-" << FinalCluster[3] << endl;
00560 if (DEBUG) cout << "Final cluster signal to noise = " << FinalCluster[6] << endl;
00561
00562 float exclusionWidth = 0.4;
00563 float TOBexclusion = 0.0;
00564 float TECexRing5 = -0.89;
00565 float TECexRing6 = -0.56;
00566 float TECexRing7 = 0.60;
00567
00568 int subdetector = ((iidd>>25) & 0x7);
00569 int ringnumber = ((iidd>>5) & 0x7);
00570
00571
00572 if((TKlayers >= 5 && TKlayers < 11)||((subdetector == 6)&&( (ringnumber >= 5)&&(ringnumber <=7) ))) {
00573
00574 float highzone = 0.0;
00575 float lowzone = 0.0;
00576 float higherr = yloc + 5.0*yErr;
00577 float lowerr = yloc - 5.0*yErr;
00578 if(TKlayers >= 5 && TKlayers < 11) {
00579
00580 highzone = TOBexclusion + exclusionWidth;
00581 lowzone = TOBexclusion - exclusionWidth;
00582 } else if (ringnumber == 5) {
00583
00584 highzone = TECexRing5 + exclusionWidth;
00585 lowzone = TECexRing5 - exclusionWidth;
00586 } else if (ringnumber == 6) {
00587
00588 highzone = TECexRing6 + exclusionWidth;
00589 lowzone = TECexRing6 - exclusionWidth;
00590 } else if (ringnumber == 7) {
00591
00592 highzone = TECexRing7 + exclusionWidth;
00593 lowzone = TECexRing7 - exclusionWidth;
00594 }
00595
00596 if((highzone <= higherr)&&(highzone >= lowerr)) withinAcceptance = false;
00597 if((lowzone >= lowerr)&&(lowzone <= higherr)) withinAcceptance = false;
00598 if((higherr <= highzone)&&(higherr >= lowzone)) withinAcceptance = false;
00599 if((lowerr >= lowzone)&&(lowerr <= highzone)) withinAcceptance = false;
00600 }
00601
00602
00603
00604 TrajGlbX = xglob;
00605 TrajGlbY = yglob;
00606 TrajGlbZ = zglob;
00607 Id = iidd;
00608 run = run_nr;
00609 event = ev_nr;
00610 bunchx = bunch_nr;
00611
00612 if ( SiStripQuality_->getBadApvs(iidd)!=0 ) {
00613 SiStripQualBad = 1;
00614 if(DEBUG) cout << "strip is bad from SiStripQuality" << endl;
00615 } else {
00616 SiStripQualBad = 0;
00617 if(DEBUG) cout << "strip is good from SiStripQuality" << endl;
00618 }
00619
00620
00621 for (unsigned int ii=0;ii< fedErrorIds->size();ii++) {
00622 if (iidd == (*fedErrorIds)[ii].rawId() )
00623 SiStripQualBad = 1;
00624 }
00625
00626 TrajLocX = xloc;
00627 TrajLocY = yloc;
00628 TrajLocErrX = xErr;
00629 TrajLocErrY = yErr;
00630 TrajLocAngleX = angleX;
00631 TrajLocAngleY = angleY;
00632 ResX = FinalCluster[0];
00633 ResXSig = FinalResSig;
00634 if (FinalResSig != FinalCluster[1]) if (DEBUG) cout << "Problem with best cluster selection because FinalResSig = " << FinalResSig << " and FinalCluster[1] = " << FinalCluster[1] << endl;
00635 ClusterLocX = FinalCluster[2];
00636 ClusterLocY = FinalCluster[4];
00637 ClusterLocErrX = FinalCluster[3];
00638 ClusterLocErrY = FinalCluster[5];
00639 ClusterStoN = FinalCluster[6];
00640
00641 if (DEBUG) cout << "before check good" << endl;
00642
00643 if ( FinalResSig < 999.0) {
00644
00645 if (DEBUG) cout << "hit being counted as good " << FinalCluster[0] << " FinalRecHit " <<
00646 iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd <<
00647 " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
00648 ModIsBad = 0;
00649 traj->Fill();
00650 }
00651 else {
00652 if (DEBUG) cout << "hit being counted as bad ######### Invalid RPhi FinalResX " << FinalCluster[0] << " FinalRecHit " <<
00653 iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd <<
00654 " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
00655 ModIsBad = 1;
00656 traj->Fill();
00657
00658 if (DEBUG) cout << " RPhi Error " << sqrt(xErr*xErr + yErr*yErr) << " ErrorX " << xErr << " yErr " << yErr << endl;
00659 } if (DEBUG) cout << "after good location check" << endl;
00660 } if (DEBUG) cout << "after list of clusters" << endl;
00661 } if (DEBUG) cout << "After layers=TKLayers if" << endl;
00662 } if (DEBUG) cout << "After looping over TrajAtValidHit list" << endl;
00663 } if (DEBUG) cout << "end TMeasurement loop" << endl;
00664 } if (DEBUG) cout << "end of trajectories loop" << endl;
00665 }
00666 }
00667
00668 void HitEff::endJob(){
00669 traj->GetDirectory()->cd();
00670 traj->Write();
00671
00672 cout << " Events Analysed " << events << endl;
00673 cout << " Number Of Tracked events " << EventTrackCKF << endl;
00674 }
00675
00676 double HitEff::checkConsistency(StripClusterParameterEstimator::LocalValues parameters, double xx, double xerr) {
00677 double error = sqrt(parameters.second.xx() + xerr*xerr);
00678 double separation = abs(parameters.first.x() - xx);
00679 double consistency = separation/error;
00680 return consistency;
00681 }
00682
00683 bool HitEff::isDoubleSided(unsigned int iidd) const {
00684 StripSubdetector strip=StripSubdetector(iidd);
00685 unsigned int subid=strip.subdetId();
00686 unsigned int layer = 0;
00687 if (subid == StripSubdetector::TIB) {
00688 TIBDetId tibid(iidd);
00689 layer = tibid.layer();
00690 if (layer == 1 || layer == 2) return true;
00691 else return false;
00692 }
00693 else if (subid == StripSubdetector::TOB) {
00694 TOBDetId tobid(iidd);
00695 layer = tobid.layer() + 4 ;
00696 if (layer == 5 || layer == 6) return true;
00697 else return false;
00698 }
00699 else if (subid == StripSubdetector::TID) {
00700 TIDDetId tidid(iidd);
00701 layer = tidid.ring() + 10;
00702 if (layer == 11 || layer == 12) return true;
00703 else return false;
00704 }
00705 else if (subid == StripSubdetector::TEC) {
00706 TECDetId tecid(iidd);
00707 layer = tecid.ring() + 13 ;
00708 if (layer == 14 || layer == 15 || layer == 18) return true;
00709 else return false;
00710 }
00711 else
00712 return false;
00713 }
00714
00715 bool HitEff::check2DPartner(unsigned int iidd, std::vector<TrajectoryMeasurement> traj) {
00716 unsigned int partner_iidd = 0;
00717 bool found2DPartner = false;
00718
00719 if ((iidd & 0x3)==1) partner_iidd = iidd+1;
00720 if ((iidd & 0x3)==2) partner_iidd = iidd-1;
00721
00722
00723 for (std::vector<TrajectoryMeasurement>::const_iterator iTM=traj.begin(); iTM!=traj.end(); ++iTM) {
00724 if (iTM->recHit()->geographicalId().rawId()==partner_iidd) {
00725 found2DPartner = true;
00726 }
00727 }
00728 return found2DPartner;
00729 }
00730
00731 unsigned int HitEff::checkLayer( unsigned int iidd) {
00732 StripSubdetector strip=StripSubdetector(iidd);
00733 unsigned int subid=strip.subdetId();
00734 if (subid == StripSubdetector::TIB) {
00735 TIBDetId tibid(iidd);
00736 return tibid.layer();
00737 }
00738 if (subid == StripSubdetector::TOB) {
00739 TOBDetId tobid(iidd);
00740 return tobid.layer() + 4 ;
00741 }
00742 if (subid == StripSubdetector::TID) {
00743 TIDDetId tidid(iidd);
00744 return tidid.wheel() + 10;
00745 }
00746 if (subid == StripSubdetector::TEC) {
00747 TECDetId tecid(iidd);
00748 return tecid.wheel() + 13 ;
00749 }
00750 return 0;
00751 }
00752
00753
00754 DEFINE_FWK_MODULE(HitEff);