CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/Validation/RecoHI/plugins/HitPixelLayersTPSelector.h

Go to the documentation of this file.
00001 #ifndef HitPixelLayersTrackSelection_h
00002 #define HitPixelLayersTrackSelection_h
00003 
00004 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00005 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
00006 
00007 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00008 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00009 
00017 class HitPixelLayersTPSelector
00018 {
00019   
00020  public:
00021   // input collection type
00022   typedef TrackingParticleCollection collection;
00023   
00024   
00025   // output collection type
00026   typedef std::vector<const TrackingParticle*> container;
00027   
00028   // iterator over result collection type.
00029   typedef container::const_iterator const_iterator;
00030   
00031   // constructor from parameter set configurability
00032   HitPixelLayersTPSelector(const edm::ParameterSet & iConfig) : 
00033     tripletSeedOnly_(iConfig.getParameter<bool>("tripletSeedOnly")),
00034     ptMin_(iConfig.getParameter<double>("ptMin")), 
00035     minRapidity_(iConfig.getParameter<double>("minRapidity")), 
00036     maxRapidity_(iConfig.getParameter<double>("maxRapidity")),
00037     tip_(iConfig.getParameter<double>("tip")), 
00038     lip_(iConfig.getParameter<double>("lip")), 
00039     minHit_(iConfig.getParameter<int>("minHit")), 
00040     signalOnly_(iConfig.getParameter<bool>("signalOnly")), 
00041     chargedOnly_(iConfig.getParameter<bool>("chargedOnly")), 
00042     primaryOnly_(iConfig.getParameter<bool>("primaryOnly")),
00043     tpStatusBased_(iConfig.getParameter<bool>("tpStatusBased")),
00044     pdgId_(iConfig.getParameter< std::vector<int> >("pdgId"))                              
00045       {};
00046     
00047     // select object from a collection and
00048     // possibly event content
00049     void select( const edm::Handle<collection> & TPCH, const edm::Event & iEvent, const edm::EventSetup & iSetup)
00050     {
00051       selected_.clear();
00052       
00053 
00054       const collection & tpc = *(TPCH.product());
00055       
00056       for (TrackingParticleCollection::size_type i=0; i<tpc.size(); i++)
00057         {
00058           TrackingParticleRef tpr(TPCH, i);
00059           
00060           // quickly reject if it is from pile-up
00061           if (signalOnly_ && !(tpr->eventId().bunchCrossing()==0 && tpr->eventId().event()==0) ) continue;
00062           if (chargedOnly_ && tpr->charge()==0) continue; //select only if charge!=0
00063           if (tpStatusBased_ && primaryOnly_ && tpr->status()!=1 ) continue; // TP status based sel primary
00064           if ((!tpStatusBased_) && primaryOnly_ && tpr->parentVertex()->nSourceTracks()!=0 ) continue; // vertex based sel for primary
00065 
00066           // loop over specified PID values 
00067           bool testId = false;
00068           unsigned int idSize = pdgId_.size();
00069           if (idSize==0) testId = true;
00070           else for (unsigned int it=0;it!=idSize;++it){
00071             if (tpr->pdgId()==pdgId_[it]) testId = true;
00072           }
00073           
00074           // selection criteria
00075           if ( tpr->matchedHit() >= minHit_ &&   
00076                sqrt(tpr->momentum().perp2()) >= ptMin_ && 
00077                tpr->momentum().eta() >= minRapidity_ && tpr->momentum().eta() <= maxRapidity_ && 
00078                sqrt(tpr->vertex().perp2()) <= tip_ &&
00079                fabs(tpr->vertex().z()) <= lip_ &&
00080                testId)
00081             {
00082               if (tripletSeedOnly_ && !goodHitPattern(pixelHitPattern(tpr)) ) continue; //findable triplet seed
00083               const TrackingParticle * trap = &(tpc[i]);
00084               selected_.push_back(trap);
00085             } 
00086           
00087         }
00088     }
00089     
00090     // return pixel layer hit pattern
00091     std::vector<bool> pixelHitPattern( const TrackingParticleRef& simTrack )
00092       {
00093         std::vector<bool> hitpattern(5,false); // PXB 0,1,2  PXF 0,1
00094         
00095         for(std::vector<PSimHit>::const_iterator simHit = simTrack->pSimHit_begin();simHit!= simTrack->pSimHit_end();simHit++){
00096           
00097           DetId id = DetId(simHit->detUnitId());
00098           uint32_t detid = id.det();
00099           uint32_t subdet = id.subdetId();
00100           
00101           if (detid == DetId::Tracker) {
00102             if (subdet == PixelSubdetector::PixelBarrel) 
00103               hitpattern[PXBDetId(id).layer()-1]=true;
00104             else if (subdet == PixelSubdetector::PixelEndcap) 
00105               hitpattern[PXFDetId(id).disk()+2]=true;
00106           }
00107           
00108         }// end simhit loop
00109         
00110         return hitpattern;
00111       }
00112     
00113     // test whether hit pattern would give a pixel triplet seed
00114     bool goodHitPattern( const std::vector<bool> hitpattern )
00115     {
00116       if( (hitpattern[0] && hitpattern[1] && hitpattern[2]) ||
00117           (hitpattern[0] && hitpattern[1] && hitpattern[3]) ||
00118           (hitpattern[0] && hitpattern[3] && hitpattern[4]) ) 
00119         return true;
00120       else 
00121         return false;
00122     }
00123     
00124     // iterators over selected objects: collection begin
00125     const_iterator begin() const
00126     {
00127       return selected_.begin();
00128     }
00129     
00130     // iterators over selected objects: collection end
00131     const_iterator end() const
00132     {
00133       return selected_.end();
00134     }
00135     
00136     // true if no object has been selected
00137     size_t size() const
00138     {
00139       return selected_.size();
00140     }
00141     
00142     //private:
00143     
00144     container selected_;
00145     bool tripletSeedOnly_;
00146     double ptMin_;
00147     double minRapidity_;
00148     double maxRapidity_;
00149     double tip_;
00150     double lip_;
00151     int    minHit_;
00152     bool signalOnly_;
00153     bool chargedOnly_;
00154     bool primaryOnly_;
00155     bool tpStatusBased_;
00156     std::vector<int> pdgId_;
00157     
00158     
00159 };
00160 
00161 
00162 #endif