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