CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
ResolutionCreator Class Reference
Inheritance diagram for ResolutionCreator:
edm::one::EDAnalyzer< edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 ResolutionCreator (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::SharedResources >
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const noexcept final
 
bool wantsGlobalRuns () const noexcept final
 
bool wantsInputProcessBlocks () const noexcept final
 
bool wantsProcessBlocks () const noexcept final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const noexcept
 
bool wantsStreamRuns () const noexcept
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 

Private Attributes

edm::EDGetTokenT< std::vector< pat::Electron > > electronsToken_
 
std::vector< double > etabinVals_
 
int etanrbins
 
TF1 * fResEtaBin [10][20]
 
TF1 * fResPtEtaBin [10][20][20]
 
edm::EDGetTokenT< TtGenEventgenEvtToken_
 
TH1F * hResEtaBin [10][20]
 
TH1F * hResPtEtaBin [10][20][20]
 
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
 
std::string labelName_
 
int matchingAlgo_
 
double maxDist_
 
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
 
double minDR_
 
edm::EDGetTokenT< std::vector< pat::Muon > > muonsToken_
 
int nrFilled
 
std::string objectType_
 
std::vector< double > pTbinVals_
 
int ptnrbins
 
edm::EDGetTokenT< std::vector< pat::Tau > > tausToken_
 
TTree * tResVar
 
bool useDeltaR_
 
bool useMaxDist_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 38 of file ResolutionCreator.cc.

Constructor & Destructor Documentation

◆ ResolutionCreator()

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

Definition at line 75 of file ResolutionCreator.cc.

References electronsToken_, etabinVals_, genEvtToken_, edm::ParameterSet::getParameter(), ProducerED_cfi::InputTag, jetsToken_, labelName_, metsToken_, minDR_, muonsToken_, nrFilled, objectType_, pTbinVals_, AlCaHLTBitMon_QueryRunRegistry::string, and tausToken_.

75  {
76  usesResource("TFileService");
77 
78  // input parameters
79  genEvtToken_ = consumes<TtGenEvent>(edm::InputTag("genEvt"));
80  objectType_ = iConfig.getParameter<std::string>("object");
81  labelName_ = iConfig.getParameter<std::string>("label");
82  if (objectType_ == "electron")
83  electronsToken_ = consumes<std::vector<pat::Electron> >(edm::InputTag(labelName_));
84  else if (objectType_ == "muon")
85  muonsToken_ = consumes<std::vector<pat::Muon> >(edm::InputTag(labelName_));
86  else if (objectType_ == "lJets" || objectType_ == "bJets")
87  jetsToken_ = consumes<std::vector<pat::Jet> >(edm::InputTag(labelName_));
88  else if (objectType_ == "met")
89  metsToken_ = consumes<std::vector<pat::MET> >(edm::InputTag(labelName_));
90  else if (objectType_ == "tau")
91  tausToken_ = consumes<std::vector<pat::Tau> >(edm::InputTag(labelName_));
92  if (objectType_ != "met") {
93  etabinVals_ = iConfig.getParameter<std::vector<double> >("etabinValues");
94  }
95  pTbinVals_ = iConfig.getParameter<std::vector<double> >("pTbinValues");
96  minDR_ = iConfig.getParameter<double>("minMatchingDR");
97 
98  nrFilled = 0;
99 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< TtGenEvent > genEvtToken_
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
edm::EDGetTokenT< std::vector< pat::Tau > > tausToken_
std::vector< double > etabinVals_
edm::EDGetTokenT< std::vector< pat::Electron > > electronsToken_
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
std::vector< double > pTbinVals_
edm::EDGetTokenT< std::vector< pat::Muon > > muonsToken_

Member Function Documentation

◆ analyze()

void ResolutionCreator::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 106 of file ResolutionCreator.cc.

References funct::abs(), b, bgen, electronAnalyzer_cfi::DeltaR, MillePedeFileConverter_cfg::e, pwdgSkimBPark_cfi::electrons, electronsToken_, reco::LeafCandidate::energy(), etabinVals_, etanrbins, TtGenEvtProducer_cfi::genEvt, genEvtToken_, hResPtEtaBin, iEvent, dqmiolumiharvest::j, PDWG_EXODelayedJetMET_cff::jets, jetsToken_, visualization-live-secondInstance_cfg::m, singleTopDQM_cfi::mets, metsToken_, minDR_, reco::CompositeRefCandidateT< D >::mother(), reco::Candidate::mother(), DiMuonV_cfg::muons, muonsToken_, nrFilled, reco::CompositeRefCandidateT< D >::numberOfMothers(), reco::Candidate::numberOfMothers(), objectType_, AlCaHLTBitMon_ParallelJobs::p, reco::LeafCandidate::p4(), reco::Candidate::pdgId(), reco::LeafCandidate::pdgId(), hitfit::phidiff(), DiDispStaMuonMonitor_cfi::pt, pTbinVals_, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, reco::LeafCandidate::status(), metsig::tau, Tau3MuMonitor_cff::taus, and tausToken_.

106  {
107  // Get the gen and cal object fourvector
108  std::vector<reco::Particle *> p4gen, p4rec;
109 
111  iEvent.getByToken(genEvtToken_, genEvt);
112 
113  if (genEvt->particles().size() < 10)
114  return;
115 
116  if (objectType_ == "electron") {
118  electrons; //to calculate the ResolutionCreator for the electrons, i used the TopElectron instead of the AOD information
119  iEvent.getByToken(electronsToken_, electrons);
120  for (size_t e = 0; e < electrons->size(); e++) {
121  for (size_t p = 0; p < genEvt->particles().size(); p++) {
122  if ((std::abs(genEvt->particles()[p].pdgId()) == 11) &&
123  (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*electrons)[e].p4()) < minDR_)) {
124  //p4gen.push_back(new reco::Particle(genEvt->particles()[p]));
125  //p4rec.push_back(new reco::Particle((pat::Electron)((*electrons)[e])));
126  }
127  }
128  }
129  } else if (objectType_ == "muon") {
131  iEvent.getByToken(muonsToken_, muons);
132  for (size_t m = 0; m < muons->size(); m++) {
133  for (size_t p = 0; p < genEvt->particles().size(); p++) {
134  if ((std::abs(genEvt->particles()[p].pdgId()) == 13) &&
135  (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*muons)[m].p4()) < minDR_)) {
136  //p4gen.push_back(new reco::Particle(genEvt->particles()[p]));
137  //p4rec.push_back(new reco::Particle((pat::Muon)((*muons)[m])));
138  }
139  }
140  }
141  } else if (objectType_ == "lJets") {
143  iEvent.getByToken(jetsToken_, jets);
144  if (jets->size() >= 4) {
145  for (unsigned int j = 0; j < 4; j++) {
146  for (size_t p = 0; p < genEvt->particles().size(); p++) {
147  if ((std::abs(genEvt->particles()[p].pdgId()) < 5) &&
148  (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*jets)[j].p4()) < minDR_)) {
149  //p4gen.push_back(new reco::Particle(genEvt->particles()[p]));
150  //p4rec.push_back(new reco::Particle((pat::Jet)(*jets)[j]));
151  }
152  }
153  }
154  }
155  } else if (objectType_ == "bJets") {
157  iEvent.getByToken(jetsToken_, jets);
158  if (jets->size() >= 4) {
159  for (unsigned int j = 0; j < 4; j++) {
160  for (size_t p = 0; p < genEvt->particles().size(); p++) {
161  if ((std::abs(genEvt->particles()[p].pdgId()) == 5) &&
162  (ROOT::Math::VectorUtil::DeltaR(genEvt->particles()[p].p4(), (*jets)[j].p4()) < minDR_)) {
163  //p4gen.push_back(new reco::Particle(genEvt->particles()[p]));
164  //p4rec.push_back(new reco::Particle((pat::Jet)(*jets)[j]));
165  }
166  }
167  }
168  }
169  } else if (objectType_ == "met") {
171  iEvent.getByToken(metsToken_, mets);
172  if (!mets->empty()) {
173  if (genEvt->isSemiLeptonic() && genEvt->singleNeutrino() != nullptr &&
174  ROOT::Math::VectorUtil::DeltaR(genEvt->singleNeutrino()->p4(), (*mets)[0].p4()) < minDR_) {
175  //p4gen.push_back(new reco::Particle(0,genEvt->singleNeutrino()->p4(),math::XYZPoint()));
176  //p4rec.push_back(new reco::Particle((pat::MET)((*mets)[0])));
177  }
178  }
179  } else if (objectType_ == "tau") {
181  iEvent.getByToken(tausToken_, taus);
182  for (std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau != taus->end(); ++tau) {
183  // find the tau (if any) that matches a MC tau from W
184  reco::GenParticle genLepton = *(tau->genLepton());
185  if (std::abs(genLepton.pdgId()) == 15 && genLepton.status() == 2 && genLepton.numberOfMothers() > 0 &&
186  std::abs(genLepton.mother(0)->pdgId()) == 15 && genLepton.mother(0)->numberOfMothers() > 0 &&
187  std::abs(genLepton.mother(0)->mother(0)->pdgId()) == 24 &&
188  ROOT::Math::VectorUtil::DeltaR(genLepton.p4(), tau->p4()) < minDR_) {
189  }
190  //p4gen.push_back(new reco::Particle(genLepton));
191  //p4rec.push_back(new reco::Particle(*tau));
192  }
193  }
194  // Fill the object's value
195  for (unsigned m = 0; m < p4gen.size(); m++) {
196  double Egen = p4gen[m]->energy();
197  double Thetagen = p4gen[m]->theta();
198  double Phigen = p4gen[m]->phi();
199  double Etgen = p4gen[m]->et();
200  double Etagen = p4gen[m]->eta();
201  double Ecal = p4rec[m]->energy();
202  double Thetacal = p4rec[m]->theta();
203  double Phical = p4rec[m]->phi();
204  double Etcal = p4rec[m]->et();
205  double Etacal = p4rec[m]->eta();
206  double phidiff = Phical - Phigen;
207  if (phidiff > 3.14159)
208  phidiff = 2. * 3.14159 - phidiff;
209  if (phidiff < -3.14159)
210  phidiff = -phidiff - 2. * 3.14159;
211 
212  // find eta and et bin
213  int etabin = 0;
214  if (etanrbins > 1) {
215  for (unsigned int b = 0; b < etabinVals_.size() - 1; b++) {
216  if (fabs(Etacal) > etabinVals_[b])
217  etabin = b;
218  }
219  }
220 
221  int ptbin = 0;
222  for (unsigned int b = 0; b < pTbinVals_.size() - 1; b++) {
223  if (p4rec[m]->pt() > pTbinVals_[b])
224  ptbin = b;
225  }
226 
227  // calculate the resolution on "a", "b", "c" & "d" according to the definition (CMS-NOTE-2006-023):
228  // p = a*|p_meas|*u_1 + b*u_2 + c*u_3
229  // E(fit) = E_meas * d
230  //
231  // with u_1 = p/|p_meas|
232  // u_3 = (u_z x u_1)/|u_z x u_1|
233  // u_2 = (u_1 x u_3)/|u_1 x u_3|
234  //
235  // The initial parameters values are chosen like (a, b, c, d) = (1., 0., 0., 1.)
236 
237  // 1/ calculate the unitary vectors of the basis u_1, u_2, u_3
238  ROOT::Math::SVector<double, 3> pcalvec(p4rec[m]->px(), p4rec[m]->py(), p4rec[m]->pz());
239  ROOT::Math::SVector<double, 3> pgenvec(p4gen[m]->px(), p4gen[m]->py(), p4gen[m]->pz());
240 
241  ROOT::Math::SVector<double, 3> u_z(0, 0, 1);
242  ROOT::Math::SVector<double, 3> u_1 = ROOT::Math::Unit(pcalvec);
243  ROOT::Math::SVector<double, 3> u_3 = ROOT::Math::Cross(u_z, u_1) / ROOT::Math::Mag(ROOT::Math::Cross(u_z, u_1));
244  ROOT::Math::SVector<double, 3> u_2 = ROOT::Math::Cross(u_1, u_3) / ROOT::Math::Mag(ROOT::Math::Cross(u_1, u_3));
245  double acal = 1.;
246  double bcal = 0.;
247  double ccal = 0.;
248  double dcal = 1.;
249  double agen = ROOT::Math::Dot(pgenvec, u_1) / ROOT::Math::Mag(pcalvec);
250  double bgen = ROOT::Math::Dot(pgenvec, u_2);
251  double cgen = ROOT::Math::Dot(pgenvec, u_3);
252  double dgen = Egen / Ecal;
253 
254  //fill histograms
255  ++nrFilled;
256  hResPtEtaBin[0][etabin][ptbin]->Fill(acal - agen);
257  hResPtEtaBin[1][etabin][ptbin]->Fill(bcal - bgen);
258  hResPtEtaBin[2][etabin][ptbin]->Fill(ccal - cgen);
259  hResPtEtaBin[3][etabin][ptbin]->Fill(dcal - dgen);
260  hResPtEtaBin[4][etabin][ptbin]->Fill(Thetacal - Thetagen);
261  hResPtEtaBin[5][etabin][ptbin]->Fill(phidiff);
262  hResPtEtaBin[6][etabin][ptbin]->Fill(Etcal - Etgen);
263  hResPtEtaBin[7][etabin][ptbin]->Fill(Etacal - Etagen);
264 
265  delete p4gen[m];
266  delete p4rec[m];
267  }
268 }
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
int status() const final
status word
double phidiff(double phi)
Normalized difference in azimuthal angles to a range between .
Definition: fourvec.cc:220
edm::EDGetTokenT< TtGenEvent > genEvtToken_
size_t numberOfMothers() const override
number of mothers
double bgen
Definition: HydjetWrapper.h:46
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
const LorentzVector & p4() const final
four-momentum Lorentz vector
int pdgId() const final
PDG identifier.
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< std::vector< pat::Tau > > tausToken_
std::vector< double > etabinVals_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const =0
PDG identifier.
TH1F * hResPtEtaBin[10][20][20]
edm::EDGetTokenT< std::vector< pat::Electron > > electronsToken_
double b
Definition: hdecay.h:120
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
std::vector< double > pTbinVals_
edm::EDGetTokenT< std::vector< pat::Muon > > muonsToken_
double energy() const final
energy

