]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCPreprocessor.cxx
651ff55761a38f93c55ba3e6eb4ae31f7e42fd14
[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 = kTRUE;
271     }
272     for(Int_t i=0; i<kNScch; i++){
273       if(  (scMap[i][2] != chMap->GetScChannel(i)) 
274         || (scMap[i][4] != chMap->GetScDetector(i)) 
275         || (scMap[i][5] != chMap->GetScSector(i)))
276          scMapUpdated = kTRUE;
277     }
278   }
279   if(adcMapUpdated || scMapUpdated) updateOCDB = kTRUE;
280   //
281   Bool_t resChMapStore = kTRUE;
282   if(updateOCDB==kTRUE){
283     Log(" A new entry ZDC/Calib/ChMap will be created");
284     //
285     // --- Initializing mapping calibration object
286     AliZDCChMap *mapCalib = new AliZDCChMap("ZDC");
287     // Writing channel map in the OCDB
288     for(Int_t k=0; k<kNModules; k++){
289       mapCalib->SetModuleMap(k, modMap[k][0], modMap[k][1], modMap[k][2]);
290     }
291     for(Int_t k=0; k<kNch; k++){
292       mapCalib->SetADCModule(k,adcMap[k][1]);
293       mapCalib->SetADCChannel(k,adcMap[k][2]);
294       mapCalib->SetDetector(k,adcMap[k][4]);
295       mapCalib->SetSector(k,adcMap[k][5]);
296     }
297     for(Int_t k=0; k<kNScch; k++){
298        mapCalib->SetScChannel(k, scMap[k][2]);
299        mapCalib->SetScDetector(k, scMap[k][4]);
300        mapCalib->SetScSector(k, scMap[k][5]);
301     }
302     //
303     mapCalib->Print("");
304     // 
305     AliCDBMetaData metaData;
306     metaData.SetBeamPeriod(0);
307     metaData.SetResponsible("Chiara Oppedisano");
308     metaData.SetComment("Filling AliZDCChMap object");  
309     //
310     resChMapStore = Store("Calib","ChMap",mapCalib, &metaData, 0, kTRUE);
311     printf("  Mapping object stored in OCDB\n");
312   }
313   else{
314     Log(" ZDC/Calib/ChMap entry in OCDB is valid and won't be updated");
315     resChMapStore = kTRUE;
316   }
317   
318   delete daqSource; daqSource=0;
319   
320   if(resChMapStore==kFALSE) return 5;
321   else return 0;
322
323 }
324
325 //______________________________________________________________________________________________
326 UInt_t AliZDCPreprocessor::ProcessppData()
327 {
328    Bool_t resEnCal=kTRUE, resTowCal=kTRUE;
329   
330    // *********** Energy calibration
331    // --- Cheking if there is already the entry in the OCDB
332    AliCDBEntry *cdbEnEntry = GetFromOCDB("Calib", "EnergyCalib");
333    if(!cdbEnEntry){   
334      Log(Form(" ZDC/Calib/EnergyCalib entry will be created"));
335      // --- Initializing calibration object
336      AliCDBMetaData metaData;
337      metaData.SetBeamPeriod(0);
338      metaData.SetResponsible("Chiara Oppedisano");
339      //
340      AliZDCEnCalib *eCalib = new AliZDCEnCalib("ZDC");
341      for(Int_t j=0; j<6; j++) eCalib->SetEnCalib(j,1.);
342      metaData.SetComment("AliZDCEnCalib object");  
343      //eCalib->Print("");
344      resEnCal = Store("Calib", "EnergyCalib", eCalib, &metaData, 0, kTRUE);
345    }
346    else{ 
347      // if entry exists it is still valid (=1 for all runs!)
348      Log(Form(" Valid ZDC/Calib/EnergyCalib object already existing in OCDB!!!"));
349      resEnCal = kTRUE;
350    }
351    
352    if(resEnCal==kFALSE) return 6;
353
354    //
355    // *********** Tower inter-calibration
356    // --- Cheking if there is already the entry in the OCDB
357    AliCDBEntry *cdbTowEntry = GetFromOCDB("Calib", "TowerCalib");
358    if(!cdbTowEntry){   
359      AliZDCTowerCalib *towCalib = new AliZDCTowerCalib("ZDC");
360      for(Int_t j=0; j<5; j++){  
361         towCalib->SetZN1EqualCoeff(j, 1.);
362         towCalib->SetZP1EqualCoeff(j, 1.);
363         towCalib->SetZN2EqualCoeff(j, 1.);
364         towCalib->SetZP2EqualCoeff(j, 1.);  
365      }
366      //towCalib->Print("");
367      // 
368      AliCDBMetaData metaData;
369      metaData.SetBeamPeriod(0);
370      metaData.SetResponsible("Chiara Oppedisano");
371      metaData.SetComment("AliZDCTowerCalib object");  
372      //
373      resTowCal = Store("Calib", "TowerCalib", towCalib, &metaData, 0, kTRUE);
374    }
375    else{ 
376      // if entry exists it is still valid (=1 for all runs!)
377      Log(Form(" Valid ZDC/Calib/TowerCalib object already existing in OCDB!!!"));
378      resTowCal = kTRUE;
379    }
380    
381    if(resTowCal==kFALSE) return 7;
382    
383    return 0;
384 }
385
386 //______________________________________________________________________________________________
387 UInt_t AliZDCPreprocessor::ProcessCalibData(Float_t beamEnergy)
388 {
389   TList* daqSources = GetFileSources(kDAQ, "EMDENERGYCALIB");
390   if(!daqSources){
391     AliError(Form("No sources for CALIBRATION_EMD run %d !", fRun));
392     return 8;
393   }
394   Log("\t List of DAQ sources for EMDENERGYCALIB id: "); daqSources->Print();
395   //
396   TIter iter2(daqSources);
397   TObjString* source = 0;
398   Int_t i=0;
399   Bool_t resEnCal=kTRUE, resTowCal=kTRUE;
400   
401   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
402     TString stringEMDFileName = GetFile(kDAQ, "EMDENERGYCALIB", source->GetName());
403     if(stringEMDFileName.Length() <= 0){
404       Log(Form("No file from source %s!", source->GetName()));
405       return 8;
406     }
407     const char* emdFileName = stringEMDFileName.Data();
408     Log(Form("\t Getting file #%d: %s from %s\n",++i,emdFileName,source->GetName()));
409     //
410     // --- Initializing energy calibration object
411     AliZDCEnCalib *eCalib = new AliZDCEnCalib("ZDC");
412     // --- Reading file with calibration data
413     if(emdFileName){
414       FILE *file;
415       if((file = fopen(emdFileName,"r")) == NULL){
416         printf("Cannot open file %s \n",emdFileName);
417         return 8;
418       }
419       Log(Form("File %s connected to process data from EM dissociation events", emdFileName));
420       //
421       Float_t fitValEMD[6];
422       for(Int_t j=0; j<6; j++){        
423         if(j<6){
424           int iread = fscanf(file,"%f",&fitValEMD[j]);
425           if(iread==0) AliDebug(3," Failing reading daa from EMD calibration data file");
426           eCalib->SetEnCalib(j, beamEnergy/fitValEMD[j]);
427         }
428       }
429       //
430       fclose(file);
431     }
432     else{
433       Log(Form("File %s not found", emdFileName));
434       return 8;
435     }
436     //eCalib->Print("");
437     // 
438     AliCDBMetaData metaData;
439     metaData.SetBeamPeriod(0);
440     metaData.SetResponsible("Chiara Oppedisano");
441     metaData.SetComment("Filling AliZDCEnCalib object");  
442     //
443     resEnCal = Store("Calib","EnergyCalib",eCalib, &metaData, 0, kTRUE);
444     if(resEnCal==kFALSE) return 6;
445   }
446   delete daqSources; daqSources = 0;
447   
448   TList* daqSourcesH = GetFileSources(kDAQ, "EMDTOWERCALIB");
449   if(!daqSourcesH){
450     AliError(Form("No sources for CALIBRATION_EMD run %d !", fRun));
451     return 9;
452   }
453   Log("\t List of DAQ sources for EMDTOWERCALIB id: "); daqSourcesH->Print();
454   //
455   TIter iter2H(daqSourcesH);
456   TObjString* sourceH = 0;
457   Int_t iH=0;
458   while((sourceH = dynamic_cast<TObjString*> (iter2H.Next()))){
459     TString stringtowEMDFileName = GetFile(kDAQ, "EMDTOWERCALIB", sourceH->GetName());
460     if(stringtowEMDFileName.Length() <= 0){
461       Log(Form("No file from source %s!", sourceH->GetName()));
462       return 9;
463     }
464     const char * towEMDFileName = stringtowEMDFileName.Data();
465     Log(Form("\t Getting file #%d: %s from source %s\n",++iH,towEMDFileName,sourceH->GetName()));
466     // --- Initializing energy calibration object
467     AliZDCTowerCalib *towCalib = new AliZDCTowerCalib("ZDC");
468     // --- Reading file with calibration data
469     if(towEMDFileName){
470       FILE *file;
471       if((file = fopen(towEMDFileName,"r")) == NULL){
472         printf("Cannot open file %s \n",towEMDFileName);
473         return 9;
474       }
475       //
476       Float_t equalCoeff[4][5];
477       for(Int_t j=0; j<4; j++){        
478          for(Int_t k=0; k<5; k++){
479            int leggi = fscanf(file,"%f",&equalCoeff[j][k]);
480            if(leggi==0) AliDebug(3," Failing reading data from EMD calibration file");
481            if(j==0)      towCalib->SetZN1EqualCoeff(k, equalCoeff[j][k]);
482            else if(j==1) towCalib->SetZP1EqualCoeff(k, equalCoeff[j][k]);
483            else if(j==2) towCalib->SetZN2EqualCoeff(k, equalCoeff[j][k]);
484            else if(j==3) towCalib->SetZP2EqualCoeff(k, equalCoeff[j][k]);  
485          }
486       }
487       //
488       fclose(file);
489     }
490     else{
491       Log(Form("File %s not found", towEMDFileName));
492       return 9;
493     }
494     //towCalib->Print("");
495     // 
496     AliCDBMetaData metaData;
497     metaData.SetBeamPeriod(0);
498     metaData.SetResponsible("Chiara Oppedisano");
499     metaData.SetComment("Filling AliZDCTowerCalib object");  
500     //
501     resTowCal = Store("Calib","TowerCalib",towCalib, &metaData, 0, kTRUE);
502     if(resTowCal==kFALSE) return 7;
503   }
504   delete daqSourcesH; daqSourcesH = 0;
505   
506      
507   return 0;
508 }
509
510 //______________________________________________________________________________________________
511 UInt_t AliZDCPreprocessor::ProcessPedestalData()
512 {
513   TList* daqSources = GetFileSources(kDAQ, "PEDESTALDATA");
514   if(!daqSources){
515     Log(Form("No source for STANDALONE_PEDESTAL run %d !", fRun));
516     return 10;
517   }
518   if(daqSources->GetEntries()==0) return 10;
519   Log("\t List of DAQ sources for PEDESTALDATA id: "); daqSources->Print();
520   //
521   TIter iter(daqSources);
522   TObjString* source;
523   Int_t i=0;
524   Bool_t resPedCal=kTRUE, resPedHist=kTRUE;
525   
526   while((source = dynamic_cast<TObjString*> (iter.Next()))){
527      TString stringPedFileName = GetFile(kDAQ, "PEDESTALDATA", source->GetName());
528      if(stringPedFileName.Length() <= 0){
529         Log(Form("No PEDESTALDATA file from source %s!", source->GetName()));
530         return 10;
531      }
532      const char* pedFileName = stringPedFileName.Data();
533      Log(Form("\t Getting file #%d: %s from %s\n",++i,pedFileName,source->GetName()));
534      //
535      // --- Initializing pedestal calibration object
536      AliZDCPedestals *pedCalib = new AliZDCPedestals("ZDC");
537      // --- Reading file with pedestal calibration data
538      // no. ADCch = (22 signal ch. + 2 reference PMs) * 2 gain chain = 48
539      const Int_t knZDCch = 48;
540      FILE *file;
541      if((file = fopen(pedFileName,"r")) == NULL){
542        printf("Cannot open file %s \n",pedFileName);
543        return 10;
544      }
545      Log(Form("File %s connected to process pedestal data", pedFileName));
546      Float_t pedVal[(3*knZDCch)][2];
547      for(Int_t k=0; k<(2*knZDCch); k++){
548         for(Int_t j=0; j<2; j++){
549            int aleggi = fscanf(file,"%f",&pedVal[k][j]);
550            if(aleggi==0) AliDebug(3," Failing reading data from pedestal file");
551            //if(j==1) printf("pedVal[%d] -> %f, %f \n",k,pedVal[k][0],pedVal[k][1]);
552         }
553         if(k<knZDCch){
554           pedCalib->SetMeanPed(k,pedVal[k][0]);
555           pedCalib->SetMeanPedWidth(k,pedVal[k][1]);
556         }
557         else if(k>=knZDCch && k<(2*knZDCch)){
558           pedCalib->SetOOTPed(k-knZDCch,pedVal[k][0]);
559           pedCalib->SetOOTPedWidth(k-knZDCch,pedVal[k][1]);
560         }
561         else if(k>=(2*knZDCch) && k<(3*knZDCch)){
562           pedCalib->SetPedCorrCoeff(k-(2*knZDCch),pedVal[k][0],pedVal[k][1]);
563         }
564      }
565      fclose(file);
566      //pedCalib->Print("");
567      // 
568      AliCDBMetaData metaData;
569      metaData.SetBeamPeriod(0);
570      metaData.SetResponsible("Chiara Oppedisano");
571      metaData.SetComment("Filling AliZDCPedestals object");  
572      //
573      resPedCal = Store("Calib","Pedestals",pedCalib, &metaData, 0, kTRUE);
574      if(resPedCal==kFALSE) return 11;
575   }
576   delete daqSources; daqSources = 0;
577   
578   TList* daqSourceH = GetFileSources(kDAQ, "PEDESTALHISTOS");
579   if(!daqSourceH){
580     Log(Form("No source for PEDESTALHISTOS id run %d !", fRun));
581     return 12;
582   }
583   Log("\t List of DAQ sources for PEDESTALHISTOS id: "); daqSourceH->Print();
584   //
585   TIter iterH(daqSourceH);
586   TObjString* sourceH = 0;
587   Int_t iH=0;
588   while((sourceH = dynamic_cast<TObjString*> (iterH.Next()))){
589      TString stringPedFileName = GetFile(kDAQ, "PEDESTALHISTOS", sourceH->GetName());
590      if(stringPedFileName.Length() <= 0){
591         Log(Form("No PEDESTALHISTOS file from source %s!", sourceH->GetName()));
592         return 12;
593      }
594      const char* pedFileName = stringPedFileName.Data();
595      Log(Form("\t Getting file #%d: %s from %s\n",++iH, pedFileName, sourceH->GetName()));
596      resPedHist = StoreReferenceFile(pedFileName, "pedestalReference.root");
597      if(resPedHist==kFALSE) return 13;
598   }
599   delete daqSourceH; daqSourceH=0;
600   
601   return 0;
602 }
603
604 //______________________________________________________________________________________________
605 UInt_t AliZDCPreprocessor::ProcessLaserData()
606 {
607   TList* daqSources = GetFileSources(kDAQ, "LASERDATA");
608   if(!daqSources){
609     AliError(Form("No sources for STANDALONE_LASER run %d !", fRun));
610     return 14;
611   }
612   if(daqSources->GetEntries()==0) return 14;
613   Log("\t List of DAQ sources for LASERDATA id: "); daqSources->Print();
614   //
615   TIter iter2(daqSources);
616   TObjString* source = 0;
617   Int_t i=0;
618   Bool_t resLaserCal=kTRUE, resLaserHist=kTRUE;
619   
620   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
621      TString stringLaserFileName = GetFile(kDAQ, "LASERDATA", source->GetName());
622      if(stringLaserFileName.Length() <= 0){
623        Log(Form("No LASER file from source %s!", source->GetName()));
624        return 14;
625      }
626      const char* laserFileName = stringLaserFileName.Data();
627      Log(Form("\t Getting file #%d: %s from %s\n",++i,laserFileName,source->GetName()));
628      //
629      // --- Initializing pedestal calibration object
630      AliZDCLaserCalib *lCalib = new AliZDCLaserCalib("ZDC");
631      // --- Reading file with pedestal calibration data
632      if(laserFileName){
633        FILE *file;
634        if((file = fopen(laserFileName,"r")) == NULL){
635          printf("Cannot open file %s \n",laserFileName);
636          return 14;
637        }
638        Log(Form("File %s connected to process data from LASER events", laserFileName));
639        //
640        Float_t ivalRead[22][4]; 
641        for(Int_t j=0; j<22; j++){
642           for(Int_t k=0; k<4; k++){
643             int aleggi = fscanf(file,"%f",&ivalRead[j][k]);
644             if(aleggi==0) AliDebug(3," Failng reading data from laser file");
645             //printf(" %d %1.0f  ",k, ivalRead[j][k]);
646           }
647           lCalib->SetDetector(j, (Int_t) ivalRead[j][0]);
648           lCalib->SetSector(j, (Int_t) ivalRead[j][1]);
649           lCalib->SetfPMValue(j, ivalRead[j][2]);
650           lCalib->SetfPMWidth(j, ivalRead[j][3]);
651        }
652        fclose(file);
653      }
654      else{
655        Log(Form("File %s not found", laserFileName));
656        return 14;
657      }
658      //lCalib->Print("");
659      // 
660      AliCDBMetaData metaData;
661      metaData.SetBeamPeriod(0);
662      metaData.SetResponsible("Chiara Oppedisano");
663      metaData.SetComment("Filling AliZDCLaserCalib object");  
664      //
665      resLaserCal = Store("Calib","LaserCalib",lCalib, &metaData, 0, kTRUE);
666      if(resLaserCal==kFALSE) return 15;
667   }
668   delete daqSources; daqSources = 0;
669
670   TList* daqSourceH = GetFileSources(kDAQ, "LASERHISTOS");
671   if(!daqSourceH){
672     AliError(Form("No sources for STANDALONE_LASER run %d !", fRun));
673     return 16;
674   }
675   Log("\t List of DAQ sources for LASERHISTOS id: "); daqSourceH->Print();
676   //
677   TIter iter2H(daqSourceH);
678   TObjString* sourceH = 0;
679   Int_t iH=0;
680   while((sourceH = dynamic_cast<TObjString*> (iter2H.Next()))){
681      Log(Form("\t Getting file #%d\n",++iH));
682      TString stringLaserFileName = GetFile(kDAQ, "LASERHISTOS", sourceH->GetName());
683      if(stringLaserFileName.Length() <= 0){
684        Log(Form("No LASER file from source %s!", sourceH->GetName()));
685        return 16;
686      }
687      resLaserHist = StoreReferenceFile(stringLaserFileName.Data(), "laserReference.root");
688      //
689      if(resLaserHist==kFALSE) return 17;
690   }
691   delete daqSourceH; daqSourceH = 0;
692   
693   return 0;
694 }
695
696
697 //______________________________________________________________________________________________
698 UInt_t AliZDCPreprocessor::ProcessMBCalibData()
699 {
700   TList* daqSources = GetFileSources(kDAQ, "MBCALIB");
701   if(!daqSources){
702     AliError(Form("No sources for CALIBRATION_MB run %d !", fRun));
703     return 18;
704   }
705   if(daqSources->GetEntries()==0) return 18;
706   Log("\t List of DAQ sources for MBCALIB id: "); daqSources->Print();
707   //
708   TIter iter2(daqSources);
709   TObjString* source = 0;
710   Int_t i=0;
711   Bool_t resMBCal=kTRUE;
712   
713   while((source = dynamic_cast<TObjString*> (iter2.Next()))){
714      TString stringMBFileName = GetFile(kDAQ, "MBCALIB", source->GetName());
715      if(stringMBFileName.Length() <= 0){
716        Log(Form("No MBCALIB file from source %s!", source->GetName()));
717        return 18;
718      }
719      const char* mbFileName = stringMBFileName.Data();
720      Log(Form("\t Getting file #%d: %s from %s\n",++i,mbFileName,source->GetName()));
721      //
722      // --- Initializing calibration object
723      AliZDCMBCalib *mbCalib = new AliZDCMBCalib("ZDC");
724      // --- Reading file with calibration data
725      if(mbFileName){
726        TFile * fileHistos = TFile::Open(mbFileName);
727        Log(Form("File %s connected to process data from CALIBRATION_MB events", mbFileName));
728        //
729        fileHistos->cd();
730        TH2F *hZDCvsZEM = (TH2F*)  fileHistos->Get("hZDCvsZEM");
731        TH2F *hZDCCvsZEM = (TH2F*) fileHistos->Get("hZDCCvsZEM");
732        TH2F *hZDCAvsZEM = (TH2F*) fileHistos->Get("hZDCAvsZEM");
733        //
734        mbCalib->SetZDCvsZEM(hZDCvsZEM);
735        mbCalib->SetZDCCvsZEM(hZDCCvsZEM);
736        mbCalib->SetZDCAvsZEM(hZDCAvsZEM);
737        //
738        //fileHistos->Close();
739      }
740      else{
741        Log(Form("File %s not found", mbFileName));
742        return 14;
743      }
744      // 
745      AliCDBMetaData metaData;
746      metaData.SetBeamPeriod(0);
747      metaData.SetResponsible("Chiara Oppedisano");
748      metaData.SetComment("Filling AliZDCMBCalib object");  
749      //
750      //mbCalib->Dump();
751      //
752      resMBCal = Store("Calib","MBCalib",mbCalib, &metaData, 0, kTRUE);
753        printf(" here 1000\n");
754      if(resMBCal==kFALSE) return 19;
755   }
756   delete daqSources; daqSources = 0;
757   
758   return 0;
759 }
760
761 //______________________________________________________________________________________________
762 UInt_t AliZDCPreprocessor::Process(TMap* dcsAliasMap)
763 {
764  UInt_t resDCS = 0;
765  UInt_t resChMap=0;
766  UInt_t resEnergyCalib=0, resPedestalCalib=0, resLaserCalib=0, resMBCalib=0;
767
768  // ************************* Process DCS data ****************************
769  if(ProcessDCS()) resDCS = ProcessDCSData(dcsAliasMap);
770   
771  // ********************************* From DAQ ************************************
772
773  const char* beamType = GetRunParameter("beamType");
774  TString runType = GetRunType();
775  Float_t beamEnergy = (Float_t)(((TString)GetRunParameter("beamEnergy")).Atof()); 
776  printf("\t **** AliZDCPreprocessor -> runType %s, beamType %s,  beamEnergy %1.0f ****\n",
777         runType.Data(),beamType,beamEnergy);
778
779  // ******************************************
780  // ADC channel mapping
781  // ******************************************
782  resChMap = ProcessChMap();
783  
784  // ******************************************
785  // Calibration param. for p-p data (all = 1)
786  // ******************************************
787  // NO ENERGY CALIBRATION -> coefficients set to 1.
788  // Temp -> also inter-calibration coefficients are set to 1.
789  if((strcmp(beamType,"p-p")==0) || (strcmp(beamType,"P-P")==0)) resEnergyCalib = ProcessppData();
790  
791  // *****************************************************
792  // CALIBRATION_EMD -> Energy calibration and equalization
793  // *****************************************************
794  if((strcmp(beamType,"A-A")==0) && (runType.CompareTo("CALIBRATION_EMD")==0)) 
795    resEnergyCalib =  ProcessCalibData(beamEnergy);
796  
797  // *****************************************************
798  // STANDALONE_PEDESTALS -> Pedestal subtraction 
799  // *****************************************************
800  if(runType.CompareTo("STANDALONE_PEDESTAL")==0) resPedestalCalib = ProcessPedestalData();
801  
802  // *****************************************************
803  // STANDALONE_LASER -> Signal stability and ageing 
804  // *****************************************************
805  else if(runType.CompareTo("STANDALONE_LASER")==0) resLaserCalib = ProcessLaserData();
806
807  // *****************************************************
808  // CALIBRATION_MB -> Signal stability and ageing 
809  // *****************************************************
810  else if(runType.CompareTo("CALIBRATION_MB")==0) resMBCalib = ProcessMBCalibData();
811  
812
813  if(resDCS!=0)                return resDCS;
814  else if(resChMap!=0)         return resChMap;
815  else if(resEnergyCalib!=0)   return resEnergyCalib;
816  else if(resPedestalCalib!=0) return resPedestalCalib;
817  else if(resLaserCalib!=0)    return resLaserCalib;
818  else if(resMBCalib!=0)       return resMBCalib;
819  
820  return 0;
821   
822 }