CMS 3D CMS Logo

PFNuclearProducer.cc
Go to the documentation of this file.
10 
12 public:
14  explicit PFNuclearProducer(const edm::ParameterSet&);
15 
17  ~PFNuclearProducer() override;
18 
19 private:
20  void beginRun(const edm::Run&, const edm::EventSetup&) override;
21  void endRun(const edm::Run&, const edm::EventSetup&) override;
22 
24  void produce(edm::Event&, const edm::EventSetup&) override;
25 
29  std::vector<edm::EDGetTokenT<reco::NuclearInteractionCollection> > nuclearContainers_;
30 
32 };
33 
36 
37 using namespace std;
38 using namespace edm;
40  : pfTransformer_(nullptr), magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()) {
41  produces<reco::PFRecTrackCollection>();
42  produces<reco::PFNuclearInteractionCollection>();
43 
44  std::vector<edm::InputTag> tags = iConfig.getParameter<vector<InputTag> >("nuclearColList");
45 
46  for (unsigned int i = 0; i < tags.size(); ++i)
47  nuclearContainers_.push_back(consumes<reco::NuclearInteractionCollection>(tags[i]));
48 
49  likelihoodCut_ = iConfig.getParameter<double>("likelihoodCut");
50 }
51 
53 
55  typedef reco::NuclearInteraction::trackRef_iterator trackRef_iterator;
56 
57  //create the empty collections
58  auto pfNuclearColl = std::make_unique<reco::PFNuclearInteractionCollection>();
59  auto pfNuclearRecTrackColl = std::make_unique<reco::PFRecTrackCollection>();
60 
61  reco::PFRecTrackRefProd pfTrackRefProd = iEvent.getRefBeforePut<reco::PFRecTrackCollection>();
62  int hid = 0;
63 
64  // loop on the nuclear interaction collections
65  for (unsigned int istr = 0; istr < nuclearContainers_.size(); istr++) {
67  iEvent.getByToken(nuclearContainers_[istr], nuclCollH);
68  const reco::NuclearInteractionCollection& nuclColl = *(nuclCollH.product());
69 
70  // loop on all NuclearInteraction
71  for (unsigned int icoll = 0; icoll < nuclColl.size(); icoll++) {
72  if (nuclColl[icoll].likelihood() < likelihoodCut_)
73  continue;
74 
75  reco::PFRecTrackRefVector pfRecTkcoll;
76 
77  // convert the secondary tracks
78  for (trackRef_iterator it = nuclColl[icoll].secondaryTracks_begin(); it != nuclColl[icoll].secondaryTracks_end();
79  it++) {
80  reco::PFRecTrack pftrack(
81  (*it)->charge(), reco::PFRecTrack::KF, it->key(), (reco::TrackRef)((*it).castTo<reco::TrackRef>()));
82  Trajectory FakeTraj;
83  bool valid = pfTransformer_->addPoints(pftrack, **it, FakeTraj);
84  if (valid) {
85  pfRecTkcoll.push_back(reco::PFRecTrackRef(pfTrackRefProd, hid++));
86  pfNuclearRecTrackColl->push_back(pftrack);
87  }
88  }
89  reco::NuclearInteractionRef niRef(nuclCollH, icoll);
90  pfNuclearColl->push_back(reco::PFNuclearInteraction(niRef, pfRecTkcoll));
91  }
92  }
93  iEvent.put(std::move(pfNuclearRecTrackColl));
94  iEvent.put(std::move(pfNuclearColl));
95 }
96 
97 // ------------ method called once each job just before starting event loop ------------
98 void PFNuclearProducer::beginRun(const edm::Run& run, const EventSetup& iSetup) {
99  auto const& magneticField = &iSetup.getData(magneticFieldToken_);
102 }
103 
104 // ------------ method called once each job just after ending the event loop ------------
105 void PFNuclearProducer::endRun(const edm::Run& run, const EventSetup& iSetup) {
106  delete pfTransformer_;
107  pfTransformer_ = nullptr;
108 }
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.cc:31
reco::PFNuclearInteraction
Definition: PFNuclearInteraction.h:15
EDProducer.h
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
edm::Ref< TrackCollection >
MakerMacros.h
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
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
reco::NuclearInteractionCollection
std::vector< NuclearInteraction > NuclearInteractionCollection
collection of NuclearInteractions
Definition: NuclearInteractionFwd.h:8
IdealMagneticFieldRecord.h
PFNuclearProducer::nuclearContainers_
std::vector< edm::EDGetTokenT< reco::NuclearInteractionCollection > > nuclearContainers_
Definition: PFNuclearProducer.cc:29
PFNuclearProducer::~PFNuclearProducer
~PFNuclearProducer() override
Destructor.
Definition: PFNuclearProducer.cc:52
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
Event.h
PFNuclearInteraction.h
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:36
MagneticField.h
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::ESGetToken< MagneticField, IdealMagneticFieldRecord >
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:54
PFNuclearProducer::likelihoodCut_
double likelihoodCut_
Definition: PFNuclearProducer.cc:28
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
PFNuclearProducer
Definition: PFNuclearProducer.cc:11
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:39
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:98
PFNuclearProducer::endRun
void endRun(const edm::Run &, const edm::EventSetup &) override
Definition: PFNuclearProducer.cc:105
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
ParameterSet.h
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.cc:27
edm::Event
Definition: Event.h:73