CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackVertexArbitrator.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <set>
3 
4 
9 
10 
16 
23 
24 
30 
36 
37 
39 
40 
41 
43  public:
45 
46 
47  virtual void produce(edm::Event &event, const edm::EventSetup &es);
48 
49  private:
50  bool trackFilter(const reco::TrackRef &track) const;
51 
56  double dRCut;
57  double distCut;
58  double sigCut;
59  double dLenFraction;
60 };
61 
63  primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
64  secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
65  trackCollection(params.getParameter<edm::InputTag>("tracks")),
66  beamSpotCollection(params.getParameter<edm::InputTag>("beamSpot")),
67  dRCut(params.getParameter<double>("dRCut")),
68  distCut(params.getParameter<double>("distCut")),
69  sigCut(params.getParameter<double>("sigCut")),
70  dLenFraction(params.getParameter<double>("dLenFraction"))
71 {
72  produces<reco::VertexCollection>();
73 }
74 
76 {
77  if (track->hitPattern().trackerLayersWithMeasurement() < 4)
78  return false;
79  if (track->pt() < 0.4 )
80  return false;
81 
82  return true;
83 }
84 
85 
87 {
88  using namespace reco;
89 
90  edm::Handle<VertexCollection> secondaryVertices;
91  event.getByLabel(secondaryVertexCollection, secondaryVertices);
92 
93  edm::Handle<VertexCollection> primaryVertices;
94  event.getByLabel(primaryVertexCollection, primaryVertices);
95 
97  event.getByLabel(trackCollection, tracks);
98 
100  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
101  trackBuilder);
102 
104  event.getByLabel(beamSpotCollection, beamSpot);
105 
106  const reco::Vertex &pv = (*primaryVertices)[0];
107 // std::cout << "PV: " << pv.position() << std::endl;
108  VertexDistance3D dist;
109 
110  double sigmacut = 3.0;
111  double Tini = 256.;
112  double ratio = 0.25;
113 
114  AdaptiveVertexFitter theAdaptiveFitter(
115  GeometricAnnealing(sigmacut, Tini, ratio),
120 
121 
122 
123  std::auto_ptr<VertexCollection> recoVertices(new VertexCollection);
124 
125  VertexDistance3D vdist;
126 
127 for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
128  sv != secondaryVertices->end(); ++sv) {
129 /* recoVertices->push_back(*sv);
130 
131 
132  for(std::vector<reco::Vertex>::iterator sv = recoVertices->begin();
133  sv != recoVertices->end(); ++sv) {
134 */
135  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
136  GlobalPoint ssv(sv->position().x(),sv->position().y(),sv->position().z());
137  GlobalVector flightDir = ssv-ppv;
138 // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl;
139  Measurement1D dlen= vdist.distance(pv,*sv);
140  std::vector<reco::TransientTrack> selTracks;
141 
142  for(TrackCollection::const_iterator track = tracks->begin();
143  track != tracks->end(); ++track) {
144 
145  TrackRef ref(tracks, track - tracks->begin());
146  if (!trackFilter(ref)) continue;
147 
148  TransientTrack tt = trackBuilder->build(ref);
149  tt.setBeamSpot(*beamSpot);
150  float w = sv->trackWeight(ref);
151  std::pair<bool,Measurement1D> ipv = IPTools::absoluteImpactParameter3D(tt,pv);
152  std::pair<bool,Measurement1D> isv = IPTools::absoluteImpactParameter3D(tt,*sv);
153  if( w > 0 || ( isv.second.significance() < sigCut && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction ) )
154  {
155  float dR = deltaR(flightDir.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
156 
157  if(isv.second.value() < ipv.second.value() && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction
158  && dR < dRCut )
159  {
160 // if(w > 0.5) std::cout << " = ";
161  // else std::cout << " + ";
162  selTracks.push_back(tt);
163  } else
164  {
165 // if(w > 0.5 && isv.second.value() > ipv.second.value() ) std::cout << " - ";
166  // else std::cout << " ";
167  //add also the tracks used in previous fitting that are still closer to Sv than Pv
168  if(w > 0.5 && isv.second.value() < ipv.second.value() && dR < dRCut) selTracks.push_back(tt);
169  }
170 
171  // std::cout << "t : " << track-tracks->begin() << " w: " << w
172  // << " svip: " << isv.second.significance() << " " << isv.second.value()
173  // << " pvip: " << ipv.second.significance() << " " << ipv.second.value() << " dr: " << dR << std::endl;
174 
175  }
176  }
177 
178  if(selTracks.size() >= 2)
179  {
180  TransientVertex singleFitVertex;
181  singleFitVertex = theAdaptiveFitter.vertex(selTracks,ssv);
182  if(singleFitVertex.isValid()) recoVertices->push_back(singleFitVertex);
183  }
184  }
185  event.put(recoVertices);
186 }
187 
TrackVertexArbitrator(const edm::ParameterSet &params)
void setBeamSpot(const reco::BeamSpot &beamSpot)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
virtual void produce(edm::Event &event, const edm::EventSetup &es)
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:138
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:140
edm::InputTag primaryVertexCollection
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
bool trackFilter(const reco::TrackRef &track) 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
edm::InputTag secondaryVertexCollection
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
const Track & track() const
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:55
double value() const
Definition: Measurement1D.h:28
bool isValid() const