ROOT logo
#include "hrecevent.h"
#include "hades.h"
#include "haddef.h"
#include "htree.h"
#include "htrack.h"
#include "hpartialevent.h"
#include "hdebug.h"
#include "TBrowser.h"
#include "TBranchClones.h"

//*-- Author : Manuel Sanchez 
//*-- Modified : 03/05/2002 by R. Holzmann
//*-- Modified : 10/01/2001 by Manuel Sanchez
//*-- Modified : 09/03/2000 by Denis Bertini
//*-- Modified : 12/05/1999 by Manuel Sanchez
//*-- Modified : 02/11/98 by Manuel Sanchez
//*-- Modified : 24/03/98 by Manuel Sanchez
//*-- Copyright : GENP (Univ. Santiago de Compostela)

//******************************************************************************
//
// HRecEvent
// 
// A HRecEvent is an event under reconstruction, and, in particular, a 
// reconstructed event. These events can be in different reconstruction levels,
// the reconstruction level of a HRecEvent is controlled by setRecLevel() and
// getRecLevel(); the reconstruction levels are identified by global constants
// (rlRaw, rlHit ...)
//
// As for the data, a HRecEvent holds reconstructed tracks, an event header
// and a list of HPartialEvent objects.
//
// An HRecEvent can be streamed in "expanded" or normal mode. This behaviour is 
// controlled with the setExpandedStreamer method.
//
// fBits assignment:
//   Bit 0: Expanded Streamer
//
//*******************************************************************************

ClassImp(HRecEvent)

HRecEvent::HRecEvent(void) : HEvent("Event","Event under Reconstruction") {
  // Default constructor
  fRecLevel=rlUndefined;
  fTracks=new TClonesArray("HTrack");
  fPartialEvs=new TObjArray(MAX_PART_EVENTS);
  fHeader = new HEventHeader(); 
  fNTracks=0;
  setExpandedStreamer(kFALSE);
}

HRecEvent::~HRecEvent(void) {
  // Destructor
  
  if (fPartialEvs) {fPartialEvs->Delete(); delete fPartialEvs;}
  if (fTracks) {fTracks->Clear(); delete fTracks;}
  if (fHeader) {delete fHeader; fHeader=0;}
}

void HRecEvent::activateBranch(TTree *tree,Int_t splitLevel) {
  // See HEvent::activateBranch()
   HPartialEvent **par;
#if DEBUG_LEVEL > 2
   gDebuger->enterFunc("HRecEvent::activateBranch");
#endif   
   if (tree==NULL) return;

   if (splitLevel==0) {
   } else if (splitLevel==1) {

     Int_t i;
     TObjArray &vector=*fPartialEvs;
     for(i=0;i<=fPartialEvs->GetLast();i++) {
       par=(HPartialEvent **)&vector[i];
       if ((*par)!=NULL) {
#if DEBUG_LEVEL>2
	 gDebuger->message("Activating %s",(*par)->GetName());
	 gDebuger->message("fPartialEvs=%p",fPartialEvs);
#endif
	 tree->SetBranchAddress((*par)->GetName(),par);
	 tree->SetBranchStatus((*par)->GetName(),1);
       }
     }

   } else if (splitLevel==2) {
         
     Int_t i;
     TObjArray &vector=*fPartialEvs;
     TString cad;
     TBranch *brHead=0;
     TBranch *brPartial=0;

     for(i=0;i<=fPartialEvs->GetLast();i++) {
       par=(HPartialEvent **)&(vector[i]);
       if ((*par)!=NULL) {
	 cad=(*par)->GetName();
	 cad+=".";
	 brPartial = tree->GetBranch(cad.Data());
	 if (brPartial) { //Post ROOT3 splitting
	   tree->SetBranchAddress(cad.Data(),par);
	   tree->SetBranchStatus(cad.Data(),1);
	   cad+="*"; //Activate subbranches
	   tree->SetBranchStatus(cad.Data(),1);
	   (*par)->activateBranch(tree,splitLevel);
	 } else {
	   cad=(*par)->GetName();
	   brPartial = tree->GetBranch(cad.Data());
	   if (brPartial) { //Pre ROOT3 Splitting
	     tree->SetBranchAddress(cad.Data(),par);
	     tree->SetBranchStatus(cad.Data(),1);
	     cad+=".*"; //Activate subbranches
	     tree->SetBranchStatus(cad.Data(),1);
	     (*par)->activateBranch(tree,splitLevel);	     
	   } else {
	     Warning("activateBranch","Partial event %s not found in tree",
		     (*par)->GetName());
	   }
	 }
       }
     }

     brHead=tree->GetBranch("EventHeader.");
     if (brHead) { //Post ROOT3 splitting
       tree->SetBranchAddress("EventHeader.",&fHeader);
       tree->SetBranchStatus("EventHeader.",1);
       tree->SetBranchStatus("EventHeader.*",1);
     } else {
       brHead=tree->GetBranch("EventHeader");
       if (brHead) { //Pre ROOT3 splitting
         tree->SetBranchAddress("EventHeader",&fHeader);
         tree->SetBranchStatus("EventHeader",1);
         tree->SetBranchStatus("EventHeader.*",1);
       }
     }
   }
#if DEBUG_LEVEL>2
   gDebuger->leaveFunc("HRecEvent::activateBranch");
#endif
}

