FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairTSBufferFunctional.cxx
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
9 
10 #include "FairLink.h" // for FairLink
11 #include "FairRootManager.h" // for FairRootManager
12 #include "FairTimeStamp.h" // for FairTimeStamp
13 
14 #include <TBranch.h> // for TBranch
15 #include <TClass.h> // for TClass
16 #include <TClonesArray.h> // for TClonesArray
17 #include <TTree.h> // for TTree
18 
20 
22  TTree* sourceTree,
23  BinaryFunctor* stopFunction,
24  BinaryFunctor* startFunction)
25  : TObject()
26  , fOutputArray(nullptr)
27  , fBufferArray(nullptr)
28  , fInputArray(nullptr)
29  , fStartFunction(startFunction)
30  , fStopFunction(stopFunction)
31  , fBranch(nullptr)
32  , fBranchIndex(-1)
33  , fTerminate(kFALSE)
34  , fVerbose(0)
35 {
36  fBranch = sourceTree->GetBranch(branchName.Data());
37  if (fBranch == 0) {
38  std::cout << "-E- FairTSBufferFunctional::FairTSBufferFunctional Branch " << branchName << " does not exist!"
39  << std::endl;
40  }
42  fInputArray = static_cast<TClonesArray*>(ioman->GetObject(branchName.Data()));
43  fBufferArray = new TClonesArray(fInputArray->GetClass()->GetName());
44  fOutputArray = new TClonesArray(fInputArray->GetClass()->GetName());
45 }
46 
47 TClonesArray* FairTSBufferFunctional::GetData(Double_t stopParameter)
48 {
49  Double_t actualTime = 0.;
50  int posBuffer = 0;
51 
52  if (fStopFunction == 0) { // no function is given ==> read in data in traditional way event by event
53  ReadInNextEntry();
54  fOutputArray->AbsorbObjects(static_cast<TClonesArray*>(fInputArray), 0, fInputArray->GetEntriesFast() - 1);
55  return fOutputArray;
56  }
57  if (fVerbose > 1) {
58  std::cout << "-I- FairTSBufferFunctional::GetData for stopParameter: " << stopParameter << std::endl;
59  }
60 
61  // if the BufferArray is empty fill it
62  if (fBufferArray->GetEntriesFast() == 0) {
63  if (fVerbose > 1) {
64  std::cout << "-I- FairTSBufferFunctional::GetData fBufferArray is empty: Read in Data" << std::endl;
65  }
66  ReadInNextFilledEntry();
67  AbsorbDataBufferArray();
68  }
69 
70  // if the BufferArray is still empty you have reached the end of your data set
71  if (fVerbose > 1) {
72  std::cout << "-I- FairTSBufferFunctional::GetData fBufferArray->GetEntriesFast(): "
73  << fBufferArray->GetEntriesFast() << std::endl;
74  }
75  FairTimeStamp* dataPoint = static_cast<FairTimeStamp*>(fBufferArray->Last());
76  if (dataPoint == 0) {
77  if (fVerbose > 0) {
78  std::cout << "-I- FairTSBufferFunctional::GetData dataPoint is empty ==> All Data read in" << std::endl;
79  }
80  return fOutputArray;
81  }
82 
83  dataPoint = static_cast<FairTimeStamp*>(fBufferArray->First());
84 
85  while (!(*fStopFunction)(dataPoint, stopParameter)) { // check if you have reached end of requested data
86  posBuffer++;
87  // if you have reached the end of the BufferArray fill it with new data from tree
88  if (posBuffer == fBufferArray->GetEntriesFast()) {
89  ReadInNextFilledEntry();
90  AbsorbDataBufferArray();
91  }
92  // if you are still at the end of the BufferArray than break (no new data in tree)
93  if (posBuffer == fBufferArray->GetEntriesFast()) {
94  break;
95  }
96  dataPoint = static_cast<FairTimeStamp*>(fBufferArray->At(posBuffer));
97  if (fVerbose > 1) {
98  std::cout << posBuffer << " TimeStampData: " << dataPoint->GetTimeStamp() << std::endl;
99  }
100  }
101 
102  if (fVerbose > 1) {
103  std::cout << "-I- FairTSBufferFunctional::GetData Index for Absorb: " << posBuffer
104  << " BufferArray size: " << fBufferArray->GetEntriesFast() << std::endl;
105  }
106  if (posBuffer < fBufferArray->GetEntriesFast() && posBuffer > 0) {
107  if (fVerbose > 1) {
108  std::cout << "-I- FairTSBufferFunctional::GetData absorb BufferArray up to posBuffer " << posBuffer
109  << " into fOutputArray" << std::endl;
110  }
111  fOutputArray->AbsorbObjects(fBufferArray, 0, posBuffer - 1);
112  posBuffer = 0;
113  return fOutputArray;
114  }
115  if (fVerbose > 1) {
116  std::cout << "Index: " << posBuffer << " BranchIndex: " << fBranchIndex << " NBranch " << fBranch->GetEntries()
117  << std::endl;
118  }
119 
120  if (posBuffer >= fBufferArray->GetEntriesFast() && posBuffer != 0 && fBranchIndex + 1 >= fBranch->GetEntries()) {
121  if (fVerbose > 1) {
122  std::cout << "-I- FairTSBufferFunctional::GetData end of data reached. Send the rest to the OutputArray!"
123  << std::endl;
124  }
125  fOutputArray->AbsorbObjects(fBufferArray, 0, fBufferArray->GetEntries() - 1);
126  }
127 
128  if (fVerbose > 1) {
129  std::cout << "-I- FairTSBufferFunctional::GetData: Read in up to entry: " << fBranchIndex << " with actualTime "
130  << actualTime << " and requested stopParameter " << stopParameter << std::endl;
131  }
132 
133  return fOutputArray;
134 }
135 
136 TClonesArray* FairTSBufferFunctional::GetData(Double_t startParameter, Double_t stopParameter)
137 {
138  if (fStartFunction != 0) {
139  fBufferArray->Delete();
140  fOutputArray->Delete();
141  Int_t startIndex = FindStartIndex(startParameter);
142  // std::cout << "StartParameter : " << startParameter << " StartIndex: " << startIndex << "/" <<
143  // GetBranchIndex() << " size BufferArray " << fBufferArray->GetEntries() <<std::endl;
144  if (startIndex > -1) {
145  ReadInEntry(fBranchIndex);
146  fBufferArray->AbsorbObjects(fInputArray, startIndex, fInputArray->GetEntries() - 1);
147  }
148  }
149  return GetData(stopParameter);
150 }
151 
152 Int_t FairTSBufferFunctional::FindStartIndex(Double_t startParameter)
153 {
154  FairTimeStamp* dataPoint;
155  Int_t tempIndex = fBranchIndex;
156  // Bool_t runBackwards = kTRUE;
157  Int_t previousIndex = -1;
158  Int_t previousBranchIndex = -1;
159 
160  ReadInEntry(tempIndex); //< Get Data out of Tree
161 
162  while (fInputArray->GetEntries() == 0 && tempIndex > 0) { //< If the entry of the tree was empty read in previous
163  // entries until you find one which is filled
164  tempIndex--;
165  ReadInEntry(tempIndex);
166  }
167 
168  if (fInputArray->GetEntries()
169  == 0) { // If the previous entries in the tree are also empty run in the forward direction
170  ReadInNextFilledEntry();
171  // runBackwards = kFALSE;
172  }
173 
174  if (fInputArray->GetEntries() == 0) { // If there is still no data the branch is empty!
175  std::cout << "-I- FairTSBufferFunctional::FindStartIndex: All entries are empty!" << std::endl;
176  return -1;
177  }
178  fBranchIndex = tempIndex;
179  // Now we have data or FindStartIndex already returned -1
180 
181  dataPoint = static_cast<FairTimeStamp*>(fInputArray->Last());
182  // std::cout << "DataPoint: " << *dataPoint << std::endl;
183  while (!(*fStartFunction)(dataPoint, startParameter)) {
184  // std::cout << "DataPoint Search Entry: " << fBranchIndex << ": " << *dataPoint << std::endl;
185  ReadInNextFilledEntry();
186  if (fInputArray->GetEntries() != 0) {
187  dataPoint = static_cast<FairTimeStamp*>(fInputArray->Last());
188  } else {
189  return -1;
190  }
191  }
192 
193  // Now you have data where the last element in the array does not fit to your request
194  Int_t startPos = fInputArray->GetEntries() - 1;
195  while ((*fStartFunction)(dataPoint, startParameter)) {
196  // std::cout << "DataPoint Search in Entry: " << fBranchIndex << ": " << *dataPoint << std::endl;
197  previousIndex = startPos;
198  previousBranchIndex = fBranchIndex;
199  startPos--;
200  if (startPos == -1) {
201  fBranchIndex = ReadInPreviousFilledEntry(fBranchIndex);
202  startPos = fInputArray->GetEntries() - 1;
203  if (startPos < 0) {
204  if (fBranchIndex == 0) {
205  return 0;
206  }
207  return -1;
208  }
209  }
210  dataPoint = static_cast<FairTimeStamp*>(fInputArray->At(startPos));
211  }
212  fBranchIndex = previousBranchIndex;
213  return previousIndex;
214 }
215 
216 void FairTSBufferFunctional::ReadInNextFilledEntry()
217 {
218  fInputArray->Delete();
219 
220  if (fVerbose > 1) {
221  std::cout << "-I- FairTSBufferFunctional::ReadInNextFilledEntry: Entries in InputArray "
222  << fInputArray->GetEntriesFast() << " Branch Entries: " << fBranch->GetEntries() << std::endl;
223  }
224  while (fInputArray->GetEntriesFast() == 0 && fBranchIndex + 1 < fBranch->GetEntries()) {
225  fBranchIndex++;
226  ReadInEntry(fBranchIndex);
227  }
228 }
229 
230 Int_t FairTSBufferFunctional::ReadInPreviousFilledEntry(Int_t startEntry)
231 {
232  Int_t tempIndex = startEntry;
233 
234  fInputArray->Delete();
235  while (fInputArray->GetEntriesFast() == 0 && tempIndex > 0) {
236  tempIndex--;
237  ReadInEntry(tempIndex);
238  }
239  return tempIndex;
240 }
241 
242 void FairTSBufferFunctional::AbsorbDataBufferArray()
243 {
244  if (fInputArray->GetEntriesFast() > 0) {
245  if (fVerbose > 1) {
246  std::cout << "-I- FairTSBufferFunctional::ReadInNextFilledEntry: Absorb InputArray into Buffer"
247  << std::endl;
248  }
249  fBufferArray->AbsorbObjects(fInputArray, 0, fInputArray->GetEntries() - 1);
250  }
251 }
252 
253 void FairTSBufferFunctional::ReadInNextEntry()
254 {
255  if (fBranchIndex + 1 < fBranch->GetEntries()) {
256  fBranchIndex++;
257  ReadInEntry(fBranchIndex);
258  }
259 }
260 
261 void FairTSBufferFunctional::ReadInEntry(Int_t number)
262 {
263  fInputArray->Delete();
264  if (number < fBranch->GetEntries()) {
265  fBranch->GetEntry(number);
266  for (int i = 0; i < fInputArray->GetEntriesFast(); i++) {
267  (static_cast<FairTimeStamp*>(fInputArray->At(i)))
268  ->SetEntryNr(FairLink(0, number, FairRootManager::Instance()->GetBranchId(fBranch->GetName()), i, 1));
269  }
270  if (fVerbose > 1)
271  std::cout << "-I- FairTSBufferFunctional::ReadInEntry BranchIndex: " << number
272  << " Entries: " << fInputArray->GetEntriesFast() << std::endl;
273  }
274 }
275 
277 {
278  if (fTerminate == kTRUE) {
279  return kTRUE;
280  }
281  if (fBranchIndex + 1 >= fBranch->GetEntries()) {
282  if (fBufferArray->GetEntriesFast() == 0) {
283  if (fOutputArray->GetEntriesFast() == 0) {
284  return kTRUE;
285  } else {
286  return kFALSE;
287  }
288  } else {
289  return kFALSE;
290  }
291  } else {
292  return kFALSE;
293  }
294 }
static FairRootManager * Instance()
ClassImp(FairEventBuilder)
Int_t FindStartIndex(Double_t startParameter)
TObject * GetObject(const char *BrName)
Double_t GetTimeStamp() const
Definition: FairTimeStamp.h:44
TClonesArray * GetData(Double_t stopParameter)
Base class for all functors which are used in the FairTSBufferFunctional.
FairTSBufferFunctional(TString branchName, TTree *sourceTree, BinaryFunctor *stopFunction, BinaryFunctor *startFunction=0)
A class to access time ordered data in a root branch.