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