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