CMS 3D CMS Logo

List of all members | Public Member Functions | Static 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
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

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 57 of file PFConversionProducer.cc.

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

58  : pfTransformer_(nullptr),
59  transientTrackToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
60  magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()) {
61  produces<reco::PFRecTrackCollection>();
62  produces<reco::PFConversionCollection>();
63 
64  pfConversionContainer_ = consumes<reco::ConversionCollection>(iConfig.getParameter<InputTag>("conversionCollection"));
65 
66  vtx_h = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("PrimaryVertexLabel"));
67 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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 69 of file PFConversionProducer.cc.

References pfTransformer_.

69 { delete pfTransformer_; }
PFTrackTransformer * pfTransformer_
PFTrackTransformer.

Member Function Documentation

◆ beginRun()

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

Definition at line 206 of file PFConversionProducer.cc.

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

206  {
207  auto const& magneticField = &iSetup.getData(magneticFieldToken_);
210 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
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 213 of file PFConversionProducer.cc.

References pfTransformer_.

213  {
214  delete pfTransformer_;
215  pfTransformer_ = nullptr;
216 }
PFTrackTransformer * pfTransformer_
PFTrackTransformer.

◆ fillDescriptions()

void PFConversionProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 46 of file PFConversionProducer.cc.

References edm::ConfigurationDescriptions::add(), and submitPVResolutionJobs::desc.

46  {
48  desc.add<edm::InputTag>("conversionCollection", {"allConversions", ""});
49  desc.add<edm::InputTag>("PrimaryVertexLabel", {"offlinePrimaryVertices"});
50  descriptions.add("pfConversions", desc);
51 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

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

Produce the PFRecTrack collection.

Definition at line 71 of file PFConversionProducer.cc.

References PFTrackTransformer::addPoints(), reco::Conversion::arbitratedMergedEcalGeneral, isoTrack_cff::chi2, ChiSquaredProbability(), MillePedeFileConverter_cfg::e, reco::PFTrajectoryPoint::ECALEntrance, DivergingColor::frac, edm::EventSetup::getData(), reco::Conversion::highPurity, heavyIonCSV_trainingSettings::idx, iEvent, reco::PFTrajectoryPoint::isValid(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, edm::Ref< C, T, F >::key(), reco::PFRecTrack::KF, eostools::move(), PixelVertexMonitor_cff::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.

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

Referenced by beginRun().

◆ pfConversionContainer_

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

Definition at line 36 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 39 of file PFConversionProducer.cc.

Referenced by produce().

◆ vtx_h

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

Definition at line 37 of file PFConversionProducer.cc.

Referenced by PFConversionProducer(), and produce().