CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CaloTowerCreatorForTauHLT.cc
Go to the documentation of this file.
1 // makes CaloTowerCandidates from CaloTowers
2 // original author: L.Lista INFN, modifyed by: F.Ratnikov UMd
3 // Author for regionality A. Nikitenko
4 // Modified by S. Gennai
5 
12 // Math
13 #include "Math/GenVector/VectorUtil.h"
14 #include <cmath>
15 
16 using namespace edm;
17 using namespace reco;
18 using namespace std;
19 using namespace l1extra;
20 
22  : mVerbose(p.getUntrackedParameter<int>("verbose", 0)),
23  mtowers_token(consumes<CaloTowerCollection>(p.getParameter<InputTag>("towers"))),
24  mCone(p.getParameter<double>("UseTowersInCone")),
25  mTauTrigger_token(consumes<L1JetParticleCollection>(p.getParameter<InputTag>("TauTrigger"))),
26  mEtThreshold(p.getParameter<double>("minimumEt")),
27  mEThreshold(p.getParameter<double>("minimumE")),
28  mTauId(p.getParameter<int>("TauId")) {
29  produces<CaloTowerCollection>();
30 }
31 
33 
36  evt.getByToken(mtowers_token, caloTowers);
37 
38  // imitate L1 seeds
40  evt.getByToken(mTauTrigger_token, jetsgen);
41 
42  std::unique_ptr<CaloTowerCollection> cands(new CaloTowerCollection);
43  cands->reserve(caloTowers->size());
44 
45  int idTau = 0;
46  L1JetParticleCollection::const_iterator myL1Jet = jetsgen->begin();
47  for (; myL1Jet != jetsgen->end(); myL1Jet++) {
48  if (idTau == mTauId) {
49  double Sum08 = 0.;
50 
51  unsigned idx = 0;
52  for (; idx < caloTowers->size(); idx++) {
53  const CaloTower* cal = &((*caloTowers)[idx]);
54  bool isAccepted = false;
55  if (mVerbose == 2) {
56  edm::LogInfo("JetDebugInfo") << "CaloTowerCreatorForTauHLT::produce-> " << idx
57  << " tower et/eta/phi/e: " << cal->et() << '/' << cal->eta() << '/' << cal->phi()
58  << '/' << cal->energy() << " is...";
59  }
60  if (cal->et() >= mEtThreshold && cal->energy() >= mEThreshold) {
61  math::PtEtaPhiELorentzVector p(cal->et(), cal->eta(), cal->phi(), cal->energy());
62  double delta = ROOT::Math::VectorUtil::DeltaR((*myL1Jet).p4().Vect(), p);
63 
64  if (delta < mCone) {
65  isAccepted = true;
66  Sum08 += cal->et();
67  cands->push_back(*cal);
68  }
69  }
70  if (mVerbose == 2) {
71  if (isAccepted)
72  edm::LogInfo("JetDebugInfo") << "accepted \n";
73  else
74  edm::LogInfo("JetDebugInfo") << "rejected \n";
75  }
76  }
77  }
78  idTau++;
79  }
80 
81  evt.put(std::move(cands));
82 }
83 
86  aDesc.add<edm::InputTag>("TauTrigger", edm::InputTag("l1extraParticles", "Tau"))
87  ->setComment("L1ExtraJet collection for seeding");
88  aDesc.add<int>("TauId", 0)->setComment("Item from L1ExtraJet collection used for seeding");
89  aDesc.add<edm::InputTag>("towers", edm::InputTag("towerMaker"))->setComment("Input tower collection");
90  aDesc.add<double>("UseTowersInCone", 0.8)->setComment("Radius of cone around seed");
91  aDesc.add<double>("minimumE", 0.8)->setComment("Minimum tower energy");
92  aDesc.add<double>("minimumEt", 0.5)->setComment("Minimum tower ET");
93  aDesc.addUntracked<int>("verbose", 0)->setComment("Verbosity level; 0=silent");
94  desc.setComment("Produce tower collection around L1ExtraJetParticle seed.");
95  desc.add("caloTowerMakerHLT", aDesc);
96 }
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
const double mCone
use only towers in cone mCone around L1 candidate for regional jet reco
~CaloTowerCreatorForTauHLT() override
destructor
void produce(edm::StreamID sid, edm::Event &evt, const edm::EventSetup &stp) const override
process one event
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< l1extra::L1JetParticleCollection > mTauTrigger_token
label of tau trigger type analysis
Log< level::Info, false > LogInfo
void setComment(std::string const &value)
PtEtaPhiELorentzVectorD PtEtaPhiELorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:27
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< CaloTowerCollection > mtowers_token
label of source collection
mVerbose(fConfig.getUntrackedParameter< bool >("verbose", false))
CaloTowerCreatorForTauHLT(const edm::ParameterSet &)
constructor from parameter set
double et(double vtxZ) const
Definition: CaloTower.h:150
const double mEThreshold
E threshold.
const double mEtThreshold
imitator of L1 seeds
double phi() const final
momentum azimuthal angle
static void fillDescriptions(edm::ConfigurationDescriptions &desc)
double energy() const final
energy
double eta() const final
momentum pseudorapidity