void HRecEvent::makeBranch(TBranch *parent, HTree* tree) {
  // See HEvent::makeBranch() 
   HPartialEvent **par;
   TBranch *b=NULL;
   Int_t sl=2;
   Char_t buffer[256];
   
#if DEBUG_LEVEL > 2
   gDebuger->enterFunc("HRecEvent::makeBranch");
#endif
   if (parent==NULL) return;

   sl=gHades->getSplitLevel();
   if (sl==0) {
   } else if (sl==1) {

     Int_t i;
     TObjArray &vector=*fPartialEvs;
     for(i=0;i<=fPartialEvs->GetLast();i++) {
       par=(HPartialEvent **)&vector[i];
       if ((*par)!=NULL) {
	 b=tree->Branch( (*par)->GetName(),
			(*par)->ClassName(),par,gHades->getTreeBufferSize(),0);
       }
     }
  
    } else if (sl==2) {

     Int_t i;
     TObjArray &vector=*fPartialEvs;
     for(i=0;i<=fPartialEvs->GetLast();i++) {
       par=(HPartialEvent **)&(vector[i]);
       if ((*par)!=NULL) {
	 sprintf(buffer,"%s.",(*par)->GetName());
	 b=tree->Branch( buffer,
			(*par)->ClassName(),par,gHades->getTreeBufferSize(),99);
	 (*par)->makeBranch(b,tree);
       }
     }
     tree->Branch("EventHeader.",fHeader->ClassName(),
		  &fHeader,gHades->getTreeBufferSize(),3);
   }
#if DEBUG_LEVEL>2
   gDebuger->leaveFunc("HRecEvent::activateBranch");
#endif
}

void HRecEvent::Clear(Option_t *) {
  // Clears the data in the event (i.e. clears the internal buffers...)
     TIter next(fPartialEvs);
     HPartialEvent *ev;
     while ( (ev=(HPartialEvent *)next())!=NULL) {
              ev->Clear();
     }
   if (fTracks) fTracks->Clear();
   fNTracks=0;
}

void HRecEvent::clearAll(Int_t level) {
  // Clears the data in the event and the event structure (list of subevents....)
  //
  //   level == 0: delete list of partial events
  //   level == 1: delete only all categories in partial events
  //   level > 1:  clear only all categories in partial events
  //
  if (level==0) {
    if (fPartialEvs) fPartialEvs->Delete();
  } else {
    TIter next(fPartialEvs);
    HPartialEvent *ev;
    while ((ev=(HPartialEvent *)next())!=NULL) {
      ev->clearAll(level);
    }
  }
  if (fTracks) fTracks->Clear();
  fNTracks=0;
}

HTrack *HRecEvent::newTrack(void) {
  // Returns a pointer to a new HTrack object.
  TClonesArray &tracks=*fTracks;
  new(tracks[fNTracks++]) HTrack;
  return ((HTrack *)tracks[fNTracks-1]);
}

