CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
PFRecoTauDiscriminationByFlightPathSignificance Class Reference
Inheritance diagram for PFRecoTauDiscriminationByFlightPathSignificance:
TauDiscriminationProducerBase< TauType, TauDiscriminator > edm::stream::EDProducer<>

Public Member Functions

void beginEvent (const edm::Event &, const edm::EventSetup &) override
 
double discriminate (const reco::PFTauRef &) const override
 
 PFRecoTauDiscriminationByFlightPathSignificance (const ParameterSet &iConfig)
 
 ~PFRecoTauDiscriminationByFlightPathSignificance () override
 
- Public Member Functions inherited from TauDiscriminationProducerBase< TauType, TauDiscriminator >
virtual double discriminate (const TauRef &tau) const =0
 
virtual void endEvent (edm::Event &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 TauDiscriminationProducerBase (const edm::ParameterSet &iConfig)
 
 TauDiscriminationProducerBase ()
 
 ~TauDiscriminationProducerBase () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from TauDiscriminationProducerBase< TauType, TauDiscriminator >
static void fillProducerDescriptions (edm::ParameterSetDescription &desc)
 

Private Member Functions

double threeProngFlightPathSig (const PFTauRef &) const
 
double vertexSignificance (reco::Vertex const &, reco::Vertex const &, GlobalVector const &) const
 

Private Attributes

bool booleanOutput
 
double flightPathSig
 
const TransientTrackBuildertransientTrackBuilder
 
reco::tau::RecoTauVertexAssociatorvertexAssociator_
 
bool withPVError
 

Additional Inherited Members

- Public Types inherited from TauDiscriminationProducerBase< TauType, TauDiscriminator >
typedef std::vector< TauType > TauCollection
 
typedef edm::Ref< TauCollectionTauRef
 
typedef edm::RefProd< TauCollectionTauRefProd
 
- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Protected Attributes inherited from TauDiscriminationProducerBase< TauType, TauDiscriminator >
std::string moduleLabel_
 
double prediscriminantFailValue_
 
edm::EDGetTokenT< TauCollectionTau_token
 
size_t tauIndex_
 
edm::InputTag TauProducer_
 

Detailed Description

Definition at line 27 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.

Constructor & Destructor Documentation

PFRecoTauDiscriminationByFlightPathSignificance::PFRecoTauDiscriminationByFlightPathSignificance ( const ParameterSet iConfig)
inlineexplicit

Definition at line 30 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.

References CaloRecoTauDiscriminationByFlightPathSignificance_cfi::flightPathSig, and edm::ParameterSet::getParameter().

32  flightPathSig = iConfig.getParameter<double>("flightPathSig");
33  withPVError = iConfig.getParameter<bool>("UsePVerror");
34  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
35  // edm::ConsumesCollector iC(consumesCollector());
36  vertexAssociator_ = new reco::tau::RecoTauVertexAssociator(iConfig.getParameter<ParameterSet>("qualityCuts"),consumesCollector());
37  }
T getParameter(std::string const &) const
TauDiscriminationProducerBase< reco::PFTau, reco::PFTauDiscriminator > PFTauDiscriminationProducerBase
PFRecoTauDiscriminationByFlightPathSignificance::~PFRecoTauDiscriminationByFlightPathSignificance ( )
inlineoverride

Definition at line 39 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.

References fillDescriptions().

39 {}

Member Function Documentation

void PFRecoTauDiscriminationByFlightPathSignificance::beginEvent ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Reimplemented from TauDiscriminationProducerBase< TauType, TauDiscriminator >.

Definition at line 58 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.

References edm::EventSetup::get(), and edm::ESHandle< T >::product().

