Fork me on GitHub

source: git/modules/TreeWriter.cc@ 4e09c3a

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

remove svn tags and fix formatting

  • Property mode set to 100644
File size: 19.5 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 * \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
53using namespace std;
54
55//------------------------------------------------------------------------------
56
57TreeWriter::TreeWriter()
58{
59}
60
61//------------------------------------------------------------------------------
62
63TreeWriter::~TreeWriter()
64{
65}
66
67//------------------------------------------------------------------------------
68
69void TreeWriter::Init()
70{
71 fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
72 fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
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;
81 fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
82 fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
83 fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit;
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 }
112
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
130void TreeWriter::Finish()
131{
132}
133
134//------------------------------------------------------------------------------
135
136void 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
171void 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;
177
178 const Double_t c_light = 2.99792458E8;
179
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
217 entry->Eta = eta;
218 entry->Phi = momentum.Phi();
219 entry->PT = pt;
220
221 entry->Rapidity = rapidity;
222
223 entry->X = position.X();
224 entry->Y = position.Y();
225 entry->Z = position.Z();
226 entry->T = position.T()*1.0E-3/c_light;
227 }
228}
229
230//------------------------------------------------------------------------------
231
232void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
233{
234 TIter iterator(array);
235 Candidate *candidate = 0;
236 Vertex *entry = 0;
237
238 const Double_t c_light = 2.99792458E8;
239
240 // loop over all vertices
241 iterator.Reset();
242 while((candidate = static_cast<Candidate*>(iterator.Next())))
243 {
244 const TLorentzVector &position = candidate->Position;
245
246 entry = static_cast<Vertex*>(branch->NewEntry());
247
248 entry->X = position.X();
249 entry->Y = position.Y();
250 entry->Z = position.Z();
251 entry->T = position.T()*1.0E-3/c_light;
252 }
253}
254
255//------------------------------------------------------------------------------
256
257void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
258{
259 TIter iterator(array);
260 Candidate *candidate = 0;
261 Candidate *particle = 0;
262 Track *entry = 0;
263 Double_t pt, signz, cosTheta, eta, rapidity;
264 const Double_t c_light = 2.99792458E8;
265
266 // loop over all tracks
267 iterator.Reset();
268 while((candidate = static_cast<Candidate*>(iterator.Next())))
269 {
270 const TLorentzVector &position = candidate->Position;
271
272 cosTheta = TMath::Abs(position.CosTheta());
273 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
274 eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
275 rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
276
277 entry = static_cast<Track*>(branch->NewEntry());
278
279 entry->SetBit(kIsReferenced);
280 entry->SetUniqueID(candidate->GetUniqueID());
281
282 entry->PID = candidate->PID;
283
284 entry->Charge = candidate->Charge;
285
286 entry->EtaOuter = eta;
287 entry->PhiOuter = position.Phi();
288
289 entry->XOuter = position.X();
290 entry->YOuter = position.Y();
291 entry->ZOuter = position.Z();
292 entry->TOuter = position.T()*1.0E-3/c_light;
293
294 entry->Dxy = candidate->Dxy;
295 entry->SDxy = candidate->SDxy ;
296 entry->Xd = candidate->Xd;
297 entry->Yd = candidate->Yd;
298 entry->Zd = candidate->Zd;
299
300 const TLorentzVector &momentum = candidate->Momentum;
301
302 pt = momentum.Pt();
303 cosTheta = TMath::Abs(momentum.CosTheta());
304 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
305 eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
306 rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
307
308 entry->Eta = eta;
309 entry->Phi = momentum.Phi();
310 entry->PT = pt;
311
312 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
313 const TLorentzVector &initialPosition = particle->Position;
314
315 entry->X = initialPosition.X();
316 entry->Y = initialPosition.Y();
317 entry->Z = initialPosition.Z();
318 entry->T = initialPosition.T()*1.0E-3/c_light;
319
320 entry->Particle = particle;
321 }
322}
323
324//------------------------------------------------------------------------------
325
326void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
327{
328 TIter iterator(array);
329 Candidate *candidate = 0;
330 Tower *entry = 0;
331 Double_t pt, signPz, cosTheta, eta, rapidity;
332 const Double_t c_light = 2.99792458E8;
333
334 // loop over all towers
335 iterator.Reset();
336 while((candidate = static_cast<Candidate*>(iterator.Next())))
337 {
338 const TLorentzVector &momentum = candidate->Momentum;
339 const TLorentzVector &position = candidate->Position;
340
341 pt = momentum.Pt();
342 cosTheta = TMath::Abs(momentum.CosTheta());
343 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
344 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
345 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
346
347 entry = static_cast<Tower*>(branch->NewEntry());
348
349 entry->SetBit(kIsReferenced);
350 entry->SetUniqueID(candidate->GetUniqueID());
351
352 entry->Eta = eta;
353 entry->Phi = momentum.Phi();
354 entry->ET = pt;
355 entry->E = momentum.E();
356 entry->Eem = candidate->Eem;
357 entry->Ehad = candidate->Ehad;
358 entry->Edges[0] = candidate->Edges[0];
359 entry->Edges[1] = candidate->Edges[1];
360 entry->Edges[2] = candidate->Edges[2];
361 entry->Edges[3] = candidate->Edges[3];
362
363 entry->T = position.T()*1.0E-3/c_light;
364
365 FillParticles(candidate, &entry->Particles);
366 }
367}
368
369//------------------------------------------------------------------------------
370
371void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
372{
373 TIter iterator(array);
374 Candidate *candidate = 0;
375 Photon *entry = 0;
376 Double_t pt, signPz, cosTheta, eta, rapidity;
377 const Double_t c_light = 2.99792458E8;
378
379 array->Sort();
380
381 // loop over all photons
382 iterator.Reset();
383 while((candidate = static_cast<Candidate*>(iterator.Next())))
384 {
385 TIter it1(candidate->GetCandidates());
386 const TLorentzVector &momentum = candidate->Momentum;
387 const TLorentzVector &position = candidate->Position;
388
389
390 pt = momentum.Pt();
391 cosTheta = TMath::Abs(momentum.CosTheta());
392 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
393 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
394 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
395
396 entry = static_cast<Photon*>(branch->NewEntry());
397
398 entry->Eta = eta;
399 entry->Phi = momentum.Phi();
400 entry->PT = pt;
401 entry->E = momentum.E();
402
403 entry->T = position.T()*1.0E-3/c_light;
404
405 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
406
407 FillParticles(candidate, &entry->Particles);
408 }
409}
410
411//------------------------------------------------------------------------------
412
413void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
414{
415 TIter iterator(array);
416 Candidate *candidate = 0;
417 Electron *entry = 0;
418 Double_t pt, signPz, cosTheta, eta, rapidity;
419 const Double_t c_light = 2.99792458E8;
420
421 array->Sort();
422
423 // loop over all electrons
424 iterator.Reset();
425 while((candidate = static_cast<Candidate*>(iterator.Next())))
426 {
427 const TLorentzVector &momentum = candidate->Momentum;
428 const TLorentzVector &position = candidate->Position;
429
430 pt = momentum.Pt();
431 cosTheta = TMath::Abs(momentum.CosTheta());
432 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
433 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
434 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
435
436 entry = static_cast<Electron*>(branch->NewEntry());
437
438 entry->Eta = eta;
439 entry->Phi = momentum.Phi();
440 entry->PT = pt;
441
442 entry->T = position.T()*1.0E-3/c_light;
443
444 entry->Charge = candidate->Charge;
445
446 entry->EhadOverEem = 0.0;
447
448 entry->Particle = candidate->GetCandidates()->At(0);
449 }
450}
451
452//------------------------------------------------------------------------------
453
454void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
455{
456 TIter iterator(array);
457 Candidate *candidate = 0;
458 Muon *entry = 0;
459 Double_t pt, signPz, cosTheta, eta, rapidity;
460
461 const Double_t c_light = 2.99792458E8;
462
463 array->Sort();
464
465 // loop over all muons
466 iterator.Reset();
467 while((candidate = static_cast<Candidate*>(iterator.Next())))
468 {
469 const TLorentzVector &momentum = candidate->Momentum;
470 const TLorentzVector &position = candidate->Position;
471
472
473 pt = momentum.Pt();
474 cosTheta = TMath::Abs(momentum.CosTheta());
475 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
476 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
477 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
478
479 entry = static_cast<Muon*>(branch->NewEntry());
480
481 entry->SetBit(kIsReferenced);
482 entry->SetUniqueID(candidate->GetUniqueID());
483
484 entry->Eta = eta;
485 entry->Phi = momentum.Phi();
486 entry->PT = pt;
487
488 entry->T = position.T()*1.0E-3/c_light;
489
490 entry->Charge = candidate->Charge;
491
492 entry->Particle = candidate->GetCandidates()->At(0);
493 }
494}
495
496//------------------------------------------------------------------------------
497
498void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
499{
500 TIter iterator(array);
501 Candidate *candidate = 0, *constituent = 0;
502 Jet *entry = 0;
503 Double_t pt, signPz, cosTheta, eta, rapidity;
504 Double_t ecalEnergy, hcalEnergy;
505 const Double_t c_light = 2.99792458E8;
506
507 array->Sort();
508
509 // loop over all jets
510 iterator.Reset();
511 while((candidate = static_cast<Candidate*>(iterator.Next())))
512 {
513 TIter itConstituents(candidate->GetCandidates());
514
515 const TLorentzVector &momentum = candidate->Momentum;
516 const TLorentzVector &position = candidate->Position;
517
518 pt = momentum.Pt();
519 cosTheta = TMath::Abs(momentum.CosTheta());
520 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
521 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
522 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
523
524 entry = static_cast<Jet*>(branch->NewEntry());
525
526 entry->Eta = eta;
527 entry->Phi = momentum.Phi();
528 entry->PT = pt;
529
530 entry->T = position.T()*1.0E-3/c_light;
531
532 entry->Mass = momentum.M();
533
534 entry->DeltaEta = candidate->DeltaEta;
535 entry->DeltaPhi = candidate->DeltaPhi;
536
537 entry->BTag = candidate->BTag;
538 entry->TauTag = candidate->TauTag;
539
540 entry->Charge = candidate->Charge;
541
542 itConstituents.Reset();
543 entry->Constituents.Clear();
544 ecalEnergy = 0.0;
545 hcalEnergy = 0.0;
546 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
547 {
548 entry->Constituents.Add(constituent);
549 ecalEnergy += constituent->Eem;
550 hcalEnergy += constituent->Ehad;
551 }
552
553 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
554
555 //--- Pile-Up Jet ID variables ----
556
557 entry->NCharged = candidate->NCharged;
558 entry->NNeutrals = candidate->NNeutrals;
559 entry->Beta = candidate->Beta;
560 entry->BetaStar = candidate->BetaStar;
561 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
562 entry->PTD = candidate->PTD;
563 entry->FracPt[0] = candidate->FracPt[0];
564 entry->FracPt[1] = candidate->FracPt[1];
565 entry->FracPt[2] = candidate->FracPt[2];
566 entry->FracPt[3] = candidate->FracPt[3];
567 entry->FracPt[4] = candidate->FracPt[4];
568
569 //--- N-subjettiness variables ----
570
571 entry->Tau1 = candidate->Tau[0];
572 entry->Tau2 = candidate->Tau[1];
573 entry->Tau3 = candidate->Tau[2];
574 entry->Tau4 = candidate->Tau[3];
575 entry->Tau5 = candidate->Tau[4];
576
577 FillParticles(candidate, &entry->Particles);
578 }
579}
580
581//------------------------------------------------------------------------------
582
583void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
584{
585 Candidate *candidate = 0;
586 MissingET *entry = 0;
587
588 // get the first entry
589 if((candidate = static_cast<Candidate*>(array->At(0))))
590 {
591 const TLorentzVector &momentum = candidate->Momentum;
592
593 entry = static_cast<MissingET*>(branch->NewEntry());
594
595 entry->Eta = (-momentum).Eta();
596 entry->Phi = (-momentum).Phi();
597 entry->MET = momentum.Pt();
598 }
599}
600
601//------------------------------------------------------------------------------
602
603void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
604{
605 Candidate *candidate = 0;
606 ScalarHT *entry = 0;
607
608 // get the first entry
609 if((candidate = static_cast<Candidate*>(array->At(0))))
610 {
611 const TLorentzVector &momentum = candidate->Momentum;
612
613 entry = static_cast<ScalarHT*>(branch->NewEntry());
614
615 entry->HT = momentum.Pt();
616 }
617}
618
619//------------------------------------------------------------------------------
620
621void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
622{
623 TIter iterator(array);
624 Candidate *candidate = 0;
625 Rho *entry = 0;
626
627 // loop over all rho
628 iterator.Reset();
629 while((candidate = static_cast<Candidate*>(iterator.Next())))
630 {
631 const TLorentzVector &momentum = candidate->Momentum;
632
633 entry = static_cast<Rho*>(branch->NewEntry());
634
635 entry->Rho = momentum.E();
636 entry->Edges[0] = candidate->Edges[0];
637 entry->Edges[1] = candidate->Edges[1];
638 }
639}
640
641//------------------------------------------------------------------------------
642
643void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
644{
645 Candidate *candidate = 0;
646 Weight *entry = 0;
647
648 // get the first entry
649 if((candidate = static_cast<Candidate*>(array->At(0))))
650 {
651 const TLorentzVector &momentum = candidate->Momentum;
652
653 entry = static_cast<Weight*>(branch->NewEntry());
654
655 entry->Weight = momentum.E();
656 }
657}
658
659//------------------------------------------------------------------------------
660
661void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
662{
663 TIter iterator(array);
664 Candidate *candidate = 0;
665 HectorHit *entry = 0;
666
667 // loop over all roman pot hits
668 iterator.Reset();
669 while((candidate = static_cast<Candidate*>(iterator.Next())))
670 {
671 const TLorentzVector &position = candidate->Position;
672 const TLorentzVector &momentum = candidate->Momentum;
673
674 entry = static_cast<HectorHit*>(branch->NewEntry());
675
676 entry->E = momentum.E();
677
678 entry->Tx = momentum.Px();
679 entry->Ty = momentum.Py();
680
681 entry->T = position.T();
682
683 entry->X = position.X();
684 entry->Y = position.Y();
685 entry->S = position.Z();
686
687 entry->Particle = candidate->GetCandidates()->At(0);
688 }
689}
690
691//------------------------------------------------------------------------------
692
693void TreeWriter::Process()
694{
695 TBranchMap::iterator itBranchMap;
696 ExRootTreeBranch *branch;
697 TProcessMethod method;
698 TObjArray *array;
699
700 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
701 {
702 branch = itBranchMap->first;
703 method = itBranchMap->second.first;
704 array = itBranchMap->second.second;
705
706 (this->*method)(branch, array);
707 }
708}
709
710//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.