void HRecEvent::addTrack(HTrack &aTrack) {
  // Adds the track aTrack to the list of reconstructed tracks.
  TClonesArray &tracks=*fTracks;
  new (tracks[fNTracks++]) HTrack(aTrack);

}

HTrack *HRecEvent::getTrack(UInt_t aId) {
  // Returns the track identified by aId (the position in the track list) 
  HTrack *track;
  track=(HTrack *)fTracks->At(aId);
  return(track);
}

void HRecEvent::clearTracks(void) {
  // Clears the track list
  fNTracks=0;
  if (fTracks) fTracks->Clear();
}

Int_t HRecEvent::getRecLevel(void) {
  // Returns the reconstruction level for this event.
 return(fRecLevel);
}

void HRecEvent::setRecLevel(Int_t aRecLevel) {
  // Sets the reconstruction level for the event.
  fRecLevel=aRecLevel;
}

HPartialEvent *HRecEvent::getPartialEvent(Cat_t idx) {
  // Returns a pointer to the partial event with number idx.
  if ((idx>>kBitCategorySize)>fPartialEvs->GetLast()) return NULL;
  return ((HPartialEvent *)fPartialEvs->At(idx>>kBitCategorySize));
}


HCategory *HRecEvent::getCategory(Cat_t aCat) {
  // Returns the category identified by aCat in the correct Partial event.
  HPartialEvent *ev;
  HCategory *cat;
  ev=getPartialEvent(aCat);
  if (!ev) return NULL;
  cat=ev->getCategory(aCat);
  return cat;
}

Bool_t HRecEvent::addCategory(Cat_t aCat,HCategory *cat,Option_t opt[]) {
  // Adds a new category to the event. The partial event it belongs to is
  // determined by aCat; if this partial event doesn't exist, one is created
  // with the name given in opt.
  HPartialEvent *event=0;
  event=getPartialEvent(aCat);
  if (!event) {
    Int_t aCatBase = aCat - aCat%kCategorySize;    // extract base category
    event=addPartialEvent(aCatBase,opt,opt);
    if (!event) return kFALSE;
  }
  return event->addCategory(aCat,cat,opt);
}


Bool_t HRecEvent::removeCategory(Cat_t aCat) {
  // Removes the category aCat from the event. If aCat is the last category in
  // the corresponding partial event, then the partial event is also removed.
  HPartialEvent *event=NULL;
  event=getPartialEvent(aCat);
  if (!event) return kFALSE;
  if (!event->removeCategory(aCat)) return kFALSE;
  if (event->isEmpty()) {
    fPartialEvs->Remove(event); // May be optimized to use the index instead of the pointer
    delete event;
  }
  return kTRUE;
}

Bool_t HRecEvent::removePartialEvent(Cat_t aCat) {
  // Remove partial event aCat from event.
  HPartialEvent *partial=getPartialEvent(aCat);
  if (!partial) return kFALSE;
  fPartialEvs->Remove(partial);
  delete partial;
  return kTRUE;
}

HPartialEvent *HRecEvent::addPartialEvent(Cat_t eventCat,
					  const Text_t *name,
					  const Text_t *title) {
  // Create and Add a new HPartialEvent to the list of HPartialEvent objects
  // in the HRecEvent.
  //
  // Input:
  //  eventCat  ---> Base category for the event (i.e. for Mdc it is catMdc)
  //  name      ---> Name of the new partial event (used to build Root trees)
  //  title     ---> Title of the new partial event
   HPartialEvent *event;
   event=new HPartialEvent(name,title,eventCat);
   if (event) fPartialEvs->AddAtAndExpand(event,eventCat>>kBitCategorySize);
   return event;
}

void HRecEvent::addPartialEvent(HPartialEvent* event) {
  // Add partial event by pointer.
   fPartialEvs->AddAtAndExpand(event,event->getBaseCat()>>kBitCategorySize);
}

