CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes
pat::PATMuonProducer Class Reference

class definition More...

Inheritance diagram for pat::PATMuonProducer:
edm::stream::EDProducer< edm::GlobalCache< PATMuonHeavyObjectCache > >

Public Member Functions

 PATMuonProducer (const edm::ParameterSet &iConfig, PATMuonHeavyObjectCache const *)
 default constructir More...
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 everything that needs to be done during the event loop More...
 
 ~PATMuonProducer () override
 default destructur More...
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< PATMuonHeavyObjectCache > >
 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)
 description of config file parameters More...
 
static void globalEndJob (PATMuonHeavyObjectCache *)
 
static std::unique_ptr< PATMuonHeavyObjectCacheinitializeGlobalCache (const edm::ParameterSet &iConfig)
 

Private Types

typedef std::vector< edm::Handle< edm::Association< reco::GenParticleCollection > > > GenAssociations
 
typedef std::vector< edm::Handle< edm::ValueMap< IsoDeposit > > > IsoDepositMaps
 
typedef std::pair< pat::IsolationKeys, edm::InputTagIsolationLabel
 
typedef std::vector< IsolationLabelIsolationLabels
 
typedef std::vector< edm::Handle< edm::ValueMap< double > > > IsolationValueMaps
 
typedef edm::RefToBase< reco::MuonMuonBaseRef
 typedefs for convenience More...
 

Private Member Functions

void embedHighLevel (pat::Muon &aMuon, reco::TrackRef track, reco::TransientTrack &tt, reco::Vertex &primaryVertex, bool primaryVertexIsValid, reco::BeamSpot &beamspot, bool beamspotIsValid)
 
void fillHltTriggerInfo (pat::Muon &muon, edm::Handle< std::vector< pat::TriggerObjectStandAlone >> &triggerObjects, const edm::TriggerNames &names, const std::vector< std::string > &collection_names)
 
void fillL1TriggerInfo (pat::Muon &muon, edm::Handle< std::vector< pat::TriggerObjectStandAlone >> &triggerObjects, const edm::TriggerNames &names, const edm::ESHandle< GlobalTrackingGeometry > &geometry)
 
void fillMuon (Muon &aMuon, const MuonBaseRef &muonRef, const reco::CandidateBaseRef &baseRef, const GenAssociations &genMatches, const IsoDepositMaps &deposits, const IsolationValueMaps &isolationValues) const
 common muon filling, for both the standard and PF2PAT case More...
 
std::optional< GlobalPointgetMuonDirection (const reco::MuonChamberMatch &chamberMatch, const edm::ESHandle< GlobalTrackingGeometry > &geometry, const DetId &chamberId)
 
double getRelMiniIsoPUCorrected (const pat::Muon &muon, double rho, const std::vector< double > &area)
 
bool isChargedHadron (long pdgid)
 
bool isNeutralHadron (long pdgid)
 
bool isPhoton (long pdgid)
 
double puppiCombinedIsolation (const pat::Muon &muon, const pat::PackedCandidateCollection *pc)
 
template<typename T >
void readIsolationLabels (const edm::ParameterSet &iConfig, const char *psetName, IsolationLabels &labels, std::vector< edm::EDGetTokenT< edm::ValueMap< T >>> &tokens)
 
double relMiniIsoPUCorrected (const pat::Muon &aMuon, double rho)
 
void setMuonMiniIso (pat::Muon &aMuon, const pat::PackedCandidateCollection *pc)
 

Private Attributes

bool addEfficiencies_
 add efficiencies to the muon (this will be data members of th muon even w/o embedding) More...
 
bool addGenMatch_
 add generator match information More...
 
bool addInverseBeta_
 add combined inverse beta measurement into the muon More...
 
bool addPuppiIsolation_
 add puppi isolation More...
 
bool addResolutions_
 add resolutions to the muon (this will be data members of th muon even w/o embedding) More...
 
bool addTriggerMatching_
 Trigger. More...
 
edm::EDGetTokenT< reco::BeamSpotbeamLineToken_
 input source of the primary vertex/beamspot More...
 
edm::EDGetTokenT< edm::ValueMap< reco::MuonMETCorrectionData > > caloMETMuonCorrsToken_
 source of caloMET muon corrections More...
 
bool computeMiniIso_
 
bool computeMuonIDMVA_
 
bool computeMuonMVA_
 standard muon selectors More...
 
bool computePuppiCombinedIso_
 
bool computeSoftMuonMVA_
 
std::vector< double > effectiveAreaVec_
 
pat::helper::EfficiencyLoader efficiencyLoader_
 helper class to add efficiencies to the muon More...
 
bool embedBestTrack_
 embed the track from best muon measurement (global pflow) More...
 
bool embedCaloMETMuonCorrs_
 embed muon MET correction info for caloMET into the muon More...
 
bool embedCombinedMuon_
 embed track of the combined fit into the muon More...
 
bool embedDytMuon_
 embed track from DYT muon fit into the muon More...
 
bool embedGenMatch_
 embed the gen match information into the muon More...
 
bool embedHighLevelSelection_
 embed high level selection variables More...
 
bool embedPFCandidate_
 embed pfCandidates into the muon More...
 
bool embedPfEcalEnergy_
 add ecal PF energy More...
 
bool embedPickyMuon_
 embed track from picky muon fit into the muon More...
 
bool embedStandAloneMuon_
 embed track from muon system into the muon More...
 
bool embedTcMETMuonCorrs_
 embed muon MET correction info for tcMET into the muon More...
 
bool embedTpfmsMuon_
 embed track from tpfms muon fit into the muon More...
 
bool embedTrack_
 embed the track from inner tracker into the muon More...
 
bool embedTunePBestTrack_
 embed the track from best muon measurement (muon only) More...
 
bool forceEmbedBestTrack_
 force separate embed of the best track even if already embedded More...
 
std::vector< edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > > genMatchTokens_
 input tags for generator match information More...
 
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecordgeometryToken_
 
std::vector< std::string > hltCollectionFilters_
 
IsolationLabels isoDepositLabels_
 input source for isoDeposits More...
 
std::vector< edm::EDGetTokenT< edm::ValueMap< IsoDeposit > > > isoDepositTokens_
 
IsolationLabels isolationValueLabels_
 input source isolation value maps More...
 
std::vector< edm::EDGetTokenT< edm::ValueMap< double > > > isolationValueTokens_
 
pat::helper::MultiIsolator isolator_
 
pat::helper::MultiIsolator::IsolationValuePairs isolatorTmpStorage_
 isolation value pair for temporary storage before being folded into the muon More...
 
std::vector< double > miniIsoParams_
 
edm::EDGetTokenT< edm::ValueMap< reco::MuonTimeExtra > > muonTimeExtraToken_
 input tag for reading inverse beta More...
 
edm::EDGetTokenT< edm::View< reco::Muon > > muonToken_
 input source More...
 
edm::EDGetTokenT< reco::JetTagCollectionmvaBTagCollectionTag_
 
edm::EDGetTokenT< reco::JetCorrectormvaL1Corrector_
 
edm::EDGetTokenT< reco::JetCorrectormvaL1L2L3ResCorrector_
 
bool mvaUseJec_
 
const edm::EDPutTokenT< std::vector< Muon > > patMuonPutToken_
 
edm::EDGetTokenT< pat::PackedCandidateCollectionpcToken_
 
edm::EDGetTokenT< reco::PFCandidateCollectionpfMuonToken_
 input source pfCandidates that will be to be transformed into pat::Muons, when using PF2PAT More...
 
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_charged_hadrons_
 
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_neutral_hadrons_
 
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_photons_
 
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_charged_hadrons_
 
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_neutral_hadrons_
 
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_photons_
 
edm::EDGetTokenT< std::vector< reco::Vertex > > pvToken_
 input source of the primary vertex More...
 
bool recomputeBasicSelectors_
 
double relMiniIsoPUCorrected_
 
pat::helper::KinResolutionsLoader resolutionLoader_
 helper class to add resolutions to the muon More...
 
edm::EDGetTokenT< double > rho_
 
edm::EDGetTokenT< edm::ValueMap< reco::MuonSimInfo > > simInfo_
 MC info. More...
 
edm::EDGetTokenT< edm::ValueMap< reco::MuonMETCorrectionData > > tcMETMuonCorrsToken_
 source of tcMET muon corrections More...
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordtransientTrackBuilderToken_
 
edm::EDGetTokenT< std::vector< pat::TriggerObjectStandAlone > > triggerObjects_
 
edm::EDGetTokenT< edm::TriggerResultstriggerResults_
 
bool useParticleFlow_
 switch to use particle flow (PF2PAT) or not More...
 
pat::PATUserDataHelper< pat::MuonuserDataHelper_
 helper class to add userData to the muon More...
 
bool useUserData_
 add user data to the muon (this will be data members of th muon even w/o embedding) More...
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< PATMuonHeavyObjectCache > >
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

class definition

Produces pat::Muon's.

The PATMuonProducer produces analysis-level pat::Muon's starting from a collection of objects of reco::Muon.

Author
Steven Lowette, Roger Wolf
Version
Id
PATMuonProducer.h,v 1.29 2012/08/22 15:02:52 bellan Exp

Definition at line 78 of file PATMuonProducer.cc.

Member Typedef Documentation

◆ GenAssociations

Definition at line 99 of file PATMuonProducer.cc.

◆ IsoDepositMaps

Definition at line 100 of file PATMuonProducer.cc.

◆ IsolationLabel

Definition at line 102 of file PATMuonProducer.cc.

◆ IsolationLabels

Definition at line 103 of file PATMuonProducer.cc.

◆ IsolationValueMaps

typedef std::vector<edm::Handle<edm::ValueMap<double> > > pat::PATMuonProducer::IsolationValueMaps
private

Definition at line 101 of file PATMuonProducer.cc.

◆ MuonBaseRef

typedefs for convenience

Definition at line 98 of file PATMuonProducer.cc.

Constructor & Destructor Documentation

◆ PATMuonProducer()

PATMuonProducer::PATMuonProducer ( const edm::ParameterSet iConfig,
PATMuonHeavyObjectCache const *   
)
explicit

default constructir

Definition at line 347 of file PATMuonProducer.cc.

References deDxTools::esConsumes().

