CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
31 
45 
49 #include "TMath.h"
50 
51 #include <memory>
52 
53 using namespace reco;
54 using namespace edm;
55 using namespace std;
56 
58  public:
59  enum Alg{useInputPV=0, useFont};
60  enum CMSSWPerigee{aCurv=0,aTheta,aPhi,aTip,aLip};
61  explicit PFTauTransverseImpactParameters(const edm::ParameterSet& iConfig);
63  virtual void produce(edm::Event&,const edm::EventSetup&);
64  private:
69 };
70 
72  PFTauToken_(consumes<std::vector<reco::PFTau> >(iConfig.getParameter<edm::InputTag>("PFTauTag"))),
73  PFTauPVAToken_(consumes<edm::AssociationVector<PFTauRefProd, std::vector<reco::VertexRef> > >(iConfig.getParameter<edm::InputTag>("PFTauPVATag"))),
74  PFTauSVAToken_(consumes<edm::AssociationVector<PFTauRefProd,std::vector<std::vector<reco::VertexRef> > > >(iConfig.getParameter<edm::InputTag>("PFTauSVATag"))),
75  useFullCalculation_(iConfig.getParameter<bool>("useFullCalculation"))
76 {
77  produces<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTauTransverseImpactParameterRef> > >();
78  produces<PFTauTransverseImpactParameterCollection>("PFTauTIP");
79 }
80 
82 
83 }
84 
86  // Obtain Collections
87  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
88  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",transTrackBuilder);
89 
91  iEvent.getByToken(PFTauToken_,Tau);
92 
94  iEvent.getByToken(PFTauPVAToken_,PFTauPVA);
95 
97  iEvent.getByToken(PFTauSVAToken_,PFTauSVA);
98 
99  // Set Association Map
100  auto_ptr<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTauTransverseImpactParameterRef> > > AVPFTauTIP(new edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTauTransverseImpactParameterRef> >(PFTauRefProd(Tau)));
101  std::auto_ptr<PFTauTransverseImpactParameterCollection> TIPCollection_out= std::auto_ptr<PFTauTransverseImpactParameterCollection>(new PFTauTransverseImpactParameterCollection());
103 
104 
105  // For each Tau Run Algorithim
106  if(Tau.isValid()) {
107  for(reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
108  reco::PFTauRef RefPFTau(Tau, iPFTau);
109  const reco::VertexRef PV=PFTauPVA->value(RefPFTau.key());
110  const std::vector<reco::VertexRef> SV=PFTauSVA->value(RefPFTau.key());
111  double dxy(-999), dxy_err(-999);
112  reco::Vertex::Point poca(0,0,0);
113  double ip3d(-999), ip3d_err(-999);
114  reco::Vertex::Point ip3d_poca(0,0,0);
115  if(RefPFTau->leadPFChargedHadrCand().isNonnull()){
116  if(RefPFTau->leadPFChargedHadrCand()->trackRef().isNonnull()){
118  reco::TransientTrack transTrk=transTrackBuilder->build(RefPFTau->leadPFChargedHadrCand()->trackRef());
119  GlobalVector direction(RefPFTau->p4().px(), RefPFTau->p4().py(), RefPFTau->p4().pz()); //To compute sign of IP
120  std::pair<bool,Measurement1D> signed_IP2D = IPTools::signedTransverseImpactParameter(transTrk, direction, (*PV));
121  dxy=signed_IP2D.second.value();
122  dxy_err=signed_IP2D.second.error();
123  std::pair<bool,Measurement1D> signed_IP3D = IPTools::signedImpactParameter3D(transTrk, direction, (*PV));
124  ip3d=signed_IP3D.second.value();
125  ip3d_err=signed_IP3D.second.error();
126  TransverseImpactPointExtrapolator extrapolator(transTrk.field());
127  GlobalPoint pos = extrapolator.extrapolate(transTrk.impactPointState(), RecoVertex::convertPos(PV->position())).globalPosition();
128  poca=reco::Vertex::Point(pos.x(),pos.y(),pos.z());
129  AnalyticalImpactPointExtrapolator extrapolator3D(transTrk.field());
130  GlobalPoint pos3d = extrapolator3D.extrapolate(transTrk.impactPointState(),RecoVertex::convertPos(PV->position())).globalPosition();
131  ip3d_poca=reco::Vertex::Point(pos3d.x(),pos3d.y(),pos3d.z());
132  }
133  else{
134  dxy_err=RefPFTau->leadPFChargedHadrCand()->trackRef()->d0Error();
135  dxy=RefPFTau->leadPFChargedHadrCand()->trackRef()->dxy(PV->position());
136  ip3d_err=RefPFTau->leadPFChargedHadrCand()->trackRef()->dzError(); //store dz, ip3d not available
137  ip3d=RefPFTau->leadPFChargedHadrCand()->trackRef()->dz(PV->position()); //store dz, ip3d not available
138  }
139  }
140  }
141  if(SV.size()>0){
143  reco::Vertex::Point v(SV.at(0)->x()-PV->x(),SV.at(0)->y()-PV->y(),SV.at(0)->z()-PV->z());
144  for(int i=0;i<reco::Vertex::dimension;i++){
145  for(int j=0;j<reco::Vertex::dimension;j++){
146  cov(i,j)=SV.at(0)->covariance(i,j)+PV->covariance(i,j);
147  }
148  }
149  GlobalVector direction(RefPFTau->px(),RefPFTau->py(),RefPFTau->pz());
150  double vSig = SecondaryVertex::computeDist3d(*PV,*SV.at(0),direction,true).significance();
151  reco::PFTauTransverseImpactParameter TIPV(poca,dxy,dxy_err,ip3d_poca,ip3d,ip3d_err,PV,v,vSig,SV.at(0));
152  reco::PFTauTransverseImpactParameterRef TIPVRef=reco::PFTauTransverseImpactParameterRef(TIPRefProd_out,TIPCollection_out->size());
153  TIPCollection_out->push_back(TIPV);
154  AVPFTauTIP->setValue(iPFTau,TIPVRef);
155  }
156  else{
157  reco::PFTauTransverseImpactParameter TIPV(poca,dxy,dxy_err,ip3d_poca,ip3d,ip3d_err,PV);
158  reco::PFTauTransverseImpactParameterRef TIPVRef=reco::PFTauTransverseImpactParameterRef(TIPRefProd_out,TIPCollection_out->size());
159  TIPCollection_out->push_back(TIPV);
160  AVPFTauTIP->setValue(iPFTau,TIPVRef);
161  }
162  }
163  }
164  iEvent.put(TIPCollection_out,"PFTauTIP");
165  iEvent.put(AVPFTauTIP);
166 }
167 
reco::Vertex::Point convertPos(const GlobalPoint &p)
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#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
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:264
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_
virtual void produce(edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
int j
Definition: DBlmapReader.cc:9
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
bool isValid() const
Definition: HandleBase.h:75
RefProd< PROD > getRefBeforePut()
Definition: Event.h:140
const T & get() const
Definition: EventSetup.h:56
edm::Ref< PFTauTransverseImpactParameterCollection > PFTauTransverseImpactParameterRef
presistent reference to a PFTauTransverseImpactParameter
PFTauTransverseImpactParameters(const edm::ParameterSet &iConfig)
TrajectoryStateOnSurface impactPointState() const
std::vector< reco::PFTauTransverseImpactParameter > PFTauTransverseImpactParameterCollection
collection of PFTauTransverseImpactParameter objects