CMS 3D CMS Logo

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  :
23  mVerbose (p.getUntrackedParameter<int> ("verbose", 0)),
24  mtowers_token (consumes<CaloTowerCollection>(p.getParameter<InputTag>("towers") ) ),
25  mCone (p.getParameter<double> ("UseTowersInCone")),
26  mTauTrigger_token (consumes<L1JetParticleCollection>(p.getParameter<InputTag>("TauTrigger") ) ),
27  mEtThreshold (p.getParameter<double> ("minimumEt")),
28  mEThreshold (p.getParameter<double> ("minimumE")),
29  mTauId (p.getParameter<int> ("TauId"))
30 {
31  produces<CaloTowerCollection>();
32 }
33 
35 }
36 
37 void CaloTowerCreatorForTauHLT::produce( StreamID sid, Event& evt, const EventSetup& stp ) const {
39  evt.getByToken( mtowers_token, caloTowers );
40 
41  // imitate L1 seeds
43  evt.getByToken( mTauTrigger_token, jetsgen);
44 
45  std::unique_ptr<CaloTowerCollection> cands( new CaloTowerCollection );
46  cands->reserve( caloTowers->size() );
47 
48  int idTau =0;
49  L1JetParticleCollection::const_iterator myL1Jet = jetsgen->begin();
50  for(;myL1Jet != jetsgen->end();myL1Jet++)
51  {
52  if(idTau == mTauId)
53  {
54  double Sum08 = 0.;
55 
56  unsigned idx = 0;
57  for (; idx < caloTowers->size(); idx++) {
58  const CaloTower* cal = &((*caloTowers) [idx]);
59  bool isAccepted = false;
60  if (mVerbose == 2) {
61  edm::LogInfo("JetDebugInfo") << "CaloTowerCreatorForTauHLT::produce-> " << idx << " tower et/eta/phi/e: "
62  << cal->et() << '/' << cal->eta() << '/' << cal->phi() << '/' << cal->energy() << " is...";
63  }
64  if (cal->et() >= mEtThreshold && cal->energy() >= mEThreshold ) {
65  math::PtEtaPhiELorentzVector p( cal->et(), cal->eta(), cal->phi(), cal->energy() );
66  double delta = ROOT::Math::VectorUtil::DeltaR((*myL1Jet).p4().Vect(), p);
67 
68  if(delta < mCone) {
69  isAccepted = true;
70  Sum08 += cal->et();
71  cands->push_back( *cal );
72  }
73  }
74  if (mVerbose == 2){
75  if (isAccepted) edm::LogInfo("JetDebugInfo") << "accepted \n";
76  else edm::LogInfo("JetDebugInfo") << "rejected \n";
77  }
78  }
79 
80  }
81  idTau++;
82  }
83 
84  evt.put(std::move(cands));
85 
86 }
87 
89 
91  aDesc.add<edm::InputTag>("TauTrigger",edm::InputTag("l1extraParticles","Tau"))->setComment("L1ExtraJet collection for seeding");
92  aDesc.add<int>("TauId",0)->setComment("Item from L1ExtraJet collection used for seeding");
93  aDesc.add<edm::InputTag>("towers",edm::InputTag("towerMaker"))->setComment("Input tower collection");
94  aDesc.add<double>("UseTowersInCone",0.8)->setComment("Radius of cone around seed");
95  aDesc.add<double>("minimumE",0.8)->setComment("Minimum tower energy");
96  aDesc.add<double>("minimumEt",0.5)->setComment("Minimum tower ET");
97  aDesc.addUntracked<int>("verbose",0)->setComment("Verbosity level; 0=silent");
98  desc.setComment("Produce tower collection around L1ExtraJetParticle seed.");
99  desc.add("caloTowerMakerHLT",aDesc);
100 
101 }
dbl * delta
Definition: mlp_gen.cc:36
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
double eta() const final
momentum pseudorapidity
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
const double mCone
use only towers in cone mCone around L1 candidate for regional jet reco
~CaloTowerCreatorForTauHLT() override
destructor
double energy() const final
energy
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< l1extra::L1JetParticleCollection > mTauTrigger_token
label of tau trigger type analysis
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)
fixed size matrix
HLT enums.
const edm::EDGetTokenT< CaloTowerCollection > mtowers_token
label of source collection
size_type size() const
CaloTowerCreatorForTauHLT(const edm::ParameterSet &)
constructor from parameter set
double et(double vtxZ) const
Definition: CaloTower.h:132
const double mEThreshold
E threshold.
const double mEtThreshold
imitator of L1 seeds
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:510
void produce(edm::StreamID sid, edm::Event &evt, const edm::EventSetup &stp) const override
process one event
static void fillDescriptions(edm::ConfigurationDescriptions &desc)