CMS 3D CMS Logo

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 
27 
33 
36 
38 
40 public:
42  typedef std::vector<const TrackingParticle*> container;
43  typedef container::const_iterator const_iterator;
44 
46 
48  double minRapidity,
49  double maxRapidity,
50  double tip,
51  double lip,
52  int minHit,
53  bool chargedOnly,
54  const std::vector<int>& pdgId = std::vector<int>())
55  : ptMin_(ptMin),
58  tip_(tip),
59  lip_(lip),
60  minHit_(minHit),
62  pdgId_(pdgId) {}
63 
65  : ptMin_(cfg.getParameter<double>("ptMin")),
66  minRapidity_(cfg.getParameter<double>("minRapidity")),
67  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
68  tip_(cfg.getParameter<double>("tip")),
69  lip_(cfg.getParameter<double>("lip")),
70  minHit_(cfg.getParameter<int>("minHit")),
71  chargedOnly_(cfg.getParameter<bool>("chargedOnly")),
72  pdgId_(cfg.getParameter<std::vector<int> >("pdgId")),
73  beamSpotToken_(iC.consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"))),
75  theMFToken_(iC.esConsumes()) {}
76 
78  selected_.clear();
80  event.getByToken(beamSpotToken_, beamSpot);
81  for (TrackingParticleCollection::const_iterator itp = c->begin(); itp != c->end(); ++itp)
82  if (operator()(TrackingParticleRef(c, itp - c->begin()), beamSpot.product(), event, setup)) {
83  selected_.push_back(&*itp);
84  }
85  }
86 
87  const_iterator begin() const { return selected_.begin(); }
88  const_iterator end() const { return selected_.end(); }
89 
91  simHitsTPAssoc = simHitsTPAssocToSet;
92  }
93 
94  // Operator() performs the selection: e.g. if (tPSelector(tp, bs, event, evtsetup)) {...
96  const reco::BeamSpot* bs,
97  const edm::Event& iEvent,
98  const edm::EventSetup& iSetup) const {
99  if (chargedOnly_ && tpr->charge() == 0)
100  return false; //select only if charge!=0
101  //bool testId = false;
102  //unsigned int idSize = pdgId_.size();
103  //if (idSize==0) testId = true;
104  //else for (unsigned int it=0;it!=idSize;++it){
105  //if (tpr->pdgId()==pdgId_[it]) testId = true;
106  //}
107 
109 
110  GlobalVector finalGV(0, 0, 0);
111  GlobalPoint finalGP(0, 0, 0);
112  GlobalVector momentum(0, 0, 0); //At the PCA
113  GlobalPoint vertex(0, 0, 0); //At the PCA
114  double radius(9999);
115  bool found(false);
116 
117  int ii = 0;
118  DetId::Detector det;
119  int subdet;
120 
121  edm::LogVerbatim("CosmicTrackingParticleSelector")
122  << "TOT Number of PSimHits = " << tpr->numberOfHits()
123  << ", Number of Tracker PSimHits = " << tpr->numberOfTrackerHits() << "\n";
124 
125  if (simHitsTPAssoc.isValid() == 0) {
126  edm::LogError("CosmicTrackingParticleSelector") << "Invalid handle!";
127  return false;
128  }
129  std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(
130  tpr, TrackPSimHitRef()); //SimHit is dummy: for simHitTPAssociationListGreater
131  // sorting only the cluster is needed
132  auto range = std::equal_range(simHitsTPAssoc->begin(),
133  simHitsTPAssoc->end(),
134  clusterTPpairWithDummyTP,
136  for (auto ip = range.first; ip != range.second; ++ip) {
137  TrackPSimHitRef it = ip->second;
138  ++ii;
139  const GeomDet* tmpDet = theGeometry->idToDet(DetId(it->detUnitId()));
140  if (!tmpDet) {
141  edm::LogVerbatim("CosmicTrackingParticleSelector")
142  << "***WARNING: PSimHit " << ii << ", no GeomDet for: " << it->detUnitId() << ". Skipping it.";
143  continue;
144  } else {
145  det = DetId(it->detUnitId()).det();
146  subdet = DetId(it->detUnitId()).subdetId();
147  }
148 
149  LocalVector lv = it->momentumAtEntry();
150  Local3DPoint lp = it->localPosition();
151  GlobalVector gv = tmpDet->surface().toGlobal(lv);
152  GlobalPoint gp = tmpDet->surface().toGlobal(lp);
153  edm::LogVerbatim("CosmicTrackingParticleSelector")
154  << "PSimHit " << ii << ", Detector = " << det << ", subdet = " << subdet << "\t Radius = " << gp.perp()
155  << ", z = " << gp.z() << "\t pt = " << gv.perp() << ", pz = " << gv.z();
156  edm::LogVerbatim("CosmicTrackingParticleSelector")
157  << "\t trackId = " << it->trackId() << ", particleType = " << it->particleType()
158  << ", processType = " << it->processType();
159 
160  // discard hits related to low energy debris from the primary particle
161  if (it->processType() != 0)
162  continue;
163 
164  if (gp.perp() < radius) {
165  found = true;
166  radius = gp.perp();
167  finalGV = gv;
168  finalGP = gp;
169  }
170  }
171  edm::LogVerbatim("CosmicTrackingParticleSelector")
172  << "\n"
173  << "FINAL State at InnerMost Hit: Radius = " << finalGP.perp() << ", z = " << finalGP.z()
174  << ", pt = " << finalGV.perp() << ", pz = " << finalGV.z();
175 
176  if (!found)
177  return false;
178  else {
179  FreeTrajectoryState ftsAtProduction(finalGP, finalGV, TrackCharge(tpr->charge()), &iSetup.getData(theMFToken_));
180  TSCBLBuilderNoMaterial tscblBuilder;
181  //as in TrackProducerAlgorithm
182  TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction, *bs);
183  if (!tsAtClosestApproach.isValid()) {
184  edm::LogVerbatim("CosmicTrackingParticleSelector")
185  << "*** WARNING in CosmicTrackingParticleSelector: tsAtClosestApproach is not valid."
186  << "\n";
187  return false;
188  } else {
189  momentum = tsAtClosestApproach.trackStateAtPCA().momentum();
190  vertex = tsAtClosestApproach.trackStateAtPCA().position();
191 
192  edm::LogVerbatim("CosmicTrackingParticleSelector")
193  << "FINAL State extrapolated at PCA: Radius = " << vertex.perp() << ", z = " << vertex.z()
194  << ", pt = " << momentum.perp() << ", pz = " << momentum.z() << "\n";
195 
196  return (tpr->numberOfTrackerLayers() >= minHit_ && sqrt(momentum.perp2()) >= ptMin_ &&
197  momentum.eta() >= minRapidity_ && momentum.eta() <= maxRapidity_ && sqrt(vertex.perp2()) <= tip_ &&
198  fabs(vertex.z()) <= lip_);
199  }
200  }
201  }
202 
203  size_t size() const { return selected_.size(); }
204 
205 private:
206  double ptMin_;
207  double minRapidity_;
208  double maxRapidity_;
209  double tip_;
210  double lip_;
211  int minHit_;
213  std::vector<int> pdgId_;
218 
220 };
221 
222 #endif
Log< level::Info, true > LogVerbatim
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T perp() const
Definition: PV3DBase.h:69
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theMFToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static bool simHitTPAssociationListGreater(SimHitTPPair i, SimHitTPPair j)
void select(const edm::Handle< collection > &c, const edm::Event &event, const edm::EventSetup &setup)
T z() const
Definition: PV3DBase.h:61
T eta() const
Definition: PV3DBase.h:73
constexpr float ptMin
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
Log< level::Error, false > LogError
const GeomDet * idToDet(DetId) const override
GlobalPoint position() const
int TrackCharge
Definition: TrackCharge.h:4
int iEvent
Definition: GenABIO.cc:224
std::vector< const TrackingParticle * > container
CosmicTrackingParticleSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
T sqrt(T t)
Definition: SSEVec.h:23
GlobalVector momentum() const
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList > simHitsTPAssoc
ii
Definition: cuy.py:589
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
Definition: DetId.h:17
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > globalTrackingGeomToken_
Detector
Definition: DetId.h:24
T perp2() const
Definition: PV3DBase.h:68
bool isValid() const
Definition: HandleBase.h:70
fixed size matrix
HLT enums.
edm::Ref< edm::PSimHitContainer > TrackPSimHitRef
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
std::vector< TrackingParticle > TrackingParticleCollection
CosmicTrackingParticleSelector(double ptMin, double minRapidity, double maxRapidity, double tip, double lip, int minHit, bool chargedOnly, const std::vector< int > &pdgId=std::vector< int >())
edm::Ref< TrackingParticleCollection > TrackingParticleRef
Definition: event.py:1