]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/AliEveTPCData.cxx
Remove trailing whitespace.
[u/mrichter/AliRoot.git] / EVE / Alieve / AliEveTPCData.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9
10 #include "AliEveTPCData.h"
11
12 #include <Alieve/AliEveTPCSectorData.h>
13
14 #include <AliSimDigits.h>
15 #include <AliTPCParam.h>
16 #include <AliTPCRawStream.h>
17 #include <TTree.h>
18
19
20 //______________________________________________________________________
21 // AliEveTPCData
22 //
23 // A central manager for TPC data of an event.  Can read digits (from
24 // a tree: LoadDigits()) and raw-data (via AliRawReader: LoadRaw()).
25 //
26 // The sector data is stored in 36 AliEveTPCSectorData objects.
27 // Sectors 0 - 17: +z side, 18 - 35: -z side.
28 // No separation of inner/outer segments, use row numbers for addressing.
29 //
30 // Threshold application and pedestal subtraction can be performed at
31 // load time: use SetLoadThreshold(thresh) and SetLoadPedestal(ped).
32 //
33 // For raw-data (loaded using LoadRaw) pedestals can be calculated
34 // automatically per pad. Use SetAutoPedestal(kTRUE) to activate it.
35 // You might still want to set load threshold (default iz zero).
36 //
37
38 ClassImp(AliEveTPCData)
39
40 AliEveTPCData::AliEveTPCData() :
41   fSectors(36), fSectorBlockSize(65536),
42   fLoadThreshold(0), fLoadPedestal(0), fAutoPedestal(kFALSE)
43 {
44   AliEveTPCSectorData::InitStatics();
45 }
46
47 AliEveTPCData::~AliEveTPCData()
48 {
49   DeleteAllSectors();
50 }
51
52 /**************************************************************************/
53
54 void AliEveTPCData::CreateSector(Int_t sector)
55 {
56   if(fSectors[sector] == 0)
57     fSectors[sector] = new AliEveTPCSectorData(sector, fSectorBlockSize);
58 }
59
60 void AliEveTPCData::CreateAllSectors()
61 {
62   for(Int_t s=0; s<36; ++s)
63     CreateSector(s);
64 }
65
66 void AliEveTPCData::DropAllSectors()
67 {
68   for(Int_t s=0; s<36; ++s) {
69     if(fSectors[s] != 0)
70       fSectors[s]->DropData();
71   }
72 }
73
74 void AliEveTPCData::DeleteAllSectors()
75 {
76   for(Int_t s=0; s<36; ++s) {
77     delete fSectors[s];
78     fSectors[s] = 0;
79   }
80 }
81
82 /**************************************************************************/
83
84 AliEveTPCSectorData* AliEveTPCData::GetSectorData(Int_t sector, Bool_t spawnSectors)
85 {
86   if(sector < 0 || sector > 35) return 0;
87   if(fSectors[sector] == 0 && spawnSectors)
88     CreateSector(sector);
89   return fSectors[sector];
90 }
91
92 /**************************************************************************/
93
94 void AliEveTPCData::LoadDigits(TTree* tree, Bool_t spawnSectors)
95 {
96   // Load data from TTree of AliSimDigits.
97   // If spawnSectors is false only sectors that have been created previously
98   // via CreateSector() are loaded.
99   // If spawnSectors is true sectors are created if data for them is encountered.
100
101   AliSimDigits digit, *digitPtr = &digit;
102   tree->GetBranch("Segment")->SetAddress(&digitPtr);
103
104   Int_t sector, row, pad, curPad;
105   Short_t time, signal;
106   Bool_t  inFill = kFALSE;
107   AliEveTPCSectorData* secData = 0;
108
109   Int_t numEnt = (Int_t) tree->GetEntries();
110   for (Int_t ent=0; ent<numEnt; ent++) {
111     tree->GetEntry(ent);
112     AliEveTPCSectorData::GetParam().AdjustSectorRow(digit.GetID(), sector, row);
113     if(sector >= 36) {
114       sector -= 36;
115       row    += AliEveTPCSectorData::GetInnSeg().GetNRows();
116     }
117     secData = GetSectorData(sector, spawnSectors);
118     if(secData == 0)
119       continue;
120
121     if(digit.First() == kFALSE)
122       continue;
123     curPad = -1;
124     do {
125       pad    = digit.CurrentColumn();
126       time   = digit.CurrentRow();
127       signal = digit.CurrentDigit();
128
129       if(pad != curPad) {
130         if(inFill)
131           secData->EndPad(fAutoPedestal, fLoadThreshold);
132         secData->BeginPad(row, pad, kFALSE);
133         curPad = pad;
134         inFill = kTRUE;
135       }
136       if(fAutoPedestal) {
137         secData->RegisterData(time, signal);
138       } else {
139         signal -= fLoadPedestal;
140         if(signal >= fLoadThreshold)
141           secData->RegisterData(time, signal);
142       }
143
144     } while (digit.Next());
145     if(inFill) {
146       secData->EndPad(fAutoPedestal, fLoadThreshold);
147       inFill = kFALSE;
148     }
149   }
150 }
151
152 void AliEveTPCData::LoadRaw(AliTPCRawStream& input, Bool_t spawnSectors, Bool_t warn)
153 {
154   // Load data from AliTPCRawStream.
155   // If spawnSectors is false only sectors that have been created previously
156   // via CreateSector() are loaded.
157   // If spawnSectors is true sectors are created if data for them is encountered.
158
159   static const TEveException eH("AliEveTPCData::LoadRaw ");
160
161   Int_t   sector = -1, row = -1, pad = -1, rowOffset = 0;
162   Short_t time,  signal;
163   Bool_t  inFill   = kFALSE;
164   Short_t lastTime = 9999;
165   Bool_t  lastTimeWarn = kFALSE;
166   AliEveTPCSectorData* secData = 0;
167
168   Short_t threshold = fLoadThreshold;
169
170   while (input.Next()) {
171     if (input.IsNewSector()) {
172       if(inFill) {
173         secData->EndPad(fAutoPedestal, threshold);
174         inFill = kFALSE;
175       }
176       sector = input.GetSector();
177       if(sector >= 36) {
178         sector -= 36;
179         rowOffset = AliEveTPCSectorData::GetInnSeg().GetNRows();
180       } else {
181         rowOffset = 0;
182       }
183       secData = GetSectorData(sector, spawnSectors);
184     }
185     if (secData == 0)
186       continue;
187
188     if (input.IsNewPad()) {
189       if(inFill) {
190         secData->EndPad(fAutoPedestal, threshold);
191         inFill = kFALSE;
192       }
193       row = input.GetRow() + rowOffset;
194       pad = input.GetPad();
195
196       if(pad >= AliEveTPCSectorData::GetNPadsInRow(row)) {
197         if(warn) {
198           Warning(eH.Data(), "pad out of range (row=%d, pad=%d, maxpad=%d).",
199                   row, pad, AliEveTPCSectorData::GetNPadsInRow(row));
200         }
201         continue;
202       }
203
204       AliEveTPCSectorData::PadRowHack* prh = secData->GetPadRowHack(row, pad);
205       if(prh != 0) {
206         threshold = prh->fThrExt + Short_t(prh->fThrFac*fLoadThreshold);
207       } else {
208         threshold = fLoadThreshold;
209       }
210
211       secData->BeginPad(row, pad, kTRUE);
212       inFill   = kTRUE;
213       lastTime = 1024;  lastTimeWarn = kFALSE;
214     }
215
216     time   = input.GetTime();
217     signal = input.GetSignal();
218     if(time >= lastTime) {
219       if(lastTimeWarn == kFALSE) {
220         if(warn)
221           Warning(eH.Data(), "time out of order (row=%d, pad=%d, time=%d, lastTime=%d).",
222                   row, pad, time, lastTime);
223         lastTimeWarn = kTRUE;
224       }
225       continue;
226     }
227     lastTime = time;
228     if(fAutoPedestal) {
229       secData->RegisterData(time, signal);
230     } else {
231       signal -= fLoadPedestal;
232       if(signal > threshold)
233         secData->RegisterData(time, signal);
234     }
235   }
236
237   if(inFill) {
238     secData->EndPad(fAutoPedestal, threshold);
239     inFill = kFALSE;
240   }
241 }