]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPreprocessor.cxx
adding some checks and (tmp) AliInfo lines to help debug nightly shuttle test failure...
[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 "AliDCSValue.h"
40 #include "AliCDBMetaData.h"
41 #include "AliEMCALTriggerTRUDCSConfig.h"
42 #include "AliEMCALTriggerDCSConfig.h"
43 #include "AliCaloCalibPedestal.h"
44 #include "AliCaloCalibSignal.h"
45 #include "AliEMCALSensorTempArray.h"
46
47 const Int_t kValCutTemp = 100;               // discard temperatures > 100 degrees
48 const Int_t kDiffCutTemp = 5;                // discard temperature differences > 5 degrees
49 const TString kPedestalRunType = "PEDESTAL";  // pedestal run identifier
50 const TString kPhysicsRunType = "PHYSICS";   // physics run identifier
51 const TString kStandAloneRunType = "STANDALONE_BC"; // standalone run identifier
52 const TString kAmandaTemp = "PT_%02d.Temperature"; // Amanda string for temperature entries
53 //const Double_t kFitFraction = 0.7;                 // Fraction of DCS sensor fits required 
54 const Double_t kFitFraction = -1.0;          // Don't require minimum number of fits during commissioning 
55
56 const TString kMetaResponsible = "David Silvermyr";
57 //legacy comments and return codes from TPC
58 const TString kMetaComment = "Preprocessor AliEMCAL data base entries.";
59 const int kReturnCodeNoInfo = 9;
60 const int kReturnCodeNoObject = 2;
61 const int kReturnCodeNoEntries = 1;
62
63 const int kNTRU = 30; // From 2011; 10 SuperModules (SM) * 3 TRU per SM
64
65 ClassImp(AliEMCALPreprocessor)
66   
67 //_______________________________________________________________________________________
68 AliEMCALPreprocessor::AliEMCALPreprocessor() :
69   AliPreprocessor("EMC",0),
70   fConfEnv(0), 
71   fTemp(0), 
72   fConfigOK(kTRUE)
73 {
74   //default constructor
75 }
76
77 //_______________________________________________________________________________________
78 AliEMCALPreprocessor::AliEMCALPreprocessor(AliShuttleInterface* shuttle):
79   AliPreprocessor("EMC",shuttle),
80   fConfEnv(0), 
81   fTemp(0), 
82   fConfigOK(kTRUE)
83 {
84   // Constructor AddRunType(kPedestalRunType);
85   
86   // define run types to be processed
87   AddRunType(kPedestalRunType);
88   AddRunType(kPhysicsRunType);
89 }
90
91 //______________________________________________________________________________________________
92 AliEMCALPreprocessor::AliEMCALPreprocessor(const AliEMCALPreprocessor&  ) :
93   AliPreprocessor("EMCAL",0),
94   fConfEnv(0), fTemp(0), fConfigOK(kTRUE)
95 {
96   Fatal("AliEMCALPreprocessor", "copy constructor not implemented");
97 }
98
99 // assignment operator; use copy ctor to make life easy.
100 //______________________________________________________________________________________________
101 AliEMCALPreprocessor& AliEMCALPreprocessor::operator = (const AliEMCALPreprocessor &source ) 
102 {
103   // assignment operator; use copy ctor
104   if (&source == this) return *this;
105   
106   new (this) AliEMCALPreprocessor(source);
107   return *this;
108 }
109
110 //____________________________________________________________________________
111 AliEMCALPreprocessor::~AliEMCALPreprocessor()
112 {
113   // destructor
114   if (fTemp) delete fTemp;
115 }
116
117 //______________________________________________________________________________________________
118 void AliEMCALPreprocessor::Initialize(Int_t run, UInt_t startTime,
119                                       UInt_t endTime)
120 {
121   // Creates AliTestDataDCS object -- start maps half an hour beforre actual run start
122   UInt_t startTimeLocal = startTime-1800;
123   AliPreprocessor::Initialize(run, startTimeLocal, endTime);
124   
125   AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
126                TTimeStamp((time_t)startTime,0).AsString(),
127                TTimeStamp((time_t)endTime,0).AsString()));
128   
129   // Preprocessor configuration
130   AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
131   if (entry) fConfEnv = (TEnv*) entry->GetObject();
132   if ( fConfEnv==0 ) {
133     Log("AliEMCALPreprocessor: Preprocessor Config OCDB entry missing.\n");
134     fConfigOK = kFALSE;
135     return;
136   }
137   
138   // Temperature sensors
139   TTree *confTree = 0;
140   
141   TString tempConf = fConfEnv->GetValue("Temperature","ON");
142   tempConf.ToUpper();
143   if (tempConf != "OFF" ) {
144     entry = GetFromOCDB("Config", "Temperature");
145     if (entry) confTree = (TTree*) entry->GetObject();
146     if ( confTree==0 ) {
147       Log("AliEMCALPreprocessor: Temperature Config OCDB entry missing.\n");
148       fConfigOK = kFALSE;
149       return;
150     }
151     fTemp = new AliEMCALSensorTempArray(startTimeLocal, fEndTime, confTree, kAmandaTemp);
152     fTemp->SetValCut(kValCutTemp);
153     fTemp->SetDiffCut(kDiffCutTemp);
154   }
155   
156   return;
157 }
158
159 //______________________________________________________________________________________________
160 UInt_t AliEMCALPreprocessor::Process(TMap* dcsAliasMap)
161 {
162   // Fills data into EMCAL calibrations objects
163   // Amanda servers provide information directly through dcsAliasMap
164   
165   if (!fConfigOK) return kReturnCodeNoInfo;
166   UInt_t result = 0;
167   TObjArray *resultArray = new TObjArray();
168   TString errorHandling = fConfEnv->GetValue("ErrorHandling","ON");
169   errorHandling.ToUpper();
170   TObject * status;
171   
172   UInt_t dcsResult=0;
173   if (errorHandling == "OFF" ) {
174     if (!dcsAliasMap) dcsResult = kReturnCodeNoEntries;
175     else if (dcsAliasMap->GetEntries() == 0 ) dcsResult = kReturnCodeNoEntries;  
176     status = new TParameter<int>("dcsResult",dcsResult);
177     resultArray->Add(status);
178   } 
179   else {
180     if (!dcsAliasMap) return kReturnCodeNoInfo;
181     else if (dcsAliasMap->GetEntries() == 0 ) return kReturnCodeNoInfo;
182   }
183   
184     
185   TString runType = GetRunType();
186   
187   // Temperature sensors are processed by AliEMCALCalTemp
188   TString tempConf = fConfEnv->GetValue("Temperature","ON");
189   tempConf.ToUpper();
190   if (tempConf != "OFF" && dcsAliasMap ) {
191     UInt_t tempResult = MapTemperature(dcsAliasMap);
192     result=tempResult;
193     status = new TParameter<int>("tempResult",tempResult);
194     resultArray->Add(status);
195   }
196   // Trigger configuration processing: only for Physics runs
197   TString triggerConf = fConfEnv->GetValue("Trigger","ON");
198   triggerConf.ToUpper();
199   if( runType == kPhysicsRunType ) {
200     //  if (triggerConf != "OFF" && dcsAliasMap ) {
201     if ( dcsAliasMap ) {
202       UInt_t triggerResult = MapTriggerConfig(dcsAliasMap);
203       AliInfo(Form("triggerConf %s\n", triggerConf.Data()));
204       result+=triggerResult;
205       status = new TParameter<int>("triggerResult",triggerResult);
206       resultArray->Add(status);
207     }
208   }
209   
210   // Other calibration information will be retrieved through FXS files
211   //  examples:
212   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
213   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
214   //
215   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
216   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
217   
218   // PEDESTAL ENTRIES:
219   
220   if ( runType == kPedestalRunType ) {
221     Int_t numSources = 1;
222     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ, AliShuttleInterface::kHLT} ;
223     TString source = fConfEnv->GetValue("Pedestal","DAQ");
224     source.ToUpper();
225     if (source != "OFF" ) { 
226       if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
227       if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
228       if (source == "HLTDAQ" ) {
229         numSources=2;
230         pedestalSource[0] = AliShuttleInterface::kHLT;
231         pedestalSource[1] = AliShuttleInterface::kDAQ;
232       }
233       if (source == "DAQHLT" ) numSources=2;
234       UInt_t pedestalResult=0;
235       for (Int_t i=0; i<numSources; i++ ) {     
236         pedestalResult = ExtractPedestals(pedestalSource[i]);
237         if ( pedestalResult == 0 ) break;
238       }
239       result += pedestalResult;
240       status = new TParameter<int>("pedestalResult",pedestalResult);
241       resultArray->Add(status);
242     }
243   }
244   
245   // SIGNAL/LED ENTRIES:
246   if( runType == kPhysicsRunType ) {
247     Int_t numSources = 1;
248     Int_t signalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
249     TString source = fConfEnv->GetValue("Signal","DAQ");
250     source.ToUpper();
251     if ( source != "OFF") { 
252       if ( source == "HLT") signalSource[0] = AliShuttleInterface::kHLT;
253       if (!GetHLTStatus()) signalSource[0] = AliShuttleInterface::kDAQ;
254       if (source == "HLTDAQ" ) {
255         numSources=2;
256         signalSource[0] = AliShuttleInterface::kHLT;
257         signalSource[1] = AliShuttleInterface::kDAQ;
258       }
259       if (source == "DAQHLT" ) numSources=2;
260       UInt_t signalResult=0;
261       for (Int_t i=0; i<numSources; i++ ) {     
262         signalResult = ExtractSignal(signalSource[i]);
263         if ( signalResult == 0 ) break;
264       }
265       result += signalResult;
266       status = new TParameter<int>("signalResult",signalResult);
267       resultArray->Add(status);
268     }
269   }
270   
271   
272   // overall status at the end
273   if (errorHandling == "OFF" ) {
274     AliCDBMetaData metaData;
275     metaData.SetBeamPeriod(0);
276     metaData.SetResponsible(kMetaResponsible);
277     metaData.SetComment("Preprocessor AliEMCAL status.");
278     Bool_t storeOK = Store("Calib", "PreprocStatus", resultArray, &metaData, 0, kFALSE);
279     resultArray->Delete();
280     result = 0;
281     if ( !storeOK )  result=1;
282     return result;
283   } 
284   else { 
285     return result;
286   }
287   
288 }
289 //______________________________________________________________________________________________
290 UInt_t AliEMCALPreprocessor::MapTemperature(TMap* dcsAliasMap)
291 { // extract DCS temperature maps. Perform fits to save space
292   UInt_t result=0;
293
294   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
295   if (map) {
296     fTemp->MakeSplineFit(map);
297     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
298     if (fitFraction > kFitFraction ) {
299       AliInfo(Form("Temperature values extracted, fits performed.\n"));
300     } 
301     else { 
302       Log ("Too few temperature maps fitted. \n");
303       result = kReturnCodeNoInfo;
304     }
305   } 
306   else {
307     Log("No temperature map extracted. \n");
308     result = kReturnCodeNoInfo;
309   }
310   delete map;
311   // Now store the final CDB file
312   
313   if ( result == 0 ) { // some info was found
314     AliCDBMetaData metaData;
315     metaData.SetBeamPeriod(0);
316     metaData.SetResponsible(kMetaResponsible);
317     metaData.SetComment(kMetaComment);
318     
319     Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
320     if ( !storeOK )  result=1;
321   }
322   
323   return result;
324 }
325
326 //______________________________________________________________________________________________
327 UInt_t AliEMCALPreprocessor::MapTriggerConfig(TMap* dcsAliasMap)
328 { // extract DCS trigger info
329   AliInfo(Form("Get TRU info from DCS DPs.\n"));
330   Int_t i, iTRU;
331   char buf[100];
332
333   AliDCSValue *dcsVal;
334   TObjArray *arrL0ALGSEL, *arrPEAKFINDER, *arrGLOBALTHRESH, *arrCOSMTHRESH;
335   TObjArray *arrMASK[6];
336
337   // overall object to hold STU and DCS config info
338   // DS comment: for now only holds TRU info, i.e. only partially filled
339   // (STU info only in raw data header; unfortunately not also picked up via DCS DPs)
340   AliEMCALTriggerDCSConfig *trigConfig = new AliEMCALTriggerDCSConfig();
341
342   // loop through all TRUs
343   bool debug = true; // debug flag for AliInfo printouts for each TRU
344   for( iTRU = 0; iTRU < kNTRU; iTRU++){
345     if (debug) AliInfo( Form("iTRU %d \n", iTRU) );
346     // get the shuttled values
347     sprintf( buf, "EMC_TRU%02d_L0ALGSEL", iTRU );
348     arrL0ALGSEL = (TObjArray*) dcsAliasMap->GetValue( buf );
349     sprintf( buf, "EMC_TRU%02d_PEAKFINDER", iTRU );
350     arrPEAKFINDER = (TObjArray*) dcsAliasMap->GetValue( buf );
351     sprintf( buf, "EMC_TRU%02d_GLOBALTHRESH", iTRU );
352     arrGLOBALTHRESH = (TObjArray*) dcsAliasMap->GetValue( buf );
353     sprintf( buf, "EMC_TRU%02d_COSMTHRESH", iTRU );
354     arrCOSMTHRESH = (TObjArray*) dcsAliasMap->GetValue( buf );
355     
356     for( i = 0; i < 6; i++ ){
357       sprintf( buf, "EMC_TRU%02d_MASK%d", iTRU, i );
358       arrMASK[i] = (TObjArray*) dcsAliasMap->GetValue( buf );
359     }
360     
361     // fill the objects
362     AliEMCALTriggerTRUDCSConfig* truConfig = trigConfig->GetTRUDCSConfig(iTRU);
363     if( ! truConfig ){
364       AliWarning( Form("EMC DCS TRU%02d config not retrieved!\n", iTRU ));
365     }
366
367     // get last entries. fill the TRU object
368     if( ! arrL0ALGSEL ){
369       AliWarning( Form("EMC DCS TRU%02d L0ALGSEL alias not found!\n", iTRU ));
370     }
371     else{
372       if (debug) AliInfo( Form("arrL0ALGSEL has %d entries \n", arrL0ALGSEL->GetEntries()) );
373       if ( arrL0ALGSEL->GetEntries() > 0 ) {
374         dcsVal = (AliDCSValue *) arrL0ALGSEL->At( arrL0ALGSEL->GetEntries() - 1 );
375         if (dcsVal) truConfig->SetL0SEL( dcsVal->GetUInt() );
376       }
377     }
378
379     if( ! arrPEAKFINDER ){
380       AliWarning( Form("EMC DCS TRU%02d PEAKFINDER alias not found!\n", iTRU ));
381     }
382     else{
383       if (debug) AliInfo( Form("arrPEAKFINDER has %d entries \n", arrPEAKFINDER->GetEntries()) );
384       if ( arrPEAKFINDER->GetEntries() > 0 ) {
385         dcsVal = (AliDCSValue *) arrPEAKFINDER->At( arrPEAKFINDER->GetEntries() - 1 );
386         if (dcsVal) truConfig->SetSELPF( dcsVal->GetUInt() );
387       }
388     }
389
390     if( ! arrGLOBALTHRESH ){
391       AliWarning( Form("EMC DCS TRU%02d GLOBALTHRESH alias not found!\n", iTRU ));
392     }
393     else{
394       if (debug) AliInfo( Form("arrGLOBALTHRESH has %d entries \n", arrGLOBALTHRESH->GetEntries()) );
395       if ( arrGLOBALTHRESH->GetEntries() > 0 ) {
396         dcsVal = (AliDCSValue *) arrGLOBALTHRESH->At( arrGLOBALTHRESH->GetEntries() - 1 );
397         if (dcsVal) truConfig->SetGTHRL0( dcsVal->GetUInt() );
398       }
399     }
400
401     if( ! arrCOSMTHRESH ){
402       AliWarning( Form("EMC DCS TRU%02d COSMTHRESH alias not found!\n", iTRU ));
403     }
404     else{
405       if (debug) AliInfo( Form("arrCOSMTHRESH has %d entries \n", arrCOSMTHRESH->GetEntries()) );
406       if ( arrCOSMTHRESH->GetEntries() > 0 ) {
407         dcsVal = (AliDCSValue *) arrCOSMTHRESH->At( arrCOSMTHRESH->GetEntries() - 1 );
408         if (dcsVal) truConfig->SetL0COSM( dcsVal->GetUInt() );
409       }
410     }
411     
412     for( i = 0; i < 6; i++ ){
413       if( ! arrMASK[i] ){
414         AliWarning( Form("EMC DCS TRU%02d MASK%d alias not found!\n", iTRU, i ));
415       }
416       else{
417         if (debug) AliInfo( Form("arrMASK[%d] has %d entries \n", i, arrMASK[i]->GetEntries()) );
418         if ( arrMASK[i]->GetEntries() > 0 ) {
419           dcsVal = (AliDCSValue *) arrMASK[i]->At( arrMASK[i]->GetEntries() - 1 );
420           if (dcsVal) truConfig->SetMaskReg( dcsVal->GetUInt(), i );
421         }
422       }
423     }
424     
425   } // TRUs
426   AliInfo(Form("TRU info retrieved.\n"));
427   // save the objects
428   AliCDBMetaData metaData;
429   metaData.SetBeamPeriod(0);
430   metaData.SetResponsible(kMetaResponsible);
431   metaData.SetComment(kMetaComment); 
432       
433   Bool_t retCode = Store("Calib", "Trigger", trigConfig, &metaData, 0, kFALSE);
434   AliInfo(Form("TRU info stored.\n"));
435   return retCode;
436 }
437
438 //______________________________________________________________________________________________
439 UInt_t AliEMCALPreprocessor::ExtractPedestals(Int_t sourceFXS)
440 {
441   //  Read pedestal file from file exchange server
442   //  Only store if new pedestal info is available
443   //
444   UInt_t result=0;
445
446   AliCaloCalibPedestal *calibPed = new AliCaloCalibPedestal(AliCaloCalibPedestal::kEmCal);
447   calibPed->Init();
448
449   TList* list = GetFileSources(sourceFXS,"pedestals");
450   if (list && list->GetEntries()>0) {
451     
452     //  loop through all files from LDCs
453
454     int changes = 0;
455     UInt_t index = 0;
456     while (list->At(index)!=NULL) {
457       TObjString* fileNameEntry = (TObjString*) list->At(index);
458       if (fileNameEntry!=NULL) {
459         TString fileName = GetFile(sourceFXS, "pedestals",
460                                    fileNameEntry->GetString().Data());
461         TFile *f = TFile::Open(fileName);
462         if (!f) {
463           Log ("Error opening pedestal file.");
464           result = kReturnCodeNoObject;
465           break;
466         }
467         AliCaloCalibPedestal *calPed;
468         f->GetObject("emcCalibPedestal",calPed);
469         if ( !calPed ) {
470           Log ("No pedestal calibration object in file.");
471           result = kReturnCodeNoObject;
472           break;
473         }
474         if ( calPed->GetNEvents()>0 && calPed->GetNChanFills()>0 ) {
475           // add info for the modules available in the present file
476           Bool_t status = calibPed->AddInfo(calPed);
477           if (status) { changes++; }
478         }
479         
480         delete calPed; 
481         f->Close();
482       }
483       index++;
484     }  // while(list)
485     
486     //
487     //  Store updated pedestal entry to OCDB
488     //
489     if (changes>0) {
490       AliCDBMetaData metaData;
491       metaData.SetBeamPeriod(0);
492       metaData.SetResponsible(kMetaResponsible);
493       metaData.SetComment(kMetaComment); 
494       
495       Bool_t storeOK = StoreReferenceData("Calib", "Pedestals", calibPed, &metaData);
496       if ( !storeOK ) result++;
497     }
498   } 
499   else {
500     Log ("Error: no entries in input file list!");
501     result = kReturnCodeNoEntries;
502   }
503   
504   return result;
505 }
506
507 //______________________________________________________________________________________________
508 UInt_t AliEMCALPreprocessor::ExtractSignal(Int_t sourceFXS)
509 { //  Read signal file from file exchange server
510   //  Only store if new signal info is available
511   //
512   UInt_t result=0;
513   AliCaloCalibSignal *calibSig = new AliCaloCalibSignal(AliCaloCalibSignal::kEmCal); 
514   
515   TList* list = GetFileSources(sourceFXS,"signal");
516   if (list && list->GetEntries()>0) {
517
518     //  loop through all files from LDCs
519     
520     int changes = 0;
521     UInt_t index = 0;
522     while (list->At(index)!=NULL) {
523       TObjString* fileNameEntry = (TObjString*) list->At(index);
524       if (fileNameEntry!=NULL) {
525         TString fileName = GetFile(sourceFXS, "signal",
526                                    fileNameEntry->GetString().Data());
527         TFile *f = TFile::Open(fileName);
528         if (!f) {
529           Log ("Error opening signal file.");
530           result = kReturnCodeNoObject;
531           break;
532         }
533         AliCaloCalibSignal *calSig;
534         f->GetObject("emcCalibSignal",calSig);
535         if ( !calSig ) {
536           Log ("No signal calibration object in file.");
537           result = kReturnCodeNoObject;
538           break;
539         }
540         if ( calSig->GetNEvents()>0 ) {
541           // add info for the modules available in the present file
542           Bool_t status = calibSig->AddInfo(calSig);
543           if (status) { changes++; }
544         }
545         
546         delete calSig; 
547         f->Close();
548       }
549       index++;
550     }  // while(list)
551     
552     //
553     //  Store updated signal entry to OCDB
554     //
555     if (changes>0) {
556       AliCDBMetaData metaData;
557       metaData.SetBeamPeriod(0);
558       metaData.SetResponsible(kMetaResponsible);
559       metaData.SetComment(kMetaComment); 
560       
561       Bool_t storeOK = Store("Calib", "LED", calibSig, &metaData, 0, kFALSE);
562       if ( !storeOK ) result++;
563     }
564   } 
565   else {
566     Log ("Error: no entries in input file list!");
567     result = kReturnCodeNoEntries;
568   }
569
570   return result;
571 }
572
573