CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
L1ECALPrefiringWeightProducer Class Reference

#include <ProdTutorial/L1ECALPrefiringWeightProducer/plugins/L1ECALPrefiringWeightProducer.cc>

Inheritance diagram for L1ECALPrefiringWeightProducer:
edm::global::EDProducer<> edm::global::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 L1ECALPrefiringWeightProducer (const edm::ParameterSet &)
 
 ~L1ECALPrefiringWeightProducer () override
 
- Public Member Functions inherited from edm::global::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) 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
 
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::vector< ModuleDescription const * > &modules, 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
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::global::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

double getPrefiringRate (double eta, double pt, TH2F *h_prefmap, fluctuations fluctuation) const
 
void produce (edm::StreamID, edm::Event &, const edm::EventSetup &) const override
 

Private Attributes

std::string dataera_
 
TH2F * h_prefmap_jet
 
TH2F * h_prefmap_photon
 
edm::EDGetTokenT< std::vector< pat::Jet > > jets_token_
 
edm::EDGetTokenT< std::vector< pat::Photon > > photons_token_
 
double prefiringRateSystUnc_
 
bool skipwarnings_
 
bool useEMpt_
 

Additional Inherited Members

- Public Types inherited from edm::global::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
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 ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
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)
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 39 of file L1ECALPrefiringWeightProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 63 of file L1ECALPrefiringWeightProducer.cc.

References gather_cfg::cout, dataera_, alignmentValidation::fname, edm::ParameterSet::getParameter(), h_prefmap_jet, h_prefmap_photon, jets_token_, photons_token_, prefiringRateSystUnc_, skipwarnings_, AlCaHLTBitMon_QueryRunRegistry::string, and useEMpt_.

64 {
65  photons_token_ = consumes<std::vector<pat::Photon> >( iConfig.getParameter<edm::InputTag>( "ThePhotons" ) );
66  jets_token_ = consumes<std::vector<pat::Jet> >( iConfig.getParameter<edm::InputTag>( "TheJets" ) );
67 
68  dataera_ = iConfig.getParameter<std::string>( "DataEra" ) ;
69  useEMpt_ = iConfig.getParameter<bool>( "UseJetEMPt" ) ;
70  prefiringRateSystUnc_ = iConfig.getParameter<double>( "PrefiringRateSystematicUncty" ) ;
71  skipwarnings_ = iConfig.getParameter<bool>( "SkipWarnings" ) ;
72 
73  TFile * file_prefiringmaps_;
74  std::string fname = iConfig.getParameter<std::string>( "L1Maps" );
75  edm::FileInPath mapsfilepath("PhysicsTools/PatUtils/data/"+fname);
76  file_prefiringmaps_ = new TFile(mapsfilepath.fullPath().c_str(),"read" );
77  if(file_prefiringmaps_ == nullptr && !skipwarnings_)
78  std::cout << "File with maps not found. All prefiring weights set to 0. " <<std::endl;
79 
80  TString mapphotonfullname= "L1prefiring_photonptvseta_"+ dataera_;
81  if(!file_prefiringmaps_->Get(mapphotonfullname)&& !skipwarnings_)
82  std::cout << "Photon map not found. All photons prefiring weights set to 0. " <<std::endl;
83  h_prefmap_photon =(TH2F*) file_prefiringmaps_->Get(mapphotonfullname);
84 
85  TString mapjetfullname= (useEMpt_) ? "L1prefiring_jetemptvseta_"+ dataera_ : "L1prefiring_jetptvseta_"+ dataera_ ;
86  if(!file_prefiringmaps_->Get(mapjetfullname)&& !skipwarnings_)
87  std::cout << "Jet map not found. All jets prefiring weights set to 0. "<< std::endl ;
88  h_prefmap_jet =(TH2F*) file_prefiringmaps_->Get(mapjetfullname);
89  file_prefiringmaps_->Close();
90  delete file_prefiringmaps_;
91  produces<double>( "nonPrefiringProb" ).setBranchAlias( "nonPrefiringProb");
92  produces<double>( "nonPrefiringProbUp" ).setBranchAlias( "nonPrefiringProbUp");
93  produces<double>( "nonPrefiringProbDown" ).setBranchAlias( "nonPrefiringProbDown");
94 
95 }
T getParameter(std::string const &) const
edm::EDGetTokenT< std::vector< pat::Photon > > photons_token_
string fname
main script
edm::EDGetTokenT< std::vector< pat::Jet > > jets_token_
L1ECALPrefiringWeightProducer::~L1ECALPrefiringWeightProducer ( )
override

