CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFConversionProducer.cc
Go to the documentation of this file.
1 #include <memory>
16 
17 typedef std::multimap<unsigned, std::vector<unsigned> > BlockMap;
18 using namespace std;
19 using namespace edm;
20 
22  produces<reco::PFRecTrackCollection>();
23  produces<reco::PFConversionCollection>();
24 
25  pfConversionContainer_ = consumes<reco::ConversionCollection>(iConfig.getParameter<InputTag>("conversionCollection"));
26 
27  vtx_h = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("PrimaryVertexLabel"));
28 }
29 
31 
33  //create the empty collections
34  auto pfConversionColl = std::make_unique<reco::PFConversionCollection>();
35  auto pfRecTrackColl = std::make_unique<reco::PFRecTrackCollection>();
36 
38  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
39  TransientTrackBuilder thebuilder = *(builder.product());
42  iEvent.getByToken(pfConversionContainer_, convCollH);
43 
44  const reco::ConversionCollection& convColl = *(convCollH.product());
45 
47  iEvent.getByToken(vtx_h, vertex);
48  //Find PV for IP calculation, if there is no PV in collection than use dummy
50  const reco::Vertex* pv = &dummy;
51  if (vertex.isValid()) {
52  pv = &*vertex->begin();
53  } else { // create a dummy PV
55  e(0, 0) = 0.0015 * 0.0015;
56  e(1, 1) = 0.0015 * 0.0015;
57  e(2, 2) = 15. * 15.;
58  reco::Vertex::Point p(0, 0, 0);
59  dummy = reco::Vertex(p, e, 0, 0, 0);
60  }
61 
62  int idx = 0; //index of track in PFRecTrack collection
63  multimap<unsigned int, unsigned int> trackmap; //Map of Collections and tracks
64  std::vector<unsigned int> conv_coll(0);
65 
66  // CLEAN CONVERSION COLLECTION FOR DUPLICATES
67  for (unsigned int icoll1 = 0; icoll1 < convColl.size(); icoll1++) {
68  if ((!convColl[icoll1].quality(reco::Conversion::arbitratedMergedEcalGeneral)) ||
69  (!convColl[icoll1].quality(reco::Conversion::highPurity)))
70  continue;
71 
72  bool greater_prob = false;
73  std::vector<edm::RefToBase<reco::Track> > tracksRefColl1 = convColl[icoll1].tracks();
74  for (unsigned it1 = 0; it1 < tracksRefColl1.size(); it1++) {
75  reco::TrackRef trackRef1 = (tracksRefColl1[it1]).castTo<reco::TrackRef>();
76 
77  for (unsigned int icoll2 = 0; icoll2 < convColl.size(); icoll2++) {
78  if (icoll1 == icoll2)
79  continue;
80  if ((!convColl[icoll2].quality(reco::Conversion::arbitratedMergedEcalGeneral)) ||
81  (!convColl[icoll2].quality(reco::Conversion::highPurity)))
82  continue;
83  std::vector<edm::RefToBase<reco::Track> > tracksRefColl2 = convColl[icoll2].tracks();
84  for (unsigned it2 = 0; it2 < tracksRefColl2.size(); it2++) {
85  reco::TrackRef trackRef2 = (tracksRefColl2[it2]).castTo<reco::TrackRef>();
86  double like1 = -999;
87  double like2 = -999;
88  //number of shared hits
89  int shared = 0;
90  for (auto const& hit1 : trackRef1->recHits())
91  if (hit1->isValid()) {
92  //count number of shared hits
93  for (auto const& hit2 : trackRef2->recHits()) {
94  if (hit2->isValid() && hit1->sharesInput(hit2, TrackingRecHit::some))
95  shared++;
96  }
97  }
98  float frac = 0;
99  //number of valid hits in tracks that are duplicates
100  float size1 = trackRef1->found();
101  float size2 = trackRef2->found();
102  //divide number of shared hits by the total number of hits for the track with less hits
103  if (size1 > size2)
104  frac = (double)shared / size2;
105  else
106  frac = (double)shared / size1;
107  if (frac > 0.9) {
108  like1 = ChiSquaredProbability(convColl[icoll1].conversionVertex().chi2(),
109  convColl[icoll1].conversionVertex().ndof());
110  like2 = ChiSquaredProbability(convColl[icoll2].conversionVertex().chi2(),
111  convColl[icoll2].conversionVertex().ndof());
112  }
113  if (like2 > like1) {
114  greater_prob = true;
115  break;
116  }
117  } //end loop over tracks in collection 2
118 
119  if (greater_prob)
120  break; //if a duplicate track is found in a collection with greater Chi^2 probability for Vertex fit then break out of comparison loop
121  } //end loop over collection 2 checking
122  if (greater_prob)
123  break; //if a duplicate track is found in a collection with greater Chi^2 probability for Vertex fit then one does not need to check the other track the collection will not be stored
124  } //end loop over tracks in collection 1
125  if (!greater_prob)
126  conv_coll.push_back(icoll1);
127  } //end loop over collection 1
128 
129  //Finally fill empty collections
130  for (unsigned iColl = 0; iColl < conv_coll.size(); iColl++) {
131  unsigned int collindex = conv_coll[iColl];
132  //std::cout<<"Filling this collection"<<collindex<<endl;
133  std::vector<reco::PFRecTrackRef> pfRecTkcoll;
134 
135  std::vector<edm::RefToBase<reco::Track> > tracksRefColl = convColl[collindex].tracks();
136  // convert the secondary tracks
137  for (unsigned it = 0; it < tracksRefColl.size(); it++) {
138  reco::TrackRef trackRef = (tracksRefColl[it]).castTo<reco::TrackRef>();
139  reco::PFRecTrack pfRecTrack(trackRef->charge(), reco::PFRecTrack::KF, trackRef.key(), trackRef);
140  //std::cout<<"Track Pt "<<trackRef->pt()<<std::endl;
141  Trajectory FakeTraj;
142  bool valid = pfTransformer_->addPoints(pfRecTrack, *trackRef, FakeTraj);
143  if (valid) {
144  double stip = -999;
145  const reco::PFTrajectoryPoint& atECAL = pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance);
146  //if extrapolation to ECAL is valid then calculate STIP
147  if (atECAL.isValid()) {
148  GlobalVector direction(pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().x(),
149  pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().y(),
150  pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().z());
151  stip = IPTools::signedTransverseImpactParameter(thebuilder.build(*trackRef), direction, *pv)
152  .second.significance();
153  }
154  pfRecTrack.setSTIP(stip);
155  pfRecTkcoll.push_back(reco::PFRecTrackRef(pfTrackRefProd, idx++));
156  pfRecTrackColl->push_back(pfRecTrack);
157  }
158  } //end loop over tracks
159  //store reference to the Conversion collection
160  reco::ConversionRef niRef(convCollH, collindex);
161  pfConversionColl->push_back(reco::PFConversion(niRef, pfRecTkcoll));
162  } //end loop over collections
163  iEvent.put(std::move(pfRecTrackColl));
164  iEvent.put(std::move(pfConversionColl));
165 }
166 
167 // ------------ method called once each job just before starting event loop ------------
170  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
171  pfTransformer_ = new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0, 0, 0))));
173 }
174 
175 // ------------ method called once each job just after ending the event loop ------------
177  delete pfTransformer_;
178  pfTransformer_ = nullptr;
179 }
T getParameter(std::string const &) const
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
edm::EDGetTokenT< reco::ConversionCollection > pfConversionContainer_
~PFConversionProducer() override
Destructor.
void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void endRun(const edm::Run &, const edm::EventSetup &) override
std::multimap< unsigned, std::vector< unsigned > > BlockMap
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
bool addPoints(reco::PFRecTrack &pftrack, const reco::Track &track, const Trajectory &traj, bool msgwarning=true) const
Add points to a PFTrack. return false if a TSOS is invalid.
#define nullptr
reco::TransientTrack build(const reco::Track *p) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
key_type key() const
Accessor for product key.
Definition: Ref.h:250
PFConversionProducer(const edm::ParameterSet &)
Constructor.
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:224
def pv(vc)
Definition: MetAnalyzer.py:7
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
float ChiSquaredProbability(double chiSquared, double nrDOF)
edm::EDGetTokenT< reco::VertexCollection > vtx_h
bool isValid() const
Definition: HandleBase.h:70
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
RefProd< PROD > getRefBeforePut()
Definition: Event.h:156
bool isValid() const
is this point valid ?
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
HLT enums.
T get() const
Definition: EventSetup.h:73
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
void beginRun(const edm::Run &, const edm::EventSetup &) override
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45