]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/TRDLoaderSingle.cxx
Fix handling of AliESDfriends.root.
[u/mrichter/AliRoot.git] / EVE / Alieve / TRDLoaderSingle.cxx
1 #include "TRDLoaderSingle.h"
2
3 //#include "AliTRDv1.h"
4
5 #include <Reve/RGTopFrame.h>
6
7 #include "TString.h"
8 #include "TObjString.h"
9 #include "TObjArray.h"
10 #include "TFile.h"
11
12 using namespace Reve;
13 using namespace Alieve;
14 using namespace std;
15
16 ClassImp(Alieve::TRDLoaderSingle)
17
18 ///////////////////////////////////////////////////////////
19 /////////////     TRDLoaderSingle     /////////////////////
20 ///////////////////////////////////////////////////////////
21
22
23 //________________________________________________________
24 TRDLoaderSingle::TRDLoaderSingle(const Text_t* n, const Text_t* t) : TRDLoader(n, t)
25 {
26 }
27
28
29 //________________________________________________________
30 Bool_t  TRDLoaderSingle::GoToEvent(const int ev)
31 {
32         fEvent = ev;
33
34         Unload();
35         
36         TTree *t = 0x0;
37         TFile *f = new TFile(Form("%s/%s", fDir.Data(), fFilename.Data()));
38         if(! f->cd(Form("Event%d", ev))){
39                 Error("GoToEvent()", Form("Could not find event %d in file %s.", ev, fFilename.Data()));
40                 return kFALSE;
41         }
42         
43         if(kLoadHits){
44                 t = (TTree*)gDirectory->Get("TreeH");
45                 if(!t) return kFALSE;
46                 if(!LoadHits(t)) return kFALSE;
47         }
48         if(kLoadDigits){
49                 t = (TTree*)gDirectory->Get("TreeD");
50                 if(!t) return kFALSE;
51                 if(!LoadDigits(t)) return kFALSE;
52         }
53         if(kLoadClusters){
54                 t = (TTree*)gDirectory->Get("TreeR");
55                 if(!t) return kFALSE;
56                 if(!LoadClusters(t)) return kFALSE;
57         }
58         if(kLoadTracks){
59                 t = (TTree*)gDirectory->Get("TreeT");
60                 if(!t) return kFALSE;
61                 if(!LoadTracklets(t)) return kFALSE;
62         }
63         f->Close();
64         
65         gReve->Redraw3D();
66         
67         return kTRUE;
68 }
69
70         
71 //________________________________________________________
72 Bool_t  TRDLoaderSingle::Open(const char *filename, const char *dir)
73 {
74         fFilename = filename;
75         fDir = dir;
76         
77         TObjArray *so = fFilename.Tokenize(".");
78
79         if(((TObjString*)(*so)[0])->GetString().CompareTo("TRD") != 0){
80                 Error("Open()", "Filename didn't fulfill naming conventions. No TRD data.");
81                 return kFALSE;
82         }
83         if(((TObjString*)(*so)[1])->GetString().CompareTo("Hits") == 0) kLoadHits = kTRUE;
84         else if(((TObjString*)(*so)[1])->GetString().CompareTo("Digits") == 0) kLoadDigits = kTRUE;
85         else if(((TObjString*)(*so)[1])->GetString().CompareTo("RecPoints") == 0) kLoadClusters = kTRUE;
86         else if(((TObjString*)(*so)[1])->GetString().CompareTo("Tracks") == 0) kLoadTracks = kTRUE;
87         else{
88                 Error("Open()", "Filename didn't fulfill naming conventions. No data type specified.");
89                 return kFALSE;
90         }
91         
92         return kTRUE;
93 }
94