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