349  useUserData_(iConfig.exists("userData")),
350  computeMuonMVA_(false),
351  computeMuonIDMVA_(false),
352  computeSoftMuonMVA_(false),
354  mvaUseJec_(false),
355  isolator_(iConfig.getParameter<edm::ParameterSet>("userIsolation"), consumesCollector(), false),
357  transientTrackBuilderToken_{esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))},
358  patMuonPutToken_{produces<std::vector<Muon>>()} {
359  // input source
360  muonToken_ = consumes<edm::View<reco::Muon>>(iConfig.getParameter<edm::InputTag>("muonSource"));
361  // embedding of tracks
362  embedBestTrack_ = iConfig.getParameter<bool>("embedMuonBestTrack");
363  embedTunePBestTrack_ = iConfig.getParameter<bool>("embedTunePMuonBestTrack");
364  forceEmbedBestTrack_ = iConfig.getParameter<bool>("forceBestTrackEmbedding");
365  embedTrack_ = iConfig.getParameter<bool>("embedTrack");
366  embedCombinedMuon_ = iConfig.getParameter<bool>("embedCombinedMuon");
367  embedStandAloneMuon_ = iConfig.getParameter<bool>("embedStandAloneMuon");
368  // embedding of muon MET correction information
369  embedCaloMETMuonCorrs_ = iConfig.getParameter<bool>("embedCaloMETMuonCorrs");
370  embedTcMETMuonCorrs_ = iConfig.getParameter<bool>("embedTcMETMuonCorrs");
372  mayConsume<edm::ValueMap<reco::MuonMETCorrectionData>>(iConfig.getParameter<edm::InputTag>("caloMETMuonCorrs"));
374  mayConsume<edm::ValueMap<reco::MuonMETCorrectionData>>(iConfig.getParameter<edm::InputTag>("tcMETMuonCorrs"));
375  // pflow specific configurables
376  useParticleFlow_ = iConfig.getParameter<bool>("useParticleFlow");
377  embedPFCandidate_ = iConfig.getParameter<bool>("embedPFCandidate");
378  pfMuonToken_ = mayConsume<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfMuonSource"));
379  embedPfEcalEnergy_ = iConfig.getParameter<bool>("embedPfEcalEnergy");
380  // embedding of tracks from TeV refit
381  embedPickyMuon_ = iConfig.getParameter<bool>("embedPickyMuon");
382  embedTpfmsMuon_ = iConfig.getParameter<bool>("embedTpfmsMuon");
383  embedDytMuon_ = iConfig.getParameter<bool>("embedDytMuon");
384  // embedding of inverse beta variable information
385  addInverseBeta_ = iConfig.getParameter<bool>("addInverseBeta");
386  if (addInverseBeta_) {
388  consumes<edm::ValueMap<reco::MuonTimeExtra>>(iConfig.getParameter<edm::InputTag>("sourceMuonTimeExtra"));
389  }
390  // Monte Carlo matching
391  addGenMatch_ = iConfig.getParameter<bool>("addGenMatch");
392  if (addGenMatch_) {
393  embedGenMatch_ = iConfig.getParameter<bool>("embedGenMatch");
394  if (iConfig.existsAs<edm::InputTag>("genParticleMatch")) {
396  iConfig.getParameter<edm::InputTag>("genParticleMatch")));
397  } else {
399  iConfig.getParameter<std::vector<edm::InputTag>>("genParticleMatch"),
400  [this](edm::InputTag const& tag) { return consumes<edm::Association<reco::GenParticleCollection>>(tag); });
401  }
402  }
403  // efficiencies
404  addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
405  if (addEfficiencies_) {
407  pat::helper::EfficiencyLoader(iConfig.getParameter<edm::ParameterSet>("efficiencies"), consumesCollector());
408  }
409  // resolutions
410  addResolutions_ = iConfig.getParameter<bool>("addResolutions");
411  if (addResolutions_) {
413  pat::helper::KinResolutionsLoader(iConfig.getParameter<edm::ParameterSet>("resolutions"), consumesCollector());
414  }
415  // puppi
416  addPuppiIsolation_ = iConfig.getParameter<bool>("addPuppiIsolation");
417  if (addPuppiIsolation_) {
419  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("puppiIsolationChargedHadrons"));
421  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("puppiIsolationNeutralHadrons"));
423  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("puppiIsolationPhotons"));
424  //puppiNoLeptons
426  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("puppiNoLeptonsIsolationChargedHadrons"));
428  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("puppiNoLeptonsIsolationNeutralHadrons"));
430  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("puppiNoLeptonsIsolationPhotons"));
431  }
432  // read isoDeposit labels, for direct embedding
433  readIsolationLabels(iConfig, "isoDeposits", isoDepositLabels_, isoDepositTokens_);
434  // read isolation value labels, for direct embedding
436  // check to see if the user wants to add user data
437  if (useUserData_) {
438  userDataHelper_ = PATUserDataHelper<Muon>(iConfig.getParameter<edm::ParameterSet>("userData"), consumesCollector());
439  }
440  // embed high level selection variables
441  embedHighLevelSelection_ = iConfig.getParameter<bool>("embedHighLevelSelection");
443  beamLineToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamLineSrc"));
444  pvToken_ = consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("pvSrc"));
445  }
446 
447  //for mini-isolation calculation
448  computeMiniIso_ = iConfig.getParameter<bool>("computeMiniIso");
449 
450  computePuppiCombinedIso_ = iConfig.getParameter<bool>("computePuppiCombinedIso");
451 
452  effectiveAreaVec_ = iConfig.getParameter<std::vector<double>>("effectiveAreaVec");
453 
454  miniIsoParams_ = iConfig.getParameter<std::vector<double>>("miniIsoParams");
455  if (computeMiniIso_ && miniIsoParams_.size() != 9) {
456  throw cms::Exception("ParameterError") << "miniIsoParams must have exactly 9 elements.\n";
457  }
459  pcToken_ = consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfCandsForMiniIso"));
460 
461  // standard selectors
462  recomputeBasicSelectors_ = iConfig.getParameter<bool>("recomputeBasicSelectors");
463  computeMuonMVA_ = iConfig.getParameter<bool>("computeMuonMVA");
464  computeMuonIDMVA_ = iConfig.getParameter<bool>("computeMuonIDMVA");
465  if (computeMuonMVA_ and not computeMiniIso_)
466  throw cms::Exception("ConfigurationError") << "MiniIso is needed for Muon MVA calculation.\n";
467 
468  if (computeMuonMVA_) {
469  // pfCombinedInclusiveSecondaryVertexV2BJetTags
470  mvaBTagCollectionTag_ = consumes<reco::JetTagCollection>(iConfig.getParameter<edm::InputTag>("mvaJetTag"));
471  mvaL1Corrector_ = consumes<reco::JetCorrector>(iConfig.getParameter<edm::InputTag>("mvaL1Corrector"));
472  mvaL1L2L3ResCorrector_ = consumes<reco::JetCorrector>(iConfig.getParameter<edm::InputTag>("mvaL1L2L3ResCorrector"));
473  rho_ = consumes<double>(iConfig.getParameter<edm::InputTag>("rho"));
474  mvaUseJec_ = iConfig.getParameter<bool>("mvaUseJec");
475  }
476 
477  computeSoftMuonMVA_ = iConfig.getParameter<bool>("computeSoftMuonMVA");
478 
479  // MC info
480  simInfo_ = consumes<edm::ValueMap<reco::MuonSimInfo>>(iConfig.getParameter<edm::InputTag>("muonSimInfo"));
481 
482  addTriggerMatching_ = iConfig.getParameter<bool>("addTriggerMatching");
483  if (addTriggerMatching_) {
485  consumes<std::vector<pat::TriggerObjectStandAlone>>(iConfig.getParameter<edm::InputTag>("triggerObjects"));
486  triggerResults_ = consumes<edm::TriggerResults>(iConfig.getParameter<edm::InputTag>("triggerResults"));
487  }
488  hltCollectionFilters_ = iConfig.getParameter<std::vector<std::string>>("hltCollectionFilters");
489 }
bool embedTpfmsMuon_
embed track from tpfms muon fit into the muon
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_charged_hadrons_
bool useUserData_
add user data to the muon (this will be data members of th muon even w/o embedding) ...
void readIsolationLabels(const edm::ParameterSet &iConfig, const char *psetName, IsolationLabels &labels, std::vector< edm::EDGetTokenT< edm::ValueMap< T >>> &tokens)
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_photons_
Assists in assimilating all pat::UserData into pat objects.
bool addPuppiIsolation_
add puppi isolation
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::ValueMap< reco::MuonMETCorrectionData > > tcMETMuonCorrsToken_
source of tcMET muon corrections
bool embedTcMETMuonCorrs_
embed muon MET correction info for tcMET into the muon
edm::EDGetTokenT< std::vector< reco::Vertex > > pvToken_
input source of the primary vertex
edm::EDGetTokenT< edm::ValueMap< reco::MuonTimeExtra > > muonTimeExtraToken_
input tag for reading inverse beta
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_photons_
edm::EDGetTokenT< double > rho_
edm::EDGetTokenT< edm::View< reco::Muon > > muonToken_
input source
bool addEfficiencies_
add efficiencies to the muon (this will be data members of th muon even w/o embedding) ...
edm::EDGetTokenT< edm::ValueMap< reco::MuonMETCorrectionData > > caloMETMuonCorrsToken_
source of caloMET muon corrections
bool embedCaloMETMuonCorrs_
embed muon MET correction info for caloMET into the muon
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< edm::ValueMap< reco::MuonSimInfo > > simInfo_
MC info.
bool addTriggerMatching_
Trigger.
bool embedBestTrack_
embed the track from best muon measurement (global pflow)
edm::EDGetTokenT< reco::JetCorrector > mvaL1Corrector_
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackBuilderToken_
bool addResolutions_
add resolutions to the muon (this will be data members of th muon even w/o embedding) ...
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_charged_hadrons_
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
std::vector< double > miniIsoParams_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< double > effectiveAreaVec_
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_neutral_hadrons_
edm::EDGetTokenT< reco::JetTagCollection > mvaBTagCollectionTag_
edm::EDGetTokenT< reco::PFCandidateCollection > pfMuonToken_
input source pfCandidates that will be to be transformed into pat::Muons, when using PF2PAT ...
pat::helper::MultiIsolator isolator_
std::vector< edm::EDGetTokenT< edm::ValueMap< IsoDeposit > > > isoDepositTokens_
edm::EDGetTokenT< reco::BeamSpot > beamLineToken_
input source of the primary vertex/beamspot
bool addInverseBeta_
add combined inverse beta measurement into the muon
bool useParticleFlow_
switch to use particle flow (PF2PAT) or not
pat::PATUserDataHelper< pat::Muon > userDataHelper_
helper class to add userData to the muon
bool embedStandAloneMuon_
embed track from muon system into the muon
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
bool embedTrack_
embed the track from inner tracker into the muon
bool addGenMatch_
add generator match information
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > geometryToken_
bool embedPfEcalEnergy_
add ecal PF energy
edm::EDGetTokenT< pat::PackedCandidateCollection > pcToken_
std::vector< std::string > hltCollectionFilters_
bool embedTunePBestTrack_
embed the track from best muon measurement (muon only)
edm::EDGetTokenT< reco::JetCorrector > mvaL1L2L3ResCorrector_
pat::helper::EfficiencyLoader efficiencyLoader_
helper class to add efficiencies to the muon
bool embedPickyMuon_
embed track from picky muon fit into the muon
bool embedDytMuon_
embed track from DYT muon fit into the muon
const edm::EDPutTokenT< std::vector< Muon > > patMuonPutToken_
edm::EDGetTokenT< std::vector< pat::TriggerObjectStandAlone > > triggerObjects_
bool embedGenMatch_
embed the gen match information into the muon
IsolationLabels isoDepositLabels_
input source for isoDeposits
std::vector< edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > > genMatchTokens_
input tags for generator match information
IsolationLabels isolationValueLabels_
input source isolation value maps
bool forceEmbedBestTrack_
force separate embed of the best track even if already embedded
bool embedPFCandidate_
embed pfCandidates into the muon
std::vector< edm::EDGetTokenT< edm::ValueMap< double > > > isolationValueTokens_
bool computeMuonMVA_
standard muon selectors
bool embedCombinedMuon_
embed track of the combined fit into the muon
pat::helper::KinResolutionsLoader resolutionLoader_
helper class to add resolutions to the muon
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_neutral_hadrons_
bool embedHighLevelSelection_
embed high level selection variables

◆ ~PATMuonProducer()

PATMuonProducer::~PATMuonProducer ( )
override

default destructur

Definition at line 491 of file PATMuonProducer.cc.

491 {}

Member Function Documentation

◆ embedHighLevel()

void PATMuonProducer::embedHighLevel ( pat::Muon aMuon,
reco::TrackRef  track,
reco::TransientTrack tt,
reco::Vertex primaryVertex,
bool  primaryVertexIsValid,
reco::BeamSpot beamspot,
bool  beamspotIsValid 
)
private

Definition at line 1390 of file PATMuonProducer.cc.

References pat::Muon::BS2D, pat::Muon::BS3D, BeamMonitor_cff::primaryVertex, pat::Muon::PV2D, pat::Muon::PV3D, pat::Muon::PVDZ, mps_fire::result, pat::Muon::setDB(), IPTools::signedImpactParameter3D(), and HLT_2022v15_cff::track.

Referenced by produce().

1396  {
1397  // Correct to PV
1398 
1399  // PV2D
1400  aMuon.setDB(track->dxy(primaryVertex.position()),
1401  track->dxyError(primaryVertex.position(), primaryVertex.covariance()),
1402  pat::Muon::PV2D);
1403 
1404  // PV3D
1405  std::pair<bool, Measurement1D> result =
1407  double d0_corr = result.second.value();
1408  double d0_err = primaryVertexIsValid ? result.second.error() : -1.0;
1409  aMuon.setDB(d0_corr, d0_err, pat::Muon::PV3D);
1410 
1411  // Correct to beam spot
1412 
1413  // BS2D
1414  aMuon.setDB(track->dxy(beamspot), track->dxyError(beamspot), pat::Muon::BS2D);
1415 
1416  // make a fake vertex out of beam spot
1417  reco::Vertex vBeamspot(beamspot.position(), beamspot.rotatedCovariance3D());
1418 
1419  // BS3D
1420  result = IPTools::signedImpactParameter3D(tt, GlobalVector(track->px(), track->py(), track->pz()), vBeamspot);
1421  d0_corr = result.second.value();
1422  d0_err = beamspotIsValid ? result.second.error() : -1.0;
1423  aMuon.setDB(d0_corr, d0_err, pat::Muon::BS3D);
1424 
1425  // PVDZ
1426  aMuon.setDB(
1427  track->dz(primaryVertex.position()), std::hypot(track->dzError(), primaryVertex.zError()), pat::Muon::PVDZ);
1428 }
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:81
Definition: TTTypes.h:54
void setDB(double dB, double edB, IPTYPE type=PV2D)
Definition: Muon.h:247
primaryVertex
hltOfflineBeamSpot for HLTMON
Global3DVector GlobalVector
Definition: GlobalVector.h:10

◆ fillDescriptions()

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

description of config file parameters

Definition at line 1245 of file PATMuonProducer.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addNode(), edm::ParameterSetDescription::addOptional(), pat::helper::KinResolutionsLoader::fillDescription(), pat::PATUserDataHelper< ObjectType >::fillDescription(), edm::ParameterSetDescription::ifValue(), HLT_2022v15_cff::InputTag, or, edm::ParameterSetDescription::setAllowAnything(), edm::ParameterSetDescription::setComment(), and edm::ParameterDescriptionNode::setComment().