◆ beginJob()

void ResolutionCreator::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 271 of file ResolutionCreator.cc.

References b, edm::errors::Configuration, defaults_cfi::etabins, etabinVals_, etanrbins, Exception, fResEtaBin, fResPtEtaBin, make_classfiles::fs, hResEtaBin, hResPtEtaBin, objectType_, LepHTMonitor_cff::ptbins, pTbinVals_, ptnrbins, and tResVar.

271  {
273  if (!fs)
274  throw edm::Exception(edm::errors::Configuration, "TFileService missing from configuration!");
275 
276  // input constants
277  TString resObsName[8] = {"_ares", "_bres", "_cres", "_dres", "_thres", "_phres", "_etres", "_etares"};
278  int resObsNrBins = 120;
279  if ((objectType_ == "muon") || (objectType_ == "electron"))
280  resObsNrBins = 80;
281  std::vector<double> resObsMin, resObsMax;
282  if (objectType_ == "electron") {
283  resObsMin.push_back(-0.15);
284  resObsMin.push_back(-0.2);
285  resObsMin.push_back(-0.1);
286  resObsMin.push_back(-0.15);
287  resObsMin.push_back(-0.0012);
288  resObsMin.push_back(-0.009);
289  resObsMin.push_back(-16);
290  resObsMin.push_back(-0.0012);
291  resObsMax.push_back(0.15);
292  resObsMax.push_back(0.2);
293  resObsMax.push_back(0.1);
294  resObsMax.push_back(0.15);
295  resObsMax.push_back(0.0012);
296  resObsMax.push_back(0.009);
297  resObsMax.push_back(16);
298  resObsMax.push_back(0.0012);
299  } else if (objectType_ == "muon") {
300  resObsMin.push_back(-0.15);
301  resObsMin.push_back(-0.1);
302  resObsMin.push_back(-0.05);
303  resObsMin.push_back(-0.15);
304  resObsMin.push_back(-0.004);
305  resObsMin.push_back(-0.003);
306  resObsMin.push_back(-8);
307  resObsMin.push_back(-0.004);
308  resObsMax.push_back(0.15);
309  resObsMax.push_back(0.1);
310  resObsMax.push_back(0.05);
311  resObsMax.push_back(0.15);
312  resObsMax.push_back(0.004);
313  resObsMax.push_back(0.003);
314  resObsMax.push_back(8);
315  resObsMax.push_back(0.004);
316  } else if (objectType_ == "tau") {
317  resObsMin.push_back(-1.);
318  resObsMin.push_back(-10.);
319  resObsMin.push_back(-10);
320  resObsMin.push_back(-1.);
321  resObsMin.push_back(-0.1);
322  resObsMin.push_back(-0.1);
323  resObsMin.push_back(-80);
324  resObsMin.push_back(-0.1);
325  resObsMax.push_back(1.);
326  resObsMax.push_back(10.);
327  resObsMax.push_back(10);
328  resObsMax.push_back(1.);
329  resObsMax.push_back(0.1);
330  resObsMax.push_back(0.1);
331  resObsMax.push_back(50);
332  resObsMax.push_back(0.1);
333  } else if (objectType_ == "lJets" || objectType_ == "bJets") {
334  resObsMin.push_back(-1.);
335  resObsMin.push_back(-10.);
336  resObsMin.push_back(-10.);
337  resObsMin.push_back(-1.);
338  resObsMin.push_back(-0.4);
339  resObsMin.push_back(-0.6);
340  resObsMin.push_back(-80);
341  resObsMin.push_back(-0.6);
342  resObsMax.push_back(1.);
343  resObsMax.push_back(10.);
344  resObsMax.push_back(10.);
345  resObsMax.push_back(1.);
346  resObsMax.push_back(0.4);
347  resObsMax.push_back(0.6);
348  resObsMax.push_back(80);
349  resObsMax.push_back(0.6);
350  } else {
351  resObsMin.push_back(-2.);
352  resObsMin.push_back(-150.);
353  resObsMin.push_back(-150.);
354  resObsMin.push_back(-2.);
355  resObsMin.push_back(-6);
356  resObsMin.push_back(-6);
357  resObsMin.push_back(-180);
358  resObsMin.push_back(-6);
359  resObsMax.push_back(3.);
360  resObsMax.push_back(150.);
361  resObsMax.push_back(150.);
362  resObsMax.push_back(3.);
363  resObsMax.push_back(6);
364  resObsMax.push_back(6);
365  resObsMax.push_back(180);
366  resObsMax.push_back(6);
367  }
368 
369  const char *resObsVsPtFit[8] = {"[0]+[1]*exp(-[2]*x)",
370  "[0]+[1]*exp(-[2]*x)",
371  "[0]+[1]*exp(-[2]*x)",
372  "[0]+[1]*exp(-[2]*x)",
373  "[0]+[1]*exp(-[2]*x)",
374  "[0]+[1]*exp(-[2]*x)",
375  "pol1",
376  "[0]+[1]*exp(-[2]*x)"};
377 
378  ptnrbins = pTbinVals_.size() - 1;
379  double *ptbins = new double[pTbinVals_.size()];
380  for (unsigned int b = 0; b < pTbinVals_.size(); b++)
381  ptbins[b] = pTbinVals_[b];
382  double *etabins;
383  if (objectType_ != "met") {
384  etanrbins = etabinVals_.size() - 1;
385  etabins = new double[etabinVals_.size()];
386  for (unsigned int b = 0; b < etabinVals_.size(); b++)
387  etabins[b] = etabinVals_[b];
388  } else {
389  etanrbins = 1;
390  etabins = new double[2];
391  etabins[0] = 0;
392  etabins[1] = 5.;
393  }
394 
395  //define the histograms booked
396  for (Int_t ro = 0; ro < 8; ro++) {
397  for (Int_t etab = 0; etab < etanrbins; etab++) {
398  for (Int_t ptb = 0; ptb < ptnrbins; ptb++) {
399  TString obsName = objectType_;
400  obsName += resObsName[ro];
401  obsName += "_etabin";
402  obsName += etab;
403  obsName += "_ptbin";
404  obsName += ptb;
405  hResPtEtaBin[ro][etab][ptb] = fs->make<TH1F>(obsName, obsName, resObsNrBins, resObsMin[ro], resObsMax[ro]);
406  fResPtEtaBin[ro][etab][ptb] = fs->make<TF1>("F_" + obsName, "gaus");
407  }
408  TString obsName2 = objectType_;
409  obsName2 += resObsName[ro];
410  obsName2 += "_etabin";
411  obsName2 += etab;
412  hResEtaBin[ro][etab] = fs->make<TH1F>(obsName2, obsName2, ptnrbins, ptbins);
413  fResEtaBin[ro][etab] =
414  fs->make<TF1>("F_" + obsName2, resObsVsPtFit[ro], pTbinVals_[0], pTbinVals_[pTbinVals_.size() - 1]);
415  }
416  }
417  tResVar = fs->make<TTree>("tResVar", "Resolution tree");
418 
419  delete[] etabins;
420  delete[] ptbins;
421 }
TF1 * fResEtaBin[10][20]
std::vector< double > etabinVals_
TH1F * hResEtaBin[10][20]
TF1 * fResPtEtaBin[10][20][20]
TH1F * hResPtEtaBin[10][20][20]
double b
Definition: hdecay.h:120
std::vector< double > pTbinVals_

