]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDBTemp.cxx
Remving, the info moved to READMEtrigger
[u/mrichter/AliRoot.git] / TPC / AliTPCDBTemp.cxx
1 /**
2 .L /afs/cern.ch/user/h/haavard/alice/tpc/temperature/AliTPCDBTemp.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 AliTPCDBTemp db
7 db->Init(run)
8 db->MakeCalib("TempSensor.txt","DCSMap.root",startTime,endTime,run)
9
10
11 **/
12 #include "AliTPCDBTemp.h"
13
14 ClassImp(AliTPCDBTemp)
15
16 const Int_t kValCut = 100;         // discard temperatures > 100 degrees
17 const Int_t kDiffCut = 5;          // discard temperature differences > 5 degrees
18
19 //______________________________________________________________________________________________
20
21 AliTPCDBTemp::AliTPCDBTemp(): 
22    fFirstRun(0),
23    fLastRun(0),
24    fTemperature(0),
25    fStorLoc(0),
26    fCalib(0),
27    fMetaData(0),
28    fConfTree(0)
29 {}
30 //______________________________________________________________________________________________
31
32 AliTPCDBTemp::AliTPCDBTemp(const AliTPCDBTemp& org):
33   TObject(org),
34   fFirstRun(org.fFirstRun),
35   fLastRun(org.fLastRun),
36   fTemperature(0),
37   fStorLoc(0),
38   fCalib(0),
39   fMetaData(0),
40   fConfTree(0)
41 {
42 //
43 //  Copy constructor
44 //
45
46  ((AliTPCDBTemp &) org).Copy(*this);
47 }
48
49 //______________________________________________________________________________________________
50 AliTPCDBTemp::~AliTPCDBTemp(){
51 //
52 // destructor
53 //
54    fCalib->Terminate();
55    delete fTemperature;
56    delete fMetaData;
57    delete fConfTree;
58 }
59
60 //______________________________________________________________________________________________
61 AliTPCDBTemp& AliTPCDBTemp::operator= (const AliTPCDBTemp& org )
62 {
63  //
64  // assignment operator
65  //
66  if (&org == this) return *this;
67
68  new (this) AliTPCDBTemp(org);
69  return *this;
70
71
72 //______________________________________________________________________________________________
73 void AliTPCDBTemp::Copy(TObject &c) const
74 {
75   //
76   // Copy function
77   //
78
79   TObject::Copy(c);
80 }
81
82
83 //______________________________________________________________________________________________
84
85 void AliTPCDBTemp::MakeCalib(const char *fList, const char *fMap,
86                              const TTimeStamp& startTime, 
87                              const TTimeStamp& endTime,
88                              Int_t run )
89 {
90    // The Terminate() function is the last function to be called during
91    // a query. It always runs on the client, it can be used to present
92    // the results graphically or save the results to file.
93
94    AliTPCSensorTempArray *fTemperature = new AliTPCSensorTempArray(fList);
95    fTemperature->SetStartTime(startTime);
96    fTemperature->SetEndTime(endTime);
97    fTemperature->SetValCut(kValCut);
98    fTemperature->SetDiffCut(kDiffCut);
99    TMap* map = SetGraphFile(fMap);
100    if (map) {
101      fTemperature->MakeSplineFit(map);
102    }
103    delete map;
104    map=0;
105    fMap=0;
106
107    SetFirstRun(run);
108    SetLastRun(run);                 
109    StoreObject("TPC/Calib/Temperature",fTemperature, fMetaData);
110 }
111
112 //______________________________________________________________________________________________
113 void AliTPCDBTemp::MakeConfig(const char *file, Int_t firstRun, Int_t lastRun )
114 {
115    //
116    // Store Configuration file to OCDB
117    //
118
119    TTree *tree = ReadListTree(file);
120    SetConfTree(tree);
121    SetFirstRun(firstRun);
122    SetLastRun(lastRun);                     
123    
124    AliCDBMetaData* metaConf=CreateMetaObject("TTree");      
125    StoreObject("TPC/Config/Temperature",fConfTree, metaConf);
126 }
127
128
129 //______________________________________________________________________________________________
130
131 AliCDBMetaData* AliTPCDBTemp::CreateMetaObject(const char* objectClassName)
132 {
133   AliCDBMetaData *md1= new AliCDBMetaData(); 
134   md1->SetObjectClassName(objectClassName);
135   md1->SetResponsible("Haavard Helstrup");
136   md1->SetBeamPeriod(2);
137   md1->SetAliRootVersion(gSystem->Getenv("ARVERSION"));
138   md1->SetComment("Temperature");
139   
140   return md1;
141 }
142 //______________________________________________________________________________________________
143
144 void AliTPCDBTemp::StoreObject(const char* cdbPath, TObject* object, AliCDBMetaData* metaData)
145 {
146
147   AliCDBId id1(cdbPath, fFirstRun, fLastRun); 
148   if (fStorLoc) fStorLoc->Put(object, id1, metaData); 
149 }
150 //______________________________________________________________________________________________
151
152 void AliTPCDBTemp::Init(Int_t run){
153
154 //   Int_t kLastRun=4000;
155    Long64_t longRun;
156    
157    SetFirstRun(run);
158    SetLastRun(run); 
159        
160    InitDB(run);
161    fCalib = AliTPCcalibDB::Instance();    
162    longRun=run;
163    fCalib->SetRun(longRun);
164    fTemperature = fCalib->GetTemperature();
165      
166 }
167 //______________________________________________________________________________________________
168
169 void AliTPCDBTemp::InitDB(Int_t run)
170
171    //   Data base generation
172    
173    char   *CDBpath="local:///afs/cern.ch/alice/tpctest/Calib/";
174
175    fMetaData = CreateMetaObject("AliTPCSensorTempArray");
176    AliCDBManager *man = AliCDBManager::Instance();
177    man->SetDefaultStorage("local:///afs/cern.ch/alice/tpctest/AliRoot/HEAD"); 
178    man->SetRun(run);
179    man->SetSpecificStorage("TPC/*/*","local:///afs/cern.ch/alice/tpctest/Calib");
180    AliCDBEntry *config = man->Get("TPC/Config/Temperature");
181    if (config) fConfTree = (TTree*)config->GetObject();
182    fStorLoc = man->GetStorage(CDBpath);
183    if (!fStorLoc)    return;
184 }
185 //_____________________________________________________________________________
186 TMap* AliTPCDBTemp::SetGraphFile(const char *fname)
187 {
188   // 
189   // Read DCS maps from file given by fname 
190   //
191   TFile file(fname);
192   TMap * map = (TMap*)file.Get("DCSMap");
193   return map;
194 }
195 //______________________________________________________________________________________________
196
197 TClonesArray * AliTPCDBTemp::ReadList(const char *fname) {
198   //
199   // read values from ascii file
200   //
201   TTree* tree = new TTree("tempConf","tempConf");
202   tree->ReadFile(fname,"");
203   TClonesArray *arr = AliTPCSensorTemp::ReadTree(tree);
204   return arr;
205 }
206
207 //______________________________________________________________________________________________
208
209 TTree * AliTPCDBTemp::ReadListTree(const char *fname) {
210   //
211   // read values from ascii file
212   //
213   TTree* tree = new TTree("tempConf","tempConf");
214   tree->ReadFile(fname,"");
215   TClonesArray *arr = AliTPCSensorTemp::ReadTree(tree);
216   arr->Delete();
217   delete arr;
218   return tree;
219 }