1245  {
1247  iDesc.setComment("PAT muon producer module");
1248 
1249  // input source
1250  iDesc.add<edm::InputTag>("muonSource", edm::InputTag("no default"))->setComment("input collection");
1251 
1252  // embedding
1253  iDesc.add<bool>("embedMuonBestTrack", true)->setComment("embed muon best track (global pflow)");
1254  iDesc.add<bool>("embedTunePMuonBestTrack", true)->setComment("embed muon best track (muon only)");
1255  iDesc.add<bool>("forceBestTrackEmbedding", true)
1256  ->setComment(
1257  "force embedding separately the best tracks even if they're already embedded e.g. as tracker or global "
1258  "tracks");
1259  iDesc.add<bool>("embedTrack", true)->setComment("embed external track");
1260  iDesc.add<bool>("embedStandAloneMuon", true)->setComment("embed external stand-alone muon");
1261  iDesc.add<bool>("embedCombinedMuon", false)->setComment("embed external combined muon");
1262  iDesc.add<bool>("embedPickyMuon", false)->setComment("embed external picky track");
1263  iDesc.add<bool>("embedTpfmsMuon", false)->setComment("embed external tpfms track");
1264  iDesc.add<bool>("embedDytMuon", false)->setComment("embed external dyt track ");
1265 
1266  // embedding of MET muon corrections
1267  iDesc.add<bool>("embedCaloMETMuonCorrs", true)->setComment("whether to add MET muon correction for caloMET or not");
1268  iDesc.add<edm::InputTag>("caloMETMuonCorrs", edm::InputTag("muonMETValueMapProducer", "muCorrData"))
1269  ->setComment("source of MET muon corrections for caloMET");
1270  iDesc.add<bool>("embedTcMETMuonCorrs", true)->setComment("whether to add MET muon correction for tcMET or not");
1271  iDesc.add<edm::InputTag>("tcMETMuonCorrs", edm::InputTag("muonTCMETValueMapProducer", "muCorrData"))
1272  ->setComment("source of MET muon corrections for tcMET");
1273 
1274  // pf specific parameters
1275  iDesc.add<edm::InputTag>("pfMuonSource", edm::InputTag("pfMuons"))->setComment("particle flow input collection");
1276  iDesc.add<bool>("useParticleFlow", false)->setComment("whether to use particle flow or not");
1277  iDesc.add<bool>("embedPFCandidate", false)->setComment("embed external particle flow object");
1278  iDesc.add<bool>("embedPfEcalEnergy", true)->setComment("add ecal energy as reconstructed by PF");
1279 
1280  // inverse beta computation
1281  iDesc.add<bool>("addInverseBeta", true)->setComment("add combined inverse beta");
1282  iDesc.add<edm::InputTag>("sourceInverseBeta", edm::InputTag("muons", "combined"))
1283  ->setComment("source of inverse beta values");
1284 
1285  // MC matching configurables
1286  iDesc.add<bool>("addGenMatch", true)->setComment("add MC matching");
1287  iDesc.add<bool>("embedGenMatch", false)->setComment("embed MC matched MC information");
1288  std::vector<edm::InputTag> emptySourceVector;
1289  iDesc
1290  .addNode(edm::ParameterDescription<edm::InputTag>("genParticleMatch", edm::InputTag(), true) xor
1291  edm::ParameterDescription<std::vector<edm::InputTag>>("genParticleMatch", emptySourceVector, true))
1292  ->setComment("input with MC match information");
1293 
1294  // mini-iso
1295  iDesc.add<bool>("computeMiniIso", false)->setComment("whether or not to compute and store electron mini-isolation");
1296  iDesc.add<bool>("computePuppiCombinedIso", false)
1297  ->setComment("whether or not to compute and store puppi combined isolation");
1298 
1299  iDesc.add<edm::InputTag>("pfCandsForMiniIso", edm::InputTag("packedPFCandidates"))
1300  ->setComment("collection to use to compute mini-iso");
1301  iDesc.add<std::vector<double>>("miniIsoParams", std::vector<double>())
1302  ->setComment("mini-iso parameters to use for muons");
1303 
1304  iDesc.add<bool>("addTriggerMatching", false)->setComment("add L1 and HLT matching to offline muon");
1305 
1307 
1308  // IsoDeposit configurables
1309  edm::ParameterSetDescription isoDepositsPSet;
1310  isoDepositsPSet.addOptional<edm::InputTag>("tracker");
1311  isoDepositsPSet.addOptional<edm::InputTag>("ecal");
1312  isoDepositsPSet.addOptional<edm::InputTag>("hcal");
1313  isoDepositsPSet.addOptional<edm::InputTag>("particle");
1314  isoDepositsPSet.addOptional<edm::InputTag>("pfChargedHadrons");
1315  isoDepositsPSet.addOptional<edm::InputTag>("pfChargedAll");
1316  isoDepositsPSet.addOptional<edm::InputTag>("pfPUChargedHadrons");
1317  isoDepositsPSet.addOptional<edm::InputTag>("pfNeutralHadrons");
1318  isoDepositsPSet.addOptional<edm::InputTag>("pfPhotons");
1319  isoDepositsPSet.addOptional<std::vector<edm::InputTag>>("user");
1320  iDesc.addOptional("isoDeposits", isoDepositsPSet);
1321 
1322  // isolation values configurables
1323  edm::ParameterSetDescription isolationValuesPSet;
1324  isolationValuesPSet.addOptional<edm::InputTag>("tracker");
1325  isolationValuesPSet.addOptional<edm::InputTag>("ecal");
1326  isolationValuesPSet.addOptional<edm::InputTag>("hcal");
1327  isolationValuesPSet.addOptional<edm::InputTag>("particle");
1328  isolationValuesPSet.addOptional<edm::InputTag>("pfChargedHadrons");
1329  isolationValuesPSet.addOptional<edm::InputTag>("pfChargedAll");
1330  isolationValuesPSet.addOptional<edm::InputTag>("pfPUChargedHadrons");
1331  isolationValuesPSet.addOptional<edm::InputTag>("pfNeutralHadrons");
1332  isolationValuesPSet.addOptional<edm::InputTag>("pfPhotons");
1333  iDesc.addOptional("isolationValues", isolationValuesPSet);
1334 
1335  iDesc.ifValue(edm::ParameterDescription<bool>("addPuppiIsolation", false, true),
1337  "puppiIsolationChargedHadrons",
1338  edm::InputTag("muonPUPPIIsolation", "h+-DR030-ThresholdVeto000-ConeVeto000"),
1339  true) and
1341  "puppiIsolationNeutralHadrons",
1342  edm::InputTag("muonPUPPIIsolation", "h0-DR030-ThresholdVeto000-ConeVeto001"),
1343  true) and
1345  "puppiIsolationPhotons",
1346  edm::InputTag("muonPUPPIIsolation", "gamma-DR030-ThresholdVeto000-ConeVeto001"),
1347  true) and
1349  "puppiNoLeptonsIsolationChargedHadrons",
1350  edm::InputTag("muonPUPPINoLeptonsIsolation", "h+-DR030-ThresholdVeto000-ConeVeto000"),
1351  true) and
1353  "puppiNoLeptonsIsolationNeutralHadrons",
1354  edm::InputTag("muonPUPPINoLeptonsIsolation", "h0-DR030-ThresholdVeto000-ConeVeto001"),
1355  true) and
1357  "puppiNoLeptonsIsolationPhotons",
1358  edm::InputTag("muonPUPPINoLeptonsIsolation", "gamma-DR030-ThresholdVeto000-ConeVeto001"),
1359  true)) or
1360  false >> edm::EmptyGroupDescription());
1361 
1362  // Efficiency configurables
1363  edm::ParameterSetDescription efficienciesPSet;
1364  efficienciesPSet.setAllowAnything(); // TODO: the pat helper needs to implement a description.
1365  iDesc.add("efficiencies", efficienciesPSet);
1366  iDesc.add<bool>("addEfficiencies", false);
1367 
1368  // Check to see if the user wants to add user data
1369  edm::ParameterSetDescription userDataPSet;
1371  iDesc.addOptional("userData", userDataPSet);
1372 
1373  edm::ParameterSetDescription isolationPSet;
1374  isolationPSet.setAllowAnything(); // TODO: the pat helper needs to implement a description.
1375  iDesc.add("userIsolation", isolationPSet);
1376 
1377  iDesc.add<bool>("embedHighLevelSelection", true)->setComment("embed high level selection");
1378  edm::ParameterSetDescription highLevelPSet;
1379  highLevelPSet.setAllowAnything();
1380  iDesc.addNode(edm::ParameterDescription<edm::InputTag>("beamLineSrc", edm::InputTag(), true))
1381  ->setComment("input with high level selection");
1383  ->setComment("input with high level selection");
1384 
1385  //descriptions.add("PATMuonProducer", iDesc);
1386 }
void setComment(std::string const &value)
ParameterDescriptionNode * ifValue(ParameterDescription< T > const &switchParameter, std::unique_ptr< ParameterDescriptionCases< T >> cases)
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
void setAllowAnything()
allow any parameter label/value pairs
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
static void fillDescription(edm::ParameterSetDescription &iDesc)
void setComment(std::string const &value)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescription(edm::ParameterSetDescription &iDesc)
Method for documentation and validation of PSet.

◆ fillHltTriggerInfo()

void PATMuonProducer::fillHltTriggerInfo ( pat::Muon muon,
edm::Handle< std::vector< pat::TriggerObjectStandAlone >> &  triggerObjects,
const edm::TriggerNames names,
const std::vector< std::string > &  collection_names 
)
private

Definition at line 556 of file PATMuonProducer.cc.

References PbPb_ZMuSkimMuonDPG_cff::deltaR, HLT_2022v15_cff::muon, Skims_PA_cff::name, names, getGTfromDQMFile::obj, trigger::TriggerMuon, and triggerMatchMonitor_cfi::triggerObjects.

Referenced by produce().

559  {
560  // WARNING: in a case of close-by muons the dR matching may select both muons.
561  // It's better to select the best match for a given collection.
562  for (const auto& triggerObject : *triggerObjects) {
563  if (triggerObject.hasTriggerObjectType(trigger::TriggerMuon)) {
564  bool keepIt = false;
565  for (const auto& name : collection_filter_names) {
566  if (triggerObject.hasCollection(name)) {
567  keepIt = true;
568  break;
569  }
570  }
571  if (not keepIt)
572  continue;
573  if (deltaR(triggerObject.p4(), muon) > 0.1)
574  continue;
575  pat::TriggerObjectStandAlone obj(triggerObject);
576  obj.unpackPathNames(names);
577  muon.addTriggerObjectMatch(obj);
578  }
579  }
580 }
const std::string names[nVars_]
Analysis-level trigger object class (stand-alone)

◆ fillL1TriggerInfo()

void PATMuonProducer::fillL1TriggerInfo ( pat::Muon muon,
edm::Handle< std::vector< pat::TriggerObjectStandAlone >> &  triggerObjects,
const edm::TriggerNames names,
const edm::ESHandle< GlobalTrackingGeometry > &  geometry 
)
private

Definition at line 504 of file PATMuonProducer.cc.

References funct::abs(), pat::PATObject< ObjectType >::addTriggerObjectMatch(), MuonSubdetId::CSC, SiPixelRawToDigiRegional_cfi::deltaPhi, PbPb_ZMuSkimMuonDPG_cff::deltaR, MuonSubdetId::DT, getMuonDirection(), reco::Muon::matches(), names, getGTfromDQMFile::obj, trigger::TriggerL1Mu, and triggerMatchMonitor_cfi::triggerObjects.

Referenced by produce().

507  {
508  // L1 trigger object parameters are defined at MB2/ME2. Use the muon
509  // chamber matching information to get the local direction of the
510  // muon trajectory and convert it to a global direction to match the
511  // trigger objects
512 
513  std::optional<GlobalPoint> muonPosition;
514  // Loop over chambers
515  // initialize muonPosition with any available match, just in case
516  // the second station is missing - it's better folling back to
517  // dR matching at IP
518  for (const auto& chamberMatch : aMuon.matches()) {
519  if (chamberMatch.id.subdetId() == MuonSubdetId::DT) {
520  DTChamberId detId(chamberMatch.id.rawId());
521  if (abs(detId.station()) > 3)
522  continue;
523  muonPosition = getMuonDirection(chamberMatch, geometry, detId);
524  if (abs(detId.station()) == 2)
525  break;
526  }
527  if (chamberMatch.id.subdetId() == MuonSubdetId::CSC) {
528  CSCDetId detId(chamberMatch.id.rawId());
529  if (abs(detId.station()) > 3)
530  continue;
531  muonPosition = getMuonDirection(chamberMatch, geometry, detId);
532  if (abs(detId.station()) == 2)
533  break;
534  }
535  }
536  if (not muonPosition)
537  return;
538  for (const auto& triggerObject : *triggerObjects) {
539  if (triggerObject.hasTriggerObjectType(trigger::TriggerL1Mu)) {
540  if (std::abs(triggerObject.eta()) < 0.001) {
541  // L1 is defined in X-Y plain
542  if (deltaPhi(triggerObject.phi(), muonPosition->phi()) > 0.1)
543  continue;
544  } else {
545  // 3D L1
546  if (deltaR(triggerObject.p4(), *muonPosition) > 0.15)
547  continue;
548  }
549  pat::TriggerObjectStandAlone obj(triggerObject);
550  obj.unpackPathNames(names);
551  aMuon.addTriggerObjectMatch(obj);
552  }
553  }
554 }
enum start value shifted to 81 so as to avoid clashes with PDG codes
const std::string names[nVars_]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr int DT
Definition: MuonSubdetId.h:11
static constexpr int CSC
Definition: MuonSubdetId.h:12
std::optional< GlobalPoint > getMuonDirection(const reco::MuonChamberMatch &chamberMatch, const edm::ESHandle< GlobalTrackingGeometry > &geometry, const DetId &chamberId)
Analysis-level trigger object class (stand-alone)

◆ fillMuon()

void PATMuonProducer::fillMuon ( Muon aMuon,
const MuonBaseRef muonRef,
const reco::CandidateBaseRef baseRef,
const GenAssociations genMatches,
const IsoDepositMaps deposits,
const IsolationValueMaps isolationValues 
) const
private

common muon filling, for both the standard and PF2PAT case

Definition at line 1081 of file PATMuonProducer.cc.

References addGenMatch_, pat::PATObject< ObjectType >::addGenParticleRef(), edm::contains(), CandIsolatorFromDeposits_cfi::deposits, reco::Muon::DYT, efficiencyLoader_, embedBestTrack_, pat::Muon::embedCombinedMuon(), embedCombinedMuon_, pat::Muon::embedDytMuon(), embedDytMuon_, embedGenMatch_, pat::PATObject< ObjectType >::embedGenParticle(), pat::Muon::embedMuonBestTrack(), pat::Muon::embedPickyMuon(), embedPickyMuon_, pat::Muon::embedStandAloneMuon(), embedStandAloneMuon_, pat::Muon::embedTpfmsMuon(), embedTpfmsMuon_, pat::Muon::embedTrack(), embedTrack_, embedTunePBestTrack_, pat::Muon::embedTunePMuonBestTrack(), pat::helper::EfficiencyLoader::enabled(), pat::helper::KinResolutionsLoader::enabled(), first, forceEmbedBestTrack_, TtFullHadEvtBuilder_cfi::genMatch, edm::RefToBase< T >::id(), reco::Muon::isAValidMuonTrack(), reco::Muon::isGlobalMuon(), isoDepositLabels_, isolationValueLabels_, displacedMuonProducer_cff::isolationValues, dqmiolumiharvest::j, pat::Muon::pfCandidateRef(), reco::Muon::Picky, resolutionLoader_, pat::helper::EfficiencyLoader::setEfficiencies(), pat::Lepton< LeptonType >::setIsoDeposit(), pat::Lepton< LeptonType >::setIsolation(), reco::LeafCandidate::setP4(), pat::helper::KinResolutionsLoader::setResolutions(), source, reco::Muon::TPFMS, and useParticleFlow_.

