ROOT logo
// #####    TO RUN THIS MACRO:
// bash$ aliroot         (due to AliInfo, AliError, ...)
// root[0] gSystem->Load("libANALYSIS");
// IN case you do not have the include path set to AliRoot includes:
// root[1] gSystem->AddIncludePath("-I\"$ALICE_ROOT/include\"");
// root[2] .L testEvent.C+;
// root[3] generate();
// root[4] filter_reco();

#include "TClonesArray.h"
#include "TChain.h"
#include "TH1.h"
#include "TMath.h"
#include "TCanvas.h"
#include "testEvent.h"   

//============= First step: generate events
void generate()
{
// Simple event generation
   AliAnalysisManager *mgr = new AliAnalysisManager("generate");
   TaskGenerate *task = new TaskGenerate("gener");
   mgr->AddTask(task);

   if (mgr->InitAnalysis()) {
      mgr->PrintStatus();
      mgr->StartAnalysis();
   }   
   delete mgr;
}   

//============= Second step: filter gammas; use TSelector functionality
void filter_reco()
{
// Filter the input events having more than 100 gammas. Reconstruct pi0's
// From gammas coming from the same vertex.
   // Get the input data as a chain
   TChain *chain = new TChain("T");
   chain->Add("event02000.root");
   chain->Add("event04000.root");
   chain->Add("event06000.root");
   chain->Add("event08000.root");
   chain->Add("event10000.root");
   // Create an analysis manager
   AliAnalysisManager *mgr = new AliAnalysisManager("testEvent");
   // Create a filter task and register it
   TaskFilter *task1 = new TaskFilter("TaskFilter");
   mgr->AddTask(task1);
   // Create a reco task and register it

   TaskRecoPi0 *task2 = new TaskRecoPi0("TaskRecoPi0");
   mgr->AddTask(task2);
   // Create containers for input/output
   AliAnalysisDataContainer *cinput = mgr->CreateContainer("input0", 
                      TTree::Class(), AliAnalysisManager::kInputContainer);
   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("output0", 
                      TTree::Class(), AliAnalysisManager::kOutputContainer);
   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("output1", 
                      TList::Class(), AliAnalysisManager::kOutputContainer);
   AliAnalysisDataContainer *coutput = mgr->CreateContainer("output", 
                      TH1::Class(), AliAnalysisManager::kOutputContainer, "output.root");
   // Connect containers to the task input/outputs
   mgr->ConnectInput(task1,0,cinput);
   mgr->ConnectOutput(task1,0,coutput1);
   mgr->ConnectOutput(task1,1,coutput2);
   mgr->ConnectInput(task2,0,coutput1);
   mgr->ConnectOutput(task2,0,coutput);

   // Connect input data
//   cinput->SetData(chain);
   // Init analysis and start event loop
   if (mgr->InitAnalysis()) {
      mgr->PrintStatus();
      mgr->StartAnalysis("local",chain);
   }
   delete mgr;   
}
   
ClassImp(TaskGenerate)

//______________________________________________________________________________
void TaskGenerate::Exec(Option_t *)
{
// Generate 10k events in 5 files
   AnaEvent::CreateEvents(2000, "event02000.root");
   AnaEvent::CreateEvents(2000, "event04000.root");
   AnaEvent::CreateEvents(2000, "event06000.root");
   AnaEvent::CreateEvents(2000, "event08000.root");
   AnaEvent::CreateEvents(2000, "event10000.root");
}

ClassImp(TaskFilter)

//______________________________________________________________________________
TaskFilter::TaskFilter(const char *name) 
           :AliAnalysisTask(name,""), fEvent(0), fOutput(0), fList(0), fHist1(0), fHist2(0)
{
// Constructor
   // Input slot #0 works with a tree of events
   DefineInput(0, TTree::Class());

   // Output slot #0 writes into a tree, #1 produces a hisogram
   DefineOutput(0, TTree::Class());
   DefineOutput(1, TList::Class());
}

//______________________________________________________________________________
void TaskFilter::ConnectInputData(Option_t *)
{
// Initialize branches.
   printf("   ConnectInputData of task %s\n", GetName());
   char ** address = (char **)GetBranchAddress(0, "event");
   if (address) {
      // One should first check if the branch address was taken by some other task
      fEvent = (AnaEvent*)(*address);
   } else {
      fEvent = new AnaEvent();
      SetBranchAddress(0, "event", &fEvent);
   } 
}
   
