CMS 3D CMS Logo

CSJetProducer.cc
Go to the documentation of this file.
5 
6 #include "fastjet/contrib/ConstituentSubtractor.hh"
7 
8 using namespace std;
9 using namespace reco;
10 using namespace edm;
11 using namespace cms;
12 
13 CSJetProducer::CSJetProducer(edm::ParameterSet const& conf):
14  VirtualJetProducer( conf ),
15  csRParam_(-1.0),
16  csAlpha_(0.)
17 {
18  //get eta range, rho and rhom map
19  etaToken_ = consumes<std::vector<double>>(conf.getParameter<edm::InputTag>( "etaMap" ));
20  rhoToken_ = consumes<std::vector<double>>(conf.getParameter<edm::InputTag>( "rho" ));
21  rhomToken_ = consumes<std::vector<double>>(conf.getParameter<edm::InputTag>( "rhom" ));
22  csRParam_ = conf.getParameter<double>("csRParam");
23  csAlpha_ = conf.getParameter<double>("csAlpha");
24 }
25 
27 {
28  // use the default production from one collection
29  VirtualJetProducer::produce( iEvent, iSetup );
30  //use runAlgorithm of this class
31 
32  //Delete allocated memory. It is allocated every time runAlgorithm is called
33  fjClusterSeq_.reset();
34 }
35 
36 //______________________________________________________________________________
38 {
39  // run algorithm
40  if ( !doAreaFastjet_ && !doRhoFastjet_) {
41  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
42  } else if (voronoiRfact_ <= 0) {
43  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
44  } else {
45  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
46  }
47 
48  fjJets_.clear();
49  std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
50 
51  //Get local rho and rhom map
55 
56  iEvent.getByToken(etaToken_, etaRanges);
57  iEvent.getByToken(rhoToken_, rhoRanges);
58  iEvent.getByToken(rhomToken_, rhomRanges);
59 
60  //Check if size of eta and background density vectors is the same
61  unsigned int bkgVecSize = etaRanges->size();
62  if(bkgVecSize<1) { throw cms::Exception("WrongBkgInput") << "Producer needs vector with background estimates\n"; }
63  if(bkgVecSize != (rhoRanges->size()+1) || bkgVecSize != (rhomRanges->size()+1) ) {
64  throw cms::Exception("WrongBkgInput") << "Size of etaRanges (" << bkgVecSize << ") and rhoRanges (" << rhoRanges->size() << ") and/or rhomRanges (" << rhomRanges->size() << ") vectors inconsistent\n";
65  }
66 
67 
68  //Allow the background densities to change within the jet
69 
70  for(fastjet::PseudoJet& ijet : tempJets ) {
71 
72  //----------------------------------------------------------------------
73  // sift ghosts and particles in the input jet
74  std::vector<fastjet::PseudoJet> particles, ghosts;
75  fastjet::SelectorIsPureGhost().sift(ijet.constituents(), ghosts, particles);
76  unsigned long nParticles=particles.size();
77  if(nParticles==0) continue; //don't subtract ghost jets
78 
79  //assign rho and rhom to ghosts according to local eta-dependent map
80  double rho = 1e-6;
81  double rhom = 1e-6;
82  for(fastjet::PseudoJet& ighost : ghosts) {
83 
84  if(ighost.eta()<=etaRanges->at(0) || bkgVecSize==1) {
85  rho = rhoRanges->at(0);
86  rhom = rhomRanges->at(0);
87  } else if(ighost.eta()>=etaRanges->at(bkgVecSize-1)) {
88  rho = rhoRanges->at(bkgVecSize-2);
89  rhom = rhomRanges->at(bkgVecSize-2);
90  } else {
91  for(unsigned int ie = 0; ie<(bkgVecSize-1); ie++) {
92  if(ighost.eta()>=etaRanges->at(ie) && ighost.eta()<etaRanges->at(ie+1)) {
93  rho = rhoRanges->at(ie);
94  rhom = rhomRanges->at(ie);
95  break;
96  }
97  }
98  }
99  double pt = rho*ighost.area();
100  double mass_squared=std::pow(rhom*ighost.area()+pt,2)-std::pow(pt,2);
101  double mass=0;
102  if (mass_squared>0) mass=sqrt(mass_squared);
103  ighost.reset_momentum_PtYPhiM(pt,ighost.rap(),ighost.phi(),mass);
104  }
105 
106  //----------------------------------------------------------------------
107  //from here use official fastjet contrib package
108 
109  fastjet::contrib::ConstituentSubtractor subtractor;
110  subtractor.set_distance_type(fastjet::contrib::ConstituentSubtractor::deltaR); // distance in eta-phi plane
111  subtractor.set_max_distance(csRParam_); // free parameter for the maximal allowed distance between particle i and ghost k
112  subtractor.set_alpha(csAlpha_); // free parameter for the distance measure (the exponent of particle pt). Note that in older versions of the package alpha was multiplied by two but in newer versions this is not the case anymore
113  subtractor.set_do_mass_subtraction(true);
114 
115  std::vector<fastjet::PseudoJet> subtracted_particles = subtractor.do_subtraction(particles,ghosts);
116 
117  //Create subtracted jets
118  fastjet::PseudoJet subtracted_jet=join(subtracted_particles);
119  if(subtracted_jet.perp()>0.)
120  fjJets_.push_back( subtracted_jet );
121  }
122  fjJets_ = fastjet::sorted_by_pt(fjJets_);
123 }
124 
126 
127  edm::ParameterSetDescription descCSJetProducer;
129  fillDescriptionsFromCSJetProducer(descCSJetProducer);
131  descCSJetProducer.add<string>("jetCollInstanceName", "" );
133  descCSJetProducer.add<bool> ("sumRecHits", false);
134 
136  descriptions.add("CSJetProducer",descCSJetProducer);
137 
138 }
139 
141 
142  desc.add<double>("csRParam",-1.);
143  desc.add<double>("csAlpha",2.);
144 
145  desc.add<edm::InputTag>("etaMap",edm::InputTag("hiFJRhoProducer","mapEtaEdges") );
146  desc.add<edm::InputTag>("rho",edm::InputTag("hiFJRhoProducer","mapToRho") );
147  desc.add<edm::InputTag>("rhom",edm::InputTag("hiFJRhoProducer","mapToRhoM") );
148 
149 }
150 
152 // define as cmssw plugin
154 
T getParameter(std::string const &) const
std::vector< fastjet::PseudoJet > fjJets_
double csAlpha_
for constituent subtraction : R parameter
Definition: CSJetProducer.h:41
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< std::vector< double > > rhoToken_
Definition: CSJetProducer.h:45
std::vector< fastjet::PseudoJet > fjInputs_
int iEvent
Definition: GenABIO.cc:230
ClusterSequencePtr fjClusterSeq_
T sqrt(T t)
Definition: SSEVec.h:18
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
edm::EDGetTokenT< std::vector< double > > etaToken_
for HI constituent subtraction : alpha (power of pt in metric)
Definition: CSJetProducer.h:44
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< std::vector< double > > rhomToken_
Definition: CSJetProducer.h:46
fixed size matrix
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
HLT enums.
AreaDefinitionPtr fjAreaDefinition_
static void fillDescriptionsFromVirtualJetProducer(edm::ParameterSetDescription &desc)
static void fillDescriptionsFromCSJetProducer(edm::ParameterSetDescription &desc)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr