CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
PFConversionProducer Class Reference
Inheritance diagram for PFConversionProducer:
edm::stream::EDProducer<>

Public Member Functions

 PFConversionProducer (const edm::ParameterSet &)
 Constructor. More...
 
 ~PFConversionProducer () override
 Destructor. More...
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

void beginRun (const edm::Run &, const edm::EventSetup &) override
 
void endRun (const edm::Run &, const edm::EventSetup &) override
 
void produce (edm::Event &, const edm::EventSetup &) override
 Produce the PFRecTrack collection. More...
 

Private Attributes

const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagneticFieldToken_
 
edm::EDGetTokenT< reco::ConversionCollectionpfConversionContainer_
 
PFTrackTransformerpfTransformer_
 PFTrackTransformer. More...
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordtransientTrackToken_
 
edm::EDGetTokenT< reco::VertexCollectionvtx_h
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 17 of file PFConversionProducer.cc.

Constructor & Destructor Documentation

◆ PFConversionProducer()

PFConversionProducer::PFConversionProducer ( const edm::ParameterSet iConfig)
explicit

Constructor.

Definition at line 48 of file PFConversionProducer.cc.

References edm::BeginRun, edm::ParameterSet::getParameter(), pfConversionContainer_, and vtx_h.

49  : pfTransformer_(nullptr),
50  transientTrackToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
51  magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()) {
52  produces<reco::PFRecTrackCollection>();
53  produces<reco::PFConversionCollection>();
54 
55  pfConversionContainer_ = consumes<reco::ConversionCollection>(iConfig.getParameter<InputTag>("conversionCollection"));
56 
57  vtx_h = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("PrimaryVertexLabel"));
58 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< reco::ConversionCollection > pfConversionContainer_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
edm::EDGetTokenT< reco::VertexCollection > vtx_h
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackToken_
PFTrackTransformer * pfTransformer_
PFTrackTransformer.

◆ ~PFConversionProducer()

PFConversionProducer::~PFConversionProducer ( )
override

Destructor.

Definition at line 60 of file PFConversionProducer.cc.

References pfTransformer_.

60 { delete pfTransformer_; }
PFTrackTransformer * pfTransformer_
PFTrackTransformer.

Member Function Documentation

◆ beginRun()