//______________________________________________________________________________
void TaskFilter::CreateOutputObjects()
{
   printf("   CreateOutputObjects of task %s\n", GetName());
   if (!fList) {
      fOutput = new TTree("TGAM", "gammas");
      TBranch *branch = fOutput->Branch("event", &fEvent, 18000,1);
      branch->SetAutoDelete(kFALSE);
      fList = new TList();
      fHist1 = new TH1I("ntracks", "Number of tracks per event", 100, 0, 1000);
      fHist1->SetLineColor(kRed);
      fHist2 = new TH1I("ngammas", "Number of gammas per event", 100, 0, 1000);
      fHist2->SetLineColor(kBlue);
      fList->Add(fHist1);
      fList->Add(fHist2);
   }   
}

//______________________________________________________________________________
void TaskFilter::Exec(Option_t *)
{
// Filtering.
   TTree *tinput = (TTree*)GetInputData(0);
   Long64_t ientry = tinput->GetReadEntry();
   // In this case fEvent address is already connected to the input
   if (!fEvent) return;
   // First check multiplicity
   Int_t ntracks = fEvent->GetNtracks();
   // Loop tracks and get rid of non-gammas
   AnaTrack *track;
   Int_t igamma = 0;
   for (Int_t i=0; i<ntracks; i++) {
      track = fEvent->GetTrack(i);
      if (track->GetMass() < 1.e-3) igamma++;
   }
   if (ientry%100 == 0) printf("TaskFilter -> Event %lld: %i tracks filtered %i gammas\n", ientry, ntracks, igamma);
   fHist1->Fill(ntracks);
   fHist2->Fill(igamma);
   if (igamma > 100) {
      fOutput->Fill();
      PostData(0, fOutput);
   }   
   PostData(1, fList);
}

//______________________________________________________________________________
void TaskFilter::Terminate(Option_t *)
{
// Draw some histogram at the end.
   if (!gROOT->IsBatch()) {
      fHist1->SetMaximum(2500);
      fHist1->Draw();
      fHist2->Draw("SAME");
   }   
}

ClassImp(TaskRecoPi0)

//______________________________________________________________________________
TaskRecoPi0::TaskRecoPi0(const char *name) 
           :AliAnalysisTask(name,""), fEvent(0), fGammas(0), fPions(0), fHist(0)
{
// Constructor
   // Input slot #0 works with a tree of events
   DefineInput(0, TTree::Class());

   // Output slot #1 produces a hisogram
   DefineOutput(0, TH1::Class());
}

//______________________________________________________________________________
TaskRecoPi0::~TaskRecoPi0()
{
// Dtor.
   if (fEvent) delete fEvent;
   if (fGammas) delete fGammas;
   if (fPions) delete fPions;
}   

//______________________________________________________________________________
void TaskRecoPi0::ConnectInputData(Option_t *)
{
// Initialize branches.
   printf("   ConnectInputData for task %s\n", GetName());
   if (!fEvent) {
      // One should first check if the branch address was taken by some other task
      char ** address = (char **)GetBranchAddress(0, "event");
      if (address) fEvent = (AnaEvent*)(*address);
      if (!fEvent) {
         fEvent = new AnaEvent();
         SetBranchAddress(0, "event", &fEvent);
      }
   } 
}

//______________________________________________________________________________
void TaskRecoPi0::CreateOutputObjects()
{
   printf("   CreateOutputObjects of task %s\n", GetName());
   if (!fHist) {
      fGammas = new TObjArray();
      fPions  = new TObjArray();
      fHist = new TH1F("Pt_pi0", "Pt distribution for pi0's", 100, 0., 10.);
      fHist->SetLineColor(kRed);
   }   
}

