Fork me on GitHub

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

Last change on this file since 1121 was 1111, checked in by Pavel Demin, 12 years ago

add pileup rotation

File size: 3.4 KB
Line 
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"
18#include "classes/DelphesFormula.h"
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() :
43 fReader(0), fItInputArray(0)
44{
45}
46
47//------------------------------------------------------------------------------
48
49PileUpMerger::~PileUpMerger()
50{
51}
52
53//------------------------------------------------------------------------------
54
55void PileUpMerger::Init()
56{
57 const char *fileName;
58
59 fMeanPileUp = GetDouble("MeanPileUp", 10);
60 fZVertexSpread = GetDouble("ZVertexSpread", 0.05)*1.0E3;
61
62 fileName = GetString("PileUpFile", "MinBias.pileup");
63 fReader = new DelphesPileUpReader(fileName);
64
65 // import input array
66 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
67 fItInputArray = fInputArray->MakeIterator();
68
69 // create output arrays
70 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
71}
72
73//------------------------------------------------------------------------------
74
75void PileUpMerger::Finish()
76{
77 if(fReader) delete fReader;
78}
79
80//------------------------------------------------------------------------------
81
82void PileUpMerger::Process()
83{
84 TDatabasePDG *pdg = TDatabasePDG::Instance();
85 TParticlePDG *pdgParticle;
86 Int_t pid;
87 Float_t x, y, z, t;
88 Float_t px, py, pz, e;
89 Double_t dz, dphi;
90 Int_t poisson, event;
91 Long64_t allEntries, entry;
92 Candidate *candidate;
93 DelphesFactory *factory;
94
95 fItInputArray->Reset();
96 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
97 {
98 fOutputArray->Add(candidate);
99 }
100
101 factory = GetFactory();
102
103 poisson = gRandom->Poisson(fMeanPileUp);
104
105 allEntries = fReader->GetEntries();
106
107 for(event = 0; event < poisson; ++event)
108 {
109 do
110 {
111 entry = TMath::Nint(gRandom->Rndm()*allEntries);
112 }
113 while(entry >= allEntries);
114
115 fReader->ReadEntry(entry);
116
117 dz = gRandom->Gaus(0.0, fZVertexSpread);
118 dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
119
120 while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e))
121 {
122 candidate = factory->NewCandidate();
123
124 candidate->PID = pid;
125
126 candidate->Status = 1;
127
128 pdgParticle = pdg->GetParticle(pid);
129 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
130 candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
131
132 candidate->IsPU = 1;
133
134 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
135 candidate->Momentum.RotateZ(dphi);
136
137 candidate->Position.SetXYZT(x, y, z + dz, t);
138 candidate->Position.RotateZ(dphi);
139
140 fOutputArray->Add(candidate);
141 }
142 }
143}
144
145//------------------------------------------------------------------------------
146
Note: See TracBrowser for help on using the repository browser.