T2LifeROOTFileManager.cc
Go to the documentation of this file.
1 
8 #include <utl/Reader.h>
9 #include <utl/TimeStamp.h>
10 #include <utl/ErrorLogger.h>
11 #include <utl/SaveCurrentTDirectory.h>
12 
13 #include <fwk/CentralConfig.h>
14 
15 #include <det/Detector.h>
16 #include <det/ValidityStamp.h>
17 
18 #include <sdet/T2LifeROOTFileManager.h>
19 
20 #include <boost/lexical_cast.hpp>
21 
22 #include <TFile.h>
23 #include <TTree.h>
24 #include <TBits.h>
25 
26 #include <sstream>
27 #include <string>
28 #include <set>
29 #include <map>
30 #include <utility>
31 #include <algorithm>
32 #include <iomanip>
33 
34 using namespace fwk;
35 using namespace det;
36 using namespace sdet;
37 using namespace utl;
38 using namespace std;
39 
40 
41 REGISTER_S_MANAGER("T2LifeROOTFileManager", T2LifeROOTFileManager);
42 
43 
44 const double T2LifeROOTFileManager::fgSearchMapBinning = 10000; //s
45 
46 
47 void
48 T2LifeROOTFileManager::Init(const std::string& configLink)
49 {
50  VManager::Init(configLink);
51 
52  // ************************ DEBUG **************************************************
53  fVerbosity = 0;
54  if (fBranch.GetChild("verbosity")) {
55  fBranch.GetChild("verbosity").GetData(fVerbosity);
56  cout << " T2LifeROOTFileManager verbosity level = " << fVerbosity << endl;
57  }
58  // **********************************************************************************
59 
60  fBranch.GetChild("file").GetData(fFileName);
61 
62  const SaveCurrentTDirectory save;
63  TFile* file = TFile::Open(fFileName.c_str(), "READ");
64  if (!file || file->IsZombie()) {
65  ostringstream err;
66  err << " T2LifeROOTFileManager input file=" << fFileName
67  << " is corrupt! Up times will not be available!" << endl;
68  INFO_TERSE(err);
69  file = 0;
70  fDataTree = 0;
71  } else
72  fDataTree = (TTree*)(file->Get("T2Life"));
73 
74  if (fDataTree) {
75 
76  // list of properties for telescope
77  fAvailableComponents.insert("station_in_acquisition");
78 
79  // create index map for fast lookup
80  TBranch* const startB = fDataTree->GetBranch("gpsStart");
81  TBranch* const endB = fDataTree->GetBranch("gpsStop");
82  UInt_t gpsStart = 0;
83  UInt_t gpsEnd = 0;
84  startB->SetAddress(&gpsStart);
85  endB->SetAddress(&gpsEnd);
86 
87  INFO_TERSE("initializing T2LifeROOTFileManager index map ... ");
88  int nT2bullshit = 0;
89  for (unsigned int i = 0; i < (unsigned int)startB->GetEntries(); ++i) {
90 
91  startB->GetEntry(i);
92  endB->GetEntry(i);
93 
94  unsigned int index_start = (unsigned int)(gpsStart / fgSearchMapBinning);
95  unsigned int index_end = (unsigned int)(gpsEnd / fgSearchMapBinning);
96 
97  if (gpsStart < 756950400 || gpsEnd < 756950400) {
98 
99  // ************************ DEBUG **************************************************
100  if (fVerbosity >= 3) {
101  ostringstream warn;
102  warn << " removing broken T2 period: start=" << gpsStart
103  << " end=" << gpsEnd;
104  WARNING(warn.str());
105  }
106  // **********************************************************************************
107 
108  ++nT2bullshit;
109  continue;
110  }
111 
112  if (gpsStart > gpsEnd) {
113  // ************************ DEBUG **************************************************
114  if (fVerbosity >= 3) {
115  cout << " strange T2 period: start=" << gpsStart
116  << " end=" << gpsEnd
117  << endl;
118  }
119  // **********************************************************************************
120  }
121  vector<unsigned int> iCurr;
122  if (index_end == index_start) {
123  iCurr.push_back(index_end);
124  } else {
125  for (unsigned int iIndex = min(index_start, index_end);
126  iIndex <= max(index_start, index_end); ++iIndex) {
127  iCurr.push_back(iIndex);
128  }
129  }
130 
131  for (unsigned int j = 0; j < iCurr.size(); ++j) {
132 
133  int iIndex = iCurr[j];
134 
135  if (fIndex.find(iIndex) == fIndex.end()) {
136 
137  fIndex[iIndex] = make_pair(i, i);
138 
139  } else { // change bounds
140 
141  if (i < fIndex[iIndex].first)
142  fIndex[iIndex].first = i;
143 
144  if (i > fIndex[iIndex].second)
145  fIndex[iIndex].second = i;
146 
147  }
148  }
149  }
150 
151  // ************************ DEBUG **************************************************
152  if (fVerbosity > 3) {
153  cout << " row of T2 files with wrong gps times: " << nT2bullshit << endl;
154  for (IndexIterator iIndex = fIndex.begin();
155  iIndex != fIndex.end(); ++iIndex) {
156  cout << " index-maxp-entry: bin=" << iIndex->first
157  << " start=" << iIndex->second.first
158  << " end=" << iIndex->second.second << endl;
159  }
160  }
161  // **********************************************************************************
162 
163  }
164 
165  ostringstream info;
166  info << "using T2 ROOT file: " << fFileName;
167  INFO(info);
168 }
169 
170 
172 T2LifeROOTFileManager::GetOkFlag(int& returnData,
173  const std::string& componentProperty,
174  const std::string& componentName,
175  const IndexMap& componentIndex)
176  const
177 {
178  if (!fDataTree ||
179  fAvailableComponents.find(componentProperty)==fAvailableComponents.end()) {
180  return VManager::eNotFound;
181  }
182 
183  if (!ReadData()) {
184  ostringstream warn;
185  warn << "*****************************************************\n"
186  " property=" << componentProperty << " component=" << componentName
187  << " index=(";
188  for (IndexMap::const_iterator imap = componentIndex.begin();
189  imap != componentIndex.end(); ++imap) {
190  warn << "\"" << imap->first << "\"=\"" << imap->second << "\" ";
191  }
192  warn << ") NOT found in data source for specified time!\n"
193  "Returning zero \n"
194  "*****************************************************";
195  WARNING(warn);
196  returnData = 0;
197  return VManager::eFound;
198  }
199 
200  // ************************ DEBUG **************************************************
201  if (fVerbosity > 2) {
202  cout << "called GetOkFlag for componentName=" << componentName
203  << " componentProperty=" << componentProperty << " and index map:\n";
204  for (IndexMap::const_iterator imap = componentIndex.begin();
205  imap != componentIndex.end(); ++imap) {
206  cout << " index=\"" << imap->first << "\" -> \"" << imap->second << "\"\n";
207  }
208  }
209  // **********************************************************************************
210 
211  if (componentProperty == "station_in_acquisition") {
212 
213  if (componentName == "is_in_acquisition") {
214 
215  int dataId = GetStationIndex(componentIndex);
216  //dataId++;
217 
218  // ************************ DEBUG **************************************************
219  if (fVerbosity > 2)
220  cout << " no stations=" << fStationList.size() << " id=" << dataId << endl;
221  // **********************************************************************************
222 
223  returnData = 0;
224  if( (unsigned int)dataId < fStationList.size())
225  returnData = fStationList[dataId];
226 
227  // ************************ DEBUG **************************************************
228  if (fVerbosity > 2)
229  cout << " IsInDAQ=" << returnData << endl;
230  // **********************************************************************************
231 
232  return VManager::eFound;
233 
234  }
235  }
236 
237  return VManager::eNotFound;
238 }
239 
240 
241 bool
242 T2LifeROOTFileManager::ReadData()
243  const
244 {
245  // ************************ DEBUG **************************************************
246  if (fVerbosity > 2)
247  INFO( " enter ReadData ");
248  // **********************************************************************************
249 
250  if (fDataValidity.IsValid()) {
251  // ************************ DEBUG **************************************************
252  if (fVerbosity > 2)
253  INFO (" T2 data still valid ");
254  // **********************************************************************************
255  return true;
256  }
257 
258  unsigned int now = det::Detector::GetInstance().GetTime().GetGPSSecond();
259 
260  // really get the data from the TTree
261  fStationList.clear();
262 
263  // using the index map to get first bounds ...
264  unsigned int iBin = (unsigned int)(now / fgSearchMapBinning);
265 
266  ConstIndexIterator iIndex = fIndex.find(iBin);
267  if (iIndex == fIndex.end()) {
268  ostringstream err;
269  err << "Index not found! now/gps=" << now << ", iBin=" << iBin;
270  ERROR(err);
271  return false;
272  }
273 
274  unsigned int lower = iIndex->second.first;
275  unsigned int upper = iIndex->second.second;
276 
277  // ************************ DEBUG **************************************************
278  if (fVerbosity > 2) {
279  ostringstream info;
280  info << " index: now=" << now
281  << " index=" << iBin
282  << " lower=" << lower
283  << " upper=" << upper;
284  INFO(info);
285  }
286  // **********************************************************************************
287 
288  // start looking into TTree
289 
290  TBranch* const startB = fDataTree->GetBranch("gpsStart");
291  TBranch* const endB = fDataTree->GetBranch("gpsStop");
292  UInt_t gpsStart = 0;
293  UInt_t gpsEnd = 0;
294  startB->SetAddress(&gpsStart);
295  endB->SetAddress(&gpsEnd);
296 
297  int nSearch = 0;
298  int index = 0 ;
299  bool foundit = false;
300  while (!foundit && upper >= lower) {
301 
302  index = lower + (upper - lower)/2; // binary search !
303 
304  startB->GetEntry(index);
305  endB->GetEntry(index);
306 
307  if (now >= gpsStart && now < gpsEnd) { // found entry ?
308 
309  foundit = true;
310 
311  } else { // not yet found right entry -> new seach boundaries
312 
313  if (now < gpsStart)
314  upper = index - 1;
315  else
316  lower = index + 1;
317  }
318 
319  ++nSearch;
320 
321  // ************************ DEBUG **************************************************
322  if (fVerbosity > 3) {
323  cout << " search loop: nSearch=" << nSearch
324  << " gpsStart=" << gpsStart
325  << " gpsEnd=" << gpsEnd
326  << " foundit=" << foundit
327  << " index=" << index
328  << endl;
329  }
330  // **********************************************************************************
331 
332  }
333 
334  // ************************ DEBUG **************************************************
335  if (fVerbosity > 2) {
336  cout << " END search loop: "
337  << " foundit=" << foundit
338  << " index=" << index
339  << endl;
340  }
341  // **********************************************************************************
342 
343  if (foundit) {
344 
345  UInt_t time_start = 0;
346  UInt_t time_stop = 0;
347  TBits* stationBits = new TBits();
348 
349  fDataTree->SetBranchAddress("gpsStart", &time_start);
350  fDataTree->SetBranchAddress("gpsStop", &time_stop);
351  fDataTree->SetBranchAddress("stationArray", &stationBits);
352 
353  fDataTree->GetEntry(index);
354 
355  const unsigned int nStations = stationBits->GetNbits();
356 
357  // station-index is the 'stationId'
358  fStationList.resize(nStations + 1);
359  for (unsigned int iStation = 1; iStation <= nStations; ++iStation) {
360 
361  if (stationBits->TestBitNumber(iStation))
362  fStationList[iStation] = 1;
363  else
364  fStationList[iStation] = 0;
365 
366  }
367 
368  fDataValidity.SetValidityInterval(TimeStamp(time_start), TimeStamp(time_stop));
369 
370  delete stationBits;
371  }
372 
373  // ************************ DEBUG **************************************************
374  if (fVerbosity > 2) {
375  cout << " ReadData result is valid from start: " << fDataValidity.GetStartTime().GetGPSSecond()
376  << " end: " << fDataValidity.GetStopTime().GetGPSSecond()
377  << endl;
378  if (fVerbosity > 4) {
379  cout << " Stations UP:\n";
380  int n = 0;
381  for (unsigned int i = 0; i<fStationList.size(); ++i) {
382  if (fStationList[i]) {
383  ++n;
384  if (n > 1)
385  cout << ", ";
386  cout << setw(4) << i;
387  }
388  if (n == 25) {
389  n = 0;
390  cout << "\n";
391  }
392  }
393  cout << "\n";
394  }
395  }
396  // **********************************************************************************
397 
398  return foundit;
399  /*
400  ostringstream query;
401  query << "SELECT * FROM upTimes WHERE "
402  << "uptime_start_time <= \"" << now << "\" AND "
403  << "uptime_end_time > \"" << now << "\"";
404 
405  TSQLResult *result = fDataTree->Query(query.str().c_str());
406  if (result->GetRowCount()) {
407  INFO_TERSE("T2LifeROOTFileManager cannot identify row");
408  }
409 
410  TSQLRow *row = result->Next();
411  row
412  */
413 }
414 
415 
416 int
417 T2LifeROOTFileManager::GetStationIndex(const IndexMap& componentIndex)
418  const
419 {
420  IndexMap::const_iterator iStationId = componentIndex.find("stationId");
421  if (iStationId == componentIndex.end()) {
422  INFO_TERSE("Called T2LifeROOTFileManager with wrong number of arguments!");
423  return -1;
424  }
425 
426  return boost::lexical_cast<int>(iStationId->second);
427 }
#define INFO_TERSE(message)
Macro for logging informational messages.
Definition: ErrorLogger.h:175
Class to manage Sd life time ROOT files.
#define INFO(message)
Macro for logging informational messages.
Definition: ErrorLogger.h:161
void Init()
Initialise the registry.
std::map< unsigned int, IndexEntry >::const_iterator ConstIndexIterator
A TimeStamp holds GPS second and nanosecond for some event.
Definition: TimeStamp.h:110
#define max(a, b)
#define REGISTER_S_MANAGER(_name_, _Type_)
const double second
Definition: GalacticUnits.h:32
#define WARNING(message)
Macro for logging warning messages.
Definition: ErrorLogger.h:163
std::map< unsigned int, IndexEntry >::iterator IndexIterator
const string file
std::map< std::string, std::string > IndexMap
Definition: VManager.h:133
#define ERROR(message)
Macro for logging error messages.
Definition: ErrorLogger.h:165
Status
Specifies success or (eventually) various possible failure modes.
Definition: VManager.h:127

, generated on Tue Sep 26 2023.