◆ endJob()

void ResolutionCreator::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 424 of file ResolutionCreator.cc.

References relativeConstraints::error, PVValHelper::eta, etabinVals_, etanrbins, fResEtaBin, fResPtEtaBin, hResEtaBin, hResPtEtaBin, createfilelist::int, nrFilled, objectType_, DiDispStaMuonMonitor_cfi::pt, pTbinVals_, ptnrbins, isotrackApplyRegressor::range, tResVar, and relativeConstraints::value.

424  {
425  TString resObsName2[8] = {"a", "b", "c", "d", "theta", "phi", "et", "eta"};
426  Int_t ro = 0;
427  Double_t pt = 0.;
428  Double_t eta = 0.;
429  Double_t value, error;
430 
431  tResVar->Branch("Pt", &pt, "Pt/D");
432  tResVar->Branch("Eta", &eta, "Eta/D");
433  tResVar->Branch("ro", &ro, "ro/I");
434  tResVar->Branch("value", &value, "value/D");
435  tResVar->Branch("error", &error, "error/D");
436 
437  for (ro = 0; ro < 8; ro++) {
438  for (int etab = 0; etab < etanrbins; etab++) {
439  //CD set eta at the center of the bin
440  eta = etanrbins > 1 ? (etabinVals_[etab] + etabinVals_[etab + 1]) / 2. : 2.5;
441  for (int ptb = 0; ptb < ptnrbins; ptb++) {
442  //CD set pt at the center of the bin
443  pt = (pTbinVals_[ptb] + pTbinVals_[ptb + 1]) / 2.;
444  double maxcontent = 0.;
445  int maxbin = 0;
446  for (int nb = 1; nb < hResPtEtaBin[ro][etab][ptb]->GetNbinsX(); nb++) {
447  if (hResPtEtaBin[ro][etab][ptb]->GetBinContent(nb) > maxcontent) {
448  maxcontent = hResPtEtaBin[ro][etab][ptb]->GetBinContent(nb);
449  maxbin = nb;
450  }
451  }
452  int range = (int)(hResPtEtaBin[ro][etab][ptb]->GetNbinsX() / 6); //in order that ~1/3 of X-axis range is fitted
453  fResPtEtaBin[ro][etab][ptb]->SetRange(hResPtEtaBin[ro][etab][ptb]->GetBinCenter(maxbin - range),
454  hResPtEtaBin[ro][etab][ptb]->GetBinCenter(maxbin + range));
455  fResPtEtaBin[ro][etab][ptb]->SetParameters(hResPtEtaBin[ro][etab][ptb]->GetMaximum(),
456  hResPtEtaBin[ro][etab][ptb]->GetMean(),
457  hResPtEtaBin[ro][etab][ptb]->GetRMS());
458  hResPtEtaBin[ro][etab][ptb]->Fit(fResPtEtaBin[ro][etab][ptb]->GetName(), "RQ SERIAL");
459  hResEtaBin[ro][etab]->SetBinContent(ptb + 1, fResPtEtaBin[ro][etab][ptb]->GetParameter(2));
460  hResEtaBin[ro][etab]->SetBinError(ptb + 1, fResPtEtaBin[ro][etab][ptb]->GetParError(2));
461  //CD: Fill the tree
462  value = fResPtEtaBin[ro][etab][ptb]->GetParameter(2); //parameter value
463  error = fResPtEtaBin[ro][etab][ptb]->GetParError(2); //parameter error
464  tResVar->Fill();
465  }
466  //CD: add a fake entry in pt=0 for the NN training
467  // for that, use a linear extrapolation.
468  pt = 0.;
469  value = ((pTbinVals_[0] + pTbinVals_[1]) / 2.) *
470  (fResPtEtaBin[ro][etab][0]->GetParameter(2) - fResPtEtaBin[ro][etab][1]->GetParameter(2)) /
471  ((pTbinVals_[2] - pTbinVals_[0]) / 2.) +
472  fResPtEtaBin[ro][etab][0]->GetParameter(2);
473  error = fResPtEtaBin[ro][etab][0]->GetParError(2) + fResPtEtaBin[ro][etab][1]->GetParError(2);
474  tResVar->Fill();
475  // standard fit
476  hResEtaBin[ro][etab]->Fit(fResEtaBin[ro][etab]->GetName(), "RQ SERIAL");
477  }
478  }
479  if (objectType_ == "lJets" && nrFilled == 0)
480  edm::LogProblem("SummaryError") << "No plots filled for light jets \n";
481  if (objectType_ == "bJets" && nrFilled == 0)
482  edm::LogProblem("SummaryError") << "No plots filled for bjets \n";
483  if (objectType_ == "muon" && nrFilled == 0)
484  edm::LogProblem("SummaryError") << "No plots filled for muons \n";
485  if (objectType_ == "electron" && nrFilled == 0)
486  edm::LogProblem("SummaryError") << "No plots filled for electrons \n";
487  if (objectType_ == "tau" && nrFilled == 0)
488  edm::LogProblem("SummaryError") << "No plots filled for taus \n";
489  if (objectType_ == "met" && nrFilled == 0)
490  edm::LogProblem("SummaryError") << "No plots filled for met \n";
491 
492  edm::LogVerbatim("MainResults") << " \n\n";
493  edm::LogVerbatim("MainResults") << " ----------------------------------------------";
494  edm::LogVerbatim("MainResults") << " ----------------------------------------------";
495  edm::LogVerbatim("MainResults") << " Resolutions on " << objectType_ << " with nrfilled: " << nrFilled;
496  edm::LogVerbatim("MainResults") << " ----------------------------------------------";
497  edm::LogVerbatim("MainResults") << " ----------------------------------------------";
498  if (nrFilled != 0 && objectType_ != "met") {
499  for (ro = 0; ro < 8; ro++) {
500  edm::LogVerbatim("MainResults") << "-------------------- ";
501  edm::LogVerbatim("MainResults") << "\n Resolutions on " << resObsName2[ro] << "\n";
502  edm::LogVerbatim("MainResults") << "-------------------- ";
503  for (int etab = 0; etab < etanrbins; etab++) {
504  if (nrFilled != 0 && ro != 6) {
505  if (etab == 0) {
506  edm::LogVerbatim("MainResults")
507  << "if(fabs(eta)<" << etabinVals_[etab + 1] << ") res = " << fResEtaBin[ro][etab]->GetParameter(0)
508  << "+" << fResEtaBin[ro][etab]->GetParameter(1) << "*exp(-(" << fResEtaBin[ro][etab]->GetParameter(2)
509  << "*pt));";
510  } else {
511  edm::LogVerbatim("MainResults")
512  << "else if(fabs(eta)<" << etabinVals_[etab + 1] << ") res = " << fResEtaBin[ro][etab]->GetParameter(0)
513  << "+" << fResEtaBin[ro][etab]->GetParameter(1) << "*exp(-(" << fResEtaBin[ro][etab]->GetParameter(2)
514  << "*pt));";
515  }
516  } else if (nrFilled != 0 && ro == 6) {
517  if (etab == 0) {
518  edm::LogVerbatim("MainResults")
519  << "if(fabs(eta)<" << etabinVals_[etab + 1] << ") res = " << fResEtaBin[ro][etab]->GetParameter(0)
520  << "+" << fResEtaBin[ro][etab]->GetParameter(1) << "*pt;";
521  } else {
522  edm::LogVerbatim("MainResults")
523  << "else if(fabs(eta)<" << etabinVals_[etab + 1] << ") res = " << fResEtaBin[ro][etab]->GetParameter(0)
524  << "+" << fResEtaBin[ro][etab]->GetParameter(1) << "*pt;";
525  }
526  }
527  }
528  }
529  } else if (nrFilled != 0 && objectType_ == "met") {
530  for (ro = 0; ro < 8; ro++) {
531  edm::LogVerbatim("MainResults") << "-------------------- ";
532  edm::LogVerbatim("MainResults") << "\n Resolutions on " << resObsName2[ro] << "\n";
533  edm::LogVerbatim("MainResults") << "-------------------- ";
534  if (nrFilled != 0 && ro != 6) {
535  edm::LogVerbatim("MainResults") << "res = " << fResEtaBin[ro][0]->GetParameter(0) << "+"
536  << fResEtaBin[ro][0]->GetParameter(1) << "*exp(-("
537  << fResEtaBin[ro][0]->GetParameter(2) << "*pt));";
538  } else if (nrFilled != 0 && ro == 6) {
539  edm::LogVerbatim("MainResults") << "res = " << fResEtaBin[ro][0]->GetParameter(0) << "+"
540  << fResEtaBin[ro][0]->GetParameter(1) << "*pt;";
541  }
542  }
543  }
544 }
Log< level::Info, true > LogVerbatim
TF1 * fResEtaBin[10][20]
std::vector< double > etabinVals_
Definition: value.py:1
TH1F * hResEtaBin[10][20]
TF1 * fResPtEtaBin[10][20][20]
TH1F * hResPtEtaBin[10][20][20]
std::vector< double > pTbinVals_
Log< level::Error, true > LogProblem

