[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"
|
---|
| 31 | #include "classes/DelphesPileUpReader.h"
|
---|
[341014c] | 32 | #include "classes/DelphesTF2.h"
|
---|
[d7d2da3] | 33 |
|
---|
| 34 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
[341014c] | 35 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 36 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
[d7d2da3] | 37 |
|
---|
| 38 | #include "TDatabasePDG.h"
|
---|
[341014c] | 39 | #include "TFormula.h"
|
---|
[d7d2da3] | 40 | #include "TLorentzVector.h"
|
---|
[341014c] | 41 | #include "TMath.h"
|
---|
| 42 | #include "TObjArray.h"
|
---|
| 43 | #include "TRandom3.h"
|
---|
| 44 | #include "TString.h"
|
---|
[d7d2da3] | 45 |
|
---|
| 46 | #include <algorithm>
|
---|
| 47 | #include <iostream>
|
---|
| 48 | #include <sstream>
|
---|
[341014c] | 49 | #include <stdexcept>
|
---|
[d7d2da3] | 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 |
|
---|
| 61 | //------------------------------------------------------------------------------
|
---|
| 62 |
|
---|
| 63 | PileUpMerger::~PileUpMerger()
|
---|
| 64 | {
|
---|
[22dc7fd] | 65 | delete fFunction;
|
---|
[d7d2da3] | 66 | }
|
---|
| 67 |
|
---|
| 68 | //------------------------------------------------------------------------------
|
---|
| 69 |
|
---|
| 70 | void PileUpMerger::Init()
|
---|
| 71 | {
|
---|
| 72 | const char *fileName;
|
---|
| 73 |
|
---|
[76d3973] | 74 | fPileUpDistribution = GetInt("PileUpDistribution", 0);
|
---|
| 75 |
|
---|
[341014c] | 76 | fMeanPileUp = GetDouble("MeanPileUp", 10);
|
---|
[0ed696f] | 77 |
|
---|
[22dc7fd] | 78 | fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
|
---|
| 79 | fTVertexSpread = GetDouble("TVertexSpread", 1.5E-09);
|
---|
| 80 |
|
---|
[2d494a6] | 81 | fInputBeamSpotX = GetDouble("InputBeamSpotX", 0.0);
|
---|
| 82 | fInputBeamSpotY = GetDouble("InputBeamSpotY", 0.0);
|
---|
| 83 | fOutputBeamSpotX = GetDouble("OutputBeamSpotX", 0.0);
|
---|
| 84 | fOutputBeamSpotY = GetDouble("OutputBeamSpotY", 0.0);
|
---|
| 85 |
|
---|
[22dc7fd] | 86 | // read vertex smearing formula
|
---|
[0ed696f] | 87 |
|
---|
[22dc7fd] | 88 | fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
|
---|
[0ed696f] | 89 | fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
|
---|
| 90 |
|
---|
[d7d2da3] | 91 | fileName = GetString("PileUpFile", "MinBias.pileup");
|
---|
| 92 | fReader = new DelphesPileUpReader(fileName);
|
---|
| 93 |
|
---|
| 94 | // import input array
|
---|
| 95 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
| 96 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 97 |
|
---|
| 98 | // create output arrays
|
---|
[d07e957] | 99 | fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
|
---|
| 100 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
[d7d2da3] | 101 | }
|
---|
| 102 |
|
---|
| 103 | //------------------------------------------------------------------------------
|
---|
| 104 |
|
---|
| 105 | void PileUpMerger::Finish()
|
---|
| 106 | {
|
---|
| 107 | if(fReader) delete fReader;
|
---|
| 108 | }
|
---|
| 109 |
|
---|
| 110 | //------------------------------------------------------------------------------
|
---|
| 111 |
|
---|
| 112 | void PileUpMerger::Process()
|
---|
| 113 | {
|
---|
| 114 | TDatabasePDG *pdg = TDatabasePDG::Instance();
|
---|
| 115 | TParticlePDG *pdgParticle;
|
---|
[2118a6a] | 116 | Int_t pid, nch, nvtx = -1;
|
---|
[8d601b7] | 117 | Float_t x, y, z, t, vx, vy;
|
---|
[2118a6a] | 118 | Float_t px, py, pz, e, pt;
|
---|
[50edcdbf] | 119 | Double_t dz, dphi, dt, sumpt2, dz0, dt0;
|
---|
[2d494a6] | 120 | Int_t numberOfEvents, event, numberOfParticles;
|
---|
[d7d2da3] | 121 | Long64_t allEntries, entry;
|
---|
[2d494a6] | 122 | Candidate *candidate, *vertex;
|
---|
[d7d2da3] | 123 | DelphesFactory *factory;
|
---|
[0ed696f] | 124 |
|
---|
[22dc7fd] | 125 | const Double_t c_light = 2.99792458E8;
|
---|
[d7d2da3] | 126 |
|
---|
| 127 | fItInputArray->Reset();
|
---|
[0ed696f] | 128 |
|
---|
[2d494a6] | 129 | // --- Deal with primary vertex first ------
|
---|
[0ed696f] | 130 |
|
---|
| 131 | fFunction->GetRandom2(dz, dt);
|
---|
| 132 |
|
---|
[50edcdbf] | 133 | dz0 = -1.0e6;
|
---|
| 134 | dt0 = -1.0e6;
|
---|
[5496767] | 135 |
|
---|
[341014c] | 136 | dt *= c_light * 1.0E3; // necessary in order to make t in mm/c
|
---|
[0ed696f] | 137 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
[5496767] | 138 |
|
---|
[8d601b7] | 139 | //cout<<dz<<","<<dt<<endl;
|
---|
[5496767] | 140 |
|
---|
[2d494a6] | 141 | vx = 0.0;
|
---|
| 142 | vy = 0.0;
|
---|
[5496767] | 143 |
|
---|
[2d494a6] | 144 | numberOfParticles = fInputArray->GetEntriesFast();
|
---|
[2118a6a] | 145 | nch = 0;
|
---|
[5496767] | 146 | sumpt2 = 0.0;
|
---|
| 147 |
|
---|
[62d3bc5] | 148 | factory = GetFactory();
|
---|
| 149 | vertex = factory->NewCandidate();
|
---|
| 150 |
|
---|
[341014c] | 151 | while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
|
---|
[d7d2da3] | 152 | {
|
---|
[2d494a6] | 153 | vx += candidate->Position.X();
|
---|
| 154 | vy += candidate->Position.Y();
|
---|
[00078bc] | 155 | z = candidate->Position.Z();
|
---|
| 156 | t = candidate->Position.T();
|
---|
[2118a6a] | 157 | pt = candidate->Momentum.Pt();
|
---|
[5496767] | 158 |
|
---|
[50edcdbf] | 159 | // take postion and time from first stable particle
|
---|
[341014c] | 160 | if(dz0 < -999999.0)
|
---|
[50edcdbf] | 161 | dz0 = z;
|
---|
[341014c] | 162 | if(dt0 < -999999.0)
|
---|
[50edcdbf] | 163 | dt0 = t;
|
---|
| 164 |
|
---|
| 165 | // cancel any possible offset in position and time the input file
|
---|
| 166 | candidate->Position.SetZ(z - dz0 + dz);
|
---|
| 167 | candidate->Position.SetT(t - dt0 + dt);
|
---|
[5496767] | 168 |
|
---|
[3e2bb2b] | 169 | candidate->IsPU = 0;
|
---|
| 170 |
|
---|
[d07e957] | 171 | fParticleOutputArray->Add(candidate);
|
---|
[5496767] | 172 |
|
---|
[341014c] | 173 | if(TMath::Abs(candidate->Charge) > 1.0E-9)
|
---|
[2118a6a] | 174 | {
|
---|
[5496767] | 175 | nch++;
|
---|
[341014c] | 176 | sumpt2 += pt * pt;
|
---|
[62d3bc5] | 177 | vertex->AddCandidate(candidate);
|
---|
[5496767] | 178 | }
|
---|
[d7d2da3] | 179 | }
|
---|
[5496767] | 180 |
|
---|
[2d494a6] | 181 | if(numberOfParticles > 0)
|
---|
| 182 | {
|
---|
[8d601b7] | 183 | vx /= sumpt2;
|
---|
| 184 | vy /= sumpt2;
|
---|
[2d494a6] | 185 | }
|
---|
[5496767] | 186 |
|
---|
[2118a6a] | 187 | nvtx++;
|
---|
[8d601b7] | 188 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
[2118a6a] | 189 | vertex->ClusterIndex = nvtx;
|
---|
| 190 | vertex->ClusterNDF = nch;
|
---|
| 191 | vertex->SumPT2 = sumpt2;
|
---|
[5496767] | 192 | vertex->GenSumPT2 = sumpt2;
|
---|
[2d494a6] | 193 | fVertexOutputArray->Add(vertex);
|
---|
[0ed696f] | 194 |
|
---|
[22dc7fd] | 195 | // --- Then with pile-up vertices ------
|
---|
[0ed696f] | 196 |
|
---|
[76d3973] | 197 | switch(fPileUpDistribution)
|
---|
| 198 | {
|
---|
[341014c] | 199 | case 0:
|
---|
| 200 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
| 201 | break;
|
---|
| 202 | case 1:
|
---|
| 203 | numberOfEvents = gRandom->Integer(2 * fMeanPileUp + 1);
|
---|
| 204 | break;
|
---|
| 205 | case 2:
|
---|
| 206 | numberOfEvents = fMeanPileUp;
|
---|
| 207 | break;
|
---|
| 208 | default:
|
---|
| 209 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
| 210 | break;
|
---|
[76d3973] | 211 | }
|
---|
[d7d2da3] | 212 |
|
---|
| 213 | allEntries = fReader->GetEntries();
|
---|
[0ed696f] | 214 |
|
---|
[76d3973] | 215 | for(event = 0; event < numberOfEvents; ++event)
|
---|
[d7d2da3] | 216 | {
|
---|
| 217 | do
|
---|
| 218 | {
|
---|
[341014c] | 219 | entry = TMath::Nint(gRandom->Rndm() * allEntries);
|
---|
| 220 | } while(entry >= allEntries);
|
---|
[d7d2da3] | 221 |
|
---|
| 222 | fReader->ReadEntry(entry);
|
---|
[0ed696f] | 223 |
|
---|
[341014c] | 224 | // --- Pile-up vertex smearing
|
---|
[0ed696f] | 225 |
|
---|
| 226 | fFunction->GetRandom2(dz, dt);
|
---|
| 227 |
|
---|
[341014c] | 228 | dt *= c_light * 1.0E3; // necessary in order to make t in mm/c
|
---|
[0ed696f] | 229 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
| 230 |
|
---|
[b68a60f] | 231 | dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
|
---|
[d7d2da3] | 232 |
|
---|
[2d494a6] | 233 | vx = 0.0;
|
---|
| 234 | vy = 0.0;
|
---|
[5496767] | 235 |
|
---|
[2d494a6] | 236 | numberOfParticles = 0;
|
---|
[2118a6a] | 237 | sumpt2 = 0.0;
|
---|
[5496767] | 238 |
|
---|
[3e2bb2b] | 239 | //factory = GetFactory();
|
---|
[5496767] | 240 | vertex = factory->NewCandidate();
|
---|
| 241 |
|
---|
[d7d2da3] | 242 | while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
|
---|
| 243 | {
|
---|
| 244 | candidate = factory->NewCandidate();
|
---|
| 245 |
|
---|
| 246 | candidate->PID = pid;
|
---|
| 247 |
|
---|
| 248 | candidate->Status = 1;
|
---|
| 249 |
|
---|
| 250 | pdgParticle = pdg->GetParticle(pid);
|
---|
[341014c] | 251 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge() / 3.0) : -999;
|
---|
[d7d2da3] | 252 | candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
|
---|
| 253 |
|
---|
| 254 | candidate->IsPU = 1;
|
---|
[0ed696f] | 255 |
|
---|
[d7d2da3] | 256 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
[b68a60f] | 257 | candidate->Momentum.RotateZ(dphi);
|
---|
[3e2bb2b] | 258 | pt = candidate->Momentum.Pt();
|
---|
[b68a60f] | 259 |
|
---|
[2d494a6] | 260 | x -= fInputBeamSpotX;
|
---|
| 261 | y -= fInputBeamSpotY;
|
---|
[00078bc] | 262 | candidate->Position.SetXYZT(x, y, z + dz, t + dt);
|
---|
[b68a60f] | 263 | candidate->Position.RotateZ(dphi);
|
---|
[2d494a6] | 264 | candidate->Position += TLorentzVector(fOutputBeamSpotX, fOutputBeamSpotY, 0.0, 0.0);
|
---|
| 265 |
|
---|
| 266 | vx += candidate->Position.X();
|
---|
| 267 | vy += candidate->Position.Y();
|
---|
[5496767] | 268 |
|
---|
[2d494a6] | 269 | ++numberOfParticles;
|
---|
[341014c] | 270 | if(TMath::Abs(candidate->Charge) > 1.0E-9)
|
---|
[2118a6a] | 271 | {
|
---|
[5496767] | 272 | nch++;
|
---|
[341014c] | 273 | sumpt2 += pt * pt;
|
---|
[5496767] | 274 | vertex->AddCandidate(candidate);
|
---|
| 275 | }
|
---|
| 276 |
|
---|
[d07e957] | 277 | fParticleOutputArray->Add(candidate);
|
---|
[d7d2da3] | 278 | }
|
---|
[2d494a6] | 279 |
|
---|
| 280 | if(numberOfParticles > 0)
|
---|
| 281 | {
|
---|
| 282 | vx /= numberOfParticles;
|
---|
| 283 | vy /= numberOfParticles;
|
---|
| 284 | }
|
---|
[5496767] | 285 |
|
---|
[2118a6a] | 286 | nvtx++;
|
---|
[2d494a6] | 287 |
|
---|
| 288 | vertex->Position.SetXYZT(vx, vy, dz, dt);
|
---|
[5496767] | 289 |
|
---|
| 290 | vertex->ClusterIndex = nvtx;
|
---|
[2118a6a] | 291 | vertex->ClusterNDF = nch;
|
---|
| 292 | vertex->SumPT2 = sumpt2;
|
---|
| 293 | vertex->GenSumPT2 = sumpt2;
|
---|
[5496767] | 294 |
|
---|
[2d494a6] | 295 | vertex->IsPU = 1;
|
---|
| 296 |
|
---|
| 297 | fVertexOutputArray->Add(vertex);
|
---|
[d7d2da3] | 298 | }
|
---|
| 299 | }
|
---|
| 300 |
|
---|
| 301 | //------------------------------------------------------------------------------
|
---|