1 | /*
|
---|
2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
4 | *
|
---|
5 | * This program is free software: you can redistribute it and/or modify
|
---|
6 | * it under the terms of the GNU General Public License as published by
|
---|
7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
8 | * (at your option) any later version.
|
---|
9 | *
|
---|
10 | * This program is distributed in the hope that it will be useful,
|
---|
11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
13 | * GNU General Public License for more details.
|
---|
14 | *
|
---|
15 | * You should have received a copy of the GNU General Public License
|
---|
16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
17 | */
|
---|
18 |
|
---|
19 | /** \class TrackPileUpSubtractor
|
---|
20 | *
|
---|
21 | * Subtract pile-up contribution from tracks.
|
---|
22 | *
|
---|
23 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
24 | *
|
---|
25 | */
|
---|
26 |
|
---|
27 | #include "modules/TrackTimingPileUpSubtractor.h"
|
---|
28 |
|
---|
29 | #include "classes/DelphesClasses.h"
|
---|
30 | #include "classes/DelphesFactory.h"
|
---|
31 | #include "classes/DelphesFormula.h"
|
---|
32 |
|
---|
33 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
34 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
35 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
36 |
|
---|
37 | #include "TDatabasePDG.h"
|
---|
38 | #include "TFormula.h"
|
---|
39 | #include "TLorentzVector.h"
|
---|
40 | #include "TMath.h"
|
---|
41 | #include "TObjArray.h"
|
---|
42 | #include "TRandom3.h"
|
---|
43 | #include "TString.h"
|
---|
44 |
|
---|
45 | #include <algorithm>
|
---|
46 | #include <iostream>
|
---|
47 | #include <sstream>
|
---|
48 | #include <stdexcept>
|
---|
49 |
|
---|
50 | using namespace std;
|
---|
51 |
|
---|
52 | //------------------------------------------------------------------------------
|
---|
53 |
|
---|
54 | TrackTimingPileUpSubtractor::TrackTimingPileUpSubtractor() :
|
---|
55 | fFormula(0)
|
---|
56 | {
|
---|
57 | fFormula = new DelphesFormula;
|
---|
58 | }
|
---|
59 |
|
---|
60 | //------------------------------------------------------------------------------
|
---|
61 |
|
---|
62 | TrackTimingPileUpSubtractor::~TrackTimingPileUpSubtractor()
|
---|
63 | {
|
---|
64 | if(fFormula) delete fFormula;
|
---|
65 | }
|
---|
66 |
|
---|
67 | //------------------------------------------------------------------------------
|
---|
68 |
|
---|
69 | void TrackTimingPileUpSubtractor::Init()
|
---|
70 | {
|
---|
71 | // import input array
|
---|
72 |
|
---|
73 | fVertexInputArray = ImportArray(GetString("VertexInputArray", "PileUpMerger/vertices"));
|
---|
74 | fItVertexInputArray = fVertexInputArray->MakeIterator();
|
---|
75 |
|
---|
76 | // read resolution formula in m
|
---|
77 | fChargedMinSignificance = GetDouble("ChargedMinSignificance", 3);
|
---|
78 | fNeutralMinSignificance = GetDouble("NeutralMinSignificance", 3);
|
---|
79 |
|
---|
80 | fPTMin = GetDouble("PTMin", 0.);
|
---|
81 |
|
---|
82 | // import arrays with output from other modules
|
---|
83 |
|
---|
84 | ExRootConfParam param = GetParam("InputArray");
|
---|
85 | Long_t i, size;
|
---|
86 | const TObjArray *array;
|
---|
87 | TIterator *iterator;
|
---|
88 |
|
---|
89 | size = param.GetSize();
|
---|
90 | for(i = 0; i < size / 2; ++i)
|
---|
91 | {
|
---|
92 | array = ImportArray(param[i * 2].GetString());
|
---|
93 | iterator = array->MakeIterator();
|
---|
94 |
|
---|
95 | fInputMap[iterator] = ExportArray(param[i * 2 + 1].GetString());
|
---|
96 | }
|
---|
97 | }
|
---|
98 |
|
---|
99 | //------------------------------------------------------------------------------
|
---|
100 |
|
---|
101 | void TrackTimingPileUpSubtractor::Finish()
|
---|
102 | {
|
---|
103 | map<TIterator *, TObjArray *>::iterator itInputMap;
|
---|
104 | TIterator *iterator;
|
---|
105 |
|
---|
106 | for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
|
---|
107 | {
|
---|
108 | iterator = itInputMap->first;
|
---|
109 |
|
---|
110 | if(iterator) delete iterator;
|
---|
111 | }
|
---|
112 |
|
---|
113 | if(fItVertexInputArray) delete fItVertexInputArray;
|
---|
114 | }
|
---|
115 |
|
---|
116 | //------------------------------------------------------------------------------
|
---|
117 |
|
---|
118 | void TrackTimingPileUpSubtractor::Process()
|
---|
119 | {
|
---|
120 | Candidate *candidate, *particle;
|
---|
121 | map<TIterator *, TObjArray *>::iterator itInputMap;
|
---|
122 | TIterator *iterator;
|
---|
123 | TObjArray *array;
|
---|
124 | Double_t z, zvtx = 0;
|
---|
125 | Double_t z_err, zvtx_err = 0;
|
---|
126 | Double_t t, tvtx = 0;
|
---|
127 | Double_t t_err, tvtx_err = 0;
|
---|
128 | Double_t sumPTSquare = 0;
|
---|
129 | Double_t tempPTSquare = 0;
|
---|
130 | Double_t pt, eta, phi, e;
|
---|
131 | Double_t distanceCharged, distanceNeutral = 0;
|
---|
132 |
|
---|
133 | // find z position of primary vertex
|
---|
134 |
|
---|
135 | fItVertexInputArray->Reset();
|
---|
136 | while((candidate = static_cast<Candidate *>(fItVertexInputArray->Next())))
|
---|
137 | {
|
---|
138 | tempPTSquare = candidate->SumPT2;
|
---|
139 | if(tempPTSquare > sumPTSquare)
|
---|
140 | {
|
---|
141 | sumPTSquare = tempPTSquare;
|
---|
142 | zvtx = candidate->Position.Z();
|
---|
143 | zvtx_err = candidate->PositionError.Z();
|
---|
144 | tvtx = candidate->Position.T();
|
---|
145 | tvtx_err = candidate->PositionError.T();
|
---|
146 | }
|
---|
147 | }
|
---|
148 |
|
---|
149 | // loop over all input arrays
|
---|
150 | for(itInputMap = fInputMap.begin(); itInputMap != fInputMap.end(); ++itInputMap)
|
---|
151 | {
|
---|
152 | iterator = itInputMap->first;
|
---|
153 | array = itInputMap->second;
|
---|
154 |
|
---|
155 | // loop over all candidates
|
---|
156 | iterator->Reset();
|
---|
157 | while((candidate = static_cast<Candidate *>(iterator->Next())))
|
---|
158 | {
|
---|
159 | particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
160 | const TLorentzVector &candidateMomentum = particle->Momentum;
|
---|
161 |
|
---|
162 | eta = candidateMomentum.Eta();
|
---|
163 | pt = candidateMomentum.Pt();
|
---|
164 | phi = candidateMomentum.Phi();
|
---|
165 | e = candidateMomentum.E();
|
---|
166 |
|
---|
167 | z = particle->Position.Z();
|
---|
168 | z_err = particle->PositionError.Z();
|
---|
169 | t = particle->InitialPosition.T();
|
---|
170 | t_err = particle->PositionError.T();
|
---|
171 |
|
---|
172 | distanceCharged = pow((zvtx - z),2)/pow((zvtx_err - z_err),2) + pow((tvtx - t),2)/pow((tvtx_err - t_err),2);
|
---|
173 | distanceNeutral = pow((tvtx - t),2)/pow((tvtx_err - t_err),2);
|
---|
174 |
|
---|
175 | if(candidate->Charge != 0 && distanceCharged < fChargedMinSignificance)
|
---|
176 | {
|
---|
177 | candidate->IsRecoPU = 1;
|
---|
178 | }
|
---|
179 | else if(candidate->Charge == 0 && distanceNeutral < fNeutralMinSignificance)
|
---|
180 | {
|
---|
181 | candidate->IsRecoPU = 1;
|
---|
182 | }
|
---|
183 | else
|
---|
184 | {
|
---|
185 | candidate->IsRecoPU = 0;
|
---|
186 | if(candidate->Momentum.Pt() > fPTMin) array->Add(candidate);
|
---|
187 | }
|
---|
188 | }
|
---|
189 | }
|
---|
190 | }
|
---|