]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPreprocessor.cxx
Test for Coverity
[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       UInt_t triggerResult = MapTriggerConfig(dcsAliasMap);
202       result+=triggerResult;
203       status = new TParameter<int>("triggerResult",triggerResult);
204       resultArray->Add(status);
205     }
206   }
207   
208   // Other calibration information will be retrieved through FXS files
209   //  examples:
210   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
211   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
212   //
213   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
214   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
215   
216   // PEDESTAL ENTRIES:
217   
218   if ( runType == kPedestalRunType ) {
219     Int_t numSources = 1;
220     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ, AliShuttleInterface::kHLT} ;
221     TString source = fConfEnv->GetValue("Pedestal","DAQ");
222     source.ToUpper();
223     if (source != "OFF" ) { 
224       if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
225       if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
226       if (source == "HLTDAQ" ) {
227         numSources=2;
228         pedestalSource[0] = AliShuttleInterface::kHLT;
229         pedestalSource[1] = AliShuttleInterface::kDAQ;
230       }
231       if (source == "DAQHLT" ) numSources=2;
232       UInt_t pedestalResult=0;
233       for (Int_t i=0; i<numSources; i++ ) {     
234         pedestalResult = ExtractPedestals(pedestalSource[i]);
235         if ( pedestalResult == 0 ) break;
236       }
237       result += pedestalResult;
238       status = new TParameter<int>("pedestalResult",pedestalResult);
239       resultArray->Add(status);
240     }
241   }
242   
243   // SIGNAL/LED ENTRIES:
244   if( runType == kPhysicsRunType ) {
245     Int_t numSources = 1;
246     Int_t signalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
247     TString source = fConfEnv->GetValue("Signal","DAQ");
248     source.ToUpper();
249     if ( source != "OFF") { 
250       if ( source == "HLT") signalSource[0] = AliShuttleInterface::kHLT;
251       if (!GetHLTStatus()) signalSource[0] = AliShuttleInterface::kDAQ;
252       if (source == "HLTDAQ" ) {
253         numSources=2;
254         signalSource[0] = AliShuttleInterface::kHLT;
255         signalSource[1] = AliShuttleInterface::kDAQ;
256       }
257       if (source == "DAQHLT" ) numSources=2;
258       UInt_t signalResult=0;
259       for (Int_t i=0; i<numSources; i++ ) {     
260         signalResult = ExtractSignal(signalSource[i]);
261         if ( signalResult == 0 ) break;
262       }
263       result += signalResult;
264       status = new TParameter<int>("signalResult",signalResult);
265       resultArray->Add(status);
266     }
267   }
268   
269   
270   // overall status at the end
271   if (errorHandling == "OFF" ) {
272     AliCDBMetaData metaData;
273     metaData.SetBeamPeriod(0);
274     metaData.SetResponsible(kMetaResponsible);
275     metaData.SetComment("Preprocessor AliEMCAL status.");
276     Bool_t storeOK = Store("Calib", "PreprocStatus", resultArray, &metaData, 0, kFALSE);
277     resultArray->Delete();
278     result = 0;
279     if ( !storeOK )  result=1;
280     return result;
281   } 
282   else { 
283     return result;
284   }
285   
286 }
287 //______________________________________________________________________________________________
288 UInt_t AliEMCALPreprocessor::MapTemperature(TMap* dcsAliasMap)
289 { // extract DCS temperature maps. Perform fits to save space
290   UInt_t result=0;
291
292   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
293   if (map) {
294     fTemp->MakeSplineFit(map);
295     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
296     if (fitFraction > kFitFraction ) {
297       AliInfo(Form("Temperature values extracted, fits performed.\n"));
298     } 
299     else { 
300       Log ("Too few temperature maps fitted. \n");
301       result = kReturnCodeNoInfo;
302     }
303   } 
304   else {
305     Log("No temperature map extracted. \n");
306     result = kReturnCodeNoInfo;
307   }
308   delete map;
309   // Now store the final CDB file
310   
311   if ( result == 0 ) { // some info was found
312     AliCDBMetaData metaData;
313     metaData.SetBeamPeriod(0);
314     metaData.SetResponsible(kMetaResponsible);
315     metaData.SetComment(kMetaComment);
316     
317     Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
318     if ( !storeOK )  result=1;
319     AliInfo(Form("Temperature info stored. result %d\n", result));
320   }
321   
322   return result;
323 }
324
325 //______________________________________________________________________________________________
326 UInt_t AliEMCALPreprocessor::MapTriggerConfig(TMap* dcsAliasMap)
327 { // extract DCS trigger info
328   AliInfo(Form("Get TRU info from DCS DPs.\n"));
329   Int_t i, iTRU;
330   const Int_t bufsize = 1000;
331   char buf[bufsize];
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   // allocate space for TRU objects
342   TClonesArray *truArr = new TClonesArray("AliEMCALTriggerTRUDCSConfig", kNTRU);
343   for( iTRU = 0; iTRU < kNTRU; iTRU++){
344     new((*truArr)[iTRU]) AliEMCALTriggerTRUDCSConfig();
345   }
346   trigConfig->SetTRUArr(truArr);
347
348   // loop through all TRUs
349   bool debug = true; // debug flag for AliInfo printouts for each TRU
350   for( iTRU = 0; iTRU < kNTRU; iTRU++){
351     if (debug) AliInfo( Form("iTRU %d \n", iTRU) );
352     // get the shuttled values
353     snprintf( buf, bufsize, "EMC_TRU%02d_L0ALGSEL", iTRU );
354     arrL0ALGSEL = (TObjArray*) dcsAliasMap->GetValue( buf );
355     snprintf( buf, bufsize, "EMC_TRU%02d_PEAKFINDER", iTRU );
356     arrPEAKFINDER = (TObjArray*) dcsAliasMap->GetValue( buf );
357     snprintf( buf, bufsize, "EMC_TRU%02d_GLOBALTHRESH", iTRU );
358     arrGLOBALTHRESH = (TObjArray*) dcsAliasMap->GetValue( buf );
359     snprintf( buf, bufsize, "EMC_TRU%02d_COSMTHRESH", iTRU );
360     arrCOSMTHRESH = (TObjArray*) dcsAliasMap->GetValue( buf );
361     
362     for( i = 0; i < 6; i++ ){
363       snprintf( buf, bufsize, "EMC_TRU%02d_MASK%d", iTRU, i );
364       arrMASK[i] = (TObjArray*) dcsAliasMap->GetValue( buf );
365     }
366     
367     // fill the objects
368     AliEMCALTriggerTRUDCSConfig* truConfig = trigConfig->GetTRUDCSConfig(iTRU);
369     if( ! truConfig ){
370       AliWarning( Form("EMC TRU%02d config not retrieved!\n", iTRU ));
371       continue;
372     }
373
374     // get last entries. fill the TRU object
375     if( ! arrL0ALGSEL ){
376       AliWarning( Form("EMC_TRU%02d_L0ALGSEL alias not found!\n", iTRU ));
377     }
378     else{
379       if (debug) AliInfo( Form("arrL0ALGSEL has %d entries \n", arrL0ALGSEL->GetEntries()) );
380       if ( arrL0ALGSEL->GetEntries() > 0 ) {
381         dcsVal = (AliDCSValue *) arrL0ALGSEL->At( arrL0ALGSEL->GetEntries() - 1 );
382         if (dcsVal) {
383                 truConfig->SetL0SEL( dcsVal->GetUInt() );
384     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
385         }
386       }
387       else
388       AliWarning( Form("EMC_TRU%02d_L0ALGSEL has no entries!\n", iTRU ));
389     }
390
391     if( ! arrPEAKFINDER ){
392       AliWarning( Form("EMC_TRU%02d_PEAKFINDER alias not found!\n", iTRU ));
393     }
394     else{
395       if (debug) AliInfo( Form("arrPEAKFINDER has %d entries \n", arrPEAKFINDER->GetEntries()) );
396       if ( arrPEAKFINDER->GetEntries() > 0 ) {
397         dcsVal = (AliDCSValue *) arrPEAKFINDER->At( arrPEAKFINDER->GetEntries() - 1 );
398         if (dcsVal){
399                 truConfig->SetSELPF( dcsVal->GetUInt() );
400     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
401         }
402       }
403       else
404       AliWarning( Form("EMC_TRU%02d_PEAKFINDER has no entries!\n", iTRU ));
405     }
406
407     if( ! arrGLOBALTHRESH ){
408       AliWarning( Form("EMC_TRU%02d_GLOBALTHRESH alias not found!\n", iTRU ));
409     }
410     else{
411       if (debug) AliInfo( Form("arrGLOBALTHRESH has %d entries \n", arrGLOBALTHRESH->GetEntries()) );
412       if ( arrGLOBALTHRESH->GetEntries() > 0 ) {
413         dcsVal = (AliDCSValue *) arrGLOBALTHRESH->At( arrGLOBALTHRESH->GetEntries() - 1 );
414         if (dcsVal){
415                 truConfig->SetGTHRL0( dcsVal->GetUInt() );
416     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
417         }
418       }
419       else
420       AliWarning( Form("EMC_TRU%02d_GLOBALTHRESH has no entries!\n", iTRU ));
421     }
422
423     if( ! arrCOSMTHRESH ){
424       AliWarning( Form("EMC_TRU%02d_COSMTHRESH alias not found!\n", iTRU ));
425     }
426     else{
427       if (debug) AliInfo( Form("arrCOSMTHRESH has %d entries \n", arrCOSMTHRESH->GetEntries()) );
428       if ( arrCOSMTHRESH->GetEntries() > 0 ) {
429         dcsVal = (AliDCSValue *) arrCOSMTHRESH->At( arrCOSMTHRESH->GetEntries() - 1 );
430         if (dcsVal){
431                 truConfig->SetL0COSM( dcsVal->GetUInt() );
432     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
433         }
434       }
435       else
436       AliWarning( Form("EMC_TRU%02d_COSMTHRESH has no entries!\n", iTRU ));
437     }
438     
439     for( i = 0; i < 6; i++ ){
440       if( ! arrMASK[i] ){
441         AliWarning( Form("EMC_TRU%02d_MASK%d alias not found!\n", iTRU, i ));
442       }
443       else{
444         if (debug) AliInfo( Form("arrMASK[%d] has %d entries \n", i, arrMASK[i]->GetEntries()) );
445         if ( arrMASK[i]->GetEntries() > 0 ) {
446           dcsVal = (AliDCSValue *) arrMASK[i]->At( arrMASK[i]->GetEntries() - 1 );
447           if (dcsVal){
448                         truConfig->SetMaskReg( dcsVal->GetUInt(), i );
449     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
450         }
451         }
452       else
453       AliWarning( Form("EMC_TRU%02d_MASK%d has no entries!\n", iTRU, i ));
454       }
455     }
456     
457   } // TRUs
458   AliInfo(Form("TRU info retrieved.\n"));
459   // save the objects
460   AliCDBMetaData metaData;
461   metaData.SetBeamPeriod(0);
462   metaData.SetResponsible(kMetaResponsible);
463   metaData.SetComment(kMetaComment); 
464       
465   UInt_t result=0;
466   Bool_t storeOK = Store("Calib", "Trigger", trigConfig, &metaData, 0, kFALSE);
467   if ( !storeOK )  result=1;
468   AliInfo(Form("TRU info stored. result %d\n", result));
469
470   return result;
471 }
472
473 //______________________________________________________________________________________________
474 UInt_t AliEMCALPreprocessor::ExtractPedestals(Int_t sourceFXS)
475 {
476   //  Read pedestal file from file exchange server
477   //  Only store if new pedestal info is available
478   //
479   UInt_t result=0;
480
481   AliCaloCalibPedestal *calibPed = new AliCaloCalibPedestal(AliCaloCalibPedestal::kEmCal);
482   calibPed->Init();
483
484   TList* list = GetFileSources(sourceFXS,"pedestals");
485   if (list && list->GetEntries()>0) {
486     
487     //  loop through all files from LDCs
488
489     int changes = 0;
490     UInt_t index = 0;
491     while (list->At(index)!=NULL) {
492       TObjString* fileNameEntry = (TObjString*) list->At(index);
493       if (fileNameEntry!=NULL) {
494         TString fileName = GetFile(sourceFXS, "pedestals",
495                                    fileNameEntry->GetString().Data());
496         TFile *f = TFile::Open(fileName);
497         if (!f) {
498           Log ("Error opening pedestal file.");
499           result = kReturnCodeNoObject;
500           break;
501         }
502         AliCaloCalibPedestal *calPed;
503         f->GetObject("emcCalibPedestal",calPed);
504         if ( !calPed ) {
505           Log ("No pedestal calibration object in file.");
506           result = kReturnCodeNoObject;
507           break;
508         }
509         if ( calPed->GetNEvents()>0 && calPed->GetNChanFills()>0 ) {
510           // add info for the modules available in the present file
511           Bool_t status = calibPed->AddInfo(calPed);
512           if (status) { changes++; }
513         }
514         
515         delete calPed; 
516         f->Close();
517       }
518       index++;
519     }  // while(list)
520     
521     //
522     //  Store updated pedestal entry to OCDB
523     //
524     if (changes>0) {
525       AliCDBMetaData metaData;
526       metaData.SetBeamPeriod(0);
527       metaData.SetResponsible(kMetaResponsible);
528       metaData.SetComment(kMetaComment); 
529       
530       Bool_t storeOK = StoreReferenceData("Calib", "Pedestals", calibPed, &metaData);
531       if ( !storeOK ) result++;
532     }
533   } 
534   else {
535     Log ("Error: no entries in input file list!");
536     result = kReturnCodeNoEntries;
537   }
538   
539   return result;
540 }
541
542 //______________________________________________________________________________________________
543 UInt_t AliEMCALPreprocessor::ExtractSignal(Int_t sourceFXS)
544 { //  Read signal file from file exchange server
545   //  Only store if new signal info is available
546   //
547   UInt_t result=0;
548   AliCaloCalibSignal *calibSig = new AliCaloCalibSignal(AliCaloCalibSignal::kEmCal); 
549   
550   TList* list = GetFileSources(sourceFXS,"signal");
551   if (list && list->GetEntries()>0) {
552
553     //  loop through all files from LDCs
554     
555     int changes = 0;
556     UInt_t index = 0;
557     while (list->At(index)!=NULL) {
558       TObjString* fileNameEntry = (TObjString*) list->At(index);
559       if (fileNameEntry!=NULL) {
560         TString fileName = GetFile(sourceFXS, "signal",
561                                    fileNameEntry->GetString().Data());
562         TFile *f = TFile::Open(fileName);
563         if (!f) {
564           Log ("Error opening signal file.");
565           result = kReturnCodeNoObject;
566           break;
567         }
568         AliCaloCalibSignal *calSig;
569         f->GetObject("emcCalibSignal",calSig);
570         if ( !calSig ) {
571           Log ("No signal calibration object in file.");
572           result = kReturnCodeNoObject;
573           break;
574         }
575         if ( calSig->GetNEvents()>0 ) {
576           // add info for the modules available in the present file
577           Bool_t status = calibSig->AddInfo(calSig);
578           if (status) { changes++; }
579         }
580         
581         delete calSig; 
582         f->Close();
583       }
584       index++;
585     }  // while(list)
586     
587     //
588     //  Store updated signal entry to OCDB
589     //
590     if (changes>0) {
591       AliCDBMetaData metaData;
592       metaData.SetBeamPeriod(0);
593       metaData.SetResponsible(kMetaResponsible);
594       metaData.SetComment(kMetaComment); 
595       
596       Bool_t storeOK = Store("Calib", "LED", calibSig, &metaData, 0, kFALSE);
597       if ( !storeOK ) result++;
598     }
599   } 
600   else {
601     Log ("Error: no entries in input file list!");
602     result = kReturnCodeNoEntries;
603   }
604
605   return result;
606 }
607
608