//______________________________________________________________________________
void TaskRecoPi0::Exec(Option_t *)
{
// Reconstruct Pi0's for one event
   AnaTrack *track = 0;
   Int_t ntracks = fEvent->GetNtracks();
   Int_t ngamma = 0;
   Int_t i,j;
   // Clear containers
   fGammas->Clear();
   fPions->Delete();
   fPions->Clear();
   // Loop tracks and move gammas to gamma container
   for (i=0; i<ntracks; i++) {
      track = fEvent->GetTrack(i);
      if (track->GetMass() < 1.e-3) {
         fGammas->Add(track);
         ngamma++;
      }
   }
   printf("TaskRecoPi0 -> Tracks %i \n", ntracks);

   // Loop gammas and check vertex position
   Double_t v1[3], v2[3];
   AnaTrack *tracko = 0;
   Double_t cutoff = 0.001;
   for (i=0; i<ngamma-1; i++) {
      track = (AnaTrack*)fGammas->At(i);
      v1[0] = track->GetVertex(0);
      v1[1] = track->GetVertex(1);
      v1[2] = track->GetVertex(2);
      for (j=i+1; j<ngamma; j++) {
         tracko = (AnaTrack*)fGammas->At(j);
         v2[0] = tracko->GetVertex(0);
         v2[1] = tracko->GetVertex(1);
         v2[2] = tracko->GetVertex(2);
         Double_t dist2 = (v2[0]-v1[0])*(v2[0]-v1[0])+
                          (v2[1]-v1[1])*(v2[1]-v1[1])+
                          (v2[2]-v1[2])*(v2[2]-v1[2]);
         if (dist2>cutoff*cutoff) continue;
         // We have found a pair candidate
         Double_t px = track->GetPx()+tracko->GetPx();
         Double_t py = track->GetPy()+tracko->GetPy();
         Double_t pz = track->GetPz()+tracko->GetPz();
         track = new AnaTrack(px,py,pz,0.135,0,v1[0],v1[1],v1[2]);
         fPions->Add(track);
         fHist->Fill(track->GetPt());
         break;
      }   
   }   
   PostData(0,fHist);         
}

//______________________________________________________________________________
void TaskRecoPi0::Terminate(Option_t *)
{
// Draw some histogram at the end.
   if (!gROOT->IsBatch()) {
      new TCanvas("pi0", "Pt for pi0's", 800,600);
      fHist->Draw();
   }   
}


ClassImp(AnaTrack)

//______________________________________________________________________________
AnaTrack::AnaTrack(Double_t random, Double_t *vertex) 
{
// Constructor
   Int_t itype;
   if (random<0.3) itype = 1;      // pi+
   else if (random<0.6) itype = 2; // pi-
   else if (random<0.9) itype = 3; // p
   else if (random<0.95) itype = 4; // pi0
   else itype = 5;                 // gamma
   gRandom->Rannor(fPx, fPy);
   Double_t vert_width = 0.1;
   
   switch (itype) {
      case 1:
         fMass = 0.13957 + gRandom->Gaus(0.,0.001);
         fCharge = 1;
         break;
      case 2:
         fMass = 0.13957 + gRandom->Gaus(0.,0.001);
         fCharge = -1;
         break;
      case 3:
         fMass = 0.938 + gRandom->Gaus(0.,0.002);
         fCharge = 1;
         fPx *= 0.15;
         fPy *= 0.15;
         break;
      case 4:
         fMass = 0.135 + gRandom->Gaus(0.,0.001);
         fCharge = 0;
         fPx *= 0.8;
         fPy *= 0.8;
         vert_width = 10.;
         break;
      case 5:
         fMass = 0.;
         fCharge = 0;
         fPx *= 0.5;
         fPy *= 1.5;
         break;
   };
   fPz = gRandom->Gaus(4., 2.);
   if (vertex) {
      fVertex[0] = vertex[0];
      fVertex[1] = vertex[1];
      fVertex[2] = vertex[2];
   } else {   
      fVertex[0] = gRandom->Gaus(0,vert_width);
      fVertex[1] = gRandom->Gaus(0,vert_width);
      fVertex[2] = gRandom->Gaus(0,0.01);
   }   
}