Referenced by produce().

1086  {
1087  // in the particle flow algorithm,
1088  // the muon momentum is recomputed.
1089  // the new value is stored as the momentum of the
1090  // resulting PFCandidate of type Muon, and choosen
1091  // as the pat::Muon momentum
1092  if (useParticleFlow_)
1093  aMuon.setP4(aMuon.pfCandidateRef()->p4());
1094  if (embedTrack_)
1095  aMuon.embedTrack();
1097  aMuon.embedStandAloneMuon();
1098  if (embedCombinedMuon_)
1099  aMuon.embedCombinedMuon();
1100 
1101  // embed the TeV refit track refs (only available for globalMuons)
1102  if (aMuon.isGlobalMuon()) {
1104  aMuon.embedPickyMuon();
1106  aMuon.embedTpfmsMuon();
1108  aMuon.embedDytMuon();
1109  }
1110 
1111  // embed best tracks (at the end, so unless forceEmbedBestTrack_ is true we can save some space not embedding them twice)
1112  if (embedBestTrack_)
1116 
1117  // store the match to the generated final state muons
1118  if (addGenMatch_) {
1119  for (auto const& genMatch : genMatches) {
1120  reco::GenParticleRef genMuon = (*genMatch)[baseRef];
1121  aMuon.addGenParticleRef(genMuon);
1122  }
1123  if (embedGenMatch_)
1124  aMuon.embedGenParticle();
1125  }
1126  if (efficiencyLoader_.enabled()) {
1127  efficiencyLoader_.setEfficiencies(aMuon, muonRef);
1128  }
1129 
1130  for (size_t j = 0, nd = deposits.size(); j < nd; ++j) {
1131  if (useParticleFlow_) {
1132  if (deposits[j]->contains(baseRef.id())) {
1133  aMuon.setIsoDeposit(isoDepositLabels_[j].first, (*deposits[j])[baseRef]);
1134  } else if (deposits[j]->contains(muonRef.id())) {
1135  aMuon.setIsoDeposit(isoDepositLabels_[j].first, (*deposits[j])[muonRef]);
1136  } else {
1137  reco::CandidatePtr source = aMuon.pfCandidateRef()->sourceCandidatePtr(0);
1139  }
1140  } else {
1141  aMuon.setIsoDeposit(isoDepositLabels_[j].first, (*deposits[j])[muonRef]);
1142  }
1143  }
1144 
1145  for (size_t j = 0; j < isolationValues.size(); ++j) {
1146  if (useParticleFlow_) {
1147  if (isolationValues[j]->contains(baseRef.id())) {
1149  } else if (isolationValues[j]->contains(muonRef.id())) {
1151  } else {
1152  reco::CandidatePtr source = aMuon.pfCandidateRef()->sourceCandidatePtr(0);
1154  }
1155  } else {
1157  }
1158  }
1159 
1160  if (resolutionLoader_.enabled()) {
1162  }
1163 }
bool embedTpfmsMuon_
embed track from tpfms muon fit into the muon
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
void embedDytMuon()
embed reference to the above dyt Track
void embedTpfmsMuon()
embed reference to the above tpfms Track
void setIsolation(IsolationKeys key, float value)
Definition: Lepton.h:115
bool contains(EventRange const &lh, EventID const &rh)
Definition: EventRange.cc:37
void embedCombinedMuon()
set reference to Track reconstructed in both tracked and muon detector (reimplemented from reco::Muon...
void embedTunePMuonBestTrack(bool force=false)
void setEfficiencies(pat::PATObject< T > &obj, const R &originalRef) const
Sets the efficiencies for this object, using the reference to the original objects.
bool isAValidMuonTrack(const MuonTrackType &type) const
void embedMuonBestTrack(bool force=false)
bool embedBestTrack_
embed the track from best muon measurement (global pflow)
void setResolutions(pat::PATObject< T > &obj) const
Sets the efficiencies for this object, using the reference to the original objects.
void embedStandAloneMuon()
set reference to Track reconstructed in the muon detector only (reimplemented from reco::Muon) ...
ProductID id() const
Definition: RefToBase.h:214
void setIsoDeposit(IsolationKeys key, const IsoDeposit &dep)
Sets the IsoDeposit associated with some key; if it is already existent, it is overwritten.
Definition: Lepton.h:191
void embedTrack()
set reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
bool useParticleFlow_
switch to use particle flow (PF2PAT) or not
bool embedStandAloneMuon_
embed track from muon system into the muon
bool embedTrack_
embed the track from inner tracker into the muon
bool addGenMatch_
add generator match information
void embedGenParticle()
Definition: PATObject.h:768
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
void addGenParticleRef(const reco::GenParticleRef &ref)
Definition: PATObject.h:751
bool embedTunePBestTrack_
embed the track from best muon measurement (muon only)
pat::helper::EfficiencyLoader efficiencyLoader_
helper class to add efficiencies to the muon
bool embedPickyMuon_
embed track from picky muon fit into the muon
bool embedDytMuon_
embed track from DYT muon fit into the muon
reco::PFCandidateRef pfCandidateRef() const
bool embedGenMatch_
embed the gen match information into the muon
IsolationLabels isoDepositLabels_
input source for isoDeposits
IsolationLabels isolationValueLabels_
input source isolation value maps
void embedPickyMuon()
embed reference to the above picky Track
bool forceEmbedBestTrack_
force separate embed of the best track even if already embedded
bool embedCombinedMuon_
embed track of the combined fit into the muon
pat::helper::KinResolutionsLoader resolutionLoader_
helper class to add resolutions to the muon
static std::string const source
Definition: EdmProvDump.cc:49
bool isGlobalMuon() const override
Definition: Muon.h:303
void setP4(const LorentzVector &p4) final
set 4-momentum
genMatch
add extra information on genMatch

◆ getMuonDirection()

std::optional< GlobalPoint > PATMuonProducer::getMuonDirection ( const reco::MuonChamberMatch chamberMatch,
const edm::ESHandle< GlobalTrackingGeometry > &  geometry,
const DetId chamberId 
)
private

Definition at line 493 of file PATMuonProducer.cc.

References GeomDet::toGlobal(), reco::MuonChamberMatch::x, and reco::MuonChamberMatch::y.

Referenced by fillL1TriggerInfo().

495  {
496  const GeomDet* chamberGeometry = geometry->idToDet(chamberId);
497  if (chamberGeometry) {
498  LocalPoint localPosition(chamberMatch.x, chamberMatch.y, 0);
499  return std::optional<GlobalPoint>(std::in_place, chamberGeometry->toGlobal(localPosition));
500  }
501  return std::optional<GlobalPoint>();
502 }
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49

◆ getRelMiniIsoPUCorrected()

double PATMuonProducer::getRelMiniIsoPUCorrected ( const pat::Muon muon,
double  rho,
const std::vector< double > &  area 
)
private

Definition at line 1180 of file PATMuonProducer.cc.

References custom_jme_cff::area, MTVHistoProducerAlgoForTrackerBlock_cfi::maxdr, MTVHistoProducerAlgoForTrackerBlock_cfi::mindr, pat::miniIsoDr(), miniIsoParams_, and pat::muonRelMiniIsoPUCorrected().

Referenced by produce().

1180  {
1181  double mindr(miniIsoParams_[0]);
1182  double maxdr(miniIsoParams_[1]);
1183  double kt_scale(miniIsoParams_[2]);
1184  double drcut = pat::miniIsoDr(muon.polarP4(), mindr, maxdr, kt_scale);
1185  return pat::muonRelMiniIsoPUCorrected(muon.miniPFIsolation(), muon.polarP4(), drcut, rho, area);
1186 }
std::vector< double > miniIsoParams_
double muonRelMiniIsoPUCorrected(const PFIsolation &iso, const reco::Candidate::PolarLorentzVector &p4, double dr, double rho, const std::vector< double > &area)
float miniIsoDr(const reco::Candidate::PolarLorentzVector &p4, float mindr, float maxdr, float kt_scale)

◆ globalEndJob()

static void pat::PATMuonProducer::globalEndJob ( PATMuonHeavyObjectCache )
inlinestatic

Definition at line 89 of file PATMuonProducer.cc.

89 {}

◆ initializeGlobalCache()

static std::unique_ptr<PATMuonHeavyObjectCache> pat::PATMuonProducer::initializeGlobalCache ( const edm::ParameterSet iConfig)
inlinestatic

Definition at line 85 of file PATMuonProducer.cc.

85  {
86  return std::make_unique<PATMuonHeavyObjectCache>(iConfig);
87  }

◆ isChargedHadron()

bool PATMuonProducer::isChargedHadron ( long  pdgid)
private

Definition at line 1240 of file PATMuonProducer.cc.

References funct::abs(), and EgammaValidation_cff::pdgid.

Referenced by puppiCombinedIsolation().

1240 { return std::abs(pdgid) == 211; }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ isNeutralHadron()

bool PATMuonProducer::isNeutralHadron ( long  pdgid)
private

Definition at line 1238 of file PATMuonProducer.cc.

References funct::abs(), and EgammaValidation_cff::pdgid.

Referenced by puppiCombinedIsolation().

1238 { return std::abs(pdgid) == 130; }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ isPhoton()

bool PATMuonProducer::isPhoton ( long  pdgid)
private

Definition at line 1242 of file PATMuonProducer.cc.

References EgammaValidation_cff::pdgid.

Referenced by puppiCombinedIsolation().

1242 { return pdgid == 22; }

◆ produce()

void PATMuonProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

everything that needs to be done during the event loop

Definition at line 582 of file PATMuonProducer.cc.

References funct::abs(), pat::PATUserDataHelper< ObjectType >::add(), addGenMatch_, addInverseBeta_, addPuppiIsolation_, addTriggerMatching_, cms::cuda::assert(), beamLineToken_, pwdgSkimBPark_cfi::beamSpot, pat::helper::MultiIsolator::beginEvent(), TransientTrackBuilder::build(), muonProducer_cfi::caloMETMuonCorrs, caloMETMuonCorrsToken_, computeMiniIso_, computeMuonIDMVA_, computeMuonMVA_, computePuppiCombinedIso_, computeSoftMuonMVA_, CandIsolatorFromDeposits_cfi::deposits, PVValHelper::dxy, PVValHelper::dz, reco::PFCandidate::ecalEnergy(), effectiveAreaVec_, efficiencyLoader_, pat::Muon::embedCaloMETMuonCorrs(), embedCaloMETMuonCorrs_, embedCombinedMuon_, embedGenMatch_, embedHighLevel(), embedHighLevelSelection_, pat::Muon::embedPFCandidate(), embedPFCandidate_, embedPfEcalEnergy_, pat::Muon::embedTcMETMuonCorrs(), embedTcMETMuonCorrs_, pat::helper::EfficiencyLoader::enabled(), pat::helper::KinResolutionsLoader::enabled(), pat::helper::MultiIsolator::enabled(), pat::helper::MultiIsolator::endEvent(), pat::helper::MultiIsolator::fill(), fillHltTriggerInfo(), fillL1TriggerInfo(), fillMuon(), genMatchTokens_, geometryToken_, edm::EventSetup::getData(), edm::EventSetup::getHandle(), getRelMiniIsoPUCorrected(), cutBasedMuonId_MuonPOG_V0_cff::globalTrack, hltCollectionFilters_, mps_fire::i, heavyIonCSV_trainingSettings::idx, iEvent, pat::Muon::initSimInfo(), susybsm::HSCParticleType::innerTrack, edm::Ref< C, T, F >::isAvailable(), edm::Ref< C, T, F >::isNonnull(), isoDepositTokens_, displacedMuonProducer_cff::isolationValues, isolationValueTokens_, isolator_, isolatorTmpStorage_, edm::HandleBase::isValid(), dqmiolumiharvest::j, reco::Muon::LowPtMvaLoose, reco::Muon::LowPtMvaMedium, muon::makeSelectorBitset(), reco::Muon::MiniIsoLoose, reco::Muon::MiniIsoMedium, reco::Muon::MiniIsoTight, reco::Muon::MiniIsoVeryTight, eostools::move(), reco::Muon::MultiIsoMedium, reco::PFCandidate::muonRef(), PDWG_BPHSkim_cff::muons, muonTimeExtraToken_, muonToken_, beam_dqm_sourceclient-live_cfg::mva, mvaBTagCollectionTag_, photons_cff::mvaID, reco::Muon::MvaIDwpMedium, reco::Muon::MvaIDwpTight, displacedMuonProducer_cff::mvaL1Corrector, mvaL1Corrector_, displacedMuonProducer_cff::mvaL1L2L3ResCorrector, mvaL1L2L3ResCorrector_, reco::Muon::MvaLoose, reco::Muon::MvaMedium, reco::Muon::MvaTight, mvaUseJec_, reco::Muon::MvaVTight, reco::Muon::MvaVVTight, pat::helper::EfficiencyLoader::newEvent(), pat::helper::KinResolutionsLoader::newEvent(), nhits, patMuonPutToken_, hemisphereProducer_cfi::patMuons, pcToken_, pfMuonToken_, BeamMonitor_cff::primaryVertex, edm::Handle< T >::product(), puppiCombinedIsolation(), PUPPIIsolation_charged_hadrons_, PUPPIIsolation_neutral_hadrons_, PUPPIIsolation_photons_, reco::Muon::PuppiIsoLoose, reco::Muon::PuppiIsoMedium, reco::Muon::PuppiIsoTight, PUPPINoLeptonsIsolation_charged_hadrons_, PUPPINoLeptonsIsolation_neutral_hadrons_, PUPPINoLeptonsIsolation_photons_, MetAnalyzer::pv(), pat::Muon::PV2D, pat::Muon::PV3D, pvToken_, pat::Muon::readTimeExtra(), recomputeBasicSelectors_, resolutionLoader_, rho_, pat::Muon::setExtSimType(), pat::Lepton< LeptonType >::setIsolation(), pat::Muon::setIsolationPUPPI(), pat::Muon::setIsolationPUPPINoLeptons(), setMuonMiniIso(), pat::Muon::setNormChi2(), pat::Muon::setNumberOfValidHits(), pat::Muon::setPFCandidateRef(), pat::Muon::setPfEcalEnergy(), pat::Muon::setSimBX(), pat::Muon::setSimEta(), pat::Muon::setSimFlavour(), pat::Muon::setSimHeaviestMotherFlavour(), pat::Muon::setSimMatchQuality(), pat::Muon::setSimMotherPdgId(), pat::Muon::setSimPdgId(), pat::Muon::setSimPhi(), pat::Muon::setSimProdRho(), pat::Muon::setSimProdZ(), pat::Muon::setSimPt(), pat::Muon::setSimTpEvent(), pat::Muon::setSimType(), simInfo_, reco::Muon::SoftMvaId, jetsAK4_CHS_cff::sort, RandomServiceHelper::t1, RandomServiceHelper::t2, muonProducer_cfi::tcMETMuonCorrs, tcMETMuonCorrsToken_, transientTrackBuilderToken_, L1TEGammaOffline_cfi::triggerNames, triggerMatchMonitor_cfi::triggerObjects, triggerObjects_, triggerResults, triggerResults_, useParticleFlow_, userDataHelper_, and useUserData_.

582  {
583  // get the tracking Geometry
584  auto geometry = iSetup.getHandle(geometryToken_);
585  if (!geometry.isValid())
586  throw cms::Exception("FatalError") << "Unable to find GlobalTrackingGeometryRecord in event!\n";
587 
588  // switch off embedding (in unschedules mode)
589  if (iEvent.isRealData()) {
590  addGenMatch_ = false;
591  embedGenMatch_ = false;
592  }
593 
595  iEvent.getByToken(muonToken_, muons);
596 
599  iEvent.getByToken(pcToken_, pc);
600 
601  // get the ESHandle for the transient track builder,
602  // if needed for high level selection embedding
603  TransientTrackBuilder const* trackBuilder = nullptr;
604 
605  if (isolator_.enabled())
606  isolator_.beginEvent(iEvent, iSetup);
611 
613  for (size_t j = 0; j < isoDepositTokens_.size(); ++j) {
614  iEvent.getByToken(isoDepositTokens_[j], deposits[j]);
615  }
616 
618  for (size_t j = 0; j < isolationValueTokens_.size(); ++j) {
620  }
621 
622  //value maps for puppi isolation
623  edm::Handle<edm::ValueMap<float>> PUPPIIsolation_charged_hadrons;
624  edm::Handle<edm::ValueMap<float>> PUPPIIsolation_neutral_hadrons;
625  edm::Handle<edm::ValueMap<float>> PUPPIIsolation_photons;
626  //value maps for puppiNoLeptons isolation
627  edm::Handle<edm::ValueMap<float>> PUPPINoLeptonsIsolation_charged_hadrons;
628  edm::Handle<edm::ValueMap<float>> PUPPINoLeptonsIsolation_neutral_hadrons;
629  edm::Handle<edm::ValueMap<float>> PUPPINoLeptonsIsolation_photons;
630  if (addPuppiIsolation_) {
631  //puppi
632  iEvent.getByToken(PUPPIIsolation_charged_hadrons_, PUPPIIsolation_charged_hadrons);
633  iEvent.getByToken(PUPPIIsolation_neutral_hadrons_, PUPPIIsolation_neutral_hadrons);
634  iEvent.getByToken(PUPPIIsolation_photons_, PUPPIIsolation_photons);
635  //puppiNoLeptons
636  iEvent.getByToken(PUPPINoLeptonsIsolation_charged_hadrons_, PUPPINoLeptonsIsolation_charged_hadrons);
637  iEvent.getByToken(PUPPINoLeptonsIsolation_neutral_hadrons_, PUPPINoLeptonsIsolation_neutral_hadrons);
638  iEvent.getByToken(PUPPINoLeptonsIsolation_photons_, PUPPINoLeptonsIsolation_photons);
639  }
640 
641  // inputs for muon mva
642  edm::Handle<reco::JetTagCollection> mvaBTagCollectionTag;
645  if (computeMuonMVA_) {
646  iEvent.getByToken(mvaBTagCollectionTag_, mvaBTagCollectionTag);
649  }
650 
651  // prepare the MC genMatchTokens_
652  GenAssociations genMatches(genMatchTokens_.size());
653  if (addGenMatch_) {
654  for (size_t j = 0, nd = genMatchTokens_.size(); j < nd; ++j) {
655  iEvent.getByToken(genMatchTokens_[j], genMatches[j]);
656  }
657  }
658 
659  // prepare the high level selection: needs beamline
660  // OR primary vertex, depending on user selection
663  bool beamSpotIsValid = false;
664  bool primaryVertexIsValid = false;
666  // get the beamspot
667  edm::Handle<reco::BeamSpot> beamSpotHandle;
668  iEvent.getByToken(beamLineToken_, beamSpotHandle);
669 
670  // get the primary vertex
672  iEvent.getByToken(pvToken_, pvHandle);
673 
674  if (beamSpotHandle.isValid()) {
675  beamSpot = *beamSpotHandle;
676  beamSpotIsValid = true;
677  } else {
678  edm::LogError("DataNotAvailable") << "No beam spot available from EventSetup, not adding high level selection \n";
679  }
680  if (pvHandle.isValid() && !pvHandle->empty()) {
681  primaryVertex = pvHandle->at(0);
682  primaryVertexIsValid = true;
683  } else {
684  edm::LogError("DataNotAvailable")
685  << "No primary vertex available from EventSetup, not adding high level selection \n";
686  }
687  // this is needed by the IPTools methods from the tracking group
688  trackBuilder = &iSetup.getData(transientTrackBuilderToken_);
689  }
690 
691  // MC info
693  bool simInfoIsAvailalbe = iEvent.getByToken(simInfo_, simInfo);
694 
695  // this will be the new object collection
696  std::vector<Muon> patMuons;
697 
699  if (useParticleFlow_) {
700  // get the PFCandidates of type muons
701  iEvent.getByToken(pfMuonToken_, pfMuons);
702 
703  unsigned index = 0;
704  for (reco::PFCandidateConstIterator i = pfMuons->begin(); i != pfMuons->end(); ++i, ++index) {
705  const reco::PFCandidate& pfmu = *i;
706  //const reco::IsolaPFCandidate& pfmu = *i;
707  const reco::MuonRef& muonRef = pfmu.muonRef();
708  assert(muonRef.isNonnull());
709 
710  MuonBaseRef muonBaseRef(muonRef);
711  Muon aMuon(muonBaseRef);
712 
713  if (useUserData_) {
714  userDataHelper_.add(aMuon, iEvent, iSetup);
715  }
716 
717  // embed high level selection
719  // get the tracks
720  reco::TrackRef innerTrack = muonBaseRef->innerTrack();
721  reco::TrackRef globalTrack = muonBaseRef->globalTrack();
722  reco::TrackRef bestTrack = muonBaseRef->muonBestTrack();
723  reco::TrackRef chosenTrack = innerTrack;
724  // Make sure the collection it points to is there
725  if (bestTrack.isNonnull() && bestTrack.isAvailable())
726  chosenTrack = bestTrack;
727 
728  if (chosenTrack.isNonnull() && chosenTrack.isAvailable()) {
729  unsigned int nhits = chosenTrack->numberOfValidHits(); // ????
730  aMuon.setNumberOfValidHits(nhits);
731 
732  reco::TransientTrack tt = trackBuilder->build(chosenTrack);
733  embedHighLevel(aMuon, chosenTrack, tt, primaryVertex, primaryVertexIsValid, beamSpot, beamSpotIsValid);
734  }
735 
736  if (globalTrack.isNonnull() && globalTrack.isAvailable() && !embedCombinedMuon_) {
737  double norm_chi2 = globalTrack->chi2() / globalTrack->ndof();
738  aMuon.setNormChi2(norm_chi2);
739  }
740  }
741  reco::PFCandidateRef pfRef(pfMuons, index);
742  //reco::PFCandidatePtr ptrToMother(pfMuons,index);
743  reco::CandidateBaseRef pfBaseRef(pfRef);
744 
745  aMuon.setPFCandidateRef(pfRef);
746  if (embedPFCandidate_)
747  aMuon.embedPFCandidate();
748  fillMuon(aMuon, muonBaseRef, pfBaseRef, genMatches, deposits, isolationValues);
749 
750  if (computeMiniIso_)
751  setMuonMiniIso(aMuon, pc.product());
752 
753  if (addPuppiIsolation_) {
754  aMuon.setIsolationPUPPI((*PUPPIIsolation_charged_hadrons)[muonBaseRef],
755  (*PUPPIIsolation_neutral_hadrons)[muonBaseRef],
756  (*PUPPIIsolation_photons)[muonBaseRef]);
757 
758  aMuon.setIsolationPUPPINoLeptons((*PUPPINoLeptonsIsolation_charged_hadrons)[muonBaseRef],
759  (*PUPPINoLeptonsIsolation_neutral_hadrons)[muonBaseRef],
760  (*PUPPINoLeptonsIsolation_photons)[muonBaseRef]);
761  } else {
762  aMuon.setIsolationPUPPI(-999., -999., -999.);
763  aMuon.setIsolationPUPPINoLeptons(-999., -999., -999.);
764  }
765 
766  if (embedPfEcalEnergy_) {
767  aMuon.setPfEcalEnergy(pfmu.ecalEnergy());
768  }
769 
770  patMuons.push_back(aMuon);
771  }
772  } else {
774  iEvent.getByToken(muonToken_, muons);
775 
776  // embedding of muon MET corrections
778  //edm::ValueMap<reco::MuonMETCorrectionData> caloMETmuCorValueMap;
781  //caloMETmuCorValueMap = *caloMETmuCorValueMap_h;
782  }
784  //edm::ValueMap<reco::MuonMETCorrectionData> tcMETmuCorValueMap;
785  if (embedTcMETMuonCorrs_) {
787  //tcMETmuCorValueMap = *tcMETmuCorValueMap_h;
788  }
789 
791  // get the PFCandidates of type muons
792  iEvent.getByToken(pfMuonToken_, pfMuons);
793  }
794 
796  if (addInverseBeta_) {
797  // get MuonTimerExtra value map
798  iEvent.getByToken(muonTimeExtraToken_, muonsTimeExtra);
799  }
800 
801  for (edm::View<reco::Muon>::const_iterator itMuon = muons->begin(); itMuon != muons->end(); ++itMuon) {
802  // construct the Muon from the ref -> save ref to original object
803  unsigned int idx = itMuon - muons->begin();
804  MuonBaseRef muonRef = muons->refAt(idx);
805  reco::CandidateBaseRef muonBaseRef(muonRef);
806 
807  Muon aMuon(muonRef);
808  fillMuon(aMuon, muonRef, muonBaseRef, genMatches, deposits, isolationValues);
809  if (computeMiniIso_)
810  setMuonMiniIso(aMuon, pc.product());
811  if (addPuppiIsolation_) {
812  aMuon.setIsolationPUPPI((*PUPPIIsolation_charged_hadrons)[muonRef],
813  (*PUPPIIsolation_neutral_hadrons)[muonRef],
814  (*PUPPIIsolation_photons)[muonRef]);
815  aMuon.setIsolationPUPPINoLeptons((*PUPPINoLeptonsIsolation_charged_hadrons)[muonRef],
816  (*PUPPINoLeptonsIsolation_neutral_hadrons)[muonRef],
817  (*PUPPINoLeptonsIsolation_photons)[muonRef]);
818  } else {
819  aMuon.setIsolationPUPPI(-999., -999., -999.);
820  aMuon.setIsolationPUPPINoLeptons(-999., -999., -999.);
821  }
822 
823  // Isolation
824  if (isolator_.enabled()) {
825  //reco::CandidatePtr mother = ptrToMother->sourceCandidatePtr(0);
827  typedef pat::helper::MultiIsolator::IsolationValuePairs IsolationValuePairs;
828  // better to loop backwards, so the vector is resized less times
829  for (IsolationValuePairs::const_reverse_iterator it = isolatorTmpStorage_.rbegin(),
830  ed = isolatorTmpStorage_.rend();
831  it != ed;
832  ++it) {
833  aMuon.setIsolation(it->first, it->second);
834  }
835  }
836 
837  // for (size_t j = 0, nd = deposits.size(); j < nd; ++j) {
838  // aMuon.setIsoDeposit(isoDepositLabels_[j].first,
839  // (*deposits[j])[muonRef]);
840  // }
841 
842  // add sel to selected
843  edm::Ptr<reco::Muon> muonsPtr = muons->ptrAt(idx);
844  if (useUserData_) {
845  userDataHelper_.add(aMuon, iEvent, iSetup);
846  }
847 
848  // embed high level selection
850  // get the tracks
851  reco::TrackRef innerTrack = itMuon->innerTrack();
852  reco::TrackRef globalTrack = itMuon->globalTrack();
853  reco::TrackRef bestTrack = itMuon->muonBestTrack();
854  reco::TrackRef chosenTrack = innerTrack;
855  // Make sure the collection it points to is there
856  if (bestTrack.isNonnull() && bestTrack.isAvailable())
857  chosenTrack = bestTrack;
858  if (chosenTrack.isNonnull() && chosenTrack.isAvailable()) {
859  unsigned int nhits = chosenTrack->numberOfValidHits(); // ????
860  aMuon.setNumberOfValidHits(nhits);
861 
862  reco::TransientTrack tt = trackBuilder->build(chosenTrack);
863  embedHighLevel(aMuon, chosenTrack, tt, primaryVertex, primaryVertexIsValid, beamSpot, beamSpotIsValid);
864  }
865 
866  if (globalTrack.isNonnull() && globalTrack.isAvailable()) {
867  double norm_chi2 = globalTrack->chi2() / globalTrack->ndof();
868  aMuon.setNormChi2(norm_chi2);
869  }
870  }
871 
872  // embed MET muon corrections
874  aMuon.embedCaloMETMuonCorrs((*caloMETMuonCorrs)[muonRef]);
876  aMuon.embedTcMETMuonCorrs((*tcMETMuonCorrs)[muonRef]);
877 
879  if (embedPfEcalEnergy_)
880  aMuon.setPfEcalEnergy(-99.0);
881  unsigned index = 0;
882  for (const reco::PFCandidate& pfmu : *pfMuons) {
883  if (pfmu.muonRef().isNonnull()) {
884  if (pfmu.muonRef().id() != muonRef.id())
885  throw cms::Exception("Configuration")
886  << "Muon reference within PF candidates does not point to the muon collection." << std::endl;
887  if (pfmu.muonRef().key() == muonRef.key()) {
888  reco::PFCandidateRef pfRef(pfMuons, index);
889  aMuon.setPFCandidateRef(pfRef);
890  if (embedPfEcalEnergy_)
891  aMuon.setPfEcalEnergy(pfmu.ecalEnergy());
892  if (embedPFCandidate_)
893  aMuon.embedPFCandidate();
894  break;
895  }
896  }
897  index++;
898  }
899  }
900 
901  if (addInverseBeta_) {
902  aMuon.readTimeExtra((*muonsTimeExtra)[muonRef]);
903  }
904  // MC info
905  aMuon.initSimInfo();
906  if (simInfoIsAvailalbe) {
907  const auto& msi = (*simInfo)[muonBaseRef];
908  aMuon.setSimType(msi.primaryClass);
909  aMuon.setExtSimType(msi.extendedClass);
910  aMuon.setSimFlavour(msi.flavour);
911  aMuon.setSimHeaviestMotherFlavour(msi.heaviestMotherFlavour);
912  aMuon.setSimPdgId(msi.pdgId);
913  aMuon.setSimMotherPdgId(msi.motherPdgId);
914  aMuon.setSimBX(msi.tpBX);
915  aMuon.setSimTpEvent(msi.tpEvent);
916  aMuon.setSimProdRho(msi.vertex.Rho());
917  aMuon.setSimProdZ(msi.vertex.Z());
918  aMuon.setSimPt(msi.p4.pt());
919  aMuon.setSimEta(msi.p4.eta());
920  aMuon.setSimPhi(msi.p4.phi());
921  aMuon.setSimMatchQuality(msi.tpAssoQuality);
922  }
923  patMuons.push_back(aMuon);
924  }
925  }
926 
927  // sort muons in pt
928  std::sort(patMuons.begin(), patMuons.end(), [](auto const& t1, auto const& t2) { return t1.pt() > t2.pt(); });
929 
930  // Store standard muon selection decisions and jet related
931  // quantaties.
932  // Need a separate loop over muons to have all inputs properly
933  // computed and stored in the object.
935  if (computeMuonMVA_)
936  iEvent.getByToken(rho_, rho);
937  const reco::Vertex* pv(nullptr);
938  if (primaryVertexIsValid)
939  pv = &primaryVertex;
940 
943  bool triggerObjectsAvailable = false;
944  bool triggerResultsAvailable = false;
945  if (addTriggerMatching_) {
946  triggerObjectsAvailable = iEvent.getByToken(triggerObjects_, triggerObjects);
947  triggerResultsAvailable = iEvent.getByToken(triggerResults_, triggerResults);
948  }
949 
950  for (auto& muon : patMuons) {
951  // trigger info
952  if (addTriggerMatching_ and triggerObjectsAvailable and triggerResultsAvailable) {
953  const edm::TriggerNames& triggerNames(iEvent.triggerNames(*triggerResults));
956  }
957 
959  muon.setSelectors(0);
960  bool isRun2016BCDEF = (272728 <= iEvent.run() && iEvent.run() <= 278808);
961  muon.setSelectors(muon::makeSelectorBitset(muon, pv, isRun2016BCDEF));
962  }
963  float miniIsoValue = -1;
964  if (computeMiniIso_) {
965  // MiniIsolation working points
966 
968 
969  muon.setSelector(reco::Muon::MiniIsoLoose, miniIsoValue < 0.40);
970  muon.setSelector(reco::Muon::MiniIsoMedium, miniIsoValue < 0.20);
971  muon.setSelector(reco::Muon::MiniIsoTight, miniIsoValue < 0.10);
972  muon.setSelector(reco::Muon::MiniIsoVeryTight, miniIsoValue < 0.05);
973  }
974 
975  double puppiCombinedIsolationPAT = -1;
977  puppiCombinedIsolationPAT = puppiCombinedIsolation(muon, pc.product());
978  muon.setSelector(reco::Muon::PuppiIsoLoose, puppiCombinedIsolationPAT < 0.27);
979  muon.setSelector(reco::Muon::PuppiIsoMedium, puppiCombinedIsolationPAT < 0.22);
980  muon.setSelector(reco::Muon::PuppiIsoTight, puppiCombinedIsolationPAT < 0.12);
981  }
982 
983  float jetPtRatio = 0.0;
984  float jetPtRel = 0.0;
985  float mva = 0.0;
986  float mva_lowpt = 0.0;
987  if (computeMuonMVA_ && primaryVertexIsValid && computeMiniIso_) {
988  if (mvaUseJec_) {
989  mva = globalCache()->muonMvaEstimator().computeMva(muon,
991  *(mvaBTagCollectionTag.product()),
992  jetPtRatio,
993  jetPtRel,
994  miniIsoValue,
995  mvaL1Corrector.product(),
996  mvaL1L2L3ResCorrector.product());
997  mva_lowpt = globalCache()->muonLowPtMvaEstimator().computeMva(muon,
999  *(mvaBTagCollectionTag.product()),
1000  jetPtRatio,
1001  jetPtRel,
1002  miniIsoValue,
1003  mvaL1Corrector.product(),
1004  mvaL1L2L3ResCorrector.product());
1005 
1006  } else {
1007  mva = globalCache()->muonMvaEstimator().computeMva(
1008  muon, primaryVertex, *mvaBTagCollectionTag, jetPtRatio, jetPtRel, miniIsoValue);
1009  mva_lowpt = globalCache()->muonLowPtMvaEstimator().computeMva(
1010  muon, primaryVertex, *mvaBTagCollectionTag, jetPtRatio, jetPtRel, miniIsoValue);
1011  }
1012 
1013  muon.setMvaValue(mva);
1014  muon.setLowPtMvaValue(mva_lowpt);
1015  muon.setJetPtRatio(jetPtRatio);
1016  muon.setJetPtRel(jetPtRel);
1017 
1018  // multi-isolation
1019  if (computeMiniIso_) {
1020  muon.setSelector(reco::Muon::MultiIsoMedium,
1021  miniIsoValue < 0.11 && (muon.jetPtRatio() > 0.74 || muon.jetPtRel() > 6.8));
1022  }
1023 
1024  // MVA working points
1025  // https://twiki.cern.ch/twiki/bin/viewauth/CMS/LeptonMVA
1026  const double dB2D = std::abs(muon.dB(pat::Muon::PV2D));
1027  const double dB3D = std::abs(muon.dB(pat::Muon::PV3D));
1028  const double edB3D = std::abs(muon.edB(pat::Muon::PV3D));
1029  const double sip3D = edB3D > 0 ? dB3D / edB3D : 0.0;
1030  const double dz = std::abs(muon.muonBestTrack()->dz(primaryVertex.position()));
1031 
1032  // muon preselection
1033  if (muon.pt() > 5 and muon.isLooseMuon() and muon.passed(reco::Muon::MiniIsoLoose) and sip3D < 8.0 and
1034  dB2D < 0.05 and dz < 0.1) {
1035  muon.setSelector(reco::Muon::MvaLoose, muon.mvaValue() > -0.60);
1036  muon.setSelector(reco::Muon::MvaMedium, muon.mvaValue() > -0.20);
1037  muon.setSelector(reco::Muon::MvaTight, muon.mvaValue() > 0.15);
1038  muon.setSelector(reco::Muon::MvaVTight, muon.mvaValue() > 0.45);
1039  muon.setSelector(reco::Muon::MvaVVTight, muon.mvaValue() > 0.9);
1040  }
1041  if (muon.pt() > 5 and muon.isLooseMuon() and sip3D < 4 and dB2D < 0.5 and dz < 1) {
1042  muon.setSelector(reco::Muon::LowPtMvaLoose, muon.lowptMvaValue() > -0.60);
1043  muon.setSelector(reco::Muon::LowPtMvaMedium, muon.lowptMvaValue() > -0.20);
1044  }
1045  }
1046 
1047  // MVA ID
1048  float mvaID = 0.0;
1049  constexpr int MVAsentinelValue = -99;
1050  constexpr float mvaIDmediumCut = 0.08;
1051  constexpr float mvaIDtightCut = 0.12;
1052  if (computeMuonIDMVA_) {
1053  const double dz = std::abs(muon.muonBestTrack()->dz(primaryVertex.position()));
1054  const double dxy = std::abs(muon.muonBestTrack()->dxy(primaryVertex.position()));
1055  if (muon.isLooseMuon()) {
1056  mvaID = globalCache()->muonMvaIDEstimator().computeMVAID(muon)[1];
1057  } else {
1058  mvaID = MVAsentinelValue;
1059  }
1060  muon.setMvaIDValue(mvaID);
1061  muon.setSelector(reco::Muon::MvaIDwpMedium, muon.mvaIDValue() > mvaIDmediumCut);
1062  muon.setSelector(reco::Muon::MvaIDwpTight, muon.mvaIDValue() > mvaIDtightCut and dz < 0.5 and dxy < 0.2);
1063  }
1064 
1065  //SOFT MVA
1066  if (computeSoftMuonMVA_) {
1067  float mva = globalCache()->softMuonMvaEstimator().computeMva(muon);
1068  muon.setSoftMvaValue(mva);
1069  //preselection in SoftMuonMvaEstimator.cc
1070  muon.setSelector(reco::Muon::SoftMvaId, muon.softMvaValue() > 0.58); //WP choose for bmm4
1071  }
1072  }
1073 
1074  // put products in Event
1076 
1077  if (isolator_.enabled())
1078  isolator_.endEvent();
1079 }
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_charged_hadrons_
bool useUserData_
add user data to the muon (this will be data members of th muon even w/o embedding) ...
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_photons_
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
void newEvent(const edm::Event &event)
To be called for each new event, reads in the ValueMaps for efficiencies.
bool addPuppiIsolation_
add puppi isolation
edm::EDGetTokenT< edm::ValueMap< reco::MuonMETCorrectionData > > tcMETMuonCorrsToken_
source of tcMET muon corrections
bool embedTcMETMuonCorrs_
embed muon MET correction info for tcMET into the muon
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
edm::EDGetTokenT< std::vector< reco::Vertex > > pvToken_
input source of the primary vertex
bool enabled() const
True if it has a non null configuration.
Definition: MultiIsolator.h:55
edm::EDGetTokenT< edm::ValueMap< reco::MuonTimeExtra > > muonTimeExtraToken_
input tag for reading inverse beta
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_photons_
edm::EDGetTokenT< double > rho_
edm::EDGetTokenT< edm::View< reco::Muon > > muonToken_
input source
edm::EDGetTokenT< edm::ValueMap< reco::MuonMETCorrectionData > > caloMETMuonCorrsToken_
source of caloMET muon corrections
bool embedCaloMETMuonCorrs_
embed muon MET correction info for caloMET into the muon
edm::EDGetTokenT< edm::ValueMap< reco::MuonSimInfo > > simInfo_
MC info.
T const * product() const
Definition: Handle.h:70
void setMuonMiniIso(pat::Muon &aMuon, const pat::PackedCandidateCollection *pc)
bool addTriggerMatching_
Trigger.
edm::EDGetTokenT< reco::JetCorrector > mvaL1Corrector_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackBuilderToken_
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_charged_hadrons_
Log< level::Error, false > LogError
PFCandidateCollection::const_iterator PFCandidateConstIterator
iterator
reco::Muon::Selector makeSelectorBitset(reco::Muon const &muon, reco::Vertex const *vertex=nullptr, bool run2016_hip_mitigation=false)
assert(be >=bs)
void embedHighLevel(pat::Muon &aMuon, reco::TrackRef track, reco::TransientTrack &tt, reco::Vertex &primaryVertex, bool primaryVertexIsValid, reco::BeamSpot &beamspot, bool beamspotIsValid)
std::vector< double > effectiveAreaVec_
key_type key() const
Accessor for product key.
Definition: Ref.h:250
edm::EDGetTokenT< edm::ValueMap< float > > PUPPINoLeptonsIsolation_neutral_hadrons_
edm::EDGetTokenT< reco::JetTagCollection > mvaBTagCollectionTag_
edm::EDGetTokenT< reco::PFCandidateCollection > pfMuonToken_
input source pfCandidates that will be to be transformed into pat::Muons, when using PF2PAT ...
pat::helper::MultiIsolator isolator_
std::vector< edm::EDGetTokenT< edm::ValueMap< IsoDeposit > > > isoDepositTokens_
void fill(const edm::View< T > &coll, int idx, IsolationValuePairs &isolations) const
Definition: MultiIsolator.h:84
reco::TransientTrack build(const reco::Track *p) const
edm::EDGetTokenT< reco::BeamSpot > beamLineToken_
input source of the primary vertex/beamspot
bool addInverseBeta_
add combined inverse beta measurement into the muon
bool useParticleFlow_
switch to use particle flow (PF2PAT) or not
pat::PATUserDataHelper< pat::Muon > userDataHelper_
helper class to add userData to the muon
std::vector< edm::Handle< edm::ValueMap< IsoDeposit > > > IsoDepositMaps
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
int iEvent
Definition: GenABIO.cc:224
Definition: TTTypes.h:54
void beginEvent(const edm::Event &event, const edm::EventSetup &eventSetup)
bool addGenMatch_
add generator match information
Definition: Muon.py:1
void fillL1TriggerInfo(pat::Muon &muon, edm::Handle< std::vector< pat::TriggerObjectStandAlone >> &triggerObjects, const edm::TriggerNames &names, const edm::ESHandle< GlobalTrackingGeometry > &geometry)
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > geometryToken_
bool embedPfEcalEnergy_
add ecal PF energy
void newEvent(const edm::Event &event, const edm::EventSetup &setup)
To be called for each new event, reads in the EventSetup object.
bool isAvailable() const
Definition: Ref.h:537
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< pat::PackedCandidateCollection > pcToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static std::string const triggerResults
Definition: EdmProvDump.cc:47
void fillMuon(Muon &aMuon, const MuonBaseRef &muonRef, const reco::CandidateBaseRef &baseRef, const GenAssociations &genMatches, const IsoDepositMaps &deposits, const IsolationValueMaps &isolationValues) const
common muon filling, for both the standard and PF2PAT case
double puppiCombinedIsolation(const pat::Muon &muon, const pat::PackedCandidateCollection *pc)
double ecalEnergy() const
return corrected Ecal energy
Definition: PFCandidate.h:221
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
reco::MuonRef muonRef() const
Definition: PFCandidate.cc:448
std::vector< std::string > hltCollectionFilters_
void add(ObjectType &patObject, edm::Event const &iEvent, edm::EventSetup const &iSetup)
void fillHltTriggerInfo(pat::Muon &muon, edm::Handle< std::vector< pat::TriggerObjectStandAlone >> &triggerObjects, const edm::TriggerNames &names, const std::vector< std::string > &collection_names)
ProductIndex id() const
Definition: ProductID.h:35
std::vector< std::pair< pat::IsolationKeys, float > > IsolationValuePairs
Definition: MultiIsolator.h:17
double getRelMiniIsoPUCorrected(const pat::Muon &muon, double rho, const std::vector< double > &area)
bool isValid() const
Definition: HandleBase.h:70
std::vector< edm::Handle< edm::Association< reco::GenParticleCollection > > > GenAssociations
edm::EDGetTokenT< reco::JetCorrector > mvaL1L2L3ResCorrector_
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
std::vector< edm::Handle< edm::ValueMap< double > > > IsolationValueMaps
pat::helper::EfficiencyLoader efficiencyLoader_
helper class to add efficiencies to the muon
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const edm::EDPutTokenT< std::vector< Muon > > patMuonPutToken_
pat::helper::MultiIsolator::IsolationValuePairs isolatorTmpStorage_
isolation value pair for temporary storage before being folded into the muon
edm::EDGetTokenT< std::vector< pat::TriggerObjectStandAlone > > triggerObjects_
bool embedGenMatch_
embed the gen match information into the muon
std::vector< edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > > genMatchTokens_
input tags for generator match information
primaryVertex
hltOfflineBeamSpot for HLTMON
edm::RefToBase< reco::Muon > MuonBaseRef
typedefs for convenience
bool embedPFCandidate_
embed pfCandidates into the muon
std::vector< edm::EDGetTokenT< edm::ValueMap< double > > > isolationValueTokens_
bool computeMuonMVA_
standard muon selectors
bool embedCombinedMuon_
embed track of the combined fit into the muon
pat::helper::KinResolutionsLoader resolutionLoader_
helper class to add resolutions to the muon
edm::EDGetTokenT< edm::ValueMap< float > > PUPPIIsolation_neutral_hadrons_
def move(src, dest)
Definition: eostools.py:511
bool embedHighLevelSelection_
embed high level selection variables

