1 |
|
---|
2 | #include "modules/MadGraphAnalysis.h"
|
---|
3 |
|
---|
4 |
|
---|
5 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
6 | #include "ExRootAnalysis/ExRootClasses.h"
|
---|
7 |
|
---|
8 | #include "ExRootAnalysis/ExRootCandidate.h"
|
---|
9 |
|
---|
10 | #include "TClonesArray.h"
|
---|
11 |
|
---|
12 | #include "TH1.h"
|
---|
13 | #include "TH2.h"
|
---|
14 | #include "TString.h"
|
---|
15 | #include "TCanvas.h"
|
---|
16 | #include "TLorentzVector.h"
|
---|
17 |
|
---|
18 | #include <iostream>
|
---|
19 |
|
---|
20 | using namespace std;
|
---|
21 |
|
---|
22 | //------------------------------------------------------------------------------
|
---|
23 |
|
---|
24 | MadGraphAnalysis::MadGraphAnalysis()
|
---|
25 | {
|
---|
26 | }
|
---|
27 |
|
---|
28 | //------------------------------------------------------------------------------
|
---|
29 |
|
---|
30 | MadGraphAnalysis::~MadGraphAnalysis()
|
---|
31 | {
|
---|
32 | }
|
---|
33 |
|
---|
34 | //------------------------------------------------------------------------------
|
---|
35 |
|
---|
36 | void MadGraphAnalysis::Init()
|
---|
37 | {
|
---|
38 | fOutputFileName = GetString("OutputFile", "pythia_plots.root");
|
---|
39 |
|
---|
40 | // import array with output from filter/classifier module
|
---|
41 |
|
---|
42 | fInputArray = ImportArray(GetString("InputArray", "merger/candidates"));
|
---|
43 |
|
---|
44 | fIsUnWeighted = GetBool("IsUnWeighted", kFALSE);
|
---|
45 |
|
---|
46 | // import ROOT tree branch
|
---|
47 |
|
---|
48 | if(fIsUnWeighted)
|
---|
49 | {
|
---|
50 | fBranchEvent = UseBranch("Event");
|
---|
51 | }
|
---|
52 | else
|
---|
53 | {
|
---|
54 | fBranchEvent = 0;
|
---|
55 | }
|
---|
56 |
|
---|
57 | }
|
---|
58 |
|
---|
59 | //------------------------------------------------------------------------------
|
---|
60 |
|
---|
61 | void MadGraphAnalysis::Finish()
|
---|
62 | {
|
---|
63 | GetPlots()->Write(fOutputFileName);
|
---|
64 |
|
---|
65 | GetPlots()->GetCanvas()->SetLogy(1);
|
---|
66 | GetPlots()->Print();
|
---|
67 | GetPlots()->GetCanvas()->SetLogy(0);
|
---|
68 | }
|
---|
69 |
|
---|
70 | //------------------------------------------------------------------------------
|
---|
71 |
|
---|
72 | void MadGraphAnalysis::Process()
|
---|
73 | {
|
---|
74 | ExRootCandidate *candidate1 = 0, *candidate2 = 0;
|
---|
75 | ParticleHistograms *histogramsParticle = 0;
|
---|
76 | PairHistograms *histogramsPair = 0;
|
---|
77 | Int_t maxEntry, entry1, entry2;
|
---|
78 |
|
---|
79 | Double_t weight = 1.0;
|
---|
80 | Double_t pt1, pt2, dr, rapidity, signPz;
|
---|
81 |
|
---|
82 | ExRootLHEFEvent *eventInfo = 0;
|
---|
83 | if(fIsUnWeighted && fBranchEvent && fBranchEvent->GetEntriesFast() == 1)
|
---|
84 | {
|
---|
85 | eventInfo = static_cast<ExRootLHEFEvent*>(fBranchEvent->At(0));
|
---|
86 |
|
---|
87 | weight = eventInfo->Weight;
|
---|
88 | }
|
---|
89 |
|
---|
90 | // fill histograms
|
---|
91 | maxEntry = fInputArray->GetEntriesFast();
|
---|
92 | for(entry1 = 0; entry1 < maxEntry; ++entry1)
|
---|
93 | {
|
---|
94 | candidate1 = static_cast<ExRootCandidate*>(fInputArray->At(entry1));
|
---|
95 |
|
---|
96 | const TLorentzVector &vector1 = candidate1->GetP4();
|
---|
97 |
|
---|
98 | pt1 = vector1.Pt();
|
---|
99 | signPz = (vector1.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
100 | rapidity = (pt1 == 0.0 ? signPz*999.9 : vector1.Rapidity());
|
---|
101 |
|
---|
102 | // fill histograms for single particles
|
---|
103 |
|
---|
104 | histogramsParticle = GetParticleHistograms(candidate1->GetName());
|
---|
105 |
|
---|
106 | histogramsParticle->fParticlePt->Fill(pt1, weight);
|
---|
107 |
|
---|
108 | histogramsParticle->fParticleRapidity->Fill(rapidity, weight);
|
---|
109 |
|
---|
110 | // skip pairs of resonances
|
---|
111 | if(candidate1->IsResonance()) continue;
|
---|
112 |
|
---|
113 | for(entry2 = entry1 + 1; entry2 < maxEntry; ++entry2)
|
---|
114 | {
|
---|
115 | candidate2 = static_cast<ExRootCandidate*>(fInputArray->At(entry2));
|
---|
116 |
|
---|
117 | // skip pairs of resonances
|
---|
118 | if(candidate2->IsResonance()) continue;
|
---|
119 |
|
---|
120 | const TLorentzVector &vector2 = candidate2->GetP4();
|
---|
121 |
|
---|
122 | pt2 = vector2.Pt();
|
---|
123 | dr = (pt1 == 0.0 || pt2 == 0.0 ? 999.9 : vector1.DeltaR(vector2));
|
---|
124 |
|
---|
125 | // fill histograms for pairs of particles
|
---|
126 |
|
---|
127 | histogramsPair = GetPairHistograms(candidate1->GetName(), candidate2->GetName());
|
---|
128 |
|
---|
129 | histogramsPair->fPairDeltaR->Fill(dr, weight);
|
---|
130 | histogramsPair->fPairMass->Fill((vector1 + vector2).M(), weight);
|
---|
131 | }
|
---|
132 | }
|
---|
133 | }
|
---|
134 |
|
---|
135 | //------------------------------------------------------------------------------
|
---|
136 |
|
---|
137 | void MadGraphAnalysis::BookParticleHistograms(MadGraphAnalysis::ParticleHistograms *histograms,
|
---|
138 | const char *name, const char *title)
|
---|
139 | {
|
---|
140 | ExRootResult *result = GetPlots();
|
---|
141 | histograms->fParticlePt = result->AddHist1D(Form("pt_%s", name),
|
---|
142 | Form("P_{T}(%s)", title),
|
---|
143 | Form("P_{T}(%s), GeV/c", title),
|
---|
144 | "pb/bin",
|
---|
145 | 60, 0.0, 300.0, 0, 1);
|
---|
146 | histograms->fParticlePt->SetStats(kTRUE);
|
---|
147 |
|
---|
148 | histograms->fParticleRapidity = result->AddHist1D(Form("y_%s", name),
|
---|
149 | Form("y(%s)", title),
|
---|
150 | Form("y(%s)", title),
|
---|
151 | "pb/bin",
|
---|
152 | 100, -5.0, 5.0, 0, 0);
|
---|
153 | histograms->fParticleRapidity->SetStats(kTRUE);
|
---|
154 |
|
---|
155 | }
|
---|
156 |
|
---|
157 | //------------------------------------------------------------------------------
|
---|
158 |
|
---|
159 | void MadGraphAnalysis::BookPairHistograms(MadGraphAnalysis::PairHistograms *histograms,
|
---|
160 | const char *name, const char *title)
|
---|
161 | {
|
---|
162 | ExRootResult *result = GetPlots();
|
---|
163 | histograms->fPairDeltaR = result->AddHist1D(Form("dr_%s", name),
|
---|
164 | Form("#DeltaR(%s)", title),
|
---|
165 | Form("#DeltaR(%s)", title),
|
---|
166 | "pb/bin",
|
---|
167 | 70, 0.0, 7.0, 0, 1);
|
---|
168 | histograms->fPairDeltaR->SetStats(kTRUE);
|
---|
169 |
|
---|
170 | histograms->fPairMass = result->AddHist1D(Form("mass_%s", name),
|
---|
171 | Form("M_{inv}(%s)", title),
|
---|
172 | Form("M_{inv}(%s), GeV/c^{2}", title),
|
---|
173 | "pb/bin",
|
---|
174 | 120, 0.0, 600.0, 0, 1);
|
---|
175 | histograms->fPairMass->SetStats(kTRUE);
|
---|
176 |
|
---|
177 | }
|
---|
178 |
|
---|
179 | //------------------------------------------------------------------------------
|
---|
180 |
|
---|
181 | MadGraphAnalysis::ParticleHistograms *
|
---|
182 | MadGraphAnalysis::GetParticleHistograms(const char *candName)
|
---|
183 | {
|
---|
184 | map<TString, ParticleHistograms *>::iterator itParticleHistogramsMap;
|
---|
185 | ParticleHistograms *histograms = 0;
|
---|
186 | TString name = Form("%s", candName);
|
---|
187 | name.ReplaceAll("{", "");
|
---|
188 | name.ReplaceAll("}", "");
|
---|
189 | name.ReplaceAll("^", "");
|
---|
190 | name.ReplaceAll("#bar", "anti_");
|
---|
191 | name.ReplaceAll("#", "");
|
---|
192 | TString title = Form("%s", candName);
|
---|
193 | itParticleHistogramsMap = fParticleHistogramsMap.find(name);
|
---|
194 | if(itParticleHistogramsMap == fParticleHistogramsMap.end())
|
---|
195 | {
|
---|
196 | histograms = new ParticleHistograms;
|
---|
197 |
|
---|
198 | BookParticleHistograms(histograms, name, title);
|
---|
199 |
|
---|
200 | fParticleHistogramsMap[name] = histograms;
|
---|
201 | }
|
---|
202 | else
|
---|
203 | {
|
---|
204 | histograms = itParticleHistogramsMap->second;
|
---|
205 | }
|
---|
206 | return histograms;
|
---|
207 | }
|
---|
208 |
|
---|
209 | //------------------------------------------------------------------------------
|
---|
210 |
|
---|
211 | MadGraphAnalysis::PairHistograms *
|
---|
212 | MadGraphAnalysis::GetPairHistograms(const char *candName1, const char *candName2)
|
---|
213 | {
|
---|
214 | map<TString, PairHistograms *>::iterator itPairHistogramsMap;
|
---|
215 | PairHistograms *histograms = 0;
|
---|
216 | TString name = Form("%s_%s", candName1, candName2);
|
---|
217 | name.ReplaceAll("{", "");
|
---|
218 | name.ReplaceAll("}", "");
|
---|
219 | name.ReplaceAll("^", "");
|
---|
220 | name.ReplaceAll("#bar", "anti_");
|
---|
221 | name.ReplaceAll("#", "");
|
---|
222 | TString title = Form("%s, %s", candName1, candName2);
|
---|
223 | itPairHistogramsMap = fPairHistogramsMap.find(name);
|
---|
224 | if(itPairHistogramsMap == fPairHistogramsMap.end())
|
---|
225 | {
|
---|
226 | histograms = new PairHistograms;
|
---|
227 |
|
---|
228 | BookPairHistograms(histograms, name, title);
|
---|
229 |
|
---|
230 | fPairHistogramsMap[name] = histograms;
|
---|
231 | }
|
---|
232 | else
|
---|
233 | {
|
---|
234 | histograms = itPairHistogramsMap->second;
|
---|
235 | }
|
---|
236 | return histograms;
|
---|
237 | }
|
---|
238 |
|
---|