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