CMS 3D CMS Logo

ImpactParameterAlgorithm.cc
Go to the documentation of this file.
2 
5 
8 
10  ip_min = -9999;
11  ip_max = 9999;
12  sip_min = 0;
13  use_sign = false;
14  use3D = false;
15 }
16 
18  ip_min = parameters.getParameter<double>("TauImpactParameterMin");
19  ip_max = parameters.getParameter<double>("TauImpactParameterMax");
20  sip_min = parameters.getParameter<double>("TauImpactParameterSignificanceMin");
21  use_sign = parameters.getParameter<bool>("UseTauImpactParameterSign");
22  use3D = parameters.getParameter<bool>("UseTau3DImpactParameter");
23 }
24 
26  transientTrackBuilder = builder;
27 }
28 
29 std::pair<float, reco::TauImpactParameterInfo> ImpactParameterAlgorithm::tag(const reco::IsolatedTauTagInfoRef& tauRef,
30  const reco::Vertex& pv) {
31  if (transientTrackBuilder == nullptr) {
32  throw cms::Exception("NullTransientTrackBuilder") << "Transient track builder is 0. ";
33  }
34 
35  reco::TauImpactParameterInfo resultExtended;
36  resultExtended.setIsolatedTauTag(tauRef);
37 
38  const reco::Jet* jet = tauRef->jet().get();
39  GlobalVector direction(jet->px(), jet->py(), jet->pz());
40 
41  const reco::TrackRefVector& tracks = tauRef->selectedTracks();
42 
44  for (iTrack = tracks.begin(); iTrack != tracks.end(); iTrack++) {
45  const reco::TransientTrack transientTrack = (transientTrackBuilder->build(&(**iTrack)));
46 
48  Measurement1D ip = stip.apply(transientTrack, direction, pv).second;
49 
50  SignedImpactParameter3D signed_ip3D;
51  Measurement1D ip3D = signed_ip3D.apply(transientTrack, direction, pv).second;
52  LogDebug("ImpactParameterAlgorithm::tag") << "check pv,ip3d " << pv.z() << " " << ip3D.value();
53  if (!use_sign) {
54  Measurement1D tmp2D(fabs(ip.value()), ip.error());
55  ip = tmp2D;
56 
57  Measurement1D tmp3D(fabs(ip3D.value()), ip3D.error());
58  ip3D = tmp3D;
59  }
60 
62 
63  theData.transverseIp = ip;
64  theData.ip3D = ip3D;
65  resultExtended.storeTrackData(*iTrack, theData);
66  }
67 
68  float discriminator = resultExtended.discriminator(ip_min, ip_max, sip_min, use_sign, use3D);
69 
70  return std::make_pair(discriminator, resultExtended);
71 }
Vector3DBase
Definition: Vector3DBase.h:8
BeamSpotPI::parameters
parameters
Definition: BeamSpotPayloadInspectorHelper.h:30
SignedTransverseImpactParameter::apply
std::pair< bool, Measurement1D > apply(const reco::TransientTrack &, const GlobalVector &, const reco::Vertex &) const
Definition: SignedTransverseImpactParameter.cc:20
Measurement1D
Definition: Measurement1D.h:11
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
MessageLogger.h
SignedImpactParameter3D
Definition: SignedImpactParameter3D.h:13
Measurement1D::value
double value() const
Definition: Measurement1D.h:25
ImpactParameterAlgorithm::use3D
bool use3D
Definition: ImpactParameterAlgorithm.h:41
reco::TauImpactParameterInfo::setIsolatedTauTag
void setIsolatedTauTag(const IsolatedTauTagInfoRef &)
Definition: TauImpactParameterInfo.cc:45
edm::Ref::get
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
ImpactParameterAlgorithm::sip_min
double sip_min
Definition: ImpactParameterAlgorithm.h:40
edm::RefVector< TrackCollection >
SignedTransverseImpactParameter
Definition: SignedTransverseImpactParameter.h:15
ImpactParameterAlgorithm::ImpactParameterAlgorithm
ImpactParameterAlgorithm()
Definition: ImpactParameterAlgorithm.cc:9
ImpactParameterAlgorithm::tag
std::pair< float, reco::TauImpactParameterInfo > tag(const reco::IsolatedTauTagInfoRef &, const reco::Vertex &)
Definition: ImpactParameterAlgorithm.cc:29
reco::TauImpactParameterTrackData::transverseIp
Measurement1D transverseIp
Definition: TauImpactParameterInfo.h:13
edm::Ref< IsolatedTauTagInfoCollection >
ImpactParameterAlgorithm::ip_min
double ip_min
Definition: ImpactParameterAlgorithm.h:40
PDWG_TauSkim_cff.discriminator
discriminator
Definition: PDWG_TauSkim_cff.py:7
reco::TauImpactParameterInfo
Definition: TauImpactParameterInfo.h:22
SignedImpactParameter3D::apply
std::pair< bool, Measurement1D > apply(const reco::TransientTrack &, const GlobalVector &direction, const reco::Vertex &vertex) const
Definition: SignedImpactParameter3D.cc:17
Measurement1D::error
double error() const
Definition: Measurement1D.h:27
TransientTrackBuilder.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:159
SignedImpactParameter3D.h
ImpactParameterAlgorithm.h
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
ImpactParameterAlgorithm::use_sign
bool use_sign
Definition: ImpactParameterAlgorithm.h:41
TransientTrackBuilder
Definition: TransientTrackBuilder.h:16
reco::TauImpactParameterInfo::storeTrackData
void storeTrackData(const reco::TrackRef &, const TauImpactParameterTrackData &)
Definition: TauImpactParameterInfo.cc:40
reco::TauImpactParameterTrackData
Definition: TauImpactParameterInfo.h:12
reco::TauImpactParameterTrackData::ip3D
Measurement1D ip3D
Definition: TauImpactParameterInfo.h:14
reco::TransientTrack
Definition: TransientTrack.h:19
ImpactParameterAlgorithm::setTransientTrackBuilder
void setTransientTrackBuilder(const TransientTrackBuilder *)
Definition: ImpactParameterAlgorithm.cc:25
metsig::jet
Definition: SignAlgoResolutions.h:47
Exception
Definition: hltDiff.cc:245
edm::RefVectorIterator
Definition: EDProductfwd.h:33
TransientTrackBuilder::build
reco::TransientTrack build(const reco::Track *p) const
Definition: TransientTrackBuilder.cc:20
reco::TauImpactParameterInfo::discriminator
float discriminator(double, double, double, bool, bool) const
Definition: TauImpactParameterInfo.cc:8
ImpactParameterAlgorithm::transientTrackBuilder
const TransientTrackBuilder * transientTrackBuilder
Definition: ImpactParameterAlgorithm.h:43
reco::Vertex
Definition: Vertex.h:35
ImpactParameterAlgorithm::ip_max
double ip_max
Definition: ImpactParameterAlgorithm.h:40
SignedTransverseImpactParameter.h