◆ puppiCombinedIsolation()

double PATMuonProducer::puppiCombinedIsolation ( const pat::Muon muon,
const pat::PackedCandidateCollection pc 
)
private

Definition at line 1188 of file PATMuonProducer.cc.

References funct::abs(), L1DTConfigBti_cff::CH, reco::deltaPhi(), reco::deltaR2(), isChargedHadron(), isNeutralHadron(), isPhoton(), LogTrace, makePlotsFromDump::NH, dqm::qstatus::OTHER, and PbPb_ZMuSkimMuonDPG_cff::particleType.

Referenced by produce().

1188  {
1189  constexpr double dR_threshold = 0.4;
1190  constexpr double dR2_threshold = dR_threshold * dR_threshold;
1191  constexpr double mix_fraction = 0.5;
1192  enum particleType { CH = 0, NH = 1, PH = 2, OTHER = 100000 };
1193  double val_PuppiWithLep = 0.0;
1194  double val_PuppiWithoutLep = 0.0;
1195 
1196  for (const auto& cand : *pc) { //pat::pat::PackedCandidate loop start
1197 
1198  const particleType pType = isChargedHadron(cand.pdgId()) ? CH
1199  : isNeutralHadron(cand.pdgId()) ? NH
1200  : isPhoton(cand.pdgId()) ? PH
1201  : OTHER;
1202  if (pType == OTHER) {
1203  if (cand.pdgId() != 1 && cand.pdgId() != 2 && abs(cand.pdgId()) != 11 && abs(cand.pdgId()) != 13) {
1204  LogTrace("PATMuonProducer") << "candidate with PDGID = " << cand.pdgId()
1205  << " is not CH/NH/PH/e/mu or 1/2 (and this is removed from isolation calculation)"
1206  << std::endl;
1207  }
1208  continue;
1209  }
1210  double d_eta = std::abs(cand.eta() - muon.eta());
1211  if (d_eta > dR_threshold)
1212  continue;
1213 
1214  double d_phi = std::abs(reco::deltaPhi(cand.phi(), muon.phi()));
1215  if (d_phi > dR_threshold)
1216  continue;
1217 
1218  double dR2 = reco::deltaR2(cand, muon);
1219  if (dR2 > dR2_threshold)
1220  continue;
1221  if (pType == CH && dR2 < 0.0001 * 0.0001)
1222  continue;
1223  if (pType == NH && dR2 < 0.01 * 0.01)
1224  continue;
1225  if (pType == PH && dR2 < 0.01 * 0.01)
1226  continue;
1227  val_PuppiWithLep += cand.pt() * cand.puppiWeight();
1228  val_PuppiWithoutLep += cand.pt() * cand.puppiWeightNoLep();
1229 
1230  } //pat::pat::PackedCandidate loop end
1231 
1232  double reliso_Puppi_withLep = val_PuppiWithLep / muon.pt();
1233  double reliso_Puppi_withoutlep = val_PuppiWithoutLep / muon.pt();
1234  double reliso_Puppi_combined = mix_fraction * reliso_Puppi_withLep + (1.0 - mix_fraction) * reliso_Puppi_withoutlep;
1235  return reliso_Puppi_combined;
1236 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
static const int OTHER
bool isChargedHadron(long pdgid)
CH
LTS and SET for low trigger suppression.
#define LogTrace(id)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
bool isPhoton(long pdgid)
bool isNeutralHadron(long pdgid)

◆ readIsolationLabels()

template<typename T >
void pat::PATMuonProducer::readIsolationLabels ( const edm::ParameterSet iConfig,
const char *  psetName,
IsolationLabels labels,
std::vector< edm::EDGetTokenT< edm::ValueMap< T >>> &  tokens 
)
private

fill label vector from the contents of the parameter set, for the embedding of isoDeposits or userIsolation values

Definition at line 275 of file PATMuonProducer.cc.

References pat::EcalIso, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), pat::HcalIso, crabWrapper::key, label, SummaryClient_cfi::labels, pat::PfAllParticleIso, pat::PfChargedAllIso, pat::PfChargedHadronIso, pat::PfGammaIso, pat::PfNeutralHadronIso, pat::PfPUChargedHadronIso, crabTemplate::psetName, pat::TrackIso, pat::UserBaseIso, and edm::vector_transform().

