]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliOADBForward.cxx
Added possiblity to use fall-back queries against the
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliOADBForward.cxx
1 #include "AliOADBForward.h"
2 #include <TBrowser.h>
3 #include <TROOT.h>
4 #include <TKey.h>
5 #include <TList.h>
6 #include <TDatime.h>
7 #include <TTree.h>
8 #include <TFile.h>
9 #include <TError.h>
10 #include <TSystem.h>
11
12 #ifndef ALIROOT_SVN_REVISION
13 # define ALIROOT_SVN_REVISION 0
14 #endif
15
16
17 ClassImp(AliOADBForward)
18 #if 0
19 ; // Do not remove - for Emacs
20 #endif
21
22 //====================================================================
23 const char* 
24 AliOADBForward::Mode2String(ERunSelectMode mode)
25 {
26   switch (mode) { 
27   case kDefault:   return "default"; 
28   case kExact:     return "exact"; 
29   case kNewest:    return "newest";
30   case kNear:      return "near";
31   case kOlder:     return "older";
32   case kNewer:     return "newer";
33   }
34   return "?"; // Should never get here 
35 }
36 AliOADBForward::ERunSelectMode
37 AliOADBForward::String2Mode(const TString& str)
38 {
39   if      (str.EqualTo("default", TString::kIgnoreCase)) return kDefault;
40   else if (str.EqualTo("exact",   TString::kIgnoreCase)) return kExact;
41   else if (str.EqualTo("newest",  TString::kIgnoreCase)) return kNewest;
42   else if (str.EqualTo("near",    TString::kIgnoreCase)) return kNear;
43   else if (str.EqualTo("older",   TString::kIgnoreCase)) return kOlder;
44   else if (str.EqualTo("newer",   TString::kIgnoreCase)) return kNewer;
45   return kDefault;
46 }
47 AliOADBForward::ERunSelectMode
48 AliOADBForward::Int2Mode(Int_t mode)
49 {
50   switch (mode) { 
51   case kDefault:   return kDefault; 
52   case kExact:     return kExact; 
53   case kNewest:    return kNewest;
54   case kNear:      return kNear;
55   case kOlder:     return kOlder;
56   case kNewer:     return kNewer;
57   }
58   return kDefault; // Should never get here 
59 }
60
61
62 //====================================================================
63 AliOADBForward::Entry::Entry(ULong_t  runNo, 
64                              UShort_t sys, 
65                              UShort_t sNN, 
66                              Short_t  field, 
67                              Bool_t   mc, 
68                              Bool_t   sat,
69                              TObject* o)
70   : fRunNo(runNo), 
71     fSys(sys), 
72     fSNN(sNN), 
73     fField(field),
74     fMC(mc), 
75     fSatellite(sat),
76     fData(o),
77     fTimestamp(0), 
78     fAliROOTRevision(0),
79     fAuthor("unknown")
80 {
81   // 
82   // Constructor 
83   // 
84 }
85
86 //____________________________________________________________________
87 AliOADBForward::Entry::Entry(const Entry& o)
88   : TObject(o), 
89     fRunNo(o.fRunNo), 
90     fSys(o.fSys), 
91     fSNN(o.fSNN), 
92     fField(o.fField),
93     fMC(o.fMC), 
94     fSatellite(o.fSatellite), 
95     fData(o.fData),
96     fTimestamp(0), 
97     fAliROOTRevision(0),
98     fAuthor(o.fAuthor)
99 {
100   // 
101   // Copy constructor 
102   //
103 }
104 //____________________________________________________________________
105 AliOADBForward::Entry&
106 AliOADBForward::Entry::operator=(const Entry& o)
107 {
108   // 
109   // Assignment operator 
110   //
111   if (this == &o) return *this;
112   fRunNo              = o.fRunNo; 
113   fSys                = o.fSys; 
114   fSNN                = o.fSNN; 
115   fField              = o.fField;
116   fMC                 = o.fMC; 
117   fSatellite          = o.fSatellite;
118   fData               = o.fData;
119   fTimestamp          = o.fTimestamp;
120   fAliROOTRevision    = o.fAliROOTRevision;
121   fAuthor             = o.fAuthor;
122
123   return *this;
124 }
125 //____________________________________________________________________
126 const char*
127 AliOADBForward::Entry::GetTitle() const
128 {
129   TDatime d(fTimestamp);
130   return Form("%09ld, %4s, %4huGeV, %+2hd, %4s, %3s, %19s: %p (%s) by %s", 
131               (fRunNo == 0xFFFFFFFF ? -1 : fRunNo), 
132               (fSys == 1 ? "pp" : fSys == 2 ? "PbPb" : fSys == 3 ? "pPb" :"?"), 
133               fSNN, fField, (fMC ? "mc" : "real"),
134               (fSatellite ? "sat" : "nom"), d.AsSQLString(), fData,
135               (fData ? fData->GetName() : "?"), fAuthor.Data());
136   
137 }
138 //____________________________________________________________________
139 void
140 AliOADBForward::Entry::Print(Option_t* /*option*/) const 
141 {
142   Printf(GetTitle());
143 }
144 //====================================================================
145 AliOADBForward::Table::Table(TTree* tree, Bool_t isNew, ERunSelectMode mode)
146   : fTree(tree), fEntry(0), fVerbose(false), fMode(mode), fFallBack(false)
147 {
148   if (!tree) return;
149
150   if (isNew) {
151     fTree->Branch("e", "AliOADBForward::Entry", &fEntry);
152     fMode = String2Mode(fTree->GetTitle());
153   }
154   else {
155     if (fMode <= kDefault || fMode > kNewer) {
156       fMode = String2Mode(fTree->GetTitle());
157       if (fMode == kDefault) fMode = kNear;
158     }
159     fTree->SetBranchAddress("e", &fEntry);
160   }
161 }
162 //____________________________________________________________________
163 AliOADBForward::Table::Table(const Table& o)
164   : TObject(o), 
165     fTree(o.fTree), 
166     fEntry(o.fEntry), 
167     fVerbose(o.fVerbose),
168     fMode(o.fMode), 
169     fFallBack(o.fFallBack)
170 {
171   //
172   // Copy constructor 
173   if (!fTree) return;
174   fTree->SetBranchAddress("e", &fEntry);
175 }
176
177 //____________________________________________________________________
178 AliOADBForward::Table::~Table()
179 {
180   // 
181   // Close this table 
182   //
183   Close();
184 }
185 //____________________________________________________________________
186 AliOADBForward::Table&
187 AliOADBForward::Table::operator=(const Table& o)
188 {
189   // 
190   // Assignment operator 
191   // 
192   if (this == &o) return *this;
193   fTree    = o.fTree;
194   fEntry   = o.fEntry;
195   fVerbose = o.fVerbose;
196   fMode    = o.fMode;
197   if (fTree) fTree->SetBranchAddress("e", &fEntry);
198
199   return *this;
200 }
201
202 //____________________________________________________________________
203 const Char_t*
204 AliOADBForward::Table::GetTableName() const
205 {
206   // 
207   // Get the table name or null
208   if (!fTree) return 0;
209   return fTree->GetName();
210 }
211 //____________________________________________________________________
212 const Char_t*
213 AliOADBForward::Table::GetName() const
214 {
215   // 
216   // Overload of TObject::GetName
217   //
218   if (!fTree) return TObject::GetName();
219   return GetTableName();
220 }
221 //____________________________________________________________________
222 Bool_t
223 AliOADBForward::Table::Update()
224 {
225   // 
226   // Flush to disk 
227   //
228   if (!IsOpen()) { 
229     Error("Update", "No tree associated");
230     return false;
231   }
232       
233   TFile* file = fTree->GetCurrentFile();
234   if (!file) { 
235     Error("Update", "No file associated with tree");
236     return false;
237   }
238   if (!file->IsWritable()) { 
239     Error("Update", "File %s not writeable", file->GetName());
240     return false;
241   }
242       
243   Int_t nBytes = file->Write();
244       
245   return (nBytes >= 0);
246 }
247 //____________________________________________________________________
248 Bool_t
249 AliOADBForward::Table::Close()
250 {
251   // 
252   // Close the connection 
253   //
254   if (!IsOpen()) { 
255     Error("Close", "No tree associated");
256     return false;
257   }
258   
259   // if (fTree)  delete fTree; 
260   // if (fEntry) delete fEntry;
261   fTree  = 0;
262   fEntry = 0;
263   return true;
264
265 //____________________________________________________________________
266 Int_t
267 AliOADBForward::Table::Query(ULong_t        runNo,
268                              ERunSelectMode mode,
269                              UShort_t       sys,
270                              UShort_t       sNN, 
271                              Short_t        fld,
272                              Bool_t         mc,
273                              Bool_t         sat) const
274 {
275   // 
276   // Query the tree 
277   //
278   return Query(runNo, mode, Conditions(sys, sNN, fld, mc, sat));
279 }
280
281 //____________________________________________________________________
282 Int_t
283 AliOADBForward::Table::Query(ULong_t        runNo,
284                              ERunSelectMode mode,
285                              const TString& q) const
286 {
287   // 
288   // Run a query against the table 
289   // 
290   
291   // Check the tree 
292   if (!IsOpen()) { 
293     Error("Close", "No tree associated");
294     return -1;
295   }
296       
297   TString query = q;
298   const char* smode = "latest";
299   if (runNo > 0) {
300     if (mode <= kDefault || mode > kNewer) mode = fMode;
301     smode = Mode2String(mode);
302     switch (mode) { 
303     case kExact:  
304       AppendToQuery(query, Form("fRunNo == %lu", runNo)); 
305       break;
306     case kNewest: 
307       break;
308     case kNear:   
309       AppendToQuery(query, Form("abs(fRunNo-%lu)<=%d",
310                                 runNo,kMaxNearDistance)); 
311       break;
312     case kOlder: 
313       AppendToQuery(query, Form("fRunNo <= %lu", runNo));
314       break;
315     case kNewer: 
316       AppendToQuery(query, Form("fRunNo >= %lu", runNo));
317       break;
318     case kDefault: 
319       Fatal("Query", "Mode should never be 'default'");
320       break;
321     }
322   }
323       
324   if (query.IsNull()) {
325     Warning("Query", "Empty query!");
326     return -1;
327   }
328
329   if (fVerbose) 
330     Printf("%s: Query is '%s'", GetName(), query.Data());
331   fTree->Draw("Entry$:fRunNo:fTimestamp", query, "goff");
332   Int_t nRows = fTree->GetSelectedRows();
333   if (nRows <= 0) return -1;
334       
335   if (fVerbose) 
336     Printf("Query: %s (%s)\n"
337            " Entry  |    Run    | Timestamp \n"
338            "--------+-----------+------------------------", 
339            query.Data(), smode);
340       
341   ULong_t  oldRun  = (mode == kNewer ? 0xFFFFFFFF : 0);
342   ULong_t  oldTim  = 0;
343   ULong_t  oldDist = 0xFFFFFFFF;
344   Int_t    entry  = -1;
345   for (Int_t row = 0; row < nRows; row++) {
346     Int_t    ent  = fTree->GetV1()[row];
347     ULong_t  run  = fTree->GetV2()[row];
348     ULong_t  tim  = fTree->GetV3()[row];
349     ULong_t  dist = (run > runNo ? run - runNo : runNo - run);
350         
351     if (fVerbose) {
352       TDatime t(tim);
353       Printf(" %6d | %9ld | %19s ", ent, run > 0x7FFFFFFF ? -1 : run, 
354              t.AsSQLString());
355     }
356
357     switch (mode) { 
358     case kExact: break; // Done in the draw `query' 
359     case kNewest: // Fall-through 
360     case kOlder: 
361       if (run < oldRun) continue;
362       break;
363     case kNewer: 
364       if (run > oldRun) continue;
365       break;
366     case kNear: 
367       if (runNo > 0 && dist > oldDist) continue;
368       break;
369     case kDefault:
370       break;
371     }
372     // If we get here, then we have the best run according to the 
373     // criteria 
374             
375     // Finally, check the timestamp 
376     if (tim < oldTim) continue;
377         
378     // Now update last values and set current best entry 
379     oldTim  = tim;
380     oldDist = dist;
381     oldRun  = run;
382     entry   = ent;
383   }
384
385   if (fVerbose) {
386     Printf("Returning entry # %d", entry);
387   }
388   return entry;
389 }
390
391 //____________________________________________________________________
392 Bool_t
393 AliOADBForward::Table::Insert(TObject* o, 
394                               ULong_t  runNo, 
395                               UShort_t sys, 
396                               UShort_t sNN, 
397                               Short_t  field, 
398                               Bool_t   mc, 
399                               Bool_t   sat,
400                               ULong_t  aliRev,
401                               const TString& author) 
402 {
403   // 
404   // Insert a new row in the table 
405   //
406
407   // Check if the file is read/write 
408   if (fVerbose) 
409     Info("Insert", "Inserting object %p for run=%lu, sys=%hu, sNN=%4hu, "
410          "field=%+2hd, mc=%d, sat=%d", o,runNo, sys, sNN, field, mc, sat);
411
412   if (!IsOpen(true)) {
413     Warning("Insert", "No tree, or not write-able");
414     return false;
415   }
416       
417   // If the entry doesn't exists 
418   if (!fEntry) fEntry = new Entry;
419
420   // Make author 
421   TString auth(author);
422   if (auth.IsNull()) { 
423     UserGroup_t* u = gSystem->GetUserInfo();
424     TInetAddress i = gSystem->GetHostByName(gSystem->HostName());
425     auth = TString::Format("%s <%s@%s>", u->fRealName.Data(), 
426                            u->fUser.Data(), i.GetHostName());
427   }
428     
429   // Set fields 
430   fEntry->fData            = o;
431   fEntry->fRunNo           = runNo; // (runNo <= 0 ? 0xFFFFFFFF : runNo);
432   fEntry->fSys             = sys;
433   fEntry->fSNN             = sNN;
434   fEntry->fField           = field;
435   fEntry->fMC              = mc;
436   fEntry->fSatellite       = sat;
437   fEntry->fAliROOTRevision = (aliRev != 0 ? aliRev : ALIROOT_SVN_REVISION);
438   fEntry->fAuthor          = auth;
439
440   TDatime now;
441   fEntry->fTimestamp       = now.Convert(true);
442
443   // Fill into tree 
444   Int_t nBytes = fTree->Fill();
445   if (nBytes <= 0) {
446     Warning("Insert", "Failed to insert new entry");
447     return false;
448   }
449     
450   // do an Auto-save and flush-baskets now 
451   fTree->AutoSave("FlushBaskets SaveSelf");
452
453   return true;
454 }
455
456 //____________________________________________________________________
457 AliOADBForward::Entry*
458 AliOADBForward::Table::Get(ULong_t        run,
459                            ERunSelectMode mode,
460                            UShort_t       sys,
461                            UShort_t       sNN, 
462                            Short_t        fld,
463                            Bool_t         mc,
464                            Bool_t         sat) const
465 {
466   // Query the tree for an object.  The strategy is as follows. 
467   // 
468   //  - First query with all fields 
469   //    - If this returns a single entry, return that. 
470   //    - If not, then ignore the run number (if given)
471   //      - If this returns a single entry, return that 
472   //      - If not, and fall-back is enabled, then 
473   //        - Ignore the collision energy (if given) 
474   //          - If this returns a single entry, return that.  
475   //          - If not, ignore all passed values 
476   //            - If this returns a single entry, return that.
477   //            - Otherwise, give up and return null
478   //      - Otherwise, give up and return null
479   //
480   // This allow us to specify default objects for a period, and for
481   // collision system, energy, and field setting.
482   //
483   if (fVerbose)
484     Printf("run=%lu mode=%s sys=%hu sNN=%hu fld=%hd mc=%d sat=%d (fall=%d)",
485            run, Mode2String(mode), sys, sNN, fld, mc, sat, fFallBack);
486   Int_t entry = Query(run, mode, sys, sNN, fld, mc, sat);
487   if (entry < 0 && run > 0) 
488     entry = Query(0, mode, sys, sNN, fld, mc, sat);
489   if (entry < 0 && fFallBack && sNN > 0) {
490     if (fVerbose)
491       Printf("Fall-back enabled, will try without sNN=%d", sNN);
492     entry = Query(run, mode, sys, 0, fld, mc, sat);
493   }
494   if (entry < 0 && fFallBack) {
495     if (fVerbose)
496       Printf("Fall-back enabled, will try without any fields");
497     entry = Query(0, mode, 0, 0, kInvalidField, false, false);  
498   }
499   if (entry < 0) {
500     Warning("Get", "No valid object could be found");
501     return 0;
502   }
503   
504   Int_t nBytes = fTree->GetEntry(entry);
505   if (nBytes <= 0) { 
506     Warning("Get", "Failed to get entry # %d\n", entry);
507     return 0;
508   }
509   if (fVerbose) fEntry->Print();
510   return fEntry;
511 }
512 //____________________________________________________________________
513 TObject*
514 AliOADBForward::Table::GetData(ULong_t        run,
515                                ERunSelectMode mode,
516                                UShort_t       sys,
517                                UShort_t       sNN, 
518                                Short_t        fld,
519                                Bool_t         mc,
520                                Bool_t         sat) const
521 {
522   // 
523   // Get data associated with a row or null. 
524   // See also AliOADBForward::Get
525   // 
526   Entry* e = Get(run, mode, sys, sNN, fld, mc, sat);
527   if (!e) return 0;
528   return e->fData;
529 }
530 //____________________________________________________________________
531 void
532 AliOADBForward::Table::Print(Option_t* option) const
533 {
534   // 
535   // Print the full table 
536   //
537   if (!IsOpen()) return;
538
539   Printf("Table %s (default mode: %s)", GetName(), Mode2String(fMode));
540   Int_t n = fTree->GetEntries();
541   for (Int_t i = 0; i < n; i++) { 
542     fTree->GetEntry(i);
543     printf("%4d/%4d: ", i, n);
544     fEntry->Print(option);
545   }
546 }
547 //____________________________________________________________________
548 void
549 AliOADBForward::Table::Browse(TBrowser* b) 
550 {
551   // Browse this table 
552   if (fTree) b->Add(fTree);
553 }
554 //____________________________________________________________________
555 Bool_t
556 AliOADBForward::Table::IsOpen(Bool_t rw) const 
557
558   if (!fTree) return false;
559   if (!rw)    return true;
560   
561   return fTree->GetCurrentFile()->IsWritable();
562 }
563 //====================================================================
564 AliOADBForward::AliOADBForward() 
565   : TObject(),
566     fTables()
567 {
568   //
569   // Constructor 
570   //
571 }
572 #if 0
573 //____________________________________________________________________
574 AliOADBForward::AliOADBForward(const AliOADBForward& other)
575   : TObject(other), 
576     fTables(other.fTables)
577 {
578   // 
579   // Copy constructor 
580   // 
581 }
582 #endif
583 //____________________________________________________________________
584 AliOADBForward::~AliOADBForward()
585 {
586   // 
587   // Destructor 
588   // 
589   Close();
590 }
591 #if 0
592 //____________________________________________________________________
593 AliOADBForward&
594 AliOADBForward::operator=(const AliOADBForward& other)
595 {
596   // 
597   // Copy constructor 
598   // 
599   if (&other == this) return *this;
600
601   fTables = other.fTables;
602
603   return *this;
604 }
605 #endif 
606
607 //____________________________________________________________________
608 Bool_t
609 AliOADBForward::Open(const  TString& fileName, 
610                      const  TString& tables, 
611                      Bool_t          rw, 
612                      Bool_t          verb,
613                      Bool_t          fallback)
614 {
615   TString  absPath(gSystem->ExpandPathName(fileName));
616   if (absPath.IsNull()) { 
617     Error("Open", "Empty path for tables %s", tables.Data());
618     return false;
619   }
620   TObject* previous = gROOT->GetListOfFiles()->FindObject(absPath);
621   TFile*   file     = 0;
622   if (previous) {
623     file = static_cast<TFile*>(previous);
624   }
625   else { 
626     file = TFile::Open(fileName, (rw ? "UPDATE" : "READ"));
627   }
628   if (!file)  { 
629     Error("Open", "Failed to open %s", GetName());
630     return false;
631   }
632   return Open(file, tables, rw, verb, fallback);
633 }
634
635 //____________________________________________________________________
636 Bool_t
637 AliOADBForward::Open(TFile*         file,
638                      const TString& tables,
639                      Bool_t         rw, 
640                      Bool_t         verb,
641                      Bool_t         fallback) 
642 {
643   // 
644   // Open database file and find or create listed tables 
645   // 
646   if (!file) return false;
647   if (rw && !file->IsWritable()) {
648     Warning("Open", "Read+write access requested, but %s opened read-only",
649             file->GetName());
650     if (file->ReOpen("UPDATE") < 0) { 
651       Error("Open", "Failed to reopen file in read+write access mode");
652       return false;
653     }
654   }
655
656   if (tables.EqualTo("*")) {
657     if (rw) { 
658       Error("Open", "Cannot open with '*' in read/write mode");
659       return false;
660     }
661     TList* l = file->GetListOfKeys();
662     TIter  next(l);
663     TKey*  key = 0;
664     while ((key = static_cast<TKey*>(next()))) { 
665       TClass* cl = gROOT->GetClass(key->GetClassName());
666       if (!cl) continue; 
667       if (!cl->InheritsFrom(TTree::Class())) continue; 
668         
669       Table* t = GetTableFromFile(file, false, key->GetName(), "DEFAULT");
670       if (!t) continue;
671
672       fTables.Add(new TObjString(key->GetName()), t);
673       t->SetVerbose(verb);
674       t->SetEnableFallBack(fallback);
675       if (verb) 
676         Printf("Found table %s. Opened with verbose=%d, fallback=%d", 
677                key->GetName(), verb, fallback);
678     }
679     // file->Close();
680     return true;
681   }
682   TObjArray*  tokens = tables.Tokenize(":,");
683   TObjString* token  = 0;
684   TIter       nextToken(tokens);
685   while ((token = static_cast<TObjString*>(nextToken()))) {
686     TString& tn = token->String();
687     if (tn.IsNull()) continue;
688      
689     TObjArray*  parts = tn.Tokenize("/");
690     TObjString* onam  = static_cast<TObjString*>(parts->At(0));
691     TString&    name  = onam->String();
692     TString     mode  = "DEFAULT";
693     if (parts->GetEntries() > 1) 
694       mode = static_cast<TObjString*>(parts->At(1))->String();
695     mode.ToUpper();
696     Table* t = GetTableFromFile(file, rw, name, mode);
697     if (!t) continue;
698
699     t->SetVerbose(verb);
700     t->SetEnableFallBack(fallback);
701     fTables.Add(onam->Clone(), t);
702       if (verb) 
703         Printf("Has table %s. Opened with verbose=%d, fallback=%d", 
704                name.Data(), verb, fallback);
705
706     delete parts;
707   }
708   delete tokens;
709
710   return true;
711 }
712
713 //____________________________________________________________________
714 Bool_t
715 AliOADBForward::Close()
716 {
717   // 
718   // Flush all tables and close all files 
719   // 
720   TList  files;
721   Int_t nFiles = GetFiles(files);
722   if (nFiles <= 0) { 
723     // Nothing to close 
724     return true;
725   }
726
727   TIter nextFile(&files);
728   TFile* file = 0;
729   while ((file = static_cast<TFile*>(nextFile()))) {
730     // if (file->IsWritable()) file->Write();
731
732     file->Close();
733   }
734     
735   fTables.DeleteAll();
736
737   return true;
738 }
739 //____________________________________________________________________
740 Bool_t
741 AliOADBForward::Update()
742 {
743   // 
744   // Flush all tables 
745   // 
746   TList  files;
747   Int_t nFiles = GetFiles(files);
748   if (nFiles <= 0) { 
749     // Nothing to close 
750     return true;
751   }
752
753   TIter nextFile(&files);
754   TFile* file = 0;
755   Int_t  nBytes = 0;
756   while ((file = static_cast<TFile*>(nextFile()))) {
757     if (!file->IsWritable()) { 
758       Error("Update", "File %s not writeable", file->GetName());
759       continue;
760     }
761
762     nBytes += file->Write();
763   }
764   return (nBytes >= 0);
765 }
766 //____________________________________________________________________
767 AliOADBForward::Entry*
768 AliOADBForward::Get(const TString& table, 
769                     ULong_t        run,
770                     ERunSelectMode mode, 
771                     UShort_t       sys,
772                     UShort_t       sNN, 
773                     Short_t        fld,
774                     Bool_t         mc,
775                     Bool_t         sat) const
776 {
777   // 
778   // Get a row from selected table 
779   // 
780   Table* t = FindTable(table);
781   if (!t) return 0;
782     
783   return t->Get(run, mode, sys, sNN, fld, mc, sat);
784 }
785 //____________________________________________________________________
786 TObject*
787 AliOADBForward::GetData(const TString& table, 
788                         ULong_t        run,
789                         ERunSelectMode mode, 
790                         UShort_t       sys,
791                         UShort_t       sNN, 
792                         Short_t        fld,
793                         Bool_t         mc,
794                         Bool_t         sat) const
795 {
796   Table* t = FindTable(table);
797   if (!t) return 0;
798
799   return t->GetData(run, mode, sys, sNN, fld, mc, sat);
800 }
801 //____________________________________________________________________
802 Bool_t
803 AliOADBForward::Insert(const TString& table, 
804                        TObject* o, 
805                        ULong_t  runNo, 
806                        UShort_t sys, 
807                        UShort_t sNN, 
808                        Short_t  field, 
809                        Bool_t   mc, 
810                        Bool_t   sat,
811                        ULong_t  aliRev,
812                        const TString& author) 
813 {
814   // 
815   // Insert a row into the selected table 
816   //
817   Table* t = FindTable(table);
818   if (!t) return false;
819
820   return t->Insert(o, runNo, sys, sNN, field, mc, sat, aliRev, author);
821 }           
822 //____________________________________________________________________
823 void
824 AliOADBForward::Print(const Option_t* option) const
825 {
826   // 
827   // Print everything 
828   //
829   TIter       nextTable(&fTables);
830   TObjString* key   = 0;
831   Table*      table  = 0;
832   while ((key = static_cast<TObjString*>(nextTable()))) {
833     Printf("Table: %p", key->GetName());
834     table = static_cast<Table*>(fTables.GetValue(key));
835     if (!table) continue;
836     table->Print(option);
837   }
838 }
839
840 //____________________________________________________________________
841 void
842 AliOADBForward::Browse(TBrowser* b) 
843 {
844   // 
845   // Browse this object
846   // 
847   TIter       nextTable(&fTables);
848   TObjString* key   = 0;
849   Table*      table  = 0;
850   while ((key = static_cast<TObjString*>(nextTable()))) {
851     table = static_cast<Table*>(fTables.GetValue(key));
852     if (!table) continue;
853     b->Add(table, key->GetName());
854   }
855 }
856 //____________________________________________________________________
857 AliOADBForward::Table*
858 AliOADBForward::FindTable(const TString& name, Bool_t quite) const
859 {
860   //
861   // Find a table by name 
862   // 
863   TPair* p = static_cast<TPair*>(fTables.FindObject(name));
864   if (!p) {
865     if (!quite)
866       Warning("FindTable", "Table %s not registered", name.Data());
867     return 0; 
868   }
869   return static_cast<Table*>(p->Value());
870 }
871 //____________________________________________________________________
872 Int_t
873 AliOADBForward::GetFiles(TList& files) const
874 {
875   // 
876   // Get all associated files 
877   // 
878   Int_t  ret    = 0;
879   TIter       nextTable(&fTables);
880   TObjString* key   = 0;
881   Table*      table  = 0;
882   while ((key = static_cast<TObjString*>(nextTable()))) {
883     table = static_cast<Table*>(fTables.GetValue(key));
884     if (!table->fTree) continue; 
885
886     TFile* f = table->fTree->GetCurrentFile();
887     if (!f) continue;
888
889     if (files.FindObject(f)) continue;
890     files.Add(f);
891     ret++;
892   }
893   return ret;
894 }
895 //____________________________________________________________________
896 AliOADBForward::Table*
897 AliOADBForward::GetTableFromFile(TFile*         file, 
898                                  Bool_t         rw, 
899                                  const TString& name,
900                                  const TString& mode) const
901 {
902   // 
903   // Get a table from a file, or make a new table 
904   // 
905   if (!file) return 0;
906   if (FindTable(name, true)) return 0;
907
908   TObject* o = file->Get(name);
909   TTree*   t = 0;
910   Bool_t   n = false;
911   if (!o) { 
912     if (!rw) { 
913       // We only fail if in read-only mode 
914       Error("Open", "No such object: %s in %s", name.Data(),
915             file->GetName());
916       return 0;
917     }
918     // Create the tree in the file 
919     t = new TTree(name, mode);
920     t->SetDirectory(file);
921     n = true;
922   }
923   else {
924     // Get the tree, and set the branch
925     t = static_cast<TTree*>(o);
926   }
927   Table* ret = new Table(t, n, String2Mode(mode));
928   return ret;
929 }
930 //____________________________________________________________________
931 void
932 AliOADBForward::AppendToQuery(TString& q, const TString& s, Bool_t andNotOr)
933 {
934   // 
935   // Helper function 
936   // 
937   if (!q.IsNull()) q.Append(andNotOr ? " && " : " || ");
938   q.Append(s);
939 }
940 //____________________________________________________________________
941 TString
942 AliOADBForward::Conditions(UShort_t       sys,
943                            UShort_t       sNN, 
944                            Short_t        fld,
945                            Bool_t         mc,
946                            Bool_t         sat)
947 {
948   // Build query string 
949   TString q;
950   if (sys   > 0)    AppendToQuery(q, Form("fSys == %hu",          sys));
951   if (sNN   > 0)    AppendToQuery(q, Form("abs(fSNN - %hu) < 11", sNN));
952   if (TMath::Abs(fld) < 10) AppendToQuery(q, Form("fField == %hd",fld));
953   // Boolean fields always queried.  In cases where these do not matter, 
954   // we always write down the false value, so we get the correct query 
955   // anyways. 
956   AppendToQuery(q, Form("%sfMC",        (mc  ? " " : "!")));
957   AppendToQuery(q, Form("%sfSatellite", (sat ? " " : "!")));
958
959   return q;
960 }
961
962 //____________________________________________________________________
963 void
964 AliOADBForward::TestGet(AliOADBForward& t, 
965                         const TString& table,
966                         ULong_t        runNo,
967                         ERunSelectMode mode,
968                         UShort_t       sys,
969                         UShort_t       sNN, 
970                         Short_t        fld,
971                         Bool_t         mc,
972                         Bool_t         sat)
973 {
974
975   Printf("=== Test query: t=%s r=%ld s=%d t=%d f=%d m=%d v=%d",
976          table.Data(), runNo, sys, sNN, fld, int(mc), int(sat));
977   AliOADBForward::Entry* e = t.Get(table, runNo, mode, sys, sNN, 
978                                    fld, mc, sat);
979   if (!e) return;
980   e->Print();
981 }
982 //____________________________________________________________________
983 void
984 AliOADBForward::TestInsert(AliOADBForward& t, 
985                            const TString&  table,
986                            ULong_t         runNo,
987                            UShort_t        sys,
988                            UShort_t        sNN, 
989                            Short_t         fld,
990                            Bool_t          mc,
991                            Bool_t          sat)
992 {
993   static Int_t cnt = 0;
994   TString what = TString::Format("%s-%03d", table.Data(), cnt++);
995   Printf("=== Insert: t=%s r=%ld s=%d t=%d f=%d m=%d v=%d w=%s",
996          table.Data(), runNo, sys, sNN, fld, int(mc), int(sat), what.Data());
997   t.Insert(table, new TObjString(what), runNo, sys, sNN, fld, mc, sat);
998   gSystem->Sleep(500);
999 }
1000
1001 //____________________________________________________________________
1002 void
1003 AliOADBForward::Test()
1004 {
1005   AliOADBForward* tt = new AliOADBForward();
1006   if (!tt->Open("db.root", "A,B", true, true)) { 
1007     ::Error("Test", "Failed to open DB");
1008     return;
1009   }
1010   AliOADBForward& t  = *tt;
1011   TestInsert(t, "A", 137161);
1012   TestInsert(t, "A", 137161);
1013   TestInsert(t, "A", 0     );
1014   TestInsert(t, "A", 999999);
1015   TestInsert(t, "A", 137166);
1016
1017
1018   TestInsert(t, "B", 137161);
1019   TestInsert(t, "B", 0     );
1020   t.Print();
1021   t.Close();
1022
1023   if (!t.Open("db.root", "A,B",false,true)) {
1024     ::Error("Test", "Failed to open DB");
1025     return;
1026   }
1027
1028   TestGet(t, "A", 137161);
1029   TestGet(t, "A", 137160);
1030   TestGet(t, "A", 0     );
1031   TestGet(t, "A", 137160, kNewest);
1032   TestGet(t, "A", 137160, kNewer);
1033   TestGet(t, "A", 137168, kOlder);
1034   TestGet(t, "A", 137161, kExact);
1035
1036   new TBrowser("b", tt);
1037 }
1038
1039 //
1040 // EOF
1041 //