[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 |
|
---|
| 20 | /** \class TreeWriter
|
---|
| 21 | *
|
---|
| 22 | * Fills ROOT tree branches.
|
---|
| 23 | *
|
---|
| 24 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 25 | *
|
---|
| 26 | */
|
---|
| 27 |
|
---|
| 28 | #include "modules/TreeWriter.h"
|
---|
| 29 |
|
---|
| 30 | #include "classes/DelphesClasses.h"
|
---|
| 31 | #include "classes/DelphesFactory.h"
|
---|
| 32 | #include "classes/DelphesFormula.h"
|
---|
| 33 |
|
---|
| 34 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 35 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 36 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 37 | #include "ExRootAnalysis/ExRootTreeBranch.h"
|
---|
| 38 |
|
---|
| 39 | #include "TROOT.h"
|
---|
| 40 | #include "TMath.h"
|
---|
| 41 | #include "TString.h"
|
---|
| 42 | #include "TFormula.h"
|
---|
| 43 | #include "TRandom3.h"
|
---|
| 44 | #include "TObjArray.h"
|
---|
| 45 | #include "TDatabasePDG.h"
|
---|
| 46 | #include "TLorentzVector.h"
|
---|
| 47 |
|
---|
| 48 | #include <algorithm>
|
---|
| 49 | #include <stdexcept>
|
---|
| 50 | #include <iostream>
|
---|
| 51 | #include <sstream>
|
---|
| 52 |
|
---|
| 53 | using namespace std;
|
---|
| 54 |
|
---|
| 55 | //------------------------------------------------------------------------------
|
---|
| 56 |
|
---|
| 57 | TreeWriter::TreeWriter()
|
---|
| 58 | {
|
---|
| 59 | }
|
---|
| 60 |
|
---|
| 61 | //------------------------------------------------------------------------------
|
---|
| 62 |
|
---|
| 63 | TreeWriter::~TreeWriter()
|
---|
| 64 | {
|
---|
| 65 | }
|
---|
| 66 |
|
---|
| 67 | //------------------------------------------------------------------------------
|
---|
| 68 |
|
---|
| 69 | void TreeWriter::Init()
|
---|
| 70 | {
|
---|
| 71 | fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
|
---|
[2de4dcd] | 72 | fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
|
---|
[d7d2da3] | 73 | fClassMap[Track::Class()] = &TreeWriter::ProcessTracks;
|
---|
| 74 | fClassMap[Tower::Class()] = &TreeWriter::ProcessTowers;
|
---|
| 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;
|
---|
| 86 | map< TClass *, TProcessMethod >::iterator itClassMap;
|
---|
| 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();
|
---|
| 99 | for(i = 0; i < size/3; ++i)
|
---|
| 100 | {
|
---|
| 101 | branchInputArray = param[i*3].GetString();
|
---|
| 102 | branchName = param[i*3 + 1].GetString();
|
---|
| 103 | branchClassName = param[i*3 + 2].GetString();
|
---|
| 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 |
|
---|
| 130 | void TreeWriter::Finish()
|
---|
| 131 | {
|
---|
| 132 | }
|
---|
| 133 |
|
---|
| 134 | //------------------------------------------------------------------------------
|
---|
| 135 |
|
---|
| 136 | void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
|
---|
| 137 | {
|
---|
| 138 | TIter it1(candidate->GetCandidates());
|
---|
| 139 | it1.Reset();
|
---|
| 140 | array->Clear();
|
---|
| 141 | while((candidate = static_cast<Candidate*>(it1.Next())))
|
---|
| 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
|
---|
| 153 | candidate = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
|
---|
| 154 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
| 155 | {
|
---|
| 156 | array->Add(candidate);
|
---|
| 157 | continue;
|
---|
| 158 | }
|
---|
| 159 |
|
---|
| 160 | // tower
|
---|
| 161 | it2.Reset();
|
---|
| 162 | while((candidate = static_cast<Candidate*>(it2.Next())))
|
---|
| 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();
|
---|
| 182 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 183 | {
|
---|
| 184 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 185 | const TLorentzVector &position = candidate->Position;
|
---|
| 186 |
|
---|
| 187 | entry = static_cast<GenParticle*>(branch->NewEntry());
|
---|
| 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;
|
---|
| 195 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 196 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 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 |
|
---|
[acd0621] | 217 | entry->D0 = candidate->D0;
|
---|
| 218 | entry->DZ = candidate->DZ;
|
---|
| 219 | entry->P = candidate->P;
|
---|
| 220 | entry->PT = candidate->PT;
|
---|
| 221 | entry->CtgTheta = candidate->CtgTheta;
|
---|
| 222 | entry->Phi = candidate->Phi;
|
---|
| 223 |
|
---|
[d7d2da3] | 224 | entry->Eta = eta;
|
---|
| 225 | entry->Phi = momentum.Phi();
|
---|
| 226 | entry->PT = pt;
|
---|
| 227 |
|
---|
| 228 | entry->Rapidity = rapidity;
|
---|
| 229 |
|
---|
| 230 | entry->X = position.X();
|
---|
| 231 | entry->Y = position.Y();
|
---|
| 232 | entry->Z = position.Z();
|
---|
[22dc7fd] | 233 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[d7d2da3] | 234 | }
|
---|
| 235 | }
|
---|
| 236 |
|
---|
| 237 | //------------------------------------------------------------------------------
|
---|
| 238 |
|
---|
[d07e957] | 239 | void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 240 | {
|
---|
| 241 | TIter iterator(array);
|
---|
[5496767] | 242 | Candidate *candidate = 0, *constituent = 0;
|
---|
[d07e957] | 243 | Vertex *entry = 0;
|
---|
[8f7db23] | 244 |
|
---|
[22dc7fd] | 245 | const Double_t c_light = 2.99792458E8;
|
---|
[5496767] | 246 |
|
---|
[332025f] | 247 | Double_t x, y, z, t, xError, yError, zError, tError, sigma, sumPT2, btvSumPT2, genDeltaZ, genSumPT2;
|
---|
[2600216] | 248 | UInt_t index, ndf;
|
---|
[5496767] | 249 |
|
---|
[3e2bb2b] | 250 | CompBase *compare = Candidate::fgCompare;
|
---|
| 251 | Candidate::fgCompare = CompSumPT2<Candidate>::Instance();
|
---|
[3c46e17] | 252 | array->Sort();
|
---|
[3e2bb2b] | 253 | Candidate::fgCompare = compare;
|
---|
[5496767] | 254 |
|
---|
[d07e957] | 255 | // loop over all vertices
|
---|
| 256 | iterator.Reset();
|
---|
| 257 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 258 | {
|
---|
[5496767] | 259 |
|
---|
[2600216] | 260 | index = candidate->ClusterIndex;
|
---|
| 261 | ndf = candidate->ClusterNDF;
|
---|
| 262 | sigma = candidate->ClusterSigma;
|
---|
| 263 | sumPT2 = candidate->SumPT2;
|
---|
| 264 | btvSumPT2 = candidate->BTVSumPT2;
|
---|
| 265 | genDeltaZ = candidate->GenDeltaZ;
|
---|
| 266 | genSumPT2 = candidate->GenSumPT2;
|
---|
[5496767] | 267 |
|
---|
[7bcca65] | 268 | x = candidate->Position.X();
|
---|
| 269 | y = candidate->Position.Y();
|
---|
| 270 | z = candidate->Position.Z();
|
---|
[2600216] | 271 | t = candidate->Position.T()*1.0E-3/c_light;
|
---|
[5496767] | 272 |
|
---|
[2600216] | 273 | xError = candidate->PositionError.X ();
|
---|
| 274 | yError = candidate->PositionError.Y ();
|
---|
| 275 | zError = candidate->PositionError.Z ();
|
---|
[332025f] | 276 | tError = candidate->PositionError.T ()*1.0E-3/c_light;
|
---|
[d07e957] | 277 |
|
---|
| 278 | entry = static_cast<Vertex*>(branch->NewEntry());
|
---|
| 279 |
|
---|
[2600216] | 280 | entry->Index = index;
|
---|
| 281 | entry->NDF = ndf;
|
---|
| 282 | entry->Sigma = sigma;
|
---|
| 283 | entry->SumPT2 = sumPT2;
|
---|
| 284 | entry->BTVSumPT2 = btvSumPT2;
|
---|
| 285 | entry->GenDeltaZ = genDeltaZ;
|
---|
| 286 | entry->GenSumPT2 = genSumPT2;
|
---|
[5496767] | 287 |
|
---|
[2600216] | 288 | entry->X = x;
|
---|
| 289 | entry->Y = y;
|
---|
| 290 | entry->Z = z;
|
---|
| 291 | entry->T = t;
|
---|
[5496767] | 292 |
|
---|
[2600216] | 293 | entry->ErrorX = xError;
|
---|
| 294 | entry->ErrorY = yError;
|
---|
| 295 | entry->ErrorZ = zError;
|
---|
[332025f] | 296 | entry->ErrorT = tError;
|
---|
[5496767] | 297 |
|
---|
| 298 |
|
---|
[3e2bb2b] | 299 | TIter itConstituents(candidate->GetCandidates());
|
---|
| 300 | itConstituents.Reset();
|
---|
| 301 | entry->Constituents.Clear();
|
---|
[5496767] | 302 | while((constituent = static_cast<Candidate*>(itConstituents.Next())))
|
---|
| 303 | {
|
---|
| 304 | entry->Constituents.Add(constituent);
|
---|
| 305 | }
|
---|
| 306 |
|
---|
[d07e957] | 307 | }
|
---|
| 308 | }
|
---|
| 309 |
|
---|
[2600216] | 310 |
|
---|
[d07e957] | 311 | //------------------------------------------------------------------------------
|
---|
| 312 |
|
---|
[d7d2da3] | 313 | void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 314 | {
|
---|
| 315 | TIter iterator(array);
|
---|
| 316 | Candidate *candidate = 0;
|
---|
| 317 | Candidate *particle = 0;
|
---|
| 318 | Track *entry = 0;
|
---|
[2a3eb22] | 319 | Double_t pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi;
|
---|
[22dc7fd] | 320 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 321 |
|
---|
[d07e957] | 322 | // loop over all tracks
|
---|
[d7d2da3] | 323 | iterator.Reset();
|
---|
| 324 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 325 | {
|
---|
| 326 | const TLorentzVector &position = candidate->Position;
|
---|
| 327 |
|
---|
| 328 | cosTheta = TMath::Abs(position.CosTheta());
|
---|
| 329 | signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 330 | eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
|
---|
| 331 | rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
|
---|
| 332 |
|
---|
| 333 | entry = static_cast<Track*>(branch->NewEntry());
|
---|
| 334 |
|
---|
| 335 | entry->SetBit(kIsReferenced);
|
---|
| 336 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 337 |
|
---|
| 338 | entry->PID = candidate->PID;
|
---|
| 339 |
|
---|
| 340 | entry->Charge = candidate->Charge;
|
---|
| 341 |
|
---|
| 342 | entry->EtaOuter = eta;
|
---|
| 343 | entry->PhiOuter = position.Phi();
|
---|
| 344 |
|
---|
| 345 | entry->XOuter = position.X();
|
---|
| 346 | entry->YOuter = position.Y();
|
---|
| 347 | entry->ZOuter = position.Z();
|
---|
[22dc7fd] | 348 | entry->TOuter = position.T()*1.0E-3/c_light;
|
---|
[5496767] | 349 |
|
---|
[acd0621] | 350 | entry->L = candidate->L;
|
---|
[5496767] | 351 |
|
---|
[acd0621] | 352 | entry->D0 = candidate->D0;
|
---|
| 353 | entry->ErrorD0 = candidate->ErrorD0;
|
---|
| 354 | entry->DZ = candidate->DZ;
|
---|
| 355 | entry->ErrorDZ = candidate->ErrorDZ;
|
---|
[2a3eb22] | 356 |
|
---|
[acd0621] | 357 | entry->ErrorP = candidate->ErrorP;
|
---|
| 358 | entry->ErrorPT = candidate->ErrorPT;
|
---|
| 359 | entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
|
---|
[5496767] | 360 | entry->ErrorPhi = candidate->ErrorPhi;
|
---|
| 361 |
|
---|
[e4c3fef] | 362 | entry->Xd = candidate->Xd;
|
---|
| 363 | entry->Yd = candidate->Yd;
|
---|
| 364 | entry->Zd = candidate->Zd;
|
---|
[9040259] | 365 |
|
---|
[d7d2da3] | 366 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 367 |
|
---|
| 368 | pt = momentum.Pt();
|
---|
[2a3eb22] | 369 | p = momentum.P();
|
---|
| 370 | phi = momentum.Phi();
|
---|
| 371 | ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1/TMath::Tan(momentum.Theta()) : 1e10;
|
---|
| 372 |
|
---|
[d7d2da3] | 373 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 374 | signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 375 | eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
|
---|
| 376 | rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
|
---|
| 377 |
|
---|
[c1ea422] | 378 | entry->P = p;
|
---|
[2a3eb22] | 379 | entry->PT = pt;
|
---|
[d7d2da3] | 380 | entry->Eta = eta;
|
---|
[2a3eb22] | 381 | entry->Phi = phi;
|
---|
| 382 | entry->CtgTheta = ctgTheta;
|
---|
[5496767] | 383 |
|
---|
[d7d2da3] | 384 | particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
|
---|
| 385 | const TLorentzVector &initialPosition = particle->Position;
|
---|
| 386 |
|
---|
| 387 | entry->X = initialPosition.X();
|
---|
| 388 | entry->Y = initialPosition.Y();
|
---|
| 389 | entry->Z = initialPosition.Z();
|
---|
[22dc7fd] | 390 | entry->T = initialPosition.T()*1.0E-3/c_light;
|
---|
[d7d2da3] | 391 |
|
---|
| 392 | entry->Particle = particle;
|
---|
[2600216] | 393 |
|
---|
| 394 | entry->VertexIndex = candidate->ClusterIndex;
|
---|
| 395 |
|
---|
[d7d2da3] | 396 | }
|
---|
| 397 | }
|
---|
| 398 |
|
---|
| 399 | //------------------------------------------------------------------------------
|
---|
| 400 |
|
---|
| 401 | void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 402 | {
|
---|
| 403 | TIter iterator(array);
|
---|
| 404 | Candidate *candidate = 0;
|
---|
| 405 | Tower *entry = 0;
|
---|
| 406 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 407 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 408 |
|
---|
[d07e957] | 409 | // loop over all towers
|
---|
[d7d2da3] | 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<Tower*>(branch->NewEntry());
|
---|
| 423 |
|
---|
| 424 | entry->SetBit(kIsReferenced);
|
---|
| 425 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 426 |
|
---|
| 427 | entry->Eta = eta;
|
---|
| 428 | entry->Phi = momentum.Phi();
|
---|
| 429 | entry->ET = pt;
|
---|
| 430 | entry->E = momentum.E();
|
---|
| 431 | entry->Eem = candidate->Eem;
|
---|
| 432 | entry->Ehad = candidate->Ehad;
|
---|
| 433 | entry->Edges[0] = candidate->Edges[0];
|
---|
| 434 | entry->Edges[1] = candidate->Edges[1];
|
---|
| 435 | entry->Edges[2] = candidate->Edges[2];
|
---|
| 436 | entry->Edges[3] = candidate->Edges[3];
|
---|
[8f7db23] | 437 |
|
---|
[22dc7fd] | 438 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[839deb7] | 439 | entry->NTimeHits = candidate->NTimeHits;
|
---|
[8f7db23] | 440 |
|
---|
[d7d2da3] | 441 | FillParticles(candidate, &entry->Particles);
|
---|
| 442 | }
|
---|
| 443 | }
|
---|
| 444 |
|
---|
| 445 | //------------------------------------------------------------------------------
|
---|
| 446 |
|
---|
| 447 | void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 448 | {
|
---|
| 449 | TIter iterator(array);
|
---|
| 450 | Candidate *candidate = 0;
|
---|
| 451 | Photon *entry = 0;
|
---|
| 452 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 453 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 454 |
|
---|
[d7d2da3] | 455 | array->Sort();
|
---|
| 456 |
|
---|
| 457 | // loop over all photons
|
---|
| 458 | iterator.Reset();
|
---|
| 459 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 460 | {
|
---|
| 461 | TIter it1(candidate->GetCandidates());
|
---|
| 462 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 463 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 464 |
|
---|
[d7d2da3] | 465 |
|
---|
| 466 | pt = momentum.Pt();
|
---|
| 467 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 468 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 469 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 470 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 471 |
|
---|
| 472 | entry = static_cast<Photon*>(branch->NewEntry());
|
---|
| 473 |
|
---|
| 474 | entry->Eta = eta;
|
---|
| 475 | entry->Phi = momentum.Phi();
|
---|
| 476 | entry->PT = pt;
|
---|
| 477 | entry->E = momentum.E();
|
---|
[8f7db23] | 478 |
|
---|
[22dc7fd] | 479 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 480 |
|
---|
[d5cae2b] | 481 | // Isolation variables
|
---|
[9040259] | 482 |
|
---|
[d5cae2b] | 483 | entry->IsolationVar = candidate->IsolationVar;
|
---|
| 484 | entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr ;
|
---|
| 485 | entry->SumPtCharged = candidate->SumPtCharged ;
|
---|
| 486 | entry->SumPtNeutral = candidate->SumPtNeutral ;
|
---|
| 487 | entry->SumPtChargedPU = candidate->SumPtChargedPU ;
|
---|
| 488 | entry->SumPt = candidate->SumPt ;
|
---|
| 489 |
|
---|
[d7d2da3] | 490 | entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
|
---|
| 491 |
|
---|
| 492 | FillParticles(candidate, &entry->Particles);
|
---|
| 493 | }
|
---|
| 494 | }
|
---|
| 495 |
|
---|
| 496 | //------------------------------------------------------------------------------
|
---|
| 497 |
|
---|
| 498 | void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 499 | {
|
---|
| 500 | TIter iterator(array);
|
---|
| 501 | Candidate *candidate = 0;
|
---|
| 502 | Electron *entry = 0;
|
---|
| 503 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[22dc7fd] | 504 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 505 |
|
---|
[d7d2da3] | 506 | array->Sort();
|
---|
| 507 |
|
---|
| 508 | // loop over all electrons
|
---|
| 509 | iterator.Reset();
|
---|
| 510 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 511 | {
|
---|
| 512 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 513 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 514 |
|
---|
[d7d2da3] | 515 | pt = momentum.Pt();
|
---|
| 516 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 517 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 518 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 519 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 520 |
|
---|
| 521 | entry = static_cast<Electron*>(branch->NewEntry());
|
---|
| 522 |
|
---|
| 523 | entry->Eta = eta;
|
---|
| 524 | entry->Phi = momentum.Phi();
|
---|
| 525 | entry->PT = pt;
|
---|
[8f7db23] | 526 |
|
---|
[22dc7fd] | 527 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 528 |
|
---|
[d5cae2b] | 529 | // Isolation variables
|
---|
[9040259] | 530 |
|
---|
[d5cae2b] | 531 | entry->IsolationVar = candidate->IsolationVar;
|
---|
| 532 | entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr ;
|
---|
| 533 | entry->SumPtCharged = candidate->SumPtCharged ;
|
---|
| 534 | entry->SumPtNeutral = candidate->SumPtNeutral ;
|
---|
| 535 | entry->SumPtChargedPU = candidate->SumPtChargedPU ;
|
---|
| 536 | entry->SumPt = candidate->SumPt ;
|
---|
| 537 |
|
---|
| 538 |
|
---|
[d7d2da3] | 539 | entry->Charge = candidate->Charge;
|
---|
| 540 |
|
---|
| 541 | entry->EhadOverEem = 0.0;
|
---|
| 542 |
|
---|
| 543 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
| 544 | }
|
---|
| 545 | }
|
---|
| 546 |
|
---|
| 547 | //------------------------------------------------------------------------------
|
---|
| 548 |
|
---|
| 549 | void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 550 | {
|
---|
| 551 | TIter iterator(array);
|
---|
| 552 | Candidate *candidate = 0;
|
---|
| 553 | Muon *entry = 0;
|
---|
| 554 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
[8f7db23] | 555 |
|
---|
[22dc7fd] | 556 | const Double_t c_light = 2.99792458E8;
|
---|
[8f7db23] | 557 |
|
---|
[d7d2da3] | 558 | array->Sort();
|
---|
| 559 |
|
---|
| 560 | // loop over all muons
|
---|
| 561 | iterator.Reset();
|
---|
| 562 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 563 | {
|
---|
| 564 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 565 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 566 |
|
---|
[d7d2da3] | 567 | pt = momentum.Pt();
|
---|
| 568 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 569 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 570 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 571 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 572 |
|
---|
| 573 | entry = static_cast<Muon*>(branch->NewEntry());
|
---|
| 574 |
|
---|
| 575 | entry->SetBit(kIsReferenced);
|
---|
| 576 | entry->SetUniqueID(candidate->GetUniqueID());
|
---|
| 577 |
|
---|
| 578 | entry->Eta = eta;
|
---|
| 579 | entry->Phi = momentum.Phi();
|
---|
| 580 | entry->PT = pt;
|
---|
| 581 |
|
---|
[22dc7fd] | 582 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[9040259] | 583 |
|
---|
[d5cae2b] | 584 | // Isolation variables
|
---|
[9040259] | 585 |
|
---|
[d5cae2b] | 586 | entry->IsolationVar = candidate->IsolationVar;
|
---|
| 587 | entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr ;
|
---|
| 588 | entry->SumPtCharged = candidate->SumPtCharged ;
|
---|
| 589 | entry->SumPtNeutral = candidate->SumPtNeutral ;
|
---|
| 590 | entry->SumPtChargedPU = candidate->SumPtChargedPU ;
|
---|
| 591 | entry->SumPt = candidate->SumPt ;
|
---|
[8f7db23] | 592 |
|
---|
[d7d2da3] | 593 | entry->Charge = candidate->Charge;
|
---|
| 594 |
|
---|
| 595 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
| 596 | }
|
---|
| 597 | }
|
---|
| 598 |
|
---|
| 599 | //------------------------------------------------------------------------------
|
---|
| 600 |
|
---|
| 601 | void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 602 | {
|
---|
| 603 | TIter iterator(array);
|
---|
| 604 | Candidate *candidate = 0, *constituent = 0;
|
---|
| 605 | Jet *entry = 0;
|
---|
| 606 | Double_t pt, signPz, cosTheta, eta, rapidity;
|
---|
| 607 | Double_t ecalEnergy, hcalEnergy;
|
---|
[22dc7fd] | 608 | const Double_t c_light = 2.99792458E8;
|
---|
[de6d698] | 609 | Int_t i;
|
---|
[8f7db23] | 610 |
|
---|
[d7d2da3] | 611 | array->Sort();
|
---|
| 612 |
|
---|
| 613 | // loop over all jets
|
---|
| 614 | iterator.Reset();
|
---|
| 615 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 616 | {
|
---|
| 617 | TIter itConstituents(candidate->GetCandidates());
|
---|
[8f7db23] | 618 |
|
---|
[d7d2da3] | 619 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
[22dc7fd] | 620 | const TLorentzVector &position = candidate->Position;
|
---|
[8f7db23] | 621 |
|
---|
[d7d2da3] | 622 | pt = momentum.Pt();
|
---|
| 623 | cosTheta = TMath::Abs(momentum.CosTheta());
|
---|
| 624 | signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
|
---|
| 625 | eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
|
---|
| 626 | rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
|
---|
| 627 |
|
---|
| 628 | entry = static_cast<Jet*>(branch->NewEntry());
|
---|
| 629 |
|
---|
| 630 | entry->Eta = eta;
|
---|
| 631 | entry->Phi = momentum.Phi();
|
---|
| 632 | entry->PT = pt;
|
---|
| 633 |
|
---|
[22dc7fd] | 634 | entry->T = position.T()*1.0E-3/c_light;
|
---|
[8f7db23] | 635 |
|
---|
[d7d2da3] | 636 | entry->Mass = momentum.M();
|
---|
| 637 |
|
---|
[ba1f1ee] | 638 | entry->Area = candidate->Area;
|
---|
| 639 |
|
---|
[d7d2da3] | 640 | entry->DeltaEta = candidate->DeltaEta;
|
---|
| 641 | entry->DeltaPhi = candidate->DeltaPhi;
|
---|
| 642 |
|
---|
[fe0273c] | 643 | entry->Flavor = candidate->Flavor;
|
---|
| 644 | entry->FlavorAlgo = candidate->FlavorAlgo;
|
---|
| 645 | entry->FlavorPhys = candidate->FlavorPhys;
|
---|
| 646 |
|
---|
[d7d2da3] | 647 | entry->BTag = candidate->BTag;
|
---|
[9040259] | 648 |
|
---|
| 649 | entry->BTagAlgo = candidate->BTagAlgo;
|
---|
[fe0273c] | 650 | entry->BTagPhys = candidate->BTagPhys;
|
---|
[9040259] | 651 |
|
---|
[d7d2da3] | 652 | entry->TauTag = candidate->TauTag;
|
---|
[7429c6a] | 653 | entry->TauWeight = candidate->TauWeight;
|
---|
[d7d2da3] | 654 |
|
---|
| 655 | entry->Charge = candidate->Charge;
|
---|
| 656 |
|
---|
| 657 | itConstituents.Reset();
|
---|
| 658 | entry->Constituents.Clear();
|
---|
| 659 | ecalEnergy = 0.0;
|
---|
| 660 | hcalEnergy = 0.0;
|
---|
| 661 | while((constituent = static_cast<Candidate*>(itConstituents.Next())))
|
---|
| 662 | {
|
---|
| 663 | entry->Constituents.Add(constituent);
|
---|
| 664 | ecalEnergy += constituent->Eem;
|
---|
| 665 | hcalEnergy += constituent->Ehad;
|
---|
| 666 | }
|
---|
| 667 |
|
---|
| 668 | entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
|
---|
[8f7db23] | 669 |
|
---|
[24d005f] | 670 | //--- Pile-Up Jet ID variables ----
|
---|
| 671 |
|
---|
| 672 | entry->NCharged = candidate->NCharged;
|
---|
| 673 | entry->NNeutrals = candidate->NNeutrals;
|
---|
| 674 | entry->Beta = candidate->Beta;
|
---|
| 675 | entry->BetaStar = candidate->BetaStar;
|
---|
| 676 | entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
|
---|
| 677 | entry->PTD = candidate->PTD;
|
---|
[9040259] | 678 |
|
---|
[de6d698] | 679 | //--- Sub-structure variables ----
|
---|
[9040259] | 680 |
|
---|
| 681 | entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
|
---|
| 682 | entry->NSubJetsPruned = candidate->NSubJetsPruned;
|
---|
[de6d698] | 683 | entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
|
---|
[9040259] | 684 |
|
---|
[ba75867] | 685 | entry->SoftDroppedJet = candidate->SoftDroppedJet ;
|
---|
| 686 | entry->SoftDroppedSubJet1 = candidate->SoftDroppedSubJet1 ;
|
---|
| 687 | entry->SoftDroppedSubJet2 = candidate->SoftDroppedSubJet2;
|
---|
| 688 |
|
---|
| 689 |
|
---|
[de6d698] | 690 | for(i = 0; i < 5; i++)
|
---|
| 691 | {
|
---|
[9040259] | 692 | entry->FracPt[i] = candidate -> FracPt[i];
|
---|
| 693 | entry->Tau[i] = candidate -> Tau[i];
|
---|
| 694 | entry->TrimmedP4[i] = candidate -> TrimmedP4[i];
|
---|
| 695 | entry->PrunedP4[i] = candidate -> PrunedP4[i];
|
---|
[de6d698] | 696 | entry->SoftDroppedP4[i] = candidate -> SoftDroppedP4[i];
|
---|
| 697 | }
|
---|
[9040259] | 698 |
|
---|
[e9c0d73] | 699 | //--- exclusive clustering variables ---
|
---|
| 700 | entry->ExclYmerge23 = candidate->ExclYmerge23;
|
---|
| 701 | entry->ExclYmerge34 = candidate->ExclYmerge34;
|
---|
| 702 | entry->ExclYmerge45 = candidate->ExclYmerge45;
|
---|
| 703 | entry->ExclYmerge56 = candidate->ExclYmerge56;
|
---|
| 704 |
|
---|
| 705 |
|
---|
[d7d2da3] | 706 | FillParticles(candidate, &entry->Particles);
|
---|
| 707 | }
|
---|
| 708 | }
|
---|
| 709 |
|
---|
| 710 | //------------------------------------------------------------------------------
|
---|
| 711 |
|
---|
| 712 | void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 713 | {
|
---|
| 714 | Candidate *candidate = 0;
|
---|
| 715 | MissingET *entry = 0;
|
---|
| 716 |
|
---|
| 717 | // get the first entry
|
---|
| 718 | if((candidate = static_cast<Candidate*>(array->At(0))))
|
---|
| 719 | {
|
---|
| 720 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 721 |
|
---|
| 722 | entry = static_cast<MissingET*>(branch->NewEntry());
|
---|
| 723 |
|
---|
[3b465ca] | 724 | entry->Eta = (-momentum).Eta();
|
---|
[d7d2da3] | 725 | entry->Phi = (-momentum).Phi();
|
---|
| 726 | entry->MET = momentum.Pt();
|
---|
| 727 | }
|
---|
| 728 | }
|
---|
| 729 |
|
---|
| 730 | //------------------------------------------------------------------------------
|
---|
| 731 |
|
---|
| 732 | void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 733 | {
|
---|
| 734 | Candidate *candidate = 0;
|
---|
| 735 | ScalarHT *entry = 0;
|
---|
| 736 |
|
---|
| 737 | // get the first entry
|
---|
| 738 | if((candidate = static_cast<Candidate*>(array->At(0))))
|
---|
| 739 | {
|
---|
| 740 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 741 |
|
---|
| 742 | entry = static_cast<ScalarHT*>(branch->NewEntry());
|
---|
| 743 |
|
---|
| 744 | entry->HT = momentum.Pt();
|
---|
| 745 | }
|
---|
| 746 | }
|
---|
| 747 |
|
---|
| 748 | //------------------------------------------------------------------------------
|
---|
| 749 |
|
---|
[71648c2] | 750 | void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 751 | {
|
---|
[3b465ca] | 752 | TIter iterator(array);
|
---|
[71648c2] | 753 | Candidate *candidate = 0;
|
---|
| 754 | Rho *entry = 0;
|
---|
| 755 |
|
---|
[3b465ca] | 756 | // loop over all rho
|
---|
| 757 | iterator.Reset();
|
---|
| 758 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
[71648c2] | 759 | {
|
---|
| 760 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 761 |
|
---|
| 762 | entry = static_cast<Rho*>(branch->NewEntry());
|
---|
| 763 |
|
---|
[8608ef8] | 764 | entry->Rho = momentum.E();
|
---|
[3b465ca] | 765 | entry->Edges[0] = candidate->Edges[0];
|
---|
| 766 | entry->Edges[1] = candidate->Edges[1];
|
---|
[71648c2] | 767 | }
|
---|
| 768 | }
|
---|
| 769 |
|
---|
| 770 | //------------------------------------------------------------------------------
|
---|
| 771 |
|
---|
[2e229c9] | 772 | void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 773 | {
|
---|
| 774 | Candidate *candidate = 0;
|
---|
| 775 | Weight *entry = 0;
|
---|
| 776 |
|
---|
| 777 | // get the first entry
|
---|
| 778 | if((candidate = static_cast<Candidate*>(array->At(0))))
|
---|
| 779 | {
|
---|
| 780 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 781 |
|
---|
| 782 | entry = static_cast<Weight*>(branch->NewEntry());
|
---|
| 783 |
|
---|
| 784 | entry->Weight = momentum.E();
|
---|
| 785 | }
|
---|
| 786 | }
|
---|
| 787 |
|
---|
| 788 | //------------------------------------------------------------------------------
|
---|
| 789 |
|
---|
[8f7db23] | 790 | void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
|
---|
| 791 | {
|
---|
| 792 | TIter iterator(array);
|
---|
| 793 | Candidate *candidate = 0;
|
---|
| 794 | HectorHit *entry = 0;
|
---|
| 795 |
|
---|
| 796 | // loop over all roman pot hits
|
---|
| 797 | iterator.Reset();
|
---|
| 798 | while((candidate = static_cast<Candidate*>(iterator.Next())))
|
---|
| 799 | {
|
---|
| 800 | const TLorentzVector &position = candidate->Position;
|
---|
| 801 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 802 |
|
---|
| 803 | entry = static_cast<HectorHit*>(branch->NewEntry());
|
---|
| 804 |
|
---|
| 805 | entry->E = momentum.E();
|
---|
| 806 |
|
---|
| 807 | entry->Tx = momentum.Px();
|
---|
| 808 | entry->Ty = momentum.Py();
|
---|
| 809 |
|
---|
| 810 | entry->T = position.T();
|
---|
| 811 |
|
---|
| 812 | entry->X = position.X();
|
---|
| 813 | entry->Y = position.Y();
|
---|
| 814 | entry->S = position.Z();
|
---|
[64a4950] | 815 |
|
---|
| 816 | entry->Particle = candidate->GetCandidates()->At(0);
|
---|
[8f7db23] | 817 | }
|
---|
| 818 | }
|
---|
| 819 |
|
---|
| 820 | //------------------------------------------------------------------------------
|
---|
| 821 |
|
---|
[d7d2da3] | 822 | void TreeWriter::Process()
|
---|
| 823 | {
|
---|
| 824 | TBranchMap::iterator itBranchMap;
|
---|
| 825 | ExRootTreeBranch *branch;
|
---|
| 826 | TProcessMethod method;
|
---|
| 827 | TObjArray *array;
|
---|
| 828 |
|
---|
| 829 | for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
|
---|
| 830 | {
|
---|
| 831 | branch = itBranchMap->first;
|
---|
| 832 | method = itBranchMap->second.first;
|
---|
| 833 | array = itBranchMap->second.second;
|
---|
| 834 |
|
---|
| 835 | (this->*method)(branch, array);
|
---|
| 836 | }
|
---|
| 837 | }
|
---|
| 838 |
|
---|
| 839 | //------------------------------------------------------------------------------
|
---|