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