Fork me on GitHub

source: git/modules/PileUpMerger.cc@ 952bbbc

Last change on this file since 952bbbc was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • Property mode set to 100644
File size: 7.6 KB
Line 
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/DelphesPileUpReader.h"
32#include "classes/DelphesTF2.h"
33
34#include "ExRootAnalysis/ExRootClassifier.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootResult.h"
37
38#include "TDatabasePDG.h"
39#include "TFormula.h"
40#include "TLorentzVector.h"
41#include "TMath.h"
42#include "TObjArray.h"
43#include "TRandom3.h"
44#include "TString.h"
45
46#include <algorithm>
47#include <iostream>
48#include <sstream>
49#include <stdexcept>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55PileUpMerger::PileUpMerger() :
56 fFunction(0), fReader(0), fItInputArray(0)
57{
58 fFunction = new DelphesTF2;
59}
60
61//------------------------------------------------------------------------------
62
63PileUpMerger::~PileUpMerger()
64{
65 delete fFunction;
66}
67
68//------------------------------------------------------------------------------
69
70void PileUpMerger::Init()
71{
72 const char *fileName;
73
74 fPileUpDistribution = GetInt("PileUpDistribution", 0);
75
76 fMeanPileUp = GetDouble("MeanPileUp", 10);
77
78 fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
79 fTVertexSpread = GetDouble("TVertexSpread", 1.5E-09);
80
81 fInputBeamSpotX = GetDouble("InputBeamSpotX", 0.0);
82 fInputBeamSpotY = GetDouble("InputBeamSpotY", 0.0);
83 fOutputBeamSpotX = GetDouble("OutputBeamSpotX", 0.0);
84 fOutputBeamSpotY = GetDouble("OutputBeamSpotY", 0.0);
85
86 // read vertex smearing formula
87
88 fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
89 fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
90
91 fileName = GetString("PileUpFile", "MinBias.pileup");
92 fReader = new DelphesPileUpReader(fileName);
93
94 // import input array
95 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
96 fItInputArray = fInputArray->MakeIterator();
97
98 // create output arrays
99 fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
100 fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
101}
102
103//------------------------------------------------------------------------------
104
105void PileUpMerger::Finish()
106{
107 if(fReader) delete fReader;
108}
109
110//------------------------------------------------------------------------------
111
112void PileUpMerger::Process()
113{
114 TDatabasePDG *pdg = TDatabasePDG::Instance();
115 TParticlePDG *pdgParticle;
116 Int_t pid, nch, nvtx = -1;
117 Float_t x, y, z, t, vx, vy;
118 Float_t px, py, pz, e, pt;
119 Double_t dz, dphi, dt, sumpt2, dz0, dt0;
120 Int_t numberOfEvents, event, numberOfParticles;
121 Long64_t allEntries, entry;
122 Candidate *candidate, *vertex;
123 DelphesFactory *factory;
124
125 const Double_t c_light = 2.99792458E8;
126
127 fItInputArray->Reset();
128
129 // --- Deal with primary vertex first ------
130
131 fFunction->GetRandom2(dz, dt);
132
133 dz0 = -1.0e6;
134 dt0 = -1.0e6;
135
136 dt *= c_light * 1.0E3; // necessary in order to make t in mm/c
137 dz *= 1.0E3; // necessary in order to make z in mm
138
139 //cout<<dz<<","<<dt<<endl;
140
141 vx = 0.0;
142 vy = 0.0;
143
144 numberOfParticles = fInputArray->GetEntriesFast();
145 nch = 0;
146 sumpt2 = 0.0;
147
148 factory = GetFactory();
149 vertex = factory->NewCandidate();
150
151 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
152 {
153 vx += candidate->Position.X();
154 vy += candidate->Position.Y();
155 z = candidate->Position.Z();
156 t = candidate->Position.T();
157 pt = candidate->Momentum.Pt();
158
159 // take postion and time from first stable particle
160 if(dz0 < -999999.0)
161 dz0 = z;
162 if(dt0 < -999999.0)
163 dt0 = t;
164
165 // cancel any possible offset in position and time the input file
166 candidate->Position.SetZ(z - dz0 + dz);
167 candidate->Position.SetT(t - dt0 + dt);
168
169 candidate->IsPU = 0;
170
171 fParticleOutputArray->Add(candidate);
172
173 if(TMath::Abs(candidate->Charge) > 1.0E-9)
174 {
175 nch++;
176 sumpt2 += pt * pt;
177 vertex->AddCandidate(candidate);
178 }
179 }
180
181 if(numberOfParticles > 0)
182 {
183 vx /= sumpt2;
184 vy /= sumpt2;
185 }
186
187 nvtx++;
188 vertex->Position.SetXYZT(vx, vy, dz, dt);
189 vertex->ClusterIndex = nvtx;
190 vertex->ClusterNDF = nch;
191 vertex->SumPT2 = sumpt2;
192 vertex->GenSumPT2 = sumpt2;
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 } while(entry >= allEntries);
221
222 fReader->ReadEntry(entry);
223
224 // --- Pile-up vertex smearing
225
226 fFunction->GetRandom2(dz, dt);
227
228 dt *= c_light * 1.0E3; // necessary in order to make t in mm/c
229 dz *= 1.0E3; // necessary in order to make z in mm
230
231 dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
232
233 vx = 0.0;
234 vy = 0.0;
235
236 numberOfParticles = 0;
237 sumpt2 = 0.0;
238
239 //factory = GetFactory();
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 pt = candidate->Momentum.Pt();
259
260 x -= fInputBeamSpotX;
261 y -= fInputBeamSpotY;
262 candidate->Position.SetXYZT(x, y, z + dz, t + dt);
263 candidate->Position.RotateZ(dphi);
264 candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
265
266 vx += candidate->Position.X();
267 vy += candidate->Position.Y();
268
269 ++numberOfParticles;
270 if(TMath::Abs(candidate->Charge) > 1.0E-9)
271 {
272 nch++;
273 sumpt2 += pt * pt;
274 vertex->AddCandidate(candidate);
275 }
276
277 fParticleOutputArray->Add(candidate);
278 }
279
280 if(numberOfParticles > 0)
281 {
282 vx /= numberOfParticles;
283 vy /= numberOfParticles;
284 }
285
286 nvtx++;
287
288 vertex->Position.SetXYZT(vx, vy, dz, dt);
289
290 vertex->ClusterIndex = nvtx;
291 vertex->ClusterNDF = nch;
292 vertex->SumPT2 = sumpt2;
293 vertex->GenSumPT2 = sumpt2;
294
295 vertex->IsPU = 1;
296
297 fVertexOutputArray->Add(vertex);
298 }
299}
300
301//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.