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