[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 |
|
---|
| 51 | using namespace std;
|
---|
| 52 |
|
---|
| 53 | //------------------------------------------------------------------------------
|
---|
| 54 |
|
---|
| 55 | PileUpMerger::PileUpMerger() :
|
---|
[22dc7fd] | 56 | fFunction(0), fReader(0), fItInputArray(0)
|
---|
[d7d2da3] | 57 | {
|
---|
[22dc7fd] | 58 | fFunction = new DelphesTF2;
|
---|
[d7d2da3] | 59 | }
|
---|
| 60 |
|
---|
[22dc7fd] | 61 |
|
---|
[d7d2da3] | 62 | //------------------------------------------------------------------------------
|
---|
| 63 |
|
---|
| 64 | PileUpMerger::~PileUpMerger()
|
---|
| 65 | {
|
---|
[22dc7fd] | 66 | delete fFunction;
|
---|
[d7d2da3] | 67 | }
|
---|
| 68 |
|
---|
| 69 | //------------------------------------------------------------------------------
|
---|
| 70 |
|
---|
| 71 | void 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 |
|
---|
| 106 | void PileUpMerger::Finish()
|
---|
| 107 | {
|
---|
| 108 | if(fReader) delete fReader;
|
---|
| 109 | }
|
---|
| 110 |
|
---|
| 111 | //------------------------------------------------------------------------------
|
---|
| 112 |
|
---|
| 113 | void PileUpMerger::Process()
|
---|
| 114 | {
|
---|
| 115 | TDatabasePDG *pdg = TDatabasePDG::Instance();
|
---|
| 116 | TParticlePDG *pdgParticle;
|
---|
[2118a6a] | 117 | Int_t pid, nch, nvtx = -1;
|
---|
[2d494a6] | 118 | Float_t x, y, z, t, vx, vy;
|
---|
[2118a6a] | 119 | Float_t px, py, pz, e, pt;
|
---|
| 120 | Double_t dz, dphi, dt, sumpt2;
|
---|
[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 |
|
---|
| 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
|
---|
[2d494a6] | 136 | vx = 0.0;
|
---|
| 137 | vy = 0.0;
|
---|
| 138 | numberOfParticles = fInputArray->GetEntriesFast();
|
---|
[2118a6a] | 139 | nch = 0;
|
---|
| 140 | sumpt2 = 0.0;
|
---|
| 141 |
|
---|
[d7d2da3] | 142 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
| 143 | {
|
---|
[2d494a6] | 144 | vx += candidate->Position.X();
|
---|
| 145 | vy += candidate->Position.Y();
|
---|
[00078bc] | 146 | z = candidate->Position.Z();
|
---|
| 147 | t = candidate->Position.T();
|
---|
[2118a6a] | 148 | pt = candidate->Momentum.Pt();
|
---|
[fb8932a] | 149 | candidate->Position.SetZ(z + dz);
|
---|
| 150 | candidate->Position.SetT(t + dt);
|
---|
[d07e957] | 151 | fParticleOutputArray->Add(candidate);
|
---|
[2118a6a] | 152 |
|
---|
| 153 | if(TMath::Abs(candidate->Charge) > 1.0E-9)
|
---|
| 154 | {
|
---|
| 155 | nch++;
|
---|
| 156 | sumpt2 += pt*pt;
|
---|
| 157 | }
|
---|
[d7d2da3] | 158 | }
|
---|
| 159 |
|
---|
[2d494a6] | 160 | if(numberOfParticles > 0)
|
---|
| 161 | {
|
---|
| 162 | vx /= numberOfParticles;
|
---|
| 163 | vy /= numberOfParticles;
|
---|
| 164 | }
|
---|
| 165 |
|
---|
[2118a6a] | 166 | nvtx++;
|
---|
[d7d2da3] | 167 | factory = GetFactory();
|
---|
[0ed696f] | 168 |
|
---|
[2d494a6] | 169 | vertex = factory->NewCandidate();
|
---|
| 170 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
[2118a6a] | 171 | vertex->ClusterIndex = nvtx;
|
---|
| 172 | vertex->ClusterNDF = nch;
|
---|
| 173 | vertex->SumPT2 = sumpt2;
|
---|
| 174 | vertex->GenSumPT2 = sumpt2;
|
---|
| 175 |
|
---|
[2d494a6] | 176 | fVertexOutputArray->Add(vertex);
|
---|
[0ed696f] | 177 |
|
---|
[22dc7fd] | 178 | // --- Then with pile-up vertices ------
|
---|
[0ed696f] | 179 |
|
---|
[76d3973] | 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;
|
---|
[76c2a3b] | 188 | case 2:
|
---|
| 189 | numberOfEvents = fMeanPileUp;
|
---|
| 190 | break;
|
---|
[76d3973] | 191 | default:
|
---|
| 192 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
| 193 | break;
|
---|
| 194 | }
|
---|
[d7d2da3] | 195 |
|
---|
| 196 | allEntries = fReader->GetEntries();
|
---|
[0ed696f] | 197 |
|
---|
[76d3973] | 198 | for(event = 0; event < numberOfEvents; ++event)
|
---|
[d7d2da3] | 199 | {
|
---|
| 200 | do
|
---|
| 201 | {
|
---|
| 202 | entry = TMath::Nint(gRandom->Rndm()*allEntries);
|
---|
| 203 | }
|
---|
| 204 | while(entry >= allEntries);
|
---|
| 205 |
|
---|
| 206 | fReader->ReadEntry(entry);
|
---|
[0ed696f] | 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 |
|
---|
[b68a60f] | 215 | dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
|
---|
[d7d2da3] | 216 |
|
---|
[2d494a6] | 217 | vx = 0.0;
|
---|
| 218 | vy = 0.0;
|
---|
[2118a6a] | 219 |
|
---|
[2d494a6] | 220 | numberOfParticles = 0;
|
---|
[2118a6a] | 221 | sumpt2 = 0.0;
|
---|
| 222 |
|
---|
[d7d2da3] | 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;
|
---|
[0ed696f] | 236 |
|
---|
[d7d2da3] | 237 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
[b68a60f] | 238 | candidate->Momentum.RotateZ(dphi);
|
---|
| 239 |
|
---|
[2d494a6] | 240 | x -= fInputBeamSpotX;
|
---|
| 241 | y -= fInputBeamSpotY;
|
---|
[00078bc] | 242 | candidate->Position.SetXYZT(x, y, z + dz, t + dt);
|
---|
[b68a60f] | 243 | candidate->Position.RotateZ(dphi);
|
---|
[2d494a6] | 244 | candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
|
---|
| 245 |
|
---|
| 246 | vx += candidate->Position.X();
|
---|
| 247 | vy += candidate->Position.Y();
|
---|
[2118a6a] | 248 |
|
---|
[2d494a6] | 249 | ++numberOfParticles;
|
---|
[2118a6a] | 250 | if(TMath::Abs(candidate->Charge) > 1.0E-9)
|
---|
| 251 | {
|
---|
| 252 | nch++;
|
---|
| 253 | sumpt2 += pt*pt;
|
---|
| 254 | }
|
---|
| 255 |
|
---|
[d07e957] | 256 | fParticleOutputArray->Add(candidate);
|
---|
[d7d2da3] | 257 | }
|
---|
[2d494a6] | 258 |
|
---|
| 259 | if(numberOfParticles > 0)
|
---|
| 260 | {
|
---|
| 261 | vx /= numberOfParticles;
|
---|
| 262 | vy /= numberOfParticles;
|
---|
| 263 | }
|
---|
[2118a6a] | 264 |
|
---|
| 265 | nvtx++;
|
---|
[2d494a6] | 266 |
|
---|
| 267 | vertex = factory->NewCandidate();
|
---|
| 268 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
[2118a6a] | 269 |
|
---|
| 270 | vertex->ClusterIndex = nvtx;
|
---|
| 271 | vertex->ClusterNDF = nch;
|
---|
| 272 | vertex->SumPT2 = sumpt2;
|
---|
| 273 | vertex->GenSumPT2 = sumpt2;
|
---|
| 274 |
|
---|
[2d494a6] | 275 | vertex->IsPU = 1;
|
---|
| 276 |
|
---|
| 277 | fVertexOutputArray->Add(vertex);
|
---|
[2118a6a] | 278 |
|
---|
[d7d2da3] | 279 | }
|
---|
| 280 | }
|
---|
| 281 |
|
---|
| 282 | //------------------------------------------------------------------------------
|
---|