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  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
20 
21 private:
22  void beginRun(const edm::Run&, const edm::EventSetup&) override;
23  void endRun(const edm::Run&, const edm::EventSetup&) override;
24 
26  void produce(edm::Event&, const edm::EventSetup&) override;
27 
31  std::vector<edm::EDGetTokenT<reco::NuclearInteractionCollection>> nuclearContainers_;
32 
34 };
35 
38 
41  // cut on the likelihood of the nuclear interaction
42  desc.add<double>("likelihoodCut", 0.1);
43  desc.add<std::vector<edm::InputTag>>("nuclearColList", {edm::InputTag("firstnuclearInteractionMaker")});
44  descriptions.add("pfNuclear", desc);
45 }
46 
47 using namespace std;
48 using namespace edm;
50  : pfTransformer_(nullptr), magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()) {
51  produces<reco::PFRecTrackCollection>();
52  produces<reco::PFNuclearInteractionCollection>();
53 
54  std::vector<edm::InputTag> tags = iConfig.getParameter<vector<InputTag>>("nuclearColList");
55 
56  for (unsigned int i = 0; i < tags.size(); ++i)
57  nuclearContainers_.push_back(consumes<reco::NuclearInteractionCollection>(tags[i]));
58 
59  likelihoodCut_ = iConfig.getParameter<double>("likelihoodCut");
60 }
61 
63 
65  typedef reco::NuclearInteraction::trackRef_iterator trackRef_iterator;
66 
67  //create the empty collections
68  auto pfNuclearColl = std::make_unique<reco::PFNuclearInteractionCollection>();
69  auto pfNuclearRecTrackColl = std::make_unique<reco::PFRecTrackCollection>();
70 
71  reco::PFRecTrackRefProd pfTrackRefProd = iEvent.getRefBeforePut<reco::PFRecTrackCollection>();
72  int hid = 0;
73 
74  // loop on the nuclear interaction collections
75  for (unsigned int istr = 0; istr < nuclearContainers_.size(); istr++) {
77  iEvent.getByToken(nuclearContainers_[istr], nuclCollH);
78  const reco::NuclearInteractionCollection& nuclColl = *(nuclCollH.product());
79 
80  // loop on all NuclearInteraction
81  for (unsigned int icoll = 0; icoll < nuclColl.size(); icoll++) {
82  if (nuclColl[icoll].likelihood() < likelihoodCut_)
83  continue;
84 
85  reco::PFRecTrackRefVector pfRecTkcoll;
86 
87  // convert the secondary tracks
88  for (trackRef_iterator it = nuclColl[icoll].secondaryTracks_begin(); it != nuclColl[icoll].secondaryTracks_end();
89  it++) {
90  reco::PFRecTrack pftrack(
91  (*it)->charge(), reco::PFRecTrack::KF, it->key(), (reco::TrackRef)((*it).castTo<reco::TrackRef>()));
92  Trajectory FakeTraj;
93  bool valid = pfTransformer_->addPoints(pftrack, **it, FakeTraj);
94  if (valid) {
95  pfRecTkcoll.push_back(reco::PFRecTrackRef(pfTrackRefProd, hid++));
96  pfNuclearRecTrackColl->push_back(pftrack);
97  }
98  }
99  reco::NuclearInteractionRef niRef(nuclCollH, icoll);
100  pfNuclearColl->push_back(reco::PFNuclearInteraction(niRef, pfRecTkcoll));
101  }
102  }
103  iEvent.put(std::move(pfNuclearRecTrackColl));
104  iEvent.put(std::move(pfNuclearColl));
105 }
106 
107 // ------------ method called once each job just before starting event loop ------------
108 void PFNuclearProducer::beginRun(const edm::Run& run, const EventSetup& iSetup) {
109  auto const& magneticField = &iSetup.getData(magneticFieldToken_);
112 }
113 
114 // ------------ method called once each job just after ending the event loop ------------
115 void PFNuclearProducer::endRun(const edm::Run& run, const EventSetup& iSetup) {
116  delete pfTransformer_;
117  pfTransformer_ = nullptr;
118 }
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void likelihood(int &npar, double *grad, double &fval, double *xval, int flag)
T const * product() const
Definition: Handle.h:70
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
std::vector< edm::EDGetTokenT< reco::NuclearInteractionCollection > > nuclearContainers_
int iEvent
Definition: GenABIO.cc:224
std::vector< NuclearInteraction > NuclearInteractionCollection
collection of NuclearInteractions
reco::Vertex::trackRef_iterator trackRef_iterator
~PFNuclearProducer() override
Destructor.
Transition
Definition: Transition.h:12
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
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.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void endRun(const edm::Run &, const edm::EventSetup &) override
HLT enums.
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
PFNuclearProducer(const edm::ParameterSet &)
Constructor.
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
void beginRun(const edm::Run &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int icoll
Definition: AMPTWrapper.h:146