]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCPreprocessor.cxx
TRD nSigma OADB related new codes and modifications and OADB root file -- Xianguo Lu
[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   delete daqSource; daqSource=0;
340   
341   TString runType = GetRunType();
342   if(runType.CompareTo("PHYSICS")==0){
343     Log(Form("RunType %s -> producing TDC calibration data",runType.Data()));
344     
345     // Reading the file for mapping from FXS
346     TList* daqSourcetdc = GetFileSources(kDAQ, "TDCDATA");
347     if(!daqSourcetdc){
348       AliError(Form("No sources for file ZDCChMappingTDCCalib.dat in run %d ", fRun));
349       return 20;
350     }
351     if(daqSourcetdc->GetEntries()==0) return 20;
352     Log("\t List of DAQ sources for TDCDATA id: "); daqSourcetdc->Print();
353     //
354     Bool_t resTDCcal = kTRUE;
355     TIter itertdc(daqSourcetdc);
356     TObjString* sourcetdc = 0;
357     Int_t isoutdc = 0;
358     //
359     while((sourcetdc = dynamic_cast<TObjString*> (itertdc.Next()))){
360      TString fileNametdc = GetFile(kDAQ, "TDCDATA", sourcetdc->GetName());
361      Log(Form("\t Getting file #%d: ZDCTDCdata.dat from %s\n",++isoutdc, sourcetdc->GetName()));
362
363      if(fileNametdc.Length() <= 0){
364        Log(Form("No file from source %s!", sourcetdc->GetName()));
365        return 20;
366      }
367      // --- Initializing TDC calibration object
368      AliZDCTDCCalib *tdcCalib = new AliZDCTDCCalib("ZDC");
369      // --- Reading file with calibration data
370      //const char* fname = fileName.Data();
371      if(fileNametdc){
372        FILE *filetdc;
373        if((filetdc = fopen(fileNametdc,"r")) == NULL){
374          printf("Cannot open file %s \n",fileNametdc.Data());
375          return 20;
376        }
377        Log(Form("File %s connected to process TDC data", fileNametdc.Data()));
378        //
379        Float_t tdcMean[6][2];
380        for(Int_t it=0; it<6; it++){
381          for(Int_t iu=0; iu<2; iu++) tdcMean[it][iu]=0.;
382        }
383        for(Int_t k=0; k<6; k++){
384         for(Int_t j=0; j<2; j++){
385            int leggi = fscanf(filetdc,"%f",&tdcMean[k][j]);
386            if(leggi==0) AliDebug(3," Failing reading data from tdc file");
387            tdcCalib->SetMeanTDC(k, tdcMean[k][0]);
388            tdcCalib->SetWidthTDC(k, tdcMean[k][1]);
389         }
390        }
391        fclose(filetdc);
392      }
393      else{
394        Log(Form("File %s not found", fileNametdc.Data()));
395        return 20;
396      }
397      //
398      AliCDBMetaData metaData;
399      metaData.SetBeamPeriod(0);
400      metaData.SetResponsible("Chiara Oppedisano");
401      metaData.SetComment("Filling AliZDCTDCCalib object");  
402      //
403      resTDCcal = Store("Calib","TDCCalib",tdcCalib, &metaData, 0, kTRUE);
404      if(resTDCcal==kFALSE) return 21;
405     }
406     delete daqSourcetdc; daqSourcetdc = 0;
407
408     Bool_t restdcHist = kTRUE;
409     TList* daqSourceH = GetFileSources(kDAQ, "TDCHISTOS");
410     if(!daqSourceH){
411       Log(Form("No source for TDCHISTOS id run %d !", fRun));
412       return 22;
413     }
414     Log("\t List of DAQ sources for TDCHISTOS id: "); daqSourceH->Print();
415     //
416     TIter iterH(daqSourceH);
417     TObjString* sourceH = 0;
418     Int_t iH=0;
419     while((sourceH = dynamic_cast<TObjString*> (iterH.Next()))){
420      TString stringTDCFileName = GetFile(kDAQ, "TDCHISTOS", sourceH->GetName());
421      if(stringTDCFileName.Length() <= 0){
422         Log(Form("No TDCHISTOS file from source %s!", sourceH->GetName()));
423         return 22;
424      }
425      const char* tdcFileName = stringTDCFileName.Data();
426      Log(Form("\t Getting file #%d: %s from %s\n",++iH, tdcFileName, sourceH->GetName()));
427      restdcHist = StoreReferenceFile(tdcFileName, "tdcReference.root");
428      if(restdcHist==kFALSE) return 23;
429     }
430     delete daqSourceH; daqSourceH=0;
431   }
432   
433     if(resChMapStore==kFALSE) return 5;
434     else return 0;
435
436 }
437
438 //______________________________________________________________________________________________
439 UInt_t AliZDCPreprocessor::ProcessppData()
440 {
441    Bool_t resEnCal=kTRUE, resTowCal=kTRUE;
442   
443    // *********** Energy calibration
444    // --- Cheking if there is already the entry in the OCDB
445    AliCDBEntry *cdbEnEntry = GetFromOCDB("Calib", "EnergyCalib");
446    if(!cdbEnEntry){   
447      Log(Form(" ZDC/Calib/EnergyCalib entry will be created"));
448      // --- Initializing calibration object
449      AliCDBMetaData metaData;
450      metaData.SetBeamPeriod(0);
451      metaData.SetResponsible("Chiara Oppedisano");
452      //
453      AliZDCEnCalib *eCalib = new AliZDCEnCalib("ZDC");
454      for(Int_t j=0; j<6; j++) eCalib->SetEnCalib(j,1.);
455      metaData.SetComment("AliZDCEnCalib object");  
456      //eCalib->Print("");
457      resEnCal = Store("Calib", "EnergyCalib", eCalib, &metaData, 0, kTRUE);
458    }
459    else{ 
460      // if entry exists it is still valid (=1 for all runs!)
461      Log(Form(" Valid ZDC/Calib/EnergyCalib object already existing in OCDB!!!"));
462      resEnCal = kTRUE;
463    }
464    
465    if(resEnCal==kFALSE) return 6;
466
467    //
468    // *********** Tower inter-calibration
469    // --- Cheking if there is already the entry in the OCDB
470    AliCDBEntry *cdbTowEntry = GetFromOCDB("Calib", "TowerCalib");
471    if(!cdbTowEntry){   
472      AliZDCTowerCalib *towCalib = new AliZDCTowerCalib("ZDC");
473      for(Int_t j=0; j<5; j++){  
474         towCalib->SetZN1EqualCoeff(j, 1.);
475         towCalib->SetZP1EqualCoeff(j, 1.);
476         towCalib->SetZN2EqualCoeff(j, 1.);
477         towCalib->SetZP2EqualCoeff(j, 1.);  
478      }
479      //towCalib->Print("");
480      // 
481      AliCDBMetaData metaData;
482      metaData.SetBeamPeriod(0);
483      metaData.SetResponsible("Chiara Oppedisano");
484      metaData.SetComment("AliZDCTowerCalib object");  
485      //
486      resTowCal = Store("Calib", "TowerCalib", towCalib, &metaData, 0, kTRUE);
487    }
488    else{ 
489      // if entry exists it is still valid (=1 for all runs!)
490      Log(Form(" Valid ZDC/Calib/TowerCalib object already existing in OCDB!!!"));
491      resTowCal = kTRUE;
492    }
493    
494    if(resTowCal==kFALSE) return 7;
495    
496    return 0;
497 }
498
499 //______________________________________________________________________________________________
500 UInt_t AliZDCPreprocessor::ProcessCalibData(Float_t beamEnergy)
501 {
502   TList* daqSources = GetFileSources(kDAQ, "EMDENERGYCALIB");
503   if(!daqSources){
504     AliError(Form("No sources for CALIBRATION_EMD run %d !", fRun));
505     return 8;
506   }
507   Log("\t List of DAQ sources for EMDENERGYCALIB id: "); daqSources->Print();
508   //
509   TIter iter2(daqSources);
510   TObjString* source = 0;
511   Int_t i=0;
512   Bool_t resEnCal=kTRUE, resTowCal=kTRUE;
513   
514   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
515     TString stringEMDFileName = GetFile(kDAQ, "EMDENERGYCALIB", source->GetName());
516     if(stringEMDFileName.Length() <= 0){
517       Log(Form("No file from source %s!", source->GetName()));
518       return 8;
519     }
520     const char* emdFileName = stringEMDFileName.Data();
521     Log(Form("\t Getting file #%d: %s from %s\n",++i,emdFileName,source->GetName()));
522     //
523     // --- Initializing energy calibration object
524     AliZDCEnCalib *eCalib = new AliZDCEnCalib("ZDC");
525     // --- Reading file with calibration data
526     if(emdFileName){
527       FILE *file;
528       if((file = fopen(emdFileName,"r")) == NULL){
529         printf("Cannot open file %s \n",emdFileName);
530         return 8;
531       }
532       Log(Form("File %s connected to process data from EM dissociation events", emdFileName));
533       //
534       Float_t fitValEMD[6];
535       for(Int_t j=0; j<6; j++){        
536         if(j<6){
537           int iread = fscanf(file,"%f",&fitValEMD[j]);
538           if(iread==0) AliDebug(3," Failing reading data from EMD calibration data file");
539           if(fitValEMD[j]!=1.) eCalib->SetEnCalib(j, beamEnergy/fitValEMD[j]);
540           else eCalib->SetEnCalib(j, fitValEMD[j]);
541         }
542       }
543       //
544       fclose(file);
545     }
546     else{
547       Log(Form("File %s not found", emdFileName));
548       return 8;
549     }
550     //eCalib->Print("");
551     // 
552     AliCDBMetaData metaData;
553     metaData.SetBeamPeriod(0);
554     metaData.SetResponsible("Chiara Oppedisano");
555     metaData.SetComment("Filling AliZDCEnCalib object");  
556     //
557     resEnCal = Store("Calib","EnergyCalib",eCalib, &metaData, 0, kTRUE);
558     if(resEnCal==kFALSE) return 6;
559   }
560   delete daqSources; daqSources = 0;
561   
562   TList* daqSourcesH = GetFileSources(kDAQ, "EMDTOWERCALIB");
563   if(!daqSourcesH){
564     AliError(Form("No sources for CALIBRATION_EMD run %d !", fRun));
565     return 9;
566   }
567   Log("\t List of DAQ sources for EMDTOWERCALIB id: "); daqSourcesH->Print();
568   //
569   TIter iter2H(daqSourcesH);
570   TObjString* sourceH = 0;
571   Int_t iH=0;
572   while((sourceH = dynamic_cast<TObjString*> (iter2H.Next()))){
573     TString stringtowEMDFileName = GetFile(kDAQ, "EMDTOWERCALIB", sourceH->GetName());
574     if(stringtowEMDFileName.Length() <= 0){
575       Log(Form("No file from source %s!", sourceH->GetName()));
576       return 9;
577     }
578     const char * towEMDFileName = stringtowEMDFileName.Data();
579     Log(Form("\t Getting file #%d: %s from source %s\n",++iH,towEMDFileName,sourceH->GetName()));
580     // --- Initializing energy calibration object
581     AliZDCTowerCalib *towCalib = new AliZDCTowerCalib("ZDC");
582     // --- Reading file with calibration data
583     if(towEMDFileName){
584       FILE *file;
585       if((file = fopen(towEMDFileName,"r")) == NULL){
586         printf("Cannot open file %s \n",towEMDFileName);
587         return 9;
588       }
589       //
590       Float_t equalCoeff[4][5];
591       for(Int_t j=0; j<4; j++){        
592          for(Int_t k=0; k<5; k++){
593            int leggi = fscanf(file,"%f",&equalCoeff[j][k]);
594            if(leggi==0) AliDebug(3," Failing reading data from EMD calibration file");
595            if(j==0)      towCalib->SetZN1EqualCoeff(k, equalCoeff[j][k]);
596            else if(j==1) towCalib->SetZP1EqualCoeff(k, equalCoeff[j][k]);
597            else if(j==2) towCalib->SetZN2EqualCoeff(k, equalCoeff[j][k]);
598            else if(j==3) towCalib->SetZP2EqualCoeff(k, equalCoeff[j][k]);  
599          }
600       }
601       //
602       fclose(file);
603     }
604     else{
605       Log(Form("File %s not found", towEMDFileName));
606       return 9;
607     }
608     //towCalib->Print("");
609     // 
610     AliCDBMetaData metaData;
611     metaData.SetBeamPeriod(0);
612     metaData.SetResponsible("Chiara Oppedisano");
613     metaData.SetComment("Filling AliZDCTowerCalib object");  
614     //
615     resTowCal = Store("Calib","TowerCalib",towCalib, &metaData, 0, kTRUE);
616     if(resTowCal==kFALSE) return 7;
617   }
618   delete daqSourcesH; daqSourcesH = 0;
619   
620      
621   return 0;
622 }
623
624 //______________________________________________________________________________________________
625 UInt_t AliZDCPreprocessor::ProcessPedestalData()
626 {
627   TList* daqSources = GetFileSources(kDAQ, "PEDESTALDATA");
628   if(!daqSources){
629     Log(Form("No source for STANDALONE_PEDESTAL run %d !", fRun));
630     return 10;
631   }
632   if(daqSources->GetEntries()==0) return 10;
633   Log("\t List of DAQ sources for PEDESTALDATA id: "); daqSources->Print();
634   //
635   TIter iter(daqSources);
636   TObjString* source;
637   Int_t i=0;
638   Bool_t resPedCal=kTRUE, resPedHist=kTRUE;
639   
640   while((source = dynamic_cast<TObjString*> (iter.Next()))){
641      TString stringPedFileName = GetFile(kDAQ, "PEDESTALDATA", source->GetName());
642      if(stringPedFileName.Length() <= 0){
643         Log(Form("No PEDESTALDATA file from source %s!", source->GetName()));
644         return 10;
645      }
646      const char* pedFileName = stringPedFileName.Data();
647      Log(Form("\t Getting file #%d: %s from %s\n",++i,pedFileName,source->GetName()));
648      //
649      // --- Initializing pedestal calibration object
650      AliZDCPedestals *pedCalib = new AliZDCPedestals("ZDC");
651      // --- Reading file with pedestal calibration data
652      // no. ADCch = (22 signal ch. + 2 reference PMs) * 2 gain chain = 48
653      const Int_t knZDCch = 48;
654      FILE *file;
655      if((file = fopen(pedFileName,"r")) == NULL){
656        printf("Cannot open file %s \n",pedFileName);
657        return 10;
658      }
659      Log(Form("File %s connected to process pedestal data", pedFileName));
660      Float_t pedVal[(3*knZDCch)][2];
661      for(Int_t k=0; k<(3*knZDCch); k++){
662         for(Int_t j=0; j<2; j++){
663            int aleggi = fscanf(file,"%f",&pedVal[k][j]);
664            if(aleggi==0) AliDebug(3," Failing reading data from pedestal file");
665            //if(j==1) printf("pedVal[%d] -> %f, %f \n",k,pedVal[k][0],pedVal[k][1]);
666         }
667         if(k<knZDCch){
668           pedCalib->SetMeanPed(k,pedVal[k][0]);
669           pedCalib->SetMeanPedWidth(k,pedVal[k][1]);
670         }
671         else if(k>=knZDCch && k<(2*knZDCch)){
672           pedCalib->SetOOTPed(k-knZDCch,pedVal[k][0]);
673           pedCalib->SetOOTPedWidth(k-knZDCch,pedVal[k][1]);
674         }
675         else if(k>=(2*knZDCch) && k<(3*knZDCch)){
676           pedCalib->SetPedCorrCoeff(k-(2*knZDCch),pedVal[k][0],pedVal[k][1]);
677         }
678      }
679      fclose(file);
680      //pedCalib->Print("");
681      // 
682      AliCDBMetaData metaData;
683      metaData.SetBeamPeriod(0);
684      metaData.SetResponsible("Chiara Oppedisano");
685      metaData.SetComment("Filling AliZDCPedestals object");  
686      //
687      resPedCal = Store("Calib","Pedestals",pedCalib, &metaData, 0, kTRUE);
688      if(resPedCal==kFALSE) return 11;
689   }
690   delete daqSources; daqSources = 0;
691   
692   TList* daqSourceH = GetFileSources(kDAQ, "PEDESTALHISTOS");
693   if(!daqSourceH){
694     Log(Form("No source for PEDESTALHISTOS id run %d !", fRun));
695     return 12;
696   }
697   Log("\t List of DAQ sources for PEDESTALHISTOS id: "); daqSourceH->Print();
698   //
699   TIter iterH(daqSourceH);
700   TObjString* sourceH = 0;
701   Int_t iH=0;
702   while((sourceH = dynamic_cast<TObjString*> (iterH.Next()))){
703      TString stringPedFileName = GetFile(kDAQ, "PEDESTALHISTOS", sourceH->GetName());
704      if(stringPedFileName.Length() <= 0){
705         Log(Form("No PEDESTALHISTOS file from source %s!", sourceH->GetName()));
706         return 12;
707      }
708      const char* pedFileName = stringPedFileName.Data();
709      Log(Form("\t Getting file #%d: %s from %s\n",++iH, pedFileName, sourceH->GetName()));
710      resPedHist = StoreReferenceFile(pedFileName, "pedestalReference.root");
711      if(resPedHist==kFALSE) return 13;
712   }
713   delete daqSourceH; daqSourceH=0;
714   
715   return 0;
716 }
717
718 //______________________________________________________________________________________________
719 UInt_t AliZDCPreprocessor::ProcessLaserData()
720 {
721   TList* daqSources = GetFileSources(kDAQ, "LASERDATA");
722   if(!daqSources){
723     AliError(Form("No sources for STANDALONE_LASER run %d !", fRun));
724     return 14;
725   }
726   if(daqSources->GetEntries()==0) return 14;
727   Log("\t List of DAQ sources for LASERDATA id: "); daqSources->Print();
728   //
729   TIter iter2(daqSources);
730   TObjString* source = 0;
731   Int_t i=0;
732   Bool_t resLaserCal=kTRUE, resLaserHist=kTRUE;
733   
734   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
735      TString stringLaserFileName = GetFile(kDAQ, "LASERDATA", source->GetName());
736      if(stringLaserFileName.Length() <= 0){
737        Log(Form("No LASER file from source %s!", source->GetName()));
738        return 14;
739      }
740      const char* laserFileName = stringLaserFileName.Data();
741      Log(Form("\t Getting file #%d: %s from %s\n",++i,laserFileName,source->GetName()));
742      //
743      // --- Initializing pedestal calibration object
744      AliZDCLaserCalib *lCalib = new AliZDCLaserCalib("ZDC");
745      // --- Reading file with pedestal calibration data
746      if(laserFileName){
747        FILE *file;
748        if((file = fopen(laserFileName,"r")) == NULL){
749          printf("Cannot open file %s \n",laserFileName);
750          return 14;
751        }
752        Log(Form("File %s connected to process data from LASER events", laserFileName));
753        //
754        Float_t ivalRead[22][4]; 
755        for(Int_t j=0; j<22; j++){
756           for(Int_t k=0; k<4; k++){
757             int aleggi = fscanf(file,"%f",&ivalRead[j][k]);
758             if(aleggi==0) AliDebug(3," Failng reading data from laser file");
759             //printf(" %d %1.0f  ",k, ivalRead[j][k]);
760           }
761           lCalib->SetDetector(j, (Int_t) ivalRead[j][0]);
762           lCalib->SetSector(j, (Int_t) ivalRead[j][1]);
763           lCalib->SetfPMValue(j, ivalRead[j][2]);
764           lCalib->SetfPMWidth(j, ivalRead[j][3]);
765        }
766        fclose(file);
767      }
768      else{
769        Log(Form("File %s not found", laserFileName));
770        return 14;
771      }
772      //lCalib->Print("");
773      // 
774      AliCDBMetaData metaData;
775      metaData.SetBeamPeriod(0);
776      metaData.SetResponsible("Chiara Oppedisano");
777      metaData.SetComment("Filling AliZDCLaserCalib object");  
778      //
779      resLaserCal = Store("Calib","LaserCalib",lCalib, &metaData, 0, kTRUE);
780      if(resLaserCal==kFALSE) return 15;
781   }
782   delete daqSources; daqSources = 0;
783
784   TList* daqSourceH = GetFileSources(kDAQ, "LASERHISTOS");
785   if(!daqSourceH){
786     AliError(Form("No sources for STANDALONE_LASER run %d !", fRun));
787     return 16;
788   }
789   Log("\t List of DAQ sources for LASERHISTOS id: "); daqSourceH->Print();
790   //
791   TIter iter2H(daqSourceH);
792   TObjString* sourceH = 0;
793   Int_t iH=0;
794   while((sourceH = dynamic_cast<TObjString*> (iter2H.Next()))){
795      Log(Form("\t Getting file #%d\n",++iH));
796      TString stringLaserFileName = GetFile(kDAQ, "LASERHISTOS", sourceH->GetName());
797      if(stringLaserFileName.Length() <= 0){
798        Log(Form("No LASER file from source %s!", sourceH->GetName()));
799        return 16;
800      }
801      resLaserHist = StoreReferenceFile(stringLaserFileName.Data(), "laserReference.root");
802      //
803      if(resLaserHist==kFALSE) return 17;
804   }
805   delete daqSourceH; daqSourceH = 0;
806   
807   return 0;
808 }
809
810
811 //______________________________________________________________________________________________
812 UInt_t AliZDCPreprocessor::ProcessMBCalibData()
813 {
814   TList* daqSources = GetFileSources(kDAQ, "MBCALIB");
815   if(!daqSources){
816     AliError(Form("No sources for CALIBRATION_MB run %d !", fRun));
817     return 18;
818   }
819   if(daqSources->GetEntries()==0) return 18;
820   Log("\t List of DAQ sources for MBCALIB id: "); daqSources->Print();
821   //
822   TIter iter2(daqSources);
823   TObjString* source = 0;
824   Int_t i=0;
825   Bool_t resMBCal=kTRUE;
826   
827   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
828      TString stringMBFileName = GetFile(kDAQ, "MBCALIB", source->GetName());
829      if(stringMBFileName.Length() <= 0){
830        Log(Form("No MBCALIB file from source %s!", source->GetName()));
831        return 18;
832      }
833      const char* mbFileName = stringMBFileName.Data();
834      Log(Form("\t Getting file #%d: %s from %s\n",++i,mbFileName,source->GetName()));
835      //
836      // --- Initializing calibration object
837      AliZDCMBCalib *mbCalib = new AliZDCMBCalib("ZDC");
838      // --- Reading file with calibration data
839      if(mbFileName){
840        TFile * fileHistos = TFile::Open(mbFileName);
841        Log(Form("File %s connected to process data from CALIBRATION_MB events", mbFileName));
842        //
843        fileHistos->cd();
844        TH2F *hZDCvsZEM = (TH2F*)  fileHistos->Get("hZDCvsZEM");
845        TH2F *hZDCCvsZEM = (TH2F*) fileHistos->Get("hZDCCvsZEM");
846        TH2F *hZDCAvsZEM = (TH2F*) fileHistos->Get("hZDCAvsZEM");
847        //
848        mbCalib->SetZDCvsZEM(hZDCvsZEM);
849        mbCalib->SetZDCCvsZEM(hZDCCvsZEM);
850        mbCalib->SetZDCAvsZEM(hZDCAvsZEM);
851        //
852        //fileHistos->Close();
853      }
854      else{
855        Log(Form("File %s not found", mbFileName));
856        return 14;
857      }
858      // 
859      AliCDBMetaData metaData;
860      metaData.SetBeamPeriod(0);
861      metaData.SetResponsible("Chiara Oppedisano");
862      metaData.SetComment("Filling AliZDCMBCalib object");  
863      //
864      //mbCalib->Dump();
865      //
866      resMBCal = Store("Calib","MBCalib",mbCalib, &metaData, 0, kTRUE);
867        printf(" here 1000\n");
868      if(resMBCal==kFALSE) return 19;
869   }
870   delete daqSources; daqSources = 0;
871   
872   return 0;
873 }
874
875 //______________________________________________________________________________________________
876 UInt_t AliZDCPreprocessor::Process(TMap* dcsAliasMap)
877 {
878  UInt_t resDCS = 0;
879  UInt_t resChMap=0;
880  UInt_t resEnergyCalib=0, resPedestalCalib=0, resLaserCalib=0, resMBCalib=0;
881
882  // ************************* Process DCS data ****************************
883  if(ProcessDCS()) resDCS = ProcessDCSData(dcsAliasMap);
884   
885  // ********************************* From DAQ ************************************
886
887  const char* beamType = GetRunParameter("beamType");
888  TString runType = GetRunType();
889  Float_t beamEnergy = (Float_t)(((TString)GetRunParameter("beamEnergy")).Atof()); 
890  printf("\t **** AliZDCPreprocessor -> runType %s, beamType %s,  beamEnergy %1.0f ****\n",
891         runType.Data(),beamType,beamEnergy);
892
893  // ******************************************
894  // ADC channel mapping
895  // ******************************************
896  resChMap = ProcessChMap();
897  
898  // ******************************************
899  // Calibration param. for p-p data (all = 1)
900  // ******************************************
901  // NO ENERGY CALIBRATION -> coefficients set to 1.
902  // Temp -> also inter-calibration coefficients are set to 1.
903  if((strcmp(beamType,"p-p")==0) || (strcmp(beamType,"P-P")==0) 
904      || (strcmp(beamType,"P-A")==0) || (strcmp(beamType,"A-P")==0)) 
905         resEnergyCalib = ProcessppData();
906  
907  // *****************************************************
908  // CALIBRATION_EMD -> Energy calibration and equalization
909  // *****************************************************
910  if((strcmp(beamType,"A-A")==0) && (runType.CompareTo("CALIBRATION_EMD")==0)) 
911    resEnergyCalib =  ProcessCalibData(beamEnergy);
912  
913  // *****************************************************
914  // STANDALONE_PEDESTALS -> Pedestal subtraction 
915  // *****************************************************
916  if(runType.CompareTo("STANDALONE_PEDESTAL")==0) resPedestalCalib = ProcessPedestalData();
917  
918  // *****************************************************
919  // STANDALONE_LASER -> Signal stability and ageing 
920  // *****************************************************
921  else if(runType.CompareTo("STANDALONE_LASER")==0) resLaserCalib = ProcessLaserData();
922
923  // *****************************************************
924  // CALIBRATION_MB -> Signal stability and ageing 
925  // *****************************************************
926  else if(runType.CompareTo("CALIBRATION_MB")==0) resMBCalib = ProcessMBCalibData();
927  
928  if(resDCS!=0)                return resDCS;
929  else if(resChMap!=0)         return resChMap;
930  else if(resEnergyCalib!=0)   return resEnergyCalib;
931  else if(resPedestalCalib!=0) return resPedestalCalib;
932  else if(resLaserCalib!=0)    return resLaserCalib;
933  else if(resMBCalib!=0)       return resMBCalib;
934  
935  return 0;
936   
937 }