CMS 3D CMS Logo

ctfseeding::HitExtractorSTRP Class Reference

#include <RecoTracker/TkSeedingLayers/src/HitExtractorSTRP.h>

Inheritance diagram for ctfseeding::HitExtractorSTRP:

ctfseeding::HitExtractor

List of all members.

Public Member Functions

virtual HitExtractorSTRPclone () const
 HitExtractorSTRP (const DetLayer *detLayer, SeedingLayer::Side &side, int idLayer)
virtual std::vector< SeedingHithits (const SeedingLayer &sl, const edm::Event &, const edm::EventSetup &) const
void useMatchedHits (const edm::InputTag &m)
void useRingSelector (int minRing, int maxRing)
void useRPhiHits (const edm::InputTag &m)
void useSimpleRphiHitsCleaner (bool use)
void useStereoHits (const edm::InputTag &m)
virtual ~HitExtractorSTRP ()

Private Member Functions

bool ringRange (int ring) const

Private Attributes

bool hasMatchedHits
bool hasRingSelector
bool hasRPhiHits
bool hasSimpleRphiHitsCleaner
bool hasStereoHits
int theIdLayer
const DetLayertheLayer
edm::InputTag theMatchedHits
int theMaxRing
int theMinRing
edm::InputTag theRPhiHits
SeedingLayer::Side theSide
edm::InputTag theStereoHits


Detailed Description

Definition at line 14 of file HitExtractorSTRP.h.


Constructor & Destructor Documentation

HitExtractorSTRP::HitExtractorSTRP ( const DetLayer detLayer,
SeedingLayer::Side side,
int  idLayer 
)

Definition at line 22 of file HitExtractorSTRP.cc.

Referenced by clone().

00024   : theLayer(detLayer), theSide(side), theIdLayer(idLayer),
00025     hasMatchedHits(false), hasRPhiHits(false), hasStereoHits(false),
00026     hasRingSelector(false), theMinRing(1), theMaxRing(0), hasSimpleRphiHitsCleaner(true)
00027 { }

virtual ctfseeding::HitExtractorSTRP::~HitExtractorSTRP (  )  [inline, virtual]

Definition at line 18 of file HitExtractorSTRP.h.

00018 {}


Member Function Documentation

virtual HitExtractorSTRP* ctfseeding::HitExtractorSTRP::clone ( void   )  const [inline, virtual]

Definition at line 21 of file HitExtractorSTRP.h.

References HitExtractorSTRP().

Referenced by SeedingLayerSetsBuilder::layers().

00021 { return new HitExtractorSTRP(*this); }

vector< SeedingHit > HitExtractorSTRP::hits ( const SeedingLayer sl,
const edm::Event ev,
const edm::EventSetup es 
) const [virtual]

Implements ctfseeding::HitExtractor.

Definition at line 43 of file HitExtractorSTRP.cc.

References edm::Event::getByLabel(), hasMatchedHits, hasRPhiHits, hasSimpleRphiHitsCleaner, hasStereoHits, it, range, HLT_VtxMuL3::result, ringRange(), TrackerLayerIdAccessor::stripTECDisk(), TrackerLayerIdAccessor::stripTIBLayer(), TrackerLayerIdAccessor::stripTIDDisk(), TrackerLayerIdAccessor::stripTOBLayer(), DetLayer::subDetector(), GeomDetEnumerators::TEC, theIdLayer, theLayer, theMatchedHits, theRPhiHits, theSide, theStereoHits, GeomDetEnumerators::TIB, GeomDetEnumerators::TID, and GeomDetEnumerators::TOB.

