Fork me on GitHub

source: git/modules/TreeWriter.cc@ cd75093

ImprovedOutputFile Timing dual_readout llp
Last change on this file since cd75093 was b443089, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix EOL characters in GPLv3 header

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