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: 2011/10/27 12:52:48 $
8  * $Revision: 1.3 $
9  *
10  */
21 
29 
31 
33 
34 public:
36  typedef std::vector<const TrackingParticle *> container;
37  typedef container::const_iterator const_iterator;
38 
40 
41  CosmicTrackingParticleSelector ( double ptMin,double minRapidity,double maxRapidity,
42  double tip,double lip,int minHit, bool chargedOnly,
43  std::vector<int> pdgId = std::vector<int>()) :
44  ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
45  tip_( tip ), lip_( lip ), minHit_( minHit ), chargedOnly_(chargedOnly), pdgId_( pdgId ) { }
46 
47 
49  ptMin_(cfg.getParameter<double>("ptMin")),
50  minRapidity_(cfg.getParameter<double>("minRapidity")),
51  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
52  tip_(cfg.getParameter<double>("tip")),
53  lip_(cfg.getParameter<double>("lip")),
54  minHit_(cfg.getParameter<int>("minHit")),
55  chargedOnly_(cfg.getParameter<bool>("chargedOnly")),
56  pdgId_(cfg.getParameter<std::vector<int> >("pdgId")) { }
57 
58 
60  selected_.clear();
62  event.getByLabel(edm::InputTag("offlineBeamSpot"), beamSpot);
63  for( TrackingParticleCollection::const_iterator itp = c->begin();
64  itp != c->end(); ++ itp )
65  if ( operator()(*itp,beamSpot.product(),event,setup) ) {
66  selected_.push_back( & * itp );
67  }
68  }
69 
70  const_iterator begin() const { return selected_.begin(); }
71  const_iterator end() const { return selected_.end(); }
72 
73  // Operator() performs the selection: e.g. if (tPSelector(tp, bs, event, evtsetup)) {...
74  bool operator()( const TrackingParticle & tp, const reco::BeamSpot* bs, const edm::Event &iEvent, const edm::EventSetup& iSetup ) const {
75  if (chargedOnly_ && tp.charge()==0) return false;//select only if charge!=0
76  //bool testId = false;
77  //unsigned int idSize = pdgId_.size();
78  //if (idSize==0) testId = true;
79  //else for (unsigned int it=0;it!=idSize;++it){
80  //if (tp.pdgId()==pdgId_[it]) testId = true;
81  //}
82 
84  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
85 
87  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
88 
89  GlobalVector finalGV (0,0,0);
90  GlobalPoint finalGP(0,0,0);
91  GlobalVector momentum(0,0,0);//At the PCA
92  GlobalPoint vertex(0,0,0);//At the PCA
93  double radius(9999);
94  bool found(0);
95 
96  const std::vector<PSimHit> & simHits = tp.trackPSimHit(DetId::Tracker);
97  for(std::vector<PSimHit>::const_iterator it=simHits.begin(); it!=simHits.end(); ++it){
98  const GeomDet* tmpDet = tracker->idToDet( DetId(it->detUnitId()) ) ;
99  LocalVector lv = it->momentumAtEntry();
100  Local3DPoint lp = it->localPosition ();
101  GlobalVector gv = tmpDet->surface().toGlobal( lv );
102  GlobalPoint gp = tmpDet->surface().toGlobal( lp );
103  if(gp.perp()<radius){
104  found=true;
105  radius = gp.perp();
106  finalGV = gv;
107  finalGP = gp;
108  }
109  }
110  if(!found) return 0;
111  else
112  {
113  FreeTrajectoryState ftsAtProduction(finalGP,finalGV,TrackCharge(tp.charge()),theMF.product());
114  TSCBLBuilderNoMaterial tscblBuilder;
115  TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,*bs);//as in TrackProducerAlgorithm
116  if(!tsAtClosestApproach.isValid()){
117  std::cout << "WARNING: tsAtClosestApproach is not valid" << std::endl;
118  return 0;
119  }
120  else
121  {
122  momentum = tsAtClosestApproach.trackStateAtPCA().momentum();
123  vertex = tsAtClosestApproach.trackStateAtPCA().position();
124  return (
125  tp.matchedHit() >= minHit_ &&
126  sqrt(momentum.perp2()) >= ptMin_ &&
127  momentum.eta() >= minRapidity_ && momentum.eta() <= maxRapidity_ &&
128  sqrt(vertex.perp2()) <= tip_ &&
129  fabs(vertex.z()) <= lip_
130  );
131  }
132  }
133  }
134 
135  size_t size() const { return selected_.size(); }
136 
137  private:
138 
139  double ptMin_;
140  double minRapidity_;
141  double maxRapidity_;
142  double tip_;
143  double lip_;
144  int minHit_;
146  std::vector<int> pdgId_;
148 
149 
150 };
151 
152 #endif
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:78
T perp() const
Definition: PV3DBase.h:71
int charge() const
electric charge
Definition: ParticleBase.h:54
std::vector< TrackingParticle > TrackingParticleCollection
void select(const edm::Handle< collection > &c, const edm::Event &event, const edm::EventSetup &setup)
const std::vector< PSimHit > & trackPSimHit() const
T perp2() const
Definition: PV3DBase.h:70
int matchedHit() const
CosmicTrackingParticleSelector(double ptMin, double minRapidity, double maxRapidity, double tip, double lip, int minHit, bool chargedOnly, std::vector< int > pdgId=std::vector< int >())
bool operator()(const TrackingParticle &tp, const reco::BeamSpot *bs, const edm::Event &iEvent, const edm::EventSetup &iSetup) const
CosmicTrackingParticleSelector(const edm::ParameterSet &cfg)
int TrackCharge
Definition: TrackCharge.h:4
int iEvent
Definition: GenABIO.cc:243
T sqrt(T t)
Definition: SSEVec.h:46
T z() const
Definition: PV3DBase.h:63
std::vector< const TrackingParticle * > container
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
GlobalVector momentum() const
Definition: DetId.h:20
GlobalPoint position() const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
tuple simHits
Definition: trackerHits.py:16
T const * product() const
Definition: Handle.h:74
T eta() const
Definition: PV3DBase.h:75
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
tuple cout
Definition: gather_cfg.py:121
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")