CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CosmicTrackingParticleSelector.h
Go to the documentation of this file.
1 #ifndef RecoSelectors_CosmicTrackingParticleSelector_h
2 #define RecoSelectors_CosmicTrackingParticleSelector_h
3 /* \class CosmicTrackingParticleSelector
4  *
5  * \author Yanyan Gao, FNAL
6  *
7  * $Date: 2013/06/24 12:25:14 $
8  * $Revision: 1.5 $
9  *
10  */
22 
30 
33 
35 
37 
38 public:
40  typedef std::vector<const TrackingParticle *> container;
41  typedef container::const_iterator const_iterator;
42 
44 
45  CosmicTrackingParticleSelector ( double ptMin,double minRapidity,double maxRapidity,
46  double tip,double lip,int minHit, bool chargedOnly,
47  const std::vector<int>& pdgId = std::vector<int>()) :
48  ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
49  tip_( tip ), lip_( lip ), minHit_( minHit ), chargedOnly_(chargedOnly), pdgId_( pdgId ) { }
50 
51 
53  ptMin_(cfg.getParameter<double>("ptMin")),
54  minRapidity_(cfg.getParameter<double>("minRapidity")),
55  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
56  tip_(cfg.getParameter<double>("tip")),
57  lip_(cfg.getParameter<double>("lip")),
58  minHit_(cfg.getParameter<int>("minHit")),
59  chargedOnly_(cfg.getParameter<bool>("chargedOnly")),
60  pdgId_(cfg.getParameter<std::vector<int> >("pdgId")),
61  beamSpotToken_(iC.consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"))) { }
62 
63 
65  selected_.clear();
67  event.getByToken(beamSpotToken_, beamSpot);
68  for( TrackingParticleCollection::const_iterator itp = c->begin();
69  itp != c->end(); ++ itp )
70  if ( operator()(TrackingParticleRef(c,itp-c->begin()),beamSpot.product(),event,setup) ) {
71  selected_.push_back( & * itp );
72  }
73  }
74 
75  const_iterator begin() const { return selected_.begin(); }
76  const_iterator end() const { return selected_.end(); }
77 
79  simHitsTPAssoc = simHitsTPAssocToSet;
80  }
81 
82  // Operator() performs the selection: e.g. if (tPSelector(tp, bs, event, evtsetup)) {...
83  bool operator()( const TrackingParticleRef tpr, const reco::BeamSpot* bs, const edm::Event &iEvent, const edm::EventSetup& iSetup ) const {
84  if (chargedOnly_ && tpr->charge()==0) return false;//select only if charge!=0
85  //bool testId = false;
86  //unsigned int idSize = pdgId_.size();
87  //if (idSize==0) testId = true;
88  //else for (unsigned int it=0;it!=idSize;++it){
89  //if (tpr->pdgId()==pdgId_[it]) testId = true;
90  //}
91 
93  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
94 
96  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
97 
98  GlobalVector finalGV (0,0,0);
99  GlobalPoint finalGP(0,0,0);
100  GlobalVector momentum(0,0,0);//At the PCA
101  GlobalPoint vertex(0,0,0);//At the PCA
102  double radius(9999);
103  bool found(0);
104 
105 
106  if (simHitsTPAssoc.isValid()==0) {
107  edm::LogError("TrackAssociation") << "Invalid handle!";
108  return false;
109  }
110  std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(tpr,TrackPSimHitRef());//SimHit is dummy: for simHitTPAssociationListGreater
111  // sorting only the cluster is needed
112  auto range = std::equal_range(simHitsTPAssoc->begin(), simHitsTPAssoc->end(),
114  for(auto ip = range.first; ip != range.second; ++ip) {
115  TrackPSimHitRef it = ip->second;
116  const GeomDet* tmpDet = tracker->idToDet( DetId(it->detUnitId()) ) ;
117  LocalVector lv = it->momentumAtEntry();
118  Local3DPoint lp = it->localPosition ();
119  GlobalVector gv = tmpDet->surface().toGlobal( lv );
120  GlobalPoint gp = tmpDet->surface().toGlobal( lp );
121  if(gp.perp()<radius){
122  found=true;
123  radius = gp.perp();
124  finalGV = gv;
125  finalGP = gp;
126  }
127  }
128  if(!found) return 0;
129  else
130  {
131  FreeTrajectoryState ftsAtProduction(finalGP,finalGV,TrackCharge(tpr->charge()),theMF.product());
132  TSCBLBuilderNoMaterial tscblBuilder;
133  TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,*bs);//as in TrackProducerAlgorithm
134  if(!tsAtClosestApproach.isValid()){
135  std::cout << "WARNING: tsAtClosestApproach is not valid" << std::endl;
136  return 0;
137  }
138  else
139  {
140  momentum = tsAtClosestApproach.trackStateAtPCA().momentum();
141  vertex = tsAtClosestApproach.trackStateAtPCA().position();
142  return (
143  tpr->numberOfTrackerLayers() >= minHit_ &&
144  sqrt(momentum.perp2()) >= ptMin_ &&
145  momentum.eta() >= minRapidity_ && momentum.eta() <= maxRapidity_ &&
146  sqrt(vertex.perp2()) <= tip_ &&
147  fabs(vertex.z()) <= lip_
148  );
149  }
150  }
151  }
152 
153  size_t size() const { return selected_.size(); }
154 
155  private:
156 
157  double ptMin_;
158  double minRapidity_;
159  double maxRapidity_;
160  double tip_;
161  double lip_;
162  int minHit_;
164  std::vector<int> pdgId_;
167 
169 
170 };
171 
172 #endif
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
T perp() const
Definition: PV3DBase.h:72
std::vector< TrackingParticle > TrackingParticleCollection
static bool simHitTPAssociationListGreater(SimHitTPPair i, SimHitTPPair j)
void select(const edm::Handle< collection > &c, const edm::Event &event, const edm::EventSetup &setup)
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
T perp2() const
Definition: PV3DBase.h:71
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
bool operator()(const TrackingParticleRef tpr, const reco::BeamSpot *bs, const edm::Event &iEvent, const edm::EventSetup &iSetup) const
int TrackCharge
Definition: TrackCharge.h:4
int iEvent
Definition: GenABIO.cc:230
CosmicTrackingParticleSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
std::vector< const TrackingParticle * > container
void initEvent(edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList > simHitsTPAssocToSet) const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList > simHitsTPAssoc
GlobalVector momentum() const
Definition: DetId.h:18
GlobalPoint position() const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:81
T eta() const
Definition: PV3DBase.h:76
edm::Ref< edm::PSimHitContainer > TrackPSimHitRef
tuple cout
Definition: gather_cfg.py:121
CosmicTrackingParticleSelector(double ptMin, double minRapidity, double maxRapidity, double tip, double lip, int minHit, bool chargedOnly, const std::vector< int > &pdgId=std::vector< int >())
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
edm::Ref< TrackingParticleCollection > TrackingParticleRef