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