void HRecEvent::merge(HRecEvent *targetEv) {
  //
  // Merge this event into target event, i.e. move all partial events
  // with all their categories to target event.
  // (needed by HRootSource in merge mode)
  
   HPartialEvent *partial=0;
   Int_t nPartialEvs=fPartialEvs->GetEntriesFast();
   for (Int_t i=0;i<nPartialEvs;i++) {
     partial=(HPartialEvent*)(fPartialEvs->At(i));
     if (partial) {
        if (targetEv->getPartialEvent(partial->getBaseCat())==0)
                                   targetEv->addPartialEvent(partial);
        fPartialEvs->RemoveAt(i);
     }
   }
   
}

void HRecEvent::Browse(TBrowser *b) {
  // Event browser.
   TIter next(fPartialEvs);
   HPartialEvent *ev;
   
   if (fTracks) b->Add(fTracks,"Tracks");
   while ( (ev=(HPartialEvent *)next())!=NULL) {
      b->Add(ev);
   }
}

void HRecEvent::setExpandedStreamer(Bool_t t) {
  // If the event is in "expanded streamer" mode each partial event and category
  // will go to a different directory.
  TIter iter(fPartialEvs);
  HPartialEvent *pev=0;

  SetBit(32,t);
  while ( (pev=(HPartialEvent *)iter()) != 0) {
    pev->setExpandedStreamer(t);
  }
}

void HRecEvent::Streamer(TBuffer &R__b)
{
  // Stream an object of class HRecEvent.
  TObjArray *partialEvNames;
  Int_t nPartialEvs;
  HPartialEvent *partial=0;
  TDirectory *dir,*top=0;
  Char_t buffer[255],b2[255];

  if (R__b.IsReading()) {
    Version_t R__v = R__b.ReadVersion(); if (R__v) { }
    HEvent::Streamer(R__b);
    R__b >> fRecLevel;
    R__b >> fNTracks;
    fHeader->Streamer(R__b);
    fTracks->Streamer(R__b);

    if (hasExpandedStreamer() || R__v == 1) {
      // Version 1 is allways expanded
      // At this point fBits is already read. 
      // And so hasExpandedStreamer makes sense.
      R__b >> partialEvNames;

      top=dir=gDirectory;
      for (Int_t i=0;i<partialEvNames->GetEntriesFast();i++) {
	strcpy(buffer,((TObjString *)partialEvNames->At(i))->GetString().Data());
	if (strcmp(buffer,"none")!=0) {
	  sprintf(b2,"dir%s",buffer);
	  top->cd(b2);
	  dir=gDirectory;
	  partial=(HPartialEvent *)dir->Get(buffer);
	  fPartialEvs->AddAtAndExpand(partial,i);
	}
      }
    } else {
      if (fPartialEvs) fPartialEvs->Delete();
      R__b >> fPartialEvs;
    } 
  } else {
//    printf(" hrecevent::streamer called \n");
    R__b.WriteVersion(HRecEvent::IsA());
    HEvent::Streamer(R__b);
    R__b << fRecLevel;
    R__b << fNTracks;
    fHeader->Streamer(R__b);
    fTracks->Streamer(R__b);
    if (hasExpandedStreamer()) {
      nPartialEvs=fPartialEvs->GetEntriesFast();
      partialEvNames=new TObjArray(nPartialEvs);
      for (Int_t i=0;i<nPartialEvs;i++) {
	partial=(HPartialEvent *)fPartialEvs->At(i);
	if (partial) {
	  partialEvNames->AddAt(new TObjString(partial->GetName()),i);
	} else {
	  partialEvNames->AddAt(new TObjString("none"),i);
	}
      }
      R__b << partialEvNames;
      top=gDirectory;
      for (Int_t i=0;i<nPartialEvs;i++) {
	partial=(HPartialEvent *)fPartialEvs->At(i);
	if (partial) {
	  sprintf(buffer,"dir%s",partial->GetName());
	  if (!(gDirectory->GetKey(buffer))) {
	    dir=gDirectory->mkdir(buffer);
	    dir->cd();
	  }
	  partial->Write();
	  top->cd();
	}
      }
    } else {
      R__b << fPartialEvs;
    }

  }
}  


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