CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/Alignment/CommonAlignmentProducer/src/AlignmentTrackSelector.cc

Go to the documentation of this file.
00001 #include "Alignment/CommonAlignmentProducer/interface/AlignmentTrackSelector.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 
00006 #include "FWCore/Framework/interface/Event.h"
00007 
00008 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00009 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00011 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00012 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00014 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerMultiRecHit.h"
00015 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00016 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00017 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00018 #include "DataFormats/Alignment/interface/AlignmentClusterFlag.h"
00019 #include "DataFormats/Alignment/interface/AliClusterValueMap.h"
00020 
00021 #include "DataFormats/DetId/interface/DetId.h"
00022 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00023 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00024 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00025 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00026 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00027 const int kBPIX = PixelSubdetector::PixelBarrel;
00028 const int kFPIX = PixelSubdetector::PixelEndcap;
00029 
00030 // constructor ----------------------------------------------------------------
00031 
00032 AlignmentTrackSelector::AlignmentTrackSelector(const edm::ParameterSet & cfg) :
00033   applyBasicCuts_( cfg.getParameter<bool>( "applyBasicCuts" ) ),
00034   applyNHighestPt_( cfg.getParameter<bool>( "applyNHighestPt" ) ),
00035   applyMultiplicityFilter_( cfg.getParameter<bool>( "applyMultiplicityFilter" ) ),
00036   seedOnlyFromAbove_( cfg.getParameter<int>( "seedOnlyFrom" ) ),
00037   applyIsolation_( cfg.getParameter<bool>( "applyIsolationCut" ) ),
00038   chargeCheck_( cfg.getParameter<bool>( "applyChargeCheck" ) ),
00039   nHighestPt_( cfg.getParameter<int>( "nHighestPt" ) ),
00040   minMultiplicity_ ( cfg.getParameter<int>( "minMultiplicity" ) ),
00041   maxMultiplicity_ ( cfg.getParameter<int>( "maxMultiplicity" ) ),
00042   multiplicityOnInput_ ( cfg.getParameter<bool>( "multiplicityOnInput" ) ),
00043   ptMin_( cfg.getParameter<double>( "ptMin" ) ),
00044   ptMax_( cfg.getParameter<double>( "ptMax" ) ),
00045   pMin_( cfg.getParameter<double>( "pMin" ) ),
00046   pMax_( cfg.getParameter<double>( "pMax" ) ),
00047   etaMin_( cfg.getParameter<double>( "etaMin" ) ),
00048   etaMax_( cfg.getParameter<double>( "etaMax" ) ),
00049   phiMin_( cfg.getParameter<double>( "phiMin" ) ),
00050   phiMax_( cfg.getParameter<double>( "phiMax" ) ),
00051   nHitMin_( cfg.getParameter<double>( "nHitMin" ) ),
00052   nHitMax_( cfg.getParameter<double>( "nHitMax" ) ),
00053   chi2nMax_( cfg.getParameter<double>( "chi2nMax" ) ),
00054   d0Min_(  cfg.getParameter<double>( "d0Min" ) ),
00055   d0Max_(  cfg.getParameter<double>( "d0Max" ) ),
00056   dzMin_(  cfg.getParameter<double>( "dzMin" ) ),
00057   dzMax_(  cfg.getParameter<double>( "dzMax" ) ),
00058   theCharge_( cfg.getParameter<int>( "theCharge" ) ),
00059   minHitChargeStrip_( cfg.getParameter<double>( "minHitChargeStrip" ) ),
00060   minHitIsolation_( cfg.getParameter<double>( "minHitIsolation" ) ),
00061   rphirecHitsTag_( cfg.getParameter<edm::InputTag>("rphirecHits") ),
00062   matchedrecHitsTag_( cfg.getParameter<edm::InputTag>("matchedrecHits") ),
00063   countStereoHitAs2D_( cfg.getParameter<bool>( "countStereoHitAs2D" ) ),
00064   nHitMin2D_( cfg.getParameter<unsigned int>( "nHitMin2D" ) ),
00065   // Ugly to use the same getParameter n times, but this allows const cut variables...
00066   minHitsinTIB_(cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTIB" ) ),
00067   minHitsinTOB_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTOB" ) ),
00068   minHitsinTID_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTID" ) ),
00069   minHitsinTEC_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTEC" ) ),
00070   minHitsinBPIX_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inBPIX" ) ),
00071   minHitsinFPIX_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inFPIX" ) ),
00072   minHitsinPIX_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inPIXEL" ) ),
00073   minHitsinTIDplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTIDplus" ) ),
00074   minHitsinTIDminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTIDminus" ) ),
00075   minHitsinTECplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTECplus" ) ),
00076   minHitsinTECminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTECminus" ) ),
00077   minHitsinFPIXplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inFPIXplus" ) ),
00078   minHitsinFPIXminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inFPIXminus" ) ),
00079   minHitsinENDCAP_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inENDCAP" ) ),
00080   minHitsinENDCAPplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inENDCAPplus" ) ),
00081   minHitsinENDCAPminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inENDCAPminus" ) ),
00082   maxHitDiffEndcaps_( cfg.getParameter<double>( "maxHitDiffEndcaps" ) ),
00083   nLostHitMax_( cfg.getParameter<double>( "nLostHitMax" ) ),
00084   clusterValueMapTag_(cfg.getParameter<edm::InputTag>("hitPrescaleMapTag")),
00085   minPrescaledHits_( cfg.getParameter<int>("minPrescaledHits")),
00086   applyPrescaledHitsFilter_(clusterValueMapTag_.encode().size() && minPrescaledHits_ > 0)
00087 {
00088 
00089   //convert track quality from string to enum
00090   std::vector<std::string> trkQualityStrings(cfg.getParameter<std::vector<std::string> >("trackQualities"));
00091   std::string qualities;
00092   if(trkQualityStrings.size()>0){
00093     applyTrkQualityCheck_=true;
00094     for (unsigned int i = 0; i < trkQualityStrings.size(); ++i) {
00095       (qualities += trkQualityStrings[i]) += ", ";
00096       trkQualities_.push_back(reco::TrackBase::qualityByName(trkQualityStrings[i]));
00097     }
00098   }
00099   else applyTrkQualityCheck_=false;
00100 
00101   std::vector<std::string> trkIterStrings(cfg.getParameter<std::vector<std::string> >("iterativeTrackingSteps"));
00102   if(trkIterStrings.size()>0){
00103     applyIterStepCheck_=true;
00104     std::string tracksteps;
00105     for (unsigned int i = 0; i < trkIterStrings.size(); ++i) {
00106       (tracksteps += trkIterStrings[i]) += ", ";
00107       trkSteps_.push_back(reco::TrackBase::algoByName(trkIterStrings[i]));
00108     }
00109   }
00110   else   applyIterStepCheck_=false;
00111 
00112   if (applyBasicCuts_){
00113       edm::LogInfo("AlignmentTrackSelector") 
00114         << "applying basic track cuts ..."
00115         << "\nptmin,ptmax:     " << ptMin_   << "," << ptMax_ 
00116         << "\npmin,pmax:       " << pMin_    << "," << pMax_ 
00117         << "\netamin,etamax:   " << etaMin_  << "," << etaMax_
00118         << "\nphimin,phimax:   " << phiMin_  << "," << phiMax_
00119         << "\nnhitmin,nhitmax: " << nHitMin_ << "," << nHitMax_
00120         << "\nnlosthitmax:     " << nLostHitMax_
00121         << "\nnhitmin2D:       " << nHitMin2D_
00122         << (countStereoHitAs2D_ ? "," : ", not") << " counting hits on SiStrip stereo modules as 2D" 
00123         << "\nchi2nmax:        " << chi2nMax_;
00124 
00125       if (applyIsolation_)
00126         edm::LogInfo("AlignmentTrackSelector") 
00127           << "only retain tracks isolated at least by " << minHitIsolation_ << 
00128           " cm from other rechits"; 
00129       
00130       if (chargeCheck_)
00131         edm::LogInfo("AlignmentTrackSelector") 
00132           << "only retain hits with at least " << minHitChargeStrip_ <<  
00133           " ADC counts of total cluster charge"; 
00134    
00135       edm::LogInfo("AlignmentTrackSelector") 
00136         << "Minimum number of hits in TIB/TID/TOB/TEC/BPIX/FPIX/PIXEL = " 
00137         << minHitsinTIB_ << "/" << minHitsinTID_ << "/" << minHitsinTOB_
00138         << "/" << minHitsinTEC_ << "/" << minHitsinBPIX_ << "/" << minHitsinFPIX_ << "/" << minHitsinPIX_;
00139 
00140       edm::LogInfo("AlignmentTrackSelector")
00141         << "Minimum number of hits in TID+/TID-/TEC+/TEC-/FPIX+/FPIX- = "
00142         << minHitsinTIDplus_ << "/" << minHitsinTIDminus_
00143         << "/" << minHitsinTECplus_ << "/" << minHitsinTECminus_
00144         << "/" << minHitsinFPIXplus_ << "/" << minHitsinFPIXminus_;
00145 
00146       edm::LogInfo("AlignmentTrackSelector")
00147         << "Minimum number of hits in EndCap (TID+TEC)/EndCap+/EndCap- = "
00148         << minHitsinENDCAP_ << "/" << minHitsinENDCAPplus_ << "/" << minHitsinENDCAPminus_;
00149 
00150       edm::LogInfo("AlignmentTrackSelector")
00151         << "Max value of |nHitsinENDCAPplus - nHitsinENDCAPminus| = "
00152         << maxHitDiffEndcaps_;
00153 
00154       if (trkQualityStrings.size()) {
00155         edm::LogInfo("AlignmentTrackSelector")
00156           << "Select tracks with these qualities: " << qualities;
00157       }    
00158   }
00159   
00160   if (applyNHighestPt_)
00161         edm::LogInfo("AlignmentTrackSelector") 
00162           << "filter N tracks with highest Pt N=" << nHighestPt_;
00163 
00164   if (applyMultiplicityFilter_)
00165         edm::LogInfo("AlignmentTrackSelector") 
00166           << "apply multiplicity filter N>= " << minMultiplicity_ << "and N<= " << maxMultiplicity_
00167           << " on " << (multiplicityOnInput_ ? "input" : "output");
00168 
00169   if (applyPrescaledHitsFilter_) {
00170     edm::LogInfo("AlignmentTrackSelector") 
00171       << "apply cut on number of prescaled hits N>= " << minPrescaledHits_
00172       << " (prescale info from " << clusterValueMapTag_ << ")";
00173     
00174   }
00175 
00176 }
00177 
00178 // destructor -----------------------------------------------------------------
00179 
00180 AlignmentTrackSelector::~AlignmentTrackSelector()
00181 {}
00182 
00183 
00184 // do selection ---------------------------------------------------------------
00185 
00186 AlignmentTrackSelector::Tracks 
00187 AlignmentTrackSelector::select(const Tracks& tracks, const edm::Event& evt) const 
00188 {
00189   
00190   if (applyMultiplicityFilter_ && multiplicityOnInput_ && 
00191       (tracks.size() < static_cast<unsigned int>(minMultiplicity_) 
00192        || tracks.size() > static_cast<unsigned int>(maxMultiplicity_))) {
00193     return Tracks(); // empty collection
00194   }
00195 
00196   Tracks result = tracks;
00197   // apply basic track cuts (if selected)
00198   if (applyBasicCuts_) result= this->basicCuts(result, evt);
00199   
00200   // filter N tracks with highest Pt (if selected)
00201   if (applyNHighestPt_) result = this->theNHighestPtTracks(result);
00202   
00203   // apply minimum multiplicity requirement (if selected)
00204   if (applyMultiplicityFilter_ && !multiplicityOnInput_) {
00205     if (result.size() < static_cast<unsigned int>(minMultiplicity_) 
00206         || result.size() > static_cast<unsigned int>(maxMultiplicity_) ) {
00207 
00208       result.clear();
00209     }
00210   }
00211 
00212   if(applyPrescaledHitsFilter_){
00213     result = this->checkPrescaledHits(result, evt);
00214   }
00215 
00216   return result;
00217 }
00218 
00220 bool AlignmentTrackSelector::useThisFilter()
00221 {
00222   return applyMultiplicityFilter_ || applyBasicCuts_ || applyNHighestPt_ || applyPrescaledHitsFilter_;
00223 
00224 }
00225 
00226 
00227 // make basic cuts ------------------------------------------------------------
00228 
00229 AlignmentTrackSelector::Tracks 
00230 AlignmentTrackSelector::basicCuts(const Tracks& tracks, const edm::Event& evt) const 
00231 {
00232   Tracks result;
00233 
00234   for (Tracks::const_iterator it=tracks.begin(); it != tracks.end(); ++it) {
00235     const reco::Track* trackp=*it;
00236     float pt=trackp->pt();
00237     float p=trackp->p();
00238     float eta=trackp->eta();
00239     float phi=trackp->phi();
00240     int nhit = trackp->numberOfValidHits(); 
00241     int nlosthit = trackp->numberOfLostHits(); 
00242     float chi2n = trackp->normalizedChi2();
00243 
00244     int q = trackp->charge();
00245     bool isChargeOk = false;
00246     if(theCharge_==-1 && q<0)         isChargeOk = true; 
00247     else if (theCharge_==1 && q>0)    isChargeOk = true; 
00248     else if (theCharge_==0)           isChargeOk = true; 
00249 
00250     float d0    = trackp->d0();
00251     float dz    = trackp->dz();
00252 
00253     // edm::LogDebug("AlignmentTrackSelector") << " pt,eta,phi,nhit: "
00254     //  <<pt<<","<<eta<<","<<phi<<","<<nhit;
00255 
00256     if (pt>ptMin_ && pt<ptMax_ 
00257        && p>pMin_ && p<pMax_ 
00258        && eta>etaMin_ && eta<etaMax_ 
00259        && phi>phiMin_ && phi<phiMax_ 
00260        && nhit>=nHitMin_ && nhit<=nHitMax_
00261        && nlosthit<=nLostHitMax_
00262        && chi2n<chi2nMax_
00263        && isChargeOk
00264        && d0>=d0Min_ && d0<=d0Max_
00265        && dz>=dzMin_ && dz<=dzMax_) {
00266       bool trkQualityOk=false ;
00267       if (!applyTrkQualityCheck_ &&!applyIterStepCheck_)trkQualityOk=true ; // nothing required
00268       else trkQualityOk = this->isOkTrkQuality(trackp);
00269 
00270       bool hitsCheckOk=this->detailedHitsCheck(trackp, evt);
00271  
00272       if (trkQualityOk && hitsCheckOk )  result.push_back(trackp);
00273     }
00274   }
00275 
00276   return result;
00277 }
00278 
00279 //-----------------------------------------------------------------------------
00280 
00281 bool AlignmentTrackSelector::detailedHitsCheck(const reco::Track *trackp, const edm::Event& evt) const
00282 {
00283   // checking hit requirements beyond simple number of valid hits
00284 
00285   if (minHitsinTIB_ || minHitsinTOB_ || minHitsinTID_ || minHitsinTEC_
00286       || minHitsinENDCAP_ || minHitsinENDCAPplus_ || minHitsinENDCAPminus_ || (maxHitDiffEndcaps_ < 999)
00287       || minHitsinTIDplus_ || minHitsinTIDminus_
00288       || minHitsinFPIXplus_ || minHitsinFPIXminus_
00289       || minHitsinTECplus_ || minHitsinTECminus_
00290       || minHitsinFPIX_ || minHitsinBPIX_ || minHitsinPIX_ ||nHitMin2D_ || chargeCheck_
00291       || applyIsolation_ || (seedOnlyFromAbove_ == 1 || seedOnlyFromAbove_ == 2)) {
00292     // any detailed hit cut is active, so have to check
00293     
00294     int nhitinTIB = 0, nhitinTOB = 0, nhitinTID = 0;
00295     int nhitinTEC = 0, nhitinBPIX = 0, nhitinFPIX = 0, nhitinPIXEL=0;
00296     int nhitinENDCAP = 0, nhitinENDCAPplus = 0, nhitinENDCAPminus = 0;
00297     int nhitinTIDplus = 0, nhitinTIDminus = 0;
00298     int nhitinFPIXplus = 0, nhitinFPIXminus = 0;
00299     int nhitinTECplus = 0, nhitinTECminus = 0;
00300     unsigned int nHit2D = 0;
00301     unsigned int thishit = 0;
00302 
00303     for (trackingRecHit_iterator iHit = trackp->recHitsBegin(); iHit != trackp->recHitsEnd(); ++iHit) {
00304       thishit++;
00305       const DetId detId((*iHit)->geographicalId());
00306       const int subdetId = detId.subdetId(); 
00307 
00308       // *** thishit == 1 means last hit in CTF *** 
00309       // (FIXME: assumption might change or not be valid for all tracking algorthms)
00310       // ==> for cosmics
00311       // seedOnlyFrom = 1 is TIB-TOB-TEC tracks only
00312       // seedOnlyFrom = 2 is TOB-TEC tracks only 
00313       if (seedOnlyFromAbove_ == 1 && thishit == 1 
00314           && (subdetId == int(SiStripDetId::TOB) || subdetId == int(SiStripDetId::TEC))){
00315         return false; 
00316       }
00317       if (seedOnlyFromAbove_ == 2 && thishit == 1 && subdetId == int(SiStripDetId::TIB)) {
00318         return false;
00319       }
00320 
00321       if (!(*iHit)->isValid()) continue; // only real hits count as in trackp->numberOfValidHits()
00322       if (detId.det() != DetId::Tracker) {
00323         edm::LogError("DetectorMismatch") << "@SUB=AlignmentTrackSelector::detailedHitsCheck"
00324                                           << "DetId.det() != DetId::Tracker (=" << DetId::Tracker
00325                                           << "), but " << detId.det() << ".";
00326       }
00327       const TrackingRecHit* therechit = (*iHit).get();
00328       if (chargeCheck_ && !(this->isOkCharge(therechit))) return false;
00329       if (applyIsolation_ && (!this->isIsolated(therechit, evt))) return false;
00330       if      (SiStripDetId::TIB == subdetId) ++nhitinTIB;
00331       else if (SiStripDetId::TOB == subdetId) ++nhitinTOB;
00332       else if (SiStripDetId::TID == subdetId) {
00333         ++nhitinTID;
00334         ++nhitinENDCAP;
00335         TIDDetId tidId(detId);
00336         if (tidId.isZMinusSide()) {
00337           ++nhitinTIDminus;
00338           ++nhitinENDCAPminus;
00339         }
00340         else if (tidId.isZPlusSide()) {
00341           ++nhitinTIDplus;
00342           ++nhitinENDCAPplus;
00343         }
00344       }
00345       else if (SiStripDetId::TEC == subdetId) {
00346         ++nhitinTEC;
00347         ++nhitinENDCAP;
00348         TECDetId tecId(detId);
00349         if (tecId.isZMinusSide()) {
00350           ++nhitinTECminus;
00351           ++nhitinENDCAPminus;
00352         }
00353         else if (tecId.isZPlusSide()) {
00354           ++nhitinTECplus;
00355           ++nhitinENDCAPplus;
00356         }
00357       }
00358       else if (            kBPIX == subdetId) {++nhitinBPIX;++nhitinPIXEL;}
00359       else if (            kFPIX == subdetId) {
00360         ++nhitinFPIX;
00361         ++nhitinPIXEL;
00362         PXFDetId fpixId(detId);
00363         if (fpixId.side()==1) ++nhitinFPIXminus;
00364         else if (fpixId.side()==2) ++nhitinFPIXplus;
00365       }
00366       // Do not call isHit2D(..) if already enough 2D hits for performance reason:
00367       if (nHit2D < nHitMin2D_ && this->isHit2D(**iHit)) ++nHit2D;
00368     } // end loop on hits
00369   
00370     return (nhitinTIB >= minHitsinTIB_ && nhitinTOB >= minHitsinTOB_ 
00371             && nhitinTID >= minHitsinTID_ && nhitinTEC >= minHitsinTEC_ 
00372             && nhitinENDCAP >= minHitsinENDCAP_ && nhitinENDCAPplus >= minHitsinENDCAPplus_ && nhitinENDCAPminus >= minHitsinENDCAPminus_
00373             && std::abs(nhitinENDCAPplus-nhitinENDCAPminus) <= maxHitDiffEndcaps_
00374             && nhitinTIDplus >= minHitsinTIDplus_ && nhitinTIDminus >= minHitsinTIDminus_
00375             && nhitinFPIXplus >= minHitsinFPIXplus_ && nhitinFPIXminus >= minHitsinFPIXminus_
00376             && nhitinTECplus >= minHitsinTECplus_ && nhitinTECminus >= minHitsinTECminus_
00377             && nhitinBPIX >= minHitsinBPIX_ 
00378             && nhitinFPIX >= minHitsinFPIX_ && nhitinPIXEL>=minHitsinPIX_ 
00379             && nHit2D >= nHitMin2D_);
00380   } else { // no cuts set, so we are just fine and can avoid loop on hits
00381     return true;
00382   }
00383   
00384 }
00385 
00386 //-----------------------------------------------------------------------------
00387 
00388 bool AlignmentTrackSelector::isHit2D(const TrackingRecHit &hit) const
00389 {
00390   // we count SiStrip stereo modules as 2D if selected via countStereoHitAs2D_
00391   // (since they provide theta information)
00392   if (!hit.isValid() ||
00393       (hit.dimension() < 2 && !countStereoHitAs2D_ && !dynamic_cast<const SiStripRecHit1D*>(&hit))){
00394     return false; // real RecHit1D - but SiStripRecHit1D depends on countStereoHitAs2D_
00395   } else {
00396     const DetId detId(hit.geographicalId());
00397     if (detId.det() == DetId::Tracker) {
00398       if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
00399         return true; // pixel is always 2D
00400       } else { // should be SiStrip now
00401         const SiStripDetId stripId(detId);
00402         if (stripId.stereo()) return countStereoHitAs2D_; // stereo modules
00403         else if (dynamic_cast<const SiStripRecHit1D*>(&hit)
00404                  || dynamic_cast<const SiStripRecHit2D*>(&hit)) return false; // rphi modules hit
00405         //the following two are not used any more since ages... 
00406         else if (dynamic_cast<const SiStripMatchedRecHit2D*>(&hit)) return true; // matched is 2D
00407         else if (dynamic_cast<const ProjectedSiStripRecHit2D*>(&hit)) {
00408           const ProjectedSiStripRecHit2D* pH = static_cast<const ProjectedSiStripRecHit2D*>(&hit);
00409           return (countStereoHitAs2D_ && this->isHit2D(pH->originalHit())); // depends on original...
00410         } else {
00411           edm::LogError("UnkownType") << "@SUB=AlignmentTrackSelector::isHit2D"
00412                                       << "Tracker hit not in pixel, neither SiStripRecHit[12]D nor "
00413                                       << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
00414           return false;
00415         }
00416       }
00417     } else { // not tracker??
00418       edm::LogWarning("DetectorMismatch") << "@SUB=AlignmentTrackSelector::isHit2D"
00419                                           << "Hit not in tracker with 'official' dimension >=2.";
00420       return true; // dimension() >= 2 so accept that...
00421     }
00422   }
00423   // never reached...
00424 }
00425 
00426 //-----------------------------------------------------------------------------
00427 
00428 bool AlignmentTrackSelector::isOkCharge(const TrackingRecHit* hit) const
00429 {
00430   if (!hit || !hit->isValid()) return true; 
00431 
00432   // check det and subdet
00433   const DetId id(hit->geographicalId());
00434   if (id.det() != DetId::Tracker) {
00435     edm::LogWarning("DetectorMismatch") << "@SUB=isOkCharge" << "Hit not in tracker!";
00436     return true;
00437   }
00438   if (id.subdetId() == kFPIX || id.subdetId() == kBPIX) {
00439     return true; // might add some requirement...
00440   }
00441 
00442   // We are in SiStrip now, so test normal hit:
00443   const std::type_info &type = typeid(*hit);
00444   
00445 
00446   if (type == typeid(SiStripRecHit2D)) {
00447     const SiStripRecHit2D *stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
00448     if (stripHit2D) {
00449       return this->isOkChargeStripHit(*stripHit2D);
00450     }
00451   }
00452   else if(type == typeid(SiStripRecHit1D)){
00453     const SiStripRecHit1D *stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
00454     if (stripHit1D) {
00455       return this->isOkChargeStripHit(*stripHit1D);
00456     } 
00457   }
00458 
00459   else if(type == typeid(SiStripMatchedRecHit2D)){  // or matched (should not occur anymore due to hit splitting since 20X)
00460     const SiStripMatchedRecHit2D *matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(hit);
00461     if (matchedHit) {  
00462       return (this->isOkChargeStripHit(matchedHit->monoHit())
00463               && this->isOkChargeStripHit(matchedHit->stereoHit()));
00464     }
00465   }
00466   else if(type == typeid(ProjectedSiStripRecHit2D)){ 
00467     // or projected (should not occur anymore due to hit splitting since 20X):
00468     const ProjectedSiStripRecHit2D *projHit = dynamic_cast<const ProjectedSiStripRecHit2D*>(hit);
00469     if (projHit) {
00470       return this->isOkChargeStripHit(projHit->originalHit());
00471     }
00472   }
00473   else{
00474     edm::LogError("AlignmentTrackSelector")<< "@SUB=isOkCharge" <<"Unknown type of a valid tracker hit in Strips "
00475                                            << " SubDet = "<<id.subdetId();
00476     return false;
00477   }
00478   
00479 
00480 
00481   // and now? SiTrackerMultiRecHit? Not here I guess!
00482   // Should we throw instead?
00483   edm::LogError("AlignmentTrackSelector") 
00484     << "@SUB=isOkCharge" << "Unknown valid tracker hit not in pixel, subdet " << id.subdetId()
00485     << ", SiTrackerMultiRecHit " << dynamic_cast<const SiTrackerMultiRecHit*>(hit)
00486     << ", BaseTrackerRecHit " << dynamic_cast<const BaseTrackerRecHit*>(hit); 
00487   
00488   return true;
00489 } 
00490 
00491 
00492 //-----------------------------------------------------------------------------
00493 
00494 bool AlignmentTrackSelector::isOkChargeStripHit(const SiStripRecHit2D & siStripRecHit2D) const
00495 {
00496   double charge = 0.;
00497 
00498   SiStripRecHit2D::ClusterRef cluster(siStripRecHit2D.cluster());
00499   const std::vector<uint8_t> &amplitudes = cluster->amplitudes();
00500 
00501   for (size_t ia = 0; ia < amplitudes.size(); ++ia) {
00502     charge += amplitudes[ia];
00503   }
00504 
00505   return (charge >= minHitChargeStrip_);
00506 }
00507 
00508 //-----------------------------------------------------------------------------
00509 
00510 bool AlignmentTrackSelector::isOkChargeStripHit(const SiStripRecHit1D & siStripRecHit1D) const
00511 {
00512   double charge = 0.;
00513 
00514   SiStripRecHit1D::ClusterRef cluster(siStripRecHit1D.cluster());
00515   const std::vector<uint8_t> &amplitudes = cluster->amplitudes();
00516 
00517   for (size_t ia = 0; ia < amplitudes.size(); ++ia) {
00518     charge += amplitudes[ia];
00519   }
00520 
00521   return (charge >= minHitChargeStrip_);
00522 }
00523 
00524 //-----------------------------------------------------------------------------
00525 
00526 bool AlignmentTrackSelector::isIsolated(const TrackingRecHit* therechit, const edm::Event& evt) const
00527 {
00528   // FIXME:
00529   // adapt to changes after introduction of SiStripRecHit1D...
00530   //
00531   // edm::ESHandle<TrackerGeometry> tracker;
00532   edm::Handle<SiStripRecHit2DCollection> rphirecHits;
00533   edm::Handle<SiStripMatchedRecHit2DCollection> matchedrecHits;
00534   // es.get<TrackerDigiGeometryRecord>().get(tracker);
00535   evt.getByLabel( rphirecHitsTag_, rphirecHits );
00536   evt.getByLabel( matchedrecHitsTag_, matchedrecHits ); 
00537 
00538   SiStripRecHit2DCollection::DataContainer::const_iterator istripSt; 
00539   SiStripMatchedRecHit2DCollection::DataContainer::const_iterator istripStm; 
00540   const SiStripRecHit2DCollection::DataContainer& stripcollSt = rphirecHits->data();
00541   const SiStripMatchedRecHit2DCollection::DataContainer& stripcollStm = matchedrecHits->data();
00542   
00543   DetId idet = therechit->geographicalId(); 
00544   
00545   // FIXME: instead of looping the full hit collection, we should explore the features of 
00546   // SiStripRecHit2DCollection::rangeRphi = rphirecHits.get(idet) and loop
00547   // only from rangeRphi.first until rangeRphi.second
00548   for( istripSt=stripcollSt.begin(); istripSt!=stripcollSt.end(); ++istripSt ) {
00549     const SiStripRecHit2D *aHit = &*(istripSt);
00550     DetId mydet1 = aHit->geographicalId(); 
00551     if (idet.rawId() != mydet1.rawId()) continue; 
00552     float theDistance = ( therechit->localPosition() - aHit->localPosition() ).mag();
00553     // std::cout << "theDistance1 = " << theDistance << "\n";
00554     if (theDistance > 0.001 && theDistance < minHitIsolation_) return false;
00555   }
00556   
00557   // FIXME: see above
00558   for( istripStm=stripcollStm.begin(); istripStm!=stripcollStm.end(); ++istripStm ) {
00559     const SiStripMatchedRecHit2D *aHit = &*(istripStm);
00560     DetId mydet2 = aHit->geographicalId(); 
00561     if (idet.rawId() != mydet2.rawId()) continue;
00562     float theDistance = (therechit->localPosition() - aHit->localPosition()).mag();
00563     // std::cout << "theDistance1 = " << theDistance << "\n";
00564     if (theDistance > 0.001 && theDistance < minHitIsolation_) return false;
00565   }
00566   return true;
00567 }
00568 
00569 //-----------------------------------------------------------------------------
00570 
00571 AlignmentTrackSelector::Tracks 
00572 AlignmentTrackSelector::theNHighestPtTracks(const Tracks& tracks) const
00573 {
00574   Tracks sortedTracks=tracks;
00575   Tracks result;
00576 
00577   // sort in pt
00578   std::sort(sortedTracks.begin(),sortedTracks.end(),ptComparator);
00579 
00580   // copy theTrackMult highest pt tracks to result vector
00581   int n=0;
00582   for (Tracks::const_iterator it=sortedTracks.begin();
00583            it!=sortedTracks.end(); ++it) {
00584         if (n<nHighestPt_) { result.push_back(*it); n++; }
00585   }
00586 
00587   return result;
00588 }
00589 
00590 //---------
00591 
00592 AlignmentTrackSelector::Tracks 
00593 AlignmentTrackSelector::checkPrescaledHits(const Tracks& tracks, const edm::Event& evt) const 
00594 {
00595   Tracks result;
00596 
00597   //take Cluster-Flag Assomap
00598   edm::Handle<AliClusterValueMap> fMap;
00599   evt.getByLabel( clusterValueMapTag_, fMap);
00600   const AliClusterValueMap &flagMap=*fMap;
00601 
00602   //for each track loop on hits and count the number of taken hits
00603   for (Tracks::const_iterator ittrk=tracks.begin(); ittrk != tracks.end(); ++ittrk) {
00604     const reco::Track* trackp=*ittrk;
00605     int ntakenhits=0;
00606     //    float pt=trackp->pt();
00607 
00608     for (trackingRecHit_iterator ith = trackp->recHitsBegin(), edh = trackp->recHitsEnd(); ith != edh; ++ith) {
00609       const TrackingRecHit *hit = ith->get(); // ith is an iterator on edm::Ref to rechit
00610       if(! hit->isValid())continue;
00611       DetId detid = hit->geographicalId();
00612       int subDet = detid.subdetId();
00613       AlignmentClusterFlag flag;
00614 
00615       bool isPixelHit=(subDet == kFPIX || subDet == kBPIX);
00616 
00617       if (!isPixelHit){
00618         const std::type_info &type = typeid(*hit);
00619 
00620         if (type == typeid(SiStripRecHit2D)) {  
00621           const SiStripRecHit2D* striphit=dynamic_cast<const  SiStripRecHit2D*>(hit); 
00622           if(striphit!=0){
00623             SiStripRecHit2D::ClusterRef stripclust(striphit->cluster());
00624             flag = flagMap[stripclust];
00625           }
00626         }
00627         else if(type == typeid(SiStripRecHit1D)){
00628           const SiStripRecHit1D* striphit = dynamic_cast<const SiStripRecHit1D*>(hit);
00629           if(striphit!=0){
00630             SiStripRecHit1D::ClusterRef stripclust(striphit->cluster());
00631             flag = flagMap[stripclust];
00632           }
00633         } 
00634         else{
00635            edm::LogError("AlignmentTrackSelector")<<"ERROR in <AlignmentTrackSelector::checkPrescaledHits>: Dynamic cast of Strip RecHit failed!"
00636                                                   <<" Skipping this hit.";
00637            continue;
00638         }
00639 
00640           
00641 
00642       }//end if hit in Strips
00643       else{ // test explicitely BPIX/FPIX
00644         const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
00645         if(pixelhit!=0){
00646           SiPixelRecHit::ClusterRef pixclust(pixelhit->cluster());
00647           flag = flagMap[pixclust];
00648         }
00649         else{
00650             edm::LogError("AlignmentTrackSelector")<<"ERROR in <AlignmentTrackSelector::checkPrescaledHits>: Dynamic cast of Pixel RecHit failed!  ";
00651         }
00652       }//end else hit is in Pixel
00653       
00654       if(flag.isTaken())ntakenhits++;
00655 
00656     }//end loop on hits
00657     if(ntakenhits >= minPrescaledHits_)result.push_back(trackp);
00658   }//end loop on tracks
00659 
00660   return result;
00661 }//end checkPrescaledHits
00662 
00663 //-----------------------------------------------------------------
00664 
00665 bool AlignmentTrackSelector::isOkTrkQuality(const reco::Track* track) const
00666 {
00667   bool qualityOk=false;
00668   bool iterStepOk=false;
00669 
00670   //check iterative step
00671   if(applyIterStepCheck_){
00672     for (unsigned int i = 0; i < trkSteps_.size(); ++i) {
00673       if (track->algo()==(trkSteps_[i])) {
00674         iterStepOk=true;
00675       }
00676     }
00677   }
00678   else iterStepOk=true;
00679 
00680   //check track quality
00681   if(applyTrkQualityCheck_){
00682     for (unsigned int i = 0; i < trkQualities_.size(); ++i) {
00683       if (track->quality(trkQualities_[i])) {
00684         qualityOk=true;
00685       }
00686     }
00687   }
00688   else  qualityOk=true;
00689 
00690   return qualityOk&&iterStepOk;
00691 }//end check on track quality