source: HiSusy/trunk/Delphes-3.0.0/modules/TauTagging.cc @ 1

Last change on this file since 1 was 1, checked in by zerwas, 11 years ago

first import of structure, PYTHIA8 and DELPHES

File size: 5.8 KB
Line 
1
2/** \class TauTagging
3 *
4 *  Determines origin of jet,
5 *  applies b-tagging efficiency (miss identification rate) formulas
6 *  and sets b-tagging flags
7 *
8 *  $Date: 2012-12-06 16:56:12 +0100 (Thu, 06 Dec 2012) $
9 *  $Revision: 871 $
10 *
11 *
12 *  \author P. Demin - UCL, Louvain-la-Neuve
13 *
14 */
15
16#include "modules/TauTagging.h"
17
18#include "classes/DelphesClasses.h"
19#include "classes/DelphesFactory.h"
20#include "classes/DelphesFormula.h"
21
22#include "ExRootAnalysis/ExRootResult.h"
23#include "ExRootAnalysis/ExRootFilter.h"
24#include "ExRootAnalysis/ExRootClassifier.h"
25
26#include "TMath.h"
27#include "TString.h"
28#include "TFormula.h"
29#include "TRandom3.h"
30#include "TObjArray.h"
31#include "TDatabasePDG.h"
32#include "TLorentzVector.h"
33
34#include <algorithm> 
35#include <stdexcept>
36#include <iostream>
37#include <sstream>
38
39using namespace std;
40
41//------------------------------------------------------------------------------
42
43class TauTaggingPartonClassifier : public ExRootClassifier
44{
45public:
46
47  TauTaggingPartonClassifier(const TObjArray *array);
48
49  Int_t GetCategory(TObject *object);
50 
51  Double_t fEtaMax, fPTMin;
52 
53  const TObjArray *fParticleInputArray;
54};
55
56//------------------------------------------------------------------------------
57TauTaggingPartonClassifier::TauTaggingPartonClassifier(const TObjArray *array) :
58  fParticleInputArray(array)
59{
60}
61
62//------------------------------------------------------------------------------
63
64Int_t TauTaggingPartonClassifier::GetCategory(TObject *object)
65{
66  Candidate *tau = static_cast<Candidate *>(object);
67  Candidate *daughter = 0;
68  const TLorentzVector &momentum = tau->Momentum;
69  Int_t pdgCode, i;
70
71  pdgCode = TMath::Abs(tau->PID);
72  if(pdgCode != 15) return -1;
73
74  if(momentum.Pt() <= fPTMin || TMath::Abs(momentum.Eta()) > fEtaMax) return -1;
75 
76  for(i = tau->D1; i <= tau->D2; ++i)
77  {
78    daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
79    pdgCode = TMath::Abs(tau->PID);
80    if(pdgCode != 11 || pdgCode != 13 || pdgCode != 15) return -1;
81  }
82
83  return 0;
84}
85
86//------------------------------------------------------------------------------
87
88TauTagging::TauTagging() :
89  fClassifier(0), fFilter(0),
90  fItPartonInputArray(0), fItJetInputArray(0)
91{
92}
93
94//------------------------------------------------------------------------------
95
96TauTagging::~TauTagging()
97{
98}
99
100//------------------------------------------------------------------------------
101
102void TauTagging::Init()
103{
104  map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
105  ExRootConfParam param;
106  DelphesFormula *formula;
107  Int_t i, size;
108
109  fDeltaR = GetDouble("DeltaR", 0.5);
110
111  // read efficiency formulas
112  param = GetParam("EfficiencyFormula");
113  size = param.GetSize();
114 
115  fEfficiencyMap.clear();
116  for(i = 0; i < size/2; ++i)
117  {
118    formula = new DelphesFormula;
119    formula->Compile(param[i*2 + 1].GetString());
120
121    fEfficiencyMap[param[i*2].GetInt()] = formula;
122  }
123
124  // set default efficiency formula
125  itEfficiencyMap = fEfficiencyMap.find(0);
126  if(itEfficiencyMap == fEfficiencyMap.end())
127  {
128    formula = new DelphesFormula;
129    formula->Compile("0.0");
130
131    fEfficiencyMap[0] = formula;
132  }
133
134  // import input array(s)
135
136  fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/particles"));
137
138  fClassifier = new TauTaggingPartonClassifier(fParticleInputArray);
139  fClassifier->fPTMin = GetDouble("TauPTMin", 1.0);
140  fClassifier->fEtaMax = GetDouble("TauEtaMax", 2.5);
141
142  fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
143  fItPartonInputArray = fPartonInputArray->MakeIterator();
144
145  fFilter = new ExRootFilter(fPartonInputArray);
146 
147  fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
148  fItJetInputArray = fJetInputArray->MakeIterator();
149
150  // create output array
151
152  fOutputArray = ExportArray(GetString("OutputArray", "taus"));
153}
154
155//------------------------------------------------------------------------------
156
157void TauTagging::Finish()
158{
159  map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
160  DelphesFormula *formula;
161
162  if(fFilter) delete fFilter;
163  if(fClassifier) delete fClassifier;
164  if(fItJetInputArray) delete fItJetInputArray;
165  if(fItPartonInputArray) delete fItPartonInputArray;
166
167  for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
168  {
169    formula = itEfficiencyMap->second;
170    if(formula) delete formula;
171  }
172}
173
174//------------------------------------------------------------------------------
175
176void TauTagging::Process()
177{
178  Candidate *jet, *tau;
179  Double_t pt, eta, phi;
180  TObjArray *tauArray;
181  map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
182  DelphesFormula *formula;
183  Int_t pdgCode, charge;
184
185  // select taus
186  fFilter->Reset();
187  tauArray = fFilter->GetSubArray(fClassifier, 0);
188 
189  if(tauArray == 0) return;
190
191  TIter itTauArray(tauArray);
192 
193  // loop over all input jets
194  fItJetInputArray->Reset();
195  while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
196  {
197    const TLorentzVector &jetMomentum = jet->Momentum;
198    pdgCode = 0;
199    charge = gRandom->Uniform() > 0.5 ? 1 : -1;
200    eta = jetMomentum.Eta();
201    phi = jetMomentum.Phi();
202    pt = jetMomentum.Pt();
203
204    // loop over all input taus
205    itTauArray.Reset();
206    while((tau = static_cast<Candidate *>(itTauArray.Next())))
207    {
208      if(jetMomentum.DeltaR(tau->Momentum) <= fDeltaR)
209      {
210        pdgCode = 15;
211        charge = tau->Charge;
212      }
213    }
214    // find an efficency formula
215    itEfficiencyMap = fEfficiencyMap.find(pdgCode);
216    if(itEfficiencyMap == fEfficiencyMap.end())
217    {
218      itEfficiencyMap = fEfficiencyMap.find(0);
219    }
220    formula = itEfficiencyMap->second;
221
222    // apply an efficency formula
223   
224    if(gRandom->Uniform() > formula->Eval(pt, eta)) continue;
225   
226    jet->Charge = charge;
227
228    fOutputArray->Add(jet);
229  }
230}
231
232//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.