source: trunk/modules/MadGraphMatchingTreeWriter.cc@ 16

Last change on this file since 16 was 16, checked in by Pavel Demin, 16 years ago

add include for gROOT

File size: 5.9 KB
Line 
1
2#include "modules/MadGraphMatchingTreeWriter.h"
3
4
5#include "ExRootAnalysis/ExRootResult.h"
6#include "ExRootAnalysis/ExRootClasses.h"
7#include "ExRootAnalysis/ExRootTreeBranch.h"
8
9#include "ExRootAnalysis/ExRootCandidate.h"
10
11#include "TROOT.h"
12#include "TClonesArray.h"
13
14#include "TH1.h"
15#include "TH2.h"
16#include "TString.h"
17#include "TCanvas.h"
18#include "TLorentzVector.h"
19
20#include <iostream>
21
22using namespace std;
23
24//------------------------------------------------------------------------------
25
26MadGraphMatchingTreeWriter::MadGraphMatchingTreeWriter()
27{
28}
29
30//------------------------------------------------------------------------------
31
32MadGraphMatchingTreeWriter::~MadGraphMatchingTreeWriter()
33{
34}
35
36//------------------------------------------------------------------------------
37
38void MadGraphMatchingTreeWriter::Init()
39{
40 fJetPTMin = GetDouble("JetPTMin", 20.0);
41 fJetEtaMax = GetDouble("JetEtaMax", 4.5);
42
43 fClassMap[ExRootGenParticle::Class()] = &MadGraphMatchingTreeWriter::ProcessPartons;
44
45 fClassMap[ExRootMatching::Class()] = &MadGraphMatchingTreeWriter::ProcessMatching;
46
47 fClassMap[ExRootGenJet::Class()] = &MadGraphMatchingTreeWriter::ProcessJets;
48
49 TBranchMap::iterator itBranchMap;
50 map< TClass *, TProcessMethod >::iterator itClassMap;
51
52 // read branch configuration and
53 // import array with output from filter/classifier/jetfinder modules
54
55 ExRootConfParam param = GetParam("Branch");
56 Long_t i, size;
57 TString branchName, branchClassName, branchInputArray;
58 TClass *branchClass;
59 const TObjArray *array;
60 ExRootTreeBranch *branch;
61
62 size = param.GetSize();
63 for(i = 0; i < size; ++i)
64 {
65 branchName = param[i][0].GetString();
66 branchClassName = param[i][1].GetString();
67 branchInputArray = param[i][2].GetString();
68
69 branchClass = gROOT->GetClass(branchClassName);
70
71 if(!branchClass)
72 {
73 cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
74 continue;
75 }
76
77 itClassMap = fClassMap.find(branchClass);
78 if(itClassMap == fClassMap.end())
79 {
80 cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
81 continue;
82 }
83
84 array = ImportArray(branchInputArray);
85 branch = NewBranch(branchName, branchClass);
86
87 fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array->MakeIterator())));
88 }
89
90}
91
92//------------------------------------------------------------------------------
93
94void MadGraphMatchingTreeWriter::Finish()
95{
96 TBranchMap::iterator itBranchMap;
97 TIterator *iterator;
98
99 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
100 {
101 iterator = itBranchMap->second.second;
102 if(iterator) delete iterator;
103 }
104}
105
106//------------------------------------------------------------------------------
107
108void MadGraphMatchingTreeWriter::ProcessPartons(ExRootTreeBranch *branch, TIterator *iterator)
109{
110 ExRootCandidate *candidate = 0;
111 ExRootGenParticle *entry = 0;
112 Double_t pt, signPz, eta, rapidity;
113
114 // loop over all partons
115 iterator->Reset();
116 while((candidate = static_cast<ExRootCandidate*>(iterator->Next())))
117 {
118 const TLorentzVector &momentum = candidate->GetP4();
119
120 entry = static_cast<ExRootGenParticle*>(branch->NewEntry());
121
122 pt = momentum.Pt();
123 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
124 eta = (pt == 0.0 ? signPz*999.9 : momentum.Eta());
125 rapidity = (pt == 0.0 ? signPz*999.9 : momentum.Rapidity());
126
127 entry->PID = candidate->GetType()->PdgCode();
128
129 entry->E = momentum.E();
130 entry->Px = momentum.Px();
131 entry->Py = momentum.Py();
132 entry->Pz = momentum.Pz();
133
134 entry->Eta = eta;
135 entry->Phi = momentum.Phi();
136 entry->PT = pt;
137
138 entry->Rapidity = rapidity;
139 }
140}
141
142//------------------------------------------------------------------------------
143
144void MadGraphMatchingTreeWriter::ProcessMatching(ExRootTreeBranch *branch, TIterator *iterator)
145{
146 ExRootMatching *matching = 0, *entry = 0;
147
148 // loop over all matching
149 iterator->Reset();
150 while((matching = static_cast<ExRootMatching*>(iterator->Next())))
151 {
152 entry = static_cast<ExRootMatching*>(branch->NewEntry());
153
154 entry->DMerge = matching->DMerge;
155 entry->YMerge = matching->YMerge;
156 }
157}
158
159//------------------------------------------------------------------------------
160
161void MadGraphMatchingTreeWriter::ProcessJets(ExRootTreeBranch *branch, TIterator *iterator)
162{
163 ExRootCandidate *candidate = 0;
164 ExRootGenJet *entry = 0;
165 Double_t pt, signPz, eta, rapidity;
166
167 // loop over all jets
168 iterator->Reset();
169 while((candidate = static_cast<ExRootCandidate*>(iterator->Next())))
170 {
171 const TLorentzVector &momentum = candidate->GetP4();
172
173 pt = momentum.Pt();
174 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
175 eta = (pt == 0.0 ? signPz*999.9 : momentum.Eta());
176 rapidity = (pt == 0.0 ? signPz*999.9 : momentum.Rapidity());
177
178 if(pt < fJetPTMin) continue;
179 if(TMath::Abs(eta) > fJetEtaMax) continue;
180
181 entry = static_cast<ExRootGenJet*>(branch->NewEntry());
182
183 entry->E = momentum.E();
184 entry->Px = momentum.Px();
185 entry->Py = momentum.Py();
186 entry->Pz = momentum.Pz();
187
188 entry->Eta = eta;
189 entry->Phi = momentum.Phi();
190 entry->PT = pt;
191
192 entry->Rapidity = rapidity;
193
194 entry->Mass = momentum.M();
195 }
196}
197
198//------------------------------------------------------------------------------
199
200void MadGraphMatchingTreeWriter::Process()
201{
202
203 TBranchMap::iterator itBranchMap;
204 ExRootTreeBranch *branch;
205 TProcessMethod method;
206 TIterator *iterator;
207
208 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
209 {
210 branch = itBranchMap->first;
211 method = itBranchMap->second.first;
212 iterator = itBranchMap->second.second;
213
214 (this->*method)(branch, iterator);
215 }
216
217}
218
219//------------------------------------------------------------------------------
220
Note: See TracBrowser for help on using the repository browser.