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