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 
21 
23  pfTransformer_(0)
24 {
25  produces<reco::PFRecTrackCollection>();
26  produces<reco::PFConversionCollection>();
27 
28  pfConversionContainer_ =consumes<reco::ConversionCollection>(iConfig.getParameter< InputTag >("conversionCollection"));
29 
30  vtx_h=consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("PrimaryVertexLabel"));
31 }
32 
34 {
35  delete pfTransformer_;
36 }
37 
38 void
40 {
41 
42  //create the empty collections
43  auto_ptr< reco::PFConversionCollection >
44  pfConversionColl (new reco::PFConversionCollection);
45  auto_ptr< reco::PFRecTrackCollection >
46  pfRecTrackColl (new reco::PFRecTrackCollection);
47 
49  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
50  TransientTrackBuilder thebuilder = *(builder.product());
53  iEvent.getByToken(pfConversionContainer_, convCollH);
54 
55  const reco::ConversionCollection& convColl = *(convCollH.product());
56 
58  iEvent.getByToken(vtx_h, vertex);
59  //Find PV for IP calculation, if there is no PV in collection than use dummy
60  reco::Vertex dummy;
61  const reco::Vertex* pv=&dummy;
62  if (vertex.isValid())
63  {
64  pv = &*vertex->begin();
65  }
66  else
67  { // create a dummy PV
69  e(0, 0) = 0.0015 * 0.0015;
70  e(1, 1) = 0.0015 * 0.0015;
71  e(2, 2) = 15. * 15.;
72  reco::Vertex::Point p(0, 0, 0);
73  dummy = reco::Vertex(p, e, 0, 0, 0);
74  }
75 
76  int idx = 0; //index of track in PFRecTrack collection
77  multimap<unsigned int, unsigned int> trackmap; //Map of Collections and tracks
78  std::vector<unsigned int> conv_coll(0);
79 
80  // CLEAN CONVERSION COLLECTION FOR DUPLICATES
81  for( unsigned int icoll1=0; icoll1 < convColl.size(); icoll1++)
82  {
83  if (( !convColl[icoll1].quality(reco::Conversion::arbitratedMergedEcalGeneral)) || (!convColl[icoll1].quality(reco::Conversion::highPurity))) continue;
84 
85  bool greater_prob=false;
86  std::vector<edm::RefToBase<reco::Track> > tracksRefColl1 = convColl[icoll1].tracks();
87  for(unsigned it1 = 0; it1 < tracksRefColl1.size(); it1++)
88  {
89  reco::TrackRef trackRef1 = (tracksRefColl1[it1]).castTo<reco::TrackRef>();
90 
91  for( unsigned int icoll2=0; icoll2 < convColl.size(); icoll2++)
92  {
93  if(icoll1==icoll2)continue;
94  if (( !convColl[icoll2].quality(reco::Conversion::arbitratedMergedEcalGeneral)) || (!convColl[icoll2].quality(reco::Conversion::highPurity))) continue;
95  std::vector<edm::RefToBase<reco::Track> > tracksRefColl2 = convColl[icoll2].tracks();
96  for(unsigned it2 = 0; it2 < tracksRefColl2.size(); it2++)
97  {
98  reco::TrackRef trackRef2 = (tracksRefColl2[it2]).castTo<reco::TrackRef>();
99  double like1=-999;
100  double like2=-999;
101  //number of shared hits
102  int shared=0;
103  for(trackingRecHit_iterator iHit1=trackRef1->recHitsBegin(); iHit1!=trackRef1->recHitsEnd(); iHit1++)
104  {
105  const TrackingRecHit *h_1=(*iHit1);
106  if(h_1->isValid()){
107  for(trackingRecHit_iterator iHit2=trackRef2->recHitsBegin(); iHit2!=trackRef2->recHitsEnd(); iHit2++)
108  {
109  const TrackingRecHit *h_2=(*iHit2);
110  if(h_2->isValid() && h_1->sharesInput(h_2, TrackingRecHit::some))shared++;//count number of shared hits
111  }
112  }
113  }
114  float frac=0;
115  //number of valid hits in tracks that are duplicates
116  float size1=trackRef1->found();
117  float size2=trackRef2->found();
118  //divide number of shared hits by the total number of hits for the track with less hits
119  if(size1>size2)frac=(double)shared/size2;
120  else frac=(double)shared/size1;
121  if(frac>0.9)
122  {
123  like1=ChiSquaredProbability(convColl[icoll1].conversionVertex().chi2(), convColl[icoll1].conversionVertex().ndof());
124  like2=ChiSquaredProbability(convColl[icoll2].conversionVertex().chi2(), convColl[icoll2].conversionVertex().ndof());
125  }
126  if(like2>like1)
127  {greater_prob=true; break;}
128  }//end loop over tracks in collection 2
129 
130  if(greater_prob)break; //if a duplicate track is found in a collection with greater Chi^2 probability for Vertex fit then break out of comparison loop
131  }//end loop over collection 2 checking
132  if(greater_prob)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
133  } //end loop over tracks in collection 1
134  if(!greater_prob)conv_coll.push_back(icoll1);
135  }//end loop over collection 1
136 
137  //Finally fill empty collections
138  for(unsigned iColl=0; iColl<conv_coll.size(); iColl++)
139  {
140  unsigned int collindex=conv_coll[iColl];
141  //std::cout<<"Filling this collection"<<collindex<<endl;
142  std::vector<reco::PFRecTrackRef> pfRecTkcoll;
143 
144  std::vector<edm::RefToBase<reco::Track> > tracksRefColl = convColl[collindex].tracks();
145  // convert the secondary tracks
146  for(unsigned it = 0; it < tracksRefColl.size(); it++)
147  {
148  reco::TrackRef trackRef = (tracksRefColl[it]).castTo<reco::TrackRef>();
149  reco::PFRecTrack pfRecTrack( trackRef->charge(),
151  trackRef.key(),
152  trackRef );
153  //std::cout<<"Track Pt "<<trackRef->pt()<<std::endl;
154  Trajectory FakeTraj;
155  bool valid = pfTransformer_->addPoints( pfRecTrack, *trackRef, FakeTraj);
156  if(valid)
157  {
158  double stip=-999;
159  const reco::PFTrajectoryPoint& atECAL=pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance);
160  //if extrapolation to ECAL is valid then calculate STIP
161  if(atECAL.isValid())
162  {
163  GlobalVector direction(pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().x(),
164  pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().y(),
165  pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().z());
166  stip = IPTools::signedTransverseImpactParameter(thebuilder.build(*trackRef), direction, *pv).second.significance();
167  }
168  pfRecTrack.setSTIP(stip);
169  pfRecTkcoll.push_back(reco::PFRecTrackRef( pfTrackRefProd, idx++));
170  pfRecTrackColl->push_back(pfRecTrack);
171  }
172  }//end loop over tracks
173  //store reference to the Conversion collection
174  reco::ConversionRef niRef(convCollH, collindex);
175  pfConversionColl->push_back( reco::PFConversion( niRef, pfRecTkcoll ));
176  }//end loop over collections
177  iEvent.put(pfRecTrackColl);
178  iEvent.put(pfConversionColl);
179 }
180 
181 // ------------ method called once each job just before starting event loop ------------
182 void
184  const EventSetup& iSetup)
185 {
187  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
188  pfTransformer_= new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0,0,0))));
190 }
191 
192 // ------------ method called once each job just after ending the event loop ------------
193 void
195  const EventSetup& iSetup) {
196  delete pfTransformer_;
197  pfTransformer_=nullptr;
198 }
T getParameter(std::string const &) const
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
edm::EDGetTokenT< reco::ConversionCollection > pfConversionContainer_
virtual void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
virtual 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:50
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.
reco::TransientTrack build(const reco::Track *p) const
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
tuple magneticField
std::vector< PFConversion > PFConversionCollection
collection of PFConversion objects
key_type key() const
Accessor for product key.
Definition: Ref.h:264
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:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
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:75
RefProd< PROD > getRefBeforePut()
Definition: Event.h:141
bool isValid() const
is this point valid ?
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool isValid() const
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
~PFConversionProducer()
Destructor.
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
Definition: Run.h:43