]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPreprocessor.cxx
LED DA now also for PHYSICS, not just STANDALONE
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 // EMCAL Preprocessor class. It runs by Shuttle at the end of the run,
20 // calculates stuff to be posted in OCDB
21 //
22 // Author: Boris Polichtchouk, 4 October 2006
23 // Adapted for EMCAL by Gustavo Conesa Balbastre, October 2006
24 // Updated by David Silvermyr May 2008, based on TPC code
25 ///////////////////////////////////////////////////////////////////////////////
26
27 //Root
28 #include "TFile.h"
29 #include "TTree.h"
30 #include "TEnv.h"
31 #include "TParameter.h"
32
33 #include <TTimeStamp.h>
34
35 //AliRoot
36 #include "AliShuttleInterface.h"
37 #include "AliEMCALPreprocessor.h"
38 #include "AliLog.h"
39 #include "AliCDBMetaData.h"
40 #include "AliCaloCalibPedestal.h"
41 #include "AliCaloCalibSignal.h"
42 #include "AliEMCALSensorTempArray.h"
43
44 const Int_t kValCutTemp = 100;               // discard temperatures > 100 degrees
45 const Int_t kDiffCutTemp = 5;                // discard temperature differences > 5 degrees
46 const TString kPedestalRunType = "PEDESTAL";  // pedestal run identifier
47 const TString kPhysicsRunType = "PHYSICS";   // physics run identifier
48 const TString kStandAloneRunType = "STANDALONE_BC"; // standalone run identifier
49 const TString kAmandaTemp = "PT_%02d.Temperature"; // Amanda string for temperature entries
50 //const Double_t kFitFraction = 0.7;                 // Fraction of DCS sensor fits required 
51 const Double_t kFitFraction = -1.0;          // Don't require minimum number of fits during commissioning 
52
53 const TString kMetaResponsible = "David Silvermyr";
54 //legacy comments and return codes from TPC
55 const TString kMetaComment = "Preprocessor AliEMCAL data base entries.";
56 const int kReturnCodeNoInfo = 9;
57 const int kReturnCodeNoObject = 2;
58 const int kReturnCodeNoEntries = 1;
59
60 ClassImp(AliEMCALPreprocessor)
61   
62 //_______________________________________________________________________________________
63 AliEMCALPreprocessor::AliEMCALPreprocessor() :
64   AliPreprocessor("EMC",0),
65   fConfEnv(0), 
66   fTemp(0), 
67   fConfigOK(kTRUE)
68 {
69   //default constructor
70 }
71
72 //_______________________________________________________________________________________
73 AliEMCALPreprocessor::AliEMCALPreprocessor(AliShuttleInterface* shuttle):
74   AliPreprocessor("EMC",shuttle),
75   fConfEnv(0), 
76   fTemp(0), 
77   fConfigOK(kTRUE)
78 {
79   // Constructor AddRunType(kPedestalRunType);
80   
81   // define run types to be processed
82   AddRunType(kPedestalRunType);
83   AddRunType(kPhysicsRunType);
84   AddRunType(kStandAloneRunType);
85 }
86
87 //______________________________________________________________________________________________
88 AliEMCALPreprocessor::AliEMCALPreprocessor(const AliEMCALPreprocessor&  ) :
89   AliPreprocessor("EMCAL",0),
90   fConfEnv(0), fTemp(0), fConfigOK(kTRUE)
91 {
92   Fatal("AliEMCALPreprocessor", "copy constructor not implemented");
93 }
94
95 // assignment operator; use copy ctor to make life easy.
96 //______________________________________________________________________________________________
97 AliEMCALPreprocessor& AliEMCALPreprocessor::operator = (const AliEMCALPreprocessor &source ) 
98 {
99   // assignment operator; use copy ctor
100   if (&source == this) return *this;
101   
102   new (this) AliEMCALPreprocessor(source);
103   return *this;
104 }
105
106 //____________________________________________________________________________
107 AliEMCALPreprocessor::~AliEMCALPreprocessor()
108 {
109   // destructor
110   if (fTemp) delete fTemp;
111 }
112
113 //______________________________________________________________________________________________
114 void AliEMCALPreprocessor::Initialize(Int_t run, UInt_t startTime,
115                                       UInt_t endTime)
116 {
117   // Creates AliTestDataDCS object -- start maps half an hour beforre actual run start
118   UInt_t startTimeLocal = startTime-1800;
119   AliPreprocessor::Initialize(run, startTimeLocal, endTime);
120   
121   AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
122                TTimeStamp((time_t)startTime,0).AsString(),
123                TTimeStamp((time_t)endTime,0).AsString()));
124   
125   // Preprocessor configuration
126   AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
127   if (entry) fConfEnv = (TEnv*) entry->GetObject();
128   if ( fConfEnv==0 ) {
129     Log("AliEMCALPreprocessor: Preprocessor Config OCDB entry missing.\n");
130     fConfigOK = kFALSE;
131     return;
132   }
133   
134   // Temperature sensors
135   TTree *confTree = 0;
136   
137   TString tempConf = fConfEnv->GetValue("Temperature","ON");
138   tempConf.ToUpper();
139   if (tempConf != "OFF" ) {
140     entry = GetFromOCDB("Config", "Temperature");
141     if (entry) confTree = (TTree*) entry->GetObject();
142     if ( confTree==0 ) {
143       Log("AliEMCALPreprocessor: Temperature Config OCDB entry missing.\n");
144       fConfigOK = kFALSE;
145       return;
146     }
147     fTemp = new AliEMCALSensorTempArray(startTimeLocal, fEndTime, confTree, kAmandaTemp);
148     fTemp->SetValCut(kValCutTemp);
149     fTemp->SetDiffCut(kDiffCutTemp);
150   }
151   
152   return;
153 }
154
155 //______________________________________________________________________________________________
156 UInt_t AliEMCALPreprocessor::Process(TMap* dcsAliasMap)
157 {
158   // Fills data into EMCAL calibrations objects
159   // Amanda servers provide information directly through dcsAliasMap
160   
161   if (!fConfigOK) return kReturnCodeNoInfo;
162   UInt_t result = 0;
163   TObjArray *resultArray = new TObjArray();
164   TString errorHandling = fConfEnv->GetValue("ErrorHandling","ON");
165   errorHandling.ToUpper();
166   TObject * status;
167   
168   UInt_t dcsResult=0;
169   if (errorHandling == "OFF" ) {
170     if (!dcsAliasMap) dcsResult = kReturnCodeNoEntries;
171     if (dcsAliasMap->GetEntries() == 0 ) dcsResult = kReturnCodeNoEntries;  
172     status = new TParameter<int>("dcsResult",dcsResult);
173     resultArray->Add(status);
174   } 
175   else {
176     if (!dcsAliasMap) return kReturnCodeNoInfo;
177     if (dcsAliasMap->GetEntries() == 0 ) return kReturnCodeNoInfo;
178   }
179   
180   TString runType = GetRunType();
181   
182   // Temperature sensors are processed by AliEMCALCalTemp
183   TString tempConf = fConfEnv->GetValue("Temperature","ON");
184   tempConf.ToUpper();
185   if (tempConf != "OFF" ) {
186     UInt_t tempResult = MapTemperature(dcsAliasMap);
187     result=tempResult;
188     status = new TParameter<int>("tempResult",tempResult);
189     resultArray->Add(status);
190   }
191   
192   // Other calibration information will be retrieved through FXS files
193   //  examples:
194   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
195   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
196   //
197   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
198   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
199   
200   // PEDESTAL ENTRIES:
201   
202   if (runType == kPedestalRunType || runType == kStandAloneRunType) {
203     Int_t numSources = 1;
204     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ, AliShuttleInterface::kHLT} ;
205     TString source = fConfEnv->GetValue("Pedestal","DAQ");
206     source.ToUpper();
207     if (source != "OFF" ) { 
208       if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
209       if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
210       if (source == "HLTDAQ" ) {
211         numSources=2;
212         pedestalSource[0] = AliShuttleInterface::kHLT;
213         pedestalSource[1] = AliShuttleInterface::kDAQ;
214       }
215       if (source == "DAQHLT" ) numSources=2;
216       UInt_t pedestalResult=0;
217       for (Int_t i=0; i<numSources; i++ ) {     
218         pedestalResult = ExtractPedestals(pedestalSource[i]);
219         if ( pedestalResult == 0 ) break;
220       }
221       result += pedestalResult;
222       status = new TParameter<int>("pedestalResult",pedestalResult);
223       resultArray->Add(status);
224     }
225   }
226   
227   // SIGNAL/LED ENTRIES:
228   if( runType == kPhysicsRunType || runType == kStandAloneRunType ) {
229     Int_t numSources = 1;
230     Int_t signalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
231     TString source = fConfEnv->GetValue("Signal","DAQ");
232     source.ToUpper();
233     if ( source != "OFF") { 
234       if ( source == "HLT") signalSource[0] = AliShuttleInterface::kHLT;
235       if (!GetHLTStatus()) signalSource[0] = AliShuttleInterface::kDAQ;
236       if (source == "HLTDAQ" ) {
237         numSources=2;
238         signalSource[0] = AliShuttleInterface::kHLT;
239         signalSource[1] = AliShuttleInterface::kDAQ;
240       }
241       if (source == "DAQHLT" ) numSources=2;
242       UInt_t signalResult=0;
243       for (Int_t i=0; i<numSources; i++ ) {     
244         signalResult = ExtractSignal(signalSource[i]);
245         if ( signalResult == 0 ) break;
246       }
247       result += signalResult;
248       status = new TParameter<int>("signalResult",signalResult);
249       resultArray->Add(status);
250     }
251   }
252   
253   
254   // overall status at the end
255   if (errorHandling == "OFF" ) {
256     AliCDBMetaData metaData;
257     metaData.SetBeamPeriod(0);
258     metaData.SetResponsible(kMetaResponsible);
259     metaData.SetComment("Preprocessor AliEMCAL status.");
260     Store("Calib", "PreprocStatus", resultArray, &metaData, 0, kFALSE);
261     resultArray->Delete();
262     return 0;
263   } 
264   else { 
265     return result;
266   }
267   
268 }
269 //______________________________________________________________________________________________
270 UInt_t AliEMCALPreprocessor::MapTemperature(TMap* dcsAliasMap)
271 {
272   // extract DCS temperature maps. Perform fits to save space
273   UInt_t result=0;
274
275   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
276   if (map) {
277     fTemp->MakeSplineFit(map);
278     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
279     if (fitFraction > kFitFraction ) {
280       AliInfo(Form("Temperature values extracted, fits performed.\n"));
281     } 
282     else { 
283       Log ("Too few temperature maps fitted. \n");
284       result = kReturnCodeNoInfo;
285     }
286   } 
287   else {
288     Log("No temperature map extracted. \n");
289     result = kReturnCodeNoInfo;
290   }
291   delete map;
292   // Now store the final CDB file
293   
294   if ( result == 0 ) { // some info was found
295     AliCDBMetaData metaData;
296     metaData.SetBeamPeriod(0);
297     metaData.SetResponsible(kMetaResponsible);
298     metaData.SetComment(kMetaComment);
299     
300     Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
301     if ( !storeOK )  result=1;
302   }
303   
304   return result;
305 }
306
307 //______________________________________________________________________________________________
308 UInt_t AliEMCALPreprocessor::ExtractPedestals(Int_t sourceFXS)
309 {
310   UInt_t result=0;
311   //
312   //  Read pedestal file from file exchange server
313   //  Only store if new pedestal info is available
314   //
315   AliCaloCalibPedestal *calibPed = new AliCaloCalibPedestal(AliCaloCalibPedestal::kEmCal);
316   
317   TList* list = GetFileSources(sourceFXS,"pedestals");
318   if (list && list->GetEntries()>0) {
319     
320     //  loop through all files from LDCs
321
322     int changes = 0;
323     UInt_t index = 0;
324     while (list->At(index)!=NULL) {
325       TObjString* fileNameEntry = (TObjString*) list->At(index);
326       if (fileNameEntry!=NULL) {
327         TString fileName = GetFile(sourceFXS, "pedestals",
328                                    fileNameEntry->GetString().Data());
329         TFile *f = TFile::Open(fileName);
330         if (!f) {
331           Log ("Error opening pedestal file.");
332           result = kReturnCodeNoObject;
333           break;
334         }
335         AliCaloCalibPedestal *calPed;
336         f->GetObject("emcCalibPedestal",calPed);
337         if ( !calPed ) {
338           Log ("No pedestal calibration object in file.");
339           result = kReturnCodeNoObject;
340           break;
341         }
342         if ( calPed->GetNEvents()>0 && calPed->GetNChanFills()>0 ) {
343           // add info for the modules available in the present file
344           Bool_t status = calibPed->AddInfo(calPed);
345           if (status) { changes++; }
346         }
347         
348         delete calPed; 
349         f->Close();
350       }
351       index++;
352     }  // while(list)
353     
354     //
355     //  Store updated pedestal entry to OCDB
356     //
357     if (changes>0) {
358       AliCDBMetaData metaData;
359       metaData.SetBeamPeriod(0);
360       metaData.SetResponsible(kMetaResponsible);
361       metaData.SetComment(kMetaComment); 
362       
363       Bool_t storeOK = StoreReferenceData("Calib", "Pedestals", calibPed, &metaData);
364       if ( !storeOK ) result++;
365     }
366   } 
367   else {
368     Log ("Error: no entries in input file list!");
369     result = kReturnCodeNoEntries;
370   }
371   
372   return result;
373 }
374
375 //______________________________________________________________________________________________
376 UInt_t AliEMCALPreprocessor::ExtractSignal(Int_t sourceFXS)
377 {
378   UInt_t result=0;
379   //
380   //  Read signal file from file exchange server
381   //  Only store if new signal info is available
382   //
383   AliCaloCalibSignal *calibSig = new AliCaloCalibSignal(AliCaloCalibSignal::kEmCal); 
384   
385   TList* list = GetFileSources(sourceFXS,"signal");
386   if (list && list->GetEntries()>0) {
387
388     //  loop through all files from LDCs
389     
390     int changes = 0;
391     UInt_t index = 0;
392     while (list->At(index)!=NULL) {
393       TObjString* fileNameEntry = (TObjString*) list->At(index);
394       if (fileNameEntry!=NULL) {
395         TString fileName = GetFile(sourceFXS, "signal",
396                                    fileNameEntry->GetString().Data());
397         TFile *f = TFile::Open(fileName);
398         if (!f) {
399           Log ("Error opening signal file.");
400           result = kReturnCodeNoObject;
401           break;
402         }
403         AliCaloCalibSignal *calSig;
404         f->GetObject("emcCalibSignal",calSig);
405         if ( !calSig ) {
406           Log ("No signal calibration object in file.");
407           result = kReturnCodeNoObject;
408           break;
409         }
410         if ( calSig->GetNEvents()>0 ) {
411           // add info for the modules available in the present file
412           Bool_t status = calibSig->AddInfo(calSig);
413           if (status) { changes++; }
414         }
415         
416         delete calSig; 
417         f->Close();
418       }
419       index++;
420     }  // while(list)
421     
422     //
423     //  Store updated signal entry to OCDB
424     //
425     if (changes>0) {
426       AliCDBMetaData metaData;
427       metaData.SetBeamPeriod(0);
428       metaData.SetResponsible(kMetaResponsible);
429       metaData.SetComment(kMetaComment); 
430       
431       Bool_t storeOK = Store("Calib", "LED", calibSig, &metaData, 0, kFALSE);
432       if ( !storeOK ) result++;
433     }
434   } 
435   else {
436     Log ("Error: no entries in input file list!");
437     result = kReturnCodeNoEntries;
438   }
439
440   return result;
441 }
442
443