//______________________________________________________________________________
Bool_t AnaTrack::Decay(Double_t &px1, Double_t &py1, Double_t &pz1, 
                       Double_t &px2, Double_t &py2, Double_t &pz2)
{
// Decay a pi0 in 2 gammas.
   if (fCharge != 0) return kFALSE;
   if (fMass<0.132 || fMass>0.138) return kFALSE;
   Double_t phi1 = 2.*TMath::Pi()*gRandom->Rndm(); // [0,2*pi]
   Double_t phi2 = phi1 + TMath::Pi();
   if (phi2 > 2.*TMath::Pi()) phi2 -= 2.*TMath::Pi();
   Double_t r2 = gRandom->Rndm();
   Double_t theta1 = TMath::ACos(1.-2.*r2);
   Double_t p0 = GetP();
   Double_t m0 = 0.135;
   Double_t p1 = 0.5*(p0*(1-2*r2)+TMath::Sqrt(p0*p0*(1-2*r2)*(1-2*r2)+2*m0*m0));
   Double_t p2 = TMath::Sqrt(p0*p0+m0*m0-p1*p1);
   Double_t theta2 = TMath::ACos((p0-p1*(1-2*r2))/p2);
   // Px, Py and Pz in the reference frame of the pion
   px1 = p1 * TMath::Sin(theta1)*TMath::Cos(phi1);
   px2 = p2 * TMath::Sin(theta2)*TMath::Cos(phi2);
   py1 = p1 * TMath::Sin(theta1)*TMath::Sin(phi1);
   py2 = p2 * TMath::Sin(theta2)*TMath::Sin(phi2);
   pz1 = p1 * TMath::Cos(theta1);
   pz2 = p2 * TMath::Cos(theta2);
   Double_t phi = TMath::ATan2(fPy, fPx) * TMath::RadToDeg();
   Double_t theta = TMath::ACos(GetPt()/p0) * TMath::RadToDeg();
   TGeoRotation r("rot", phi+90., theta, 0.);
   Double_t loc[3], vect[3];
   p0 = TMath::Sqrt(px1*px1+py1*py1+pz1*pz1);
   loc[0] = px1/p0;
   loc[1] = py1/p0;
   loc[2] = pz1/p0;
   r.LocalToMasterVect(loc, vect);
   px1 = vect[0]*p0;
   py1 = vect[1]*p0;
   pz1 = vect[2]*p0;
//   t1 = new AnaTrack(1., fVertex);
   p0 = TMath::Sqrt(px2*px2+py2*py2+pz2*pz2);
   loc[0] = px2/p0;
   loc[1] = py2/p0;
   loc[2] = pz2/p0;
   r.LocalToMasterVect(loc, vect);
   px2 = vect[0]*p0;
   py2 = vect[1]*p0;
   pz2 = vect[2]*p0;
//   t2 = new AnaTrack(1., fVertex);
   return kTRUE;
}

//______________________________________________________________________________
ClassImp(AnaEvent)

TClonesArray *AnaEvent::fgTracks = 0;

//______________________________________________________________________________
AnaEvent::AnaEvent()
{
// Ctor
   if (!fgTracks) fgTracks = new TClonesArray("AnaTrack", 2000);
   fTracks = fgTracks;
   fEventNumber = 0;
   fNtracks = 0;
}

//______________________________________________________________________________
AnaTrack *AnaEvent::AddTrack(Double_t rnd, Double_t *vert)
{
// Add a random track
//   printf("track %d\n", fNtracks);
   TClonesArray &tracks = *fTracks;
   AnaTrack *track = new(tracks[fNtracks++]) AnaTrack(rnd, vert);
   return track;
}   

//______________________________________________________________________________
void AnaEvent::Clear(Option_t *)
{
// Clears current event.
   fTracks->Clear();
   fNtracks = 0;
   fEventNumber = 0;
}   

//______________________________________________________________________________
Int_t AnaEvent::Build(Int_t ev)
{
// Create a random event
   Clear();
   fEventNumber = ev;
   Int_t ntracks = Int_t(gRandom->Gaus(500., 100.));
   if (ntracks < 1) ntracks = 1;
   if (ntracks>1000) ntracks = 1000;
   AnaTrack *track, *track0;
   for (Int_t i=0; i<ntracks; i++) {
      Double_t rnd = gRandom->Rndm();
      if (rnd>=0.90 && rnd<0.95) {
      // Pi0 -> decay the track in 2 gammas
         track0 = new AnaTrack(rnd);
         Double_t vert[3];
         vert[0] = track0->GetVertex(0);
         vert[1] = track0->GetVertex(1);
         vert[2] = track0->GetVertex(2);
         Double_t px1,py1,pz1,px2,py2,pz2;
         if (track0->Decay(px1,py1,pz1,px2,py2,pz2)) {
            track = AddTrack(1.,vert);
            track->SetPx(px1);
            track->SetPy(py1);
            track->SetPz(pz1);
            track = AddTrack(1.,vert);
            track->SetPx(px2);
            track->SetPy(py2);
            track->SetPz(pz2);
         }
         delete track0;
      } else {   
         track = AddTrack(rnd);
      }   
   }
   return fNtracks;
}   

