CMS 3D CMS Logo

PFTauTransverseImpactParameters.cc
Go to the documentation of this file.
1 /* class PFTauTransverseImpactParamters
2  * EDProducer of the
3  * authors: Ian M. Nugent
4  * This work is based on the impact parameter work by Rosamaria Venditti and reconstructing the 3 prong taus.
5  * The idea of the fully reconstructing the tau using a kinematic fit comes from
6  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
7  * work was continued by Ian M. Nugent and Vladimir Cherepanov.
8  * Thanks goes to Christian Veelken and Evan Klose Friis for their help and suggestions.
9  */
10 
11 
21 
30 
39 
43 #include "TMath.h"
44 
45 #include <memory>
46 
47 using namespace reco;
48 using namespace edm;
49 using namespace std;
50 
52  public:
53  enum CMSSWPerigee{aCurv=0,aTheta,aPhi,aTip,aLip};
54  explicit PFTauTransverseImpactParameters(const edm::ParameterSet& iConfig);
56  void produce(edm::Event&,const edm::EventSetup&) override;
57  private:
62 };
63 
65  PFTauToken_(consumes<std::vector<reco::PFTau> >(iConfig.getParameter<edm::InputTag>("PFTauTag"))),
66  PFTauPVAToken_(consumes<edm::AssociationVector<PFTauRefProd, std::vector<reco::VertexRef> > >(iConfig.getParameter<edm::InputTag>("PFTauPVATag"))),
67  PFTauSVAToken_(consumes<edm::AssociationVector<PFTauRefProd,std::vector<std::vector<reco::VertexRef> > > >(iConfig.getParameter<edm::InputTag>("PFTauSVATag"))),
68  useFullCalculation_(iConfig.getParameter<bool>("useFullCalculation"))
69 {
70  produces<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTauTransverseImpactParameterRef> > >();
71  produces<PFTauTransverseImpactParameterCollection>("PFTauTIP");
72 }
73 
75 
76 }
77 
79  // Obtain Collections
80  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
81  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",transTrackBuilder);
82 
84  iEvent.getByToken(PFTauToken_,Tau);
85 
87  iEvent.getByToken(PFTauPVAToken_,PFTauPVA);
88 
90  iEvent.getByToken(PFTauSVAToken_,PFTauSVA);
91 
92  // Set Association Map
93  auto AVPFTauTIP = std::make_unique< edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTauTransverseImpactParameterRef>>>(PFTauRefProd(Tau));
94  auto TIPCollection_out = std::make_unique<PFTauTransverseImpactParameterCollection>();
96 
97  // For each Tau Run Algorithim
98  if(Tau.isValid()) {
99  for(reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
100  reco::PFTauRef RefPFTau(Tau, iPFTau);
101  const reco::VertexRef PV=PFTauPVA->value(RefPFTau.key());
102  const std::vector<reco::VertexRef> SV=PFTauSVA->value(RefPFTau.key());
103  double dxy(-999), dxy_err(-999);
104  reco::Vertex::Point poca(0,0,0);
105  double ip3d(-999), ip3d_err(-999);
106  reco::Vertex::Point ip3d_poca(0,0,0);
107  if(RefPFTau->leadPFChargedHadrCand().isNonnull()){
108  const reco::Track* track = nullptr;
109  if(RefPFTau->leadPFChargedHadrCand()->trackRef().isNonnull())
110  track = RefPFTau->leadPFChargedHadrCand()->trackRef().get();
111  else if(RefPFTau->leadPFChargedHadrCand()->gsfTrackRef().isNonnull())
112  track = RefPFTau->leadPFChargedHadrCand()->gsfTrackRef().get();
113  if(track != nullptr){
115  reco::TransientTrack transTrk=transTrackBuilder->build(*track);
116  GlobalVector direction(RefPFTau->p4().px(), RefPFTau->p4().py(), RefPFTau->p4().pz()); //To compute sign of IP
117  std::pair<bool,Measurement1D> signed_IP2D = IPTools::signedTransverseImpactParameter(transTrk, direction, (*PV));
118  dxy=signed_IP2D.second.value();
119  dxy_err=signed_IP2D.second.error();
120  std::pair<bool,Measurement1D> signed_IP3D = IPTools::signedImpactParameter3D(transTrk, direction, (*PV));
121  ip3d=signed_IP3D.second.value();
122  ip3d_err=signed_IP3D.second.error();
123  TransverseImpactPointExtrapolator extrapolator(transTrk.field());
124  GlobalPoint pos = extrapolator.extrapolate(transTrk.impactPointState(), RecoVertex::convertPos(PV->position())).globalPosition();
125  poca=reco::Vertex::Point(pos.x(),pos.y(),pos.z());
126  AnalyticalImpactPointExtrapolator extrapolator3D(transTrk.field());
127  GlobalPoint pos3d = extrapolator3D.extrapolate(transTrk.impactPointState(),RecoVertex::convertPos(PV->position())).globalPosition();
128  ip3d_poca=reco::Vertex::Point(pos3d.x(),pos3d.y(),pos3d.z());
129  }
130  else{
131  dxy_err=track->d0Error();
132  dxy=track->dxy(PV->position());
133  ip3d_err=track->dzError(); //store dz, ip3d not available
134  ip3d=track->dz(PV->position()); //store dz, ip3d not available
135  }
136  }
137  }
138  if(!SV.empty()){
140  reco::Vertex::Point v(SV.at(0)->x()-PV->x(),SV.at(0)->y()-PV->y(),SV.at(0)->z()-PV->z());
141  for(int i=0;i<reco::Vertex::dimension;i++){
142  for(int j=0;j<reco::Vertex::dimension;j++){
143  cov(i,j)=SV.at(0)->covariance(i,j)+PV->covariance(i,j);
144  }
145  }
146  GlobalVector direction(RefPFTau->px(),RefPFTau->py(),RefPFTau->pz());
147  double vSig = SecondaryVertex::computeDist3d(*PV,*SV.at(0),direction,true).significance();
148  reco::PFTauTransverseImpactParameter TIPV(poca,dxy,dxy_err,ip3d_poca,ip3d,ip3d_err,PV,v,vSig,SV.at(0));
149  reco::PFTauTransverseImpactParameterRef TIPVRef=reco::PFTauTransverseImpactParameterRef(TIPRefProd_out,TIPCollection_out->size());
150  TIPCollection_out->push_back(TIPV);
151  AVPFTauTIP->setValue(iPFTau,TIPVRef);
152  }
153  else{
154  reco::PFTauTransverseImpactParameter TIPV(poca,dxy,dxy_err,ip3d_poca,ip3d,ip3d_err,PV);
155  reco::PFTauTransverseImpactParameterRef TIPVRef=reco::PFTauTransverseImpactParameterRef(TIPRefProd_out,TIPCollection_out->size());
156  TIPCollection_out->push_back(TIPV);
157  AVPFTauTIP->setValue(iPFTau,TIPVRef);
158  }
159  }
160  }
161  iEvent.put(std::move(TIPCollection_out),"PFTauTIP");
162  iEvent.put(std::move(AVPFTauTIP));
163 }
164 
reco::Vertex::Point convertPos(const GlobalPoint &p)
double d0Error() const
error on d0
Definition: TrackBase.h:802
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:50
reco::TransientTrack build(const reco::Track *p) const
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:71
const MagneticField * field() const
key_type key() const
Accessor for product key.
Definition: Ref.h:265
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
Definition: PFTauFwd.h:15
uint16_t size_type
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
edm::EDGetTokenT< std::vector< reco::PFTau > > PFTauToken_
edm::EDGetTokenT< edm::AssociationVector< PFTauRefProd, std::vector< std::vector< reco::VertexRef > > > > PFTauSVAToken_
edm::EDGetTokenT< edm::AssociationVector< PFTauRefProd, std::vector< reco::VertexRef > > > PFTauPVAToken_
int iEvent
Definition: GenABIO.cc:230
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:245
bool isValid() const
Definition: HandleBase.h:74
RefProd< PROD > getRefBeforePut()
Definition: Event.h:167
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:609
double dzError() const
error on dz
Definition: TrackBase.h:814
edm::Ref< PFTauTransverseImpactParameterCollection > PFTauTransverseImpactParameterRef
presistent reference to a PFTauTransverseImpactParameter
PFTauTransverseImpactParameters(const edm::ParameterSet &iConfig)
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:63
TrajectoryStateOnSurface impactPointState() const
std::vector< reco::PFTauTransverseImpactParameter > PFTauTransverseImpactParameterCollection
collection of PFTauTransverseImpactParameter objects
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:591
def move(src, dest)
Definition: eostools.py:510
void produce(edm::Event &, const edm::EventSetup &) override