CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFTau3ProngReco.cc
Go to the documentation of this file.
1 /* class PFTau3ProngReco
2  * EDProducer of the
3  * author: Ian M. Nugent
4  * The idea of the fully reconstructing the tau using a kinematic fit comes from
5  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
6  * code is a result of the continuation of this work by Ian M. Nugent and Vladimir Cherepanov.
7  */
8 
9 
19 
24 
39 
45 
48 #include <memory>
49 #include <boost/foreach.hpp>
50 #include <TFormula.h>
51 
52 #include <memory>
53 
64 #include "TLorentzVector.h"
65 
67 
68 using namespace reco;
69 using namespace edm;
70 using namespace std;
71 using namespace tauImpactParameter;
72 
73 class PFTau3ProngReco : public EDProducer {
74  public:
75  enum Alg{useKalmanFit=0, useTrackHelix};
76 
77  struct DiscCutPair{
78  DiscCutPair():cutFormula_(0){}
79  ~DiscCutPair(){delete cutFormula_;}
82  double cut_;
83  TFormula* cutFormula_;
84  };
85  typedef std::vector<DiscCutPair*> DiscCutPairVec;
86 
87  explicit PFTau3ProngReco(const edm::ParameterSet& iConfig);
88  ~PFTau3ProngReco();
89  virtual void produce(edm::Event&,const edm::EventSetup&);
90  private:
95  std::auto_ptr<StringCutObjectSelector<reco::PFTau> > cut_;
96  int ndfPVT_;
98 };
99 
101  PFTauTag_(iConfig.getParameter<edm::InputTag>("PFTauTag")),
102  PFTauTIPTag_(iConfig.getParameter<edm::InputTag>("PFTauTIPTag")),
103  Algorithm_(iConfig.getParameter<int>("Algorithm")),
104  ndfPVT_(iConfig.getUntrackedParameter("ndfPVT",(int)5))
105 {
107  std::vector<edm::ParameterSet> discriminators =iConfig.getParameter<std::vector<edm::ParameterSet> >("discriminators");
108  // Build each of our cuts
109  BOOST_FOREACH(const edm::ParameterSet &pset, discriminators) {
110  DiscCutPair* newCut = new DiscCutPair();
111  newCut->inputTag_ = pset.getParameter<edm::InputTag>("discriminator");
112  if ( pset.existsAs<std::string>("selectionCut") ) newCut->cutFormula_ = new TFormula("selectionCut", pset.getParameter<std::string>("selectionCut").data());
113  else newCut->cut_ = pset.getParameter<double>("selectionCut");
114  discriminators_.push_back(newCut);
115  }
116  // Build a string cut if desired
117  if (iConfig.exists("cut")) cut_.reset(new StringCutObjectSelector<reco::PFTau>(iConfig.getParameter<std::string>( "cut" )));
119  produces<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef> > >();
120  produces<PFTau3ProngSummaryCollection>("PFTau3ProngSummary");
121 }
122 
124 
125 }
126 
128  // Obtain Collections
129  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
130  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",transTrackBuilder);
131 
133  iEvent.getByLabel(PFTauTag_,Tau);
134 
136  iEvent.getByLabel(PFTauTIPTag_,TIPAV);
137 
138  auto_ptr<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef> > > AVPFTau3PS(new edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef> >(PFTauRefProd(Tau)));
139  std::auto_ptr<PFTau3ProngSummaryCollection> PFTau3PSCollection_out= std::auto_ptr<PFTau3ProngSummaryCollection>(new PFTau3ProngSummaryCollection());
140  reco::PFTau3ProngSummaryRefProd PFTau3RefProd_out = iEvent.getRefBeforePut<reco::PFTau3ProngSummaryCollection>("PFTau3ProngSummary");
141 
142  // Load each discriminator
143  BOOST_FOREACH(DiscCutPair *disc, discriminators_) {iEvent.getByLabel(disc->inputTag_, disc->handle_);}
144 
145  // For each Tau Run Algorithim
146  if(Tau.isValid()){
147  for(reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
148  reco::PFTauRef tau(Tau, iPFTau);
149  reco::PFTau3ProngSummary PFTau3PS;
151  // Check if it passed all the discrimiantors
152  bool passed(true);
153  BOOST_FOREACH(const DiscCutPair* disc, discriminators_) {
154  // Check this discriminator passes
155  bool passedDisc = true;
156  if ( disc->cutFormula_ )passedDisc = (disc->cutFormula_->Eval((*disc->handle_)[tau]) > 0.5);
157  else passedDisc = ((*disc->handle_)[tau] > disc->cut_);
158  if ( !passedDisc ){passed = false; break;}
159  }
160  if (passed && cut_.get()){passed = (*cut_)(*tau);}
161  if (passed){
162  PDGInfo pdgInfo;
163  const reco::PFTauTransverseImpactParameterRef theTIP=TIPAV->value(tau.key());
164  const reco::VertexRef primaryVertex=theTIP->primaryVertex();
166  // Now compute the 3 prong Tau
167  bool SecondaryVtxOK(false);
169  if(theTIP->hasSecondaryVertex() && primaryVertex->ndof()>ndfPVT_){
170  const VertexRef secVtx=theTIP->secondaryVertex();
171  GlobalPoint sv(secVtx->position().x(),secVtx->position().y(),secVtx->position().z());
172  double vtxchi2(0), vtxndf(1);
174  vtxchi2=secVtx->chi2();
175  vtxndf=secVtx->ndof();
176  const std::vector<reco::Track>& selectedTracks=secVtx->refittedTracks();
177  std::vector<reco::TransientTrack> transTrkVect;
178  for(unsigned int i = 0; i!=selectedTracks.size();i++) transTrkVect.push_back(transTrackBuilder->build(selectedTracks[i]));
180  float piMassSigma(1.e-6), piChi(0.0), piNdf(0.0);
181  std::vector<RefCountedKinematicParticle> pions;
182  for(unsigned int i = 0; i<transTrkVect.size();i++) pions.push_back(kinFactory.particle(transTrkVect[i],pdgInfo.pi_mass(),piChi,piNdf,sv,piMassSigma));
183  RefCountedKinematicTree jpTree = kpvFitter_.fit(pions);
184  jpTree->movePointerToTheTop();
185  const KinematicParameters parameters = jpTree->currentParticle()->currentState().kinematicParameters();
186  AlgebraicSymMatrix77 cov=jpTree->currentParticle()->currentState().kinematicParametersError().matrix();
187  // get pions
188  double c(0);
189  std::vector<reco::Track> Tracks;
190  std::vector<LorentzVectorParticle> ReFitPions;
191  for(unsigned int i=0;i<transTrkVect.size();i++){
192  c+=transTrkVect[i].charge();
193  ReFitPions.push_back(ParticleBuilder::createLorentzVectorParticle(transTrkVect[i],*secVtx,true,true));
194  }
195  // now covert a1 into LorentzVectorParticle
196  TVectorT<double> a1_par(LorentzVectorParticle::NLorentzandVertexPar);
197  TMatrixTSym<double> a1_cov(LorentzVectorParticle::NLorentzandVertexPar);
198  for(int i = 0; i<LorentzVectorParticle::NLorentzandVertexPar; i++){
199  a1_par(i)=parameters(i);
200  for(int j = 0; j<LorentzVectorParticle::NLorentzandVertexPar; j++){
201  a1_cov(i,j)=cov(i,j);
202  }
203  }
204  a1=LorentzVectorParticle(a1_par,a1_cov,abs(PdtPdgMini::a_1_plus)*c,c,transTrackBuilder->field()->inInverseGeV(sv).z());
205  SecondaryVtxOK=true;
206  PFTau3PS=reco::PFTau3ProngSummary(theTIP,a1.p4(),vtxchi2,vtxndf);
207  }
208  else if(useTrackHelix==Algorithm_){
209  // use Track Helix
210  std::vector<TrackParticle> pions;
211  GlobalPoint pvpoint(primaryVertex->position().x(),primaryVertex->position().y(),primaryVertex->position().z());
212  const std::vector<edm::Ptr<reco::PFCandidate> > cands = tau->signalPFChargedHadrCands();
213  for (std::vector<edm::Ptr<reco::PFCandidate> >::const_iterator iter = cands.begin(); iter!=cands.end(); ++iter) {
214  if(iter->get()->trackRef().isNonnull()){
215  reco::TransientTrack transTrk=transTrackBuilder->build(iter->get()->trackRef());
216  pions.push_back(ParticleBuilder::createTrackParticle(transTrk,pvpoint,true,true));
217  }
218  else if(iter->get()->gsfTrackRef().isNonnull()){
219  //reco::TransientTrack transTrk=transTrackBuilder->build(iter->get()->gsfTrackRef());
220  //pions.push_back(ParticleBuilder::CreateTrackParticle(transTrk,pvpoint,true,true));
221  }
222  }
223  TVector3 pv(secVtx->position().x(),secVtx->position().y(),secVtx->position().z());
224  Chi2VertexFitter chi2v(pions,pv);
225  SecondaryVtxOK=chi2v.fit();
226  double c(0); for(unsigned int i=0;i<pions.size();i++){c+=pions[i].charge();}
227  int pdgid=abs(PdtPdgMini::a_1_plus)*c;
228  a1=chi2v.getMother(pdgid);
229  PFTau3PS =reco::PFTau3ProngSummary(theTIP,a1.p4(),vtxchi2,vtxndf);
230  }
231  }
232  if(SecondaryVtxOK){
233  // Tau Solver
234  TVector3 pv(primaryVertex->position().x(),primaryVertex->position().y(),primaryVertex->position().z());
235  TMatrixTSym<double> pvcov(LorentzVectorParticle::NVertex);
237  primaryVertex->fill(pvCov);
238  for(int i = 0; i<LorentzVectorParticle::NVertex; i++){
239  for(int j = 0; j<LorentzVectorParticle::NVertex; j++){
240  pvcov(i,j)=pvCov(i,j);
241  }
242  }
243  for(unsigned int i=0; i<PFTau3ProngSummary::nsolutions;i++){
244  TauA1NuConstrainedFitter TauA1NU(i,a1,pv,pvcov);
245  bool isFitOK=TauA1NU.fit();
246  if(isFitOK){
247  LorentzVectorParticle theTau=TauA1NU.getMother();
248  std::vector<LorentzVectorParticle> daughter=TauA1NU.getRefitDaughters();
249  std::vector<TLorentzVector> daughter_p4;
250  std::vector<int> daughter_charge,daughter_PDGID;
251  for(unsigned int d=0;d<daughter.size();d++){
252  daughter_p4.push_back(daughter[d].p4());
253  daughter_charge.push_back((int)daughter[d].charge());
254  daughter_PDGID.push_back(daughter[d].pdgId());
255  }
256  PFTau3PS.AddSolution(i,theTau.p4(),daughter_p4,daughter_charge,daughter_PDGID,(isFitOK),0.0,-999);
257  }
258  }
259  }
260  }
261  reco::PFTau3ProngSummaryRef PFTau3PSRef=reco::PFTau3ProngSummaryRef(PFTau3RefProd_out,PFTau3PSCollection_out->size());
262  PFTau3PSCollection_out->push_back(PFTau3PS);
263  AVPFTau3PS->setValue(iPFTau,PFTau3PSRef);
264  }
265  }
266  iEvent.put(PFTau3PSCollection_out,"PFTau3ProngSummary");
267  iEvent.put(AVPFTau3PS);
268 }
269 
std::vector< reco::PFTau3ProngSummary > PFTau3ProngSummaryCollection
collection of PFTau3ProngSummary objects
KinematicParticleVertexFitter kpvFitter_
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
PFTau3ProngReco(const edm::ParameterSet &iConfig)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::InputTag PFTauTIPTag_
ErrorD< N >::type type
Definition: Error.h:39
key_type key() const
Accessor for product key.
Definition: Ref.h:266
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
Definition: PFTauFwd.h:15
virtual void produce(edm::Event &, const edm::EventSetup &)
uint16_t size_type
std::auto_ptr< StringCutObjectSelector< reco::PFTau > > cut_
ROOT::Math::SMatrix< double, 7, 7, ROOT::Math::MatRepSym< double, 7 > > AlgebraicSymMatrix77
Definition: Matrices.h:16
std::vector< LorentzVectorParticle > getRefitDaughters()
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
edm::Handle< reco::PFTauDiscriminator > handle_
static double pi_mass()
Definition: PDGInfo.h:13
DiscCutPairVec discriminators_
bool isValid() const
Definition: HandleBase.h:76
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:402
RefProd< PROD > getRefBeforePut()
Definition: Event.h:133
edm::Ref< PFTau3ProngSummaryCollection > PFTau3ProngSummaryRef
presistent reference to a PFTau3ProngSummary
virtual bool AddSolution(unsigned int solution, const TLorentzVector &tau, const std::vector< TLorentzVector > &daughter_p4, const std::vector< int > &daughter_charge, const std::vector< int > &daughter_PDGID, bool has3ProngSolution, double solutionChi2, double thetaGJsig)
const T & get() const
Definition: EventSetup.h:55
std::vector< DiscCutPair * > DiscCutPairVec
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
edm::InputTag PFTauTag_