59  {
60 
61  vertexAssociator_->setEvent(iEvent);
62 
63  // Transient Tracks
65  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
66  transientTrackBuilder = builder.product();
67 
68 }
void setEvent(const edm::Event &evt)
Load the vertices from the event.
T get() const
Definition: EventSetup.h:71
T const * product() const
Definition: ESHandle.h:86
double PFRecoTauDiscriminationByFlightPathSignificance::discriminate ( const reco::PFTauRef tau) const
override
void PFRecoTauDiscriminationByFlightPathSignificance::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 123 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.

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

123  {
124  // pfRecoTauDiscriminationByFlightPathSignificance
126  desc.add<double>("flightPathSig", 1.5);
127  desc.add<edm::InputTag>("PVProducer", edm::InputTag("offlinePrimaryVertices"));
128  desc.add<bool>("BooleanOutput", true);
129  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
130  {
132  psd0.add<std::string>("BooleanOperator", "and");
133  {
135  psd1.add<double>("cut");
136  psd1.add<edm::InputTag>("Producer");
137  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
138  }
139  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
140  }
141  desc.add<bool>("UsePVerror", true);
142  descriptions.add("pfRecoTauDiscriminationByFlightPathSignificance", desc);
143 }
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double PFRecoTauDiscriminationByFlightPathSignificance::threeProngFlightPathSig ( const PFTauRef tau) const
private

Definition at line 76 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.

References edm::Ref< C, T, F >::isNull(), TransientVertex::isValid(), impactParameterTagInfos_cfi::primaryVertex, and KalmanVertexFitter::vertex().

77  {
78  double flightPathSignificance = 0;
79 
81 
82  if (primaryVertex.isNull()) {
83  edm::LogError("FlightPathSignficance") << "Could not get vertex associated"
84  << " to tau, returning -999!" << std::endl;
85  return -999;
86  }
87 
88  //Secondary vertex
89  vector<TransientTrack> transientTracks;
90  for(const auto& pfSignalCand : tau->signalPFChargedHadrCands()){
91  if(pfSignalCand->trackRef().isNonnull()){
92  const TransientTrack transientTrack = transientTrackBuilder->build(pfSignalCand->trackRef());
93  transientTracks.push_back(transientTrack);
94  }
95  else if(pfSignalCand->gsfTrackRef().isNonnull()){
96  const TransientTrack transientTrack = transientTrackBuilder->build(pfSignalCand->gsfTrackRef());
97  transientTracks.push_back(transientTrack);
98  }
99  }
100  if(transientTracks.size() > 1){
101  KalmanVertexFitter kvf(true);
102  TransientVertex tv = kvf.vertex(transientTracks);
103 
104  if(tv.isValid()){
105  GlobalVector tauDir(tau->px(),
106  tau->py(),
107  tau->pz());
108  Vertex secVer = tv;
109  // We have to un-const the PV for some reason
110  reco::Vertex primaryVertexNonConst = *primaryVertex;
111  flightPathSignificance = vertexSignificance(primaryVertexNonConst,secVer,tauDir);
112  }
113  }
114  return flightPathSignificance;
115 }
double vertexSignificance(reco::Vertex const &, reco::Vertex const &, GlobalVector const &) const
reco::TransientTrack build(const reco::Track *p) const
bool isNull() const
Checks for null.
Definition: Ref.h:248
reco::VertexRef associatedVertex(const Jet &jet) const
bool isValid() const
double PFRecoTauDiscriminationByFlightPathSignificance::vertexSignificance ( reco::Vertex const &  pv,
reco::Vertex const &  sv,
GlobalVector const &  direction 
) const
private

Member Data Documentation

bool PFRecoTauDiscriminationByFlightPathSignificance::booleanOutput
private
double PFRecoTauDiscriminationByFlightPathSignificance::flightPathSig
private
const TransientTrackBuilder* PFRecoTauDiscriminationByFlightPathSignificance::transientTrackBuilder
private
reco::tau::RecoTauVertexAssociator* PFRecoTauDiscriminationByFlightPathSignificance::vertexAssociator_
private
bool PFRecoTauDiscriminationByFlightPathSignificance::withPVError
private