StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
MapTableTest.C
1  {
2  // Copyright(c) 2001 [BNL] Brookhaven National Laboratory, Valeri Fine (fine@bnl.gov). All right reserved",
3  //
4  // TGenericTable and TTableMap classes test macro
5  // ----------------------------------------------
6  // This example creates a simple event containing one "track" table and one "hit" table
7  // The "track" table is in possession of the "column" with variable length array referencing
8  // the appropriated "hit table rows"
9  //
10  // Macro:
11  // 1. Create an event
12  // 2. Write it out
13  // 3. Close ROOT file
14  // 4. Re-open ROOT file
15  // 5. Read the event back
16  // 6. Plot a histogram for all hits of the first track (for instance)
17  //
18  // --------------------------------------------------
19  // One needs NO extra share library to be loaded !!!!
20  // --------------------------------------------------
21 
22  //---------------------------------------------------
23  // 1. Load "libRootKernel" share library
24  //---------------------------------------------------
25  gSystem->Load("libTable");
26 
27  //---------------------------------------------------
28  // 2. Define the C-structure to describe "track" and "hit"
29  //---------------------------------------------------
30  struct hit {
31  float energy; // energy
32  int detectorId; // geometry id
33  };
34 
35  // Pay attention each track element has a pointer to the container
36  // of hit-references (watch "hitList" data-member)
37 
38  struct track {
39  float curvature; // curvature
40  Ptr_t hitList; // the list of this track hits from hit table
41  };
42 
43  //---------------------------------------------------
44  // 3. Create and fill the hit table object first
45  //---------------------------------------------------
46 
47  TGenericTable *allHits = new TGenericTable("hit","hits",1000);
48  hit a;
49  int i = 0;
50  for (i=0; i<100; i++) {
51  a.energy = sin(i*0.1);
52  a.detectorId++;
53  allHits->AddAt(&a);
54  }
55  allHits->Print();
56  //---------------------------------------------------
57  // 4. Create and fill the track table object
58  //---------------------------------------------------
59  TGenericTable *allTracks = new TGenericTable("track","tracks",21);
60  allTracks->Print();
61  track t;
62  int i = 0;
63  for (i=0; i<20; i++) {
64  t.curvature = gRandom->Rndm();
65  t.hitList = new TTableMap(allHits);
66  // TGenericTable::iterator h = allHits->begin();
67  // TGenericTable::iterator hLast = allHits->end();
68  //---------------------------------------------------
69  // 5. Simulate track->hit association
70  //---------------------------------------------------
71  const char *h = allHits->GetArray();
72  Int_t len = allHits->GetNRows();
73  UInt_t indx;
74  for (indx = 0;indx < len; indx++ ) {
75  if (gRandom->Rndm() < 0.2 )
76  t.hitList.Push_back(indx);
77  // In compiled code one can use the "regular" vector::push_back(indx)
78  // as follows:
79  // t.hitList.push_back(indx);
80  }
81  allTracks->AddAt(&t);
82  }
83  allTracks->Print(0,5);
84  //---------------------------------------------------
85  // 6. Create the full event
86  //---------------------------------------------------
87  TDataSet *event = new TDataSet("event");
88  event->Add(allHits);
89  event->Add(allTracks);
90  event->ls(9);
91  //---------------------------------------------------
92  // 6. Write the full event out
93  //---------------------------------------------------
94  TFile ff("generic.root","RECREATE");
95  event->Write();
96  ff.Write();
97  ff.Close();
98  printf(" One event has been written out\n");
99  //-----------------------------------------------------
100  // Read the event back now
101  //-----------------------------------------------------
102  printf("\n ----- \n Now we will try to read it back\n");
103  printf(" ----- \n");
104 
105  // delete event;
106  TFile newFile("generic.root");
107  event = (TDataSet *)newFile.Get("event");
108  event->ls(4);
109  TGenericTable *rdTracks = (TGenericTable *)event->FindByName("tracks");
110  rdTracks->Print(0,5);
111  //-----------------------------------------------------
112  // Getting the list of the hits for some 1-st track
113  //-----------------------------------------------------
114  TH1F *ehist = new TH1F("ehist","Energy deposit",30,-1,1);
115  track &firstTrack = *(track *)rdTracks->GetTable();
116  TTable::iterator fHit = firstTrack.hitList.Begin();
117  TTable::iterator lHit = firstTrack.hitList.End();
118  for (;fHit != lHit; fHit++) {
119  hit &nextHit = *(hit *)fHit.rowPtr();
120  // in compiled code one can apply (operator *())
121  // hit &nextHit = *(hit *)*fHit;
122  ehist->Fill(nextHit.energy);
123  }
124  ehist->Draw();
125  printf("\n ------------ >> finish !!! << ---------------- \n");
126  }
virtual Long_t GetNRows() const
Returns the number of the used rows for the wrapped table.
Definition: TTable.cxx:1388
virtual Int_t AddAt(const void *c)
Definition: TTable.cxx:1122
virtual Char_t * Print(Char_t *buf, Int_t n) const
Create IDL table defintion (to be used for XDF I/O)
Definition: TTable.cxx:1548