Fork me on GitHub

source: git/modules/PileUpMerger.cc@ 5a989eb

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5a989eb was 50edcdbf, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

reset eventual z,t shift

  • Property mode set to 100644
File size: 7.6 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[d7d2da3]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"
[22dc7fd]31#include "classes/DelphesTF2.h"
[d7d2da3]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() :
[22dc7fd]56 fFunction(0), fReader(0), fItInputArray(0)
[d7d2da3]57{
[22dc7fd]58 fFunction = new DelphesTF2;
[d7d2da3]59}
60
[22dc7fd]61
[d7d2da3]62//------------------------------------------------------------------------------
63
64PileUpMerger::~PileUpMerger()
65{
[22dc7fd]66 delete fFunction;
[d7d2da3]67}
68
69//------------------------------------------------------------------------------
70
71void PileUpMerger::Init()
72{
73 const char *fileName;
74
[76d3973]75 fPileUpDistribution = GetInt("PileUpDistribution", 0);
76
[d7d2da3]77 fMeanPileUp = GetDouble("MeanPileUp", 10);
[0ed696f]78
[22dc7fd]79 fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
80 fTVertexSpread = GetDouble("TVertexSpread", 1.5E-09);
81
[2d494a6]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
[22dc7fd]87 // read vertex smearing formula
[0ed696f]88
[22dc7fd]89 fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
[0ed696f]90 fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
91
[d7d2da3]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
[d07e957]100 fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
101 fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
[d7d2da3]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;
[2118a6a]117 Int_t pid, nch, nvtx = -1;
[6fbd089]118 Float_t x, y, z, t, vx, vy, vz, vt;
[2118a6a]119 Float_t px, py, pz, e, pt;
[50edcdbf]120 Double_t dz, dphi, dt, sumpt2, dz0, dt0;
[2d494a6]121 Int_t numberOfEvents, event, numberOfParticles;
[d7d2da3]122 Long64_t allEntries, entry;
[2d494a6]123 Candidate *candidate, *vertex;
[d7d2da3]124 DelphesFactory *factory;
[0ed696f]125
[22dc7fd]126 const Double_t c_light = 2.99792458E8;
[d7d2da3]127
128 fItInputArray->Reset();
[0ed696f]129
[2d494a6]130 // --- Deal with primary vertex first ------
[0ed696f]131
132 fFunction->GetRandom2(dz, dt);
133
[50edcdbf]134 dz0 = -1.0e6;
135 dt0 = -1.0e6;
136
[0ed696f]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
[50edcdbf]139
[2d494a6]140 vx = 0.0;
141 vy = 0.0;
[6fbd089]142 vz = 0.0;
143 vt = 0.0;
144
[2d494a6]145 numberOfParticles = fInputArray->GetEntriesFast();
[2118a6a]146 nch = 0;
147 sumpt2 = 0.0;
148
[d7d2da3]149 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
150 {
[2d494a6]151 vx += candidate->Position.X();
152 vy += candidate->Position.Y();
[00078bc]153 z = candidate->Position.Z();
154 t = candidate->Position.T();
[2118a6a]155 pt = candidate->Momentum.Pt();
[6fbd089]156
[50edcdbf]157 // take postion and time from first stable particle
158 if (dz0 < -999999.0)
159 dz0 = z;
160 if (dt0 < -999999.0)
161 dt0 = t;
162
163 // cancel any possible offset in position and time the input file
164 candidate->Position.SetZ(z - dz0 + dz);
165 candidate->Position.SetT(t - dt0 + dt);
166
167 vz += z - dz0 + dz;
168 vt += t - dt0 + dt;
[6fbd089]169
[d07e957]170 fParticleOutputArray->Add(candidate);
[2118a6a]171
172 if(TMath::Abs(candidate->Charge) > 1.0E-9)
173 {
174 nch++;
175 sumpt2 += pt*pt;
176 }
[d7d2da3]177 }
178
[2d494a6]179 if(numberOfParticles > 0)
180 {
181 vx /= numberOfParticles;
182 vy /= numberOfParticles;
[6fbd089]183 vz /= numberOfParticles;
[50edcdbf]184 vt /= numberOfParticles;
[2d494a6]185 }
186
[2118a6a]187 nvtx++;
[d7d2da3]188 factory = GetFactory();
[0ed696f]189
[2d494a6]190 vertex = factory->NewCandidate();
[6fbd089]191 vertex->Position.SetXYZT(vx, vy, vz, vt);
[2118a6a]192 vertex->ClusterIndex = nvtx;
193 vertex->ClusterNDF = nch;
194 vertex->SumPT2 = sumpt2;
[6fbd089]195 vertex->GenSumPT2 = sumpt2;
[2118a6a]196
[2d494a6]197 fVertexOutputArray->Add(vertex);
[0ed696f]198
[22dc7fd]199 // --- Then with pile-up vertices ------
[0ed696f]200
[76d3973]201 switch(fPileUpDistribution)
202 {
203 case 0:
204 numberOfEvents = gRandom->Poisson(fMeanPileUp);
205 break;
206 case 1:
207 numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
208 break;
[76c2a3b]209 case 2:
210 numberOfEvents = fMeanPileUp;
211 break;
[76d3973]212 default:
213 numberOfEvents = gRandom->Poisson(fMeanPileUp);
214 break;
215 }
[d7d2da3]216
217 allEntries = fReader->GetEntries();
[0ed696f]218
[76d3973]219 for(event = 0; event < numberOfEvents; ++event)
[d7d2da3]220 {
221 do
222 {
223 entry = TMath::Nint(gRandom->Rndm()*allEntries);
224 }
225 while(entry >= allEntries);
226
227 fReader->ReadEntry(entry);
[0ed696f]228
229 // --- Pile-up vertex smearing
230
231 fFunction->GetRandom2(dz, dt);
232
233 dt *= c_light*1.0E3; // necessary in order to make t in mm/c
234 dz *= 1.0E3; // necessary in order to make z in mm
235
[b68a60f]236 dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
[d7d2da3]237
[2d494a6]238 vx = 0.0;
239 vy = 0.0;
[2118a6a]240
[2d494a6]241 numberOfParticles = 0;
[2118a6a]242 sumpt2 = 0.0;
243
[d7d2da3]244 while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
245 {
246 candidate = factory->NewCandidate();
247
248 candidate->PID = pid;
249
250 candidate->Status = 1;
251
252 pdgParticle = pdg->GetParticle(pid);
253 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
254 candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
255
256 candidate->IsPU = 1;
[0ed696f]257
[d7d2da3]258 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
[b68a60f]259 candidate->Momentum.RotateZ(dphi);
260
[2d494a6]261 x -= fInputBeamSpotX;
262 y -= fInputBeamSpotY;
[00078bc]263 candidate->Position.SetXYZT(x, y, z + dz, t + dt);
[b68a60f]264 candidate->Position.RotateZ(dphi);
[2d494a6]265 candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
266
267 vx += candidate->Position.X();
268 vy += candidate->Position.Y();
[50edcdbf]269 vz += z+dz;
270 vt += t+dt;
[2118a6a]271
[2d494a6]272 ++numberOfParticles;
[2118a6a]273 if(TMath::Abs(candidate->Charge) > 1.0E-9)
274 {
275 nch++;
276 sumpt2 += pt*pt;
277 }
278
[d07e957]279 fParticleOutputArray->Add(candidate);
[d7d2da3]280 }
[2d494a6]281
282 if(numberOfParticles > 0)
283 {
284 vx /= numberOfParticles;
285 vy /= numberOfParticles;
286 }
[2118a6a]287
288 nvtx++;
[2d494a6]289
290 vertex = factory->NewCandidate();
291 vertex->Position.SetXYZT(vx, vy, dz, dt);
[2118a6a]292
293 vertex->ClusterIndex = nvtx;
294 vertex->ClusterNDF = nch;
295 vertex->SumPT2 = sumpt2;
296 vertex->GenSumPT2 = sumpt2;
297
[2d494a6]298 vertex->IsPU = 1;
299
300 fVertexOutputArray->Add(vertex);
[2118a6a]301
[d7d2da3]302 }
303}
304
305//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.