void PFConversionProducer::beginRun ( const edm::Run run,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 197 of file PFConversionProducer.cc.

References edm::EventSetup::getData(), HLT_2022v15_cff::magneticField, magneticFieldToken_, PFTrackTransformer::OnlyProp(), and pfTransformer_.

197  {
198  auto const& magneticField = &iSetup.getData(magneticFieldToken_);
201 }
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
PFTrackTransformer * pfTransformer_
PFTrackTransformer.

◆ endRun()

void PFConversionProducer::endRun ( const edm::Run run,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 204 of file PFConversionProducer.cc.

References pfTransformer_.

204  {
205  delete pfTransformer_;
206  pfTransformer_ = nullptr;
207 }
PFTrackTransformer * pfTransformer_
PFTrackTransformer.

◆ produce()

void PFConversionProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Produce the PFRecTrack collection.

Definition at line 62 of file PFConversionProducer.cc.

References PFTrackTransformer::addPoints(), reco::Conversion::arbitratedMergedEcalGeneral, hltPixelTracks_cff::chi2, ChiSquaredProbability(), MillePedeFileConverter_cfg::e, reco::PFTrajectoryPoint::ECALEntrance, DivergingColor::frac, edm::EventSetup::getData(), reco::Conversion::highPurity, heavyIonCSV_trainingSettings::idx, iEvent, reco::PFTrajectoryPoint::isValid(), edm::Ref< C, T, F >::key(), reco::PFRecTrack::KF, eostools::move(), ndof, AlCaHLTBitMon_ParallelJobs::p, pfConversionContainer_, pfTransformer_, MetAnalyzer::pv(), quality, edm::second(), IPTools::signedTransverseImpactParameter(), TrackingRecHit::some, transientTrackToken_, validateGeometry_cfg::valid, bphysicsOniaDQM_cfi::vertex, HltBtagValidation_cff::Vertex, and vtx_h.

62  {
63  //create the empty collections
64  auto pfConversionColl = std::make_unique<reco::PFConversionCollection>();
65  auto pfRecTrackColl = std::make_unique<reco::PFRecTrackCollection>();
66 
67  TransientTrackBuilder const& thebuilder = iSetup.getData(transientTrackToken_);
68 
69  reco::PFRecTrackRefProd pfTrackRefProd = iEvent.getRefBeforePut<reco::PFRecTrackCollection>();
71  iEvent.getByToken(pfConversionContainer_, convCollH);
72 
73  const reco::ConversionCollection& convColl = *(convCollH.product());
74 
76  iEvent.getByToken(vtx_h, vertex);
77  //Find PV for IP calculation, if there is no PV in collection than use dummy
79  const reco::Vertex* pv = &dummy;
80  if (vertex.isValid()) {
81  pv = &*vertex->begin();
82  } else { // create a dummy PV
84  e(0, 0) = 0.0015 * 0.0015;
85  e(1, 1) = 0.0015 * 0.0015;
86  e(2, 2) = 15. * 15.;
87  reco::Vertex::Point p(0, 0, 0);
88  dummy = reco::Vertex(p, e, 0, 0, 0);
89  }
90 
91  int idx = 0; //index of track in PFRecTrack collection
92  multimap<unsigned int, unsigned int> trackmap; //Map of Collections and tracks
93  std::vector<unsigned int> conv_coll(0);
94 
95  // CLEAN CONVERSION COLLECTION FOR DUPLICATES
96  for (unsigned int icoll1 = 0; icoll1 < convColl.size(); icoll1++) {
97  if ((!convColl[icoll1].quality(reco::Conversion::arbitratedMergedEcalGeneral)) ||
98  (!convColl[icoll1].quality(reco::Conversion::highPurity)))
99  continue;
100 
101  bool greater_prob = false;
102  std::vector<edm::RefToBase<reco::Track> > tracksRefColl1 = convColl[icoll1].tracks();
103  for (unsigned it1 = 0; it1 < tracksRefColl1.size(); it1++) {
104  reco::TrackRef trackRef1 = (tracksRefColl1[it1]).castTo<reco::TrackRef>();
105 
106  for (unsigned int icoll2 = 0; icoll2 < convColl.size(); icoll2++) {
107  if (icoll1 == icoll2)
108  continue;
109  if ((!convColl[icoll2].quality(reco::Conversion::arbitratedMergedEcalGeneral)) ||
110  (!convColl[icoll2].quality(reco::Conversion::highPurity)))
111  continue;
112  std::vector<edm::RefToBase<reco::Track> > tracksRefColl2 = convColl[icoll2].tracks();
113  for (unsigned it2 = 0; it2 < tracksRefColl2.size(); it2++) {
114  reco::TrackRef trackRef2 = (tracksRefColl2[it2]).castTo<reco::TrackRef>();
115  double like1 = -999;
116  double like2 = -999;
117  //number of shared hits
118  int shared = 0;
119  for (auto const& hit1 : trackRef1->recHits())
120  if (hit1->isValid()) {
121  //count number of shared hits
122  for (auto const& hit2 : trackRef2->recHits()) {
123  if (hit2->isValid() && hit1->sharesInput(hit2, TrackingRecHit::some))
124  shared++;
125  }
126  }
127  float frac = 0;
128  //number of valid hits in tracks that are duplicates
129  float size1 = trackRef1->found();
130  float size2 = trackRef2->found();
131  //divide number of shared hits by the total number of hits for the track with less hits
132  if (size1 > size2)
133  frac = (double)shared / size2;
134  else
135  frac = (double)shared / size1;
136  if (frac > 0.9) {
137  like1 = ChiSquaredProbability(convColl[icoll1].conversionVertex().chi2(),
138  convColl[icoll1].conversionVertex().ndof());
139  like2 = ChiSquaredProbability(convColl[icoll2].conversionVertex().chi2(),
140  convColl[icoll2].conversionVertex().ndof());
141  }
142  if (like2 > like1) {
143  greater_prob = true;
144  break;
145  }
146  } //end loop over tracks in collection 2
147 
148  if (greater_prob)
149  break; //if a duplicate track is found in a collection with greater Chi^2 probability for Vertex fit then break out of comparison loop
150  } //end loop over collection 2 checking
151  if (greater_prob)
152  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
153  } //end loop over tracks in collection 1
154  if (!greater_prob)
155  conv_coll.push_back(icoll1);
156  } //end loop over collection 1
157 
158  //Finally fill empty collections
159  for (unsigned iColl = 0; iColl < conv_coll.size(); iColl++) {
160  unsigned int collindex = conv_coll[iColl];
161  //std::cout<<"Filling this collection"<<collindex<<endl;
162  std::vector<reco::PFRecTrackRef> pfRecTkcoll;
163 
164  std::vector<edm::RefToBase<reco::Track> > tracksRefColl = convColl[collindex].tracks();
165  // convert the secondary tracks
166  for (unsigned it = 0; it < tracksRefColl.size(); it++) {
167  reco::TrackRef trackRef = (tracksRefColl[it]).castTo<reco::TrackRef>();
168  reco::PFRecTrack pfRecTrack(trackRef->charge(), reco::PFRecTrack::KF, trackRef.key(), trackRef);
169  //std::cout<<"Track Pt "<<trackRef->pt()<<std::endl;
170  Trajectory FakeTraj;
171  bool valid = pfTransformer_->addPoints(pfRecTrack, *trackRef, FakeTraj);
172  if (valid) {
173  double stip = -999;
174  const reco::PFTrajectoryPoint& atECAL = pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance);
175  //if extrapolation to ECAL is valid then calculate STIP
176  if (atECAL.isValid()) {
177  GlobalVector direction(pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().x(),
178  pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().y(),
179  pfRecTrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().z());
180  stip = IPTools::signedTransverseImpactParameter(thebuilder.build(*trackRef), direction, *pv)
181  .second.significance();
182  }
183  pfRecTrack.setSTIP(stip);
184  pfRecTkcoll.push_back(reco::PFRecTrackRef(pfTrackRefProd, idx++));
185  pfRecTrackColl->push_back(pfRecTrack);
186  }
187  } //end loop over tracks
188  //store reference to the Conversion collection
189  reco::ConversionRef niRef(convCollH, collindex);
190  pfConversionColl->push_back(reco::PFConversion(niRef, pfRecTkcoll));
191  } //end loop over collections
192  iEvent.put(std::move(pfRecTrackColl));
193  iEvent.put(std::move(pfConversionColl));
194 }
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
edm::EDGetTokenT< reco::ConversionCollection > pfConversionContainer_
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
key_type key() const
Accessor for product key.
Definition: Ref.h:250
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
U second(std::pair< T, U > const &p)
bool isValid() const
is this point valid ?
int iEvent
Definition: GenABIO.cc:224
def pv(vc)
Definition: MetAnalyzer.py:7
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
float ChiSquaredProbability(double chiSquared, double nrDOF)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
edm::EDGetTokenT< reco::VertexCollection > vtx_h
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackToken_
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.
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
string quality
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ magneticFieldToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> PFConversionProducer::magneticFieldToken_
private

Definition at line 38 of file PFConversionProducer.cc.

Referenced by beginRun().

◆ pfConversionContainer_

edm::EDGetTokenT<reco::ConversionCollection> PFConversionProducer::pfConversionContainer_
private

Definition at line 34 of file PFConversionProducer.cc.

Referenced by PFConversionProducer(), and produce().

◆ pfTransformer_

PFTrackTransformer* PFConversionProducer::pfTransformer_
private

◆ transientTrackToken_

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> PFConversionProducer::transientTrackToken_
private

Definition at line 37 of file PFConversionProducer.cc.

Referenced by produce().

◆ vtx_h

edm::EDGetTokenT<reco::VertexCollection> PFConversionProducer::vtx_h
private

Definition at line 35 of file PFConversionProducer.cc.

Referenced by PFConversionProducer(), and produce().