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