CMS 3D CMS Logo

PFNuclearProducer.cc
Go to the documentation of this file.
1 #include <memory>
6 
7 using namespace std;
8 using namespace edm;
10  : pfTransformer_(nullptr), magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()) {
11  produces<reco::PFRecTrackCollection>();
12  produces<reco::PFNuclearInteractionCollection>();
13 
14  std::vector<edm::InputTag> tags = iConfig.getParameter<vector<InputTag> >("nuclearColList");
15 
16  for (unsigned int i = 0; i < tags.size(); ++i)
17  nuclearContainers_.push_back(consumes<reco::NuclearInteractionCollection>(tags[i]));
18 
19  likelihoodCut_ = iConfig.getParameter<double>("likelihoodCut");
20 }
21 
23 
25  typedef reco::NuclearInteraction::trackRef_iterator trackRef_iterator;
26 
27  //create the empty collections
28  auto pfNuclearColl = std::make_unique<reco::PFNuclearInteractionCollection>();
29  auto pfNuclearRecTrackColl = std::make_unique<reco::PFRecTrackCollection>();
30 
31  reco::PFRecTrackRefProd pfTrackRefProd = iEvent.getRefBeforePut<reco::PFRecTrackCollection>();
32  int hid = 0;
33 
34  // loop on the nuclear interaction collections
35  for (unsigned int istr = 0; istr < nuclearContainers_.size(); istr++) {
37  iEvent.getByToken(nuclearContainers_[istr], nuclCollH);
38  const reco::NuclearInteractionCollection& nuclColl = *(nuclCollH.product());
39 
40  // loop on all NuclearInteraction
41  for (unsigned int icoll = 0; icoll < nuclColl.size(); icoll++) {
42  if (nuclColl[icoll].likelihood() < likelihoodCut_)
43  continue;
44 
45  reco::PFRecTrackRefVector pfRecTkcoll;
46 
47  // convert the secondary tracks
48  for (trackRef_iterator it = nuclColl[icoll].secondaryTracks_begin(); it != nuclColl[icoll].secondaryTracks_end();
49  it++) {
50  reco::PFRecTrack pftrack(
51  (*it)->charge(), reco::PFRecTrack::KF, it->key(), (reco::TrackRef)((*it).castTo<reco::TrackRef>()));
52  Trajectory FakeTraj;
53  bool valid = pfTransformer_->addPoints(pftrack, **it, FakeTraj);
54  if (valid) {
55  pfRecTkcoll.push_back(reco::PFRecTrackRef(pfTrackRefProd, hid++));
56  pfNuclearRecTrackColl->push_back(pftrack);
57  }
58  }
59  reco::NuclearInteractionRef niRef(nuclCollH, icoll);
60  pfNuclearColl->push_back(reco::PFNuclearInteraction(niRef, pfRecTkcoll));
61  }
62  }
63  iEvent.put(std::move(pfNuclearRecTrackColl));
64  iEvent.put(std::move(pfNuclearColl));
65 }
66 
67 // ------------ method called once each job just before starting event loop ------------
68 void PFNuclearProducer::beginRun(const edm::Run& run, const EventSetup& iSetup) {
69  auto const& magneticField = &iSetup.getData(magneticFieldToken_);
72 }
73 
74 // ------------ method called once each job just after ending the event loop ------------
75 void PFNuclearProducer::endRun(const edm::Run& run, const EventSetup& iSetup) {
76  delete pfTransformer_;
77  pfTransformer_ = nullptr;
78 }
edm::RefProd
Definition: EDProductfwd.h:25
mps_fire.i
i
Definition: mps_fire.py:428
edm::Handle::product
T const * product() const
Definition: Handle.h:70
ESHandle.h
edm::Run
Definition: Run.h:45
edm
HLT enums.
Definition: AlignableModifier.h:19
PFNuclearProducer::magneticFieldToken_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
Definition: PFNuclearProducer.h:34
reco::PFNuclearInteraction
Definition: PFNuclearInteraction.h:15
reco::PFRecTrack::KF
Definition: PFRecTrack.h:25
edm::RefVector< PFRecTrackCollection >
HLT_FULL_cff.magneticField
magneticField
Definition: HLT_FULL_cff.py:348
PFTrackTransformer.h
edm::Handle
Definition: AssociativeIterator.h:50
likelihood
void likelihood(int &npar, double *grad, double &fval, double *xval, int flag)
Definition: MuScleFitUtils.cc:1784
PFNuclearProducer.h
edm::Ref< TrackCollection >
PFTrackTransformer::addPoints
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.
Definition: PFTrackTransformer.cc:40
reco::NuclearInteractionCollection
std::vector< NuclearInteraction > NuclearInteractionCollection
collection of NuclearInteractions
Definition: NuclearInteractionFwd.h:8
PFNuclearProducer::nuclearContainers_
std::vector< edm::EDGetTokenT< reco::NuclearInteractionCollection > > nuclearContainers_
Definition: PFNuclearProducer.h:32
PFNuclearProducer::~PFNuclearProducer
~PFNuclearProducer() override
Destructor.
Definition: PFNuclearProducer.cc:22
Point3DBase< float, GlobalTag >
reco::TrackRef
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
edm::ParameterSet
Definition: ParameterSet.h:47
edm::Transition
Transition
Definition: Transition.h:12
reco::NuclearInteraction::trackRef_iterator
reco::Vertex::trackRef_iterator trackRef_iterator
Definition: NuclearInteraction.h:15
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:58
edm::RefVector::push_back
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
Trajectory.h
PFNuclearProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
Definition: PFNuclearProducer.cc:24
PFNuclearProducer::likelihoodCut_
double likelihoodCut_
Definition: PFNuclearProducer.h:31
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
reco::PFRecTrack
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
Trajectory
Definition: Trajectory.h:38
PFNuclearProducer::PFNuclearProducer
PFNuclearProducer(const edm::ParameterSet &)
Constructor.
Definition: PFNuclearProducer.cc:9
PFTrackTransformer::OnlyProp
void OnlyProp()
Definition: PFTrackTransformer.h:54
edm::Transition::BeginRun
PFNuclearProducer::beginRun
void beginRun(const edm::Run &, const edm::EventSetup &) override
Definition: PFNuclearProducer.cc:68
PFNuclearProducer::endRun
void endRun(const edm::Run &, const edm::EventSetup &) override
Definition: PFNuclearProducer.cc:75
icoll
int icoll
Definition: AMPTWrapper.h:146
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
triggerMatcherToHLTDebug_cfi.tags
tags
Definition: triggerMatcherToHLTDebug_cfi.py:9
reco::PFRecTrackCollection
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
RunInfoPI::valid
Definition: RunInfoPayloadInspectoHelper.h:16
PFTrackTransformer
Definition: PFTrackTransformer.h:34
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
PFNuclearProducer::pfTransformer_
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
Definition: PFNuclearProducer.h:30
edm::Event
Definition: Event.h:73