CMS 3D CMS Logo

QjetsAdder.cc
Go to the documentation of this file.
1 #include <numeric>
2 
4 #include "fastjet/PseudoJet.hh"
5 #include "fastjet/ClusterSequence.hh"
6 
8 
9 using namespace std;
10 
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  const 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),
68  QjetsPreclustering_);
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 }
106 
virtual double energy() const =0
energy
double pt() const final
transverse momentum
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.
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: QjetsAdder.cc:11
virtual Constituents getJetConstituents() const
list of constituents
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
def move(src, dest)
Definition: eostools.py:511