]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCGenDBTemp.cxx
Detectors can have more than one AMANDA server. SHUTTLE queries the servers sequentially,
[u/mrichter/AliRoot.git] / TPC / AliTPCGenDBTemp.cxx
1
2 // .L /afs/cern.ch/user/h/haavard/alice/tpc/temperature/AliTPCGenDBTemp.C+
3 // TTimeStamp startTime(2006,10,18,0,0,0,0,kFALSE)
4 // TTimeStamp endTime(2006,10,19,0,0,0,0,kFALSE)
5 // Int_t run=2546
6 // AliTPCGenDBTemp db
7 // db->Init(run,"TPC/Config/Temperature","TPC/*/*")
8 // db->MakeCalib("TempSensor.txt","DCSMap.root",startTime,endTime,run)
9
10
11 #include "AliTPCGenDBTemp.h"
12
13 ClassImp(AliTPCGenDBTemp)
14
15 const Int_t kValCut = 100;         // discard temperatures > 100 degrees
16 const Int_t kDiffCut = 5;          // discard temperature differences > 5 degrees
17
18 //______________________________________________________________________________________________
19
20 AliTPCGenDBTemp::AliTPCGenDBTemp():
21    AliDCSGenDB()
22 {
23 }
24
25 //______________________________________________________________________________________________
26
27 AliTPCGenDBTemp::AliTPCGenDBTemp(const AliTPCGenDBTemp& org):
28   AliDCSGenDB(org)
29 {
30
31 //
32 //  Copy constructor
33 //
34
35  ((AliTPCGenDBTemp &) org).Copy(*this);
36 }
37
38 //______________________________________________________________________________________________
39 AliTPCGenDBTemp::~AliTPCGenDBTemp(){
40 //
41 // destructor
42 //
43
44 }
45 //______________________________________________________________________________________________
46 AliTPCGenDBTemp& AliTPCGenDBTemp::operator= (const AliTPCGenDBTemp& org )
47 {
48  //
49  // assignment operator
50  //
51  if (&org == this) return *this;
52
53  new (this) AliTPCGenDBTemp(org);
54  return *this;
55 }
56
57
58 //______________________________________________________________________________________________
59
60 void AliTPCGenDBTemp::MakeCalib(const char *fList, const char *fMap,
61                              const TTimeStamp& startTime,
62                              const TTimeStamp& endTime,
63                              Int_t run )
64 {
65    // The Terminate() function is the last function to be called during
66    // a query. It always runs on the client, it can be used to present
67    // the results graphically or save the results to file.
68
69    AliTPCSensorTempArray *temperature = new AliTPCSensorTempArray(fList);
70    temperature->SetStartTime(startTime);
71    temperature->SetEndTime(endTime);
72    temperature->SetValCut(kValCut);
73    temperature->SetDiffCut(kDiffCut);
74    TMap* map = SetGraphFile(fMap);
75    if (map) {
76      temperature->MakeSplineFit(map);
77    }
78    delete map;
79    map=0;
80    fMap=0;
81
82    SetFirstRun(run);
83    SetLastRun(run);
84    SetSensorArray(temperature);
85    StoreObject("TPC/Calib/Temperature",temperature, fMetaData);
86 }
87
88 //______________________________________________________________________________________________
89
90 TClonesArray * AliTPCGenDBTemp::ReadList(const char *fname) {
91   //
92   // read values from ascii file
93   //
94   TTree* tree = new TTree("tempConf","tempConf");
95   tree->ReadFile(fname,"");
96   TClonesArray *arr = AliTPCSensorTemp::ReadTree(tree);
97   delete tree;
98   return arr;
99 }
100
101 //______________________________________________________________________________________________
102
103 TTree * AliTPCGenDBTemp::ReadListTree(const char *fname) {
104   //
105   // read values from ascii file
106   //
107   TTree* tree = new TTree("tempConf","tempConf");
108   tree->ReadFile(fname,"");
109   TClonesArray *arr = AliTPCSensorTemp::ReadTree(tree);
110   arr->Delete();
111   delete arr;
112   return tree;
113 }