Member Data Documentation

◆ electronsToken_

edm::EDGetTokenT<std::vector<pat::Electron> > ResolutionCreator::electronsToken_
private

Definition at line 50 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ etabinVals_

std::vector<double> ResolutionCreator::etabinVals_
private

Definition at line 55 of file ResolutionCreator.cc.

Referenced by analyze(), beginJob(), endJob(), and ResolutionCreator().

◆ etanrbins

int ResolutionCreator::etanrbins
private

Definition at line 61 of file ResolutionCreator.cc.

Referenced by analyze(), beginJob(), and endJob().

◆ fResEtaBin

TF1* ResolutionCreator::fResEtaBin[10][20]
private

Definition at line 66 of file ResolutionCreator.cc.

Referenced by beginJob(), and endJob().

◆ fResPtEtaBin

TF1* ResolutionCreator::fResPtEtaBin[10][20][20]
private

Definition at line 65 of file ResolutionCreator.cc.

Referenced by beginJob(), and endJob().

◆ genEvtToken_

edm::EDGetTokenT<TtGenEvent> ResolutionCreator::genEvtToken_
private

Definition at line 48 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ hResEtaBin

TH1F* ResolutionCreator::hResEtaBin[10][20]
private

Definition at line 68 of file ResolutionCreator.cc.

