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