CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QjetsAdder.cc
Go to the documentation of this file.
2 #include "fastjet/PseudoJet.hh"
3 #include "fastjet/ClusterSequence.hh"
4 
6 
7 using namespace std;
8 
10  // read input collection
11  //edm::Handle<edm::View<pat::Jet> > jets;
13  iEvent.getByToken(src_token_, jets);
14 
15  // prepare room for output
16  std::vector<float> QjetsVolatility; QjetsVolatility.reserve(jets->size());
17 
18  for ( typename edm::View<reco::Jet>::const_iterator jetIt = jets->begin() ; jetIt != jets->end() ; ++jetIt ) {
19  reco::Jet newCand(*jetIt);
20 
21  if(newCand.pt()<cutoff_)
22  {
23  QjetsVolatility.push_back(-1);
24  continue;
25  }
26 
27  //refill and recluster
28  vector<fastjet::PseudoJet> allconstits;
29  //for (unsigned k=0; k < newCand.getPFConstituents().size(); k++){
30  for (unsigned k=0; k < newCand.getJetConstituents().size(); k++){
31  const edm::Ptr<reco::Candidate> thisParticle = newCand.getJetConstituents().at(k);
32  allconstits.push_back( fastjet::PseudoJet( thisParticle->px(), thisParticle->py(), thisParticle->pz(), thisParticle->energy() ) );
33  }
34 
35  fastjet::JetDefinition jetDef(fastjet::cambridge_algorithm, jetRad_);
36  if (mJetAlgo_== "AK") jetDef.set_jet_algorithm( fastjet::antikt_algorithm );
37  else if (mJetAlgo_ == "CA") jetDef.set_jet_algorithm( fastjet::cambridge_algorithm );
38  else throw cms::Exception("GroomedJetFiller") << " unknown jet algorithm " << std::endl;
39 
40  fastjet::ClusterSequence thisClustering_basic(allconstits, jetDef);
41  std::vector<fastjet::PseudoJet> out_jets_basic = sorted_by_pt(thisClustering_basic.inclusive_jets(cutoff_));
42  //std::cout << newCand.pt() << " " << out_jets_basic.size() <<std::endl;
43  if(out_jets_basic.size()==0){ // jet reclustering failed, most likely due to the higher cutoff. Use a recognizeable default value for this jet
44  QjetsVolatility.push_back(-1);
45  continue;
46  }
47  assert( out_jets_basic.size()==1 ); // jet reclustering of one jet should yield exactly one jet at this stage
48 
49  // setup objects for qjets computation
50  fastjet::JetDefinition qjet_def(&qjetsAlgo_);
51 
52  std::vector<double> qjetmass;
53 
54  vector<fastjet::PseudoJet> constits;
55  unsigned int nqjetconstits = out_jets_basic.at(0).constituents().size(); // there should always be exactly one reclsutered jet => always "at(0)"
56  if (nqjetconstits < (unsigned int) QjetsPreclustering_) constits = out_jets_basic.at(0).constituents();
57  else constits = out_jets_basic.at(0).associated_cluster_sequence()->exclusive_subjets_up_to(out_jets_basic.at(0),QjetsPreclustering_);
58 
60  CLHEP::HepRandomEngine* engine = &rng->getEngine(iEvent.streamID());
61  qjetsAlgo_.SetRNEngine(engine);
62  // create probabilistic recusterings
63  for(unsigned int ii = 0 ; ii < (unsigned int) ntrial_ ; ii++){
64  //qjetsAlgo_.SetRandSeed(iEvent.id().event()*100 + (jetIt - jets->begin())*ntrial_ + ii );// set random seed for reprudcibility. We need a smarted scheme
65  fastjet::ClusterSequence qjet_seq(constits, qjet_def);
66  vector<fastjet::PseudoJet> inclusive_jets2 = sorted_by_pt(qjet_seq.inclusive_jets(cutoff_));
67  if (inclusive_jets2.size()>0){ // fill the massvalue only if the reclustering was successfull
68  qjetmass.push_back(inclusive_jets2[0].m());
69  }
70 
71  }//end loop on trials
72 
73  if (qjetmass.size()<1) {//protection against dummy case
74  QjetsVolatility.push_back(-1);
75  continue;
76  }
77 
78  double mean = std::accumulate( qjetmass.begin( ) , qjetmass.end( ) , 0 ) /qjetmass.size() ;
79  float totalsquared = 0.;
80  for (unsigned int i = 0; i < qjetmass.size(); i++){
81  totalsquared += (qjetmass[i] - mean)*(qjetmass[i] - mean) ;
82  }
83  float variance = sqrt( totalsquared/qjetmass.size() );
84 
85  QjetsVolatility.push_back(variance/mean);
86  }//end loop on jets
87 
88  std::auto_ptr<edm::ValueMap<float> > outQJV(new edm::ValueMap<float>());
89  edm::ValueMap<float>::Filler fillerQJV(*outQJV);
90  fillerQJV.insert(jets, QjetsVolatility.begin(), QjetsVolatility.end());
91  fillerQJV.fill();
92 
93  iEvent.put(outQJV,"QjetsVolatility");
94 }
95 
96 
97 
98 
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#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:9
assert(m_qm.get())
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
virtual Constituents getJetConstituents() const
list of constituents
int ii
Definition: cuy.py:588
virtual double pt() const
transverse momentum
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &) const =0
Use this engine in event methods.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
StreamID streamID() const
Definition: Event.h:74