[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 |
|
---|
| 82 | // read vertex smearing formula
|
---|
[0ed696f] | 83 |
|
---|
[22dc7fd] | 84 | fFunction->Compile(GetString("VertexDistributionFormula", "0.0"));
|
---|
[0ed696f] | 85 | fFunction->SetRange(-fZVertexSpread, -fTVertexSpread, fZVertexSpread, fTVertexSpread);
|
---|
| 86 |
|
---|
[d7d2da3] | 87 | fileName = GetString("PileUpFile", "MinBias.pileup");
|
---|
| 88 | fReader = new DelphesPileUpReader(fileName);
|
---|
| 89 |
|
---|
| 90 | // import input array
|
---|
| 91 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
| 92 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 93 |
|
---|
| 94 | // create output arrays
|
---|
[d07e957] | 95 | fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles"));
|
---|
| 96 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
[d7d2da3] | 97 | }
|
---|
| 98 |
|
---|
| 99 | //------------------------------------------------------------------------------
|
---|
| 100 |
|
---|
| 101 | void PileUpMerger::Finish()
|
---|
| 102 | {
|
---|
| 103 | if(fReader) delete fReader;
|
---|
| 104 | }
|
---|
| 105 |
|
---|
| 106 | //------------------------------------------------------------------------------
|
---|
| 107 |
|
---|
| 108 | void PileUpMerger::Process()
|
---|
| 109 | {
|
---|
| 110 | TDatabasePDG *pdg = TDatabasePDG::Instance();
|
---|
| 111 | TParticlePDG *pdgParticle;
|
---|
| 112 | Int_t pid;
|
---|
| 113 | Float_t x, y, z, t;
|
---|
| 114 | Float_t px, py, pz, e;
|
---|
[22dc7fd] | 115 | Double_t dz, dphi, dt;
|
---|
[76d3973] | 116 | Int_t numberOfEvents, event;
|
---|
[d7d2da3] | 117 | Long64_t allEntries, entry;
|
---|
[22dc7fd] | 118 | Candidate *candidate, *vertexcandidate;
|
---|
[d7d2da3] | 119 | DelphesFactory *factory;
|
---|
[0ed696f] | 120 |
|
---|
[22dc7fd] | 121 | const Double_t c_light = 2.99792458E8;
|
---|
[d7d2da3] | 122 |
|
---|
| 123 | fItInputArray->Reset();
|
---|
[0ed696f] | 124 |
|
---|
| 125 | // --- Deal with Primary vertex first ------
|
---|
| 126 |
|
---|
| 127 | fFunction->GetRandom2(dz, dt);
|
---|
| 128 |
|
---|
| 129 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
| 130 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
| 131 |
|
---|
[d7d2da3] | 132 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
| 133 | {
|
---|
[00078bc] | 134 | z = candidate->Position.Z();
|
---|
| 135 | t = candidate->Position.T();
|
---|
[fb8932a] | 136 | candidate->Position.SetZ(z + dz);
|
---|
| 137 | candidate->Position.SetT(t + dt);
|
---|
[d07e957] | 138 | fParticleOutputArray->Add(candidate);
|
---|
[d7d2da3] | 139 | }
|
---|
| 140 |
|
---|
| 141 | factory = GetFactory();
|
---|
[0ed696f] | 142 |
|
---|
[22dc7fd] | 143 | vertexcandidate = factory->NewCandidate();
|
---|
| 144 | vertexcandidate->Position.SetXYZT(0.0, 0.0, dz, dt);
|
---|
| 145 | fVertexOutputArray->Add(vertexcandidate);
|
---|
[0ed696f] | 146 |
|
---|
[22dc7fd] | 147 | // --- Then with pile-up vertices ------
|
---|
[0ed696f] | 148 |
|
---|
[76d3973] | 149 | switch(fPileUpDistribution)
|
---|
| 150 | {
|
---|
| 151 | case 0:
|
---|
| 152 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
| 153 | break;
|
---|
| 154 | case 1:
|
---|
| 155 | numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
|
---|
| 156 | break;
|
---|
| 157 | default:
|
---|
| 158 | numberOfEvents = gRandom->Poisson(fMeanPileUp);
|
---|
| 159 | break;
|
---|
| 160 | }
|
---|
[d7d2da3] | 161 |
|
---|
| 162 | allEntries = fReader->GetEntries();
|
---|
[0ed696f] | 163 |
|
---|
[76d3973] | 164 | for(event = 0; event < numberOfEvents; ++event)
|
---|
[d7d2da3] | 165 | {
|
---|
| 166 | do
|
---|
| 167 | {
|
---|
| 168 | entry = TMath::Nint(gRandom->Rndm()*allEntries);
|
---|
| 169 | }
|
---|
| 170 | while(entry >= allEntries);
|
---|
| 171 |
|
---|
| 172 | fReader->ReadEntry(entry);
|
---|
[0ed696f] | 173 |
|
---|
| 174 | // --- Pile-up vertex smearing
|
---|
| 175 |
|
---|
| 176 | fFunction->GetRandom2(dz, dt);
|
---|
| 177 |
|
---|
| 178 | dt *= c_light*1.0E3; // necessary in order to make t in mm/c
|
---|
| 179 | dz *= 1.0E3; // necessary in order to make z in mm
|
---|
| 180 |
|
---|
[b68a60f] | 181 | dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
|
---|
[d7d2da3] | 182 |
|
---|
[22dc7fd] | 183 | vertexcandidate = factory->NewCandidate();
|
---|
| 184 | vertexcandidate->Position.SetXYZT(0.0, 0.0, dz, dt);
|
---|
| 185 | vertexcandidate->IsPU = 1;
|
---|
[0ed696f] | 186 |
|
---|
[22dc7fd] | 187 | fVertexOutputArray->Add(vertexcandidate);
|
---|
[d07e957] | 188 |
|
---|
[d7d2da3] | 189 | while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
|
---|
| 190 | {
|
---|
| 191 | candidate = factory->NewCandidate();
|
---|
| 192 |
|
---|
| 193 | candidate->PID = pid;
|
---|
| 194 |
|
---|
| 195 | candidate->Status = 1;
|
---|
| 196 |
|
---|
| 197 | pdgParticle = pdg->GetParticle(pid);
|
---|
| 198 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
|
---|
| 199 | candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
|
---|
| 200 |
|
---|
| 201 | candidate->IsPU = 1;
|
---|
[0ed696f] | 202 |
|
---|
[d7d2da3] | 203 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
[b68a60f] | 204 | candidate->Momentum.RotateZ(dphi);
|
---|
| 205 |
|
---|
[00078bc] | 206 | candidate->Position.SetXYZT(x, y, z + dz, t + dt);
|
---|
[b68a60f] | 207 | candidate->Position.RotateZ(dphi);
|
---|
[0ed696f] | 208 |
|
---|
[d07e957] | 209 | fParticleOutputArray->Add(candidate);
|
---|
[d7d2da3] | 210 | }
|
---|
| 211 | }
|
---|
| 212 | }
|
---|
| 213 |
|
---|
| 214 | //------------------------------------------------------------------------------
|
---|
| 215 |
|
---|