Fork me on GitHub

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

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

add vertex block

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