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; QjetsVolatility.reserve(jets->size());
19 
20  for ( typename edm::View<reco::Jet>::const_iterator jetIt = jets->begin() ; jetIt != jets->end() ; ++jetIt ) {
21  reco::Jet newCand(*jetIt);
22 
23  if(newCand.pt()<cutoff_)
24  {
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( fastjet::PseudoJet( thisParticle->px(), thisParticle->py(), thisParticle->pz(), thisParticle->energy() ) );
35  }
36 
37  fastjet::JetDefinition jetDef(fastjet::cambridge_algorithm, jetRad_);
38  if (mJetAlgo_== "AK") jetDef.set_jet_algorithm( fastjet::antikt_algorithm );
39  else if (mJetAlgo_ == "CA") jetDef.set_jet_algorithm( fastjet::cambridge_algorithm );
40  else throw cms::Exception("GroomedJetFiller") << " unknown jet algorithm " << std::endl;
41 
42  fastjet::ClusterSequence thisClustering_basic(allconstits, jetDef);
43  std::vector<fastjet::PseudoJet> out_jets_basic = sorted_by_pt(thisClustering_basic.inclusive_jets(cutoff_));
44  //std::cout << newCand.pt() << " " << out_jets_basic.size() <<std::endl;
45  if(out_jets_basic.size()!=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
46  QjetsVolatility.push_back(-1);
47  continue;
48  }
49 
50  // setup objects for qjets computation
51  fastjet::JetDefinition qjet_def(&qjetsAlgo_);
52 
53  std::vector<double> qjetmass;
54 
55  vector<fastjet::PseudoJet> constits;
56  unsigned int nqjetconstits = out_jets_basic.at(0).constituents().size(); // there should always be exactly one reclsutered jet => always "at(0)"
57  if (nqjetconstits < (unsigned int) QjetsPreclustering_) constits = out_jets_basic.at(0).constituents();
58  else constits = out_jets_basic.at(0).associated_cluster_sequence()->exclusive_subjets_up_to(out_jets_basic.at(0),QjetsPreclustering_);
59 
61  CLHEP::HepRandomEngine* engine = &rng->getEngine(iEvent.streamID());
62  qjetsAlgo_.SetRNEngine(engine);
63  // create probabilistic recusterings
64  for(unsigned int ii = 0 ; ii < (unsigned int) ntrial_ ; ii++){
65  //qjetsAlgo_.SetRandSeed(iEvent.id().event()*100 + (jetIt - jets->begin())*ntrial_ + ii );// set random seed for reprudcibility. We need a smarted scheme
66  fastjet::ClusterSequence qjet_seq(constits, qjet_def);
67  vector<fastjet::PseudoJet> inclusive_jets2 = sorted_by_pt(qjet_seq.inclusive_jets(cutoff_));
68  if (inclusive_jets2.size()>0){ // fill the massvalue only if the reclustering was successfull
69  qjetmass.push_back(inclusive_jets2[0].m());
70  }
71 
72  }//end loop on trials
73 
74  if (qjetmass.size()<1) {//protection against dummy case
75  QjetsVolatility.push_back(-1);
76  continue;
77  }
78 
79  double mean = std::accumulate( qjetmass.begin( ) , qjetmass.end( ) , 0 ) /qjetmass.size() ;
80  float totalsquared = 0.;
81  for (unsigned int i = 0; i < qjetmass.size(); i++){
82  totalsquared += (qjetmass[i] - mean)*(qjetmass[i] - mean) ;
83  }
84  float variance = sqrt( totalsquared/qjetmass.size() );
85 
86  QjetsVolatility.push_back(variance/mean);
87  }//end loop on jets
88 
89  auto outQJV = std::make_unique<edm::ValueMap<float>>();
90  edm::ValueMap<float>::Filler fillerQJV(*outQJV);
91  fillerQJV.insert(jets, QjetsVolatility.begin(), QjetsVolatility.end());
92  fillerQJV.fill();
93 
94  iEvent.put(std::move(outQJV),"QjetsVolatility");
95 }
96 
97 
98 
99 
virtual double pt() const final
transverse momentum
virtual double pz() const =0
z coordinate of momentum vector
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Base class for all types of Jets.
Definition: Jet.h:20
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
Definition: QjetsAdder.cc:11
virtual Constituents getJetConstituents() const
list of constituents
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
virtual double energy() const =0
energy
virtual double py() const =0
y coordinate of momentum vector
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
ii
Definition: cuy.py:588
int k[5][pyjets_maxn]
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
StreamID streamID() const
Definition: Event.h:81
virtual double px() const =0
x coordinate of momentum vector
def move(src, dest)
Definition: eostools.py:510