1 | /** \class PileUpMerger
|
---|
2 | *
|
---|
3 | * Merges particles from pile-up sample into event
|
---|
4 | *
|
---|
5 | *
|
---|
6 | * $Date: 2013-02-12 15:13:59 +0100 (Tue, 12 Feb 2013) $
|
---|
7 | * $Revision: 907 $
|
---|
8 | *
|
---|
9 | *
|
---|
10 | * \author M. Selvaggi - UCL, Louvain-la-Neuve
|
---|
11 | *
|
---|
12 | */
|
---|
13 |
|
---|
14 | #include "modules/PileUpMerger.h"
|
---|
15 |
|
---|
16 | #include "classes/DelphesClasses.h"
|
---|
17 | #include "classes/DelphesFactory.h"
|
---|
18 | #include "classes/DelphesTF2.h"
|
---|
19 | #include "classes/DelphesPileUpReader.h"
|
---|
20 |
|
---|
21 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
22 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
23 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
24 |
|
---|
25 | #include "TMath.h"
|
---|
26 | #include "TString.h"
|
---|
27 | #include "TFormula.h"
|
---|
28 | #include "TRandom3.h"
|
---|
29 | #include "TObjArray.h"
|
---|
30 | #include "TDatabasePDG.h"
|
---|
31 | #include "TLorentzVector.h"
|
---|
32 |
|
---|
33 | #include <algorithm>
|
---|
34 | #include <stdexcept>
|
---|
35 | #include <iostream>
|
---|
36 | #include <sstream>
|
---|
37 |
|
---|
38 | using namespace std;
|
---|
39 |
|
---|
40 | //------------------------------------------------------------------------------
|
---|
41 |
|
---|
42 | PileUpMerger::PileUpMerger() :
|
---|
43 | fFunction(0), fReader(0), fItInputArray(0)
|
---|
44 | {
|
---|
45 | fFunction = new DelphesTF2;
|
---|
46 | }
|
---|
47 |
|
---|
48 |
|
---|
49 | //------------------------------------------------------------------------------
|
---|
50 |
|
---|
51 | PileUpMerger::~PileUpMerger()
|
---|
52 | {
|
---|
53 | delete fFunction;
|
---|
54 | }
|
---|
55 |
|
---|
56 | //------------------------------------------------------------------------------
|
---|
57 |
|
---|
58 | void PileUpMerger::Init()
|
---|
59 | {
|
---|
60 | const char *fileName;
|
---|
61 |
|
---|
62 | fPileUpDistribution = GetInt("PileUpDistribution", 0);
|
---|
63 |
|
---|
64 | fMeanPileUp = GetDouble("MeanPileUp", 10);
|
---|
65 |
|
---|
66 | fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
|
---|
67 | fTVertexSpread = GetDouble("TVertexSpread", 1.5E-09);
|
---|
68 |
|
---|
69 | // read vertex smearing formula
|
---|
70 |
|
---|
71 | fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
|
---|
72 | fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
|
---|
73 |
|
---|
74 | fileName = GetString("PileUpFile", "MinBias.pileup");
|
---|
75 | fReader = new DelphesPileUpReader(fileName);
|
---|
76 |
|
---|
77 | // import input array
|
---|
78 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
79 | fItInputArray = fInputArray->MakeIterator();
|
---|
80 |
|
---|
81 | // create output arrays
|
---|
82 | fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
|
---|
83 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
84 | }
|
---|
85 |
|
---|
86 | //------------------------------------------------------------------------------
|
---|
87 |
|
---|
88 | void PileUpMerger::Finish()
|
---|
89 | {
|
---|
90 | if(fReader) delete fReader;
|
---|
91 | }
|
---|
92 |
|
---|
93 | //------------------------------------------------------------------------------
|
---|
94 |
|
---|
95 | void PileUpMerger::Process()
|
---|
96 | {
|
---|
97 | TDatabasePDG *pdg = TDatabasePDG::Instance();
|
---|
98 | TParticlePDG *pdgParticle;
|
---|
99 | Int_t pid;
|
---|
100 | Float_t x, y, z, t;
|
---|
101 | Float_t px, py, pz, e;
|
---|
102 | Double_t dz, dphi, dt;
|
---|
103 | Int_t numberOfEvents, event;
|
---|
104 | Long64_t allEntries, entry;
|
---|
105 | Candidate *candidate, *vertexcandidate;
|
---|
106 | DelphesFactory *factory;
|
---|
107 |
|
---|
108 | const Double_t c_light = 2.99792458E8;
|
---|
109 |
|
---|
110 | fItInputArray->Reset();
|
---|
111 |
|
---|
112 | // --- Deal with Primary vertex first ------
|
---|
113 |
|
---|
114 | fFunction->GetRandom2(dz, dt);
|
---|
115 |
|
---|
116 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
117 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
118 |
|
---|
119 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
120 | {
|
---|
121 | z = candidate->Position.Z();
|
---|
122 | t = candidate->Position.T();
|
---|
123 | candidate->Position.SetZ(z + dz);
|
---|
124 | candidate->Position.SetT(t + dt);
|
---|
125 | fParticleOutputArray->Add(candidate);
|
---|
126 | }
|
---|
127 |
|
---|
128 | factory = GetFactory();
|
---|
129 |
|
---|
130 | vertexcandidate = factory->NewCandidate();
|
---|
131 | vertexcandidate->Position.SetXYZT(0.0, 0.0, dz, dt);
|
---|
132 | fVertexOutputArray->Add(vertexcandidate);
|
---|
133 |
|
---|
134 | // --- Then with pile-up vertices ------
|
---|
135 |
|
---|
136 | switch(fPileUpDistribution)
|
---|
137 | {
|
---|
138 | case 0:
|
---|
139 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
140 | break;
|
---|
141 | case 1:
|
---|
142 | numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
|
---|
143 | break;
|
---|
144 | default:
|
---|
145 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
146 | break;
|
---|
147 | }
|
---|
148 |
|
---|
149 | allEntries = fReader->GetEntries();
|
---|
150 |
|
---|
151 | for(event = 0; event < numberOfEvents; ++event)
|
---|
152 | {
|
---|
153 | do
|
---|
154 | {
|
---|
155 | entry = TMath::Nint(gRandom->Rndm()*allEntries);
|
---|
156 | }
|
---|
157 | while(entry >= allEntries);
|
---|
158 |
|
---|
159 | fReader->ReadEntry(entry);
|
---|
160 |
|
---|
161 | // --- Pile-up vertex smearing
|
---|
162 |
|
---|
163 | fFunction->GetRandom2(dz, dt);
|
---|
164 |
|
---|
165 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
166 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
167 |
|
---|
168 | dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
|
---|
169 |
|
---|
170 | vertexcandidate = factory->NewCandidate();
|
---|
171 | vertexcandidate->Position.SetXYZT(0.0, 0.0, dz, dt);
|
---|
172 | vertexcandidate->IsPU = 1;
|
---|
173 |
|
---|
174 | fVertexOutputArray->Add(vertexcandidate);
|
---|
175 |
|
---|
176 | while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
|
---|
177 | {
|
---|
178 | candidate = factory->NewCandidate();
|
---|
179 |
|
---|
180 | candidate->PID = pid;
|
---|
181 |
|
---|
182 | candidate->Status = 1;
|
---|
183 |
|
---|
184 | pdgParticle = pdg->GetParticle(pid);
|
---|
185 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
|
---|
186 | candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
|
---|
187 |
|
---|
188 | candidate->IsPU = 1;
|
---|
189 |
|
---|
190 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
191 | candidate->Momentum.RotateZ(dphi);
|
---|
192 |
|
---|
193 | candidate->Position.SetXYZT(x, y, z + dz, t + dt);
|
---|
194 | candidate->Position.RotateZ(dphi);
|
---|
195 |
|
---|
196 | fParticleOutputArray->Add(candidate);
|
---|
197 | }
|
---|
198 | }
|
---|
199 | }
|
---|
200 |
|
---|
201 | //------------------------------------------------------------------------------
|
---|
202 |
|
---|