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 PileUpMergerPythia8
|
---|
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/PileUpMergerPythia8.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 "Pythia.h"
|
---|
39 |
|
---|
40 | #include "TMath.h"
|
---|
41 | #include "TString.h"
|
---|
42 | #include "TFormula.h"
|
---|
43 | #include "TRandom3.h"
|
---|
44 | #include "TObjArray.h"
|
---|
45 | #include "TDatabasePDG.h"
|
---|
46 | #include "TLorentzVector.h"
|
---|
47 |
|
---|
48 | #include <algorithm>
|
---|
49 | #include <stdexcept>
|
---|
50 | #include <iostream>
|
---|
51 | #include <sstream>
|
---|
52 |
|
---|
53 | using namespace std;
|
---|
54 |
|
---|
55 | //------------------------------------------------------------------------------
|
---|
56 |
|
---|
57 | PileUpMergerPythia8::PileUpMergerPythia8() :
|
---|
58 | fFunction(0), fPythia(0), fItInputArray(0)
|
---|
59 | {
|
---|
60 | fFunction = new DelphesTF2;
|
---|
61 | }
|
---|
62 |
|
---|
63 | //------------------------------------------------------------------------------
|
---|
64 |
|
---|
65 | PileUpMergerPythia8::~PileUpMergerPythia8()
|
---|
66 | {
|
---|
67 | delete fFunction;
|
---|
68 | }
|
---|
69 |
|
---|
70 | //------------------------------------------------------------------------------
|
---|
71 |
|
---|
72 | void PileUpMergerPythia8::Init()
|
---|
73 | {
|
---|
74 | const char *fileName;
|
---|
75 |
|
---|
76 | fPileUpDistribution = GetInt("PileUpDistribution", 0);
|
---|
77 |
|
---|
78 | fMeanPileUp = GetDouble("MeanPileUp", 10);
|
---|
79 |
|
---|
80 | fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
|
---|
81 | fTVertexSpread = GetDouble("TVertexSpread", 1.5E-09);
|
---|
82 |
|
---|
83 | fInputBeamSpotX = GetDouble("InputBeamSpotX", 0.0);
|
---|
84 | fInputBeamSpotY = GetDouble("InputBeamSpotY", 0.0);
|
---|
85 | fOutputBeamSpotX = GetDouble("OutputBeamSpotX", 0.0);
|
---|
86 | fOutputBeamSpotY = GetDouble("OutputBeamSpotY", 0.0);
|
---|
87 |
|
---|
88 | fPTMin = GetDouble("PTMin", 0.0);
|
---|
89 |
|
---|
90 | fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
|
---|
91 | fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
|
---|
92 |
|
---|
93 | fileName = GetString("ConfigFile", "MinBias.cmnd");
|
---|
94 | fPythia = new Pythia8::Pythia();
|
---|
95 | fPythia->readFile(fileName);
|
---|
96 |
|
---|
97 | // import input array
|
---|
98 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
99 | fItInputArray = fInputArray->MakeIterator();
|
---|
100 |
|
---|
101 | // create output arrays
|
---|
102 | fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
|
---|
103 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
104 | }
|
---|
105 |
|
---|
106 | //------------------------------------------------------------------------------
|
---|
107 |
|
---|
108 | void PileUpMergerPythia8::Finish()
|
---|
109 | {
|
---|
110 | if(fPythia) delete fPythia;
|
---|
111 | }
|
---|
112 |
|
---|
113 | //------------------------------------------------------------------------------
|
---|
114 |
|
---|
115 | void PileUpMergerPythia8::Process()
|
---|
116 | {
|
---|
117 | TDatabasePDG *pdg = TDatabasePDG::Instance();
|
---|
118 | TParticlePDG *pdgParticle;
|
---|
119 | Int_t pid, status;
|
---|
120 | Float_t x, y, z, t, vx, vy;
|
---|
121 | Float_t px, py, pz, e;
|
---|
122 | Double_t dz, dphi, dt;
|
---|
123 | Int_t numberOfEvents, event, numberOfParticles, i;
|
---|
124 | Candidate *candidate, *vertex;
|
---|
125 | DelphesFactory *factory;
|
---|
126 |
|
---|
127 | const Double_t c_light = 2.99792458E8;
|
---|
128 |
|
---|
129 | fItInputArray->Reset();
|
---|
130 |
|
---|
131 | // --- Deal with primary vertex first ------
|
---|
132 |
|
---|
133 | fFunction->GetRandom2(dz, dt);
|
---|
134 |
|
---|
135 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
136 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
137 | vx = 0.0;
|
---|
138 | vy = 0.0;
|
---|
139 | numberOfParticles = fInputArray->GetEntriesFast();
|
---|
140 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
141 | {
|
---|
142 | vx += candidate->Position.X();
|
---|
143 | vy += candidate->Position.Y();
|
---|
144 | z = candidate->Position.Z();
|
---|
145 | t = candidate->Position.T();
|
---|
146 | candidate->Position.SetZ(z + dz);
|
---|
147 | candidate->Position.SetT(t + dt);
|
---|
148 | fParticleOutputArray->Add(candidate);
|
---|
149 | }
|
---|
150 |
|
---|
151 | if(numberOfParticles > 0)
|
---|
152 | {
|
---|
153 | vx /= numberOfParticles;
|
---|
154 | vy /= numberOfParticles;
|
---|
155 | }
|
---|
156 |
|
---|
157 | factory = GetFactory();
|
---|
158 |
|
---|
159 | vertex = factory->NewCandidate();
|
---|
160 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
161 | fVertexOutputArray->Add(vertex);
|
---|
162 |
|
---|
163 | // --- Then with pile-up vertices ------
|
---|
164 |
|
---|
165 | switch(fPileUpDistribution)
|
---|
166 | {
|
---|
167 | case 0:
|
---|
168 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
169 | break;
|
---|
170 | case 1:
|
---|
171 | numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
|
---|
172 | break;
|
---|
173 | default:
|
---|
174 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
175 | break;
|
---|
176 | }
|
---|
177 |
|
---|
178 | for(event = 0; event < numberOfEvents; ++event)
|
---|
179 | {
|
---|
180 | while(!fPythia->next());
|
---|
181 |
|
---|
182 | // --- Pile-up vertex smearing
|
---|
183 |
|
---|
184 | fFunction->GetRandom2(dz, dt);
|
---|
185 |
|
---|
186 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
187 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
188 |
|
---|
189 | dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
|
---|
190 |
|
---|
191 | vx = 0.0;
|
---|
192 | vy = 0.0;
|
---|
193 | numberOfParticles = fPythia->event.size();
|
---|
194 | for(i = 1; i < numberOfParticles; ++i)
|
---|
195 | {
|
---|
196 | Pythia8::Particle &particle = fPythia->event[i];
|
---|
197 |
|
---|
198 | status = particle.statusHepMC();
|
---|
199 |
|
---|
200 | if(status != 1 || !particle.isVisible() || particle.pT() <= fPTMin) continue;
|
---|
201 |
|
---|
202 | pid = particle.id();
|
---|
203 | px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.e();
|
---|
204 | x = particle.xProd(); y = particle.yProd(); z = particle.zProd(); t = particle.tProd();
|
---|
205 |
|
---|
206 | candidate = factory->NewCandidate();
|
---|
207 |
|
---|
208 | candidate->PID = pid;
|
---|
209 |
|
---|
210 | candidate->Status = 1;
|
---|
211 |
|
---|
212 | pdgParticle = pdg->GetParticle(pid);
|
---|
213 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
|
---|
214 | candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
|
---|
215 |
|
---|
216 | candidate->IsPU = 1;
|
---|
217 |
|
---|
218 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
219 | candidate->Momentum.RotateZ(dphi);
|
---|
220 |
|
---|
221 | x -= fInputBeamSpotX;
|
---|
222 | y -= fInputBeamSpotY;
|
---|
223 | candidate->Position.SetXYZT(x, y, z + dz, t + dt);
|
---|
224 | candidate->Position.RotateZ(dphi);
|
---|
225 | candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
|
---|
226 |
|
---|
227 | vx += candidate->Position.X();
|
---|
228 | vy += candidate->Position.Y();
|
---|
229 |
|
---|
230 | fParticleOutputArray->Add(candidate);
|
---|
231 | }
|
---|
232 |
|
---|
233 | if(numberOfParticles > 0)
|
---|
234 | {
|
---|
235 | vx /= numberOfParticles;
|
---|
236 | vy /= numberOfParticles;
|
---|
237 | }
|
---|
238 |
|
---|
239 | vertex = factory->NewCandidate();
|
---|
240 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
241 | vertex->IsPU = 1;
|
---|
242 |
|
---|
243 | fVertexOutputArray->Add(vertex);
|
---|
244 | }
|
---|
245 | }
|
---|
246 |
|
---|
247 | //------------------------------------------------------------------------------
|
---|