]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPreprocessor.cxx
Writing back ITS alignement data (Ruben).
[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         
330         
331         AliInfo("Print DCS alias map content");
332         dcsAliasMap->Print();
333         
334   AliInfo(Form("Get TRU info from DCS DPs.\n"));
335   Int_t i, iTRU;
336   const Int_t bufsize = 1000;
337   char buf[bufsize];
338
339   AliDCSValue *dcsVal;
340   TObjArray *arrL0ALGSEL, *arrPEAKFINDER, *arrGLOBALTHRESH, *arrCOSMTHRESH;
341   TObjArray *arrMASK[6];
342         
343   TObjArray *arrSTUG[3][2], *arrSTUJ[3][2];
344   TObjArray *arrSTUD, *arrSTUR, *arrSTUF;
345         
346   // overall object to hold STU and DCS config info
347   // DS comment: for now only holds TRU info, i.e. only partially filled
348   // (STU info only in raw data header; unfortunately not also picked up via DCS DPs)
349   AliEMCALTriggerDCSConfig *trigConfig = new AliEMCALTriggerDCSConfig();
350   // allocate space for TRU objects
351   TClonesArray *truArr = new TClonesArray("AliEMCALTriggerTRUDCSConfig", kNTRU);
352   for( iTRU = 0; iTRU < kNTRU; iTRU++){
353     new((*truArr)[iTRU]) AliEMCALTriggerTRUDCSConfig();
354   }
355   trigConfig->SetTRUArr(truArr);
356
357   AliEMCALTriggerSTUDCSConfig *stuConfig = new AliEMCALTriggerSTUDCSConfig();
358         
359   // loop through all TRUs
360   bool debug = true; // debug flag for AliInfo printouts for each TRU
361   for( iTRU = 0; iTRU < kNTRU; iTRU++){
362     if (debug) AliInfo( Form("iTRU %d \n", iTRU) );
363     // get the shuttled values
364     snprintf( buf, bufsize, "EMC_TRU%02d_L0ALGSEL", iTRU );
365     arrL0ALGSEL = (TObjArray*) dcsAliasMap->GetValue( buf );
366     snprintf( buf, bufsize, "EMC_TRU%02d_PEAKFINDER", iTRU );
367     arrPEAKFINDER = (TObjArray*) dcsAliasMap->GetValue( buf );
368     snprintf( buf, bufsize, "EMC_TRU%02d_GLOBALTHRESH", iTRU );
369     arrGLOBALTHRESH = (TObjArray*) dcsAliasMap->GetValue( buf );
370     snprintf( buf, bufsize, "EMC_TRU%02d_COSMTHRESH", iTRU );
371     arrCOSMTHRESH = (TObjArray*) dcsAliasMap->GetValue( buf );
372     
373     for( i = 0; i < 6; i++ ){
374       snprintf( buf, bufsize, "EMC_TRU%02d_MASK%d", iTRU, i );
375       arrMASK[i] = (TObjArray*) dcsAliasMap->GetValue( buf );
376     }
377     
378     // fill the objects
379     AliEMCALTriggerTRUDCSConfig* truConfig = trigConfig->GetTRUDCSConfig(iTRU);
380     if( ! truConfig ){
381       AliWarning( Form("EMC TRU%02d config not retrieved!\n", iTRU ));
382       continue;
383     }
384
385     // get last entries. fill the TRU object
386     if( ! arrL0ALGSEL ){
387       AliWarning( Form("EMC_TRU%02d_L0ALGSEL alias not found!\n", iTRU ));
388     }
389     else{
390       if (debug) AliInfo( Form("arrL0ALGSEL has %d entries \n", arrL0ALGSEL->GetEntries()) );
391       if ( arrL0ALGSEL->GetEntries() > 0 ) {
392         dcsVal = (AliDCSValue *) arrL0ALGSEL->At( arrL0ALGSEL->GetEntries() - 1 );
393         if (dcsVal) {
394                 truConfig->SetL0SEL( dcsVal->GetUInt() );
395     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
396         }
397       }
398       else
399       AliWarning( Form("EMC_TRU%02d_L0ALGSEL has no entries!\n", iTRU ));
400     }
401
402     if( ! arrPEAKFINDER ){
403       AliWarning( Form("EMC_TRU%02d_PEAKFINDER alias not found!\n", iTRU ));
404     }
405     else{
406       if (debug) AliInfo( Form("arrPEAKFINDER has %d entries \n", arrPEAKFINDER->GetEntries()) );
407       if ( arrPEAKFINDER->GetEntries() > 0 ) {
408         dcsVal = (AliDCSValue *) arrPEAKFINDER->At( arrPEAKFINDER->GetEntries() - 1 );
409         if (dcsVal){
410                 truConfig->SetSELPF( dcsVal->GetUInt() );
411     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
412         }
413       }
414       else
415       AliWarning( Form("EMC_TRU%02d_PEAKFINDER has no entries!\n", iTRU ));
416     }
417
418     if( ! arrGLOBALTHRESH ){
419       AliWarning( Form("EMC_TRU%02d_GLOBALTHRESH alias not found!\n", iTRU ));
420     }
421     else{
422       if (debug) AliInfo( Form("arrGLOBALTHRESH has %d entries \n", arrGLOBALTHRESH->GetEntries()) );
423       if ( arrGLOBALTHRESH->GetEntries() > 0 ) {
424         dcsVal = (AliDCSValue *) arrGLOBALTHRESH->At( arrGLOBALTHRESH->GetEntries() - 1 );
425         if (dcsVal){
426                 truConfig->SetGTHRL0( dcsVal->GetUInt() );
427     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
428         }
429       }
430       else
431       AliWarning( Form("EMC_TRU%02d_GLOBALTHRESH has no entries!\n", iTRU ));
432     }
433
434     if( ! arrCOSMTHRESH ){
435       AliWarning( Form("EMC_TRU%02d_COSMTHRESH alias not found!\n", iTRU ));
436     }
437     else{
438       if (debug) AliInfo( Form("arrCOSMTHRESH has %d entries \n", arrCOSMTHRESH->GetEntries()) );
439       if ( arrCOSMTHRESH->GetEntries() > 0 ) {
440         dcsVal = (AliDCSValue *) arrCOSMTHRESH->At( arrCOSMTHRESH->GetEntries() - 1 );
441         if (dcsVal){
442                 truConfig->SetL0COSM( dcsVal->GetUInt() );
443     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
444         }
445       }
446       else
447       AliWarning( Form("EMC_TRU%02d_COSMTHRESH has no entries!\n", iTRU ));
448     }
449     
450     for( i = 0; i < 6; i++ ){
451       if( ! arrMASK[i] ){
452         AliWarning( Form("EMC_TRU%02d_MASK%d alias not found!\n", iTRU, i ));
453       }
454       else{
455         if (debug) AliInfo( Form("arrMASK[%d] has %d entries \n", i, arrMASK[i]->GetEntries()) );
456         if ( arrMASK[i]->GetEntries() > 0 ) {
457           dcsVal = (AliDCSValue *) arrMASK[i]->At( arrMASK[i]->GetEntries() - 1 );
458           if (dcsVal){
459                         truConfig->SetMaskReg( dcsVal->GetUInt(), i );
460     if (debug) AliInfo( Form("saving value: %u\n", dcsVal->GetUInt()) );
461         }
462         }
463       else
464       AliWarning( Form("EMC_TRU%02d_MASK%d has no entries!\n", iTRU, i ));
465       }
466     }
467     
468   } // TRUs
469   AliInfo(Form("TRU info retrieved.\n"));
470                 
471   // STU
472   for (i = 0; i < 3; i++) {
473                 for (int j = 0; j < 2; j++) {
474                         arrSTUG[i][j] = (TObjArray*)dcsAliasMap->GetValue(Form("EMC_STU_G%c%d", i + 65, j));
475                         arrSTUJ[i][j] = (TObjArray*)dcsAliasMap->GetValue(Form("EMC_STU_J%c%d", i + 65, j));    
476                         
477                         if (!arrSTUG[i][j]) {
478                                 AliWarning(Form("EMC_STU_G%c%d alias not found!", i + 65, j));
479                         } else {
480                                 if (debug) AliInfo( Form("EMC_STU_G%c%d has %d entries", i + 65, j, arrSTUG[i][j]->GetEntries()));
481                                 if (arrSTUG[i][j]->GetEntries() > 0) {
482                                         dcsVal = (AliDCSValue*)arrSTUG[i][j]->At(arrSTUG[i][j]->GetEntries() - 1);
483                                         if (dcsVal) {
484                                                 stuConfig->SetG(i, j, dcsVal->GetInt());
485                                                 if (debug) AliInfo(Form("saving value: %u", dcsVal->GetInt()));
486                                         }
487                                 } else
488                                         AliWarning(Form("EMC_STU_G%c%d has no entry!", i + 65, j));
489                         }
490                         
491                         if (!arrSTUJ[i][j]) {
492                                 AliWarning(Form("EMC_STU_J%c%d alias not found!", i + 65, j));
493                         } else {
494                                 if (debug) AliInfo( Form("EMC_STU_J%c%d has %d entries", i + 65, j, arrSTUJ[i][j]->GetEntries()));
495                                 if (arrSTUJ[i][j]->GetEntries() > 0) {
496                                         dcsVal = (AliDCSValue*)arrSTUJ[i][j]->At(arrSTUJ[i][j]->GetEntries() - 1);
497                                         if (dcsVal) {
498                                                 stuConfig->SetJ(i, j, dcsVal->GetInt());
499                                                 if (debug) AliInfo(Form("saving value: %u", dcsVal->GetInt()));
500                                         }
501                                 } else
502                                         AliWarning(Form("EMC_STU_J%c%d has no entry!", i + 65, j));
503                         }
504                 }
505         }
506         
507         arrSTUD = (TObjArray*)dcsAliasMap->GetValue("EMC_STU_GETRAW");
508         arrSTUR = (TObjArray*)dcsAliasMap->GetValue("EMC_STU_REGION");  
509         arrSTUF = (TObjArray*)dcsAliasMap->GetValue("EMC_STU_FWVERS");  
510
511         if (!arrSTUD) {
512                 AliWarning("EMC_STU_GETRAW alias not found!");
513         } else {
514                 if (debug) AliInfo(Form("EMC_STU_GETRAW has %d entries", arrSTUD->GetEntries()));
515                 if (arrSTUD->GetEntries() > 0) {
516                         dcsVal = (AliDCSValue*)arrSTUD->At(arrSTUD->GetEntries() - 1);
517                         if (dcsVal) {
518                                 stuConfig->SetRawData(dcsVal->GetInt());
519                                 if (debug) AliInfo(Form("saving value: %u", dcsVal->GetInt()));
520                         }
521                 } else
522                         AliWarning("EMC_STU_GETRAW has no entry!");
523         }
524         
525         if (!arrSTUR) {
526                 AliWarning("EMC_STU_REGION");
527         } else {
528                 if (debug) AliInfo( Form("EMC_STU_REGION has %d entries", arrSTUR->GetEntries()));
529                 if (arrSTUR->GetEntries() > 0) {
530                         dcsVal = (AliDCSValue*)arrSTUR->At(arrSTUR->GetEntries() - 1);
531                         if (dcsVal) {
532                                 stuConfig->SetRegion(dcsVal->GetInt());
533                                 if (debug) AliInfo(Form("saving value: %u", dcsVal->GetInt()));
534                         }
535                 } else
536                         AliWarning("EMC_STU_REGION has no entry!");
537         }
538         
539         if (!arrSTUF) {
540                 AliWarning("EMC_STU_FWVERS");
541         } else {
542                 if (debug) AliInfo(Form("EMC_STU_FWVERS has %d entries", arrSTUF->GetEntries()));
543                 if (arrSTUF->GetEntries() > 0) {
544                         dcsVal = (AliDCSValue*)arrSTUF->At(arrSTUF->GetEntries() - 1);
545                         if (dcsVal) {
546                                 stuConfig->SetFw(dcsVal->GetInt());
547                                 if (debug) AliInfo(Form("saving value: %u", dcsVal->GetInt()));
548                         }
549                 } else
550                         AliWarning("EMC_STU_FWVERS has no entry!");
551         }
552         
553         trigConfig->SetSTUObj(stuConfig);
554         
555         AliInfo(Form("STU info retrieved."));
556
557         
558   // save the objects
559   AliCDBMetaData metaData;
560   metaData.SetBeamPeriod(0);
561   metaData.SetResponsible(kMetaResponsible);
562   metaData.SetComment(kMetaComment); 
563       
564   UInt_t result=0;
565   Bool_t storeOK = Store("Calib", "Trigger", trigConfig, &metaData, 0, kFALSE);
566   if ( !storeOK )  result=1;
567   AliInfo(Form("TRU info stored. result %d\n", result));
568
569   return result;
570 }
571
572 //______________________________________________________________________________________________
573 UInt_t AliEMCALPreprocessor::ExtractPedestals(Int_t sourceFXS)
574 {
575   //  Read pedestal file from file exchange server
576   //  Only store if new pedestal info is available
577   //
578   UInt_t result=0;
579
580   AliCaloCalibPedestal *calibPed = new AliCaloCalibPedestal(AliCaloCalibPedestal::kEmCal);
581   calibPed->Init();
582
583   TList* list = GetFileSources(sourceFXS,"pedestals");
584   if (list && list->GetEntries()>0) {
585     
586     //  loop through all files from LDCs
587
588     int changes = 0;
589     UInt_t index = 0;
590     while (list->At(index)!=NULL) {
591       TObjString* fileNameEntry = (TObjString*) list->At(index);
592       if (fileNameEntry!=NULL) {
593         TString fileName = GetFile(sourceFXS, "pedestals",
594                                    fileNameEntry->GetString().Data());
595         TFile *f = TFile::Open(fileName);
596         if (!f) {
597           Log ("Error opening pedestal file.");
598           result = kReturnCodeNoObject;
599           break;
600         }
601         AliCaloCalibPedestal *calPed;
602         f->GetObject("emcCalibPedestal",calPed);
603         if ( !calPed ) {
604           Log ("No pedestal calibration object in file.");
605           result = kReturnCodeNoObject;
606           break;
607         }
608         if ( calPed->GetNEvents()>0 && calPed->GetNChanFills()>0 ) {
609           // add info for the modules available in the present file
610           Bool_t status = calibPed->AddInfo(calPed);
611           if (status) { changes++; }
612         }
613         
614         delete calPed; 
615         f->Close();
616       }
617       index++;
618     }  // while(list)
619     
620     //
621     //  Store updated pedestal entry to OCDB
622     //
623     if (changes>0) {
624       AliCDBMetaData metaData;
625       metaData.SetBeamPeriod(0);
626       metaData.SetResponsible(kMetaResponsible);
627       metaData.SetComment(kMetaComment); 
628       
629       Bool_t storeOK = StoreReferenceData("Calib", "Pedestals", calibPed, &metaData);
630       if ( !storeOK ) result++;
631     }
632   } 
633   else {
634     Log ("Error: no entries in input file list!");
635     result = kReturnCodeNoEntries;
636   }
637   
638   return result;
639 }
640
641 //______________________________________________________________________________________________
642 UInt_t AliEMCALPreprocessor::ExtractSignal(Int_t sourceFXS)
643 { //  Read signal file from file exchange server
644   //  Only store if new signal info is available
645   //
646   UInt_t result=0;
647   AliCaloCalibSignal *calibSig = new AliCaloCalibSignal(AliCaloCalibSignal::kEmCal); 
648   
649   TList* list = GetFileSources(sourceFXS,"signal");
650   if (list && list->GetEntries()>0) {
651
652     //  loop through all files from LDCs
653     
654     int changes = 0;
655     UInt_t index = 0;
656     while (list->At(index)!=NULL) {
657       TObjString* fileNameEntry = (TObjString*) list->At(index);
658       if (fileNameEntry!=NULL) {
659         TString fileName = GetFile(sourceFXS, "signal",
660                                    fileNameEntry->GetString().Data());
661         TFile *f = TFile::Open(fileName);
662         if (!f) {
663           Log ("Error opening signal file.");
664           result = kReturnCodeNoObject;
665           break;
666         }
667         AliCaloCalibSignal *calSig;
668         f->GetObject("emcCalibSignal",calSig);
669         if ( !calSig ) {
670           Log ("No signal calibration object in file.");
671           result = kReturnCodeNoObject;
672           break;
673         }
674         if ( calSig->GetNEvents()>0 ) {
675           // add info for the modules available in the present file
676           Bool_t status = calibSig->AddInfo(calSig);
677           if (status) { changes++; }
678         }
679         
680         delete calSig; 
681         f->Close();
682       }
683       index++;
684     }  // while(list)
685     
686     //
687     //  Store updated signal entry to OCDB
688     //
689     if (changes>0) {
690       AliCDBMetaData metaData;
691       metaData.SetBeamPeriod(0);
692       metaData.SetResponsible(kMetaResponsible);
693       metaData.SetComment(kMetaComment); 
694       
695       Bool_t storeOK = Store("Calib", "LED", calibSig, &metaData, 0, kFALSE);
696       if ( !storeOK ) result++;
697     }
698   } 
699   else {
700     Log ("Error: no entries in input file list!");
701     result = kReturnCodeNoEntries;
702   }
703
704   return result;
705 }
706
707