Member Function Documentation

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

Definition at line 209 of file L1ECALPrefiringWeightProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), DEFINE_FWK_MODULE, and AlCaHLTBitMon_QueryRunRegistry::string.

209  {
210 
212 
213  desc.add<edm::InputTag>("ThePhotons", edm::InputTag("slimmedPhotons"));
214  desc.add<edm::InputTag>("TheJets", edm::InputTag("slimmedJets"));
215  desc.add<std::string>("L1Maps","L1PrefiringMaps.root");
216  desc.add<std::string>("DataEra", "2017BtoF");
217  desc.add<bool>("UseJetEMPt", false);
218  desc.add<double>("PrefiringRateSystematicUncty",0.2);
219  desc.add<bool>("SkipWarnings", true);
220  descriptions.add("l1ECALPrefiringWeightProducer", desc);
221 
222 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double L1ECALPrefiringWeightProducer::getPrefiringRate ( double  eta,
double  pt,
TH2F *  h_prefmap,
fluctuations  fluctuation 
) const
private

Definition at line 186 of file L1ECALPrefiringWeightProducer.cc.

References gather_cfg::cout, down, SiStripPI::max, min(), funct::pow(), prefiringRateSystUnc_, skipwarnings_, mathSSE::sqrt(), and up.

Referenced by produce().

186  {
187 
188  if(h_prefmap==nullptr&& !skipwarnings_)std::cout << "Prefiring map not found, setting prefiring rate to 0 " <<std::endl;
189  if(h_prefmap==nullptr ) return 0.;
190  //Check pt is not above map overflow
191  int nbinsy = h_prefmap->GetNbinsY();
192  double maxy= h_prefmap->GetYaxis()->GetBinLowEdge(nbinsy+1);
193  if(pt>=maxy) pt = maxy-0.01;
194  int thebin= h_prefmap->FindBin(eta,pt);
195 
196  double prefrate = h_prefmap->GetBinContent(thebin);
197 
198  double statuncty = h_prefmap->GetBinError(thebin) ;
199  double systuncty = prefiringRateSystUnc_*prefrate;
200 
201  if(fluctuation == up) prefrate = std::min(1.,prefrate+sqrt(pow(statuncty,2)+pow(systuncty,2)));
202  if(fluctuation == down) prefrate = std::max(0.,prefrate-sqrt(pow(statuncty,2)+pow(systuncty,2)));
203  return prefrate;
204 
205 }
T sqrt(T t)
Definition: SSEVec.h:18
T min(T a, T b)
Definition: MathUtil.h:58
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void L1ECALPrefiringWeightProducer::produce ( edm::StreamID  ,
edm::Event iEvent,
const edm::EventSetup iSetup 
) const
overrideprivate

Definition at line 104 of file L1ECALPrefiringWeightProducer.cc.

References central, reco::deltaR(), down, PATTauDiscriminationAgainstElectronDeadECAL_cfi::dR, edm::Event::getByToken(), getPrefiringRate(), h_prefmap_jet, h_prefmap_photon, metsig::jet, jets_token_, eostools::move(), muons2muons_cfi::photon, photons_token_, edm::Event::put(), up, and useEMpt_.

105 {
106  using namespace edm;
107 
108  //Photons
110  iEvent.getByToken(photons_token_,thePhotons);
111 
112  //Jets
114  iEvent.getByToken(jets_token_,theJets );
115 
116  //Probability for the event NOT to prefire, computed with the prefiring maps per object.
117  //Up and down values correspond to the resulting value when shifting up/down all prefiring rates in prefiring maps.
118  double nonPrefiringProba[3]={1.,1.,1.};//0: central, 1: up, 2: down
119 
120  for ( const auto fluct : { fluctuations::central, fluctuations::up, fluctuations::down } ) {
121  for ( const auto & photon : *thePhotons ) {
122  double pt_gam= photon.pt();
123  double eta_gam= photon.eta();
124  if( pt_gam < 20.) continue;
125  if( fabs(eta_gam) <2.) continue;
126  if( fabs(eta_gam) >3.) continue;
127  double prefiringprob_gam= getPrefiringRate( eta_gam, pt_gam, h_prefmap_photon, fluct);
128  nonPrefiringProba[fluct] *= (1.-prefiringprob_gam);
129  }
130 
131  //Now applying the prefiring maps to jets in the affected regions.
132  for ( const auto &jet: *theJets ) {
133 
134  double pt_jet= jet.pt();
135  double eta_jet= jet.eta();
136  double phi_jet= jet.phi();
137  if( pt_jet < 20.) continue;
138  if( fabs(eta_jet) <2.) continue;
139  if( fabs(eta_jet) >3.) continue;
140 
141  //Loop over photons to remove overlap
142  double nonprefiringprobfromoverlappingphotons =1.;
143  for ( const auto & photon : *thePhotons ) {
144  double pt_gam= photon.pt();
145  double eta_gam= photon.eta();
146  double phi_gam= photon.phi();
147  if( pt_gam < 20.) continue;
148  if( fabs(eta_gam) <2.) continue;
149  if( fabs(eta_gam) >3.) continue;
150  double dR = reco::deltaR( eta_jet,phi_jet,eta_gam,phi_gam );
151  if(dR>0.4)continue;
152  double prefiringprob_gam = getPrefiringRate( eta_gam, pt_gam, h_prefmap_photon , fluct);
153  nonprefiringprobfromoverlappingphotons *= (1.-prefiringprob_gam) ;
154  }
155  //useEMpt =true if one wants to use maps parametrized vs Jet EM pt instead of pt.
156  if (useEMpt_)
157  pt_jet *= (jet.neutralEmEnergyFraction()+jet.chargedEmEnergyFraction());
158  double nonprefiringprobfromoverlappingjet = 1. - getPrefiringRate( eta_jet, pt_jet , h_prefmap_jet , fluct);
159 
160  if(nonprefiringprobfromoverlappingphotons ==1.){
161  nonPrefiringProba[fluct] *= nonprefiringprobfromoverlappingjet ;
162  }
163  //If overlapping photons have a non prefiring rate larger than the jet, then replace these weights by the jet one
164  else if(nonprefiringprobfromoverlappingphotons > nonprefiringprobfromoverlappingjet ) {
165  if(nonprefiringprobfromoverlappingphotons !=0.) {
166  nonPrefiringProba[fluct] *= nonprefiringprobfromoverlappingjet / nonprefiringprobfromoverlappingphotons;
167  }
168  else{
169  nonPrefiringProba[fluct] = 0.;
170  }
171  }
172  //Last case: if overlapping photons have a non prefiring rate smaller than the jet, don't consider the jet in the event weight, and do nothing.
173  }
174  }
175 
176  auto nonPrefiringProb =std::make_unique <double> ( nonPrefiringProba[0]);
177  auto nonPrefiringProbUp =std::make_unique <double> ( nonPrefiringProba[1]);
178  auto nonPrefiringProbDown =std::make_unique <double> ( nonPrefiringProba[2]);
179  iEvent.put( std::move(nonPrefiringProb), "nonPrefiringProb" );
180  iEvent.put( std::move(nonPrefiringProbUp), "nonPrefiringProbUp" );
181  iEvent.put( std::move(nonPrefiringProbDown), "nonPrefiringProbDown" );
182 
183 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
double getPrefiringRate(double eta, double pt, TH2F *h_prefmap, fluctuations fluctuation) const
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
edm::EDGetTokenT< std::vector< pat::Photon > > photons_token_
HLT enums.
edm::EDGetTokenT< std::vector< pat::Jet > > jets_token_
def move(src, dest)
Definition: eostools.py:510

Member Data Documentation

std::string L1ECALPrefiringWeightProducer::dataera_
private

Definition at line 56 of file L1ECALPrefiringWeightProducer.cc.

Referenced by L1ECALPrefiringWeightProducer().

TH2F* L1ECALPrefiringWeightProducer::h_prefmap_jet
private
TH2F* L1ECALPrefiringWeightProducer::h_prefmap_photon
private
edm::EDGetTokenT<std::vector< pat::Jet> > L1ECALPrefiringWeightProducer::jets_token_
private

Definition at line 52 of file L1ECALPrefiringWeightProducer.cc.

Referenced by L1ECALPrefiringWeightProducer(), and produce().

edm::EDGetTokenT<std::vector< pat::Photon> > L1ECALPrefiringWeightProducer::photons_token_
private

Definition at line 51 of file L1ECALPrefiringWeightProducer.cc.

Referenced by L1ECALPrefiringWeightProducer(), and produce().

double L1ECALPrefiringWeightProducer::prefiringRateSystUnc_
private
bool L1ECALPrefiringWeightProducer::skipwarnings_
private
bool L1ECALPrefiringWeightProducer::useEMpt_
private

Definition at line 57 of file L1ECALPrefiringWeightProducer.cc.

Referenced by L1ECALPrefiringWeightProducer(), and produce().