Referenced by beginJob(), and endJob().

◆ hResPtEtaBin

TH1F* ResolutionCreator::hResPtEtaBin[10][20][20]
private

Definition at line 67 of file ResolutionCreator.cc.

Referenced by analyze(), beginJob(), and endJob().

◆ jetsToken_

edm::EDGetTokenT<std::vector<pat::Jet> > ResolutionCreator::jetsToken_
private

Definition at line 52 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ labelName_

std::string ResolutionCreator::labelName_
private

Definition at line 49 of file ResolutionCreator.cc.

Referenced by ResolutionCreator().

◆ matchingAlgo_

int ResolutionCreator::matchingAlgo_
private

Definition at line 57 of file ResolutionCreator.cc.

◆ maxDist_

double ResolutionCreator::maxDist_
private

Definition at line 60 of file ResolutionCreator.cc.

◆ metsToken_

edm::EDGetTokenT<std::vector<pat::MET> > ResolutionCreator::metsToken_
private

Definition at line 53 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ minDR_

double ResolutionCreator::minDR_
private

Definition at line 56 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ muonsToken_

edm::EDGetTokenT<std::vector<pat::Muon> > ResolutionCreator::muonsToken_
private

Definition at line 51 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ nrFilled

int ResolutionCreator::nrFilled
private

