1 | /*
|
---|
2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
4 | *
|
---|
5 | * This program is free software: you can redistribute it and/or modify
|
---|
6 | * it under the terms of the GNU General Public License as published by
|
---|
7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
8 | * (at your option) any later version.
|
---|
9 | *
|
---|
10 | * This program is distributed in the hope that it will be useful,
|
---|
11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
13 | * GNU General Public License for more details.
|
---|
14 | *
|
---|
15 | * You should have received a copy of the GNU General Public License
|
---|
16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
17 | */
|
---|
18 |
|
---|
19 | /** \class TreeWriter
|
---|
20 | *
|
---|
21 | * Fills ROOT tree branches.
|
---|
22 | *
|
---|
23 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
24 | *
|
---|
25 | */
|
---|
26 |
|
---|
27 | #include "modules/TreeWriter.h"
|
---|
28 |
|
---|
29 | #include "classes/DelphesClasses.h"
|
---|
30 | #include "classes/DelphesFactory.h"
|
---|
31 | #include "classes/DelphesFormula.h"
|
---|
32 |
|
---|
33 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
34 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
35 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
36 | #include "ExRootAnalysis/ExRootTreeBranch.h"
|
---|
37 |
|
---|
38 | #include "TDatabasePDG.h"
|
---|
39 | #include "TFormula.h"
|
---|
40 | #include "TLorentzVector.h"
|
---|
41 | #include "TMath.h"
|
---|
42 | #include "TObjArray.h"
|
---|
43 | #include "TROOT.h"
|
---|
44 | #include "TRandom3.h"
|
---|
45 | #include "TString.h"
|
---|
46 |
|
---|
47 | #include <algorithm>
|
---|
48 | #include <iostream>
|
---|
49 | #include <sstream>
|
---|
50 | #include <stdexcept>
|
---|
51 |
|
---|
52 | using namespace std;
|
---|
53 |
|
---|
54 | //------------------------------------------------------------------------------
|
---|
55 |
|
---|
56 | TreeWriter::TreeWriter()
|
---|
57 | {
|
---|
58 | }
|
---|
59 |
|
---|
60 | //------------------------------------------------------------------------------
|
---|
61 |
|
---|
62 | TreeWriter::~TreeWriter()
|
---|
63 | {
|
---|
64 | }
|
---|
65 |
|
---|
66 | //------------------------------------------------------------------------------
|
---|
67 |
|
---|
68 | void TreeWriter::Init()
|
---|
69 | {
|
---|
70 | fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
|
---|
71 | fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
|
---|
72 | fClassMap[Track::Class()] = &TreeWriter::ProcessTracks;
|
---|
73 | fClassMap[Tower::Class()] = &TreeWriter::ProcessTowers;
|
---|
74 | fClassMap[ParticleFlowCandidate::Class()] = &TreeWriter::ProcessParticleFlowCandidates;
|
---|
75 | fClassMap[Photon::Class()] = &TreeWriter::ProcessPhotons;
|
---|
76 | fClassMap[Electron::Class()] = &TreeWriter::ProcessElectrons;
|
---|
77 | fClassMap[Muon::Class()] = &TreeWriter::ProcessMuons;
|
---|
78 | fClassMap[CscCluster::Class()] = &TreeWriter::ProcessCscCluster;
|
---|
79 | fClassMap[Jet::Class()] = &TreeWriter::ProcessJets;
|
---|
80 | fClassMap[MissingET::Class()] = &TreeWriter::ProcessMissingET;
|
---|
81 | fClassMap[ScalarHT::Class()] = &TreeWriter::ProcessScalarHT;
|
---|
82 | fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
|
---|
83 | fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
|
---|
84 | fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit;
|
---|
85 |
|
---|
86 | TBranchMap::iterator itBranchMap;
|
---|
87 | map<TClass *, TProcessMethod>::iterator itClassMap;
|
---|
88 |
|
---|
89 | // read branch configuration and
|
---|
90 | // import array with output from filter/classifier/jetfinder modules
|
---|
91 |
|
---|
92 | ExRootConfParam param = GetParam("Branch");
|
---|
93 | Long_t i, size;
|
---|
94 | TString branchName, branchClassName, branchInputArray;
|
---|
95 | TClass *branchClass;
|
---|
96 | TObjArray *array;
|
---|
97 | ExRootTreeBranch *branch;
|
---|
98 |
|
---|
99 | size = param.GetSize();
|
---|
100 | for(i = 0; i < size / 3; ++i)
|
---|
101 | {
|
---|
102 | branchInputArray = param[i * 3].GetString();
|
---|
103 | branchName = param[i * 3 + 1].GetString();
|
---|
104 | branchClassName = param[i * 3 + 2].GetString();
|
---|
105 |
|
---|
106 | branchClass = gROOT->GetClass(branchClassName);
|
---|
107 |
|
---|
108 | if(!branchClass)
|
---|
109 | {
|
---|
110 | cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
|
---|
111 | continue;
|
---|
112 | }
|
---|
113 |
|
---|
114 | itClassMap = fClassMap.find(branchClass);
|
---|
115 | if(itClassMap == fClassMap.end())
|
---|
116 | {
|
---|
117 | cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
|
---|
118 | continue;
|
---|
119 | }
|
---|
120 |
|
---|
121 | array = ImportArray(branchInputArray);
|
---|
122 | branch = NewBranch(branchName, branchClass);
|
---|
123 |
|
---|
124 | fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array)));
|
---|
125 | }
|
---|
126 |
|
---|
127 | param = GetParam("Info");
|
---|
128 | TString infoName;
|
---|
129 | Double_t infoValue;
|
---|
130 |
|
---|
131 | size = param.GetSize();
|
---|
132 | for(i = 0; i < size / 2; ++i)
|
---|
133 | {
|
---|
134 | infoName = param[i * 2].GetString();
|
---|
135 | infoValue = param[i * 2 + 1].GetDouble();
|
---|
136 |
|
---|
137 | AddInfo(infoName, infoValue);
|
---|
138 | }
|
---|
139 | }
|
---|
140 |
|
---|
141 | //------------------------------------------------------------------------------
|
---|
142 |
|
---|
143 | void TreeWriter::Finish()
|
---|
144 | {
|
---|
145 | }
|
---|
146 |
|
---|
147 | //------------------------------------------------------------------------------
|
---|
148 |
|
---|
149 | void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
|
---|
150 | {
|
---|
151 | TIter it1(candidate->GetCandidates());
|
---|
152 | it1.Reset();
|
---|
153 | array->Clear();
|
---|
154 |
|
---|
155 | while((candidate = static_cast<Candidate *>(it1.Next())))
|
---|
156 | {
|
---|
157 | TIter it2(candidate->GetCandidates());
|
---|
158 |
|
---|
159 | // particle
|
---|
160 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
161 | {
|
---|
162 | array->Add(candidate);
|
---|
163 | continue;
|
---|
164 | }
|
---|
165 |
|
---|
166 | // track
|
---|
167 | candidate = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
168 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
169 | {
|
---|
170 | array->Add(candidate);
|
---|
171 | continue;
|
---|
172 | }
|
---|
173 |
|
---|
174 | // tower
|
---|
175 | it2.Reset();
|
---|
176 | while((candidate = static_cast<Candidate *>(it2.Next())))
|
---|
177 | {
|
---|
178 | array->Add(candidate->GetCandidates()->At(0));
|
---|
179 | }
|
---|
180 | }
|
---|
181 | }
|
---|
182 |
|
---|
183 | //------------------------------------------------------------------------------
|
---|
184 |
|
---|
185 | void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
|
---|
186 | {
|
---|
187 | TIter iterator(array);
|
---|
188 | Candidate *candidate = 0;
|
---|
189 | GenParticle *entry = 0;
|
---|
190 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
191 |
|
---|
192 | const Double_t c_light = 2.99792458E8;
|
---|
193 |
|
---|
194 | // loop over all particles
|
---|
195 | iterator.Reset();
|
---|
196 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
197 | {
|
---|
198 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
199 | const TLorentzVector &position = candidate->Position;
|
---|
200 | const TLorentzVector &DecayPosition = candidate->DecayPosition;
|
---|
201 |
|
---|
202 | entry = static_cast<GenParticle *>(branch->NewEntry());
|
---|
203 |
|
---|
204 | entry->SetBit(kIsReferenced);
|
---|
205 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
206 |
|
---|
207 | pt = momentum.Pt();
|
---|
208 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
209 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
210 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
211 | rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
|
---|
212 |
|
---|
213 | entry->PID = candidate->PID;
|
---|
214 |
|
---|
215 | entry->Status = candidate->Status;
|
---|
216 | entry->IsPU = candidate->IsPU;
|
---|
217 |
|
---|
218 | entry->M1 = candidate->M1;
|
---|
219 | entry->M2 = candidate->M2;
|
---|
220 |
|
---|
221 | entry->D1 = candidate->D1;
|
---|
222 | entry->D2 = candidate->D2;
|
---|
223 |
|
---|
224 | entry->Charge = candidate->Charge;
|
---|
225 | entry->Mass = candidate->Mass;
|
---|
226 |
|
---|
227 | entry->E = momentum.E();
|
---|
228 | entry->Px = momentum.Px();
|
---|
229 | entry->Py = momentum.Py();
|
---|
230 | entry->Pz = momentum.Pz();
|
---|
231 |
|
---|
232 | entry->Eta = eta;
|
---|
233 | entry->Phi = momentum.Phi();
|
---|
234 | entry->PT = pt;
|
---|
235 |
|
---|
236 | entry->Rapidity = rapidity;
|
---|
237 |
|
---|
238 | entry->X = position.X();
|
---|
239 | entry->Y = position.Y();
|
---|
240 | entry->Z = position.Z();
|
---|
241 | entry->T = position.T() * 1.0E-3 / c_light;
|
---|
242 |
|
---|
243 | entry->decayX = DecayPosition.X();
|
---|
244 | entry->decayY = DecayPosition.Y();
|
---|
245 | entry->decayZ = DecayPosition.Z();
|
---|
246 | entry->decayT = DecayPosition.T()* 1.0E-3 / c_light;
|
---|
247 | }
|
---|
248 | }
|
---|
249 |
|
---|
250 | //------------------------------------------------------------------------------
|
---|
251 |
|
---|
252 | void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
|
---|
253 | {
|
---|
254 | TIter iterator(array);
|
---|
255 | Candidate *candidate = 0, *constituent = 0;
|
---|
256 | Vertex *entry = 0;
|
---|
257 |
|
---|
258 | const Double_t c_light = 2.99792458E8;
|
---|
259 |
|
---|
260 | Double_t x, y, z, t, xError, yError, zError, tError, sigma, sumPT2, btvSumPT2, genDeltaZ, genSumPT2;
|
---|
261 | UInt_t index, ndf;
|
---|
262 |
|
---|
263 | CompBase *compare = Candidate::fgCompare;
|
---|
264 | Candidate::fgCompare = CompSumPT2<Candidate>::Instance();
|
---|
265 | array->Sort();
|
---|
266 | Candidate::fgCompare = compare;
|
---|
267 |
|
---|
268 | // loop over all vertices
|
---|
269 | iterator.Reset();
|
---|
270 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
271 | {
|
---|
272 |
|
---|
273 | index = candidate->ClusterIndex;
|
---|
274 | ndf = candidate->ClusterNDF;
|
---|
275 | sigma = candidate->ClusterSigma;
|
---|
276 | sumPT2 = candidate->SumPT2;
|
---|
277 | btvSumPT2 = candidate->BTVSumPT2;
|
---|
278 | genDeltaZ = candidate->GenDeltaZ;
|
---|
279 | genSumPT2 = candidate->GenSumPT2;
|
---|
280 |
|
---|
281 | x = candidate->Position.X();
|
---|
282 | y = candidate->Position.Y();
|
---|
283 | z = candidate->Position.Z();
|
---|
284 | t = candidate->Position.T() * 1.0E-3 / c_light;
|
---|
285 |
|
---|
286 | xError = candidate->PositionError.X();
|
---|
287 | yError = candidate->PositionError.Y();
|
---|
288 | zError = candidate->PositionError.Z();
|
---|
289 | tError = candidate->PositionError.T() * 1.0E-3 / c_light;
|
---|
290 |
|
---|
291 | entry = static_cast<Vertex *>(branch->NewEntry());
|
---|
292 |
|
---|
293 | entry->Index = index;
|
---|
294 | entry->NDF = ndf;
|
---|
295 | entry->Sigma = sigma;
|
---|
296 | entry->SumPT2 = sumPT2;
|
---|
297 | entry->BTVSumPT2 = btvSumPT2;
|
---|
298 | entry->GenDeltaZ = genDeltaZ;
|
---|
299 | entry->GenSumPT2 = genSumPT2;
|
---|
300 |
|
---|
301 | entry->X = x;
|
---|
302 | entry->Y = y;
|
---|
303 | entry->Z = z;
|
---|
304 | entry->T = t;
|
---|
305 |
|
---|
306 | entry->ErrorX = xError;
|
---|
307 | entry->ErrorY = yError;
|
---|
308 | entry->ErrorZ = zError;
|
---|
309 | entry->ErrorT = tError;
|
---|
310 |
|
---|
311 | TIter itConstituents(candidate->GetCandidates());
|
---|
312 | itConstituents.Reset();
|
---|
313 | entry->Constituents.Clear();
|
---|
314 | while((constituent = static_cast<Candidate *>(itConstituents.Next())))
|
---|
315 | {
|
---|
316 | entry->Constituents.Add(constituent);
|
---|
317 | }
|
---|
318 | }
|
---|
319 | }
|
---|
320 |
|
---|
321 | //------------------------------------------------------------------------------
|
---|
322 |
|
---|
323 | void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
|
---|
324 | {
|
---|
325 | TIter iterator(array);
|
---|
326 | Candidate *candidate = 0;
|
---|
327 | Candidate *particle = 0;
|
---|
328 | Track *entry = 0;
|
---|
329 | Double_t pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi, m;
|
---|
330 | const Double_t c_light = 2.99792458E8;
|
---|
331 |
|
---|
332 | // loop over all tracks
|
---|
333 | iterator.Reset();
|
---|
334 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
335 | {
|
---|
336 | const TLorentzVector &position = candidate->Position;
|
---|
337 |
|
---|
338 | cosTheta = TMath::Abs(position.CosTheta());
|
---|
339 | signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
340 | eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
|
---|
341 | rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
|
---|
342 |
|
---|
343 | entry = static_cast<Track *>(branch->NewEntry());
|
---|
344 |
|
---|
345 | entry->SetBit(kIsReferenced);
|
---|
346 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
347 |
|
---|
348 | entry->PID = candidate->PID;
|
---|
349 |
|
---|
350 | entry->Charge = candidate->Charge;
|
---|
351 |
|
---|
352 | entry->EtaOuter = eta;
|
---|
353 | entry->PhiOuter = position.Phi();
|
---|
354 |
|
---|
355 | entry->XOuter = position.X();
|
---|
356 | entry->YOuter = position.Y();
|
---|
357 | entry->ZOuter = position.Z();
|
---|
358 | entry->TOuter = position.T() * 1.0E-3 / c_light;
|
---|
359 |
|
---|
360 | entry->L = candidate->L;
|
---|
361 |
|
---|
362 | entry->D0 = candidate->D0;
|
---|
363 | entry->DZ = candidate->DZ;
|
---|
364 | entry->Nclusters = candidate->Nclusters;
|
---|
365 | entry->dNdx = candidate->dNdx;
|
---|
366 |
|
---|
367 | entry->ErrorP = candidate->ErrorP;
|
---|
368 | entry->ErrorPT = candidate->ErrorPT;
|
---|
369 |
|
---|
370 | // diagonal covariance matrix terms
|
---|
371 | entry->ErrorD0 = candidate->ErrorD0;
|
---|
372 | entry->ErrorC = candidate->ErrorC;
|
---|
373 | entry->ErrorPhi = candidate->ErrorPhi;
|
---|
374 | entry->ErrorDZ = candidate->ErrorDZ;
|
---|
375 | entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
|
---|
376 |
|
---|
377 | // add some offdiagonal covariance matrix elements
|
---|
378 | entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
|
---|
379 | entry->ErrorD0C = candidate->TrackCovariance(0,2);
|
---|
380 | entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
|
---|
381 | entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
|
---|
382 | entry->ErrorPhiC = candidate->TrackCovariance(1,2);
|
---|
383 | entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
|
---|
384 | entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
|
---|
385 | entry->ErrorCDZ = candidate->TrackCovariance(2,3);
|
---|
386 | entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
|
---|
387 | entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
|
---|
388 |
|
---|
389 | entry->Xd = candidate->Xd;
|
---|
390 | entry->Yd = candidate->Yd;
|
---|
391 | entry->Zd = candidate->Zd;
|
---|
392 |
|
---|
393 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
394 |
|
---|
395 | pt = momentum.Pt();
|
---|
396 | p = momentum.P();
|
---|
397 | phi = momentum.Phi();
|
---|
398 | m = momentum.M();
|
---|
399 | ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
|
---|
400 |
|
---|
401 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
402 | signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
403 | eta = (cosTheta == 1.0 ? signz * 999.9 : momentum.Eta());
|
---|
404 | rapidity = (cosTheta == 1.0 ? signz * 999.9 : momentum.Rapidity());
|
---|
405 |
|
---|
406 | entry->P = p;
|
---|
407 | entry->PT = pt;
|
---|
408 | entry->Eta = eta;
|
---|
409 | entry->Phi = phi;
|
---|
410 | entry->CtgTheta = ctgTheta;
|
---|
411 | entry->C = candidate->C;
|
---|
412 | entry->Mass = m;
|
---|
413 |
|
---|
414 | particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
415 | //const TLorentzVector &initialPosition = particle->Position;
|
---|
416 | const TLorentzVector &initialPosition = candidate->InitialPosition;
|
---|
417 |
|
---|
418 | entry->X = initialPosition.X();
|
---|
419 | entry->Y = initialPosition.Y();
|
---|
420 | entry->Z = initialPosition.Z();
|
---|
421 | entry->T = initialPosition.T() * 1.0E-3 / c_light;
|
---|
422 | entry->ErrorT =candidate-> ErrorT * 1.0E-3 / c_light;
|
---|
423 |
|
---|
424 | entry->Particle = particle;
|
---|
425 |
|
---|
426 | entry->VertexIndex = candidate->ClusterIndex;
|
---|
427 | }
|
---|
428 | }
|
---|
429 |
|
---|
430 | //------------------------------------------------------------------------------
|
---|
431 |
|
---|
432 | void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
|
---|
433 | {
|
---|
434 | TIter iterator(array);
|
---|
435 | Candidate *candidate = 0;
|
---|
436 | Tower *entry = 0;
|
---|
437 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
438 | const Double_t c_light = 2.99792458E8;
|
---|
439 |
|
---|
440 | // loop over all towers
|
---|
441 | iterator.Reset();
|
---|
442 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
443 | {
|
---|
444 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
445 | const TLorentzVector &position = candidate->Position;
|
---|
446 |
|
---|
447 | pt = momentum.Pt();
|
---|
448 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
449 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
450 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
451 | rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
|
---|
452 |
|
---|
453 | entry = static_cast<Tower *>(branch->NewEntry());
|
---|
454 |
|
---|
455 | entry->SetBit(kIsReferenced);
|
---|
456 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
457 |
|
---|
458 | entry->Eta = eta;
|
---|
459 | entry->Phi = momentum.Phi();
|
---|
460 | entry->ET = pt;
|
---|
461 | entry->E = momentum.E();
|
---|
462 | entry->Eem = candidate->Eem;
|
---|
463 | entry->Ehad = candidate->Ehad;
|
---|
464 | entry->Etrk = candidate->Etrk;
|
---|
465 | entry->Edges[0] = candidate->Edges[0];
|
---|
466 | entry->Edges[1] = candidate->Edges[1];
|
---|
467 | entry->Edges[2] = candidate->Edges[2];
|
---|
468 | entry->Edges[3] = candidate->Edges[3];
|
---|
469 |
|
---|
470 | entry->T = position.T() * 1.0E-3 / c_light;
|
---|
471 | entry->NTimeHits = candidate->NTimeHits;
|
---|
472 |
|
---|
473 | FillParticles(candidate, &entry->Particles);
|
---|
474 | }
|
---|
475 | }
|
---|
476 |
|
---|
477 | //------------------------------------------------------------------------------
|
---|
478 |
|
---|
479 | void TreeWriter::ProcessParticleFlowCandidates(ExRootTreeBranch *branch, TObjArray *array)
|
---|
480 | {
|
---|
481 |
|
---|
482 | TIter iterator(array);
|
---|
483 | Candidate *candidate = 0;
|
---|
484 | Candidate *particle = 0;
|
---|
485 | ParticleFlowCandidate *entry = 0;
|
---|
486 | Double_t e, pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi, m;
|
---|
487 | const Double_t c_light = 2.99792458E8;
|
---|
488 |
|
---|
489 | // loop over all tracks
|
---|
490 | iterator.Reset();
|
---|
491 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
492 | {
|
---|
493 | const TLorentzVector &position = candidate->Position;
|
---|
494 |
|
---|
495 | cosTheta = TMath::Abs(position.CosTheta());
|
---|
496 | signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
497 | eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
|
---|
498 | rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
|
---|
499 |
|
---|
500 | entry = static_cast<ParticleFlowCandidate *>(branch->NewEntry());
|
---|
501 |
|
---|
502 | entry->SetBit(kIsReferenced);
|
---|
503 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
504 |
|
---|
505 | entry->PID = candidate->PID;
|
---|
506 |
|
---|
507 | entry->Charge = candidate->Charge;
|
---|
508 |
|
---|
509 | entry->EtaOuter = eta;
|
---|
510 | entry->PhiOuter = position.Phi();
|
---|
511 |
|
---|
512 | entry->XOuter = position.X();
|
---|
513 | entry->YOuter = position.Y();
|
---|
514 | entry->ZOuter = position.Z();
|
---|
515 | entry->TOuter = position.T() * 1.0E-3 / c_light;
|
---|
516 |
|
---|
517 | entry->L = candidate->L;
|
---|
518 |
|
---|
519 | entry->D0 = candidate->D0;
|
---|
520 | entry->DZ = candidate->DZ;
|
---|
521 | entry->Nclusters = candidate->Nclusters;
|
---|
522 | entry->dNdx = candidate->dNdx;
|
---|
523 |
|
---|
524 | entry->ErrorP = candidate->ErrorP;
|
---|
525 | entry->ErrorPT = candidate->ErrorPT;
|
---|
526 | entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
|
---|
527 |
|
---|
528 |
|
---|
529 | // diagonal covariance matrix terms
|
---|
530 |
|
---|
531 | entry->ErrorD0 = candidate->ErrorD0;
|
---|
532 | entry->ErrorC = candidate->ErrorC;
|
---|
533 | entry->ErrorPhi = candidate->ErrorPhi;
|
---|
534 | entry->ErrorDZ = candidate->ErrorDZ;
|
---|
535 | entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
|
---|
536 |
|
---|
537 | // add some offdiagonal covariance matrix elements
|
---|
538 | entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
|
---|
539 | entry->ErrorD0C = candidate->TrackCovariance(0,2);
|
---|
540 | entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
|
---|
541 | entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
|
---|
542 | entry->ErrorPhiC = candidate->TrackCovariance(1,2);
|
---|
543 | entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
|
---|
544 | entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
|
---|
545 | entry->ErrorCDZ = candidate->TrackCovariance(2,3);
|
---|
546 | entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
|
---|
547 | entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
|
---|
548 |
|
---|
549 | entry->Xd = candidate->Xd;
|
---|
550 | entry->Yd = candidate->Yd;
|
---|
551 | entry->Zd = candidate->Zd;
|
---|
552 |
|
---|
553 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
554 |
|
---|
555 | e = momentum.E();
|
---|
556 | pt = momentum.Pt();
|
---|
557 | p = momentum.P();
|
---|
558 | phi = momentum.Phi();
|
---|
559 | m = momentum.M();
|
---|
560 | ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
|
---|
561 |
|
---|
562 | entry->E = e;
|
---|
563 | entry->P = p;
|
---|
564 | entry->PT = pt;
|
---|
565 | entry->Eta = eta;
|
---|
566 | entry->Phi = phi;
|
---|
567 | entry->CtgTheta = ctgTheta;
|
---|
568 | entry->C = candidate->C;
|
---|
569 | entry->Mass = m;
|
---|
570 |
|
---|
571 | particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
572 | //const TLorentzVector &initialPosition = particle->Position;
|
---|
573 | const TLorentzVector &initialPosition = candidate->InitialPosition;
|
---|
574 |
|
---|
575 | entry->X = initialPosition.X();
|
---|
576 | entry->Y = initialPosition.Y();
|
---|
577 | entry->Z = initialPosition.Z();
|
---|
578 | entry->T = initialPosition.T() * 1.0E-3 / c_light;
|
---|
579 | entry->ErrorT = candidate-> ErrorT * 1.0E-3 / c_light;
|
---|
580 |
|
---|
581 | entry->VertexIndex = candidate->ClusterIndex;
|
---|
582 |
|
---|
583 | entry->Eem = candidate->Eem;
|
---|
584 | entry->Ehad = candidate->Ehad;
|
---|
585 | entry->Etrk = candidate->Etrk;
|
---|
586 | entry->Edges[0] = candidate->Edges[0];
|
---|
587 | entry->Edges[1] = candidate->Edges[1];
|
---|
588 | entry->Edges[2] = candidate->Edges[2];
|
---|
589 | entry->Edges[3] = candidate->Edges[3];
|
---|
590 |
|
---|
591 | //entry->T = position.T() * 1.0E-3 / c_light;
|
---|
592 | entry->NTimeHits = candidate->NTimeHits;
|
---|
593 |
|
---|
594 | FillParticles(candidate, &entry->Particles);
|
---|
595 |
|
---|
596 | }
|
---|
597 | }
|
---|
598 |
|
---|
599 | //------------------------------------------------------------------------------
|
---|
600 |
|
---|
601 | void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
602 | {
|
---|
603 | TIter iterator(array);
|
---|
604 | Candidate *candidate = 0;
|
---|
605 | Photon *entry = 0;
|
---|
606 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
607 | const Double_t c_light = 2.99792458E8;
|
---|
608 |
|
---|
609 | array->Sort();
|
---|
610 |
|
---|
611 | // loop over all photons
|
---|
612 | iterator.Reset();
|
---|
613 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
614 | {
|
---|
615 | TIter it1(candidate->GetCandidates());
|
---|
616 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
617 | const TLorentzVector &position = candidate->Position;
|
---|
618 |
|
---|
619 | pt = momentum.Pt();
|
---|
620 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
621 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
622 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
623 | rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
|
---|
624 |
|
---|
625 | entry = static_cast<Photon *>(branch->NewEntry());
|
---|
626 |
|
---|
627 | entry->Eta = eta;
|
---|
628 | entry->Phi = momentum.Phi();
|
---|
629 | entry->PT = pt;
|
---|
630 | entry->E = momentum.E();
|
---|
631 | entry->T = position.T() * 1.0E-3 / c_light;
|
---|
632 |
|
---|
633 | // Isolation variables
|
---|
634 |
|
---|
635 | entry->IsolationVar = candidate->IsolationVar;
|
---|
636 | entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
|
---|
637 | entry->SumPtCharged = candidate->SumPtCharged;
|
---|
638 | entry->SumPtNeutral = candidate->SumPtNeutral;
|
---|
639 | entry->SumPtChargedPU = candidate->SumPtChargedPU;
|
---|
640 | entry->SumPt = candidate->SumPt;
|
---|
641 |
|
---|
642 | entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad / candidate->Eem : 999.9;
|
---|
643 |
|
---|
644 | // 1: prompt -- 2: non prompt -- 3: fake
|
---|
645 | entry->Status = candidate->Status;
|
---|
646 |
|
---|
647 | FillParticles(candidate, &entry->Particles);
|
---|
648 | }
|
---|
649 | }
|
---|
650 |
|
---|
651 | //------------------------------------------------------------------------------
|
---|
652 |
|
---|
653 | void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
654 | {
|
---|
655 | TIter iterator(array);
|
---|
656 | Candidate *candidate = 0;
|
---|
657 | Electron *entry = 0;
|
---|
658 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
659 | const Double_t c_light = 2.99792458E8;
|
---|
660 |
|
---|
661 | array->Sort();
|
---|
662 |
|
---|
663 | // loop over all electrons
|
---|
664 | iterator.Reset();
|
---|
665 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
666 | {
|
---|
667 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
668 | const TLorentzVector &position = candidate->Position;
|
---|
669 |
|
---|
670 | pt = momentum.Pt();
|
---|
671 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
672 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
673 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
674 | rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
|
---|
675 |
|
---|
676 | entry = static_cast<Electron *>(branch->NewEntry());
|
---|
677 |
|
---|
678 | entry->Eta = eta;
|
---|
679 | entry->Phi = momentum.Phi();
|
---|
680 | entry->PT = pt;
|
---|
681 |
|
---|
682 | entry->T = position.T() * 1.0E-3 / c_light;
|
---|
683 |
|
---|
684 | // displacement
|
---|
685 | entry->D0 = candidate->D0;
|
---|
686 | entry->ErrorD0 = candidate->ErrorD0;
|
---|
687 | entry->DZ = candidate->DZ;
|
---|
688 | entry->ErrorDZ = candidate->ErrorDZ;
|
---|
689 |
|
---|
690 | // Isolation variables
|
---|
691 | entry->IsolationVar = candidate->IsolationVar;
|
---|
692 | entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
|
---|
693 | entry->SumPtCharged = candidate->SumPtCharged;
|
---|
694 | entry->SumPtNeutral = candidate->SumPtNeutral;
|
---|
695 | entry->SumPtChargedPU = candidate->SumPtChargedPU;
|
---|
696 | entry->SumPt = candidate->SumPt;
|
---|
697 |
|
---|
698 | entry->Charge = candidate->Charge;
|
---|
699 |
|
---|
700 | entry->EhadOverEem = 0.0;
|
---|
701 |
|
---|
702 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
703 | }
|
---|
704 | }
|
---|
705 |
|
---|
706 | //------------------------------------------------------------------------------
|
---|
707 |
|
---|
708 | void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
709 | {
|
---|
710 | TIter iterator(array);
|
---|
711 | Candidate *candidate = 0;
|
---|
712 | Muon *entry = 0;
|
---|
713 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
714 |
|
---|
715 | const Double_t c_light = 2.99792458E8;
|
---|
716 |
|
---|
717 | array->Sort();
|
---|
718 |
|
---|
719 | // loop over all muons
|
---|
720 | iterator.Reset();
|
---|
721 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
722 | {
|
---|
723 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
724 | const TLorentzVector &position = candidate->Position;
|
---|
725 |
|
---|
726 | pt = momentum.Pt();
|
---|
727 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
728 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
729 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
730 | rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
|
---|
731 |
|
---|
732 | entry = static_cast<Muon *>(branch->NewEntry());
|
---|
733 |
|
---|
734 | entry->SetBit(kIsReferenced);
|
---|
735 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
736 |
|
---|
737 | entry->Eta = eta;
|
---|
738 | entry->Phi = momentum.Phi();
|
---|
739 | entry->PT = pt;
|
---|
740 |
|
---|
741 | entry->T = position.T() * 1.0E-3 / c_light;
|
---|
742 |
|
---|
743 | // displacement
|
---|
744 | entry->D0 = candidate->D0;
|
---|
745 | entry->ErrorD0 = candidate->ErrorD0;
|
---|
746 | entry->DZ = candidate->DZ;
|
---|
747 | entry->ErrorDZ = candidate->ErrorDZ;
|
---|
748 |
|
---|
749 | // Isolation variables
|
---|
750 |
|
---|
751 | entry->IsolationVar = candidate->IsolationVar;
|
---|
752 | entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
|
---|
753 | entry->SumPtCharged = candidate->SumPtCharged;
|
---|
754 | entry->SumPtNeutral = candidate->SumPtNeutral;
|
---|
755 | entry->SumPtChargedPU = candidate->SumPtChargedPU;
|
---|
756 | entry->SumPt = candidate->SumPt;
|
---|
757 |
|
---|
758 | entry->Charge = candidate->Charge;
|
---|
759 |
|
---|
760 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
761 | }
|
---|
762 | }
|
---|
763 |
|
---|
764 | //------------------------------------------------------------------------------
|
---|
765 |
|
---|
766 | void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
|
---|
767 | {
|
---|
768 | TIter iterator(array);
|
---|
769 | Candidate *candidate = 0, *constituent = 0;
|
---|
770 | Jet *entry = 0;
|
---|
771 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
772 | Double_t ecalEnergy, hcalEnergy;
|
---|
773 | const Double_t c_light = 2.99792458E8;
|
---|
774 | Int_t i;
|
---|
775 |
|
---|
776 | array->Sort();
|
---|
777 |
|
---|
778 | // loop over all jets
|
---|
779 | iterator.Reset();
|
---|
780 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
781 | {
|
---|
782 | TIter itConstituents(candidate->GetCandidates());
|
---|
783 |
|
---|
784 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
785 | const TLorentzVector &position = candidate->Position;
|
---|
786 |
|
---|
787 | pt = momentum.Pt();
|
---|
788 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
789 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
790 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
791 | rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
|
---|
792 |
|
---|
793 | entry = static_cast<Jet *>(branch->NewEntry());
|
---|
794 |
|
---|
795 | entry->Eta = eta;
|
---|
796 | entry->Phi = momentum.Phi();
|
---|
797 | entry->PT = pt;
|
---|
798 |
|
---|
799 | entry->T = position.T() * 1.0E-3 / c_light;
|
---|
800 |
|
---|
801 | entry->Mass = momentum.M();
|
---|
802 |
|
---|
803 | entry->Area = candidate->Area;
|
---|
804 |
|
---|
805 | entry->DeltaEta = candidate->DeltaEta;
|
---|
806 | entry->DeltaPhi = candidate->DeltaPhi;
|
---|
807 |
|
---|
808 | entry->Flavor = candidate->Flavor;
|
---|
809 | entry->FlavorAlgo = candidate->FlavorAlgo;
|
---|
810 | entry->FlavorPhys = candidate->FlavorPhys;
|
---|
811 |
|
---|
812 | entry->BTag = candidate->BTag;
|
---|
813 |
|
---|
814 | entry->BTagAlgo = candidate->BTagAlgo;
|
---|
815 | entry->BTagPhys = candidate->BTagPhys;
|
---|
816 |
|
---|
817 | entry->TauTag = candidate->TauTag;
|
---|
818 | entry->TauWeight = candidate->TauWeight;
|
---|
819 |
|
---|
820 | entry->Charge = candidate->Charge;
|
---|
821 |
|
---|
822 | itConstituents.Reset();
|
---|
823 | entry->Constituents.Clear();
|
---|
824 | ecalEnergy = 0.0;
|
---|
825 | hcalEnergy = 0.0;
|
---|
826 | while((constituent = static_cast<Candidate *>(itConstituents.Next())))
|
---|
827 | {
|
---|
828 | entry->Constituents.Add(constituent);
|
---|
829 | ecalEnergy += constituent->Eem;
|
---|
830 | hcalEnergy += constituent->Ehad;
|
---|
831 | }
|
---|
832 |
|
---|
833 | entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy / ecalEnergy : 999.9;
|
---|
834 |
|
---|
835 | //--- Pile-Up Jet ID variables ----
|
---|
836 |
|
---|
837 | entry->NCharged = candidate->NCharged;
|
---|
838 | entry->NNeutrals = candidate->NNeutrals;
|
---|
839 |
|
---|
840 | entry->NeutralEnergyFraction = candidate->NeutralEnergyFraction;
|
---|
841 | entry->ChargedEnergyFraction = candidate->ChargedEnergyFraction;
|
---|
842 | entry->Beta = candidate->Beta;
|
---|
843 | entry->BetaStar = candidate->BetaStar;
|
---|
844 | entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
|
---|
845 | entry->PTD = candidate->PTD;
|
---|
846 |
|
---|
847 | //--- Sub-structure variables ----
|
---|
848 |
|
---|
849 | entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
|
---|
850 | entry->NSubJetsPruned = candidate->NSubJetsPruned;
|
---|
851 | entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
|
---|
852 |
|
---|
853 | entry->SoftDroppedJet = candidate->SoftDroppedJet;
|
---|
854 | entry->SoftDroppedSubJet1 = candidate->SoftDroppedSubJet1;
|
---|
855 | entry->SoftDroppedSubJet2 = candidate->SoftDroppedSubJet2;
|
---|
856 |
|
---|
857 | for(i = 0; i < 5; i++)
|
---|
858 | {
|
---|
859 | entry->FracPt[i] = candidate->FracPt[i];
|
---|
860 | entry->Tau[i] = candidate->Tau[i];
|
---|
861 | entry->TrimmedP4[i] = candidate->TrimmedP4[i];
|
---|
862 | entry->PrunedP4[i] = candidate->PrunedP4[i];
|
---|
863 | entry->SoftDroppedP4[i] = candidate->SoftDroppedP4[i];
|
---|
864 | }
|
---|
865 |
|
---|
866 | //--- exclusive clustering variables ---
|
---|
867 | entry->ExclYmerge23 = candidate->ExclYmerge23;
|
---|
868 | entry->ExclYmerge34 = candidate->ExclYmerge34;
|
---|
869 | entry->ExclYmerge45 = candidate->ExclYmerge45;
|
---|
870 | entry->ExclYmerge56 = candidate->ExclYmerge56;
|
---|
871 |
|
---|
872 | FillParticles(candidate, &entry->Particles);
|
---|
873 | }
|
---|
874 | }
|
---|
875 |
|
---|
876 | //------------------------------------------------------------------------------
|
---|
877 |
|
---|
878 | void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
|
---|
879 | {
|
---|
880 | Candidate *candidate = 0;
|
---|
881 | MissingET *entry = 0;
|
---|
882 |
|
---|
883 | // get the first entry
|
---|
884 | if((candidate = static_cast<Candidate *>(array->At(0))))
|
---|
885 | {
|
---|
886 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
887 |
|
---|
888 | entry = static_cast<MissingET *>(branch->NewEntry());
|
---|
889 |
|
---|
890 | entry->Eta = (-momentum).Eta();
|
---|
891 | entry->Phi = (-momentum).Phi();
|
---|
892 | entry->MET = momentum.Pt();
|
---|
893 | }
|
---|
894 | }
|
---|
895 | //------------------------------------------------------------------------------
|
---|
896 |
|
---|
897 | void TreeWriter::ProcessCscCluster(ExRootTreeBranch *branch, TObjArray *array)
|
---|
898 | {
|
---|
899 | TIter iterator(array);
|
---|
900 | Candidate *candidate = 0;
|
---|
901 | CscCluster *entry = 0;
|
---|
902 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
903 |
|
---|
904 | const Double_t c_light = 2.99792458E8; // in unit of m/s
|
---|
905 |
|
---|
906 | array->Sort();
|
---|
907 |
|
---|
908 |
|
---|
909 | // loop over all clusters
|
---|
910 | iterator.Reset();
|
---|
911 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
912 | {
|
---|
913 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
914 | const TLorentzVector &position = candidate->DecayPosition;
|
---|
915 |
|
---|
916 | pt = momentum.Pt();
|
---|
917 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
918 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
919 | eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
|
---|
920 |
|
---|
921 | entry = static_cast<CscCluster *>(branch->NewEntry());
|
---|
922 |
|
---|
923 | entry->SetBit(kIsReferenced);
|
---|
924 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
925 |
|
---|
926 | entry->Eta = eta;
|
---|
927 | entry->Phi = momentum.Phi();
|
---|
928 |
|
---|
929 | entry->PT = momentum.Pt(); // pt of LLP
|
---|
930 | entry->Px = momentum.Px();// px of LLP
|
---|
931 | entry->Py = momentum.Py();// py of LLP
|
---|
932 | entry->Pz = momentum.Pz();// pz of LLP
|
---|
933 | entry->E = momentum.E(); // E of LLP
|
---|
934 | entry->pid = candidate->PID; // LLP pid
|
---|
935 | entry->Eem = candidate->Eem; // LLP Eem
|
---|
936 | entry->Ehad = candidate->Ehad; // LLP Ehad
|
---|
937 | Double_t beta = momentum.P()/momentum.E();
|
---|
938 | Double_t gamma = 1.0/sqrt(1-beta*beta);
|
---|
939 | Double_t decayDistance = sqrt(pow(position.X(),2)+pow(position.Y(),2)+pow(position.Z(),2)); // mm
|
---|
940 | entry->beta = beta; // LLP pid
|
---|
941 | entry->ctau = decayDistance/(beta * gamma); // LLP travel time in its rest frame
|
---|
942 | entry->T = decayDistance*(1./beta-1)* 1.0E-3/c_light*1e9; // ns
|
---|
943 | entry->X = position.X(); // LLP decay x
|
---|
944 | entry->Y = position.Y(); // LLP decay y
|
---|
945 | entry->Z = position.Z(); // LLP decay z
|
---|
946 | }
|
---|
947 | }
|
---|
948 |
|
---|
949 | //------------------------------------------------------------------------------
|
---|
950 |
|
---|
951 | void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
|
---|
952 | {
|
---|
953 | Candidate *candidate = 0;
|
---|
954 | ScalarHT *entry = 0;
|
---|
955 |
|
---|
956 | // get the first entry
|
---|
957 | if((candidate = static_cast<Candidate *>(array->At(0))))
|
---|
958 | {
|
---|
959 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
960 |
|
---|
961 | entry = static_cast<ScalarHT *>(branch->NewEntry());
|
---|
962 |
|
---|
963 | entry->HT = momentum.Pt();
|
---|
964 | }
|
---|
965 | }
|
---|
966 |
|
---|
967 | //------------------------------------------------------------------------------
|
---|
968 |
|
---|
969 | void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
|
---|
970 | {
|
---|
971 | TIter iterator(array);
|
---|
972 | Candidate *candidate = 0;
|
---|
973 | Rho *entry = 0;
|
---|
974 |
|
---|
975 | // loop over all rho
|
---|
976 | iterator.Reset();
|
---|
977 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
978 | {
|
---|
979 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
980 |
|
---|
981 | entry = static_cast<Rho *>(branch->NewEntry());
|
---|
982 |
|
---|
983 | entry->Rho = momentum.E();
|
---|
984 | entry->Edges[0] = candidate->Edges[0];
|
---|
985 | entry->Edges[1] = candidate->Edges[1];
|
---|
986 | }
|
---|
987 | }
|
---|
988 |
|
---|
989 | //------------------------------------------------------------------------------
|
---|
990 |
|
---|
991 | void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
|
---|
992 | {
|
---|
993 | Candidate *candidate = 0;
|
---|
994 | Weight *entry = 0;
|
---|
995 |
|
---|
996 | // get the first entry
|
---|
997 | if((candidate = static_cast<Candidate *>(array->At(0))))
|
---|
998 | {
|
---|
999 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
1000 |
|
---|
1001 | entry = static_cast<Weight *>(branch->NewEntry());
|
---|
1002 |
|
---|
1003 | entry->Weight = momentum.E();
|
---|
1004 | }
|
---|
1005 | }
|
---|
1006 |
|
---|
1007 | //------------------------------------------------------------------------------
|
---|
1008 |
|
---|
1009 | void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
|
---|
1010 | {
|
---|
1011 | TIter iterator(array);
|
---|
1012 | Candidate *candidate = 0;
|
---|
1013 | HectorHit *entry = 0;
|
---|
1014 |
|
---|
1015 | // loop over all roman pot hits
|
---|
1016 | iterator.Reset();
|
---|
1017 | while((candidate = static_cast<Candidate *>(iterator.Next())))
|
---|
1018 | {
|
---|
1019 | const TLorentzVector &position = candidate->Position;
|
---|
1020 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
1021 |
|
---|
1022 | entry = static_cast<HectorHit *>(branch->NewEntry());
|
---|
1023 |
|
---|
1024 | entry->E = momentum.E();
|
---|
1025 |
|
---|
1026 | entry->Tx = momentum.Px();
|
---|
1027 | entry->Ty = momentum.Py();
|
---|
1028 |
|
---|
1029 | entry->T = position.T();
|
---|
1030 |
|
---|
1031 | entry->X = position.X();
|
---|
1032 | entry->Y = position.Y();
|
---|
1033 | entry->S = position.Z();
|
---|
1034 |
|
---|
1035 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
1036 | }
|
---|
1037 | }
|
---|
1038 |
|
---|
1039 | //------------------------------------------------------------------------------
|
---|
1040 |
|
---|
1041 | void TreeWriter::Process()
|
---|
1042 | {
|
---|
1043 | TBranchMap::iterator itBranchMap;
|
---|
1044 | ExRootTreeBranch *branch;
|
---|
1045 | TProcessMethod method;
|
---|
1046 | TObjArray *array;
|
---|
1047 |
|
---|
1048 | for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
|
---|
1049 | {
|
---|
1050 | branch = itBranchMap->first;
|
---|
1051 | method = itBranchMap->second.first;
|
---|
1052 | array = itBranchMap->second.second;
|
---|
1053 |
|
---|
1054 | (this->*method)(branch, array);
|
---|
1055 | }
|
---|
1056 | }
|
---|
1057 |
|
---|
1058 | //------------------------------------------------------------------------------
|
---|