Fork me on GitHub

source: git/modules/PileUpMerger.cc@ b291b0b

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b291b0b was 2118a6a, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

fixed missing vertexing variables

  • Property mode set to 100644
File size: 7.2 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/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
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55PileUpMerger::PileUpMerger() :
56 fFunction(0), fReader(0), fItInputArray(0)
57{
58 fFunction = new DelphesTF2;
59}
60
61
62//------------------------------------------------------------------------------
63
64PileUpMerger::~PileUpMerger()
65{
66 delete fFunction;
67}
68
69//------------------------------------------------------------------------------
70
71void 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
106void PileUpMerger::Finish()
107{
108 if(fReader) delete fReader;
109}
110
111//------------------------------------------------------------------------------
112
113void 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;
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 dt *= c_light*1.0E3; // necessary in order to make t in mm/c
135 dz *= 1.0E3; // necessary in order to make z in mm
136 vx = 0.0;
137 vy = 0.0;
138 numberOfParticles = fInputArray->GetEntriesFast();
139 nch = 0;
140 sumpt2 = 0.0;
141
142 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
143 {
144 vx += candidate->Position.X();
145 vy += candidate->Position.Y();
146 z = candidate->Position.Z();
147 t = candidate->Position.T();
148 pt = candidate->Momentum.Pt();
149 candidate->Position.SetZ(z + dz);
150 candidate->Position.SetT(t + dt);
151 fParticleOutputArray->Add(candidate);
152
153 if(TMath::Abs(candidate->Charge) > 1.0E-9)
154 {
155 nch++;
156 sumpt2 += pt*pt;
157 }
158 }
159
160 if(numberOfParticles > 0)
161 {
162 vx /= numberOfParticles;
163 vy /= numberOfParticles;
164 }
165
166 nvtx++;
167 factory = GetFactory();
168
169 vertex = factory->NewCandidate();
170 vertex->Position.SetXYZT(vx, vy, dz, dt);
171 vertex->ClusterIndex = nvtx;
172 vertex->ClusterNDF = nch;
173 vertex->SumPT2 = sumpt2;
174 vertex->GenSumPT2 = sumpt2;
175
176 fVertexOutputArray->Add(vertex);
177
178 // --- Then with pile-up vertices ------
179
180 switch(fPileUpDistribution)
181 {
182 case 0:
183 numberOfEvents = gRandom->Poisson(fMeanPileUp);
184 break;
185 case 1:
186 numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
187 break;
188 case 2:
189 numberOfEvents = fMeanPileUp;
190 break;
191 default:
192 numberOfEvents = gRandom->Poisson(fMeanPileUp);
193 break;
194 }
195
196 allEntries = fReader->GetEntries();
197
198 for(event = 0; event < numberOfEvents; ++event)
199 {
200 do
201 {
202 entry = TMath::Nint(gRandom->Rndm()*allEntries);
203 }
204 while(entry >= allEntries);
205
206 fReader->ReadEntry(entry);
207
208 // --- Pile-up vertex smearing
209
210 fFunction->GetRandom2(dz, dt);
211
212 dt *= c_light*1.0E3; // necessary in order to make t in mm/c
213 dz *= 1.0E3; // necessary in order to make z in mm
214
215 dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
216
217 vx = 0.0;
218 vy = 0.0;
219
220 numberOfParticles = 0;
221 sumpt2 = 0.0;
222
223 while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
224 {
225 candidate = factory->NewCandidate();
226
227 candidate->PID = pid;
228
229 candidate->Status = 1;
230
231 pdgParticle = pdg->GetParticle(pid);
232 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
233 candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
234
235 candidate->IsPU = 1;
236
237 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
238 candidate->Momentum.RotateZ(dphi);
239
240 x -= fInputBeamSpotX;
241 y -= fInputBeamSpotY;
242 candidate->Position.SetXYZT(x, y, z + dz, t + dt);
243 candidate->Position.RotateZ(dphi);
244 candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
245
246 vx += candidate->Position.X();
247 vy += candidate->Position.Y();
248
249 ++numberOfParticles;
250 if(TMath::Abs(candidate->Charge) > 1.0E-9)
251 {
252 nch++;
253 sumpt2 += pt*pt;
254 }
255
256 fParticleOutputArray->Add(candidate);
257 }
258
259 if(numberOfParticles > 0)
260 {
261 vx /= numberOfParticles;
262 vy /= numberOfParticles;
263 }
264
265 nvtx++;
266
267 vertex = factory->NewCandidate();
268 vertex->Position.SetXYZT(vx, vy, dz, dt);
269
270 vertex->ClusterIndex = nvtx;
271 vertex->ClusterNDF = nch;
272 vertex->SumPT2 = sumpt2;
273 vertex->GenSumPT2 = sumpt2;
274
275 vertex->IsPU = 1;
276
277 fVertexOutputArray->Add(vertex);
278
279 }
280}
281
282//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.