Definition at line 62 of file ResolutionCreator.cc.

Referenced by analyze(), endJob(), and ResolutionCreator().

◆ objectType_

std::string ResolutionCreator::objectType_
private

Definition at line 49 of file ResolutionCreator.cc.

Referenced by analyze(), beginJob(), endJob(), and ResolutionCreator().

◆ pTbinVals_

std::vector<double> ResolutionCreator::pTbinVals_
private

Definition at line 55 of file ResolutionCreator.cc.

Referenced by analyze(), beginJob(), endJob(), and ResolutionCreator().

◆ ptnrbins

int ResolutionCreator::ptnrbins
private

Definition at line 61 of file ResolutionCreator.cc.

Referenced by beginJob(), and endJob().

◆ tausToken_

edm::EDGetTokenT<std::vector<pat::Tau> > ResolutionCreator::tausToken_
private

Definition at line 54 of file ResolutionCreator.cc.

Referenced by analyze(), and ResolutionCreator().

◆ tResVar

TTree* ResolutionCreator::tResVar
private

Definition at line 69 of file ResolutionCreator.cc.

Referenced by beginJob(), and endJob().

◆ useDeltaR_

bool ResolutionCreator::useDeltaR_
private

Definition at line 59 of file ResolutionCreator.cc.

◆ useMaxDist_

bool ResolutionCreator::useMaxDist_
private

Definition at line 58 of file ResolutionCreator.cc.