CMS 3D CMS Logo

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 <TFormula.h>
49 
50 #include <memory>
51 
62 #include "TLorentzVector.h"
63 
65 
66 using namespace reco;
67 using namespace edm;
68 using namespace std;
69 using namespace tauImpactParameter;
70 
71 class PFTau3ProngReco : public EDProducer {
72  public:
73  enum Alg{useKalmanFit=0, useTrackHelix};
74 
75  struct DiscCutPair{
76  DiscCutPair():cutFormula_(nullptr){}
77  ~DiscCutPair(){delete cutFormula_;}
80  double cut_;
81  TFormula* cutFormula_;
82  };
83  typedef std::vector<DiscCutPair*> DiscCutPairVec;
84 
85  explicit PFTau3ProngReco(const edm::ParameterSet& iConfig);
86  ~PFTau3ProngReco() override;
87  void produce(edm::Event&,const edm::EventSetup&) override;
88  private:
92  DiscCutPairVec discriminators_;
93  std::unique_ptr<StringCutObjectSelector<reco::PFTau> > cut_;
94  int ndfPVT_;
96 };
97 
99  PFTauTag_(iConfig.getParameter<edm::InputTag>("PFTauTag")),
100  PFTauTIPTag_(iConfig.getParameter<edm::InputTag>("PFTauTIPTag")),
101  Algorithm_(iConfig.getParameter<int>("Algorithm")),
102  ndfPVT_(iConfig.getUntrackedParameter("ndfPVT",(int)5))
103 {
105  std::vector<edm::ParameterSet> discriminators =iConfig.getParameter<std::vector<edm::ParameterSet> >("discriminators");
106  // Build each of our cuts
107  for(auto const& pset : discriminators) {
108  DiscCutPair* newCut = new DiscCutPair();
109  newCut->inputTag_ = pset.getParameter<edm::InputTag>("discriminator");
110  if ( pset.existsAs<std::string>("selectionCut") ) newCut->cutFormula_ = new TFormula("selectionCut", pset.getParameter<std::string>("selectionCut").data());
111  else newCut->cut_ = pset.getParameter<double>("selectionCut");
112  discriminators_.push_back(newCut);
113  }
114  // Build a string cut if desired
115  if (iConfig.exists("cut")) cut_ = std::make_unique<StringCutObjectSelector<reco::PFTau>>(iConfig.getParameter<std::string>("cut"));
117  produces<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef> > >();
118  produces<PFTau3ProngSummaryCollection>("PFTau3ProngSummary");
119 }
120 
122 
123 }
124 
126  // Obtain Collections
127  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
128  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",transTrackBuilder);
129 
131  iEvent.getByLabel(PFTauTag_,Tau);
132 
134  iEvent.getByLabel(PFTauTIPTag_,TIPAV);
135 
136  auto AVPFTau3PS = std::make_unique<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef>>>(PFTauRefProd(Tau));
137  auto PFTau3PSCollection_out = std::make_unique<PFTau3ProngSummaryCollection>();
138  reco::PFTau3ProngSummaryRefProd PFTau3RefProd_out = iEvent.getRefBeforePut<reco::PFTau3ProngSummaryCollection>("PFTau3ProngSummary");
139 
140  // Load each discriminator
141  for(auto& disc : discriminators_) {iEvent.getByLabel(disc->inputTag_, disc->handle_);}
142 
143  // For each Tau Run Algorithim
144  if(Tau.isValid()){
145  for(reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
146  reco::PFTauRef tau(Tau, iPFTau);
147  reco::PFTau3ProngSummary PFTau3PS;
149  // Check if it passed all the discrimiantors
150  bool passed(true);
151  for(auto const& disc : discriminators_) {
152  // Check this discriminator passes
153  bool passedDisc = true;
154  if ( disc->cutFormula_ )passedDisc = (disc->cutFormula_->Eval((*disc->handle_)[tau]) > 0.5);
155  else passedDisc = ((*disc->handle_)[tau] > disc->cut_);
156  if ( !passedDisc ){passed = false; break;}
157  }
158  if (passed && cut_.get()){passed = (*cut_)(*tau);}
159  if (passed){
160  PDGInfo pdgInfo;
161  const reco::PFTauTransverseImpactParameterRef theTIP=TIPAV->value(tau.key());
162  const reco::VertexRef primaryVertex=theTIP->primaryVertex();
164  // Now compute the 3 prong Tau
165  bool SecondaryVtxOK(false);
167  if(theTIP->hasSecondaryVertex() && primaryVertex->ndof()>ndfPVT_){
168  const VertexRef secVtx=theTIP->secondaryVertex();
169  GlobalPoint sv(secVtx->position().x(),secVtx->position().y(),secVtx->position().z());
170  double vtxchi2(0), vtxndf(1);
172  vtxchi2=secVtx->chi2();
173  vtxndf=secVtx->ndof();
174  const std::vector<reco::Track>& selectedTracks=secVtx->refittedTracks();
175  std::vector<reco::TransientTrack> transTrkVect;
176  for(unsigned int i = 0; i!=selectedTracks.size();i++) transTrkVect.push_back(transTrackBuilder->build(selectedTracks[i]));
178  float piMassSigma(1.e-6), piChi(0.0), piNdf(0.0);
179  std::vector<RefCountedKinematicParticle> pions;
180  for(unsigned int i = 0; i<transTrkVect.size();i++) pions.push_back(kinFactory.particle(transTrkVect[i],pdgInfo.pi_mass(),piChi,piNdf,sv,piMassSigma));
181  RefCountedKinematicTree jpTree = kpvFitter_.fit(pions);
182  jpTree->movePointerToTheTop();
183  const KinematicParameters parameters = jpTree->currentParticle()->currentState().kinematicParameters();
184  AlgebraicSymMatrix77 cov=jpTree->currentParticle()->currentState().kinematicParametersError().matrix();
185  // get pions
186  double c(0);
187  std::vector<reco::Track> Tracks;
188  std::vector<LorentzVectorParticle> ReFitPions;
189  for(unsigned int i=0;i<transTrkVect.size();i++){
190  c+=transTrkVect[i].charge();
191  ReFitPions.push_back(ParticleBuilder::createLorentzVectorParticle(transTrkVect[i],*secVtx,true,true));
192  }
193  // now covert a1 into LorentzVectorParticle
194  TVectorT<double> a1_par(LorentzVectorParticle::NLorentzandVertexPar);
195  TMatrixTSym<double> a1_cov(LorentzVectorParticle::NLorentzandVertexPar);
196  for(int i = 0; i<LorentzVectorParticle::NLorentzandVertexPar; i++){
197  a1_par(i)=parameters(i);
198  for(int j = 0; j<LorentzVectorParticle::NLorentzandVertexPar; j++){
199  a1_cov(i,j)=cov(i,j);
200  }
201  }
202  a1=LorentzVectorParticle(a1_par,a1_cov,abs(PdtPdgMini::a_1_plus)*c,c,transTrackBuilder->field()->inInverseGeV(sv).z());
203  SecondaryVtxOK=true;
204  PFTau3PS=reco::PFTau3ProngSummary(theTIP,a1.p4(),vtxchi2,vtxndf);
205  }
206  else if(useTrackHelix==Algorithm_){
207  // use Track Helix
208  std::vector<TrackParticle> pions;
209  GlobalPoint pvpoint(primaryVertex->position().x(),primaryVertex->position().y(),primaryVertex->position().z());
210  const std::vector<edm::Ptr<reco::PFCandidate> > cands = tau->signalPFChargedHadrCands();
211  for (std::vector<edm::Ptr<reco::PFCandidate> >::const_iterator iter = cands.begin(); iter!=cands.end(); ++iter) {
212  if(iter->get()->trackRef().isNonnull()){
213  reco::TransientTrack transTrk=transTrackBuilder->build(iter->get()->trackRef());
214  pions.push_back(ParticleBuilder::createTrackParticle(transTrk,pvpoint,true,true));
215  }
216  else if(iter->get()->gsfTrackRef().isNonnull()){
217  //reco::TransientTrack transTrk=transTrackBuilder->build(iter->get()->gsfTrackRef());
218  //pions.push_back(ParticleBuilder::CreateTrackParticle(transTrk,pvpoint,true,true));
219  }
220  }
221  TVector3 pv(secVtx->position().x(),secVtx->position().y(),secVtx->position().z());
222  Chi2VertexFitter chi2v(pions,pv);
223  SecondaryVtxOK=chi2v.fit();
224  double c(0); for(unsigned int i=0;i<pions.size();i++){c+=pions[i].charge();}
226  a1=chi2v.getMother(pdgid);
227  PFTau3PS =reco::PFTau3ProngSummary(theTIP,a1.p4(),vtxchi2,vtxndf);
228  }
229  }
230  if(SecondaryVtxOK){
231  // Tau Solver
232  TVector3 pv(primaryVertex->position().x(),primaryVertex->position().y(),primaryVertex->position().z());
233  TMatrixTSym<double> pvcov(LorentzVectorParticle::NVertex);
235  primaryVertex->fill(pvCov);
236  for(int i = 0; i<LorentzVectorParticle::NVertex; i++){
237  for(int j = 0; j<LorentzVectorParticle::NVertex; j++){
238  pvcov(i,j)=pvCov(i,j);
239  }
240  }
241  for(unsigned int i=0; i<PFTau3ProngSummary::nsolutions;i++){
242  TauA1NuConstrainedFitter TauA1NU(i,a1,pv,pvcov);
243  bool isFitOK=TauA1NU.fit();
244  if(isFitOK){
245  LorentzVectorParticle theTau=TauA1NU.getMother();
246  std::vector<LorentzVectorParticle> daughter=TauA1NU.getRefitDaughters();
247  std::vector<TLorentzVector> daughter_p4;
248  std::vector<int> daughter_charge,daughter_PDGID;
249  for(unsigned int d=0;d<daughter.size();d++){
250  daughter_p4.push_back(daughter[d].p4());
251  daughter_charge.push_back((int)daughter[d].charge());
252  daughter_PDGID.push_back(daughter[d].pdgId());
253  }
254  PFTau3PS.AddSolution(i,theTau.p4(),daughter_p4,daughter_charge,daughter_PDGID,(isFitOK),0.0,-999);
255  }
256  }
257  }
258  }
259  reco::PFTau3ProngSummaryRef PFTau3PSRef=reco::PFTau3ProngSummaryRef(PFTau3RefProd_out,PFTau3PSCollection_out->size());
260  PFTau3PSCollection_out->push_back(PFTau3PS);
261  AVPFTau3PS->setValue(iPFTau,PFTau3PSRef);
262  }
263  }
264  iEvent.put(std::move(PFTau3PSCollection_out),"PFTau3ProngSummary");
265  iEvent.put(std::move(AVPFTau3PS));
266 }
267 
std::vector< reco::PFTau3ProngSummary > PFTau3ProngSummaryCollection
collection of PFTau3ProngSummary objects
KinematicParticleVertexFitter kpvFitter_
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
PFTau3ProngReco(const edm::ParameterSet &iConfig)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::unique_ptr< StringCutObjectSelector< reco::PFTau > > cut_
reco::TransientTrack build(const reco::Track *p) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::InputTag PFTauTIPTag_
ErrorD< N >::type type
Definition: Error.h:33
key_type key() const
Accessor for product key.
Definition: Ref.h:263
#define nullptr
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
Definition: PFTauFwd.h:15
uint16_t size_type
ROOT::Math::SMatrix< double, 7, 7, ROOT::Math::MatRepSym< double, 7 > > AlgebraicSymMatrix77
Definition: Matrices.h:9
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:41
std::vector< LorentzVectorParticle > getRefitDaughters()
int iEvent
Definition: GenABIO.cc:230
double p4[4]
Definition: TauolaWrapper.h:92
T z() const
Definition: PV3DBase.h:64
const MagneticField * field() const
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Handle< reco::PFTauDiscriminator > handle_
static double pi_mass()
Definition: PDGInfo.h:13
DiscCutPairVec discriminators_
bool isValid() const
Definition: HandleBase.h:74
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
RefProd< PROD > getRefBeforePut()
Definition: Event.h:167
edm::Ref< PFTau3ProngSummaryCollection > PFTau3ProngSummaryRef
presistent reference to a PFTau3ProngSummary
~PFTau3ProngReco() override
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)
std::vector< DiscCutPair * > DiscCutPairVec
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:68
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
void produce(edm::Event &, const edm::EventSetup &) override
edm::InputTag PFTauTag_
def move(src, dest)
Definition: eostools.py:511