CMS 3D CMS Logo

L1TkEmParticleProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 //
4 // dummy producer for a TkEm
5 //
6 
7 // system include files
8 #include <memory>
9 
10 // user include files
13 
19 
21 
25 
28 
29 // for L1Tracks:
31 
32 #include <string>
33 #include <cmath>
34 
35 static constexpr float EtaECal = 1.479;
36 static constexpr float REcal = 129.;
37 static constexpr float ZEcal = 315.4;
38 using namespace l1t;
39 //
40 // class declaration
41 //
42 
44 public:
46  typedef std::vector<L1TTTrackType> L1TTTrackCollectionType;
47 
49  ~L1TkEmParticleProducer() override;
50 
51  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
52 
53  float CorrectedEta(float eta, float zv) const;
54 
55 private:
56  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
57 
58  // ----------member data ---------------------------
59 
63 
65 
66  const float etMin_; // min ET in GeV of L1EG objects
67  const float zMax_; // |z_track| < zMax_ in cm
68  const float chi2Max_;
69  const float pTMinTra_;
70  const float dR2Min_;
71  const float dR2Max_;
72  const bool primaryVtxConstrain_; // use the primary vertex (default = false)
73  const float deltaZMax_; // | z_track - z_primaryvtx | < deltaZMax_ in cm.
74  // Used only when primaryVtxConstrain_ = True.
75  const float isoCut_;
76  const bool relativeIsolation_;
77 };
78 
79 //
80 // constructors and destructor
81 //
83  : egToken_(consumes<EGammaBxCollection>(iConfig.getParameter<edm::InputTag>("L1EGammaInputTag"))),
84  trackToken_(consumes<L1TTTrackCollectionType>(iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
85  vertexToken_(consumes<TkPrimaryVertexCollection>(iConfig.getParameter<edm::InputTag>("L1VertexInputTag"))),
86  label_(iConfig.getParameter<std::string>("label")),
87  etMin_((float)iConfig.getParameter<double>("ETmin")),
88  zMax_((float)iConfig.getParameter<double>("ZMAX")),
89  chi2Max_((float)iConfig.getParameter<double>("CHI2MAX")),
90  pTMinTra_((float)iConfig.getParameter<double>("PTMINTRA")),
91  dR2Min_((float)iConfig.getParameter<double>("DRmin") * (float)iConfig.getParameter<double>("DRmin")),
92  dR2Max_((float)iConfig.getParameter<double>("DRmax") * (float)iConfig.getParameter<double>("DRmax")),
93  primaryVtxConstrain_(iConfig.getParameter<bool>("PrimaryVtxConstrain")),
94  deltaZMax_((float)iConfig.getParameter<double>("DeltaZMax")),
95  isoCut_((float)iConfig.getParameter<double>("IsoCut")),
96  relativeIsolation_(iConfig.getParameter<bool>("RelativeIsolation")) {
97  produces<TkEmCollection>(label_);
98 }
99 
101 
102 // ------------ method called to produce the data ------------
104  using namespace edm;
105 
106  auto result = std::make_unique<TkEmCollection>();
107 
108  // the L1Tracks
109  edm::Handle<L1TTTrackCollectionType> L1TTTrackHandle;
110  iEvent.getByToken(trackToken_, L1TTTrackHandle);
111  if (!L1TTTrackHandle.isValid()) {
112  LogError("L1TkEmParticleProducer") << "\nWarning: L1TTTrackCollectionType not found in the event. Exit."
113  << std::endl;
114  return;
115  }
116 
117  // the primary vertex (used only if primaryVtxConstrain_ = true)
118  float zvtxL1tk = -999;
121  iEvent.getByToken(vertexToken_, L1VertexHandle);
122  if (!L1VertexHandle.isValid()) {
123  LogWarning("L1TkEmParticleProducer")
124  << "Warning: TkPrimaryVertexCollection not found in the event. Won't use any PrimaryVertex constraint."
125  << std::endl;
126  primaryVtxConstrain = false;
127  } else {
128  std::vector<TkPrimaryVertex>::const_iterator vtxIter = L1VertexHandle->begin();
129  // by convention, the first vertex in the collection is the one that should
130  // be used by default
131  zvtxL1tk = vtxIter->zvertex();
132  }
133 
134  if (!L1TTTrackHandle.isValid()) {
135  LogError("L1TkEmParticleProducer") << "\nWarning: L1TTTrackCollectionType not found in the event. Exit."
136  << std::endl;
137  return;
138  }
139 
140  // Now loop over the L1EGamma objects
141 
142  edm::Handle<EGammaBxCollection> eGammaHandle;
143  iEvent.getByToken(egToken_, eGammaHandle);
144  if (!eGammaHandle.isValid()) {
145  LogError("L1TkEmParticleProducer") << "\nWarning: L1EmCollection not found in the event. Exit." << std::endl;
146  return;
147  }
148  EGammaBxCollection eGammaCollection = (*eGammaHandle.product());
150 
151  int ieg = 0;
152  for (egIter = eGammaCollection.begin(0); egIter != eGammaCollection.end(0); ++egIter) // considering BX = only
153  {
154  edm::Ref<EGammaBxCollection> EGammaRef(eGammaHandle, ieg);
155  ieg++;
156 
157  float et = egIter->et();
158  if (et < etMin_)
159  continue;
160 
161  float eta = egIter->eta();
162  // The eta of the L1EG object is seen from (0,0,0).
163  // if primaryVtxConstrain_ = true, and for the PV constrained iso, use the zvtxL1tk to correct the eta(L1EG)
164  // that is used in the calculation of DeltaR.
165  float etaPV = CorrectedEta(eta, zvtxL1tk);
166 
167  float phi = egIter->phi();
168 
169  // calculate the isolation of the L1EG object with
170  // respect to L1Tracks.
171 
172  float sumPt = 0;
173  float sumPtPV = 0;
174 
175  for (const auto& track : *L1TTTrackHandle) {
176  float Pt = track.momentum().perp();
177  float Eta = track.momentum().eta();
178  float Phi = track.momentum().phi();
179  float z = track.POCA().z();
180  if (fabs(z) > zMax_)
181  continue;
182  if (Pt < pTMinTra_)
183  continue;
184  float chi2 = track.chi2();
185  if (chi2 > chi2Max_)
186  continue;
187 
188  float dr2 = reco::deltaR2(Eta, Phi, eta, phi);
189  if (dr2 < dR2Max_ && dr2 >= dR2Min_) {
190  sumPt += Pt;
191  }
192 
193  if (zvtxL1tk > -999 && std::abs(z - zvtxL1tk) >= deltaZMax_)
194  continue; // Now, PV constrained trackSum:
195 
196  dr2 = reco::deltaR2(Eta, Phi, etaPV, phi); // recompute using the corrected eta
197  if (dr2 < dR2Max_ && dr2 >= dR2Min_) {
198  sumPtPV += Pt;
199  }
200  } // end loop over tracks
201 
202  float trkisol = -999;
203  float trkisolPV = -999;
204  if (relativeIsolation_) {
205  if (et > 0) {
206  trkisol = sumPt / et;
207  trkisolPV = sumPtPV / et;
208  } // relative isolation
209  } else { // absolute isolation
210  trkisol = sumPt;
211  trkisolPV = sumPtPV;
212  }
213 
214  float isolation = primaryVtxConstrain ? trkisolPV : trkisol;
215 
216  const math::XYZTLorentzVector P4 = egIter->p4();
217  TkEm trkEm(P4, EGammaRef, trkisol, trkisolPV);
218 
219  if (isoCut_ <= 0) {
220  // write the L1TkEm particle to the collection,
221  // irrespective of its relative isolation
222  result->push_back(trkEm);
223  } else {
224  // the object is written to the collection only
225  // if it passes the isolation cut
226  if (isolation <= isoCut_)
227  result->push_back(trkEm);
228  }
229 
230  } // end loop over EGamma objects
231 
232  iEvent.put(std::move(result), label_);
233 }
234 
235 // --------------------------------------------------------------------------------------
236 
237 float L1TkEmParticleProducer::CorrectedEta(float eta, float zv) const {
238  // Correct the eta of the L1EG object once we know the zvertex
239 
240  if (zv == 0.)
241  return eta;
242 
243  if (eta == 0) {
244  float thetaprime = atan(-REcal / zv);
245  if (thetaprime < 0)
246  thetaprime = thetaprime + M_PI;
247  float etaprime = -log(tan(0.5 * thetaprime));
248  return etaprime;
249  }
250 
251  bool IsBarrel = (std::abs(eta) < EtaECal);
252 
253  float tanhalftheta = exp(-eta);
254  float tantheta = 2. * tanhalftheta / (1. - tanhalftheta * tanhalftheta);
255 
256  float delta;
257  if (IsBarrel)
258  delta = REcal / tantheta;
259  else
260  delta = eta > 0 ? ZEcal : -ZEcal;
261 
262  float etaprime;
263  if (delta == zv) {
264  etaprime = 0.;
265  } else {
266  float tanthetaprime = delta * tantheta / (delta - zv);
267  float thetaprime = atan(tanthetaprime);
268  if (thetaprime < 0)
269  thetaprime = thetaprime + M_PI;
270  etaprime = -log(tan(0.5 * thetaprime));
271  }
272 
273  return etaprime;
274 }
275 
276 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
278  //The following says we do not know what parameters are allowed so do no validation
279  // Please change this to state exactly what you do use, even if it is no parameters
281  desc.add<string>("label", "EG");
282  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("simCaloStage2Digis"));
283  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
284  desc.add<edm::InputTag>("L1VertexInputTag", edm::InputTag("L1TkPrimaryVertex"));
285  desc.add<double>("ETmin", -1.);
286  desc.add<bool>("RelativeIsolation", true);
287  desc.add<double>("IsoCut", 0.23);
288  desc.add<double>("ZMAX", 25.);
289  desc.add<double>("CHI2MAX", 100.);
290  desc.add<double>("PTMINTRA", 2.);
291  desc.add<double>("DRmin", 0.07);
292  desc.add<double>("DRmax", 0.30);
293  desc.add<bool>("PrimaryVtxConstrain", false);
294  desc.add<double>("DeltaZMax", 0.6);
295  descriptions.add("l1TkEmParticleProducer", desc);
296 }
297 
298 //define this as a plug-in
static constexpr float EtaECal
edm::Ref< EGammaBxCollection > EGammaRef
Definition: TkEGTau.h:35
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float *__restrict__ zv
T const * product() const
Definition: Handle.h:70
const edm::EDGetTokenT< EGammaBxCollection > egToken_
delete x;
Definition: CaloConfig.h:22
Log< level::Error, false > LogError
static constexpr float REcal
const edm::EDGetTokenT< TkPrimaryVertexCollection > vertexToken_
const_iterator begin(int bx) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const edm::EDGetTokenT< L1TTTrackCollectionType > trackToken_
std::vector< T >::const_iterator const_iterator
Definition: BXVector.h:18
int iEvent
Definition: GenABIO.cc:224
Definition: TkEm.h:19
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define M_PI
L1TkEmParticleProducer(const edm::ParameterSet &)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
std::vector< TkPrimaryVertex > TkPrimaryVertexCollection
float CorrectedEta(float eta, float zv) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
const_iterator end(int bx) const
HLT enums.
static constexpr float ZEcal
std::vector< L1TTTrackType > L1TTTrackCollectionType
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511