CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
QjetsAdder Class Reference

#include <QjetsAdder.h>

Inheritance diagram for QjetsAdder:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 QjetsAdder (const edm::ParameterSet &iConfig)
 
virtual ~QjetsAdder ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

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

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 13 of file QjetsAdder.h.

Constructor & Destructor Documentation

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

Definition at line 15 of file QjetsAdder.h.

15  :
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  {
29  produces<edm::ValueMap<float> >("QjetsVolatility");
30  }
T getParameter(std::string const &) const
int ntrial_
Definition: QjetsAdder.h:40
double jetRad_
Definition: QjetsAdder.h:42
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
QjetsPlugin qjetsAlgo_
Definition: QjetsAdder.h:39
int QjetsPreclustering_
Definition: QjetsAdder.h:44
double cutoff_
Definition: QjetsAdder.h:41
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
Definition: QjetsAdder.h:38
std::string mJetAlgo_
Definition: QjetsAdder.h:43
edm::InputTag src_
Definition: QjetsAdder.h:37
virtual QjetsAdder::~QjetsAdder ( )
inlinevirtual

Definition at line 32 of file QjetsAdder.h.

32 {}

Member Function Documentation

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

Implements edm::EDProducer.

Definition at line 9 of file QjetsAdder.cc.

References assert(), edm::hlt::Exception, edm::helper::Filler< Map >::fill(), edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), reco::Jet::getJetConstituents(), i, cuy::ii, edm::helper::Filler< Map >::insert(), fwrapper::jets, relval_steps::k, visualization-live-secondInstance_cfg::m, timingPdfMaker::mean, reco::LeafCandidate::pt(), edm::Event::put(), mathSSE::sqrt(), and edm::Event::streamID().

9  {
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 }
int i
Definition: DBlmapReader.cc:9
void SetRNEngine(CLHEP::HepRandomEngine *rnEngine)
Definition: QjetsPlugin.h:26
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
int ntrial_
Definition: QjetsAdder.h:40
Base class for all types of Jets.
Definition: Jet.h:20
assert(m_qm.get())
int ii
Definition: cuy.py:588
double jetRad_
Definition: QjetsAdder.h:42
QjetsPlugin qjetsAlgo_
Definition: QjetsAdder.h:39
int QjetsPreclustering_
Definition: QjetsAdder.h:44
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
T sqrt(T t)
Definition: SSEVec.h:48
double cutoff_
Definition: QjetsAdder.h:41
vector< PseudoJet > jets
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
Definition: QjetsAdder.h:38
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &) const =0
Use this engine in event methods.
std::string mJetAlgo_
Definition: QjetsAdder.h:43
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
StreamID streamID() const
Definition: Event.h:72

Member Data Documentation

double QjetsAdder::cutoff_
private

Definition at line 41 of file QjetsAdder.h.

double QjetsAdder::jetRad_
private

Definition at line 42 of file QjetsAdder.h.

std::string QjetsAdder::mJetAlgo_
private

Definition at line 43 of file QjetsAdder.h.

int QjetsAdder::ntrial_
private

Definition at line 40 of file QjetsAdder.h.

QjetsPlugin QjetsAdder::qjetsAlgo_
private

Definition at line 39 of file QjetsAdder.h.

int QjetsAdder::QjetsPreclustering_
private

Definition at line 44 of file QjetsAdder.h.

edm::Service<edm::RandomNumberGenerator> QjetsAdder::rng_
private

Definition at line 45 of file QjetsAdder.h.

edm::InputTag QjetsAdder::src_
private

Definition at line 37 of file QjetsAdder.h.

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

Definition at line 38 of file QjetsAdder.h.