CMS 3D CMS Logo

HiFJRhoProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HiJetBackground/HiFJRhoProducer
4 // Class: HiFJRhoProducer
5 //
13 //
14 // Original Author: Marta Verweij
15 // Created: Thu, 16 Jul 2015 10:57:12 GMT
16 //
17 //
18 
19 #include <memory>
20 #include <string>
21 
23 
28 
34 
35 //using namespace edm;
36 using namespace reco;
37 //using namespace pat;
38 
39 //
40 // constants, enums and typedefs
41 //
42 
43 //
44 // static data member definitions
45 //
46 
47 //
48 // constructors and destructor
49 //
51  : src_(iConfig.getParameter<edm::InputTag>("jetSource")),
52  nExcl_(iConfig.getParameter<int>("nExcl")),
53  etaMaxExcl_(iConfig.getParameter<double>("etaMaxExcl")),
54  ptMinExcl_(iConfig.getParameter<double>("ptMinExcl")),
55  nExcl2_(iConfig.getParameter<int>("nExcl2")),
56  etaMaxExcl2_(iConfig.getParameter<double>("etaMaxExcl2")),
57  ptMinExcl2_(iConfig.getParameter<double>("ptMinExcl2")),
58  checkJetCand(true),
59  usingPackedCand(false) {
60  jetsToken_ = consumes<edm::View<reco::Jet>>(src_);
61 
62  //register your products
63  produces<std::vector<double>>("mapEtaEdges");
64  produces<std::vector<double>>("mapToRho");
65  produces<std::vector<double>>("mapToRhoM");
66  produces<std::vector<double>>("ptJets");
67  produces<std::vector<double>>("areaJets");
68  produces<std::vector<double>>("etaJets");
69  etaRanges = iConfig.getParameter<std::vector<double>>("etaRanges");
70 }
71 
73  // do anything here that needs to be done at desctruction time
74  // (e.g. close files, deallocate resources etc.)
75 }
76 
77 // ------------ method called to produce the data ------------
79  // Get the vector of jets
81  iEvent.getByToken(jetsToken_, jets);
82 
83  int neta = (int)etaRanges.size();
84  auto mapEtaRangesOut = std::make_unique<std::vector<double>>(neta, -999.);
85 
86  for (int ieta = 0; ieta < neta; ieta++) {
87  mapEtaRangesOut->at(ieta) = etaRanges[ieta];
88  }
89  auto mapToRhoOut = std::make_unique<std::vector<double>>(neta - 1, 1e-6);
90  auto mapToRhoMOut = std::make_unique<std::vector<double>>(neta - 1, 1e-6);
91 
92  int njets = jets->size();
93 
94  auto ptJetsOut = std::make_unique<std::vector<double>>(njets, 1e-6);
95  auto areaJetsOut = std::make_unique<std::vector<double>>(njets, 1e-6);
96  auto etaJetsOut = std::make_unique<std::vector<double>>(njets, 1e-6);
97 
98  double rhoVec[999];
99  double rhomVec[999];
100  double etaVec[999];
101 
102  // int neta = (int)mapEtaRangesOut->size();
103  int nacc = 0;
104  unsigned int njetsEx = 0;
105  unsigned int njetsEx2 = 0;
106  for (auto jet = jets->begin(); jet != jets->end(); ++jet) {
107  if (njetsEx < nExcl_ && fabs(jet->eta()) < etaMaxExcl_ && jet->pt() > ptMinExcl_) {
108  njetsEx++;
109  continue;
110  }
111  if (njetsEx2 < nExcl2_ && fabs(jet->eta()) < etaMaxExcl2_ && fabs(jet->eta()) > etaMaxExcl_ &&
112  jet->pt() > ptMinExcl2_) {
113  njetsEx2++;
114  continue;
115  }
116  float pt = jet->pt();
117  float area = jet->jetArea();
118  float eta = jet->eta();
119 
120  if (eta < mapEtaRangesOut->at(0) || eta > mapEtaRangesOut->at(neta - 1))
121  continue;
122  if (area > 0.) {
123  rhoVec[nacc] = pt / area;
124  rhomVec[nacc] = calcMd(&*jet) / area;
125  etaVec[nacc] = eta;
126  ptJetsOut->at(nacc) = pt;
127  areaJetsOut->at(nacc) = area;
128  etaJetsOut->at(nacc) = eta;
129  ++nacc;
130  }
131  }
132 
133  ptJetsOut->resize(nacc);
134  areaJetsOut->resize(nacc);
135  etaJetsOut->resize(nacc);
136  //calculate rho and rhom in eta ranges
137  double radius = 0.2; //distance kt clusters needs to be from edge
138  for (int ieta = 0; ieta < (neta - 1); ieta++) {
139  std::vector<double> rhoVecCur;
140  std::vector<double> rhomVecCur;
141 
142  double etaMin = mapEtaRangesOut->at(ieta) + radius;
143  double etaMax = mapEtaRangesOut->at(ieta + 1) - radius;
144 
145  int naccCur = 0;
146  double rhoCurSum = 0.;
147  double rhomCurSum = 0.;
148  for (int i = 0; i < nacc; i++) {
149  if (etaVec[i] >= etaMin && etaVec[i] < etaMax) {
150  rhoVecCur.push_back(rhoVec[i]);
151  rhomVecCur.push_back(rhomVec[i]);
152 
153  rhoCurSum += rhoVec[i];
154  rhomCurSum += rhomVec[i];
155  ++naccCur;
156  } //eta selection
157  } //accepted jet loop
158 
159  if (naccCur > 0) {
160  double rhoCur = calcMedian(rhoVecCur);
161  double rhomCur = calcMedian(rhomVecCur);
162  mapToRhoOut->at(ieta) = rhoCur;
163  mapToRhoMOut->at(ieta) = rhomCur;
164  }
165  } //eta ranges
166 
167  iEvent.put(std::move(mapEtaRangesOut), "mapEtaEdges");
168  iEvent.put(std::move(mapToRhoOut), "mapToRho");
169  iEvent.put(std::move(mapToRhoMOut), "mapToRhoM");
170 
171  iEvent.put(std::move(ptJetsOut), "ptJets");
172  iEvent.put(std::move(areaJetsOut), "areaJets");
173  iEvent.put(std::move(etaJetsOut), "etaJets");
174 
175  return;
176 }
177 
178 // ------------ method called once each stream before processing any runs, lumis or events ------------
180 
181 // ------------ method called once each stream after processing all runs, lumis and events ------------
183 
185  //
186  //get md as defined in http://arxiv.org/pdf/1211.2811.pdf
187  //
188 
189  //Loop over the jet constituents
190  double sum = 0.;
191  for (auto daughter : jet->getJetConstituentsQuick()) {
192  if (isPackedCandidate(daughter)) { //packed candidate situation
193  auto part = static_cast<const pat::PackedCandidate*>(daughter);
194  sum += sqrt(part->mass() * part->mass() + part->pt() * part->pt()) - part->pt();
195  } else {
196  auto part = static_cast<const reco::PFCandidate*>(daughter);
197  sum += sqrt(part->mass() * part->mass() + part->pt() * part->pt()) - part->pt();
198  }
199  }
200 
201  return sum;
202 }
203 
205  if (checkJetCand) {
206  if (typeid(pat::PackedCandidate) == typeid(*candidate))
207  usingPackedCand = true;
208  else if (typeid(reco::PFCandidate) == typeid(*candidate))
209  usingPackedCand = false;
210  else
211  throw cms::Exception("WrongJetCollection", "Jet constituents are not particle flow candidates");
212  checkJetCand = false;
213  }
214  return usingPackedCand;
215 }
216 
219  desc.add<edm::InputTag>("jetSource", edm::InputTag("kt4PFJets"));
220  desc.add<int>("nExcl", 2);
221  desc.add<double>("etaMaxExcl", 2.);
222  desc.add<double>("ptMinExcl", 20.);
223  desc.add<int>("nExcl2", 2);
224  desc.add<double>("etaMaxExcl2", 2.);
225  desc.add<double>("ptMinExcl2", 20.);
226  desc.add<std::vector<double>>("etaRanges", {});
227  descriptions.add("hiFJRhoProducer", desc);
228 }
229 
230 //--------- method to calculate median ------------------
231 double HiFJRhoProducer::calcMedian(std::vector<double>& v) {
232  //post-condition: After returning, the elements in v may be reordered and the resulting order is implementation defined.
233  //works for even and odd collections
234  if (v.empty()) {
235  return 0.0;
236  }
237  auto n = v.size() / 2;
238  std::nth_element(v.begin(), v.begin() + n, v.end());
239  auto med = v[n];
240  if (!(v.size() & 1)) { //If the set size is even
241  auto max_it = std::max_element(v.begin(), v.begin() + n);
242  med = (*max_it + med) / 2.0;
243  }
244  return med;
245 }
246 
247 //define this as a plug-in
edm::StreamID
Definition: StreamID.h:30
neta
const int neta
Definition: CMTRawAnalyzer.h:423
Handle.h
ALCARECOTkAlBeamHalo_cff.etaMin
etaMin
GeV.
Definition: ALCARECOTkAlBeamHalo_cff.py:32
mps_fire.i
i
Definition: mps_fire.py:428
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
HiFJRhoProducer.h
HiFJRhoProducer::ptMinExcl2_
double ptMinExcl2_
Definition: HiFJRhoProducer.h:53
funct::false
false
Definition: Factorize.h:29
HiFJRhoProducer::HiFJRhoProducer
HiFJRhoProducer(const edm::ParameterSet &)
Definition: HiFJRhoProducer.cc:50
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
ESHandle.h
PFCandidate.h
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89281
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
HiFJRhoProducer::calcMd
double calcMd(const reco::Jet *jet)
Definition: HiFJRhoProducer.cc:184
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
HiFJRhoProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: HiFJRhoProducer.cc:78
findQualityFiles.v
v
Definition: findQualityFiles.py:179
edm::Handle
Definition: AssociativeIterator.h:50
HiFJRhoProducer::src_
edm::InputTag src_
Definition: HiFJRhoProducer.h:47
HiFJRhoProducer::endStream
void endStream() override
Definition: HiFJRhoProducer.cc:182
MakerMacros.h
part
part
Definition: HCALResponse.h:20
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
HiFJRhoProducer::usingPackedCand
bool usingPackedCand
Definition: HiFJRhoProducer.h:55
PVValHelper::eta
Definition: PVValidationHelpers.h:70
HiFJRhoProducer::nExcl2_
unsigned int nExcl2_
Definition: HiFJRhoProducer.h:51
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
HiFJRhoProducer::calcMedian
double calcMedian(std::vector< double > &v)
Definition: HiFJRhoProducer.cc:231
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
funct::true
true
Definition: Factorize.h:173
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
HiFJRhoProducer::etaRanges
std::vector< double > etaRanges
Definition: HiFJRhoProducer.h:54
edm::ParameterSet
Definition: ParameterSet.h:47
HiFJRhoProducer::etaMaxExcl2_
double etaMaxExcl2_
Definition: HiFJRhoProducer.h:52
HiFJRhoProducer::isPackedCandidate
bool isPackedCandidate(const reco::Candidate *candidate)
Definition: HiFJRhoProducer.cc:204
pat::PackedCandidate
Definition: PackedCandidate.h:22
PackedCandidate.h
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
HiFJRhoProducer::etaMaxExcl_
double etaMaxExcl_
Definition: HiFJRhoProducer.h:49
edm::EventSetup
Definition: EventSetup.h:58
Jet.h
HiFJRhoProducer::ptMinExcl_
double ptMinExcl_
Definition: HiFJRhoProducer.h:50
HiFJRhoProducer::~HiFJRhoProducer
~HiFJRhoProducer() override
Definition: HiFJRhoProducer.cc:72
reco::Candidate
Definition: Candidate.h:27
HiFJRhoProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HiFJRhoProducer.cc:217
HiFJRhoProducer::jetsToken_
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken_
Definition: HiFJRhoProducer.h:44
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
HiFJRhoProducer::nExcl_
unsigned int nExcl_
Definition: HiFJRhoProducer.h:48
eostools.move
def move(src, dest)
Definition: eostools.py:511
HiFJRhoProducer::checkJetCand
bool checkJetCand
Definition: HiFJRhoProducer.h:55
BTaggingMonitoring_cff.njets
njets
Definition: BTaggingMonitoring_cff.py:10
metsig::jet
Definition: SignAlgoResolutions.h:47
Exception
Definition: hltDiff.cc:245
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
ALCARECOTkAlBeamHalo_cff.etaMax
etaMax
Definition: ALCARECOTkAlBeamHalo_cff.py:33
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HiFJRhoProducer
Definition: HiFJRhoProducer.h:26
HiFJRhoProducer::beginStream
void beginStream(edm::StreamID) override
Definition: HiFJRhoProducer.cc:179
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
ConsumesCollector.h
View.h
ParameterSet.h
edm::Event
Definition: Event.h:73
custom_jme_cff.area
area
Definition: custom_jme_cff.py:140
edm::InputTag
Definition: InputTag.h:15
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37