]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
DAs updated, including documentation lines
[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 kStandAloneRunType = "STANDALONE"; // standalone run identifier
42 const TString kDaqRunType = "DAQ"; // DAQ run identifier
43 const TString kAmandaTemp = "tpc_PT_%d.Temperature"; // Amanda string for temperature entries
44 const Double_t kFitFraction = 0.7;                 // Fraction of DCS sensor fits required              
45
46 //
47 // This class is the SHUTTLE preprocessor for the TPC detector.
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 || runType == kStandAloneRunType || 
208       runType == kDaqRunType ) {    
209     Int_t ceSource = AliShuttleInterface::kDAQ;
210     TString source = fConfEnv->GetValue("CE","DAQ");
211     source.ToUpper();
212     if ( source == "HLT" ) ceSource = AliShuttleInterface::kHLT;
213     if (!GetHLTStatus()) ceSource = AliShuttleInterface::kDAQ;
214     UInt_t ceResult = ExtractCE(ceSource);
215     result += ceResult;
216
217   }
218
219   return result;
220 }
221 //______________________________________________________________________________________________
222 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
223 {
224
225    // extract DCS temperature maps. Perform fits to save space
226
227   UInt_t result=0;
228   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
229   if (map) {
230     fTemp->MakeSplineFit(map);
231     Double_t fitFraction = fTemp->NumFits()/fTemp->NumSensors(); 
232     if (fitFraction > kFitFraction ) {
233       AliInfo(Form("Temperature values extracted, fits performed.\n"));
234     } else { 
235       Log ("Too few temperature maps fitted. \n");
236       result = 9;
237     }
238   } else {
239     Log("No temperature map extracted. \n");
240     result=9;
241   }
242   delete map;
243   // Now store the final CDB file
244
245   if ( result == 0 ) {
246         AliCDBMetaData metaData;
247         metaData.SetBeamPeriod(0);
248         metaData.SetResponsible("Haavard Helstrup");
249         metaData.SetComment("Preprocessor AliTPC data base entries.");
250
251         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
252         if ( !storeOK )  result=1;
253
254    }
255
256    return result;
257
258 }
259
260 //______________________________________________________________________________________________
261 UInt_t AliTPCPreprocessor::MapHighVoltage(TMap* dcsAliasMap)
262 {
263
264    // extract DCS HV maps. Perform fits to save space
265
266   UInt_t result=0;
267   TMap *map = fHighVoltage->ExtractDCS(dcsAliasMap);
268   if (map) {
269     fHighVoltage->MakeSplineFit(map);
270     Double_t fitFraction = fHighVoltage->NumFits()/fHighVoltage->NumSensors(); 
271     if (fitFraction > kFitFraction ) {
272       AliInfo(Form("High volatge recordings extracted, fits performed.\n"));
273     } else { 
274       Log ("Too few high voltage recordings fitted. \n");
275       result = 9;
276     }
277   } else {
278     Log("No high voltage recordings extracted. \n");
279     result=9;
280   }
281   delete map;
282   // Now store the final CDB file
283
284   if ( result == 0 ) {
285         AliCDBMetaData metaData;
286         metaData.SetBeamPeriod(0);
287         metaData.SetResponsible("Haavard Helstrup");
288         metaData.SetComment("Preprocessor AliTPC data base entries.");
289
290         Bool_t storeOK = Store("Calib", "HighVoltage", fHighVoltage, &metaData, 0, kFALSE);
291         if ( !storeOK )  result=1;
292
293    }
294
295    return result;
296
297 }
298
299
300 //______________________________________________________________________________________________
301
302 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
303 {
304  //
305  //  Read pedestal file from file exchage server
306  //  Keep original entry from OCDB in case no new pedestals are available
307  //
308  AliTPCCalPad *calPadPed=0;
309  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
310  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
311  if ( calPadPed==NULL ) {
312      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
313      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
314  }
315
316  AliTPCCalPad *calPadRMS=0;
317  entry = GetFromOCDB("Calib", "Noise");
318  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
319  if ( calPadRMS==NULL ) {
320      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
321      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
322  }
323
324
325  UInt_t result=0;
326
327  Int_t nSectors = fROC->GetNSectors();
328  TList* list = GetFileSources(sourceFXS,"pedestals");
329  
330  if (list && list->GetEntries()>0) {
331
332 //  loop through all files from LDCs
333
334     Bool_t changed=false;
335     UInt_t index = 0;
336     while (list->At(index)!=NULL) {
337      TObjString* fileNameEntry = (TObjString*) list->At(index);
338      if (fileNameEntry!=NULL) {
339         TString fileName = GetFile(sourceFXS, "pedestals",
340                                          fileNameEntry->GetString().Data());
341         TFile *f = TFile::Open(fileName);
342         if (!f) {
343           Log ("Error opening pedestal file.");
344           result =2;
345           break;
346         }
347         AliTPCCalibPedestal *calPed;
348         f->GetObject("AliTPCCalibPedestal",calPed);
349         if ( !calPed ) {
350           Log ("No pedestal calibration object in file.");
351           result = 2;
352           break;
353         }
354
355         //  replace entries for the sectors available in the present file
356
357         changed=true;
358         for (Int_t sector=0; sector<nSectors; sector++) {
359            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
360            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
361            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
362            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
363         }
364       }
365      ++index;
366     }  // while(list)
367 //
368 //  Store updated pedestal entry to OCDB
369 //
370     if (changed) {
371      AliCDBMetaData metaData;
372      metaData.SetBeamPeriod(0);
373      metaData.SetResponsible("Haavard Helstrup");
374      metaData.SetComment("Preprocessor AliTPC data base entries."); 
375  
376      Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
377      if ( !storeOK ) ++result;
378      storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
379      if ( !storeOK ) ++result;
380     }
381   } else {
382     Log ("Error: no entries in input file list!");
383     result = 1;
384   }
385
386   return result;
387 }
388
389 //______________________________________________________________________________________________
390
391
392 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
393 {
394  //
395  //  Read pulser calibration file from file exchage server
396  //  Keep original entry from OCDB in case no new pulser calibration is available
397  //
398  TObjArray    *pulserObjects=0;
399  AliTPCCalPad *pulserTmean=0;
400  AliTPCCalPad *pulserTrms=0;
401  AliTPCCalPad *pulserQmean=0;
402  AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
403  if (entry) pulserObjects = (TObjArray*)entry->GetObject();
404  if ( pulserObjects==NULL ) {
405      Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
406      pulserObjects = new TObjArray;    
407  }
408
409  pulserTmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserTmean");
410  if ( !pulserTmean ) {
411     pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
412     pulserObjects->Add(pulserTmean);
413  }
414  pulserTrms = (AliTPCCalPad*)pulserObjects->FindObject("PulserTrms");
415  if ( !pulserTrms )  { 
416     pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
417     pulserObjects->Add(pulserTrms);
418  }
419  pulserQmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserQmean");
420  if ( !pulserQmean )  { 
421     pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
422     pulserObjects->Add(pulserQmean);
423  }
424
425
426  UInt_t result=0;
427
428  Int_t nSectors = fROC->GetNSectors();
429  TList* list = GetFileSources(sourceFXS,"pulser");
430  
431  if (list && list->GetEntries()>0) {
432
433 //  loop through all files from LDCs
434
435     Bool_t changed=false;
436     UInt_t index = 0;
437     while (list->At(index)!=NULL) {
438      TObjString* fileNameEntry = (TObjString*) list->At(index);
439      if (fileNameEntry!=NULL) {
440         TString fileName = GetFile(sourceFXS, "pulser",
441                                          fileNameEntry->GetString().Data());
442         TFile *f = TFile::Open(fileName);
443         if (!f) {
444           Log ("Error opening pulser file.");
445           result =2;
446           break;
447         }
448         AliTPCCalibPulser *calPulser;
449         f->GetObject("AliTPCCalibPulser",calPulser);
450         if ( !calPulser ) {
451           Log ("No pulser calibration object in file.");
452           result = 2;
453           break;
454         }
455
456         //  replace entries for the sectors available in the present file
457
458         changed=true;
459         for (Int_t sector=0; sector<nSectors; sector++) {
460            AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
461            if ( rocTmean )  pulserTmean->SetCalROC(rocTmean,sector);
462            AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
463            if ( rocTrms )  pulserTrms->SetCalROC(rocTrms,sector);
464            AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
465            if ( rocQmean )  pulserQmean->SetCalROC(rocQmean,sector);
466         }
467       }
468      ++index;
469     }  // while(list)
470 //
471 //  Store updated pedestal entry to OCDB
472 //
473     if (changed) {
474      AliCDBMetaData metaData;
475      metaData.SetBeamPeriod(0);
476      metaData.SetResponsible("Haavard Helstrup");
477      metaData.SetComment("Preprocessor AliTPC data base entries.");
478
479      Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
480      if ( !storeOK ) ++result;
481     }  
482   } else {
483     Log ("Error: no entries in input file list!");
484     result = 1;
485   }
486
487   return result;
488 }
489
490 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
491 {
492  //
493  //  Read Central Electrode file from file exchage server
494  //  Keep original entry from OCDB in case no new CE calibration is available
495  //
496  TObjArray    *ceObjects=0;
497  AliTPCCalPad *ceTmean=0;
498  AliTPCCalPad *ceTrms=0;
499  AliTPCCalPad *ceQmean=0;
500  AliCDBEntry* entry = GetFromOCDB("Calib", "CE");
501  if (entry) ceObjects = (TObjArray*)entry->GetObject();
502  if ( ceObjects==NULL ) {
503      Log("AliTPCPreprocsessor: No previous TPC central electrode entry available.\n");
504      ceObjects = new TObjArray;    
505  }
506
507  ceTmean = (AliTPCCalPad*)ceObjects->FindObject("CETmean");
508  if ( !ceTmean ) {
509     ceTmean = new AliTPCCalPad("CETmean","CETmean");
510     ceObjects->Add(ceTmean);
511  }
512  ceTrms = (AliTPCCalPad*)ceObjects->FindObject("CETrms");
513  if ( !ceTrms )  { 
514     ceTrms = new AliTPCCalPad("CETrms","CETrms");
515     ceObjects->Add(ceTrms);
516  }
517  ceQmean = (AliTPCCalPad*)ceObjects->FindObject("CEQmean");
518  if ( !ceQmean )  { 
519     ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
520     ceObjects->Add(ceQmean);
521  }
522
523
524  UInt_t result=0;
525
526  Int_t nSectors = fROC->GetNSectors();
527  TList* list = GetFileSources(sourceFXS,"CE");
528  
529  if (list && list->GetEntries()>0) {
530
531 //  loop through all files from LDCs
532
533     UInt_t index = 0;
534     while (list->At(index)!=NULL) {
535      TObjString* fileNameEntry = (TObjString*) list->At(index);
536      if (fileNameEntry!=NULL) {
537         TString fileName = GetFile(sourceFXS, "CE",
538                                          fileNameEntry->GetString().Data());
539         TFile *f = TFile::Open(fileName);
540         if (!f) {
541           Log ("Error opening central electrode file.");
542           result =2;
543           break;
544         }
545         AliTPCCalibCE *calCE;
546         f->GetObject("AliTPCCalibCE",calCE);
547
548         //  replace entries for the sectors available in the present file
549
550         for (Int_t sector=0; sector<nSectors; sector++) {
551            AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
552            if ( rocTmean )  ceTmean->SetCalROC(rocTmean,sector);
553            AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
554            if ( rocTrms )  ceTrms->SetCalROC(rocTrms,sector);
555            AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
556            if ( rocQmean )  ceQmean->SetCalROC(rocQmean,sector);
557         }
558       }
559      ++index;
560     }  // while(list)
561 //
562 //  Store updated pedestal entry to OCDB
563 //
564     AliCDBMetaData metaData;
565     metaData.SetBeamPeriod(0);
566     metaData.SetResponsible("Haavard Helstrup");
567     metaData.SetComment("Preprocessor AliTPC data base entries.");
568
569     Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
570     if ( !storeOK ) ++result;
571     
572   } else {
573     Log ("Error: no entries!");
574     result = 1;
575   }
576
577   return result;
578 }