Fork me on GitHub

source: git/modules/PileUpMergerPythia8.cc@ a0b6d15

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a0b6d15 was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

delete duplicate license file and prepend GPLv3 header to all source code files

  • Property mode set to 100644
File size: 4.4 KB
Line 
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 *
24 * $Date$
25 * $Revision$
26 *
27 *
28 * \author P. Selvaggi - UCL, Louvain-la-Neuve
29 *
30 */
31
32#include "modules/PileUpMergerPythia8.h"
33
34#include "classes/DelphesClasses.h"
35#include "classes/DelphesFactory.h"
36#include "classes/DelphesFormula.h"
37#include "classes/DelphesPileUpReader.h"
38
39#include "ExRootAnalysis/ExRootResult.h"
40#include "ExRootAnalysis/ExRootFilter.h"
41#include "ExRootAnalysis/ExRootClassifier.h"
42
43#include "Pythia8/Pythia.h"
44
45#include "TMath.h"
46#include "TString.h"
47#include "TFormula.h"
48#include "TRandom3.h"
49#include "TObjArray.h"
50#include "TDatabasePDG.h"
51#include "TLorentzVector.h"
52
53#include <algorithm>
54#include <stdexcept>
55#include <iostream>
56#include <sstream>
57
58using namespace std;
59
60//------------------------------------------------------------------------------
61
62PileUpMergerPythia8::PileUpMergerPythia8() :
63 fPythia(0), fItInputArray(0)
64{
65}
66
67//------------------------------------------------------------------------------
68
69PileUpMergerPythia8::~PileUpMergerPythia8()
70{
71}
72
73//------------------------------------------------------------------------------
74
75void PileUpMergerPythia8::Init()
76{
77 const char *fileName;
78
79 fMeanPileUp = GetDouble("MeanPileUp", 10);
80 fZVertexSpread = GetDouble("ZVertexSpread", 0.05)*1.0E3;
81
82 fPTMin = GetDouble("PTMin", 0.0);
83
84 fileName = GetString("ConfigFile", "MinBias.cmnd");
85 fPythia = new Pythia8::Pythia();
86 fPythia->readFile(fileName);
87
88 // import input array
89 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
90 fItInputArray = fInputArray->MakeIterator();
91
92 // create output arrays
93 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
94}
95
96//------------------------------------------------------------------------------
97
98void PileUpMergerPythia8::Finish()
99{
100 if(fPythia) delete fPythia;
101}
102
103//------------------------------------------------------------------------------
104
105void PileUpMergerPythia8::Process()
106{
107 TDatabasePDG *pdg = TDatabasePDG::Instance();
108 TParticlePDG *pdgParticle;
109 Int_t pid;
110 Float_t x, y, z, t;
111 Float_t px, py, pz, e;
112 Double_t dz, dphi;
113 Int_t poisson, event, i;
114 Candidate *candidate;
115 DelphesFactory *factory;
116
117 fItInputArray->Reset();
118 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
119 {
120 fOutputArray->Add(candidate);
121 }
122
123 factory = GetFactory();
124
125 poisson = gRandom->Poisson(fMeanPileUp);
126
127 for(event = 0; event < poisson; ++event)
128 {
129 while(!fPythia->next());
130
131 dz = gRandom->Gaus(0.0, fZVertexSpread);
132 dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
133
134 for(i = 0; i < fPythia->event.size(); ++i)
135 {
136 Pythia8::Particle &particle = fPythia->event[i];
137
138 if(particle.status() != 1 || !particle.isVisible() || particle.pT() <= fPTMin) continue;
139
140 pid = particle.id();
141 px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.e();
142 x = particle.xProd(); y = particle.yProd(); z = particle.zProd(); t = particle.tProd();
143
144 candidate = factory->NewCandidate();
145
146 candidate->PID = pid;
147
148 candidate->Status = 1;
149
150 pdgParticle = pdg->GetParticle(pid);
151 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
152 candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
153
154 candidate->IsPU = 1;
155
156 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
157 candidate->Momentum.RotateZ(dphi);
158
159 candidate->Position.SetXYZT(x, y, z + dz, t);
160 candidate->Position.RotateZ(dphi);
161
162 fOutputArray->Add(candidate);
163 }
164 }
165}
166
167//------------------------------------------------------------------------------
168
Note: See TracBrowser for help on using the repository browser.