[d7d2da3] | 1 |
|
---|
| 2 | /** \class TreeWriter
|
---|
| 3 | *
|
---|
| 4 | * Fills ROOT tree branches.
|
---|
| 5 | *
|
---|
| 6 | * $Date$
|
---|
| 7 | * $Revision$
|
---|
| 8 | *
|
---|
| 9 | *
|
---|
| 10 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 11 | *
|
---|
| 12 | */
|
---|
| 13 |
|
---|
| 14 | #include "modules/TreeWriter.h"
|
---|
| 15 |
|
---|
| 16 | #include "classes/DelphesClasses.h"
|
---|
| 17 | #include "classes/DelphesFactory.h"
|
---|
| 18 | #include "classes/DelphesFormula.h"
|
---|
| 19 |
|
---|
| 20 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 21 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 22 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 23 | #include "ExRootAnalysis/ExRootTreeBranch.h"
|
---|
| 24 |
|
---|
| 25 | #include "TROOT.h"
|
---|
| 26 | #include "TMath.h"
|
---|
| 27 | #include "TString.h"
|
---|
| 28 | #include "TFormula.h"
|
---|
| 29 | #include "TRandom3.h"
|
---|
| 30 | #include "TObjArray.h"
|
---|
| 31 | #include "TDatabasePDG.h"
|
---|
| 32 | #include "TLorentzVector.h"
|
---|
| 33 |
|
---|
| 34 | #include <algorithm>
|
---|
| 35 | #include <stdexcept>
|
---|
| 36 | #include <iostream>
|
---|
| 37 | #include <sstream>
|
---|
| 38 |
|
---|
| 39 | using namespace std;
|
---|
| 40 |
|
---|
| 41 | //------------------------------------------------------------------------------
|
---|
| 42 |
|
---|
| 43 | TreeWriter::TreeWriter()
|
---|
| 44 | {
|
---|
| 45 | }
|
---|
| 46 |
|
---|
| 47 | //------------------------------------------------------------------------------
|
---|
| 48 |
|
---|
| 49 | TreeWriter::~TreeWriter()
|
---|
| 50 | {
|
---|
| 51 | }
|
---|
| 52 |
|
---|
| 53 | //------------------------------------------------------------------------------
|
---|
| 54 |
|
---|
| 55 | void TreeWriter::Init()
|
---|
| 56 | {
|
---|
| 57 | fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
|
---|
[2de4dcd] | 58 | fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
|
---|
[d7d2da3] | 59 | fClassMap[Track::Class()] = &TreeWriter::ProcessTracks;
|
---|
| 60 | fClassMap[Tower::Class()] = &TreeWriter::ProcessTowers;
|
---|
| 61 | fClassMap[Photon::Class()] = &TreeWriter::ProcessPhotons;
|
---|
| 62 | fClassMap[Electron::Class()] = &TreeWriter::ProcessElectrons;
|
---|
| 63 | fClassMap[Muon::Class()] = &TreeWriter::ProcessMuons;
|
---|
| 64 | fClassMap[Jet::Class()] = &TreeWriter::ProcessJets;
|
---|
| 65 | fClassMap[MissingET::Class()] = &TreeWriter::ProcessMissingET;
|
---|
| 66 | fClassMap[ScalarHT::Class()] = &TreeWriter::ProcessScalarHT;
|
---|
[71648c2] | 67 | fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
|
---|
[2e229c9] | 68 | fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
|
---|
[8f7db23] | 69 | fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit;
|
---|
[d7d2da3] | 70 |
|
---|
| 71 | TBranchMap::iterator itBranchMap;
|
---|
| 72 | map< TClass *, TProcessMethod >::iterator itClassMap;
|
---|
| 73 |
|
---|
| 74 | // read branch configuration and
|
---|
| 75 | // import array with output from filter/classifier/jetfinder modules
|
---|
| 76 |
|
---|
| 77 | ExRootConfParam param = GetParam("Branch");
|
---|
| 78 | Long_t i, size;
|
---|
| 79 | TString branchName, branchClassName, branchInputArray;
|
---|
| 80 | TClass *branchClass;
|
---|
| 81 | TObjArray *array;
|
---|
| 82 | ExRootTreeBranch *branch;
|
---|
| 83 |
|
---|
| 84 | size = param.GetSize();
|
---|
| 85 | for(i = 0; i < size/3; ++i)
|
---|
| 86 | {
|
---|
| 87 | branchInputArray = param[i*3].GetString();
|
---|
| 88 | branchName = param[i*3 + 1].GetString();
|
---|
| 89 | branchClassName = param[i*3 + 2].GetString();
|
---|
| 90 |
|
---|
| 91 | branchClass = gROOT->GetClass(branchClassName);
|
---|
| 92 |
|
---|
| 93 | if(!branchClass)
|
---|
| 94 | {
|
---|
| 95 | cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
|
---|
| 96 | continue;
|
---|
| 97 | }
|
---|
[8f7db23] | 98 |
|
---|
[d7d2da3] | 99 | itClassMap = fClassMap.find(branchClass);
|
---|
| 100 | if(itClassMap == fClassMap.end())
|
---|
| 101 | {
|
---|
| 102 | cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
|
---|
| 103 | continue;
|
---|
| 104 | }
|
---|
| 105 |
|
---|
| 106 | array = ImportArray(branchInputArray);
|
---|
| 107 | branch = NewBranch(branchName, branchClass);
|
---|
| 108 |
|
---|
| 109 | fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array)));
|
---|
| 110 | }
|
---|
| 111 |
|
---|
| 112 | }
|
---|
| 113 |
|
---|
| 114 | //------------------------------------------------------------------------------
|
---|
| 115 |
|
---|
| 116 | void TreeWriter::Finish()
|
---|
| 117 | {
|
---|
| 118 | }
|
---|
| 119 |
|
---|
| 120 | //------------------------------------------------------------------------------
|
---|
| 121 |
|
---|
| 122 | void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
|
---|
| 123 | {
|
---|
| 124 | TIter it1(candidate->GetCandidates());
|
---|
| 125 | it1.Reset();
|
---|
| 126 | array->Clear();
|
---|
| 127 | while((candidate = static_cast<Candidate*>(it1.Next())))
|
---|
| 128 | {
|
---|
| 129 | TIter it2(candidate->GetCandidates());
|
---|
| 130 |
|
---|
| 131 | // particle
|
---|
| 132 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
| 133 | {
|
---|
| 134 | array->Add(candidate);
|
---|
| 135 | continue;
|
---|
| 136 | }
|
---|
| 137 |
|
---|
| 138 | // track
|
---|
| 139 | candidate = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
|
---|
| 140 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
| 141 | {
|
---|
| 142 | array->Add(candidate);
|
---|
| 143 | continue;
|
---|
| 144 | }
|
---|
| 145 |
|
---|
| 146 | // tower
|
---|
| 147 | it2.Reset();
|
---|
| 148 | while((candidate = static_cast<Candidate*>(it2.Next())))
|
---|
| 149 | {
|
---|
| 150 | array->Add(candidate->GetCandidates()->At(0));
|
---|
| 151 | }
|
---|
| 152 | }
|
---|
| 153 | }
|
---|
| 154 |
|
---|
| 155 | //------------------------------------------------------------------------------
|
---|
| 156 |
|
---|
| 157 | void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 158 | {
|
---|
| 159 | TIter iterator(array);
|
---|
| 160 | Candidate *candidate = 0;
|
---|
| 161 | GenParticle *entry = 0;
|
---|
| 162 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[8f7db23] | 163 |
|
---|
[22dc7fd] | 164 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 165 |
|
---|
[d7d2da3] | 166 | // loop over all particles
|
---|
| 167 | iterator.Reset();
|
---|
| 168 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 169 | {
|
---|
| 170 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 171 | const TLorentzVector &position = candidate->Position;
|
---|
| 172 |
|
---|
| 173 | entry = static_cast<GenParticle*>(branch->NewEntry());
|
---|
| 174 |
|
---|
| 175 | entry->SetBit(kIsReferenced);
|
---|
| 176 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 177 |
|
---|
| 178 | pt = momentum.Pt();
|
---|
| 179 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 180 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 181 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 182 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 183 |
|
---|
| 184 | entry->PID = candidate->PID;
|
---|
| 185 |
|
---|
| 186 | entry->Status = candidate->Status;
|
---|
| 187 | entry->IsPU = candidate->IsPU;
|
---|
| 188 |
|
---|
| 189 | entry->M1 = candidate->M1;
|
---|
| 190 | entry->M2 = candidate->M2;
|
---|
| 191 |
|
---|
| 192 | entry->D1 = candidate->D1;
|
---|
| 193 | entry->D2 = candidate->D2;
|
---|
| 194 |
|
---|
| 195 | entry->Charge = candidate->Charge;
|
---|
| 196 | entry->Mass = candidate->Mass;
|
---|
| 197 |
|
---|
| 198 | entry->E = momentum.E();
|
---|
| 199 | entry->Px = momentum.Px();
|
---|
| 200 | entry->Py = momentum.Py();
|
---|
| 201 | entry->Pz = momentum.Pz();
|
---|
| 202 |
|
---|
| 203 | entry->Eta = eta;
|
---|
| 204 | entry->Phi = momentum.Phi();
|
---|
| 205 | entry->PT = pt;
|
---|
| 206 |
|
---|
| 207 | entry->Rapidity = rapidity;
|
---|
| 208 |
|
---|
| 209 | entry->X = position.X();
|
---|
| 210 | entry->Y = position.Y();
|
---|
| 211 | entry->Z = position.Z();
|
---|
[22dc7fd] | 212 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[d7d2da3] | 213 | }
|
---|
| 214 | }
|
---|
| 215 |
|
---|
| 216 | //------------------------------------------------------------------------------
|
---|
| 217 |
|
---|
[d07e957] | 218 | void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 219 | {
|
---|
| 220 | TIter iterator(array);
|
---|
| 221 | Candidate *candidate = 0;
|
---|
| 222 | Vertex *entry = 0;
|
---|
[8f7db23] | 223 |
|
---|
[22dc7fd] | 224 | const Double_t c_light = 2.99792458E8;
|
---|
[d07e957] | 225 |
|
---|
| 226 | // loop over all vertices
|
---|
| 227 | iterator.Reset();
|
---|
| 228 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 229 | {
|
---|
| 230 | const TLorentzVector &position = candidate->Position;
|
---|
| 231 |
|
---|
| 232 | entry = static_cast<Vertex*>(branch->NewEntry());
|
---|
| 233 |
|
---|
| 234 | entry->X = position.X();
|
---|
| 235 | entry->Y = position.Y();
|
---|
| 236 | entry->Z = position.Z();
|
---|
[22dc7fd] | 237 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[d07e957] | 238 | }
|
---|
| 239 | }
|
---|
| 240 |
|
---|
| 241 | //------------------------------------------------------------------------------
|
---|
| 242 |
|
---|
[d7d2da3] | 243 | void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 244 | {
|
---|
| 245 | TIter iterator(array);
|
---|
| 246 | Candidate *candidate = 0;
|
---|
| 247 | Candidate *particle = 0;
|
---|
| 248 | Track *entry = 0;
|
---|
| 249 | Double_t pt, signz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 250 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 251 |
|
---|
[d07e957] | 252 | // loop over all tracks
|
---|
[d7d2da3] | 253 | iterator.Reset();
|
---|
| 254 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 255 | {
|
---|
| 256 | const TLorentzVector &position = candidate->Position;
|
---|
| 257 |
|
---|
| 258 | cosTheta = TMath::Abs(position.CosTheta());
|
---|
| 259 | signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 260 | eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
|
---|
| 261 | rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
|
---|
| 262 |
|
---|
| 263 | entry = static_cast<Track*>(branch->NewEntry());
|
---|
| 264 |
|
---|
| 265 | entry->SetBit(kIsReferenced);
|
---|
| 266 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 267 |
|
---|
| 268 | entry->PID = candidate->PID;
|
---|
| 269 |
|
---|
| 270 | entry->Charge = candidate->Charge;
|
---|
| 271 |
|
---|
| 272 | entry->EtaOuter = eta;
|
---|
| 273 | entry->PhiOuter = position.Phi();
|
---|
| 274 |
|
---|
| 275 | entry->XOuter = position.X();
|
---|
| 276 | entry->YOuter = position.Y();
|
---|
| 277 | entry->ZOuter = position.Z();
|
---|
[22dc7fd] | 278 | entry->TOuter = position.T()*1.0E-3/c_light;
|
---|
[a0431dc] | 279 |
|
---|
| 280 | entry->Dxy = candidate->Dxy;
|
---|
| 281 | entry->SDxy = candidate->SDxy ;
|
---|
| 282 | entry->Xd = candidate -> Xd;
|
---|
| 283 | entry->Yd = candidate -> Yd;
|
---|
| 284 | entry->Zd = candidate -> Zd;
|
---|
| 285 |
|
---|
[d7d2da3] | 286 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 287 |
|
---|
| 288 | pt = momentum.Pt();
|
---|
| 289 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 290 | signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 291 | eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
|
---|
| 292 | rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
|
---|
| 293 |
|
---|
| 294 | entry->Eta = eta;
|
---|
| 295 | entry->Phi = momentum.Phi();
|
---|
| 296 | entry->PT = pt;
|
---|
| 297 |
|
---|
| 298 | particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
|
---|
| 299 | const TLorentzVector &initialPosition = particle->Position;
|
---|
| 300 |
|
---|
| 301 | entry->X = initialPosition.X();
|
---|
| 302 | entry->Y = initialPosition.Y();
|
---|
| 303 | entry->Z = initialPosition.Z();
|
---|
[22dc7fd] | 304 | entry->T = initialPosition.T()*1.0E-3/c_light;
|
---|
[d7d2da3] | 305 |
|
---|
| 306 | entry->Particle = particle;
|
---|
| 307 | }
|
---|
| 308 | }
|
---|
| 309 |
|
---|
| 310 | //------------------------------------------------------------------------------
|
---|
| 311 |
|
---|
| 312 | void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 313 | {
|
---|
| 314 | TIter iterator(array);
|
---|
| 315 | Candidate *candidate = 0;
|
---|
| 316 | Tower *entry = 0;
|
---|
| 317 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 318 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 319 |
|
---|
[d07e957] | 320 | // loop over all towers
|
---|
[d7d2da3] | 321 | iterator.Reset();
|
---|
| 322 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 323 | {
|
---|
| 324 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 325 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 326 |
|
---|
[d7d2da3] | 327 | pt = momentum.Pt();
|
---|
| 328 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 329 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 330 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 331 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 332 |
|
---|
| 333 | entry = static_cast<Tower*>(branch->NewEntry());
|
---|
| 334 |
|
---|
| 335 | entry->SetBit(kIsReferenced);
|
---|
| 336 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 337 |
|
---|
| 338 | entry->Eta = eta;
|
---|
| 339 | entry->Phi = momentum.Phi();
|
---|
| 340 | entry->ET = pt;
|
---|
| 341 | entry->E = momentum.E();
|
---|
| 342 | entry->Eem = candidate->Eem;
|
---|
| 343 | entry->Ehad = candidate->Ehad;
|
---|
| 344 | entry->Edges[0] = candidate->Edges[0];
|
---|
| 345 | entry->Edges[1] = candidate->Edges[1];
|
---|
| 346 | entry->Edges[2] = candidate->Edges[2];
|
---|
| 347 | entry->Edges[3] = candidate->Edges[3];
|
---|
[8f7db23] | 348 |
|
---|
[22dc7fd] | 349 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 350 |
|
---|
[d7d2da3] | 351 | FillParticles(candidate, &entry->Particles);
|
---|
| 352 | }
|
---|
| 353 | }
|
---|
| 354 |
|
---|
| 355 | //------------------------------------------------------------------------------
|
---|
| 356 |
|
---|
| 357 | void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 358 | {
|
---|
| 359 | TIter iterator(array);
|
---|
| 360 | Candidate *candidate = 0;
|
---|
| 361 | Photon *entry = 0;
|
---|
| 362 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 363 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 364 |
|
---|
[d7d2da3] | 365 | array->Sort();
|
---|
| 366 |
|
---|
| 367 | // loop over all photons
|
---|
| 368 | iterator.Reset();
|
---|
| 369 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 370 | {
|
---|
| 371 | TIter it1(candidate->GetCandidates());
|
---|
| 372 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 373 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 374 |
|
---|
[d7d2da3] | 375 |
|
---|
| 376 | pt = momentum.Pt();
|
---|
| 377 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 378 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 379 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 380 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 381 |
|
---|
| 382 | entry = static_cast<Photon*>(branch->NewEntry());
|
---|
| 383 |
|
---|
| 384 | entry->Eta = eta;
|
---|
| 385 | entry->Phi = momentum.Phi();
|
---|
| 386 | entry->PT = pt;
|
---|
| 387 | entry->E = momentum.E();
|
---|
[8f7db23] | 388 |
|
---|
[22dc7fd] | 389 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 390 |
|
---|
[d7d2da3] | 391 | entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
|
---|
| 392 |
|
---|
| 393 | FillParticles(candidate, &entry->Particles);
|
---|
| 394 | }
|
---|
| 395 | }
|
---|
| 396 |
|
---|
| 397 | //------------------------------------------------------------------------------
|
---|
| 398 |
|
---|
| 399 | void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 400 | {
|
---|
| 401 | TIter iterator(array);
|
---|
| 402 | Candidate *candidate = 0;
|
---|
| 403 | Electron *entry = 0;
|
---|
| 404 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 405 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 406 |
|
---|
[d7d2da3] | 407 | array->Sort();
|
---|
| 408 |
|
---|
| 409 | // loop over all electrons
|
---|
| 410 | iterator.Reset();
|
---|
| 411 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 412 | {
|
---|
| 413 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 414 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 415 |
|
---|
[d7d2da3] | 416 | pt = momentum.Pt();
|
---|
| 417 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 418 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 419 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 420 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 421 |
|
---|
| 422 | entry = static_cast<Electron*>(branch->NewEntry());
|
---|
| 423 |
|
---|
| 424 | entry->Eta = eta;
|
---|
| 425 | entry->Phi = momentum.Phi();
|
---|
| 426 | entry->PT = pt;
|
---|
[8f7db23] | 427 |
|
---|
[22dc7fd] | 428 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 429 |
|
---|
[d7d2da3] | 430 | entry->Charge = candidate->Charge;
|
---|
| 431 |
|
---|
| 432 | entry->EhadOverEem = 0.0;
|
---|
| 433 |
|
---|
| 434 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
| 435 | }
|
---|
| 436 | }
|
---|
| 437 |
|
---|
| 438 | //------------------------------------------------------------------------------
|
---|
| 439 |
|
---|
| 440 | void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 441 | {
|
---|
| 442 | TIter iterator(array);
|
---|
| 443 | Candidate *candidate = 0;
|
---|
| 444 | Muon *entry = 0;
|
---|
| 445 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[8f7db23] | 446 |
|
---|
[22dc7fd] | 447 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 448 |
|
---|
[d7d2da3] | 449 | array->Sort();
|
---|
| 450 |
|
---|
| 451 | // loop over all muons
|
---|
| 452 | iterator.Reset();
|
---|
| 453 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 454 | {
|
---|
| 455 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 456 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 457 |
|
---|
[d7d2da3] | 458 |
|
---|
| 459 | pt = momentum.Pt();
|
---|
| 460 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 461 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 462 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 463 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 464 |
|
---|
| 465 | entry = static_cast<Muon*>(branch->NewEntry());
|
---|
| 466 |
|
---|
| 467 | entry->SetBit(kIsReferenced);
|
---|
| 468 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 469 |
|
---|
| 470 | entry->Eta = eta;
|
---|
| 471 | entry->Phi = momentum.Phi();
|
---|
| 472 | entry->PT = pt;
|
---|
| 473 |
|
---|
[22dc7fd] | 474 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 475 |
|
---|
[22dc7fd] | 476 | // cout<<entry->PT<<","<<entry->T<<endl;
|
---|
| 477 |
|
---|
[d7d2da3] | 478 | entry->Charge = candidate->Charge;
|
---|
| 479 |
|
---|
| 480 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
| 481 | }
|
---|
| 482 | }
|
---|
| 483 |
|
---|
| 484 | //------------------------------------------------------------------------------
|
---|
| 485 |
|
---|
| 486 | void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 487 | {
|
---|
| 488 | TIter iterator(array);
|
---|
| 489 | Candidate *candidate = 0, *constituent = 0;
|
---|
| 490 | Jet *entry = 0;
|
---|
| 491 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
| 492 | Double_t ecalEnergy, hcalEnergy;
|
---|
[22dc7fd] | 493 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 494 |
|
---|
[d7d2da3] | 495 | array->Sort();
|
---|
| 496 |
|
---|
| 497 | // loop over all jets
|
---|
| 498 | iterator.Reset();
|
---|
| 499 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 500 | {
|
---|
| 501 | TIter itConstituents(candidate->GetCandidates());
|
---|
[8f7db23] | 502 |
|
---|
[d7d2da3] | 503 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 504 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 505 |
|
---|
[d7d2da3] | 506 | pt = momentum.Pt();
|
---|
| 507 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 508 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 509 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 510 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 511 |
|
---|
| 512 | entry = static_cast<Jet*>(branch->NewEntry());
|
---|
| 513 |
|
---|
| 514 | entry->Eta = eta;
|
---|
| 515 | entry->Phi = momentum.Phi();
|
---|
| 516 | entry->PT = pt;
|
---|
| 517 |
|
---|
[22dc7fd] | 518 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 519 |
|
---|
[d7d2da3] | 520 | entry->Mass = momentum.M();
|
---|
| 521 |
|
---|
| 522 | entry->DeltaEta = candidate->DeltaEta;
|
---|
| 523 | entry->DeltaPhi = candidate->DeltaPhi;
|
---|
| 524 |
|
---|
| 525 | entry->BTag = candidate->BTag;
|
---|
| 526 | entry->TauTag = candidate->TauTag;
|
---|
| 527 |
|
---|
| 528 | entry->Charge = candidate->Charge;
|
---|
| 529 |
|
---|
| 530 | itConstituents.Reset();
|
---|
| 531 | entry->Constituents.Clear();
|
---|
| 532 | ecalEnergy = 0.0;
|
---|
| 533 | hcalEnergy = 0.0;
|
---|
| 534 | while((constituent = static_cast<Candidate*>(itConstituents.Next())))
|
---|
| 535 | {
|
---|
| 536 | entry->Constituents.Add(constituent);
|
---|
| 537 | ecalEnergy += constituent->Eem;
|
---|
| 538 | hcalEnergy += constituent->Ehad;
|
---|
| 539 | }
|
---|
| 540 |
|
---|
| 541 | entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
|
---|
[8f7db23] | 542 |
|
---|
[24d005f] | 543 | //--- Pile-Up Jet ID variables ----
|
---|
| 544 |
|
---|
| 545 | entry->NCharged = candidate->NCharged;
|
---|
| 546 | entry->NNeutrals = candidate->NNeutrals;
|
---|
| 547 | entry->Beta = candidate->Beta;
|
---|
| 548 | entry->BetaStar = candidate->BetaStar;
|
---|
| 549 | entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
|
---|
| 550 | entry->PTD = candidate->PTD;
|
---|
| 551 | entry->FracPt[0] = candidate->FracPt[0];
|
---|
| 552 | entry->FracPt[1] = candidate->FracPt[1];
|
---|
| 553 | entry->FracPt[2] = candidate->FracPt[2];
|
---|
| 554 | entry->FracPt[3] = candidate->FracPt[3];
|
---|
| 555 | entry->FracPt[4] = candidate->FracPt[4];
|
---|
[9687203] | 556 |
|
---|
| 557 | //--- N-subjettiness variables ----
|
---|
| 558 |
|
---|
| 559 | entry->Tau1 = candidate->Tau1;
|
---|
| 560 | entry->Tau2 = candidate->Tau2;
|
---|
| 561 | entry->Tau3 = candidate->Tau3;
|
---|
| 562 | entry->Tau4 = candidate->Tau4;
|
---|
| 563 | entry->Tau5 = candidate->Tau5;
|
---|
| 564 |
|
---|
[d7d2da3] | 565 | FillParticles(candidate, &entry->Particles);
|
---|
| 566 | }
|
---|
| 567 | }
|
---|
| 568 |
|
---|
| 569 | //------------------------------------------------------------------------------
|
---|
| 570 |
|
---|
| 571 | void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 572 | {
|
---|
| 573 | Candidate *candidate = 0;
|
---|
| 574 | MissingET *entry = 0;
|
---|
| 575 |
|
---|
| 576 | // get the first entry
|
---|
| 577 | if((candidate = static_cast<Candidate*>(array->At(0))))
|
---|
| 578 | {
|
---|
| 579 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 580 |
|
---|
| 581 | entry = static_cast<MissingET*>(branch->NewEntry());
|
---|
| 582 |
|
---|
[3b465ca] | 583 | entry->Eta = (-momentum).Eta();
|
---|
[d7d2da3] | 584 | entry->Phi = (-momentum).Phi();
|
---|
| 585 | entry->MET = momentum.Pt();
|
---|
| 586 | }
|
---|
| 587 | }
|
---|
| 588 |
|
---|
| 589 | //------------------------------------------------------------------------------
|
---|
| 590 |
|
---|
| 591 | void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 592 | {
|
---|
| 593 | Candidate *candidate = 0;
|
---|
| 594 | ScalarHT *entry = 0;
|
---|
| 595 |
|
---|
| 596 | // get the first entry
|
---|
| 597 | if((candidate = static_cast<Candidate*>(array->At(0))))
|
---|
| 598 | {
|
---|
| 599 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 600 |
|
---|
| 601 | entry = static_cast<ScalarHT*>(branch->NewEntry());
|
---|
| 602 |
|
---|
| 603 | entry->HT = momentum.Pt();
|
---|
| 604 | }
|
---|
| 605 | }
|
---|
| 606 |
|
---|
| 607 | //------------------------------------------------------------------------------
|
---|
| 608 |
|
---|
[71648c2] | 609 | void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 610 | {
|
---|
[3b465ca] | 611 | TIter iterator(array);
|
---|
[71648c2] | 612 | Candidate *candidate = 0;
|
---|
| 613 | Rho *entry = 0;
|
---|
| 614 |
|
---|
[3b465ca] | 615 | // loop over all rho
|
---|
| 616 | iterator.Reset();
|
---|
| 617 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
[71648c2] | 618 | {
|
---|
| 619 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 620 |
|
---|
| 621 | entry = static_cast<Rho*>(branch->NewEntry());
|
---|
| 622 |
|
---|
[8608ef8] | 623 | entry->Rho = momentum.E();
|
---|
[3b465ca] | 624 | entry->Edges[0] = candidate->Edges[0];
|
---|
| 625 | entry->Edges[1] = candidate->Edges[1];
|
---|
[71648c2] | 626 | }
|
---|
| 627 | }
|
---|
| 628 |
|
---|
| 629 | //------------------------------------------------------------------------------
|
---|
| 630 |
|
---|
[2e229c9] | 631 | void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 632 | {
|
---|
| 633 | Candidate *candidate = 0;
|
---|
| 634 | Weight *entry = 0;
|
---|
| 635 |
|
---|
| 636 | // get the first entry
|
---|
| 637 | if((candidate = static_cast<Candidate*>(array->At(0))))
|
---|
| 638 | {
|
---|
| 639 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 640 |
|
---|
| 641 | entry = static_cast<Weight*>(branch->NewEntry());
|
---|
| 642 |
|
---|
| 643 | entry->Weight = momentum.E();
|
---|
| 644 | }
|
---|
| 645 | }
|
---|
| 646 |
|
---|
| 647 | //------------------------------------------------------------------------------
|
---|
| 648 |
|
---|
[8f7db23] | 649 | void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 650 | {
|
---|
| 651 | TIter iterator(array);
|
---|
| 652 | Candidate *candidate = 0;
|
---|
| 653 | HectorHit *entry = 0;
|
---|
| 654 |
|
---|
| 655 | // loop over all roman pot hits
|
---|
| 656 | iterator.Reset();
|
---|
| 657 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 658 | {
|
---|
| 659 | const TLorentzVector &position = candidate->Position;
|
---|
| 660 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 661 |
|
---|
| 662 | entry = static_cast<HectorHit*>(branch->NewEntry());
|
---|
| 663 |
|
---|
| 664 | entry->E = momentum.E();
|
---|
| 665 |
|
---|
| 666 | entry->Tx = momentum.Px();
|
---|
| 667 | entry->Ty = momentum.Py();
|
---|
| 668 |
|
---|
| 669 | entry->T = position.T();
|
---|
| 670 |
|
---|
| 671 | entry->X = position.X();
|
---|
| 672 | entry->Y = position.Y();
|
---|
| 673 | entry->S = position.Z();
|
---|
[64a4950] | 674 |
|
---|
| 675 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
[8f7db23] | 676 | }
|
---|
| 677 | }
|
---|
| 678 |
|
---|
| 679 | //------------------------------------------------------------------------------
|
---|
| 680 |
|
---|
[d7d2da3] | 681 | void TreeWriter::Process()
|
---|
| 682 | {
|
---|
| 683 | TBranchMap::iterator itBranchMap;
|
---|
| 684 | ExRootTreeBranch *branch;
|
---|
| 685 | TProcessMethod method;
|
---|
| 686 | TObjArray *array;
|
---|
| 687 |
|
---|
| 688 | for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
|
---|
| 689 | {
|
---|
| 690 | branch = itBranchMap->first;
|
---|
| 691 | method = itBranchMap->second.first;
|
---|
| 692 | array = itBranchMap->second.second;
|
---|
| 693 |
|
---|
| 694 | (this->*method)(branch, array);
|
---|
| 695 | }
|
---|
| 696 | }
|
---|
| 697 |
|
---|
| 698 | //------------------------------------------------------------------------------
|
---|