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
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 minHitChargeStrip_( cfg.getParameter<double>( "minHitChargeStrip" ) ),
00059 minHitIsolation_( cfg.getParameter<double>( "minHitIsolation" ) ),
00060 rphirecHitsTag_( cfg.getParameter<edm::InputTag>("rphirecHits") ),
00061 matchedrecHitsTag_( cfg.getParameter<edm::InputTag>("matchedrecHits") ),
00062 countStereoHitAs2D_( cfg.getParameter<bool>( "countStereoHitAs2D" ) ),
00063 nHitMin2D_( cfg.getParameter<unsigned int>( "nHitMin2D" ) ),
00064
00065 minHitsinTIB_(cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTIB" ) ),
00066 minHitsinTOB_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTOB" ) ),
00067 minHitsinTID_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTID" ) ),
00068 minHitsinTEC_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTEC" ) ),
00069 minHitsinBPIX_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inBPIX" ) ),
00070 minHitsinFPIX_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inFPIX" ) ),
00071 minHitsinPIX_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inPIXEL" ) ),
00072 minHitsinTIDplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTIDplus" ) ),
00073 minHitsinTIDminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTIDminus" ) ),
00074 minHitsinTECplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTECplus" ) ),
00075 minHitsinTECminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inTECminus" ) ),
00076 minHitsinFPIXplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inFPIXplus" ) ),
00077 minHitsinFPIXminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inFPIXminus" ) ),
00078 minHitsinENDCAP_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inENDCAP" ) ),
00079 minHitsinENDCAPplus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inENDCAPplus" ) ),
00080 minHitsinENDCAPminus_ (cfg.getParameter<edm::ParameterSet>( "minHitsPerSubDet" ).getParameter<int>( "inENDCAPminus" ) ),
00081 maxHitDiffEndcaps_( cfg.getParameter<double>( "maxHitDiffEndcaps" ) ),
00082 nLostHitMax_( cfg.getParameter<double>( "nLostHitMax" ) ),
00083 clusterValueMapTag_(cfg.getParameter<edm::InputTag>("hitPrescaleMapTag")),
00084 minPrescaledHits_( cfg.getParameter<int>("minPrescaledHits")),
00085 applyPrescaledHitsFilter_(clusterValueMapTag_.encode().size() && minPrescaledHits_ > 0)
00086 {
00087
00088
00089 std::vector<std::string> trkQualityStrings(cfg.getParameter<std::vector<std::string> >("trackQualities"));
00090 std::string qualities;
00091 if(trkQualityStrings.size()>0){
00092 applyTrkQualityCheck_=true;
00093 for (unsigned int i = 0; i < trkQualityStrings.size(); ++i) {
00094 (qualities += trkQualityStrings[i]) += ", ";
00095 trkQualities_.push_back(reco::TrackBase::qualityByName(trkQualityStrings[i]));
00096 }
00097 }
00098 else applyTrkQualityCheck_=false;
00099
00100 std::vector<std::string> trkIterStrings(cfg.getParameter<std::vector<std::string> >("iterativeTrackingSteps"));
00101 if(trkIterStrings.size()>0){
00102 applyIterStepCheck_=true;
00103 std::string tracksteps;
00104 for (unsigned int i = 0; i < trkIterStrings.size(); ++i) {
00105 (tracksteps += trkIterStrings[i]) += ", ";
00106 trkSteps_.push_back(reco::TrackBase::algoByName(trkIterStrings[i]));
00107 }
00108 }
00109 else applyIterStepCheck_=false;
00110
00111 if (applyBasicCuts_){
00112 edm::LogInfo("AlignmentTrackSelector")
00113 << "applying basic track cuts ..."
00114 << "\nptmin,ptmax: " << ptMin_ << "," << ptMax_
00115 << "\npmin,pmax: " << pMin_ << "," << pMax_
00116 << "\netamin,etamax: " << etaMin_ << "," << etaMax_
00117 << "\nphimin,phimax: " << phiMin_ << "," << phiMax_
00118 << "\nnhitmin,nhitmax: " << nHitMin_ << "," << nHitMax_
00119 << "\nnlosthitmax: " << nLostHitMax_
00120 << "\nnhitmin2D: " << nHitMin2D_
00121 << (countStereoHitAs2D_ ? "," : ", not") << " counting hits on SiStrip stereo modules as 2D"
00122 << "\nchi2nmax: " << chi2nMax_;
00123
00124 if (applyIsolation_)
00125 edm::LogInfo("AlignmentTrackSelector")
00126 << "only retain tracks isolated at least by " << minHitIsolation_ <<
00127 " cm from other rechits";
00128
00129 if (chargeCheck_)
00130 edm::LogInfo("AlignmentTrackSelector")
00131 << "only retain hits with at least " << minHitChargeStrip_ <<
00132 " ADC counts of total cluster charge";
00133
00134 edm::LogInfo("AlignmentTrackSelector")
00135 << "Minimum number of hits in TIB/TID/TOB/TEC/BPIX/FPIX/PIXEL = "
00136 << minHitsinTIB_ << "/" << minHitsinTID_ << "/" << minHitsinTOB_
00137 << "/" << minHitsinTEC_ << "/" << minHitsinBPIX_ << "/" << minHitsinFPIX_ << "/" << minHitsinPIX_;
00138
00139 edm::LogInfo("AlignmentTrackSelector")
00140 << "Minimum number of hits in TID+/TID-/TEC+/TEC-/FPIX+/FPIX- = "
00141 << minHitsinTIDplus_ << "/" << minHitsinTIDminus_
00142 << "/" << minHitsinTECplus_ << "/" << minHitsinTECminus_
00143 << "/" << minHitsinFPIXplus_ << "/" << minHitsinFPIXminus_;
00144
00145 edm::LogInfo("AlignmentTrackSelector")
00146 << "Minimum number of hits in EndCap (TID+TEC)/EndCap+/EndCap- = "
00147 << minHitsinENDCAP_ << "/" << minHitsinENDCAPplus_ << "/" << minHitsinENDCAPminus_;
00148
00149 edm::LogInfo("AlignmentTrackSelector")
00150 << "Max value of |nHitsinENDCAPplus - nHitsinENDCAPminus| = "
00151 << maxHitDiffEndcaps_;
00152
00153 if (trkQualityStrings.size()) {
00154 edm::LogInfo("AlignmentTrackSelector")
00155 << "Select tracks with these qualities: " << qualities;
00156 }
00157 }
00158
00159 if (applyNHighestPt_)
00160 edm::LogInfo("AlignmentTrackSelector")
00161 << "filter N tracks with highest Pt N=" << nHighestPt_;
00162
00163 if (applyMultiplicityFilter_)
00164 edm::LogInfo("AlignmentTrackSelector")
00165 << "apply multiplicity filter N>= " << minMultiplicity_ << "and N<= " << maxMultiplicity_
00166 << " on " << (multiplicityOnInput_ ? "input" : "output");
00167
00168 if (applyPrescaledHitsFilter_) {
00169 edm::LogInfo("AlignmentTrackSelector")
00170 << "apply cut on number of prescaled hits N>= " << minPrescaledHits_
00171 << " (prescale info from " << clusterValueMapTag_ << ")";
00172
00173 }
00174
00175 }
00176
00177
00178
00179 AlignmentTrackSelector::~AlignmentTrackSelector()
00180 {}
00181
00182
00183
00184
00185 AlignmentTrackSelector::Tracks
00186 AlignmentTrackSelector::select(const Tracks& tracks, const edm::Event& evt) const
00187 {
00188
00189 if (applyMultiplicityFilter_ && multiplicityOnInput_ &&
00190 (tracks.size() < static_cast<unsigned int>(minMultiplicity_)
00191 || tracks.size() > static_cast<unsigned int>(maxMultiplicity_))) {
00192 return Tracks();
00193 }
00194
00195 Tracks result = tracks;
00196
00197 if (applyBasicCuts_) result= this->basicCuts(result, evt);
00198
00199
00200 if (applyNHighestPt_) result = this->theNHighestPtTracks(result);
00201
00202
00203 if (applyMultiplicityFilter_ && !multiplicityOnInput_) {
00204 if (result.size() < static_cast<unsigned int>(minMultiplicity_)
00205 || result.size() > static_cast<unsigned int>(maxMultiplicity_) ) {
00206
00207 result.clear();
00208 }
00209 }
00210
00211 if(applyPrescaledHitsFilter_){
00212 result = this->checkPrescaledHits(result, evt);
00213 }
00214
00215 return result;
00216 }
00217
00219 bool AlignmentTrackSelector::useThisFilter()
00220 {
00221 return applyMultiplicityFilter_ || applyBasicCuts_ || applyNHighestPt_ || applyPrescaledHitsFilter_;
00222
00223 }
00224
00225
00226
00227
00228 AlignmentTrackSelector::Tracks
00229 AlignmentTrackSelector::basicCuts(const Tracks& tracks, const edm::Event& evt) const
00230 {
00231 Tracks result;
00232
00233 for (Tracks::const_iterator it=tracks.begin(); it != tracks.end(); ++it) {
00234 const reco::Track* trackp=*it;
00235 float pt=trackp->pt();
00236 float p=trackp->p();
00237 float eta=trackp->eta();
00238 float phi=trackp->phi();
00239 int nhit = trackp->numberOfValidHits();
00240 int nlosthit = trackp->numberOfLostHits();
00241 float chi2n = trackp->normalizedChi2();
00242 float d0 = trackp->d0();
00243 float dz = trackp->dz();
00244
00245
00246
00247
00248 if (pt>ptMin_ && pt<ptMax_
00249 && p>pMin_ && p<pMax_
00250 && eta>etaMin_ && eta<etaMax_
00251 && phi>phiMin_ && phi<phiMax_
00252 && nhit>=nHitMin_ && nhit<=nHitMax_
00253 && nlosthit<=nLostHitMax_
00254 && chi2n<chi2nMax_
00255 && d0>=d0Min_ && d0<=d0Max_
00256 && dz>=dzMin_ && dz<=dzMax_) {
00257 bool trkQualityOk=false ;
00258 if (!applyTrkQualityCheck_ &&!applyIterStepCheck_)trkQualityOk=true ;
00259 else trkQualityOk = this->isOkTrkQuality(trackp);
00260
00261 bool hitsCheckOk=this->detailedHitsCheck(trackp, evt);
00262
00263 if (trkQualityOk && hitsCheckOk ) result.push_back(trackp);
00264 }
00265 }
00266
00267 return result;
00268 }
00269
00270
00271
00272 bool AlignmentTrackSelector::detailedHitsCheck(const reco::Track *trackp, const edm::Event& evt) const
00273 {
00274
00275
00276 if (minHitsinTIB_ || minHitsinTOB_ || minHitsinTID_ || minHitsinTEC_
00277 || minHitsinENDCAP_ || minHitsinENDCAPplus_ || minHitsinENDCAPminus_ || (maxHitDiffEndcaps_ < 999)
00278 || minHitsinTIDplus_ || minHitsinTIDminus_
00279 || minHitsinFPIXplus_ || minHitsinFPIXminus_
00280 || minHitsinTECplus_ || minHitsinTECminus_
00281 || minHitsinFPIX_ || minHitsinBPIX_ || minHitsinPIX_ ||nHitMin2D_ || chargeCheck_
00282 || applyIsolation_ || (seedOnlyFromAbove_ == 1 || seedOnlyFromAbove_ == 2)) {
00283
00284
00285 int nhitinTIB = 0, nhitinTOB = 0, nhitinTID = 0;
00286 int nhitinTEC = 0, nhitinBPIX = 0, nhitinFPIX = 0, nhitinPIXEL=0;
00287 int nhitinENDCAP = 0, nhitinENDCAPplus = 0, nhitinENDCAPminus = 0;
00288 int nhitinTIDplus = 0, nhitinTIDminus = 0;
00289 int nhitinFPIXplus = 0, nhitinFPIXminus = 0;
00290 int nhitinTECplus = 0, nhitinTECminus = 0;
00291 unsigned int nHit2D = 0;
00292 unsigned int thishit = 0;
00293
00294 for (trackingRecHit_iterator iHit = trackp->recHitsBegin(); iHit != trackp->recHitsEnd(); ++iHit) {
00295 thishit++;
00296 const DetId detId((*iHit)->geographicalId());
00297 const int subdetId = detId.subdetId();
00298
00299
00300
00301
00302
00303
00304 if (seedOnlyFromAbove_ == 1 && thishit == 1
00305 && (subdetId == int(SiStripDetId::TOB) || subdetId == int(SiStripDetId::TEC))){
00306 return false;
00307 }
00308 if (seedOnlyFromAbove_ == 2 && thishit == 1 && subdetId == int(SiStripDetId::TIB)) {
00309 return false;
00310 }
00311
00312 if (!(*iHit)->isValid()) continue;
00313 if (detId.det() != DetId::Tracker) {
00314 edm::LogError("DetectorMismatch") << "@SUB=AlignmentTrackSelector::detailedHitsCheck"
00315 << "DetId.det() != DetId::Tracker (=" << DetId::Tracker
00316 << "), but " << detId.det() << ".";
00317 }
00318 const TrackingRecHit* therechit = (*iHit).get();
00319 if (chargeCheck_ && !(this->isOkCharge(therechit))) return false;
00320 if (applyIsolation_ && (!this->isIsolated(therechit, evt))) return false;
00321 if (SiStripDetId::TIB == subdetId) ++nhitinTIB;
00322 else if (SiStripDetId::TOB == subdetId) ++nhitinTOB;
00323 else if (SiStripDetId::TID == subdetId) {
00324 ++nhitinTID;
00325 ++nhitinENDCAP;
00326 TIDDetId tidId(detId);
00327 if (tidId.isZMinusSide()) {
00328 ++nhitinTIDminus;
00329 ++nhitinENDCAPminus;
00330 }
00331 else if (tidId.isZPlusSide()) {
00332 ++nhitinTIDplus;
00333 ++nhitinENDCAPplus;
00334 }
00335 }
00336 else if (SiStripDetId::TEC == subdetId) {
00337 ++nhitinTEC;
00338 ++nhitinENDCAP;
00339 TECDetId tecId(detId);
00340 if (tecId.isZMinusSide()) {
00341 ++nhitinTECminus;
00342 ++nhitinENDCAPminus;
00343 }
00344 else if (tecId.isZPlusSide()) {
00345 ++nhitinTECplus;
00346 ++nhitinENDCAPplus;
00347 }
00348 }
00349 else if ( kBPIX == subdetId) {++nhitinBPIX;++nhitinPIXEL;}
00350 else if ( kFPIX == subdetId) {
00351 ++nhitinFPIX;
00352 ++nhitinPIXEL;
00353 PXFDetId fpixId(detId);
00354 if (fpixId.side()==1) ++nhitinFPIXminus;
00355 else if (fpixId.side()==2) ++nhitinFPIXplus;
00356 }
00357
00358 if (nHit2D < nHitMin2D_ && this->isHit2D(**iHit)) ++nHit2D;
00359 }
00360
00361 return (nhitinTIB >= minHitsinTIB_ && nhitinTOB >= minHitsinTOB_
00362 && nhitinTID >= minHitsinTID_ && nhitinTEC >= minHitsinTEC_
00363 && nhitinENDCAP >= minHitsinENDCAP_ && nhitinENDCAPplus >= minHitsinENDCAPplus_ && nhitinENDCAPminus >= minHitsinENDCAPminus_
00364 && std::abs(nhitinENDCAPplus-nhitinENDCAPminus) <= maxHitDiffEndcaps_
00365 && nhitinTIDplus >= minHitsinTIDplus_ && nhitinTIDminus >= minHitsinTIDminus_
00366 && nhitinFPIXplus >= minHitsinFPIXplus_ && nhitinFPIXminus >= minHitsinFPIXminus_
00367 && nhitinTECplus >= minHitsinTECplus_ && nhitinTECminus >= minHitsinTECminus_
00368 && nhitinBPIX >= minHitsinBPIX_
00369 && nhitinFPIX >= minHitsinFPIX_ && nhitinPIXEL>=minHitsinPIX_
00370 && nHit2D >= nHitMin2D_);
00371 } else {
00372 return true;
00373 }
00374
00375 }
00376
00377
00378
00379 bool AlignmentTrackSelector::isHit2D(const TrackingRecHit &hit) const
00380 {
00381
00382
00383 if (!hit.isValid() ||
00384 (hit.dimension() < 2 && !countStereoHitAs2D_ && !dynamic_cast<const SiStripRecHit1D*>(&hit))){
00385 return false;
00386 } else {
00387 const DetId detId(hit.geographicalId());
00388 if (detId.det() == DetId::Tracker) {
00389 if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
00390 return true;
00391 } else {
00392 const SiStripDetId stripId(detId);
00393 if (stripId.stereo()) return countStereoHitAs2D_;
00394 else if (dynamic_cast<const SiStripRecHit1D*>(&hit)
00395 || dynamic_cast<const SiStripRecHit2D*>(&hit)) return false;
00396
00397 else if (dynamic_cast<const SiStripMatchedRecHit2D*>(&hit)) return true;
00398 else if (dynamic_cast<const ProjectedSiStripRecHit2D*>(&hit)) {
00399 const ProjectedSiStripRecHit2D* pH = static_cast<const ProjectedSiStripRecHit2D*>(&hit);
00400 return (countStereoHitAs2D_ && this->isHit2D(pH->originalHit()));
00401 } else {
00402 edm::LogError("UnkownType") << "@SUB=AlignmentTrackSelector::isHit2D"
00403 << "Tracker hit not in pixel, neither SiStripRecHit[12]D nor "
00404 << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
00405 return false;
00406 }
00407 }
00408 } else {
00409 edm::LogWarning("DetectorMismatch") << "@SUB=AlignmentTrackSelector::isHit2D"
00410 << "Hit not in tracker with 'official' dimension >=2.";
00411 return true;
00412 }
00413 }
00414
00415 }
00416
00417
00418
00419 bool AlignmentTrackSelector::isOkCharge(const TrackingRecHit* hit) const
00420 {
00421 if (!hit || !hit->isValid()) return true;
00422
00423
00424 const DetId id(hit->geographicalId());
00425 if (id.det() != DetId::Tracker) {
00426 edm::LogWarning("DetectorMismatch") << "@SUB=isOkCharge" << "Hit not in tracker!";
00427 return true;
00428 }
00429 if (id.subdetId() == kFPIX || id.subdetId() == kBPIX) {
00430 return true;
00431 }
00432
00433
00434 const std::type_info &type = typeid(*hit);
00435
00436
00437 if (type == typeid(SiStripRecHit2D)) {
00438 const SiStripRecHit2D *stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
00439 if (stripHit2D) {
00440 return this->isOkChargeStripHit(stripHit2D);
00441 }
00442 }
00443 else if(type == typeid(SiStripRecHit1D)){
00444 const SiStripRecHit1D *stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
00445 if (stripHit1D) {
00446 return this->isOkChargeStripHit(stripHit1D);
00447 }
00448 }
00449
00450 else if(type == typeid(SiStripMatchedRecHit2D)){
00451 const SiStripMatchedRecHit2D *matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(hit);
00452 if (matchedHit) {
00453 return (this->isOkChargeStripHit(matchedHit->monoHit())
00454 && this->isOkChargeStripHit(matchedHit->stereoHit()));
00455 }
00456 }
00457 else if(type == typeid(ProjectedSiStripRecHit2D)){
00458
00459 const ProjectedSiStripRecHit2D *projHit = dynamic_cast<const ProjectedSiStripRecHit2D*>(hit);
00460 if (projHit) {
00461 return this->isOkChargeStripHit(&projHit->originalHit());
00462 }
00463 }
00464 else{
00465 edm::LogError("AlignmentTrackSelector")<< "@SUB=isOkCharge" <<"Unknown type of a valid tracker hit in Strips "
00466 << " SubDet = "<<id.subdetId();
00467 return false;
00468 }
00469
00470
00471
00472
00473
00474 edm::LogError("AlignmentTrackSelector")
00475 << "@SUB=isOkCharge" << "Unknown valid tracker hit not in pixel, subdet " << id.subdetId()
00476 << ", SiTrackerMultiRecHit " << dynamic_cast<const SiTrackerMultiRecHit*>(hit)
00477 << ", BaseSiTrackerRecHit2DLocalPos " << dynamic_cast<const BaseSiTrackerRecHit2DLocalPos*>(hit);
00478
00479 return true;
00480 }
00481
00482
00483
00484
00485 bool AlignmentTrackSelector::isOkChargeStripHit(const SiStripRecHit2D *siStripRecHit2D) const
00486 {
00487 double charge = 0.;
00488
00489 SiStripRecHit2D::ClusterRef cluster(siStripRecHit2D->cluster());
00490 const std::vector<uint8_t> &litudes = cluster->amplitudes();
00491
00492 for (size_t ia = 0; ia < amplitudes.size(); ++ia) {
00493 charge += amplitudes[ia];
00494 }
00495
00496 return (charge >= minHitChargeStrip_);
00497 }
00498
00499
00500
00501 bool AlignmentTrackSelector::isOkChargeStripHit(const SiStripRecHit1D *siStripRecHit1D) const
00502 {
00503 double charge = 0.;
00504
00505 SiStripRecHit1D::ClusterRef cluster(siStripRecHit1D->cluster());
00506 const std::vector<uint8_t> &litudes = cluster->amplitudes();
00507
00508 for (size_t ia = 0; ia < amplitudes.size(); ++ia) {
00509 charge += amplitudes[ia];
00510 }
00511
00512 return (charge >= minHitChargeStrip_);
00513 }
00514
00515
00516
00517 bool AlignmentTrackSelector::isIsolated(const TrackingRecHit* therechit, const edm::Event& evt) const
00518 {
00519
00520
00521
00522
00523 edm::Handle<SiStripRecHit2DCollection> rphirecHits;
00524 edm::Handle<SiStripMatchedRecHit2DCollection> matchedrecHits;
00525
00526 evt.getByLabel( rphirecHitsTag_, rphirecHits );
00527 evt.getByLabel( matchedrecHitsTag_, matchedrecHits );
00528
00529 SiStripRecHit2DCollection::DataContainer::const_iterator istripSt;
00530 SiStripMatchedRecHit2DCollection::DataContainer::const_iterator istripStm;
00531 const SiStripRecHit2DCollection::DataContainer& stripcollSt = rphirecHits->data();
00532 const SiStripMatchedRecHit2DCollection::DataContainer& stripcollStm = matchedrecHits->data();
00533
00534 DetId idet = therechit->geographicalId();
00535
00536
00537
00538
00539 for( istripSt=stripcollSt.begin(); istripSt!=stripcollSt.end(); ++istripSt ) {
00540 const SiStripRecHit2D *aHit = &*(istripSt);
00541 DetId mydet1 = aHit->geographicalId();
00542 if (idet.rawId() != mydet1.rawId()) continue;
00543 float theDistance = ( therechit->localPosition() - aHit->localPosition() ).mag();
00544
00545 if (theDistance > 0.001 && theDistance < minHitIsolation_) return false;
00546 }
00547
00548
00549 for( istripStm=stripcollStm.begin(); istripStm!=stripcollStm.end(); ++istripStm ) {
00550 const SiStripMatchedRecHit2D *aHit = &*(istripStm);
00551 DetId mydet2 = aHit->geographicalId();
00552 if (idet.rawId() != mydet2.rawId()) continue;
00553 float theDistance = (therechit->localPosition() - aHit->localPosition()).mag();
00554
00555 if (theDistance > 0.001 && theDistance < minHitIsolation_) return false;
00556 }
00557 return true;
00558 }
00559
00560
00561
00562 AlignmentTrackSelector::Tracks
00563 AlignmentTrackSelector::theNHighestPtTracks(const Tracks& tracks) const
00564 {
00565 Tracks sortedTracks=tracks;
00566 Tracks result;
00567
00568
00569 std::sort(sortedTracks.begin(),sortedTracks.end(),ptComparator);
00570
00571
00572 int n=0;
00573 for (Tracks::const_iterator it=sortedTracks.begin();
00574 it!=sortedTracks.end(); ++it) {
00575 if (n<nHighestPt_) { result.push_back(*it); n++; }
00576 }
00577
00578 return result;
00579 }
00580
00581
00582
00583 AlignmentTrackSelector::Tracks
00584 AlignmentTrackSelector::checkPrescaledHits(const Tracks& tracks, const edm::Event& evt) const
00585 {
00586 Tracks result;
00587
00588
00589 edm::Handle<AliClusterValueMap> fMap;
00590 evt.getByLabel( clusterValueMapTag_, fMap);
00591 const AliClusterValueMap &flagMap=*fMap;
00592
00593
00594 for (Tracks::const_iterator ittrk=tracks.begin(); ittrk != tracks.end(); ++ittrk) {
00595 const reco::Track* trackp=*ittrk;
00596 int ntakenhits=0;
00597
00598
00599 for (trackingRecHit_iterator ith = trackp->recHitsBegin(), edh = trackp->recHitsEnd(); ith != edh; ++ith) {
00600 const TrackingRecHit *hit = ith->get();
00601 if(! hit->isValid())continue;
00602 DetId detid = hit->geographicalId();
00603 int subDet = detid.subdetId();
00604 AlignmentClusterFlag flag;
00605
00606 bool isPixelHit=(subDet == kFPIX || subDet == kBPIX);
00607
00608 if (!isPixelHit){
00609 const std::type_info &type = typeid(*hit);
00610
00611 if (type == typeid(SiStripRecHit2D)) {
00612 const SiStripRecHit2D* striphit=dynamic_cast<const SiStripRecHit2D*>(hit);
00613 if(striphit!=0){
00614 SiStripRecHit2D::ClusterRef stripclust(striphit->cluster());
00615 flag = flagMap[stripclust];
00616 }
00617 }
00618 else if(type == typeid(SiStripRecHit1D)){
00619 const SiStripRecHit1D* striphit = dynamic_cast<const SiStripRecHit1D*>(hit);
00620 if(striphit!=0){
00621 SiStripRecHit1D::ClusterRef stripclust(striphit->cluster());
00622 flag = flagMap[stripclust];
00623 }
00624 }
00625 else{
00626 edm::LogError("AlignmentTrackSelector")<<"ERROR in <AlignmentTrackSelector::checkPrescaledHits>: Dynamic cast of Strip RecHit failed!"
00627 <<" Skipping this hit.";
00628 continue;
00629 }
00630
00631
00632
00633 }
00634 else{
00635 const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
00636 if(pixelhit!=0){
00637 SiPixelRecHit::ClusterRef pixclust(pixelhit->cluster());
00638 flag = flagMap[pixclust];
00639 }
00640 else{
00641 edm::LogError("AlignmentTrackSelector")<<"ERROR in <AlignmentTrackSelector::checkPrescaledHits>: Dynamic cast of Pixel RecHit failed! ";
00642 }
00643 }
00644
00645 if(flag.isTaken())ntakenhits++;
00646
00647 }
00648 if(ntakenhits >= minPrescaledHits_)result.push_back(trackp);
00649 }
00650
00651 return result;
00652 }
00653
00654
00655
00656 bool AlignmentTrackSelector::isOkTrkQuality(const reco::Track* track) const
00657 {
00658 bool qualityOk=false;
00659 bool iterStepOk=false;
00660
00661
00662 if(applyIterStepCheck_){
00663 for (unsigned int i = 0; i < trkSteps_.size(); ++i) {
00664 if (track->algo()==(trkSteps_[i])) {
00665 iterStepOk=true;
00666 }
00667 }
00668 }
00669 else iterStepOk=true;
00670
00671
00672 if(applyTrkQualityCheck_){
00673 for (unsigned int i = 0; i < trkQualities_.size(); ++i) {
00674 if (track->quality(trkQualities_[i])) {
00675 qualityOk=true;
00676 }
00677 }
00678 }
00679 else qualityOk=true;
00680
00681 return qualityOk&&iterStepOk;
00682 }