]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCPreprocessor.cxx
Detector Algorithm for pedestal runs.
[u/mrichter/AliRoot.git] / ZDC / AliZDCPreprocessor.cxx
1 // --- ROOT system
2 #include <TFile.h>
3 #include <TClonesArray.h>
4 #include <TList.h>
5 #include <TObjString.h>
6 #include <TTimeStamp.h>
7
8 #include "AliZDCPreprocessor.h"
9 #include "AliCDBManager.h"
10 #include "AliCDBEntry.h"
11 #include "AliCDBMetaData.h"
12 #include "AliDCSValue.h"
13 #include "AliAlignObj.h"
14 #include "AliAlignObjParams.h"
15 #include "AliLog.h"
16 #include "AliZDCDataDCS.h"
17 #include "AliZDCChMap.h"
18 #include "AliZDCPedestals.h"
19 #include "AliZDCLaserCalib.h"
20 #include "AliZDCEnCalib.h"
21 #include "AliZDCTowerCalib.h"
22 #include "AliZDCMBCalib.h"
23
24 /////////////////////////////////////////////////////////////////////
25 //                                                                 //
26 // Class implementing Shuttle ZDC pre-processor.                   //
27 // It takes data from DCS and DAQ and writes calibration objects   //
28 // in the OCDB and reference values/histos in the ReferenceData.   //
29 //                                                                 //
30 /////////////////////////////////////////////////////////////////////
31
32 // ******************************************************************
33 //    RETURN CODES:
34 // return 0 : everything OK
35 // return 1 : no DCS input data Map
36 // return 2 : error storing DCS data in RefData 
37 // return 3 : error storing alignment object in OCDB
38 // return 4 : error in ZDCMapping.dat file retrieved from DAQ FXS (not existing|empty|corrupted)
39 // return 5 : error storing mapping obj. in OCDB
40 // return 6 : error storing energy calibration obj. in OCDB
41 // return 7 : error storing tower inter-calibration obj. in OCDB
42 // return 8 : error in ZDCEnergyCalib.dat file retrieved from DAQ FXS 
43 // return 9 : error in ZDCTowerCalib.dat file retrieved from DAQ FXS 
44 // return 10: error in ZDCPedestal.dat file retrieved from DAQ FXS 
45 // return 11: error storing pedestal calibration obj. in OCDB
46 // return 12: error in ZDCPedHisto.root file retrieved from DAQ FXS 
47 // return 13: error storing pedestal histos in RefData 
48 // return 14: error in ZDCLaserCalib.dat file retrieved from DAQ FXS 
49 // return 15: error storing laser calibration obj. in OCDB
50 // return 16: error in ZDCLaserHisto.root file retrieved from DAQ FXS 
51 // return 17: error storing laser histos in RefData 
52 // return 18: error in ZDCMBCalib.root file retrieved from DAQ FXS 
53 // return 19: error storing MB calibration obj. in OCDB
54 // ******************************************************************
55
56 ClassImp(AliZDCPreprocessor)
57
58 //______________________________________________________________________________________________
59 AliZDCPreprocessor::AliZDCPreprocessor(AliShuttleInterface* shuttle) :
60   AliPreprocessor("ZDC", shuttle),
61   fData(0)
62 {
63   // constructor
64   // May 2009 - run types updated according to
65   // http://alice-ecs.web.cern.ch/alice-ecs/runtypes_3.36.html
66   AddRunType("STANDALONE_PEDESTAL");
67   AddRunType("STANDALONE_LASER");
68   AddRunType("STANDALONE_COSMIC");
69   AddRunType("CALIBRATION_EMD");
70   AddRunType("CALIBRATION_MB");
71   AddRunType("CALIBRATION_CENTRAL");
72   AddRunType("CALIBRATION_SEMICENTRAL");
73   AddRunType("CALIBRATION_BC");
74   AddRunType("PHYSICS");
75 }
76
77
78 //______________________________________________________________________________________________
79 AliZDCPreprocessor::~AliZDCPreprocessor()
80 {
81   // destructor
82 }
83
84
85 //______________________________________________________________________________________________
86 void AliZDCPreprocessor::Initialize(Int_t run, UInt_t startTime, UInt_t endTime)
87 {
88   // Creates AliZDCDataDCS object
89
90   AliPreprocessor::Initialize(run, startTime, endTime);
91
92   AliDebug(2,Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s \n\tStartTime DCS Query %s \n\tEndTime DCS Query %s", run,
93                 TTimeStamp(startTime).AsString(),
94                 TTimeStamp(endTime).AsString(), ((TTimeStamp)GetStartTimeDCSQuery()).AsString(), ((TTimeStamp)GetEndTimeDCSQuery()).AsString()));
95
96   fRun = run;
97   fStartTime = startTime;
98   fEndTime = endTime;
99
100   fData = new AliZDCDataDCS(fRun, fStartTime, fEndTime, GetStartTimeDCSQuery(), GetEndTimeDCSQuery());
101 }
102
103 //_____________________________________________________________________________
104 Bool_t AliZDCPreprocessor::ProcessDCS(){
105
106   // tells whether DCS should be processed or not
107
108   TString runType = GetRunType();
109   Log(Form("RunType %s",runType.Data()));
110
111   if (runType=="STANDALONE_COSMIC" || runType=="STANDALONE_PEDESTAL"){
112     return kFALSE;
113   }
114
115   return kTRUE;
116 }
117
118 //______________________________________________________________________________________________
119 UInt_t AliZDCPreprocessor::ProcessDCSData(TMap* dcsAliasMap)
120 {
121   
122   // Fills data into a AliZDCDataDCS object
123   if(!dcsAliasMap){
124     Log(" No DCS map found: ZDC exiting from Shuttle");
125     if(fData){
126       delete fData;
127       fData = 0;
128     }
129     return 1;
130   }
131
132   Log(Form("Processing data from DCS"));
133    
134   // The processing of the DCS input data is forwarded to AliZDCDataDCS
135   //dcsAliasMap->Print(""); 
136   Bool_t resDCSProcess = fData->ProcessData(*dcsAliasMap);
137   if(resDCSProcess==kFALSE){
138     Log(" Problems in processing DCS DP");
139     return 1;
140   }  
141   
142   // ------------------------------------------------------
143   // Change introduced 26/9/09 in order NOT to process the
144   // HV DP since some of them are never found in amanda DB
145   // ------------------------------------------------------
146   // Store DCS data as reference
147   AliCDBMetaData metadata;
148   metadata.SetResponsible("Chiara Oppedisano");
149   metadata.SetComment("DCS DP TMap for ZDC");
150   Bool_t resDCSRef = kTRUE;
151   resDCSRef = StoreReferenceData("DCS","Data", dcsAliasMap, &metadata);
152   
153   if(resDCSRef==kFALSE) return 2;
154
155   // --- Writing ZDC table positions into alignment object
156   TClonesArray *array = new TClonesArray("AliAlignObjParams",10);
157   TClonesArray &alobj = *array;
158   AliAlignObjParams a;
159   Double_t dx=0., dz=0., dpsi=0., dtheta=0., dphi=0.;
160   // Vertical table position in mm from DCS
161   Double_t dyZN1 = (Double_t) (fData->GetAlignData(0)/10.);
162   Double_t dyZP1 = (Double_t) (fData->GetAlignData(1)/10.);
163   Double_t dyZN2 = (Double_t) (fData->GetAlignData(2)/10.);
164   Double_t dyZP2 = (Double_t) (fData->GetAlignData(3)/10.);
165   //
166   const char *n1ZDC="ZDC/NeutronZDC_C";  
167   const char *p1ZDC="ZDC/ProtonZDC_C";
168   const char *n2ZDC="ZDC/NeutronZDC_A";
169   const char *p2ZDC="ZDC/ProtonZDC_A";
170   //
171   UShort_t iIndex=0;
172   AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
173   UShort_t volid = AliGeomManager::LayerToVolUID(iLayer,iIndex);
174   //
175   new(alobj[0]) AliAlignObjParams(n1ZDC, volid, dx, dyZN1, dz, dpsi, dtheta, dphi, kTRUE);
176   new(alobj[1]) AliAlignObjParams(p1ZDC, volid, dx, dyZP1, dz, dpsi, dtheta, dphi, kTRUE);
177   new(alobj[2]) AliAlignObjParams(n2ZDC, volid, dx, dyZN2, dz, dpsi, dtheta, dphi, kTRUE);
178   new(alobj[3]) AliAlignObjParams(p2ZDC, volid, dx, dyZP2, dz, dpsi, dtheta, dphi, kTRUE);
179  
180   // save in CDB storage
181   AliCDBMetaData mdDCS;
182   mdDCS.SetResponsible("Chiara Oppedisano");
183   mdDCS.SetComment("Alignment object for ZDC");
184   Bool_t resultAl = Store("Align","Data", array, &mdDCS, 0, kFALSE);
185   if(resultAl==kFALSE)  return 3;
186   
187   return 0;
188 }
189
190 //______________________________________________________________________________________________
191 UInt_t AliZDCPreprocessor::ProcessChMap()
192
193   const int kNModules=10, kNch=48, kNScch=32;
194   
195   // Reading the file for mapping from FXS
196   TList* daqSource = GetFileSources(kDAQ, "MAPPING");
197   if(!daqSource){
198     AliError(Form("No sources for file ZDCChMapping.dat in run %d ", fRun));
199     return 4;
200   }
201   if(daqSource->GetEntries()==0) return 4;
202   Log("\t List of DAQ sources for MAPPING id: "); daqSource->Print();
203   //
204   TIter iter(daqSource);
205   TObjString* source = 0;
206   Int_t isou = 0;
207   Int_t modMap[kNModules][3], adcMap[kNch][6], scMap[kNScch][6]; 
208   //
209   while((source = dynamic_cast<TObjString*> (iter.Next()))){
210      TString fileName = GetFile(kDAQ, "MAPPING", source->GetName());
211      Log(Form("\t Getting file #%d: ZDCChMapping.dat from %s\n",++isou, source->GetName()));
212
213      if(fileName.Length() <= 0){
214        Log(Form("No file from source %s!", source->GetName()));
215        return 4;
216      }
217      // --- Reading file with calibration data
218      //const char* fname = fileName.Data();
219      if(fileName){
220        FILE *file;
221        if((file = fopen(fileName,"r")) == NULL){
222          printf("Cannot open file %s \n",fileName.Data());
223          return 4;
224        }
225        Log(Form("File %s connected to process data for ADC mapping", fileName.Data()));
226        //
227        for(Int_t j=0; j<kNch; j++){       
228            for(Int_t k=0; k<6; k++){
229              int read = fscanf(file,"%d",&adcMap[j][k]);
230              if(read == 0) AliDebug(3," Failing in reading data from mapping file");
231            }
232        }
233        for(Int_t j=kNch; j<kNch+kNScch; j++){     
234            for(Int_t k=0; k<6; k++){
235              int read = fscanf(file,"%d",&scMap[j-kNch][k]);
236              if(read == 0) AliDebug(3," Failing in reading data from mapping file");
237            }
238        }
239        for(Int_t j=kNch+kNScch; j<kNModules+kNch+kNScch; j++){    
240            for(Int_t k=0; k<3; k++){
241              int read = fscanf(file,"%d",&modMap[j-kNch-kNScch][k]);
242              if(read == 0) AliDebug(3," Failing in reading data from mapping file");
243            }
244        }
245        fclose(file);
246      }
247      else{
248        Log(Form("File %s not found", fileName.Data()));
249        return 4;
250      }
251   }
252   
253   // Store the currently read map ONLY IF it is different
254   // from the entry in the OCDB
255   Bool_t adcMapUpdated=kFALSE, scMapUpdated=kFALSE;
256   Bool_t updateOCDB = kFALSE;
257   
258   AliCDBEntry *cdbEntry = GetFromOCDB("Calib","ChMap");
259   if(!cdbEntry){
260     Log(" No existing CDB entry for ADC mapping");
261     updateOCDB = kTRUE;
262   }
263   else{
264     AliZDCChMap *chMap = (AliZDCChMap*) cdbEntry->GetObject();
265     for(Int_t i=0; i<kNch; i++){
266       if(  (adcMap[i][1] == chMap->GetADCModule(i)) 
267         && (adcMap[i][2] == chMap->GetADCChannel(i)) 
268         && (adcMap[i][4] == chMap->GetDetector(i)) 
269         && (adcMap[i][5] == chMap->GetSector(i))){
270          adcMapUpdated = kFALSE;
271       }
272       else adcMapUpdated = kTRUE;
273     }
274     for(Int_t i=0; i<kNScch; i++){
275       if(  (scMap[i][2] == chMap->GetScChannel(i)) 
276         && (scMap[i][4] == chMap->GetScDetector(i)) 
277         && (scMap[i][5] == chMap->GetScSector(i))){
278          scMapUpdated = kFALSE;
279       }
280       else scMapUpdated = kTRUE;
281     }
282   }
283   if(adcMapUpdated || scMapUpdated) updateOCDB = kTRUE;
284   //
285   Bool_t resChMapStore = kTRUE;
286   if(updateOCDB==kTRUE){
287     Log(" A new entry ZDC/Calib/ChMap will be created");
288     //
289     // --- Initializing mapping calibration object
290     AliZDCChMap *mapCalib = new AliZDCChMap("ZDC");
291     // Writing channel map in the OCDB
292     for(Int_t k=0; k<kNModules; k++){
293       mapCalib->SetModuleMap(k, modMap[k][0], modMap[k][1], modMap[k][2]);
294     }
295     for(Int_t k=0; k<kNch; k++){
296       mapCalib->SetADCModule(k,adcMap[k][1]);
297       mapCalib->SetADCChannel(k,adcMap[k][2]);
298       mapCalib->SetDetector(k,adcMap[k][4]);
299       mapCalib->SetSector(k,adcMap[k][5]);
300     }
301     for(Int_t k=0; k<kNScch; k++){
302        mapCalib->SetScChannel(k, scMap[k][2]);
303        mapCalib->SetScDetector(k, scMap[k][4]);
304        mapCalib->SetScSector(k, scMap[k][5]);
305     }
306     //
307     mapCalib->Print("");
308     // 
309     AliCDBMetaData metaData;
310     metaData.SetBeamPeriod(0);
311     metaData.SetResponsible("Chiara Oppedisano");
312     metaData.SetComment("Filling AliZDCChMap object");  
313     //
314     resChMapStore = Store("Calib","ChMap",mapCalib, &metaData, 0, kTRUE);
315     printf("  Mapping object stored in OCDB\n");
316   }
317   else{
318     Log(" ZDC/Calib/ChMap entry in OCDB is valid and won't be updated");
319     resChMapStore = kTRUE;
320   }
321   
322   delete daqSource; daqSource=0;
323   
324   if(resChMapStore==kFALSE) return 5;
325   else return 0;
326
327 }
328
329 //______________________________________________________________________________________________
330 UInt_t AliZDCPreprocessor::ProcessppData()
331 {
332    Bool_t resEnCal=kTRUE, resTowCal=kTRUE;
333   
334    // *********** Energy calibration
335    // --- Cheking if there is already the entry in the OCDB
336    AliCDBEntry *cdbEnEntry = GetFromOCDB("Calib", "EnergyCalib");
337    if(!cdbEnEntry){   
338      Log(Form(" ZDC/Calib/EnergyCalib entry will be created"));
339      // --- Initializing calibration object
340      AliCDBMetaData metaData;
341      metaData.SetBeamPeriod(0);
342      metaData.SetResponsible("Chiara Oppedisano");
343      //
344      AliZDCEnCalib *eCalib = new AliZDCEnCalib("ZDC");
345      for(Int_t j=0; j<6; j++) eCalib->SetEnCalib(j,1.);
346      metaData.SetComment("AliZDCEnCalib object");  
347      //eCalib->Print("");
348      resEnCal = Store("Calib", "EnergyCalib", eCalib, &metaData, 0, kTRUE);
349    }
350    else{ 
351      // if entry exists it is still valid (=1 for all runs!)
352      Log(Form(" Valid ZDC/Calib/EnergyCalib object already existing in OCDB!!!"));
353      resEnCal = kTRUE;
354    }
355    
356    if(resEnCal==kFALSE) return 6;
357
358    //
359    // *********** Tower inter-calibration
360    // --- Cheking if there is already the entry in the OCDB
361    AliCDBEntry *cdbTowEntry = GetFromOCDB("Calib", "TowerCalib");
362    if(!cdbTowEntry){   
363      AliZDCTowerCalib *towCalib = new AliZDCTowerCalib("ZDC");
364      for(Int_t j=0; j<5; j++){  
365         towCalib->SetZN1EqualCoeff(j, 1.);
366         towCalib->SetZP1EqualCoeff(j, 1.);
367         towCalib->SetZN2EqualCoeff(j, 1.);
368         towCalib->SetZP2EqualCoeff(j, 1.);  
369      }
370      //towCalib->Print("");
371      // 
372      AliCDBMetaData metaData;
373      metaData.SetBeamPeriod(0);
374      metaData.SetResponsible("Chiara Oppedisano");
375      metaData.SetComment("AliZDCTowerCalib object");  
376      //
377      resTowCal = Store("Calib", "TowerCalib", towCalib, &metaData, 0, kTRUE);
378    }
379    else{ 
380      // if entry exists it is still valid (=1 for all runs!)
381      Log(Form(" Valid ZDC/Calib/TowerCalib object already existing in OCDB!!!"));
382      resTowCal = kTRUE;
383    }
384    
385    if(resTowCal==kFALSE) return 7;
386    
387    return 0;
388 }
389
390 //______________________________________________________________________________________________
391 UInt_t AliZDCPreprocessor::ProcessCalibData()
392 {
393   TList* daqSources = GetFileSources(kDAQ, "EMDENERGYCALIB");
394   if(!daqSources){
395     AliError(Form("No sources for CALIBRATION_EMD run %d !", fRun));
396     return 8;
397   }
398   Log("\t List of DAQ sources for EMDENERGYCALIB id: "); daqSources->Print();
399   //
400   TIter iter2(daqSources);
401   TObjString* source = 0;
402   Int_t i=0;
403   Bool_t resEnCal=kTRUE, resTowCal=kTRUE;
404   
405   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
406     TString stringEMDFileName = GetFile(kDAQ, "EMDENERGYCALIB", source->GetName());
407     if(stringEMDFileName.Length() <= 0){
408       Log(Form("No file from source %s!", source->GetName()));
409       return 8;
410     }
411     const char* emdFileName = stringEMDFileName.Data();
412     Log(Form("\t Getting file #%d: %s from %s\n",++i,emdFileName,source->GetName()));
413     //
414     // --- Initializing energy calibration object
415     AliZDCEnCalib *eCalib = new AliZDCEnCalib("ZDC");
416     // --- Reading file with calibration data
417     if(emdFileName){
418       FILE *file;
419       if((file = fopen(emdFileName,"r")) == NULL){
420         printf("Cannot open file %s \n",emdFileName);
421         return 8;
422       }
423       Log(Form("File %s connected to process data from EM dissociation events", emdFileName));
424       //
425       Float_t fitValEMD[6];
426       for(Int_t j=0; j<6; j++){        
427         if(j<6){
428           int iread = fscanf(file,"%f",&fitValEMD[j]);
429           if(iread==0) AliDebug(3," Failing reading daa from EMD calibration data file");
430           eCalib->SetEnCalib(j,fitValEMD[j]);
431         }
432       }
433       //
434       fclose(file);
435     }
436     else{
437       Log(Form("File %s not found", emdFileName));
438       return 8;
439     }
440     //eCalib->Print("");
441     // 
442     AliCDBMetaData metaData;
443     metaData.SetBeamPeriod(0);
444     metaData.SetResponsible("Chiara Oppedisano");
445     metaData.SetComment("Filling AliZDCEnCalib object");  
446     //
447     resEnCal = Store("Calib","EnergyCalib",eCalib, &metaData, 0, kTRUE);
448     if(resEnCal==kFALSE) return 6;
449   }
450   delete daqSources; daqSources = 0;
451   
452   TList* daqSourcesH = GetFileSources(kDAQ, "EMDTOWERCALIB");
453   if(!daqSourcesH){
454     AliError(Form("No sources for CALIBRATION_EMD run %d !", fRun));
455     return 9;
456   }
457   Log("\t List of DAQ sources for EMDTOWERCALIB id: "); daqSourcesH->Print();
458   //
459   TIter iter2H(daqSourcesH);
460   TObjString* sourceH = 0;
461   Int_t iH=0;
462   while((sourceH = dynamic_cast<TObjString*> (iter2H.Next()))){
463     TString stringtowEMDFileName = GetFile(kDAQ, "EMDTOWERCALIB", sourceH->GetName());
464     if(stringtowEMDFileName.Length() <= 0){
465       Log(Form("No file from source %s!", sourceH->GetName()));
466       return 9;
467     }
468     const char * towEMDFileName = stringtowEMDFileName.Data();
469     Log(Form("\t Getting file #%d: %s from source %s\n",++iH,towEMDFileName,sourceH->GetName()));
470     // --- Initializing energy calibration object
471     AliZDCTowerCalib *towCalib = new AliZDCTowerCalib("ZDC");
472     // --- Reading file with calibration data
473     if(towEMDFileName){
474       FILE *file;
475       if((file = fopen(towEMDFileName,"r")) == NULL){
476         printf("Cannot open file %s \n",towEMDFileName);
477         return 9;
478       }
479       //
480       Float_t equalCoeff[4][5];
481       for(Int_t j=0; j<4; j++){        
482          for(Int_t k=0; k<5; k++){
483            int leggi = fscanf(file,"%f",&equalCoeff[j][k]);
484            if(leggi==0) AliDebug(3," Failing reading data from EMD calibration file");
485            if(j==0)      towCalib->SetZN1EqualCoeff(k, equalCoeff[j][k]);
486            else if(j==1) towCalib->SetZP1EqualCoeff(k, equalCoeff[j][k]);
487            else if(j==2) towCalib->SetZN2EqualCoeff(k, equalCoeff[j][k]);
488            else if(j==3) towCalib->SetZP2EqualCoeff(k, equalCoeff[j][k]);  
489          }
490       }
491       //
492       fclose(file);
493     }
494     else{
495       Log(Form("File %s not found", towEMDFileName));
496       return 9;
497     }
498     //towCalib->Print("");
499     // 
500     AliCDBMetaData metaData;
501     metaData.SetBeamPeriod(0);
502     metaData.SetResponsible("Chiara Oppedisano");
503     metaData.SetComment("Filling AliZDCTowerCalib object");  
504     //
505     resTowCal = Store("Calib","TowerCalib",towCalib, &metaData, 0, kTRUE);
506     if(resTowCal==kFALSE) return 7;
507   }
508   delete daqSourcesH; daqSourcesH = 0;
509   
510      
511   return 0;
512 }
513
514 //______________________________________________________________________________________________
515 UInt_t AliZDCPreprocessor::ProcessPedestalData()
516 {
517   TList* daqSources = GetFileSources(kDAQ, "PEDESTALDATA");
518   if(!daqSources){
519     Log(Form("No source for STANDALONE_PEDESTAL run %d !", fRun));
520     return 10;
521   }
522   if(daqSources->GetEntries()==0) return 10;
523   Log("\t List of DAQ sources for PEDESTALDATA id: "); daqSources->Print();
524   //
525   TIter iter(daqSources);
526   TObjString* source;
527   Int_t i=0;
528   Bool_t resPedCal=kTRUE, resPedHist=kTRUE;
529   
530   while((source = dynamic_cast<TObjString*> (iter.Next()))){
531      TString stringPedFileName = GetFile(kDAQ, "PEDESTALDATA", source->GetName());
532      if(stringPedFileName.Length() <= 0){
533         Log(Form("No PEDESTALDATA file from source %s!", source->GetName()));
534         return 10;
535      }
536      const char* pedFileName = stringPedFileName.Data();
537      Log(Form("\t Getting file #%d: %s from %s\n",++i,pedFileName,source->GetName()));
538      //
539      // --- Initializing pedestal calibration object
540      AliZDCPedestals *pedCalib = new AliZDCPedestals("ZDC");
541      // --- Reading file with pedestal calibration data
542      // no. ADCch = (22 signal ch. + 2 reference PMs) * 2 gain chain = 48
543      const Int_t knZDCch = 48;
544      FILE *file;
545      if((file = fopen(pedFileName,"r")) == NULL){
546        printf("Cannot open file %s \n",pedFileName);
547        return 10;
548      }
549      Log(Form("File %s connected to process pedestal data", pedFileName));
550      Float_t pedVal[(2*knZDCch)][2];
551      for(Int_t k=0; k<(2*knZDCch); k++){
552         for(Int_t j=0; j<2; j++){
553            int aleggi = fscanf(file,"%f",&pedVal[k][j]);
554            if(aleggi==0) AliDebug(3," Failing reading data from pedestal file");
555            //if(j==1) printf("pedVal[%d] -> %f, %f \n",k,pedVal[k][0],pedVal[k][1]);
556         }
557         if(k<knZDCch){
558           pedCalib->SetMeanPed(k,pedVal[k][0]);
559           pedCalib->SetMeanPedWidth(k,pedVal[k][1]);
560         }
561         else if(k>=knZDCch && k<(2*knZDCch)){
562           pedCalib->SetOOTPed(k-knZDCch,pedVal[k][0]);
563           pedCalib->SetOOTPedWidth(k-knZDCch,pedVal[k][1]);
564         }
565         else if(k>=(2*knZDCch) && k<(3*knZDCch)){
566           pedCalib->SetPedCorrCoeff(k-(2*knZDCch),pedVal[k][0],pedVal[k][1]);
567         }
568      }
569      fclose(file);
570      //pedCalib->Print("");
571      // 
572      AliCDBMetaData metaData;
573      metaData.SetBeamPeriod(0);
574      metaData.SetResponsible("Chiara Oppedisano");
575      metaData.SetComment("Filling AliZDCPedestals object");  
576      //
577      resPedCal = Store("Calib","Pedestals",pedCalib, &metaData, 0, kTRUE);
578      if(resPedCal==kFALSE) return 11;
579   }
580   delete daqSources; daqSources = 0;
581   
582   TList* daqSourceH = GetFileSources(kDAQ, "PEDESTALHISTOS");
583   if(!daqSourceH){
584     Log(Form("No source for PEDESTALHISTOS id run %d !", fRun));
585     return 12;
586   }
587   Log("\t List of DAQ sources for PEDESTALHISTOS id: "); daqSourceH->Print();
588   //
589   TIter iterH(daqSourceH);
590   TObjString* sourceH = 0;
591   Int_t iH=0;
592   while((sourceH = dynamic_cast<TObjString*> (iterH.Next()))){
593      TString stringPedFileName = GetFile(kDAQ, "PEDESTALHISTOS", sourceH->GetName());
594      if(stringPedFileName.Length() <= 0){
595         Log(Form("No PEDESTALHISTOS file from source %s!", sourceH->GetName()));
596         return 12;
597      }
598      const char* pedFileName = stringPedFileName.Data();
599      Log(Form("\t Getting file #%d: %s from %s\n",++iH, pedFileName, sourceH->GetName()));
600      resPedHist = StoreReferenceFile(pedFileName, "pedestalReference.root");
601      if(resPedHist==kFALSE) return 13;
602   }
603   delete daqSourceH; daqSourceH=0;
604   
605   return 0;
606 }
607
608 //______________________________________________________________________________________________
609 UInt_t AliZDCPreprocessor::ProcessLaserData()
610 {
611   TList* daqSources = GetFileSources(kDAQ, "LASERDATA");
612   if(!daqSources){
613     AliError(Form("No sources for STANDALONE_LASER run %d !", fRun));
614     return 14;
615   }
616   if(daqSources->GetEntries()==0) return 14;
617   Log("\t List of DAQ sources for LASERDATA id: "); daqSources->Print();
618   //
619   TIter iter2(daqSources);
620   TObjString* source = 0;
621   Int_t i=0;
622   Bool_t resLaserCal=kTRUE, resLaserHist=kTRUE;
623   
624   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
625      TString stringLaserFileName = GetFile(kDAQ, "LASERDATA", source->GetName());
626      if(stringLaserFileName.Length() <= 0){
627        Log(Form("No LASER file from source %s!", source->GetName()));
628        return 14;
629      }
630      const char* laserFileName = stringLaserFileName.Data();
631      Log(Form("\t Getting file #%d: %s from %s\n",++i,laserFileName,source->GetName()));
632      //
633      // --- Initializing pedestal calibration object
634      AliZDCLaserCalib *lCalib = new AliZDCLaserCalib("ZDC");
635      // --- Reading file with pedestal calibration data
636      if(laserFileName){
637        FILE *file;
638        if((file = fopen(laserFileName,"r")) == NULL){
639          printf("Cannot open file %s \n",laserFileName);
640          return 14;
641        }
642        Log(Form("File %s connected to process data from LASER events", laserFileName));
643        //
644        Float_t ivalRead[22][4]; 
645        for(Int_t j=0; j<22; j++){
646           for(Int_t k=0; k<4; k++){
647             int aleggi = fscanf(file,"%f",&ivalRead[j][k]);
648             if(aleggi==0) AliDebug(3," Failng reading data from laser file");
649             //printf(" %d %1.0f  ",k, ivalRead[j][k]);
650           }
651           lCalib->SetDetector(j, (Int_t) ivalRead[j][0]);
652           lCalib->SetSector(j, (Int_t) ivalRead[j][1]);
653           lCalib->SetfPMValue(j, ivalRead[j][2]);
654           lCalib->SetfPMWidth(j, ivalRead[j][3]);
655        }
656        fclose(file);
657      }
658      else{
659        Log(Form("File %s not found", laserFileName));
660        return 14;
661      }
662      //lCalib->Print("");
663      // 
664      AliCDBMetaData metaData;
665      metaData.SetBeamPeriod(0);
666      metaData.SetResponsible("Chiara Oppedisano");
667      metaData.SetComment("Filling AliZDCLaserCalib object");  
668      //
669      resLaserCal = Store("Calib","LaserCalib",lCalib, &metaData, 0, kTRUE);
670      if(resLaserCal==kFALSE) return 15;
671   }
672   delete daqSources; daqSources = 0;
673
674   TList* daqSourceH = GetFileSources(kDAQ, "LASERHISTOS");
675   if(!daqSourceH){
676     AliError(Form("No sources for STANDALONE_LASER run %d !", fRun));
677     return 16;
678   }
679   Log("\t List of DAQ sources for LASERHISTOS id: "); daqSourceH->Print();
680   //
681   TIter iter2H(daqSourceH);
682   TObjString* sourceH = 0;
683   Int_t iH=0;
684   while((sourceH = dynamic_cast<TObjString*> (iter2H.Next()))){
685      Log(Form("\t Getting file #%d\n",++iH));
686      TString stringLaserFileName = GetFile(kDAQ, "LASERHISTOS", sourceH->GetName());
687      if(stringLaserFileName.Length() <= 0){
688        Log(Form("No LASER file from source %s!", sourceH->GetName()));
689        return 16;
690      }
691      resLaserHist = StoreReferenceFile(stringLaserFileName.Data(), "laserReference.root");
692      //
693      if(resLaserHist==kFALSE) return 17;
694   }
695   delete daqSourceH; daqSourceH = 0;
696   
697   return 0;
698 }
699
700
701 //______________________________________________________________________________________________
702 UInt_t AliZDCPreprocessor::ProcessMBCalibData()
703 {
704   TList* daqSources = GetFileSources(kDAQ, "MBCALIB");
705   if(!daqSources){
706     AliError(Form("No sources for CALIBRATION_MB run %d !", fRun));
707     return 18;
708   }
709   if(daqSources->GetEntries()==0) return 18;
710   Log("\t List of DAQ sources for MBCALIB id: "); daqSources->Print();
711   //
712   TIter iter2(daqSources);
713   TObjString* source = 0;
714   Int_t i=0;
715   Bool_t resMBCal=kTRUE;
716   
717   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
718      TString stringMBFileName = GetFile(kDAQ, "MBCALIB", source->GetName());
719      if(stringMBFileName.Length() <= 0){
720        Log(Form("No MBCALIB file from source %s!", source->GetName()));
721        return 18;
722      }
723      const char* mbFileName = stringMBFileName.Data();
724      Log(Form("\t Getting file #%d: %s from %s\n",++i,mbFileName,source->GetName()));
725      //
726      // --- Initializing calibration object
727      AliZDCMBCalib *mbCalib = new AliZDCMBCalib("ZDC");
728      // --- Reading file with calibration data
729      if(mbFileName){
730        TFile * fileHistos = TFile::Open(mbFileName);
731        Log(Form("File %s connected to process data from CALIBRATION_MB events", mbFileName));
732        //
733        fileHistos->cd();
734        TH2F *hZDCvsZEM = (TH2F*)  fileHistos->Get("hZDCvsZEM");
735        TH2F *hZDCCvsZEM = (TH2F*) fileHistos->Get("hZDCCvsZEM");
736        TH2F *hZDCAvsZEM = (TH2F*) fileHistos->Get("hZDCAvsZEM");
737        //
738        mbCalib->SetZDCvsZEM(hZDCvsZEM);
739        mbCalib->SetZDCCvsZEM(hZDCCvsZEM);
740        mbCalib->SetZDCAvsZEM(hZDCAvsZEM);
741        //
742        //fileHistos->Close();
743      }
744      else{
745        Log(Form("File %s not found", mbFileName));
746        return 14;
747      }
748      // 
749      AliCDBMetaData metaData;
750      metaData.SetBeamPeriod(0);
751      metaData.SetResponsible("Chiara Oppedisano");
752      metaData.SetComment("Filling AliZDCMBCalib object");  
753      //
754      //mbCalib->Dump();
755      //
756      resMBCal = Store("Calib","MBCalib",mbCalib, &metaData, 0, kTRUE);
757        printf(" here 1000\n");
758      if(resMBCal==kFALSE) return 19;
759   }
760   delete daqSources; daqSources = 0;
761   
762   return 0;
763 }
764
765 //______________________________________________________________________________________________
766 UInt_t AliZDCPreprocessor::Process(TMap* dcsAliasMap)
767 {
768  UInt_t resDCS = 0;
769  UInt_t resChMap=0;
770  UInt_t resEnergyCalib=0, resPedestalCalib=0, resLaserCalib=0, resMBCalib=0;
771
772  // ************************* Process DCS data ****************************
773  if(ProcessDCS()) resDCS = ProcessDCSData(dcsAliasMap);
774   
775  // ********************************* From DAQ ************************************
776
777  const char* beamType = GetRunParameter("beamType");
778  TString runType = GetRunType();
779  printf("\t **** AliZDCPreprocessor -> beamType %s, runType %s ****\n",beamType,runType.Data());
780
781  // ******************************************
782  // ADC channel mapping
783  // ******************************************
784  resChMap = ProcessChMap();
785  
786  // ******************************************
787  // Calibration param. for p-p data (all = 1)
788  // ******************************************
789  // NO ENERGY CALIBRATION -> coefficients set to 1.
790  // Temp -> also inter-calibration coefficients are set to 1.
791  if((strcmp(beamType,"p-p")==0) || (strcmp(beamType,"P-P")==0)) resEnergyCalib = ProcessppData();
792  
793  // *****************************************************
794  // CALIBRATION_EMD -> Energy calibration and equalization
795  // *****************************************************
796  if((strcmp(beamType,"A-A")==0) && (runType.CompareTo("CALIBRATION_EMD")==0)) 
797    resEnergyCalib =  ProcessCalibData();
798  
799  // *****************************************************
800  // STANDALONE_PEDESTALS -> Pedestal subtraction 
801  // *****************************************************
802  if(runType.CompareTo("STANDALONE_PEDESTAL")==0) resPedestalCalib = ProcessPedestalData();
803  
804  // *****************************************************
805  // STANDALONE_LASER -> Signal stability and ageing 
806  // *****************************************************
807  else if(runType.CompareTo("STANDALONE_LASER")==0) resLaserCalib = ProcessLaserData();
808
809  // *****************************************************
810  // CALIBRATION_MB -> Signal stability and ageing 
811  // *****************************************************
812  else if(runType.CompareTo("CALIBRATION_MB")==0) resMBCalib = ProcessMBCalibData();
813  
814
815  if(resDCS!=0)                return resDCS;
816  else if(resChMap!=0)         return resChMap;
817  else if(resEnergyCalib!=0)   return resEnergyCalib;
818  else if(resPedestalCalib!=0) return resPedestalCalib;
819  else if(resLaserCalib!=0)    return resLaserCalib;
820  else if(resMBCalib!=0)       return resMBCalib;
821  
822  return 0;
823   
824 }