Fork me on GitHub

source: git/modules/PileUpMerger.cc@ fb8932a

ImprovedOutputFile Timing dual_readout llp
Last change on this file since fb8932a was fb8932a, checked in by pavel <pavel@…>, 11 years ago

fixed candidate position in PileUpMerger

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