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