00044 {
00045   TrackerLayerIdAccessor accessor;
00046   std::vector<SeedingHit> result;
00047 
00048   //
00049   // TIB
00050   //
00051   if (theLayer->subDetector() == GeomDetEnumerators::TIB) {
00052     if (hasMatchedHits) {
00053       edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00054       ev.getByLabel( theMatchedHits, matchedHits);
00055       const SiStripMatchedRecHit2DCollection::range range =
00056         matchedHits->get(accessor.stripTIBLayer(theIdLayer) );
00057       for(SiStripMatchedRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00058         result.push_back( SeedingHit(&(*it), sl, es) );
00059       }
00060     }
00061     if (hasRPhiHits) {
00062       edm::Handle<SiStripRecHit2DCollection> rphiHits;
00063       ev.getByLabel( theRPhiHits, rphiHits);
00064       if (hasMatchedHits){ 
00065         if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
00066           const SiStripRecHit2DCollection::range range =
00067             rphiHits->get(accessor.stripTIBLayer(theIdLayer) );
00068           for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00069             result.push_back( SeedingHit(&(*it), sl, es) );
00070           }
00071         }
00072       } else {
00073         const SiStripRecHit2DCollection::range range =
00074           rphiHits->get(accessor.stripTIBLayer(theIdLayer) );
00075         for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00076           result.push_back( SeedingHit(&(*it), sl, es) );
00077         }
00078       }
00079     }
00080     if (hasStereoHits) {
00081       edm::Handle<SiStripRecHit2DCollection> stereoHits;
00082       ev.getByLabel( theStereoHits, stereoHits);
00083       const SiStripRecHit2DCollection::range range =
00084         stereoHits->get(accessor.stripTIBLayer(theIdLayer) );
00085       for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00086         result.push_back( SeedingHit(&(*it), sl, es) );
00087       }
00088     }
00089   }
00090   
00091   //
00092   // TID
00093   //
00094   else if (theLayer->subDetector() == GeomDetEnumerators::TID) {
00095     if (hasMatchedHits) {
00096       edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00097       ev.getByLabel( theMatchedHits, matchedHits);
00098       const SiStripMatchedRecHit2DCollection::range range =
00099         matchedHits->get(accessor.stripTIDDisk(theSide,theIdLayer) );
00100       for(SiStripMatchedRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00101         int ring = TIDDetId( it->geographicalId() ).ring();
00102         if (ringRange(ring))result.push_back( SeedingHit(&(*it), sl, es) );
00103       }
00104     }
00105     if (hasRPhiHits) {
00106       edm::Handle<SiStripRecHit2DCollection> rphiHits;
00107       ev.getByLabel( theRPhiHits, rphiHits);
00108       const SiStripRecHit2DCollection::range range =
00109         rphiHits->get(accessor.stripTIDDisk(theSide,theIdLayer) );
00110       for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00111         int ring = TIDDetId( it->geographicalId() ).ring();
00112         if (ringRange(ring)){
00113           bool  isInGlued = SiStripDetId(it->geographicalId() ).partnerDetId();
00114           if (hasMatchedHits){
00115             if (isInGlued){
00116               // this is a brutal "cleaning". Add something smarter in the future
00117               if(!hasSimpleRphiHitsCleaner) result.push_back( SeedingHit(&(*it), sl, es) );
00118             }else{
00119               result.push_back( SeedingHit(&(*it), sl, es) );
00120             }
00121           } else {
00122             result.push_back( SeedingHit(&(*it), sl, es) );
00123           }
00124         }
00125       }
00126     }
00127     if (hasStereoHits) {
00128       edm::Handle<SiStripRecHit2DCollection> stereoHits;
00129         ev.getByLabel( theStereoHits, stereoHits);
00130         const SiStripRecHit2DCollection::range range =
00131                         stereoHits->get(accessor.stripTIDDisk(theSide,theIdLayer) );
00132         for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00133          int ring = TIDDetId( it->geographicalId() ).ring();
00134            if (ringRange(ring))result.push_back( SeedingHit(&(*it), sl, es) );
00135         }  
00136     }
00137   }
00138   //
00139   // TOB
00140   //
00141   else if (theLayer->subDetector() == GeomDetEnumerators::TOB) {
00142     if (hasMatchedHits) {
00143       edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00144       ev.getByLabel( theMatchedHits, matchedHits);
00145       const SiStripMatchedRecHit2DCollection::range range =
00146         matchedHits->get(accessor.stripTOBLayer(theIdLayer) );
00147       for(SiStripMatchedRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00148         result.push_back( SeedingHit(&(*it), sl, es) );
00149       }
00150     }
00151     if (hasRPhiHits) {
00152       edm::Handle<SiStripRecHit2DCollection> rphiHits;
00153       ev.getByLabel( theRPhiHits, rphiHits);
00154       if (hasMatchedHits){ 
00155         if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
00156           const SiStripRecHit2DCollection::range range =
00157             rphiHits->get(accessor.stripTOBLayer(theIdLayer) );
00158           for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00159             result.push_back( SeedingHit(&(*it), sl, es) );
00160           }
00161         }
00162       } else {
00163         const SiStripRecHit2DCollection::range range =
00164           rphiHits->get(accessor.stripTOBLayer(theIdLayer) );
00165         for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00166           result.push_back( SeedingHit(&(*it), sl, es) );
00167         }
00168       }
00169     }
00170     if (hasStereoHits) {
00171       edm::Handle<SiStripRecHit2DCollection> stereoHits;
00172       ev.getByLabel( theStereoHits, stereoHits);
00173       const SiStripRecHit2DCollection::range range =
00174         stereoHits->get(accessor.stripTOBLayer(theIdLayer) );
00175       for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00176         result.push_back( SeedingHit(&(*it), sl, es) );
00177       }
00178     }
00179   }
00180 
00181   //
00182   // TEC
00183   //
00184   else if (theLayer->subDetector() == GeomDetEnumerators::TEC) {
00185     //std::cout << "=== in TEC, hasSimpleCleaner: " << hasSimpleRphiHitsCleaner << std::endl;    
00186     if (hasMatchedHits) {
00187       edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00188       ev.getByLabel( theMatchedHits, matchedHits);
00189       const SiStripMatchedRecHit2DCollection::range range =
00190         matchedHits->get(accessor.stripTECDisk(theSide,theIdLayer) );
00191       for(SiStripMatchedRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00192         int ring = TECDetId( it->geographicalId() ).ring();
00193         if (ringRange(ring))result.push_back( SeedingHit(&(*it), sl, es) );
00194       }
00195     }
00196     if (hasRPhiHits) {
00197       edm::Handle<SiStripRecHit2DCollection> rphiHits;
00198       ev.getByLabel( theRPhiHits, rphiHits);
00199       const SiStripRecHit2DCollection::range range =
00200         rphiHits->get(accessor.stripTECDisk(theSide,theIdLayer) );
00201       for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00202         int ring = TECDetId( it->geographicalId() ).ring();
00203         //      std::cout << "layer " << theIdLayer << " ring " << ring << std::endl;
00204         if (ringRange(ring)){
00205           bool  isInGlued = SiStripDetId(it->geographicalId() ).partnerDetId();
00206           if (hasMatchedHits){
00207             if(isInGlued){
00208               // this is a brutal "cleaning". Add something smarter in the future
00209               if (!hasSimpleRphiHitsCleaner) result.push_back( SeedingHit(&(*it), sl, es) );       
00210             }else{
00211               result.push_back( SeedingHit(&(*it), sl, es) );
00212             }
00213           }else {
00214             result.push_back( SeedingHit(&(*it), sl, es) );
00215           }
00216         }
00217       }
00218     }
00219     if (hasStereoHits) {
00220       edm::Handle<SiStripRecHit2DCollection> stereoHits;
00221       ev.getByLabel( theStereoHits, stereoHits);
00222       const SiStripRecHit2DCollection::range range =
00223         stereoHits->get(accessor.stripTECDisk(theSide,theIdLayer) );
00224       for(SiStripRecHit2DCollection::const_iterator it=range.first; it!=range.second; it++){
00225         int ring = TECDetId( it->geographicalId() ).ring();
00226         if (ringRange(ring))result.push_back( SeedingHit(&(*it), sl, es) );
00227       }
00228     }
00229   }
00230   
00231 
00232   return result;
00233 }

