CMS 3D CMS Logo

PFNuclearProducer.cc
Go to the documentation of this file.
1 #include <memory>
8 using namespace std;
9 using namespace edm;
10 PFNuclearProducer::PFNuclearProducer(const ParameterSet& iConfig) : pfTransformer_(nullptr) {
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 
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) {
70  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
71  pfTransformer_ = new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0, 0, 0))));
73 }
74 
75 // ------------ method called once each job just after ending the event loop ------------
76 void PFNuclearProducer::endRun(const edm::Run& run, const EventSetup& iSetup) {
77  delete pfTransformer_;
78  pfTransformer_ = nullptr;
79 }
T getParameter(std::string const &) const
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
void likelihood(int &npar, double *grad, double &fval, double *xval, int flag)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
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.
#define nullptr
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.
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
RefProd< PROD > getRefBeforePut()
Definition: Event.h:156
T const * product() const
Definition: Handle.h:69
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
void endRun(const edm::Run &, const edm::EventSetup &) override
HLT enums.
T get() const
Definition: EventSetup.h:73
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
int icoll
Definition: AMPTWrapper.h:146