//______________________________________________________________________________
void AnaEvent::CreateEvents(Int_t nevents, const char *filename)
{
// Create nevents in one tree and put them in filename.
   TFile *hfile = new TFile(filename, "RECREATE", "Some AnaEvents...");
   TTree *tree  = new TTree("T", "Tree of AnaEvents");
   tree->SetAutoSave(1000000000);  // autosave when 1 Gbyte written
   Int_t bufsize = 16000;
   AnaEvent *event = new AnaEvent();
   TBranch *branch = tree->Branch("event", &event, bufsize,1);
   branch->SetAutoDelete(kFALSE);
   
   for (Int_t ev=0; ev<nevents; ev++) {
      Int_t ntracks = event->Build(ev);
      if (ev%100 == 0) printf("event: %d  ntracks=%d\n", ev, ntracks);
      tree->Fill();
   }   
   hfile->Write();
   tree->Print();
   hfile->Close();
}  
   
 testEvent.C:1
 testEvent.C:2
 testEvent.C:3
 testEvent.C:4
 testEvent.C:5
 testEvent.C:6
 testEvent.C:7
 testEvent.C:8
 testEvent.C:9
 testEvent.C:10
 testEvent.C:11
 testEvent.C:12
 testEvent.C:13
 testEvent.C:14
 testEvent.C:15
 testEvent.C:16
 testEvent.C:17
 testEvent.C:18
 testEvent.C:19
 testEvent.C:20
 testEvent.C:21
 testEvent.C:22
 testEvent.C:23
 testEvent.C:24
 testEvent.C:25
 testEvent.C:26
 testEvent.C:27
 testEvent.C:28
 testEvent.C:29
 testEvent.C:30
 testEvent.C:31
 testEvent.C:32
 testEvent.C:33
 testEvent.C:34
 testEvent.C:35
 testEvent.C:36
 testEvent.C:37
 testEvent.C:38
 testEvent.C:39
 testEvent.C:40
 testEvent.C:41
 testEvent.C:42
 testEvent.C:43
 testEvent.C:44
 testEvent.C:45
 testEvent.C:46
 testEvent.C:47
 testEvent.C:48
 testEvent.C:49
 testEvent.C:50
 testEvent.C:51
 testEvent.C:52
 testEvent.C:53
 testEvent.C:54
 testEvent.C:55
 testEvent.C:56
 testEvent.C:57
 testEvent.C:58
 testEvent.C:59
 testEvent.C:60
 testEvent.C:61
 testEvent.C:62
 testEvent.C:63
 testEvent.C:64
 testEvent.C:65
 testEvent.C:66
 testEvent.C:67
 testEvent.C:68
 testEvent.C:69
 testEvent.C:70
 testEvent.C:71
 testEvent.C:72
 testEvent.C:73
 testEvent.C:74
 testEvent.C:75
 testEvent.C:76
 testEvent.C:77
 testEvent.C:78
 testEvent.C:79
 testEvent.C:80
 testEvent.C:81
 testEvent.C:82
 testEvent.C:83
 testEvent.C:84
 testEvent.C:85
 testEvent.C:86
 testEvent.C:87
 testEvent.C:88
 testEvent.C:89
 testEvent.C:90
 testEvent.C:91
 testEvent.C:92
 testEvent.C:93
 testEvent.C:94
 testEvent.C:95
 testEvent.C:96
 testEvent.C:97
 testEvent.C:98
 testEvent.C:99
 testEvent.C:100
 testEvent.C:101
 testEvent.C:102
 testEvent.C:103
 testEvent.C:104
 testEvent.C:105
 testEvent.C:106
 testEvent.C:107
 testEvent.C:108
 testEvent.C:109
 testEvent.C:110
 testEvent.C:111
 testEvent.C:112
 testEvent.C:113
 testEvent.C:114
 testEvent.C:115
 testEvent.C:116
 testEvent.C:117
 testEvent.C:118
 testEvent.C:119
 testEvent.C:120
 testEvent.C:121
 testEvent.C:122
 testEvent.C:123
 testEvent.C:124
 testEvent.C:125
 testEvent.C:126
 testEvent.C:127
 testEvent.C:128
 testEvent.C:129
 testEvent.C:130
 testEvent.C:131
 testEvent.C:132
 testEvent.C:133
 testEvent.C:134
 testEvent.C:135
 testEvent.C:136
 testEvent.C:137
 testEvent.C:138
 testEvent.C:139
 testEvent.C:140
 testEvent.C:141
 testEvent.C:142
 testEvent.C:143
 testEvent.C:144
 testEvent.C:145
 testEvent.C:146
 testEvent.C:147
 testEvent.C:148
 testEvent.C:149
 testEvent.C:150
 testEvent.C:151
 testEvent.C:152
 testEvent.C:153
 testEvent.C:154
 testEvent.C:155
 testEvent.C:156
 testEvent.C:157
 testEvent.C:158
 testEvent.C:159
 testEvent.C:160
 testEvent.C:161
 testEvent.C:162
 testEvent.C:163
 testEvent.C:164
 testEvent.C:165
 testEvent.C:166
 testEvent.C:167
 testEvent.C:168
 testEvent.C:169
 testEvent.C:170
 testEvent.C:171
 testEvent.C:172
 testEvent.C:173
 testEvent.C:174
 testEvent.C:175
 testEvent.C:176
 testEvent.C:177
 testEvent.C:178
 testEvent.C:179
 testEvent.C:180
 testEvent.C:181
 testEvent.C:182
 testEvent.C:183
 testEvent.C:184
 testEvent.C:185
 testEvent.C:186
 testEvent.C:187
 testEvent.C:188
 testEvent.C:189
 testEvent.C:190
 testEvent.C:191
 testEvent.C:192
 testEvent.C:193
 testEvent.C:194
 testEvent.C:195
 testEvent.C:196
 testEvent.C:197
 testEvent.C:198
 testEvent.C:199
 testEvent.C:200
 testEvent.C:201
 testEvent.C:202
 testEvent.C:203
 testEvent.C:204
 testEvent.C:205
 testEvent.C:206
 testEvent.C:207
 testEvent.C:208
 testEvent.C:209
 testEvent.C:210
 testEvent.C:211
 testEvent.C:212
 testEvent.C:213
 testEvent.C:214
 testEvent.C:215
 testEvent.C:216
 testEvent.C:217
 testEvent.C:218
 testEvent.C:219
 testEvent.C:220
 testEvent.C:221
 testEvent.C:222
 testEvent.C:223
 testEvent.C:224
 testEvent.C:225
 testEvent.C:226
 testEvent.C:227
 testEvent.C:228
 testEvent.C:229
 testEvent.C:230
 testEvent.C:231
 testEvent.C:232
 testEvent.C:233
 testEvent.C:234
 testEvent.C:235
 testEvent.C:236
 testEvent.C:237
 testEvent.C:238
 testEvent.C:239
 testEvent.C:240
 testEvent.C:241
 testEvent.C:242
 testEvent.C:243
 testEvent.C:244
 testEvent.C:245
 testEvent.C:246
 testEvent.C:247
 testEvent.C:248
 testEvent.C:249
 testEvent.C:250
 testEvent.C:251
 testEvent.C:252
 testEvent.C:253
 testEvent.C:254
 testEvent.C:255
 testEvent.C:256
 testEvent.C:257
 testEvent.C:258
 testEvent.C:259
 testEvent.C:260
 testEvent.C:261
 testEvent.C:262
 testEvent.C:263
 testEvent.C:264
 testEvent.C:265
 testEvent.C:266
 testEvent.C:267
 testEvent.C:268
 testEvent.C:269
 testEvent.C:270
 testEvent.C:271
 testEvent.C:272
 testEvent.C:273
 testEvent.C:274
 testEvent.C:275
 testEvent.C:276
 testEvent.C:277
 testEvent.C:278
 testEvent.C:279
 testEvent.C:280
 testEvent.C:281
 testEvent.C:282
 testEvent.C:283
 testEvent.C:284
 testEvent.C:285
 testEvent.C:286
 testEvent.C:287
 testEvent.C:288
 testEvent.C:289
 testEvent.C:290
 testEvent.C:291
 testEvent.C:292
 testEvent.C:293
 testEvent.C:294
 testEvent.C:295
 testEvent.C:296
 testEvent.C:297
 testEvent.C:298
 testEvent.C:299
 testEvent.C:300
 testEvent.C:301
 testEvent.C:302
 testEvent.C:303
 testEvent.C:304
 testEvent.C:305
 testEvent.C:306
 testEvent.C:307
 testEvent.C:308
 testEvent.C:309
 testEvent.C:310
 testEvent.C:311
 testEvent.C:312
 testEvent.C:313
 testEvent.C:314
 testEvent.C:315
 testEvent.C:316
 testEvent.C:317
 testEvent.C:318
 testEvent.C:319
 testEvent.C:320
 testEvent.C:321
 testEvent.C:322
 testEvent.C:323
 testEvent.C:324
 testEvent.C:325
 testEvent.C:326
 testEvent.C:327
 testEvent.C:328
 testEvent.C:329
 testEvent.C:330
 testEvent.C:331
 testEvent.C:332
 testEvent.C:333
 testEvent.C:334
 testEvent.C:335
 testEvent.C:336
 testEvent.C:337
 testEvent.C:338
 testEvent.C:339
 testEvent.C:340
 testEvent.C:341
 testEvent.C:342
 testEvent.C:343
 testEvent.C:344
 testEvent.C:345
 testEvent.C:346
 testEvent.C:347
 testEvent.C:348
 testEvent.C:349
 testEvent.C:350
 testEvent.C:351
 testEvent.C:352
 testEvent.C:353
 testEvent.C:354
 testEvent.C:355
 testEvent.C:356
 testEvent.C:357
 testEvent.C:358
 testEvent.C:359
 testEvent.C:360
 testEvent.C:361
 testEvent.C:362
 testEvent.C:363
 testEvent.C:364
 testEvent.C:365
 testEvent.C:366
 testEvent.C:367
 testEvent.C:368
 testEvent.C:369
 testEvent.C:370
 testEvent.C:371
 testEvent.C:372
 testEvent.C:373
 testEvent.C:374
 testEvent.C:375
 testEvent.C:376
 testEvent.C:377
 testEvent.C:378
 testEvent.C:379
 testEvent.C:380
 testEvent.C:381
 testEvent.C:382
 testEvent.C:383
 testEvent.C:384
 testEvent.C:385
 testEvent.C:386
 testEvent.C:387
 testEvent.C:388
 testEvent.C:389
 testEvent.C:390
 testEvent.C:391
 testEvent.C:392
 testEvent.C:393
 testEvent.C:394
 testEvent.C:395
 testEvent.C:396
 testEvent.C:397
 testEvent.C:398
 testEvent.C:399
 testEvent.C:400
 testEvent.C:401
 testEvent.C:402
 testEvent.C:403
 testEvent.C:404
 testEvent.C:405
 testEvent.C:406
 testEvent.C:407
 testEvent.C:408
 testEvent.C:409
 testEvent.C:410
 testEvent.C:411
 testEvent.C:412
 testEvent.C:413
 testEvent.C:414
 testEvent.C:415
 testEvent.C:416
 testEvent.C:417
 testEvent.C:418
 testEvent.C:419
 testEvent.C:420
 testEvent.C:421
 testEvent.C:422
 testEvent.C:423
 testEvent.C:424
 testEvent.C:425
 testEvent.C:426
 testEvent.C:427
 testEvent.C:428
 testEvent.C:429
 testEvent.C:430
 testEvent.C:431
 testEvent.C:432
 testEvent.C:433
 testEvent.C:434
 testEvent.C:435
 testEvent.C:436
 testEvent.C:437
 testEvent.C:438
 testEvent.C:439
 testEvent.C:440
 testEvent.C:441
 testEvent.C:442
 testEvent.C:443
 testEvent.C:444
 testEvent.C:445
 testEvent.C:446
 testEvent.C:447
 testEvent.C:448
 testEvent.C:449
 testEvent.C:450
 testEvent.C:451
 testEvent.C:452
 testEvent.C:453
 testEvent.C:454
 testEvent.C:455
 testEvent.C:456
 testEvent.C:457
 testEvent.C:458
 testEvent.C:459
 testEvent.C:460
 testEvent.C:461
 testEvent.C:462
 testEvent.C:463
 testEvent.C:464
 testEvent.C:465
 testEvent.C:466
 testEvent.C:467
 testEvent.C:468
 testEvent.C:469
 testEvent.C:470
 testEvent.C:471
 testEvent.C:472
 testEvent.C:473
 testEvent.C:474
 testEvent.C:475
 testEvent.C:476
 testEvent.C:477
 testEvent.C:478
 testEvent.C:479
 testEvent.C:480
 testEvent.C:481
 testEvent.C:482
 testEvent.C:483
 testEvent.C:484
 testEvent.C:485
 testEvent.C:486
 testEvent.C:487
 testEvent.C:488
 testEvent.C:489
 testEvent.C:490
 testEvent.C:491
 testEvent.C:492