CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HiFJRhoFlowModulationProducer Class Reference
Inheritance diagram for HiFJRhoFlowModulationProducer:
edm::stream::EDProducer<>

Public Member Functions

 HiFJRhoFlowModulationProducer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const bool doEvtPlane_
 
const bool doFreePlaneFit_
 
const bool doJettyExclusion_
 
const int evtPlaneLevel_
 
const edm::EDGetTokenT< reco::EvtPlaneCollectionevtPlaneToken_
 
std::unique_ptr< TF1 > flowFit_p_
 
const edm::EDGetTokenT< reco::JetViewjetToken_
 
std::unique_ptr< TF1 > lineFit_p_
 
const edm::EDGetTokenT< reco::PFCandidateCollectionpfCandsToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 45 of file HiFJRhoFlowModulationProducer.cc.

Constructor & Destructor Documentation

◆ HiFJRhoFlowModulationProducer()

HiFJRhoFlowModulationProducer::HiFJRhoFlowModulationProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 63 of file HiFJRhoFlowModulationProducer.cc.

References flowFit_p_, lineFit_p_, and Pi.

64  : doEvtPlane_(iConfig.getParameter<bool>("doEvtPlane")),
65  doFreePlaneFit_(iConfig.getParameter<bool>("doFreePlaneFit")),
66  doJettyExclusion_(iConfig.getParameter<bool>("doJettyExclusion")),
67  evtPlaneLevel_(iConfig.getParameter<int>("evtPlaneLevel")),
68  jetToken_(doJettyExclusion_ ? consumes<reco::JetView>(iConfig.getParameter<edm::InputTag>("jetTag"))
70  pfCandsToken_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfCandSource"))),
71  evtPlaneToken_(consumes<reco::EvtPlaneCollection>(iConfig.getParameter<edm::InputTag>("EvtPlane"))) {
72  produces<std::vector<double>>("rhoFlowFitParams");
73  TMinuitMinimizer::UseStaticMinuit(false);
74  lineFit_p_ = std::make_unique<TF1>("lineFit", lineFunction, -TMath::Pi(), TMath::Pi());
75  flowFit_p_ = std::make_unique<TF1>("flowFit", flowFunction, -TMath::Pi(), TMath::Pi());
76 }
const double Pi
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::EDGetTokenT< reco::EvtPlaneCollection > evtPlaneToken_
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandsToken_
const edm::EDGetTokenT< reco::JetView > jetToken_

Member Function Documentation

◆ fillDescriptions()

void HiFJRhoFlowModulationProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 205 of file HiFJRhoFlowModulationProducer.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, and HLT_2022v15_cff::InputTag.

205  {
207  desc.add<bool>("doEvtPlane", false);
208  desc.add<edm::InputTag>("EvtPlane", edm::InputTag("hiEvtPlane"));
209  desc.add<bool>("doJettyExclusion", false);
210  desc.add<bool>("doFreePlaneFit", false);
211  desc.add<bool>("doFlatTest", false);
212  desc.add<edm::InputTag>("jetTag", edm::InputTag("ak4PFJets"));
213  desc.add<edm::InputTag>("pfCandSource", edm::InputTag("particleFlow"));
214  desc.add<int>("evtPlaneLevel", 0);
215  descriptions.add("hiFJRhoFlowModulationProducer", desc);
216 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

void HiFJRhoFlowModulationProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 79 of file HiFJRhoFlowModulationProducer.cc.

References funct::abs(), cms::cuda::assert(), funct::cos(), HLTMuonOfflineAnalyzer_cfi::deltaR2, doEvtPlane_, doFreePlaneFit_, doJettyExclusion_, MillePedeFileConverter_cfg::e, evtPlaneLevel_, evtPlaneToken_, flowFit_p_, hi::HF2, hi::HF3, HiEvtPlane_cfi::hiEvtPlane, iEvent, metsig::jet, PDWG_EXODelayedJetMET_cff::jets, jetToken_, lineFit_p_, SiStripPI::max, eostools::move(), Skims_PA_cff::name, ecaldqm::binning::nPhiBins, hi::NumEPNames, pfCandsToken_, Pi, funct::sin(), AlCaHLTBitMon_QueryRunRegistry::string, cond::impl::to_string(), and HcalDetIdTransform::transform().

79  {
80  int nPhiBins = 10;
81  auto const& pfCands = iEvent.get(pfCandsToken_);
82  static constexpr int kMaxEvtPlane = 29;
83  std::array<float, kMaxEvtPlane> hiEvtPlane;
84 
85  if (doEvtPlane_) {
86  auto const& evtPlanes = iEvent.get(evtPlaneToken_);
87  assert(evtPlanes.size() == hi::NumEPNames);
88  std::transform(evtPlanes.begin(), evtPlanes.end(), hiEvtPlane.begin(), [this](auto const& ePlane) -> float {
89  return ePlane.angle(evtPlaneLevel_);
90  });
91  }
92 
95  iEvent.getByToken(jetToken_, jets);
96 
97  int nFill = 0;
98 
99  double eventPlane2 = -100;
100  double eventPlane3 = -100;
101 
102  constexpr int nParamVals = 9;
103  auto rhoFlowFitParamsOut = std::make_unique<std::vector<double>>(nParamVals, 1e-6);
104 
105  rhoFlowFitParamsOut->at(0) = 0;
106  rhoFlowFitParamsOut->at(1) = 0;
107  rhoFlowFitParamsOut->at(2) = 0;
108  rhoFlowFitParamsOut->at(3) = 0;
109  rhoFlowFitParamsOut->at(4) = 0;
110 
111  double eventPlane2Cos = 0;
112  double eventPlane2Sin = 0;
113 
114  double eventPlane3Cos = 0;
115  double eventPlane3Sin = 0;
116 
117  std::vector<bool> pfcuts(pfCands.size(), false);
118  int iCand = -1;
119  for (auto const& pfCandidate : pfCands) {
120  iCand++;
121  if (pfCandidate.particleId() != 1 || std::abs(pfCandidate.eta()) > 1.0 || pfCandidate.pt() < .3 ||
122  pfCandidate.pt() > 3.) {
123  continue;
124  }
125 
126  if (doJettyExclusion_) {
127  bool isGood = true;
128  for (auto const& jet : *jets) {
129  if (deltaR2(jet, pfCandidate) < .16) {
130  isGood = false;
131  break;
132  }
133  }
134  if (!isGood) {
135  continue;
136  }
137  }
138 
139  nFill++;
140  pfcuts[iCand] = true;
141 
142  if (!doEvtPlane_) {
143  eventPlane2Cos += std::cos(2 * pfCandidate.phi());
144  eventPlane2Sin += std::sin(2 * pfCandidate.phi());
145 
146  eventPlane3Cos += std::cos(3 * pfCandidate.phi());
147  eventPlane3Sin += std::sin(3 * pfCandidate.phi());
148  }
149  }
150 
151  if (!doEvtPlane_) {
152  eventPlane2 = std::atan2(eventPlane2Sin, eventPlane2Cos) / 2.;
153  eventPlane3 = std::atan2(eventPlane3Sin, eventPlane3Cos) / 3.;
154  } else {
155  eventPlane2 = hiEvtPlane[hi::HF2];
156  eventPlane3 = hiEvtPlane[hi::HF3];
157  }
158  int pfcuts_count = 0;
159  if (nFill >= 100 && eventPlane2 > -99) {
160  nPhiBins = std::max(10, nFill / 30);
161 
162  std::string name = "phiTestIEta4_" + std::to_string(iEvent.id().event()) + "_h";
163  std::string nameFlat = "phiTestIEta4_Flat_" + std::to_string(iEvent.id().event()) + "_h";
164  std::unique_ptr<TH1F> phi_h = std::make_unique<TH1F>(name.data(), "", nPhiBins, -TMath::Pi(), TMath::Pi());
165  phi_h->SetDirectory(nullptr);
166  for (auto const& pfCandidate : pfCands) {
167  if (pfcuts.at(pfcuts_count))
168  phi_h->Fill(pfCandidate.phi());
169  pfcuts_count++;
170  }
171  flowFit_p_->SetParameter(0, 10);
172  flowFit_p_->SetParameter(1, 0);
173  flowFit_p_->SetParameter(2, eventPlane2);
174  flowFit_p_->SetParameter(3, 0);
175  flowFit_p_->SetParameter(4, eventPlane3);
176  if (!doFreePlaneFit_) {
177  flowFit_p_->FixParameter(2, eventPlane2);
178  flowFit_p_->FixParameter(4, eventPlane3);
179  }
180 
181  lineFit_p_->SetParameter(0, 10);
182 
183  phi_h->Fit(flowFit_p_.get(), "Q SERIAL", "", -TMath::Pi(), TMath::Pi());
184  phi_h->Fit(lineFit_p_.get(), "Q SERIAL", "", -TMath::Pi(), TMath::Pi());
185  rhoFlowFitParamsOut->at(0) = flowFit_p_->GetParameter(0);
186  rhoFlowFitParamsOut->at(1) = flowFit_p_->GetParameter(1);
187  rhoFlowFitParamsOut->at(2) = flowFit_p_->GetParameter(2);
188  rhoFlowFitParamsOut->at(3) = flowFit_p_->GetParameter(3);
189  rhoFlowFitParamsOut->at(4) = flowFit_p_->GetParameter(4);
190 
191  rhoFlowFitParamsOut->at(5) = flowFit_p_->GetChisquare();
192  rhoFlowFitParamsOut->at(6) = flowFit_p_->GetNDF();
193 
194  rhoFlowFitParamsOut->at(7) = lineFit_p_->GetChisquare();
195  rhoFlowFitParamsOut->at(8) = lineFit_p_->GetNDF();
196 
197  phi_h.reset();
198  pfcuts.clear();
199  }
200 
201  iEvent.put(std::move(rhoFlowFitParamsOut), "rhoFlowFitParams");
202 }
const double Pi
const edm::EDGetTokenT< reco::EvtPlaneCollection > evtPlaneToken_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::string to_string(const V &value)
Definition: OMSAccess.h:71
assert(be >=bs)
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandsToken_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::EDGetTokenT< reco::JetView > jetToken_
static const int NumEPNames
def move(src, dest)
Definition: eostools.py:511
unsigned transform(const HcalDetId &id, unsigned transformCode)

Member Data Documentation

◆ doEvtPlane_

const bool HiFJRhoFlowModulationProducer::doEvtPlane_
private

Definition at line 53 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().

◆ doFreePlaneFit_

const bool HiFJRhoFlowModulationProducer::doFreePlaneFit_
private

Definition at line 54 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().

◆ doJettyExclusion_

const bool HiFJRhoFlowModulationProducer::doJettyExclusion_
private

Definition at line 55 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().

◆ evtPlaneLevel_

const int HiFJRhoFlowModulationProducer::evtPlaneLevel_
private

Definition at line 56 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().

◆ evtPlaneToken_

const edm::EDGetTokenT<reco::EvtPlaneCollection> HiFJRhoFlowModulationProducer::evtPlaneToken_
private

Definition at line 59 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().

◆ flowFit_p_

std::unique_ptr<TF1> HiFJRhoFlowModulationProducer::flowFit_p_
private

Definition at line 61 of file HiFJRhoFlowModulationProducer.cc.

Referenced by HiFJRhoFlowModulationProducer(), and produce().

◆ jetToken_

const edm::EDGetTokenT<reco::JetView> HiFJRhoFlowModulationProducer::jetToken_
private

Definition at line 57 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().

◆ lineFit_p_

std::unique_ptr<TF1> HiFJRhoFlowModulationProducer::lineFit_p_
private

Definition at line 60 of file HiFJRhoFlowModulationProducer.cc.

Referenced by HiFJRhoFlowModulationProducer(), and produce().

◆ pfCandsToken_

const edm::EDGetTokenT<reco::PFCandidateCollection> HiFJRhoFlowModulationProducer::pfCandsToken_
private

Definition at line 58 of file HiFJRhoFlowModulationProducer.cc.

Referenced by produce().