bool HitExtractorSTRP::ringRange ( int  ring  )  const [private]

Definition at line 36 of file HitExtractorSTRP.cc.

References hasRingSelector, theMaxRing, and theMinRing.

Referenced by hits().

00037 {
00038   if (!hasRingSelector) return true;
00039   else if ( ring >= theMinRing && ring <= theMaxRing) return true;
00040   else return false;
00041 }

void ctfseeding::HitExtractorSTRP::useMatchedHits ( const edm::InputTag m  )  [inline]

Definition at line 23 of file HitExtractorSTRP.h.

References hasMatchedHits, and theMatchedHits.

Referenced by SeedingLayerSetsBuilder::layers().

00023 { hasMatchedHits = true; theMatchedHits = m; }

void HitExtractorSTRP::useRingSelector ( int  minRing,
int  maxRing 
)

Definition at line 29 of file HitExtractorSTRP.cc.

References hasRingSelector, theMaxRing, and theMinRing.

Referenced by SeedingLayerSetsBuilder::layers().

00030 {
00031   hasRingSelector=true;
00032   theMinRing=minRing;
00033   theMaxRing=maxRing; 
00034 }

void ctfseeding::HitExtractorSTRP::useRPhiHits ( const edm::InputTag m  )  [inline]

Definition at line 24 of file HitExtractorSTRP.h.