278  {
279  labels.clear();
280 
281  if (iConfig.exists(psetName)) {
283 
284  if (depconf.exists("tracker"))
285  labels.emplace_back(pat::TrackIso, depconf.getParameter<edm::InputTag>("tracker"));
286  if (depconf.exists("ecal"))
287  labels.emplace_back(pat::EcalIso, depconf.getParameter<edm::InputTag>("ecal"));
288  if (depconf.exists("hcal"))
289  labels.emplace_back(pat::HcalIso, depconf.getParameter<edm::InputTag>("hcal"));
290  if (depconf.exists("pfAllParticles")) {
291  labels.emplace_back(pat::PfAllParticleIso, depconf.getParameter<edm::InputTag>("pfAllParticles"));
292  }
293  if (depconf.exists("pfChargedHadrons")) {
294  labels.emplace_back(pat::PfChargedHadronIso, depconf.getParameter<edm::InputTag>("pfChargedHadrons"));
295  }
296  if (depconf.exists("pfChargedAll")) {
297  labels.emplace_back(pat::PfChargedAllIso, depconf.getParameter<edm::InputTag>("pfChargedAll"));
298  }
299  if (depconf.exists("pfPUChargedHadrons")) {
300  labels.emplace_back(pat::PfPUChargedHadronIso, depconf.getParameter<edm::InputTag>("pfPUChargedHadrons"));
301  }
302  if (depconf.exists("pfNeutralHadrons")) {
303  labels.emplace_back(pat::PfNeutralHadronIso, depconf.getParameter<edm::InputTag>("pfNeutralHadrons"));
304  }
305  if (depconf.exists("pfPhotons")) {
306  labels.emplace_back(pat::PfGammaIso, depconf.getParameter<edm::InputTag>("pfPhotons"));
307  }
308  if (depconf.exists("user")) {
309  std::vector<edm::InputTag> userdeps = depconf.getParameter<std::vector<edm::InputTag>>("user");
310  std::vector<edm::InputTag>::const_iterator it = userdeps.begin(), ed = userdeps.end();
312  for (; it != ed; ++it, ++key) {
313  labels.push_back(std::make_pair(pat::IsolationKeys(key), *it));
314  }
315  tokens = edm::vector_transform(
316  labels, [this](IsolationLabel const& label) { return consumes<edm::ValueMap<T>>(label.second); });
317  }
318  }
320  return consumes<edm::ValueMap<T>>(label.second);
321  });
322 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool exists(std::string const &parameterName) const
checks if a parameter exists
IsolationKeys
Enum defining isolation keys.
Definition: Isolation.h:9
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
char const * label
std::pair< pat::IsolationKeys, edm::InputTag > IsolationLabel

