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 PileUpMerger
|
---|
20 | *
|
---|
21 | * Merges particles from pile-up sample into event
|
---|
22 | *
|
---|
23 | * \author M. Selvaggi - UCL, Louvain-la-Neuve
|
---|
24 | *
|
---|
25 | */
|
---|
26 |
|
---|
27 | #include "modules/PileUpMerger.h"
|
---|
28 |
|
---|
29 | #include "classes/DelphesClasses.h"
|
---|
30 | #include "classes/DelphesFactory.h"
|
---|
31 | #include "classes/DelphesTF2.h"
|
---|
32 | #include "classes/DelphesPileUpReader.h"
|
---|
33 |
|
---|
34 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
35 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
36 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
37 |
|
---|
38 | #include "TMath.h"
|
---|
39 | #include "TString.h"
|
---|
40 | #include "TFormula.h"
|
---|
41 | #include "TRandom3.h"
|
---|
42 | #include "TObjArray.h"
|
---|
43 | #include "TDatabasePDG.h"
|
---|
44 | #include "TLorentzVector.h"
|
---|
45 |
|
---|
46 | #include <algorithm>
|
---|
47 | #include <stdexcept>
|
---|
48 | #include <iostream>
|
---|
49 | #include <sstream>
|
---|
50 |
|
---|
51 | using namespace std;
|
---|
52 |
|
---|
53 | //------------------------------------------------------------------------------
|
---|
54 |
|
---|
55 | PileUpMerger::PileUpMerger() :
|
---|
56 | fFunction(0), fReader(0), fItInputArray(0)
|
---|
57 | {
|
---|
58 | fFunction = new DelphesTF2;
|
---|
59 | }
|
---|
60 |
|
---|
61 |
|
---|
62 | //------------------------------------------------------------------------------
|
---|
63 |
|
---|
64 | PileUpMerger::~PileUpMerger()
|
---|
65 | {
|
---|
66 | delete fFunction;
|
---|
67 | }
|
---|
68 |
|
---|
69 | //------------------------------------------------------------------------------
|
---|
70 |
|
---|
71 | void PileUpMerger::Init()
|
---|
72 | {
|
---|
73 | const char *fileName;
|
---|
74 |
|
---|
75 | fPileUpDistribution = GetInt("PileUpDistribution", 0);
|
---|
76 |
|
---|
77 | fMeanPileUp = GetDouble("MeanPileUp", 10);
|
---|
78 |
|
---|
79 | fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
|
---|
80 | fTVertexSpread = GetDouble("TVertexSpread", 1.5E-09);
|
---|
81 |
|
---|
82 | fInputBeamSpotX = GetDouble("InputBeamSpotX", 0.0);
|
---|
83 | fInputBeamSpotY = GetDouble("InputBeamSpotY", 0.0);
|
---|
84 | fOutputBeamSpotX = GetDouble("OutputBeamSpotX", 0.0);
|
---|
85 | fOutputBeamSpotY = GetDouble("OutputBeamSpotY", 0.0);
|
---|
86 |
|
---|
87 | // read vertex smearing formula
|
---|
88 |
|
---|
89 | fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
|
---|
90 | fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
|
---|
91 |
|
---|
92 | fileName = GetString("PileUpFile", "MinBias.pileup");
|
---|
93 | fReader = new DelphesPileUpReader(fileName);
|
---|
94 |
|
---|
95 | // import input array
|
---|
96 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
97 | fItInputArray = fInputArray->MakeIterator();
|
---|
98 |
|
---|
99 | // create output arrays
|
---|
100 | fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
|
---|
101 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
102 | }
|
---|
103 |
|
---|
104 | //------------------------------------------------------------------------------
|
---|
105 |
|
---|
106 | void PileUpMerger::Finish()
|
---|
107 | {
|
---|
108 | if(fReader) delete fReader;
|
---|
109 | }
|
---|
110 |
|
---|
111 | //------------------------------------------------------------------------------
|
---|
112 |
|
---|
113 | void PileUpMerger::Process()
|
---|
114 | {
|
---|
115 | TDatabasePDG *pdg = TDatabasePDG::Instance();
|
---|
116 | TParticlePDG *pdgParticle;
|
---|
117 | Int_t pid, nch, nvtx = -1;
|
---|
118 | Float_t x, y, z, t, vx, vy;
|
---|
119 | Float_t px, py, pz, e, pt;
|
---|
120 | Double_t dz, dphi, dt, sumpt2, dz0, dt0;
|
---|
121 | Int_t numberOfEvents, event, numberOfParticles;
|
---|
122 | Long64_t allEntries, entry;
|
---|
123 | Candidate *candidate, *vertex;
|
---|
124 | DelphesFactory *factory;
|
---|
125 |
|
---|
126 | const Double_t c_light = 2.99792458E8;
|
---|
127 |
|
---|
128 | fItInputArray->Reset();
|
---|
129 |
|
---|
130 | // --- Deal with primary vertex first ------
|
---|
131 |
|
---|
132 | fFunction->GetRandom2(dz, dt);
|
---|
133 |
|
---|
134 | dz0 = -1.0e6;
|
---|
135 | dt0 = -1.0e6;
|
---|
136 |
|
---|
137 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
138 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
139 |
|
---|
140 | //cout<<dz<<","<<dt<<endl;
|
---|
141 |
|
---|
142 | vx = 0.0;
|
---|
143 | vy = 0.0;
|
---|
144 |
|
---|
145 | numberOfParticles = fInputArray->GetEntriesFast();
|
---|
146 | nch = 0;
|
---|
147 | sumpt2 = 0.0;
|
---|
148 |
|
---|
149 | factory = GetFactory();
|
---|
150 | vertex = factory->NewCandidate();
|
---|
151 |
|
---|
152 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
153 | {
|
---|
154 | vx += candidate->Position.X();
|
---|
155 | vy += candidate->Position.Y();
|
---|
156 | z = candidate->Position.Z();
|
---|
157 | t = candidate->Position.T();
|
---|
158 | pt = candidate->Momentum.Pt();
|
---|
159 |
|
---|
160 | // take postion and time from first stable particle
|
---|
161 | if (dz0 < -999999.0)
|
---|
162 | dz0 = z;
|
---|
163 | if (dt0 < -999999.0)
|
---|
164 | dt0 = t;
|
---|
165 |
|
---|
166 | // cancel any possible offset in position and time the input file
|
---|
167 | candidate->Position.SetZ(z - dz0 + dz);
|
---|
168 | candidate->Position.SetT(t - dt0 + dt);
|
---|
169 |
|
---|
170 | fParticleOutputArray->Add(candidate);
|
---|
171 |
|
---|
172 | if(TMath::Abs(candidate->Charge) > 1.0E-9)
|
---|
173 | {
|
---|
174 | nch++;
|
---|
175 | sumpt2 += pt*pt;
|
---|
176 | vertex->AddCandidate(candidate);
|
---|
177 | }
|
---|
178 | }
|
---|
179 |
|
---|
180 | if(numberOfParticles > 0)
|
---|
181 | {
|
---|
182 | vx /= sumpt2;
|
---|
183 | vy /= sumpt2;
|
---|
184 | }
|
---|
185 |
|
---|
186 | nvtx++;
|
---|
187 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
188 | vertex->ClusterIndex = nvtx;
|
---|
189 | vertex->ClusterNDF = nch;
|
---|
190 | vertex->SumPT2 = sumpt2;
|
---|
191 | vertex->GenSumPT2 = sumpt2;
|
---|
192 |
|
---|
193 | fVertexOutputArray->Add(vertex);
|
---|
194 |
|
---|
195 | // --- Then with pile-up vertices ------
|
---|
196 |
|
---|
197 | switch(fPileUpDistribution)
|
---|
198 | {
|
---|
199 | case 0:
|
---|
200 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
201 | break;
|
---|
202 | case 1:
|
---|
203 | numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
|
---|
204 | break;
|
---|
205 | case 2:
|
---|
206 | numberOfEvents = fMeanPileUp;
|
---|
207 | break;
|
---|
208 | default:
|
---|
209 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
210 | break;
|
---|
211 | }
|
---|
212 |
|
---|
213 | allEntries = fReader->GetEntries();
|
---|
214 |
|
---|
215 | for(event = 0; event < numberOfEvents; ++event)
|
---|
216 | {
|
---|
217 | do
|
---|
218 | {
|
---|
219 | entry = TMath::Nint(gRandom->Rndm()*allEntries);
|
---|
220 | }
|
---|
221 | while(entry >= allEntries);
|
---|
222 |
|
---|
223 | fReader->ReadEntry(entry);
|
---|
224 |
|
---|
225 | // --- Pile-up vertex smearing
|
---|
226 |
|
---|
227 | fFunction->GetRandom2(dz, dt);
|
---|
228 |
|
---|
229 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
230 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
231 |
|
---|
232 | dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
|
---|
233 |
|
---|
234 | vx = 0.0;
|
---|
235 | vy = 0.0;
|
---|
236 |
|
---|
237 | numberOfParticles = 0;
|
---|
238 | sumpt2 = 0.0;
|
---|
239 |
|
---|
240 | vertex = factory->NewCandidate();
|
---|
241 |
|
---|
242 | while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
|
---|
243 | {
|
---|
244 | candidate = factory->NewCandidate();
|
---|
245 |
|
---|
246 | candidate->PID = pid;
|
---|
247 |
|
---|
248 | candidate->Status = 1;
|
---|
249 |
|
---|
250 | pdgParticle = pdg->GetParticle(pid);
|
---|
251 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
|
---|
252 | candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
|
---|
253 |
|
---|
254 | candidate->IsPU = 1;
|
---|
255 |
|
---|
256 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
257 | candidate->Momentum.RotateZ(dphi);
|
---|
258 |
|
---|
259 | x -= fInputBeamSpotX;
|
---|
260 | y -= fInputBeamSpotY;
|
---|
261 | candidate->Position.SetXYZT(x, y, z + dz, t + dt);
|
---|
262 | candidate->Position.RotateZ(dphi);
|
---|
263 | candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
|
---|
264 |
|
---|
265 | vx += candidate->Position.X();
|
---|
266 | vy += candidate->Position.Y();
|
---|
267 |
|
---|
268 | ++numberOfParticles;
|
---|
269 | if(TMath::Abs(candidate->Charge) > 1.0E-9)
|
---|
270 | {
|
---|
271 | nch++;
|
---|
272 | sumpt2 += pt*pt;
|
---|
273 | vertex->AddCandidate(candidate);
|
---|
274 | }
|
---|
275 |
|
---|
276 | fParticleOutputArray->Add(candidate);
|
---|
277 | }
|
---|
278 |
|
---|
279 | if(numberOfParticles > 0)
|
---|
280 | {
|
---|
281 | vx /= numberOfParticles;
|
---|
282 | vy /= numberOfParticles;
|
---|
283 | }
|
---|
284 |
|
---|
285 | nvtx++;
|
---|
286 |
|
---|
287 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
288 |
|
---|
289 | vertex->ClusterIndex = nvtx;
|
---|
290 | vertex->ClusterNDF = nch;
|
---|
291 | vertex->SumPT2 = sumpt2;
|
---|
292 | vertex->GenSumPT2 = sumpt2;
|
---|
293 |
|
---|
294 | vertex->IsPU = 1;
|
---|
295 |
|
---|
296 | fVertexOutputArray->Add(vertex);
|
---|
297 |
|
---|
298 | }
|
---|
299 | }
|
---|
300 |
|
---|
301 | //------------------------------------------------------------------------------
|
---|