CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
QjetsAdder Class Reference

#include <QjetsAdder.h>

Inheritance diagram for QjetsAdder:
edm::stream::EDProducer<>

Public Member Functions

void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 QjetsAdder (const edm::ParameterSet &iConfig)
 
 ~QjetsAdder () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Attributes

double cutoff_
 
double jetRad_
 
std::string mJetAlgo_
 
int ntrial_
 
QjetsPlugin qjetsAlgo_
 
int QjetsPreclustering_
 
edm::InputTag src_
 
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 13 of file QjetsAdder.h.

Constructor & Destructor Documentation

◆ QjetsAdder()

QjetsAdder::QjetsAdder ( const edm::ParameterSet iConfig)
inlineexplicit

Definition at line 15 of file QjetsAdder.h.

16  : src_(iConfig.getParameter<edm::InputTag>("src")),
18  qjetsAlgo_(iConfig.getParameter<double>("zcut"),
19  iConfig.getParameter<double>("dcutfctr"),
20  iConfig.getParameter<double>("expmin"),
21  iConfig.getParameter<double>("expmax"),
22  iConfig.getParameter<double>("rigidity")),
23  ntrial_(iConfig.getParameter<int>("ntrial")),
24  cutoff_(iConfig.getParameter<double>("cutoff")),
25  jetRad_(iConfig.getParameter<double>("jetRad")),
26  mJetAlgo_(iConfig.getParameter<std::string>("jetAlgo")),
27  QjetsPreclustering_(iConfig.getParameter<int>("preclustering")) {
28  produces<edm::ValueMap<float>>("QjetsVolatility");
29  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int ntrial_
Definition: QjetsAdder.h:39
double jetRad_
Definition: QjetsAdder.h:41
QjetsPlugin qjetsAlgo_
Definition: QjetsAdder.h:38
int QjetsPreclustering_
Definition: QjetsAdder.h:43
double cutoff_
Definition: QjetsAdder.h:40
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
Definition: QjetsAdder.h:37
std::string mJetAlgo_
Definition: QjetsAdder.h:42
edm::InputTag src_
Definition: QjetsAdder.h:36

◆ ~QjetsAdder()

QjetsAdder::~QjetsAdder ( )
inlineoverride

Definition at line 31 of file QjetsAdder.h.

31 {}

Member Function Documentation

◆ produce()

void QjetsAdder::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 11 of file QjetsAdder.cc.

References reco::Candidate::energy(), Exception, edm::RandomNumberGenerator::getEngine(), reco::Jet::getJetConstituents(), mps_fire::i, iEvent, cuy::ii, createfilelist::int, singleTopDQM_cfi::jets, dqmdumpme::k, visualization-live-secondInstance_cfg::m, SiStripPI::mean, eostools::move(), reco::LeafCandidate::pt(), reco::Candidate::px(), reco::Candidate::py(), reco::Candidate::pz(), and mathSSE::sqrt().

11  {
12  // read input collection
13  //edm::Handle<edm::View<pat::Jet> > jets;
15  iEvent.getByToken(src_token_, jets);
16 
17  // prepare room for output
18  std::vector<float> QjetsVolatility;
19  QjetsVolatility.reserve(jets->size());
20 
21  for (typename edm::View<reco::Jet>::const_iterator jetIt = jets->begin(); jetIt != jets->end(); ++jetIt) {
22  reco::Jet newCand(*jetIt);
23 
24  if (newCand.pt() < cutoff_) {
25  QjetsVolatility.push_back(-1);
26  continue;
27  }
28 
29  //refill and recluster
30  vector<fastjet::PseudoJet> allconstits;
31  //for (unsigned k=0; k < newCand.getPFConstituents().size(); k++){
32  for (unsigned k = 0; k < newCand.getJetConstituents().size(); k++) {
33  const edm::Ptr<reco::Candidate> thisParticle = newCand.getJetConstituents().at(k);
34  allconstits.push_back(
35  fastjet::PseudoJet(thisParticle->px(), thisParticle->py(), thisParticle->pz(), thisParticle->energy()));
36  }
37 
38  fastjet::JetDefinition jetDef(fastjet::cambridge_algorithm, jetRad_);
39  if (mJetAlgo_ == "AK")
40  jetDef.set_jet_algorithm(fastjet::antikt_algorithm);
41  else if (mJetAlgo_ == "CA")
42  jetDef.set_jet_algorithm(fastjet::cambridge_algorithm);
43  else
44  throw cms::Exception("GroomedJetFiller") << " unknown jet algorithm " << std::endl;
45 
46  fastjet::ClusterSequence thisClustering_basic(allconstits, jetDef);
47  std::vector<fastjet::PseudoJet> out_jets_basic = sorted_by_pt(thisClustering_basic.inclusive_jets(cutoff_));
48  //std::cout << newCand.pt() << " " << out_jets_basic.size() <<std::endl;
49  if (out_jets_basic.size() !=
50  1) { // jet reclustering did not return exactly 1 jet, most likely due to the higher cutoff or large cone size. Use a recognizeable default value for this jet
51  QjetsVolatility.push_back(-1);
52  continue;
53  }
54 
55  // setup objects for qjets computation
56  fastjet::JetDefinition qjet_def(&qjetsAlgo_);
57 
58  std::vector<double> qjetmass;
59 
60  vector<fastjet::PseudoJet> constits;
61  unsigned int nqjetconstits = out_jets_basic.at(0)
62  .constituents()
63  .size(); // there should always be exactly one reclsutered jet => always "at(0)"
64  if (nqjetconstits < (unsigned int)QjetsPreclustering_)
65  constits = out_jets_basic.at(0).constituents();
66  else
67  constits = out_jets_basic.at(0).associated_cluster_sequence()->exclusive_subjets_up_to(out_jets_basic.at(0),
69 
71  CLHEP::HepRandomEngine* engine = &rng->getEngine(iEvent.streamID());
72  qjetsAlgo_.SetRNEngine(engine);
73  // create probabilistic recusterings
74  for (unsigned int ii = 0; ii < (unsigned int)ntrial_; ii++) {
75  //qjetsAlgo_.SetRandSeed(iEvent.id().event()*100 + (jetIt - jets->begin())*ntrial_ + ii );// set random seed for reprudcibility. We need a smarted scheme
76  fastjet::ClusterSequence qjet_seq(constits, qjet_def);
77  vector<fastjet::PseudoJet> inclusive_jets2 = sorted_by_pt(qjet_seq.inclusive_jets(cutoff_));
78  if (!inclusive_jets2.empty()) { // fill the massvalue only if the reclustering was successfull
79  qjetmass.push_back(inclusive_jets2[0].m());
80  }
81 
82  } //end loop on trials
83 
84  if (qjetmass.empty()) { //protection against dummy case
85  QjetsVolatility.push_back(-1);
86  continue;
87  }
88 
89  double mean = std::accumulate(qjetmass.begin(), qjetmass.end(), 0) / qjetmass.size();
90  float totalsquared = 0.;
91  for (unsigned int i = 0; i < qjetmass.size(); i++) {
92  totalsquared += (qjetmass[i] - mean) * (qjetmass[i] - mean);
93  }
94  float variance = sqrt(totalsquared / qjetmass.size());
95 
96  QjetsVolatility.push_back(variance / mean);
97  } //end loop on jets
98 
99  auto outQJV = std::make_unique<edm::ValueMap<float>>();
100  edm::ValueMap<float>::Filler fillerQJV(*outQJV);
101  fillerQJV.insert(jets, QjetsVolatility.begin(), QjetsVolatility.end());
102  fillerQJV.fill();
103 
104  iEvent.put(std::move(outQJV), "QjetsVolatility");
105 }
virtual double energy() const =0
energy
void SetRNEngine(CLHEP::HepRandomEngine *rnEngine)
Definition: QjetsPlugin.h:28
int ntrial_
Definition: QjetsAdder.h:39
Base class for all types of Jets.
Definition: Jet.h:20
virtual double pz() const =0
z coordinate of momentum vector
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
double jetRad_
Definition: QjetsAdder.h:41
QjetsPlugin qjetsAlgo_
Definition: QjetsAdder.h:38
int QjetsPreclustering_
Definition: QjetsAdder.h:43
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
double cutoff_
Definition: QjetsAdder.h:40
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
Definition: QjetsAdder.h:37
virtual double py() const =0
y coordinate of momentum vector
ii
Definition: cuy.py:589
virtual double px() const =0
x coordinate of momentum vector
std::string mJetAlgo_
Definition: QjetsAdder.h:42
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ cutoff_

double QjetsAdder::cutoff_
private

Definition at line 40 of file QjetsAdder.h.

◆ jetRad_

double QjetsAdder::jetRad_
private

Definition at line 41 of file QjetsAdder.h.

◆ mJetAlgo_

std::string QjetsAdder::mJetAlgo_
private

Definition at line 42 of file QjetsAdder.h.

◆ ntrial_

int QjetsAdder::ntrial_
private

Definition at line 39 of file QjetsAdder.h.

◆ qjetsAlgo_

QjetsPlugin QjetsAdder::qjetsAlgo_
private

Definition at line 38 of file QjetsAdder.h.

◆ QjetsPreclustering_

int QjetsAdder::QjetsPreclustering_
private

Definition at line 43 of file QjetsAdder.h.

◆ src_

edm::InputTag QjetsAdder::src_
private

Definition at line 36 of file QjetsAdder.h.

◆ src_token_

edm::EDGetTokenT<edm::View<reco::Jet> > QjetsAdder::src_token_
private

Definition at line 37 of file QjetsAdder.h.