◆ relMiniIsoPUCorrected()

double pat::PATMuonProducer::relMiniIsoPUCorrected ( const pat::Muon aMuon,
double  rho 
)
private

◆ setMuonMiniIso()

void PATMuonProducer::setMuonMiniIso ( pat::Muon aMuon,
const pat::PackedCandidateCollection pc 
)
private

Definition at line 1165 of file PATMuonProducer.cc.

References pat::getMiniPFIsolation(), miniIsoParams_, reco::LeafCandidate::polarP4(), and pat::Lepton< LeptonType >::setMiniPFIsolation().

Referenced by produce().

1165  {
1167  aMuon.polarP4(),
1168  miniIsoParams_[0],
1169  miniIsoParams_[1],
1170  miniIsoParams_[2],
1171  miniIsoParams_[3],
1172  miniIsoParams_[4],
1173  miniIsoParams_[5],
1174  miniIsoParams_[6],
1175  miniIsoParams_[7],
1176  miniIsoParams_[8]);
1177  aMuon.setMiniPFIsolation(miniiso);
1178 }
std::vector< double > miniIsoParams_
const PolarLorentzVector & polarP4() const final
four-momentum Lorentz vector
void setMiniPFIsolation(PFIsolation const &iso)
Definition: Lepton.h:217
PFIsolation getMiniPFIsolation(const pat::PackedCandidateCollection *pfcands, const reco::Candidate::PolarLorentzVector &p4, float mindr=0.05, float maxdr=0.2, float kt_scale=10.0, float ptthresh=0.5, float deadcone_ch=0.0001, float deadcone_pu=0.01, float deadcone_ph=0.01, float deadcone_nh=0.01, float dZ_cut=0.0)

Member Data Documentation

◆ addEfficiencies_

bool pat::PATMuonProducer::addEfficiencies_
private

add efficiencies to the muon (this will be data members of th muon even w/o embedding)

Definition at line 221 of file PATMuonProducer.cc.

◆ addGenMatch_

bool pat::PATMuonProducer::addGenMatch_
private

add generator match information

Definition at line 193 of file PATMuonProducer.cc.

Referenced by fillMuon(), and produce().

◆ addInverseBeta_

bool pat::PATMuonProducer::addInverseBeta_
private

add combined inverse beta measurement into the muon

Definition at line 189 of file PATMuonProducer.cc.

Referenced by produce().

◆ addPuppiIsolation_

bool pat::PATMuonProducer::addPuppiIsolation_
private

add puppi isolation

Definition at line 227 of file PATMuonProducer.cc.

Referenced by produce().

◆ addResolutions_

bool pat::PATMuonProducer::addResolutions_
private

add resolutions to the muon (this will be data members of th muon even w/o embedding)

Definition at line 199 of file PATMuonProducer.cc.

◆ addTriggerMatching_

bool pat::PATMuonProducer::addTriggerMatching_
private

Trigger.

Definition at line 261 of file PATMuonProducer.cc.

Referenced by produce().

◆ beamLineToken_

