source: HiSusy/trunk/Delphes-3.0.0/modules/Isolation.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: 4.2 KB
Line 
1
2/** \class Isolation
3 *
4 *  Sums transverse momenta of isolation objects (tracks, calorimeter towers, etc)
5 *  within a DeltaR cone around a candidate and calculates fraction of this sum
6 *  to the candidate's transverse momentum. outputs candidates that have
7 *  the transverse momenta fraction within (PTRatioMin, PTRatioMax].
8 *
9 *  $Date: 2013-01-11 11:08:39 +0100 (Fri, 11 Jan 2013) $
10 *  $Revision: 878 $
11 *
12 *
13 *  \author P. Demin - UCL, Louvain-la-Neuve
14 *
15 */
16
17#include "modules/Isolation.h"
18
19#include "classes/DelphesClasses.h"
20#include "classes/DelphesFactory.h"
21#include "classes/DelphesFormula.h"
22
23#include "ExRootAnalysis/ExRootResult.h"
24#include "ExRootAnalysis/ExRootFilter.h"
25#include "ExRootAnalysis/ExRootClassifier.h"
26
27#include "TMath.h"
28#include "TString.h"
29#include "TFormula.h"
30#include "TRandom3.h"
31#include "TObjArray.h"
32#include "TDatabasePDG.h"
33#include "TLorentzVector.h"
34
35#include <algorithm> 
36#include <stdexcept>
37#include <iostream>
38#include <sstream>
39
40using namespace std;
41
42//------------------------------------------------------------------------------
43
44class IsolationClassifier : public ExRootClassifier
45{
46public:
47
48  IsolationClassifier() {}
49
50  Int_t GetCategory(TObject *object);
51 
52  Double_t fPTMin;
53};
54
55//------------------------------------------------------------------------------
56
57Int_t IsolationClassifier::GetCategory(TObject *object)
58{
59  Candidate *track = static_cast<Candidate*>(object);
60  const TLorentzVector &momentum = track->Momentum;
61
62  if(momentum.Pt() < fPTMin) return -1;
63
64  return 0;
65}
66
67//------------------------------------------------------------------------------
68
69Isolation::Isolation() :
70  fClassifier(0), fFilter(0),
71  fItIsolationInputArray(0), fItCandidateInputArray(0)
72{
73  fClassifier = new IsolationClassifier;
74}
75
76//------------------------------------------------------------------------------
77
78Isolation::~Isolation()
79{
80}
81
82//------------------------------------------------------------------------------
83
84void Isolation::Init()
85{
86  fDeltaRMin = GetDouble("DeltaRMin", 0.0);
87  fDeltaRMax = GetDouble("DeltaRMax", 0.5);
88
89  fPTRatioMin = GetDouble("PTRatioMin", 0.0);
90  fPTRatioMax = GetDouble("PTRatioMax", 0.1);
91
92  fClassifier->fPTMin = GetDouble("PTMin", 0.5);
93
94  // import input array(s)
95
96  fIsolationInputArray = ImportArray(GetString("IsolationInputArray", "Delphes/partons"));
97  fItIsolationInputArray = fIsolationInputArray->MakeIterator();
98
99  fFilter = new ExRootFilter(fIsolationInputArray);
100 
101  fCandidateInputArray = ImportArray(GetString("CandidateInputArray", "Calorimeter/electrons"));
102  fItCandidateInputArray = fCandidateInputArray->MakeIterator();
103 
104  // create output array
105
106  fOutputArray = ExportArray(GetString("OutputArray", "electrons"));
107}
108
109//------------------------------------------------------------------------------
110
111void Isolation::Finish()
112{
113  if(fFilter) delete fFilter;
114  if(fItCandidateInputArray) delete fItCandidateInputArray;
115  if(fItIsolationInputArray) delete fItIsolationInputArray;
116}
117
118//------------------------------------------------------------------------------
119
120void Isolation::Process()
121{
122  Candidate *candidate, *isolation;
123  TObjArray *isolationArray;
124  Double_t sumPT, ratio;
125
126  // select isolation objects
127  fFilter->Reset();
128  isolationArray = fFilter->GetSubArray(fClassifier, 0);
129 
130  if(isolationArray == 0) return;
131
132  TIter itIsolationArray(isolationArray);
133 
134  // loop over all input jets
135  fItCandidateInputArray->Reset();
136  while((candidate = static_cast<Candidate*>(fItCandidateInputArray->Next())))
137  {
138    const TLorentzVector &candidateMomentum = candidate->Momentum;
139
140    // loop over all input tracks
141    sumPT = 0.0;
142    itIsolationArray.Reset();
143    while((isolation = static_cast<Candidate*>(itIsolationArray.Next())))
144    {
145      const TLorentzVector &isolationMomentum = isolation->Momentum;
146      if(candidateMomentum.DeltaR(isolationMomentum) > fDeltaRMin &&
147         candidateMomentum.DeltaR(isolationMomentum) <= fDeltaRMax)
148      {
149        sumPT += isolationMomentum.Pt();
150      }
151    }
152
153    ratio = sumPT/candidateMomentum.Pt();
154    if(ratio <= fPTRatioMin || ratio > fPTRatioMax) continue;
155
156    fOutputArray->Add(candidate);
157
158  }
159}
160
161//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.