CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/CalibTracker/SiStripHitEfficiency/src/HitEff.cc

Go to the documentation of this file.
00001 
00002 // Package:          CalibTracker/SiStripHitEfficiency
00003 // Class:            HitEff
00004 // Original Author:  Keith Ulmer--University of Colorado
00005 //                   keith.ulmer@colorado.edu
00006 //
00008 
00009 // system include files
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 // constructors and destructor
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 // Virtual destructor needed.
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   //  bool DEBUG = false;
00137 
00138   if (DEBUG)  cout << "beginning analyze from HitEff" << endl;
00139 
00140   using namespace edm;
00141   using namespace reco;
00142   // Step A: Get Inputs 
00143 
00144   int run_nr = e.id().run();
00145   int ev_nr = e.id().event();
00146   int bunch_nr = e.bunchCrossing();
00147 
00148   //CombinatoriaTrack
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   // Clusters
00158   // get the SiStripClusters from the event
00159   edm::Handle< edmNew::DetSetVector<SiStripCluster> > theClusters;
00160   e.getByLabel("siStripClusters", theClusters);
00161 
00162   //get tracker geometry
00163   edm::ESHandle<TrackerGeometry> tracker;
00164   es.get<TrackerDigiGeometryRecord>().get(tracker);
00165   const TrackerGeometry * tkgeom=&(* tracker);
00166 
00167   //get Cluster Parameter Estimator
00168   //std::string cpe = conf_.getParameter<std::string>("StripCPE");
00169   edm::ESHandle<StripClusterParameterEstimator> parameterestimator;
00170   es.get<TkStripCPERecord>().get("StripCPEfromTrackAngle", parameterestimator); 
00171   const StripClusterParameterEstimator &stripcpe(*parameterestimator);
00172 
00173   // get the SiStripQuality records
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   // get the list of module IDs with FED-detected errors
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   // *************** SiStripCluster Collection
00203   const edmNew::DetSetVector<SiStripCluster>& input = *theClusters;
00204 
00205   //go through clusters to write out global position of good clusters for the layer understudy for comparison
00206   // Loop through clusters just to print out locations
00207   // Commented out to avoid discussion, should really be deleted.
00208   /*
00209   for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
00210     // DSViter is a vector of SiStripClusters located on a single module
00211     unsigned int ClusterId = DSViter->id();
00212     DetId ClusterDetId(ClusterId);
00213     const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
00214     
00215     edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
00216     edmNew::DetSet<SiStripCluster>::const_iterator end  =DSViter->end();
00217     for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
00218       //iter is a single SiStripCluster
00219       StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
00220       
00221       const Surface* surface;
00222       surface = &(tracker->idToDet(ClusterDetId)->surface());
00223       LocalPoint lp = parameters.first;
00224       GlobalPoint gp = surface->toGlobal(lp);
00225       unsigned int layer = checkLayer(ClusterId);
00226             if(DEBUG) cout << "Found hit in cluster collection layer = " << layer << " with id = " << ClusterId << "   local X position = " << lp.x() << " +- " << sqrt(parameters.second.xx()) << "   matched/stereo/rphi = " << ((ClusterId & 0x3)==0) << "/" << ((ClusterId & 0x3)==1) << "/" << ((ClusterId & 0x3)==2) << endl;
00227     }
00228   }
00229   */
00230   
00231   // Tracking 
00232   const   reco::TrackCollection *tracksCKF=trackCollectionCKF.product();
00233   if (DEBUG)  cout << "number ckf tracks found = " << tracksCKF->size() << endl;
00234   //if (tracksCKF->size() == 1 ){
00235   if (tracksCKF->size() > 0 && tracksCKF->size()<100) {
00236     if (DEBUG)    cout << "starting checking good event with < 100 tracks" << endl;
00237 
00238     EventTrackCKF++;  
00239 
00240     /*
00241 
00242     //get dEdx info if available
00243     Handle<ValueMap<DeDxData> >          dEdxUncalibHandle;
00244     if (e.getByLabel("dedxMedianCTF", dEdxUncalibHandle)) {
00245       const ValueMap<DeDxData> dEdxTrackUncalib = *dEdxUncalibHandle.product();
00246       
00247       reco::TrackRef itTrack  = reco::TrackRef( trackCollectionCKF, 0 );
00248       dedx = dEdxTrackUncalib[itTrack].dEdx();
00249       dedxNOM  = dEdxTrackUncalib[itTrack].numberOfMeasurements();
00250     } else {
00251       dedx = -999.0; dedxNOM = -999;
00252     }
00253 
00254     */
00255 
00256     //get muon and ecal timing info if available
00257     Handle<MuonCollection> muH;
00258     if(e.getByLabel("muonsWitht0Correction",muH)){
00259       const MuonCollection & muonsT0  =  *muH.product();
00260       if(muonsT0.size()!=0) {
00261         MuonTime mt0 = muonsT0[0].time();
00262         timeDT = mt0.timeAtIpInOut; 
00263         timeDTErr = mt0.timeAtIpInOutErr;
00264         timeDTDOF = mt0.nDof;
00265         
00266         bool hasCaloEnergyInfo = muonsT0[0].isEnergyValid();
00267         if (hasCaloEnergyInfo) timeECAL = muonsT0[0].calEnergy().ecal_time;
00268       }
00269     } else {
00270       timeDT = -999.0; timeDTErr = -999.0; timeDTDOF = -999; timeECAL = -999.0;
00271     }
00272 
00273     // actually should do a loop over all the tracks in the event here
00274 
00275     for (vector<Trajectory>::const_iterator itraj = TrajectoryCollectionCKF.product()->begin();
00276          itraj != TrajectoryCollectionCKF.product()->end();
00277          itraj++) {
00278 
00279       // for each track, fill some variables such as number of hits and momentum
00280       nHits = itraj->foundHits();
00281       nLostHits = itraj->lostHits();
00282       chi2 = (itraj->chiSquared()/itraj->ndof());
00283       pT = sqrt( ( itraj->lastMeasurement().updatedState().globalMomentum().x() *
00284                    itraj->lastMeasurement().updatedState().globalMomentum().x()) +
00285                  ( itraj->lastMeasurement().updatedState().globalMomentum().y() *
00286                    itraj->lastMeasurement().updatedState().globalMomentum().y()) );
00287       p = itraj->lastMeasurement().updatedState().globalMomentum().mag();
00288       
00289       //Put in code to check track quality
00290       
00291       
00292       std::vector<TrajectoryMeasurement> TMeas=itraj->measurements();
00293       vector<TrajectoryMeasurement>::iterator itm;
00294       double xloc = 0.;
00295       double yloc = 0.;
00296       double xErr = 0.;
00297       double yErr = 0.;
00298       double angleX = -999.;
00299       double angleY = -999.;
00300       double xglob,yglob,zglob;
00301       
00302       for (itm=TMeas.begin();itm!=TMeas.end();itm++){
00303         ConstReferenceCountingPointer<TransientTrackingRecHit> theInHit;
00304         theInHit = (*itm).recHit();
00305         
00306         if(DEBUG) cout << "theInHit is valid = " << theInHit->isValid() << endl;
00307         
00308         unsigned int iidd = theInHit->geographicalId().rawId();
00309         
00310         unsigned int TKlayers = checkLayer(iidd);
00311         if (DEBUG)      cout << "TKlayer from trajectory: " << TKlayers << "  from module = " << iidd <<  "   matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
00312 
00313         // If Trajectory measurement from TOB 6 or TEC 9, skip it because it's always valid they are filled later
00314         if ( TKlayers == 10 || TKlayers == 22 ) {
00315           if (DEBUG) cout << "skipping original TM for TOB 6 or TEC 9" << endl;
00316           continue;
00317         }
00318 
00319         // Make vector of TrajectoryAtInvalidHits to hold the trajectories
00320         std::vector<TrajectoryAtInvalidHit> TMs;
00321         
00322         // Make AnalyticalPropagator to use in TAVH constructor
00323         AnalyticalPropagator propagator(magField_,anyDirection); 
00324         
00325         // for double sided layers check both sensors--if no hit was found on either sensor surface,
00326         // the trajectory measurements only have one invalid hit entry on the matched surface
00327         // so get the TrajectoryAtInvalidHit for both surfaces and include them in the study
00328         if (isDoubleSided(iidd) &&  ((iidd & 0x3)==0) ) {
00329           // do hit eff check twice--once for each sensor
00330           //add a TM for each surface
00331           TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 1));
00332           TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 2));
00333         } else if ( isDoubleSided(iidd) && (!check2DPartner(iidd, TMeas)) ) {
00334           // if only one hit was found the trajectory measurement is on that sensor surface, and the other surface from
00335           // the matched layer should be added to the study as well
00336           TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 1));
00337           TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator, 2));
00338           if (DEBUG) cout << " found a hit with a missing partner" << endl;
00339         } else {
00340           //only add one TM for the single surface and the other will be added in the next iteration
00341           TMs.push_back(TrajectoryAtInvalidHit(*itm,tkgeom, propagator));
00342         }
00343 
00345         //Now check for tracks at TOB6 and TEC9
00346 
00347         // to make sure we only propagate on the last TOB5 hit check the next entry isn't also in TOB5
00348         // to avoid bias, make sure the TOB5 hit is valid (an invalid hit on TOB5 could only exist with a valid hit on TOB6)
00349 
00350         bool isValid = theInHit->isValid();
00351         bool isLast = (itm==(TMeas.end()-1));
00352         bool isLastTOB5 = true;
00353         if (!isLast) {
00354           if ( checkLayer((++itm)->recHit()->geographicalId().rawId()) == 9 ) isLastTOB5 = false;
00355           else isLastTOB5 = true;
00356           --itm;
00357         }
00358         
00359         if ( TKlayers==9 && isValid && isLastTOB5 ) {
00360           //      if ( TKlayers==9 && itm==TMeas.rbegin()) {
00361         //        if ( TKlayers==9 && (itm==TMeas.back()) ) {     // to check for only the last entry in the trajectory for propagation
00362           std::vector< BarrelDetLayer*> barrelTOBLayers = measurementTrackerHandle->geometricSearchTracker()->tobLayers() ;
00363           const DetLayer* tob6 = barrelTOBLayers[barrelTOBLayers.size()-1];
00364           const MeasurementEstimator* estimator = est.product();
00365           const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(&*measurementTrackerHandle);
00366           const TrajectoryStateOnSurface tsosTOB5 = itm->updatedState();
00367           vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements(*tob6, tsosTOB5, *thePropagator, *estimator);
00368           
00369           if ( !tmp.empty()) {
00370             if (DEBUG) cout << "size of TM from propagation = " << tmp.size() << endl;
00371 
00372             // take the last of the TMs, which is always an invalid hit
00373             // if no detId is available, ie detId==0, then no compatible layer was crossed
00374             // otherwise, use that TM for the efficiency measurement
00375             TrajectoryMeasurement tob6TM(tmp.back());
00376             ConstReferenceCountingPointer<TransientTrackingRecHit> tob6Hit;
00377             tob6Hit = tob6TM.recHit();
00378             
00379             if (tob6Hit->geographicalId().rawId()!=0) {
00380               if (DEBUG) cout << "tob6 hit actually being added to TM vector" << endl;
00381               TMs.push_back(TrajectoryAtInvalidHit(tob6TM,tkgeom, propagator));
00382             }
00383           }
00384         }
00385 
00386         bool isLastTEC8 = true;
00387         if (!isLast) {
00388           if ( checkLayer((++itm)->recHit()->geographicalId().rawId()) == 21 ) isLastTEC8 = false;
00389           else isLastTEC8 = true;
00390           --itm;
00391         }
00392         
00393         if ( TKlayers==21 && isValid && isLastTEC8 ) {
00394           
00395           std::vector< ForwardDetLayer*> posTecLayers = measurementTrackerHandle->geometricSearchTracker()->posTecLayers() ;
00396           const DetLayer* tec9pos = posTecLayers[posTecLayers.size()-1];
00397           std::vector< ForwardDetLayer*> negTecLayers = measurementTrackerHandle->geometricSearchTracker()->negTecLayers() ;
00398           const DetLayer* tec9neg = negTecLayers[negTecLayers.size()-1];
00399           
00400           const MeasurementEstimator* estimator = est.product();
00401           const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(&*measurementTrackerHandle);
00402           const TrajectoryStateOnSurface tsosTEC9 = itm->updatedState();
00403           
00404           // check if track on positive or negative z
00405           if (!iidd ==  StripSubdetector::TEC) cout << "there is a problem with TEC 9 extrapolation" << endl;
00406           TECDetId tecdetid(iidd);
00407           //cout << " tec9 id = " << iidd << " and side = " << tecdetid.side() << endl;
00408           vector<TrajectoryMeasurement> tmp;
00409           if ( tecdetid.side() == 1 ) {
00410             tmp = theLayerMeasurements->measurements(*tec9neg, tsosTEC9, *thePropagator, *estimator);
00411             //cout << "on negative side" << endl;
00412           }
00413           if ( tecdetid.side() == 2 ) {
00414             tmp = theLayerMeasurements->measurements(*tec9pos, tsosTEC9, *thePropagator, *estimator);
00415             //cout << "on positive side" << endl;
00416           }
00417 
00418           if ( !tmp.empty()) {
00419             // take the last of the TMs, which is always an invalid hit
00420             // if no detId is available, ie detId==0, then no compatible layer was crossed
00421             // otherwise, use that TM for the efficiency measurement
00422             TrajectoryMeasurement tec9TM(tmp.back());
00423             ConstReferenceCountingPointer<TransientTrackingRecHit> tec9Hit;
00424             tec9Hit = tec9TM.recHit();
00425             
00426             unsigned int tec9id = tec9Hit->geographicalId().rawId();
00427             if (DEBUG) cout << "tec9id = " << tec9id << " is Double sided = " <<  isDoubleSided(tec9id) << "  and 0x3 = " << (tec9id & 0x3) << endl;
00428             
00429             if (tec9Hit->geographicalId().rawId()!=0) {
00430               if (DEBUG) cout << "tec9 hit actually being added to TM vector" << endl;
00431               // in tec the hit can be single or doubled sided. whenever the invalid hit at the end of vector of TMs is
00432               // double sided it is always on the matched surface, so we need to split it into the true sensor surfaces
00433               if (isDoubleSided(tec9id)) {
00434                 TMs.push_back(TrajectoryAtInvalidHit(tec9TM,tkgeom, propagator, 1));
00435                 TMs.push_back(TrajectoryAtInvalidHit(tec9TM,tkgeom, propagator, 2));
00436               } else 
00437                 TMs.push_back(TrajectoryAtInvalidHit(tec9TM,tkgeom, propagator));
00438             }
00439           } //else cout << "tec9 tmp empty" << endl;
00440         }
00441         
00443         
00444         // Modules Constraints
00445 
00446         for(std::vector<TrajectoryAtInvalidHit>::const_iterator TM=TMs.begin();TM!=TMs.end();++TM) {
00447           
00448           // --> Get trajectory from combinatedState 
00449           iidd = TM->monodet_id();
00450           if (DEBUG) cout << "setting iidd = " << iidd << " before checking efficiency and ";
00451           
00452           xloc = TM->localX();
00453           yloc = TM->localY();
00454           
00455           xErr =  TM->localErrorX();
00456           yErr =  TM->localErrorY();
00457           
00458           angleX = atan( TM->localDxDz() );
00459           angleY = atan( TM->localDyDz() );
00460           
00461           xglob = TM->globalX();
00462           yglob = TM->globalY();
00463           zglob = TM->globalZ();
00464           withinAcceptance = TM->withinAcceptance();
00465           
00466           trajHitValid = TM->validHit();
00467 
00468           // reget layer from iidd here, to account for TOB 6 and TEC 9 TKlayers being off
00469           TKlayers = checkLayer(iidd);
00470 
00471           if ( (layers == TKlayers) || (layers == 0) ) {   // Look at the layer not used to reconstruct the track
00472             whatlayer = TKlayers;
00473             if (DEBUG)    cout << "Looking at layer under study" << endl;
00474             TrajGlbX = 0.0; TrajGlbY = 0.0; TrajGlbZ = 0.0; ModIsBad = 2; Id = 0; SiStripQualBad = 0; 
00475             run = 0; event = 0; TrajLocX = 0.0; TrajLocY = 0.0; TrajLocErrX = 0.0; TrajLocErrY = 0.0; 
00476             TrajLocAngleX = -999.0; TrajLocAngleY = -999.0;     ResX = 0.0; ResXSig = 0.0;
00477             ClusterLocX = 0.0; ClusterLocY = 0.0; ClusterLocErrX = 0.0; ClusterLocErrY = 0.0; ClusterStoN = 0.0;
00478             bunchx = 0;
00479             
00480             // RPhi RecHit Efficiency 
00481             
00482             if (input.size() > 0 ) {  
00483               if (DEBUG) cout << "Checking clusters with size = " << input.size() << endl;
00484               int nClusters = 0;
00485               std::vector< std::vector<float> > VCluster_info; //fill with X residual, X residual pull, local X, sig(X), local Y, sig(Y), StoN
00486               for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
00487                 // DSViter is a vector of SiStripClusters located on a single module
00488                 //if (DEBUG)      cout << "the ID from the DSViter = " << DSViter->id() << endl; 
00489                 unsigned int ClusterId = DSViter->id();
00490                 if (ClusterId == iidd) {
00491                   if (DEBUG) cout << "found  (ClusterId == iidd) with ClusterId = " << ClusterId << " and iidd = " << iidd << endl;
00492                   DetId ClusterDetId(ClusterId);
00493                   const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
00494                   
00495                   for(edmNew::DetSet<SiStripCluster>::const_iterator iter=DSViter->begin();iter!=DSViter->end();++iter) {
00496                     //iter is a single SiStripCluster
00497                     StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
00498                     float res = (parameters.first.x() - xloc);
00499                     float sigma = checkConsistency(parameters , xloc, xErr);
00500                     // The consistency is probably more accurately measured with the Chi2MeasurementEstimator. To use it
00501                     // you need a TransientTrackingRecHit instead of the cluster
00502                     //theEstimator=       new Chi2MeasurementEstimator(30);
00503                     //const Chi2MeasurementEstimator *theEstimator(100);
00504                     //theEstimator->estimate(TM->tsos(), TransientTrackingRecHit);
00505                     
00506                     SiStripClusterInfo clusterInfo = SiStripClusterInfo(*iter, es, ClusterId);  
00507                     // signal to noise from SiStripClusterInfo not working in 225. I'll fix this after the interface
00508                     // redesign in 300 -ku
00509                     //float cluster_info[7] = {res, sigma, parameters.first.x(), sqrt(parameters.second.xx()), parameters.first.y(), sqrt(parameters.second.yy()), signal_to_noise};
00510                     std::vector< float > cluster_info;
00511                     cluster_info.push_back(res); 
00512                     cluster_info.push_back(sigma);
00513                     cluster_info.push_back(parameters.first.x()); 
00514                     cluster_info.push_back(sqrt(parameters.second.xx()));
00515                     cluster_info.push_back(parameters.first.y());
00516                     cluster_info.push_back(sqrt(parameters.second.yy()));
00517                     cluster_info.push_back( clusterInfo.signalOverNoise() );
00518                     //cluster_info.push_back( clusterInfo.getSignalOverNoise() );
00519                     VCluster_info.push_back(cluster_info);
00520                     nClusters++;
00521                     if (DEBUG) cout << "Have ID match. residual = " << VCluster_info.back()[0] << "  res sigma = " << VCluster_info.back()[1] << endl;
00522                     if (DEBUG) cout << "trajectory measurement compatability estimate = " << (*itm).estimate() << endl;
00523                     if (DEBUG) cout << "hit position = " << parameters.first.x() << "  hit error = " << sqrt(parameters.second.xx()) << "  trajectory position = " << xloc << "  traj error = " << xErr << endl;
00524                   }
00525                 }
00526               }
00527               float FinalResSig = 1000.0;
00528               float FinalCluster[7]= {1000.0, 1000.0, 0.0, 0.0, 0.0, 0.0, 0.0};
00529               if (nClusters > 0) {
00530                 if (DEBUG) cout << "found clusters > 0" << endl;
00531                 if (nClusters > 1) {
00532                   //get the smallest one
00533                   vector< vector<float> >::iterator ires;
00534                   for (ires=VCluster_info.begin(); ires!=VCluster_info.end(); ires++){
00535                     if ( abs((*ires)[1]) < abs(FinalResSig)) {
00536                       FinalResSig = (*ires)[1];
00537                       for (unsigned int i = 0; i<ires->size(); i++) {
00538                         if (DEBUG) cout << "filling final cluster. i = " << i << " before fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
00539                         FinalCluster[i] = (*ires)[i];
00540                         if (DEBUG) cout << "filling final cluster. i = " << i << " after fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
00541                       }
00542                     }
00543                     if (DEBUG) cout << "iresidual = " << (*ires)[0] << "  isigma = " << (*ires)[1] << "  and FinalRes = " << FinalCluster[0] << endl;
00544                   }
00545                 }
00546                 else {
00547                   FinalResSig = VCluster_info.at(0)[1];
00548                   for (unsigned int i = 0; i<VCluster_info.at(0).size(); i++) {
00549                     FinalCluster[i] = VCluster_info.at(0)[i];
00550                   }
00551                 }
00552                 nClusters=0;
00553                 VCluster_info.clear();
00554               }
00555               
00556               if (DEBUG) cout << "Final residual in X = " << FinalCluster[0] << "+-" << (FinalCluster[0]/FinalResSig) << endl;
00557               if (DEBUG) cout << "Checking location of trajectory: abs(yloc) = " << abs(yloc) << "  abs(xloc) = " << abs(xloc) << endl;
00558               if (DEBUG) cout << "Checking location of cluster hit: yloc = " << FinalCluster[4] << "+-" << FinalCluster[5] << "  xloc = " << FinalCluster[2] << "+-" << FinalCluster[3] << endl;
00559               if (DEBUG) cout << "Final cluster signal to noise = " << FinalCluster[6] << endl;
00560               
00561               float exclusionWidth = 0.4;
00562               float TOBexclusion = 0.0;
00563               float TECexRing5 = -0.89;
00564               float TECexRing6 = -0.56;
00565               float TECexRing7 = 0.60;
00566               //Added by Chris Edelmaier to do TEC bonding exclusion
00567               int subdetector = ((iidd>>25) & 0x7);            
00568               int ringnumber = ((iidd>>5) & 0x7);
00569               
00570               //New TOB and TEC bonding region exclusion zone
00571               if((TKlayers >= 5 && TKlayers < 11)||((subdetector == 6)&&( (ringnumber >= 5)&&(ringnumber <=7) ))) {
00572                 //There are only 2 cases that we need to exclude for
00573                 float highzone = 0.0;
00574                 float lowzone = 0.0;
00575                 float higherr = yloc + 5.0*yErr;
00576                 float lowerr = yloc - 5.0*yErr;
00577                 if(TKlayers >= 5 && TKlayers < 11) {
00578                   //TOB zone
00579                   highzone = TOBexclusion + exclusionWidth;
00580                   lowzone = TOBexclusion - exclusionWidth;
00581                 } else if (ringnumber == 5) {
00582                   //TEC ring 5
00583                   highzone = TECexRing5 + exclusionWidth;
00584                   lowzone = TECexRing5 - exclusionWidth;
00585                 } else if (ringnumber == 6) {
00586                   //TEC ring 6
00587                   highzone = TECexRing6 + exclusionWidth;
00588                   lowzone = TECexRing6 - exclusionWidth;
00589                 } else if (ringnumber == 7) {
00590                   //TEC ring 7
00591                   highzone = TECexRing7 + exclusionWidth;
00592                   lowzone = TECexRing7 - exclusionWidth;
00593                 }
00594                 //Now that we have our exclusion region, we just have to properly identify it
00595                 if((highzone <= higherr)&&(highzone >= lowerr)) withinAcceptance = false;
00596                 if((lowzone >= lowerr)&&(lowzone <= higherr)) withinAcceptance = false;
00597                 if((higherr <= highzone)&&(higherr >= lowzone)) withinAcceptance = false;
00598                 if((lowerr >= lowzone)&&(lowerr <= highzone)) withinAcceptance = false;
00599               }
00600               
00601               // fill ntuple varibles
00602               //get global position from module id number iidd
00603               TrajGlbX = xglob;
00604               TrajGlbY = yglob;
00605               TrajGlbZ = zglob;   
00606               Id = iidd;
00607               run = run_nr;
00608               event = ev_nr;
00609               bunchx = bunch_nr;
00610               //if ( SiStripQuality_->IsModuleBad(iidd) ) {
00611               if ( SiStripQuality_->getBadApvs(iidd)!=0 ) {
00612                 SiStripQualBad = 1; 
00613                 if(DEBUG) cout << "strip is bad from SiStripQuality" << endl;
00614               } else {
00615                 SiStripQualBad = 0; 
00616                 if(DEBUG) cout << "strip is good from SiStripQuality" << endl;
00617               }
00618 
00619               //check for FED-detected errors and include those in SiStripQualBad
00620               for (unsigned int ii=0;ii< fedErrorIds->size();ii++) {
00621                 if (iidd == (*fedErrorIds)[ii].rawId() )
00622                   SiStripQualBad = 1;
00623               }
00624               
00625               TrajLocX = xloc;
00626               TrajLocY = yloc;
00627               TrajLocErrX = xErr;
00628               TrajLocErrY = yErr;
00629               TrajLocAngleX = angleX;
00630               TrajLocAngleY = angleY;
00631               ResX = FinalCluster[0];
00632               ResXSig = FinalResSig;
00633               if (FinalResSig != FinalCluster[1]) if (DEBUG) cout << "Problem with best cluster selection because FinalResSig = " << FinalResSig << " and FinalCluster[1] = " << FinalCluster[1] << endl;
00634               ClusterLocX = FinalCluster[2];
00635               ClusterLocY = FinalCluster[4];
00636               ClusterLocErrX = FinalCluster[3];
00637               ClusterLocErrY = FinalCluster[5];
00638               ClusterStoN = FinalCluster[6];
00639               
00640               if (DEBUG)              cout << "before check good" << endl;
00641               
00642               if ( FinalResSig < 999.0) {  //could make requirement on track/hit consistency, but for
00643                 //now take anything with a hit on the module
00644                 if (DEBUG) cout << "hit being counted as good " << FinalCluster[0] << " FinalRecHit " << 
00645                   iidd << "   TKlayers  "  <<  TKlayers  << " xloc " <<  xloc << " yloc  " << yloc << " module " << iidd << 
00646                   "   matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
00647                 ModIsBad = 0;
00648                 traj->Fill();
00649               }
00650               else {
00651                 if (DEBUG)  cout << "hit being counted as bad   ######### Invalid RPhi FinalResX " << FinalCluster[0] << " FinalRecHit " << 
00652                   iidd << "   TKlayers  "  <<  TKlayers  << " xloc " <<  xloc << " yloc  " << yloc << " module " << iidd << 
00653                   "   matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
00654                 ModIsBad = 1;
00655                 traj->Fill();
00656                 
00657                 if (DEBUG) cout << " RPhi Error " << sqrt(xErr*xErr + yErr*yErr) << " ErrorX " << xErr << " yErr " <<  yErr <<  endl;
00658               } if (DEBUG) cout << "after good location check" << endl;
00659             } if (DEBUG) cout << "after list of clusters" << endl;
00660           } if (DEBUG) cout << "After layers=TKLayers if" << endl;
00661         } if (DEBUG) cout << "After looping over TrajAtValidHit list" << endl;
00662       } if (DEBUG) cout << "end TMeasurement loop" << endl;
00663     } if (DEBUG) cout << "end of trajectories loop" << endl;
00664   }
00665 }
00666 
00667 void HitEff::endJob(){
00668   traj->GetDirectory()->cd();
00669   traj->Write();  
00670   
00671   if (DEBUG) cout << " Events Analysed             " <<  events          << endl;
00672   if (DEBUG) cout << " Number Of Tracked events    " <<  EventTrackCKF   << endl;
00673 }
00674 
00675 double HitEff::checkConsistency(StripClusterParameterEstimator::LocalValues parameters, double xx, double xerr) {
00676   double error = sqrt(parameters.second.xx() + xerr*xerr);
00677   double separation = abs(parameters.first.x() - xx);
00678   double consistency = separation/error;
00679   return consistency;
00680 }
00681 
00682 bool HitEff::isDoubleSided(unsigned int iidd) const {
00683   StripSubdetector strip=StripSubdetector(iidd);
00684   unsigned int subid=strip.subdetId();
00685   unsigned int layer = 0;
00686   if (subid ==  StripSubdetector::TIB) { 
00687     TIBDetId tibid(iidd);
00688     layer = tibid.layer();
00689     if (layer == 1 || layer == 2) return true;
00690     else return false;
00691   }
00692   else if (subid ==  StripSubdetector::TOB) { 
00693     TOBDetId tobid(iidd);
00694     layer = tobid.layer() + 4 ; 
00695     if (layer == 5 || layer == 6) return true;
00696     else return false;
00697   }
00698   else if (subid ==  StripSubdetector::TID) { 
00699     TIDDetId tidid(iidd);
00700     layer = tidid.ring() + 10;
00701     if (layer == 11 || layer == 12) return true;
00702     else return false;
00703   }
00704   else if (subid ==  StripSubdetector::TEC) { 
00705     TECDetId tecid(iidd);
00706     layer = tecid.ring() + 13 ; 
00707     if (layer == 14 || layer == 15 || layer == 18) return true;
00708     else return false;
00709   }
00710   else
00711     return false;
00712 }
00713 
00714 bool HitEff::check2DPartner(unsigned int iidd, std::vector<TrajectoryMeasurement> traj) {
00715   unsigned int partner_iidd = 0;
00716   bool found2DPartner = false;
00717   // first get the id of the other detector
00718   if ((iidd & 0x3)==1) partner_iidd = iidd+1;
00719   if ((iidd & 0x3)==2) partner_iidd = iidd-1;
00720   // next look in the trajectory measurements for a measurement from that detector
00721   // loop through trajectory measurements to find the partner_iidd
00722   for (std::vector<TrajectoryMeasurement>::const_iterator iTM=traj.begin(); iTM!=traj.end(); ++iTM) {
00723     if (iTM->recHit()->geographicalId().rawId()==partner_iidd) {
00724       found2DPartner = true;
00725     }
00726   }
00727   return found2DPartner;
00728 }
00729 
00730 unsigned int HitEff::checkLayer( unsigned int iidd) {
00731   StripSubdetector strip=StripSubdetector(iidd);
00732   unsigned int subid=strip.subdetId();
00733   if (subid ==  StripSubdetector::TIB) { 
00734     TIBDetId tibid(iidd);
00735     return tibid.layer();
00736   }
00737   if (subid ==  StripSubdetector::TOB) { 
00738     TOBDetId tobid(iidd);
00739     return tobid.layer() + 4 ; 
00740   }
00741   if (subid ==  StripSubdetector::TID) { 
00742     TIDDetId tidid(iidd);
00743     return tidid.wheel() + 10;
00744   }
00745   if (subid ==  StripSubdetector::TEC) { 
00746     TECDetId tecid(iidd);
00747     return tecid.wheel() + 13 ; 
00748   }
00749   return 0;
00750 }
00751 
00752 //define this as a plug-in
00753 DEFINE_FWK_MODULE(HitEff);