edm::EDGetTokenT<reco::BeamSpot> pat::PATMuonProducer::beamLineToken_
private

input source of the primary vertex/beamspot

Definition at line 211 of file PATMuonProducer.cc.

Referenced by produce().

◆ caloMETMuonCorrsToken_

edm::EDGetTokenT<edm::ValueMap<reco::MuonMETCorrectionData> > pat::PATMuonProducer::caloMETMuonCorrsToken_
private

source of caloMET muon corrections

Definition at line 177 of file PATMuonProducer.cc.

Referenced by produce().

◆ computeMiniIso_

bool pat::PATMuonProducer::computeMiniIso_
private

Definition at line 156 of file PATMuonProducer.cc.

Referenced by produce().

◆ computeMuonIDMVA_

bool pat::PATMuonProducer::computeMuonIDMVA_
private

Definition at line 238 of file PATMuonProducer.cc.

Referenced by produce().

◆ computeMuonMVA_

bool pat::PATMuonProducer::computeMuonMVA_
private

standard muon selectors

Definition at line 237 of file PATMuonProducer.cc.

Referenced by produce().

◆ computePuppiCombinedIso_

bool pat::PATMuonProducer::computePuppiCombinedIso_
private

Definition at line 157 of file PATMuonProducer.cc.

Referenced by produce().

◆ computeSoftMuonMVA_

bool pat::PATMuonProducer::computeSoftMuonMVA_
private

Definition at line 239 of file PATMuonProducer.cc.

Referenced by produce().

◆ effectiveAreaVec_

std::vector<double> pat::PATMuonProducer::effectiveAreaVec_
private

Definition at line 158 of file PATMuonProducer.cc.

Referenced by produce().

◆ efficiencyLoader_

pat::helper::EfficiencyLoader pat::PATMuonProducer::efficiencyLoader_
private

helper class to add efficiencies to the muon

Definition at line 253 of file PATMuonProducer.cc.

Referenced by fillMuon(), and produce().

◆ embedBestTrack_

bool pat::PATMuonProducer::embedBestTrack_
private

embed the track from best muon measurement (global pflow)

Definition at line 163 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ embedCaloMETMuonCorrs_

bool pat::PATMuonProducer::embedCaloMETMuonCorrs_
private

embed muon MET correction info for caloMET into the muon

Definition at line 175 of file PATMuonProducer.cc.

Referenced by produce().

◆ embedCombinedMuon_

bool pat::PATMuonProducer::embedCombinedMuon_
private

embed track of the combined fit into the muon

Definition at line 173 of file PATMuonProducer.cc.

Referenced by fillMuon(), and produce().

◆ embedDytMuon_

bool pat::PATMuonProducer::embedDytMuon_
private

embed track from DYT muon fit into the muon

Definition at line 187 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ embedGenMatch_

bool pat::PATMuonProducer::embedGenMatch_
private

embed the gen match information into the muon

Definition at line 197 of file PATMuonProducer.cc.

Referenced by fillMuon(), and produce().

◆ embedHighLevelSelection_

bool pat::PATMuonProducer::embedHighLevelSelection_
private

embed high level selection variables

Definition at line 209 of file PATMuonProducer.cc.

Referenced by produce().

◆ embedPFCandidate_

bool pat::PATMuonProducer::embedPFCandidate_
private

embed pfCandidates into the muon

Definition at line 207 of file PATMuonProducer.cc.

Referenced by produce().

◆ embedPfEcalEnergy_

bool pat::PATMuonProducer::embedPfEcalEnergy_
private

add ecal PF energy

Definition at line 225 of file PATMuonProducer.cc.

Referenced by produce().

◆ embedPickyMuon_

bool pat::PATMuonProducer::embedPickyMuon_
private

embed track from picky muon fit into the muon

Definition at line 183 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ embedStandAloneMuon_

bool pat::PATMuonProducer::embedStandAloneMuon_
private

embed track from muon system into the muon

Definition at line 171 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ embedTcMETMuonCorrs_

bool pat::PATMuonProducer::embedTcMETMuonCorrs_
private

embed muon MET correction info for tcMET into the muon

Definition at line 179 of file PATMuonProducer.cc.

Referenced by produce().

◆ embedTpfmsMuon_

bool pat::PATMuonProducer::embedTpfmsMuon_
private

embed track from tpfms muon fit into the muon

Definition at line 185 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ embedTrack_

bool pat::PATMuonProducer::embedTrack_
private

embed the track from inner tracker into the muon

Definition at line 169 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ embedTunePBestTrack_

bool pat::PATMuonProducer::embedTunePBestTrack_
private

embed the track from best muon measurement (muon only)

Definition at line 165 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ forceEmbedBestTrack_

bool pat::PATMuonProducer::forceEmbedBestTrack_
private

force separate embed of the best track even if already embedded

Definition at line 167 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ genMatchTokens_

std::vector<edm::EDGetTokenT<edm::Association<reco::GenParticleCollection> > > pat::PATMuonProducer::genMatchTokens_
private

input tags for generator match information

Definition at line 195 of file PATMuonProducer.cc.

Referenced by produce().

◆ geometryToken_

const edm::ESGetToken<GlobalTrackingGeometry, GlobalTrackingGeometryRecord> pat::PATMuonProducer::geometryToken_
private

Definition at line 266 of file PATMuonProducer.cc.

Referenced by produce().

◆ hltCollectionFilters_

std::vector<std::string> pat::PATMuonProducer::hltCollectionFilters_
private

Definition at line 264 of file PATMuonProducer.cc.

Referenced by produce().

◆ isoDepositLabels_

IsolationLabels pat::PATMuonProducer::isoDepositLabels_
private

input source for isoDeposits

Definition at line 215 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ isoDepositTokens_

std::vector<edm::EDGetTokenT<edm::ValueMap<IsoDeposit> > > pat::PATMuonProducer::isoDepositTokens_
private

Definition at line 216 of file PATMuonProducer.cc.

Referenced by produce().

◆ isolationValueLabels_

IsolationLabels pat::PATMuonProducer::isolationValueLabels_
private

input source isolation value maps

Definition at line 218 of file PATMuonProducer.cc.

Referenced by fillMuon().

◆ isolationValueTokens_

std::vector<edm::EDGetTokenT<edm::ValueMap<double> > > pat::PATMuonProducer::isolationValueTokens_
private

Definition at line 219 of file PATMuonProducer.cc.

Referenced by produce().

◆ isolator_

pat::helper::MultiIsolator pat::PATMuonProducer::isolator_
private

— tools — helper class to add userdefined isolation values to the muon

Definition at line 249 of file PATMuonProducer.cc.

Referenced by produce().

◆ isolatorTmpStorage_

pat::helper::MultiIsolator::IsolationValuePairs pat::PATMuonProducer::isolatorTmpStorage_
private

isolation value pair for temporary storage before being folded into the muon

Definition at line 251 of file PATMuonProducer.cc.

Referenced by produce().

◆ miniIsoParams_

std::vector<double> pat::PATMuonProducer::miniIsoParams_
private

Definition at line 159 of file PATMuonProducer.cc.

Referenced by getRelMiniIsoPUCorrected(), and setMuonMiniIso().

◆ muonTimeExtraToken_

edm::EDGetTokenT<edm::ValueMap<reco::MuonTimeExtra> > pat::PATMuonProducer::muonTimeExtraToken_
private

input tag for reading inverse beta

Definition at line 191 of file PATMuonProducer.cc.

Referenced by produce().

◆ muonToken_

edm::EDGetTokenT<edm::View<reco::Muon> > pat::PATMuonProducer::muonToken_
private

input source

Definition at line 152 of file PATMuonProducer.cc.

Referenced by produce().

◆ mvaBTagCollectionTag_

edm::EDGetTokenT<reco::JetTagCollection> pat::PATMuonProducer::mvaBTagCollectionTag_
private

Definition at line 242 of file PATMuonProducer.cc.

Referenced by produce().

◆ mvaL1Corrector_

edm::EDGetTokenT<reco::JetCorrector> pat::PATMuonProducer::mvaL1Corrector_
private

Definition at line 243 of file PATMuonProducer.cc.

Referenced by produce().

◆ mvaL1L2L3ResCorrector_

edm::EDGetTokenT<reco::JetCorrector> pat::PATMuonProducer::mvaL1L2L3ResCorrector_
private

Definition at line 244 of file PATMuonProducer.cc.

Referenced by produce().

◆ mvaUseJec_

bool pat::PATMuonProducer::mvaUseJec_
private

Definition at line 241 of file PATMuonProducer.cc.

Referenced by produce().

◆ patMuonPutToken_

const edm::EDPutTokenT<std::vector<Muon> > pat::PATMuonProducer::patMuonPutToken_
private

Definition at line 269 of file PATMuonProducer.cc.

Referenced by produce().

◆ pcToken_

edm::EDGetTokenT<pat::PackedCandidateCollection> pat::PATMuonProducer::pcToken_
private

Definition at line 155 of file PATMuonProducer.cc.

Referenced by produce().

◆ pfMuonToken_

edm::EDGetTokenT<reco::PFCandidateCollection> pat::PATMuonProducer::pfMuonToken_
private

input source pfCandidates that will be to be transformed into pat::Muons, when using PF2PAT

Definition at line 205 of file PATMuonProducer.cc.

Referenced by produce().

◆ PUPPIIsolation_charged_hadrons_

edm::EDGetTokenT<edm::ValueMap<float> > pat::PATMuonProducer::PUPPIIsolation_charged_hadrons_
private

Definition at line 229 of file PATMuonProducer.cc.

Referenced by produce().

◆ PUPPIIsolation_neutral_hadrons_

edm::EDGetTokenT<edm::ValueMap<float> > pat::PATMuonProducer::PUPPIIsolation_neutral_hadrons_
private

Definition at line 230 of file PATMuonProducer.cc.

Referenced by produce().

◆ PUPPIIsolation_photons_

edm::EDGetTokenT<edm::ValueMap<float> > pat::PATMuonProducer::PUPPIIsolation_photons_
private

Definition at line 231 of file PATMuonProducer.cc.

Referenced by produce().

◆ PUPPINoLeptonsIsolation_charged_hadrons_

edm::EDGetTokenT<edm::ValueMap<float> > pat::PATMuonProducer::PUPPINoLeptonsIsolation_charged_hadrons_
private

Definition at line 233 of file PATMuonProducer.cc.

Referenced by produce().

◆ PUPPINoLeptonsIsolation_neutral_hadrons_

edm::EDGetTokenT<edm::ValueMap<float> > pat::PATMuonProducer::PUPPINoLeptonsIsolation_neutral_hadrons_
private

Definition at line 234 of file PATMuonProducer.cc.

Referenced by produce().

◆ PUPPINoLeptonsIsolation_photons_

edm::EDGetTokenT<edm::ValueMap<float> > pat::PATMuonProducer::PUPPINoLeptonsIsolation_photons_
private

Definition at line 235 of file PATMuonProducer.cc.

Referenced by produce().

◆ pvToken_

edm::EDGetTokenT<std::vector<reco::Vertex> > pat::PATMuonProducer::pvToken_
private

input source of the primary vertex

Definition at line 213 of file PATMuonProducer.cc.

Referenced by produce().

◆ recomputeBasicSelectors_

bool pat::PATMuonProducer::recomputeBasicSelectors_
private

Definition at line 240 of file PATMuonProducer.cc.

Referenced by produce().

◆ relMiniIsoPUCorrected_

double pat::PATMuonProducer::relMiniIsoPUCorrected_
private

Definition at line 160 of file PATMuonProducer.cc.

◆ resolutionLoader_

pat::helper::KinResolutionsLoader pat::PATMuonProducer::resolutionLoader_
private

helper class to add resolutions to the muon

Definition at line 201 of file PATMuonProducer.cc.

Referenced by fillMuon(), and produce().

◆ rho_

edm::EDGetTokenT<double> pat::PATMuonProducer::rho_
private

Definition at line 245 of file PATMuonProducer.cc.

Referenced by produce().

◆ simInfo_

edm::EDGetTokenT<edm::ValueMap<reco::MuonSimInfo> > pat::PATMuonProducer::simInfo_
private

MC info.

Definition at line 258 of file PATMuonProducer.cc.

Referenced by produce().

◆ tcMETMuonCorrsToken_

edm::EDGetTokenT<edm::ValueMap<reco::MuonMETCorrectionData> > pat::PATMuonProducer::tcMETMuonCorrsToken_
private

source of tcMET muon corrections

Definition at line 181 of file PATMuonProducer.cc.

Referenced by produce().

◆ transientTrackBuilderToken_

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> pat::PATMuonProducer::transientTrackBuilderToken_
private

Definition at line 267 of file PATMuonProducer.cc.

Referenced by produce().

◆ triggerObjects_

edm::EDGetTokenT<std::vector<pat::TriggerObjectStandAlone> > pat::PATMuonProducer::triggerObjects_
private

Definition at line 262 of file PATMuonProducer.cc.

Referenced by produce().

◆ triggerResults_

edm::EDGetTokenT<edm::TriggerResults> pat::PATMuonProducer::triggerResults_
private

Definition at line 263 of file PATMuonProducer.cc.

Referenced by produce().

◆ useParticleFlow_

bool pat::PATMuonProducer::useParticleFlow_
private

switch to use particle flow (PF2PAT) or not

Definition at line 203 of file PATMuonProducer.cc.

Referenced by fillMuon(), and produce().

◆ userDataHelper_

pat::PATUserDataHelper<pat::Muon> pat::PATMuonProducer::userDataHelper_
private

helper class to add userData to the muon

Definition at line 255 of file PATMuonProducer.cc.

Referenced by produce().

◆ useUserData_

bool pat::PATMuonProducer::useUserData_
private

add user data to the muon (this will be data members of th muon even w/o embedding)

Definition at line 223 of file PATMuonProducer.cc.

Referenced by produce().