References hasRPhiHits, and theRPhiHits.

Referenced by SeedingLayerSetsBuilder::layers().

00024 { hasRPhiHits    = true; theRPhiHits = m; }

void ctfseeding::HitExtractorSTRP::useSimpleRphiHitsCleaner ( bool  use  )  [inline]

Definition at line 27 of file HitExtractorSTRP.h.

References hasSimpleRphiHitsCleaner.

Referenced by SeedingLayerSetsBuilder::layers().

00027 {hasSimpleRphiHitsCleaner = use;}

void ctfseeding::HitExtractorSTRP::useStereoHits ( const edm::InputTag m  )  [inline]

Definition at line 25 of file HitExtractorSTRP.h.

References hasStereoHits, and theStereoHits.

Referenced by SeedingLayerSetsBuilder::layers().

00025 { hasStereoHits = true; theStereoHits = m; }


Member Data Documentation

bool ctfseeding::HitExtractorSTRP::hasMatchedHits [private]

Definition at line 35 of file HitExtractorSTRP.h.

Referenced by hits(), and useMatchedHits().

bool ctfseeding::HitExtractorSTRP::hasRingSelector [private]

Definition at line 38 of file HitExtractorSTRP.h.

Referenced by ringRange(), and useRingSelector().

bool ctfseeding::HitExtractorSTRP::hasRPhiHits [private]

Definition at line 36 of file HitExtractorSTRP.h.

Referenced by hits(), and useRPhiHits().

bool ctfseeding::HitExtractorSTRP::hasSimpleRphiHitsCleaner [private]

Definition at line 39 of file HitExtractorSTRP.h.

Referenced by hits(), and useSimpleRphiHitsCleaner().

bool ctfseeding::HitExtractorSTRP::hasStereoHits [private]

Definition at line 37 of file HitExtractorSTRP.h.

Referenced by hits(), and useStereoHits().

int ctfseeding::HitExtractorSTRP::theIdLayer [private]

Definition at line 34 of file HitExtractorSTRP.h.

Referenced by hits().

const DetLayer* ctfseeding::HitExtractorSTRP::theLayer [private]

Definition at line 32 of file HitExtractorSTRP.h.

Referenced by hits().

edm::InputTag ctfseeding::HitExtractorSTRP::theMatchedHits [private]

Definition at line 35 of file HitExtractorSTRP.h.

Referenced by hits(), and useMatchedHits().

int ctfseeding::HitExtractorSTRP::theMaxRing [private]

Definition at line 38 of file HitExtractorSTRP.h.

Referenced by ringRange(), and useRingSelector().

int ctfseeding::HitExtractorSTRP::theMinRing [private]

Definition at line 38 of file HitExtractorSTRP.h.

Referenced by ringRange(), and useRingSelector().

edm::InputTag ctfseeding::HitExtractorSTRP::theRPhiHits [private]

Definition at line 36 of file HitExtractorSTRP.h.

Referenced by hits(), and useRPhiHits().

SeedingLayer::Side ctfseeding::HitExtractorSTRP::theSide [private]

Definition at line 33 of file HitExtractorSTRP.h.

Referenced by hits().

edm::InputTag ctfseeding::HitExtractorSTRP::theStereoHits [private]

Definition at line 37 of file HitExtractorSTRP.h.

Referenced by hits(), and useStereoHits().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:37:09 2009 for CMSSW by  doxygen 1.5.4