CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

CosmicTrackingParticleSelector Class Reference

#include <CosmicTrackingParticleSelector.h>

List of all members.

Public Types

typedef TrackingParticleCollection collection
typedef container::const_iterator const_iterator
typedef std::vector< const
TrackingParticle * > 
container

Public Member Functions

const_iterator begin () const
 CosmicTrackingParticleSelector ()
 CosmicTrackingParticleSelector (double ptMin, double minRapidity, double maxRapidity, double tip, double lip, int minHit, bool chargedOnly, std::vector< int > pdgId=std::vector< int >())
 CosmicTrackingParticleSelector (const edm::ParameterSet &cfg)
const_iterator end () const
void initEvent (edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList > simHitsTPAssocToSet) const
bool operator() (const TrackingParticleRef tpr, const reco::BeamSpot *bs, const edm::Event &iEvent, const edm::EventSetup &iSetup) const
void select (const edm::Handle< collection > &c, const edm::Event &event, const edm::EventSetup &setup)
size_t size () const

Private Attributes

bool chargedOnly_
double lip_
double maxRapidity_
int minHit_
double minRapidity_
std::vector< int > pdgId_
double ptMin_
container selected_
edm::Handle
< SimHitTPAssociationProducer::SimHitTPAssociationList
simHitsTPAssoc
double tip_

Detailed Description

Definition at line 35 of file CosmicTrackingParticleSelector.h.


Member Typedef Documentation

Definition at line 38 of file CosmicTrackingParticleSelector.h.

typedef container::const_iterator CosmicTrackingParticleSelector::const_iterator

Definition at line 40 of file CosmicTrackingParticleSelector.h.

Definition at line 39 of file CosmicTrackingParticleSelector.h.


Constructor & Destructor Documentation

CosmicTrackingParticleSelector::CosmicTrackingParticleSelector ( ) [inline]

Definition at line 42 of file CosmicTrackingParticleSelector.h.

{}
CosmicTrackingParticleSelector::CosmicTrackingParticleSelector ( double  ptMin,
double  minRapidity,
double  maxRapidity,
double  tip,
double  lip,
int  minHit,
bool  chargedOnly,
std::vector< int >  pdgId = std::vector<int>() 
) [inline]

Definition at line 44 of file CosmicTrackingParticleSelector.h.

                                                                            :
    ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
    tip_( tip ), lip_( lip ), minHit_( minHit ), chargedOnly_(chargedOnly), pdgId_( pdgId ) { }
CosmicTrackingParticleSelector::CosmicTrackingParticleSelector ( const edm::ParameterSet cfg) [inline]

Definition at line 51 of file CosmicTrackingParticleSelector.h.

                                                                    :
      ptMin_(cfg.getParameter<double>("ptMin")),
      minRapidity_(cfg.getParameter<double>("minRapidity")),
      maxRapidity_(cfg.getParameter<double>("maxRapidity")),
      tip_(cfg.getParameter<double>("tip")),
      lip_(cfg.getParameter<double>("lip")),
      minHit_(cfg.getParameter<int>("minHit")),
      chargedOnly_(cfg.getParameter<bool>("chargedOnly")),
      pdgId_(cfg.getParameter<std::vector<int> >("pdgId")) { }

Member Function Documentation

const_iterator CosmicTrackingParticleSelector::begin ( void  ) const [inline]

Definition at line 73 of file CosmicTrackingParticleSelector.h.

References selected_.

{ return selected_.begin(); }
const_iterator CosmicTrackingParticleSelector::end ( void  ) const [inline]

Definition at line 74 of file CosmicTrackingParticleSelector.h.

References selected_.

{ return selected_.end(); }
void CosmicTrackingParticleSelector::initEvent ( edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList simHitsTPAssocToSet) const [inline]

Definition at line 76 of file CosmicTrackingParticleSelector.h.

References simHitsTPAssoc.

Referenced by MultiTrackValidator::analyze().

                                                                                                            {
      simHitsTPAssoc = simHitsTPAssocToSet;
    }
bool CosmicTrackingParticleSelector::operator() ( const TrackingParticleRef  tpr,
const reco::BeamSpot bs,
const edm::Event iEvent,
const edm::EventSetup iSetup 
) const [inline]

Definition at line 81 of file CosmicTrackingParticleSelector.h.

