CMS 3D CMS Logo

PFTICLProducer.cc
Go to the documentation of this file.
1 // This producer converts a list of TICLCandidates to a list of PFCandidates.
2 
5 
8 
10 
13 
15 
17 
19 public:
21  ~PFTICLProducer() override {}
22 
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25  void produce(edm::Event&, const edm::EventSetup&) override;
26 
27 private:
28  // parameters
29  const bool useMTDTiming_;
30  const bool useTimingAverage_;
33  // inputs
37  // For PFMuonAlgo
38  std::unique_ptr<PFMuonAlgo> pfmu_;
39 };
40 
42 
44  : useMTDTiming_(conf.getParameter<bool>("useMTDTiming")),
45  useTimingAverage_(conf.getParameter<bool>("useTimingAverage")),
46  timingQualityThreshold_(conf.getParameter<double>("timingQualityThreshold")),
47  energy_from_regression_(conf.getParameter<bool>("energyFromRegression")),
48  ticl_candidates_(consumes<edm::View<TICLCandidate>>(conf.getParameter<edm::InputTag>("ticlCandidateSrc"))),
49  srcTrackTime_(consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("trackTimeValueMap"))),
50  srcTrackTimeError_(consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("trackTimeErrorMap"))),
51  srcTrackTimeQuality_(consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("trackTimeQualityMap"))),
52  muons_(consumes<reco::MuonCollection>(conf.getParameter<edm::InputTag>("muonSrc"))),
53  pfmu_(std::make_unique<PFMuonAlgo>(conf.getParameterSet("pfMuonAlgoParameters"),
54  false)) { // postMuonCleaning = false
55  produces<reco::PFCandidateCollection>();
56 }
57 
60  desc.add<edm::InputTag>("ticlCandidateSrc", edm::InputTag("ticlTrackstersMerge"));
61  desc.add<edm::InputTag>("trackTimeValueMap", edm::InputTag("tofPID:t0"));
62  desc.add<edm::InputTag>("trackTimeErrorMap", edm::InputTag("tofPID:sigmat0"));
63  desc.add<edm::InputTag>("trackTimeQualityMap", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
64  desc.add<bool>("energyFromRegression", true);
65  desc.add<double>("timingQualityThreshold", 0.5);
66  desc.add<bool>("useMTDTiming", true);
67  desc.add<bool>("useTimingAverage", false);
68  // For PFMuonAlgo
69  desc.add<edm::InputTag>("muonSrc", edm::InputTag("muons1stStep"));
70  edm::ParameterSetDescription psd_PFMuonAlgo;
71  PFMuonAlgo::fillPSetDescription(psd_PFMuonAlgo);
72  desc.add<edm::ParameterSetDescription>("pfMuonAlgoParameters", psd_PFMuonAlgo);
73  //
74  descriptions.add("pfTICLProducer", desc);
75 }
76 
78  //get TICLCandidates
80  evt.getByToken(ticl_candidates_, ticl_cand_h);
81  const auto ticl_candidates = *ticl_cand_h;
82  edm::Handle<edm::ValueMap<float>> trackTimeH, trackTimeErrH, trackTimeQualH;
83  evt.getByToken(srcTrackTime_, trackTimeH);
84  evt.getByToken(srcTrackTimeError_, trackTimeErrH);
85  evt.getByToken(srcTrackTimeQuality_, trackTimeQualH);
86  const auto muonH = evt.getHandle(muons_);
87  const auto& muons = *muonH;
88 
89  auto candidates = std::make_unique<reco::PFCandidateCollection>();
90 
91  for (const auto& ticl_cand : ticl_candidates) {
92  const auto abs_pdg_id = std::abs(ticl_cand.pdgId());
93  const auto charge = ticl_cand.charge();
94  const auto& four_mom = ticl_cand.p4();
95  float total_raw_energy = 0.f;
96  float total_em_raw_energy = 0.f;
97  for (const auto& t : ticl_cand.tracksters()) {
98  total_raw_energy += t->raw_energy();
99  total_em_raw_energy += t->raw_em_energy();
100  }
101  float ecal_energy_fraction = total_em_raw_energy / total_raw_energy;
102  float ecal_energy = energy_from_regression_ ? ticl_cand.p4().energy() * ecal_energy_fraction
103  : ticl_cand.rawEnergy() * ecal_energy_fraction;
104  float hcal_energy =
105  energy_from_regression_ ? ticl_cand.p4().energy() - ecal_energy : ticl_cand.rawEnergy() - ecal_energy;
106  // fix for floating point rounding could go slightly below 0
107  hcal_energy = std::max(0.f, hcal_energy);
109  switch (abs_pdg_id) {
110  case 11:
111  part_type = reco::PFCandidate::e;
112  break;
113  case 13:
114  part_type = reco::PFCandidate::mu;
115  break;
116  case 22:
117  part_type = reco::PFCandidate::gamma;
118  break;
119  case 130:
120  part_type = reco::PFCandidate::h0;
121  break;
122  case 211:
123  part_type = reco::PFCandidate::h;
124  break;
125  // default also handles neutral pions (111) for the time being (not yet foreseen in PFCandidate)
126  default:
127  part_type = reco::PFCandidate::X;
128  }
129 
130  candidates->emplace_back(charge, four_mom, part_type);
131 
132  auto& candidate = candidates->back();
133  candidate.setEcalEnergy(ecal_energy, ecal_energy);
134  candidate.setHcalEnergy(hcal_energy, hcal_energy);
135  if (candidate.charge()) { // otherwise PFCandidate throws
136  // Construct edm::Ref from edm::Ptr. As of now, assumes type to be reco::Track. To be extended (either via
137  // dynamic type checking or configuration) if additional track types are needed.
138  reco::TrackRef trackref(ticl_cand.trackPtr().id(), int(ticl_cand.trackPtr().key()), &evt.productGetter());
139  candidate.setTrackRef(trackref);
140  // Utilize PFMuonAlgo
141  const int muId = PFMuonAlgo::muAssocToTrack(trackref, muons);
142  if (muId != -1) {
143  const reco::MuonRef muonref = reco::MuonRef(muonH, muId);
144  const bool allowLoose = (part_type == reco::PFCandidate::mu);
145  // Redefine pfmuon candidate kinematics and add muonref
146  pfmu_->reconstructMuon(candidate, muonref, allowLoose);
147  }
148  }
149 
150  // HGCAL timing as default values
151  auto time = ticl_cand.time();
152  auto timeE = ticl_cand.timeError();
153 
154  if (useMTDTiming_ and candidate.charge()) {
155  // Ignore HGCAL timing until it will be TOF corrected
156  time = -99.;
157  timeE = -1.;
158  // Check MTD timing availability
159  const bool assocQuality = (*trackTimeQualH)[candidate.trackRef()] > timingQualityThreshold_;
160  if (assocQuality) {
161  const auto timeHGC = time;
162  const auto timeEHGC = timeE;
163  const auto timeMTD = (*trackTimeH)[candidate.trackRef()];
164  const auto timeEMTD = (*trackTimeErrH)[candidate.trackRef()];
165 
166  if (useTimingAverage_ && (timeEMTD > 0 && timeEHGC > 0)) {
167  // Compute weighted average between HGCAL and MTD timing
168  const auto invTimeESqHGC = pow(timeEHGC, -2);
169  const auto invTimeESqMTD = pow(timeEMTD, -2);
170  timeE = (invTimeESqHGC * invTimeESqMTD) / (invTimeESqHGC + invTimeESqMTD);
171  time = (timeHGC * invTimeESqHGC + timeMTD * invTimeESqMTD) * timeE;
172  timeE = sqrt(timeE);
173  } else if (timeEMTD > 0) { // Ignore HGCal timing until it will be TOF corrected
174  time = timeMTD;
175  timeE = timeEMTD;
176  }
177  }
178  }
179  candidate.setTime(time, timeE);
180  }
181 
182  evt.put(std::move(candidates));
183 }
static int muAssocToTrack(const reco::TrackRef &trackref, const reco::MuonCollection &muons)
Definition: PFMuonAlgo.cc:479
const bool energy_from_regression_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
ParticleType
particle types
Definition: PFCandidate.h:44
~PFTICLProducer() override
const bool useTimingAverage_
const edm::EDGetTokenT< edm::ValueMap< float > > srcTrackTime_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
PFTICLProducer(const edm::ParameterSet &)
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
const edm::EDGetTokenT< edm::ValueMap< float > > srcTrackTimeError_
EDProductGetter const & productGetter() const
Definition: Event.cc:106
const edm::EDGetTokenT< edm::ValueMap< float > > srcTrackTimeQuality_
const edm::EDGetTokenT< edm::View< TICLCandidate > > ticl_candidates_
T sqrt(T t)
Definition: SSEVec.h:19
std::unique_ptr< PFMuonAlgo > pfmu_
const bool useMTDTiming_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
edm::Ref< MuonCollection > MuonRef
presistent reference to a Muon
Definition: MuonFwd.h:13
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const float timingQualityThreshold_
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
Definition: PFMuonAlgo.cc:1043
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void produce(edm::Event &, const edm::EventSetup &) override
ParameterSet const & getParameterSet(ParameterSetID const &id)
fixed size matrix
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:564
HLT enums.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511
const edm::EDGetTokenT< reco::MuonCollection > muons_