Fork me on GitHub

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

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

add pileup rotation

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