]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
AliDCSGenDB/AliTPCGenDBConf/AliTPCGenDBTemp: Hard coded reference to specific storage...
[u/mrichter/AliRoot.git] / TPC / AliTPCPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved.      *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 #include "AliTPCPreprocessor.h"
18 #include "AliShuttleInterface.h"
19
20 #include "AliCDBMetaData.h"
21 #include "AliDCSValue.h"
22 #include "AliLog.h"
23 #include "AliTPCSensorTempArray.h"
24 #include "AliTPCROC.h"
25 #include "AliTPCCalROC.h"
26 #include "AliTPCCalPad.h"
27 #include "AliTPCCalibPedestal.h"
28 #include "AliTPCCalibPulser.h"
29 #include "AliTPCCalibCE.h"
30 #include "TFile.h"
31 #include "TTree.h"
32 #include "TEnv.h"
33
34 #include <TTimeStamp.h>
35
36 const Int_t kValCutTemp = 100;               // discard temperatures > 100 degrees
37 const Int_t kDiffCutTemp = 5;                // discard temperature differences > 5 degrees
38 const TString kPedestalRunType = "PEDESTAL_RUN";  // pedestal run identifier
39 const TString kPulserRunType = "PULSER_RUN";   // pulser run identifier
40
41 //
42 // This class is the SHUTTLE preprocessor for the TPC detector.
43 // It contains several components, this far the part containing
44 // temperatures is implemented
45 //
46
47 ClassImp(AliTPCPreprocessor)
48
49 //______________________________________________________________________________________________
50 AliTPCPreprocessor::AliTPCPreprocessor(AliShuttleInterface* shuttle) :
51   AliPreprocessor("TPC",shuttle),
52   fConfEnv(0), fTemp(0), fPressure(0), fConfigOK(kTRUE), fROC(0)
53 {
54   // constructor
55   fROC = AliTPCROC::Instance();
56 }
57 //______________________________________________________________________________________________
58 // AliTPCPreprocessor::AliTPCPreprocessor(const AliTPCPreprocessor& org) :
59 //   AliPreprocessor(org),
60 //   fConfEnv(0), fTemp(0), fPressure(0), fConfigOK(kTRUE)
61 // {
62 //   // copy constructor not implemented
63 //   //   -- missing underlying copy constructor in AliPreprocessor
64 //
65 //   Fatal("AliTPCPreprocessor", "copy constructor not implemented");
66 //
67 // //  fTemp = new AliTPCSensorTempArray(*(org.fTemp));
68 // }
69
70 //______________________________________________________________________________________________
71 AliTPCPreprocessor::~AliTPCPreprocessor()
72 {
73   // destructor
74
75   delete fTemp;
76   delete fPressure;
77 }
78 //______________________________________________________________________________________________
79 AliTPCPreprocessor& AliTPCPreprocessor::operator = (const AliTPCPreprocessor& )
80 {
81   Fatal("operator =", "assignment operator not implemented");
82   return *this;
83 }
84
85
86 //______________________________________________________________________________________________
87 void AliTPCPreprocessor::Initialize(Int_t run, UInt_t startTime,
88         UInt_t endTime)
89 {
90   // Creates AliTestDataDCS object
91
92   AliPreprocessor::Initialize(run, startTime, endTime);
93
94         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
95                 TTimeStamp(startTime).AsString(),
96                 TTimeStamp(endTime).AsString()));
97
98   // Preprocessor configuration
99
100         AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
101         if (entry) fConfEnv = (TEnv*) entry->GetObject();
102         if ( fConfEnv==0 ) {
103            Log("AliTPCPreprocsessor: Preprocessor Config OCDB entry missing.\n");
104            fConfigOK = kFALSE;
105            return;
106         }
107
108   // Temperature sensors
109
110         TTree *confTree = 0;
111         entry = GetFromOCDB("Config", "Temperature");
112         if (entry) confTree = (TTree*) entry->GetObject();
113         if ( confTree==0 ) {
114            Log("AliTPCPreprocsessor: Temperature Config OCDB entry missing.\n");
115            fConfigOK = kFALSE;
116            return;
117         }
118         fTemp = new AliTPCSensorTempArray(fStartTime, fEndTime, confTree);
119         fTemp->SetValCut(kValCutTemp);
120         fTemp->SetDiffCut(kDiffCutTemp);
121
122   // Pressure sensors
123
124         confTree=0;
125         entry=0;
126         entry = GetFromOCDB("Config", "Pressure");
127         if (entry) confTree = (TTree*) entry->GetObject();
128         if ( confTree==0 ) {
129            Log("AliTPCPreprocsessor: Pressure Config OCDB entry missing.\n");
130            fConfigOK = kFALSE;
131            return;
132         }
133         fPressure = new AliDCSSensorArray(fStartTime, fEndTime, confTree);
134
135 }
136
137 //______________________________________________________________________________________________
138 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
139 {
140   // Fills data into TPC calibrations objects
141
142   // Amanda servers provide information directly through dcsAliasMap
143
144   if (!dcsAliasMap) return 9;
145   if (dcsAliasMap->GetEntries() == 0 ) return 9;
146   if (!fConfigOK) return 9;
147
148   TString runType = GetRunType();
149
150   // Temperature sensors are processed by AliTPCCalTemp
151
152
153   UInt_t tempResult = MapTemperature(dcsAliasMap);
154   UInt_t result=tempResult;
155
156   // Pressure sensors
157
158   UInt_t pressureResult = MapPressure(dcsAliasMap);
159   result += pressureResult;
160
161   // Other calibration information will be retrieved through FXS files
162   //  examples:
163   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
164   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
165   //
166   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
167   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
168
169   // pedestal entries
170
171   if(runType == kPedestalRunType) {
172     Int_t pedestalSource = AliShuttleInterface::kDAQ;
173     TString source = fConfEnv->GetValue("Pedestal","DAQ");
174     source.ToUpper();
175     if ( source == "HLT" ) pedestalSource = AliShuttleInterface::kHLT;
176     if (!GetHLTStatus()) pedestalSource = AliShuttleInterface::kDAQ;
177     UInt_t pedestalResult = ExtractPedestals(pedestalSource);
178     result += pedestalResult;
179
180   }
181
182   // pulser trigger processing
183
184   if( runType == kPulserRunType ) {
185     Int_t pulserSource = AliShuttleInterface::kDAQ;
186     TString source = fConfEnv->GetValue("Pulser","DAQ");
187     source.ToUpper();
188     if ( source == "HLT" ) pulserSource = AliShuttleInterface::kHLT;
189     if (!GetHLTStatus()) pulserSource = AliShuttleInterface::kDAQ;
190     UInt_t pulserResult = ExtractPulser(pulserSource);
191     result += pulserResult;
192
193   }
194
195   // Central Electrode processing
196
197   if( false ) {    // CE input file not generated yet
198     Int_t ceSource = AliShuttleInterface::kDAQ;
199     TString source = fConfEnv->GetValue("CE","DAQ");
200     source.ToUpper();
201     if ( source == "HLT" ) ceSource = AliShuttleInterface::kHLT;
202     if (!GetHLTStatus()) ceSource = AliShuttleInterface::kDAQ;
203     UInt_t ceResult = ExtractCE(ceSource);
204     result += ceResult;
205
206   }
207
208   return result;
209 }
210 //______________________________________________________________________________________________
211 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
212 {
213
214    // extract DCS temperature maps. Perform fits to save space
215
216   UInt_t result=0;
217   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
218   if (map) {
219     fTemp->MakeSplineFit(map);
220     AliInfo(Form("Temperature values extracted, fits performed.\n"));
221   } else {
222     Log("AliTPCPreprocsessor: no temperature map extracted. \n");
223     result=9;
224   }
225   delete map;
226   // Now store the final CDB file
227
228   if ( result == 0 ) {
229         AliCDBMetaData metaData;
230         metaData.SetBeamPeriod(0);
231         metaData.SetResponsible("Haavard Helstrup");
232         metaData.SetComment("Preprocessor AliTPC data base entries.");
233
234         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
235         if ( !storeOK )  result=1;
236
237    }
238
239    return result;
240
241 }
242 //______________________________________________________________________________________________
243
244 UInt_t AliTPCPreprocessor::MapPressure(TMap* dcsAliasMap)
245 {
246
247    // extract DCS temperature maps. Perform fits to save space
248
249   UInt_t result=0;
250   TMap *map = fPressure->ExtractDCS(dcsAliasMap);
251   if (map) {
252     fPressure->MakeSplineFit(map);
253     AliInfo(Form("Pressure values extracted, fits performed.\n"));
254   } else {
255     Log("AliTPCPreprocsessor: no atmospheric pressure map extracted. \n");
256     result=9;
257   }
258   delete map;
259   // Now store the final CDB file
260
261   if ( result == 0 ) {
262         AliCDBMetaData metaData;
263         metaData.SetBeamPeriod(0);
264         metaData.SetResponsible("Haavard Helstrup");
265         metaData.SetComment("Preprocessor AliTPC data base entries.");
266
267         Bool_t storeOK = Store("Calib", "Pressure", fPressure, &metaData, 0, 0);
268         if ( !storeOK ) result=1;
269
270    }
271
272    return result;
273
274 }
275
276
277 //______________________________________________________________________________________________
278
279 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
280 {
281  //
282  //  Read pedestal file from file exchage server
283  //  Keep original entry from OCDB in case no new pedestals are available
284  //
285  AliTPCCalPad *calPadPed=0;
286  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
287  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
288  if ( calPadPed==NULL ) {
289      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
290      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
291  }
292
293  AliTPCCalPad *calPadRMS=0;
294  entry = GetFromOCDB("Calib", "Noise");
295  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
296  if ( calPadRMS==NULL ) {
297      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
298      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
299  }
300
301
302  UInt_t result=0;
303
304  Int_t nSectors = fROC->GetNSectors();
305  TList* list = GetFileSources(sourceFXS,"pedestals");
306  
307  if (list && list->GetEntries()>0) {
308
309 //  loop through all files from LDCs
310
311     UInt_t index = 0;
312     while (list->At(index)!=NULL) {
313      TObjString* fileNameEntry = (TObjString*) list->At(index);
314      if (fileNameEntry!=NULL) {
315         TString fileName = GetFile(sourceFXS, "pedestals",
316                                          fileNameEntry->GetString().Data());
317         TFile *f = TFile::Open(fileName);
318         if (!f) {
319           Log ("Error opening pedestal file.");
320           result =2;
321           break;
322         }
323         AliTPCCalibPedestal *calPed;
324         f->GetObject("AliTPCCalibPedestal",calPed);
325
326         //  replace entries for the sectors available in the present file
327
328         for (Int_t sector=0; sector<nSectors; sector++) {
329            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
330            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
331            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
332            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
333         }
334       }
335      ++index;
336     }  // while(list)
337 //
338 //  Store updated pedestal entry to OCDB
339 //
340     AliCDBMetaData metaData;
341     metaData.SetBeamPeriod(0);
342     metaData.SetResponsible("Haavard Helstrup");
343     metaData.SetComment("Preprocessor AliTPC data base entries.");
344
345     Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
346     if ( !storeOK ) ++result;
347     storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
348     if ( !storeOK ) ++result;
349
350   } else {
351     Log ("Error: no entries!");
352     result = 1;
353   }
354
355   return result;
356 }
357
358 //______________________________________________________________________________________________
359
360
361 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
362 {
363  //
364  //  Read pulser calibration file from file exchage server
365  //  Keep original entry from OCDB in case no new pulser calibration is available
366  //
367  TObjArray    *pulserObjects=0;
368  AliTPCCalPad *pulserTmean=0;
369  AliTPCCalPad *pulserTrms=0;
370  AliTPCCalPad *pulserQmean=0;
371  AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
372  if (entry) pulserObjects = (TObjArray*)entry->GetObject();
373  if ( pulserObjects==NULL ) {
374      Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
375      pulserObjects = new TObjArray;    
376  }
377
378  pulserTmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserTmean");
379  if ( !pulserTmean ) {
380     pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
381     pulserObjects->Add(pulserTmean);
382  }
383  pulserTrms = (AliTPCCalPad*)pulserObjects->FindObject("PulserTrms");
384  if ( !pulserTrms )  { 
385     pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
386     pulserObjects->Add(pulserTrms);
387  }
388  pulserQmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserQmean");
389  if ( !pulserQmean )  { 
390     pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
391     pulserObjects->Add(pulserQmean);
392  }
393
394
395  UInt_t result=0;
396
397  Int_t nSectors = fROC->GetNSectors();
398  TList* list = GetFileSources(sourceFXS,"pulser");
399  
400  if (list && list->GetEntries()>0) {
401
402 //  loop through all files from LDCs
403
404     UInt_t index = 0;
405     while (list->At(index)!=NULL) {
406      TObjString* fileNameEntry = (TObjString*) list->At(index);
407      if (fileNameEntry!=NULL) {
408         TString fileName = GetFile(sourceFXS, "pulser",
409                                          fileNameEntry->GetString().Data());
410         TFile *f = TFile::Open(fileName);
411         if (!f) {
412           Log ("Error opening pulser file.");
413           result =2;
414           break;
415         }
416         AliTPCCalibPulser *calPulser;
417         f->GetObject("AliTPCCalibPulser",calPulser);
418
419         //  replace entries for the sectors available in the present file
420
421         for (Int_t sector=0; sector<nSectors; sector++) {
422            AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
423            if ( rocTmean )  pulserTmean->SetCalROC(rocTmean,sector);
424            AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
425            if ( rocTrms )  pulserTrms->SetCalROC(rocTrms,sector);
426            AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
427            if ( rocQmean )  pulserQmean->SetCalROC(rocQmean,sector);
428         }
429       }
430      ++index;
431     }  // while(list)
432 //
433 //  Store updated pedestal entry to OCDB
434 //
435     AliCDBMetaData metaData;
436     metaData.SetBeamPeriod(0);
437     metaData.SetResponsible("Haavard Helstrup");
438     metaData.SetComment("Preprocessor AliTPC data base entries.");
439
440     Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
441     if ( !storeOK ) ++result;
442     
443   } else {
444     Log ("Error: no entries!");
445     result = 1;
446   }
447
448   return result;
449 }
450
451 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
452 {
453  //
454  //  Read Central Electrode file from file exchage server
455  //  Keep original entry from OCDB in case no new CE calibration is available
456  //
457  TObjArray    *ceObjects=0;
458  AliTPCCalPad *ceTmean=0;
459  AliTPCCalPad *ceTrms=0;
460  AliTPCCalPad *ceQmean=0;
461  AliCDBEntry* entry = GetFromOCDB("Calib", "CE");
462  if (entry) ceObjects = (TObjArray*)entry->GetObject();
463  if ( ceObjects==NULL ) {
464      Log("AliTPCPreprocsessor: No previous TPC central electrode entry available.\n");
465      ceObjects = new TObjArray;    
466  }
467
468  ceTmean = (AliTPCCalPad*)ceObjects->FindObject("CETmean");
469  if ( !ceTmean ) {
470     ceTmean = new AliTPCCalPad("CETmean","CETmean");
471     ceObjects->Add(ceTmean);
472  }
473  ceTrms = (AliTPCCalPad*)ceObjects->FindObject("CETrms");
474  if ( !ceTrms )  { 
475     ceTrms = new AliTPCCalPad("CETrms","CETrms");
476     ceObjects->Add(ceTrms);
477  }
478  ceQmean = (AliTPCCalPad*)ceObjects->FindObject("CEQmean");
479  if ( !ceQmean )  { 
480     ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
481     ceObjects->Add(ceQmean);
482  }
483
484
485  UInt_t result=0;
486
487  Int_t nSectors = fROC->GetNSectors();
488  TList* list = GetFileSources(sourceFXS,"CE");
489  
490  if (list && list->GetEntries()>0) {
491
492 //  loop through all files from LDCs
493
494     UInt_t index = 0;
495     while (list->At(index)!=NULL) {
496      TObjString* fileNameEntry = (TObjString*) list->At(index);
497      if (fileNameEntry!=NULL) {
498         TString fileName = GetFile(sourceFXS, "CE",
499                                          fileNameEntry->GetString().Data());
500         TFile *f = TFile::Open(fileName);
501         if (!f) {
502           Log ("Error opening central electrode file.");
503           result =2;
504           break;
505         }
506         AliTPCCalibCE *calCE;
507         f->GetObject("AliTPCCalibCE",calCE);
508
509         //  replace entries for the sectors available in the present file
510
511         for (Int_t sector=0; sector<nSectors; sector++) {
512            AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
513            if ( rocTmean )  ceTmean->SetCalROC(rocTmean,sector);
514            AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
515            if ( rocTrms )  ceTrms->SetCalROC(rocTrms,sector);
516            AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
517            if ( rocQmean )  ceQmean->SetCalROC(rocQmean,sector);
518         }
519       }
520      ++index;
521     }  // while(list)
522 //
523 //  Store updated pedestal entry to OCDB
524 //
525     AliCDBMetaData metaData;
526     metaData.SetBeamPeriod(0);
527     metaData.SetResponsible("Haavard Helstrup");
528     metaData.SetComment("Preprocessor AliTPC data base entries.");
529
530     Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
531     if ( !storeOK ) ++result;
532     
533   } else {
534     Log ("Error: no entries!");
535     result = 1;
536   }
537
538   return result;
539 }