References chargedOnly_, gather_cfg::cout, PV3DBase< T, PVType, FrameType >::eta(), newFWLiteAna::found, edm::EventSetup::get(), edm::HandleBase::isValid(), TrajectoryStateClosestToBeamLine::isValid(), lip_, maxRapidity_, minHit_, minRapidity_, FreeTrajectoryState::momentum(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::perp2(), FreeTrajectoryState::position(), edm::ESHandle< T >::product(), ptMin_, CosmicsPD_Skims::radius, simHitsTPAssoc, SimHitTPAssociationProducer::simHitTPAssociationListGreater(), mathSSE::sqrt(), GeomDet::surface(), tip_, Surface::toGlobal(), patCandidatesForDimuonsSequences_cff::tracker, TrajectoryStateClosestToBeamLine::trackStateAtPCA(), and PV3DBase< T, PVType, FrameType >::z().

                                                                                                                                        { 
      if (chargedOnly_ && tpr->charge()==0) return false;//select only if charge!=0
      //bool testId = false;
      //unsigned int idSize = pdgId_.size();
      //if (idSize==0) testId = true;
      //else for (unsigned int it=0;it!=idSize;++it){
        //if (tpr->pdgId()==pdgId_[it]) testId = true;
      //}
      
      edm::ESHandle<TrackerGeometry> tracker;
      iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
      
      edm::ESHandle<MagneticField> theMF;
      iSetup.get<IdealMagneticFieldRecord>().get(theMF);
      
      GlobalVector finalGV (0,0,0);
      GlobalPoint finalGP(0,0,0);      
      GlobalVector momentum(0,0,0);//At the PCA
      GlobalPoint vertex(0,0,0);//At the PCA
      double radius(9999);
      bool found(0);


      if (simHitsTPAssoc.isValid()==0) {
        edm::LogError("TrackAssociation") << "Invalid handle!";
        return false;
      }
      std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(tpr,TrackPSimHitRef());//SimHit is dummy: for simHitTPAssociationListGreater 
                                                                                                      // sorting only the cluster is needed
      auto range = std::equal_range(simHitsTPAssoc->begin(), simHitsTPAssoc->end(), 
                                    clusterTPpairWithDummyTP, SimHitTPAssociationProducer::simHitTPAssociationListGreater);
      for(auto ip = range.first; ip != range.second; ++ip) {
        TrackPSimHitRef it = ip->second;
        const GeomDet* tmpDet  = tracker->idToDet( DetId(it->detUnitId()) ) ;
        LocalVector  lv = it->momentumAtEntry();
        Local3DPoint lp = it->localPosition ();
        GlobalVector gv = tmpDet->surface().toGlobal( lv );
        GlobalPoint  gp = tmpDet->surface().toGlobal( lp );
        if(gp.perp()<radius){
          found=true;
          radius = gp.perp();
          finalGV = gv;
          finalGP = gp;
        }
      }
      if(!found) return 0;
      else
        {
          FreeTrajectoryState ftsAtProduction(finalGP,finalGV,TrackCharge(tpr->charge()),theMF.product());
          TSCBLBuilderNoMaterial tscblBuilder;
          TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,*bs);//as in TrackProducerAlgorithm
          if(!tsAtClosestApproach.isValid()){
            std::cout << "WARNING: tsAtClosestApproach is not valid" << std::endl;
            return 0;
          }
          else
            {
              momentum = tsAtClosestApproach.trackStateAtPCA().momentum();
              vertex = tsAtClosestApproach.trackStateAtPCA().position();
              return (
                      tpr->numberOfTrackerLayers() >= minHit_ &&
                      sqrt(momentum.perp2()) >= ptMin_ && 
                      momentum.eta() >= minRapidity_ && momentum.eta() <= maxRapidity_ && 
                      sqrt(vertex.perp2()) <= tip_ &&
                      fabs(vertex.z()) <= lip_ 
                      );
            }
        }
    }
void CosmicTrackingParticleSelector::select ( const edm::Handle< collection > &  c,
const edm::Event event,
const edm::EventSetup setup 
) [inline]

Definition at line 62 of file CosmicTrackingParticleSelector.h.

References SiPixelRawToDigiRegional_cfi::beamSpot, event(), edm::Handle< T >::product(), selected_, and HcalObjRepresent::setup().

                                                                                                       {
        selected_.clear();
        edm::Handle<reco::BeamSpot> beamSpot;
        event.getByLabel(edm::InputTag("offlineBeamSpot"), beamSpot); 
        for( TrackingParticleCollection::const_iterator itp = c->begin(); 
             itp != c->end(); ++ itp )
          if ( operator()(TrackingParticleRef(c,itp-c->begin()),beamSpot.product(),event,setup) ) {
            selected_.push_back( & * itp );
          }
      }
size_t CosmicTrackingParticleSelector::size ( void  ) const [inline]

Definition at line 151 of file CosmicTrackingParticleSelector.h.

References selected_.

{ return selected_.size(); }

Member Data Documentation

Definition at line 161 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().

Definition at line 159 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().

Definition at line 157 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().

Definition at line 160 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().

Definition at line 156 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().

std::vector<int> CosmicTrackingParticleSelector::pdgId_ [private]

Definition at line 162 of file CosmicTrackingParticleSelector.h.

Definition at line 155 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().

Definition at line 163 of file CosmicTrackingParticleSelector.h.

Referenced by begin(), end(), select(), and size().

Definition at line 165 of file CosmicTrackingParticleSelector.h.

Referenced by initEvent(), and operator()().

Definition at line 158 of file CosmicTrackingParticleSelector.h.

Referenced by operator()().