Fork me on GitHub

source: svn/trunk/modules/PileUpMerger.cc@ 1390

Last change on this file since 1390 was 1352, checked in by Pavel Demin, 11 years ago

fix formatting in PileUpMerger

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