]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFPreprocessor.cxx
Added protection against null list of file source
[u/mrichter/AliRoot.git] / TOF / AliTOFPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //#include <Riostream.h>
19 //#include <stdio.h>
20 //#include <stdlib.h>
21
22 #include <TFile.h>
23 //#include <TH1.h>
24 //#include <TH1F.h>
25 //#include <TH1S.h>
26 #include <TH2S.h>
27 #include <TMath.h>
28 #include <TObjArray.h>
29 #include <TObjString.h>
30 #include <TTimeStamp.h>
31
32 #include "AliCDBMetaData.h"
33 #include "AliCDBEntry.h"
34 #include "AliLog.h"
35 #include "AliTOFChannelOnline.h"
36 #include "AliTOFChannelOnlineStatus.h"
37 #include "AliTOFDataDCS.h"
38 #include "AliTOFGeometry.h"
39 #include "AliTOFPreprocessor.h"
40
41 //class TF1;
42 //class AliDCSValue;
43 //class AliTOFGeometry;
44
45 // TOF preprocessor class.
46 // It takes data from DCS and passes them to the class AliTOFDataDCS, which
47 // processes them. The result is then written to the CDB.
48 // analogously, it takes data form DAQ (both at Run level and inclusive - 
49 // of all the runs - level, processes them, and stores both Reference Data
50 // and Online Calibration files in the CDB. 
51 // Processing of Pulser/Noise Run data and of TOF FEE DCS map
52
53 // return codes:
54 // return=0 : all ok
55 // return=1 : no DCS input data Map
56 // return=2 : no DCS input data processing
57 // return=3 : no DCS processed data was stored in Ref Data
58 // return=4 : no DAQ input for Ref Data
59 // return=5 : failed to store DAQ Ref Data
60 // return=6 : failed to retrieve DAQ data for calibration 
61 // return=7 : problems in processing histos in the input DAQ file 
62 // return=8 : failed to store Online Delays
63 // return=9 : failed to store Reference Data for Pulser
64 // return=10: failed to retrieve Pulser data 
65 // return=11: failed to store Pulser map in OCDB
66 // return=12: failed to store Reference Data for Noise
67 // return=13: failed to retrieve Noise data 
68 // return=14: failed to store Noise map in OCDB
69
70 ClassImp(AliTOFPreprocessor)
71
72 const Int_t    AliTOFPreprocessor::fgkBinRangeAve = 13;    // number of bins where to calculate the mean 
73 const Double_t AliTOFPreprocessor::fgkIntegralThr = 100;   // min number of entries to perform computation of delay per channel 
74 const Double_t AliTOFPreprocessor::fgkThrPar      = 0.013; // parameter used to trigger the calculation of the delay
75
76 //_____________________________________________________________________________
77
78 AliTOFPreprocessor::AliTOFPreprocessor(AliShuttleInterface* shuttle) :
79   AliPreprocessor("TOF", shuttle),
80   fData(0),
81   fh2(0),
82   fCal(0),
83   fCalStatus(0),
84   fNChannels(0),
85   fStoreRefData(kTRUE),
86   fFDRFlag(kTRUE)
87 {
88   // constructor
89
90 }
91
92 //_____________________________________________________________________________
93
94 AliTOFPreprocessor::~AliTOFPreprocessor()
95 {
96   // destructor
97   if (fData){
98     delete fData;
99     fData = 0;
100   }
101   if (fh2){
102     delete fh2;
103     fh2 = 0;
104   }
105   if (fCal){
106     //    fCal->Clear();
107     delete fCal;
108     fCal = 0;
109   }
110   if (fCalStatus){
111     delete fCalStatus;
112     fCalStatus = 0;
113   }
114 }
115
116 //______________________________________________________________________________
117 void AliTOFPreprocessor::Initialize(Int_t run, UInt_t startTime,
118         UInt_t endTime)
119 {
120   // Creates AliTOFDataDCS object
121
122   AliPreprocessor::Initialize(run, startTime, endTime);
123
124         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
125                 TTimeStamp(startTime).AsString(),
126                 TTimeStamp(endTime).AsString()));
127
128         fData = new AliTOFDataDCS(fRun, fStartTime, fEndTime);
129         fData->SetFDRFlag(fFDRFlag);
130         fh2 = 0x0;
131         fNChannels = AliTOFGeometry::NSectors()*(2*(AliTOFGeometry::NStripC()+AliTOFGeometry::NStripB())+AliTOFGeometry::NStripA())*AliTOFGeometry::NpadZ()*AliTOFGeometry::NpadX();
132         fCal = new TObjArray(fNChannels);
133         fCal->SetOwner();
134         for (Int_t ich = 0; ich<fNChannels; ich ++){
135           AliTOFChannelOnline * calChOnline = new AliTOFChannelOnline();
136           fCal->AddAt(calChOnline,ich);
137         }
138         fCalStatus = new TObjArray(fNChannels);
139         fCalStatus->SetOwner();
140         for (Int_t ich = 0; ich<fNChannels; ich ++){
141           AliTOFChannelOnlineStatus * calChOnlineStatus = new AliTOFChannelOnlineStatus();
142           fCalStatus->AddAt(calChOnlineStatus,ich);
143         }
144 }
145 //_____________________________________________________________________________
146 Bool_t AliTOFPreprocessor::ProcessDCS(){
147
148   // check whether DCS should be processed or not...
149
150   TString runType = GetRunType();
151   Log(Form("RunType %s",runType.Data()));
152
153   if (runType != "PHYSICS"){
154     return kFALSE;
155   }
156
157   return kTRUE;
158 }
159 //_____________________________________________________________________________
160
161 UInt_t AliTOFPreprocessor::ProcessDCSDataPoints(TMap* dcsAliasMap)
162 {
163   // Fills data into a AliTOFDataDCS object
164
165
166   Log("Processing DCS DP");
167   TH1::AddDirectory(0);
168
169   Bool_t resultDCSMap=kFALSE;
170   Bool_t resultDCSStore=kFALSE;
171
172   // processing DCS
173
174   if (!dcsAliasMap){
175     Log("No DCS map found: TOF exiting from Shuttle");
176     return 1;// return error Code for DCS input data not found 
177   }
178   else {
179   // The processing of the DCS input data is forwarded to AliTOFDataDCS
180     resultDCSMap=fData->ProcessData(*dcsAliasMap);
181     if(!resultDCSMap){
182       Log("Some problems occurred while processing DCS data, TOF exiting from Shuttle");
183       return 2;// return error Code for processed DCS data not stored 
184     }
185     else{
186       AliCDBMetaData metaDataDCS;
187       metaDataDCS.SetBeamPeriod(0);
188       metaDataDCS.SetResponsible("Chiara Zampolli");
189       metaDataDCS.SetComment("This preprocessor fills an AliTOFDataDCS object.");
190       AliInfo("Storing DCS Data");
191       resultDCSStore = StoreReferenceData("Calib","DCSData",fData, &metaDataDCS);
192       if (!resultDCSStore){
193         Log("Some problems occurred while storing DCS data results in Reference Data, TOF exiting from Shuttle");
194         return 3;// return error Code for processed DCS data not stored 
195                  // in reference data
196         
197       }
198     }
199   }
200   return 0;
201 }
202 //_____________________________________________________________________________
203
204 UInt_t AliTOFPreprocessor::ProcessOnlineDelays()
205 {
206   // Processing data from DAQ for online calibration 
207
208   Log("Processing DAQ delays");
209
210   TH1::AddDirectory(0);
211
212   Bool_t resultDAQRef=kFALSE;
213   Bool_t resultTOFPP=kFALSE;
214   // processing DAQ
215   
216   TFile * daqFile=0x0;
217   
218   if(fStoreRefData){
219     //retrieving data at Run level
220     TList* list = GetFileSources(kDAQ, "RUNLevel");
221     if (list !=0x0 && list->GetEntries()!=0)
222       {
223         AliInfo("The following sources produced files with the id RUNLevel");
224         list->Print();
225         for (Int_t jj=0;jj<list->GetEntries();jj++){
226           TObjString * str = dynamic_cast<TObjString*> (list->At(jj));
227           AliInfo(Form("found source %s", str->String().Data()));
228           // file to be stored run per run
229           TString fileNameRun = GetFile(kDAQ, "RUNLevel", str->GetName());
230           if (fileNameRun.Length()>0){
231             AliInfo(Form("Got the file %s, now we can store the Reference Data for the current Run.", fileNameRun.Data()));
232             daqFile = new TFile(fileNameRun.Data(),"READ");
233             fh2 = (TH2S*) daqFile->Get("htof");
234             AliCDBMetaData metaDataHisto;
235             metaDataHisto.SetBeamPeriod(0);
236             metaDataHisto.SetResponsible("Chiara Zampolli");
237             metaDataHisto.SetComment("This preprocessor stores the array of histos object as Reference Data.");
238             AliInfo("Storing Reference Data");
239             resultDAQRef = StoreReferenceData("Calib","DAQData",fh2, &metaDataHisto);
240             if (!resultDAQRef){
241               Log("some problems occurred::No Reference Data stored, TOF exiting from Shuttle");
242               return 5;//return error code for failure in storing Ref Data 
243             }
244             daqFile->Close();
245             delete daqFile;
246           }
247           
248           else{
249             Log("The input data file from DAQ (run-level) was not found, TOF exiting from Shuttle "); 
250             return 4;//return error code for failure in retrieving Ref Data 
251           }
252         }
253         delete list;
254       }
255     else{
256       Log("The input data file list from DAQ (run-level) was not found, TOF exiting from Shuttle "); 
257       return 4;//return error code for failure in retrieving Ref Data 
258     }   
259   }
260
261
262 //Total files, with cumulative histos
263   
264   TList* listTot = GetFileSources(kDAQ, "DELAYS");
265   if (listTot !=0x0 && listTot->GetEntries()!=0)
266     {
267       AliInfo("The following sources produced files with the id DELAYS");
268       listTot->Print();
269       for (Int_t jj=0;jj<listTot->GetEntries();jj++){
270         TObjString * str = dynamic_cast<TObjString*> (listTot->At(jj));
271         AliInfo(Form("found source %s", str->String().Data()));
272
273         // file with summed histos, to extract calib params
274         TString fileName = GetFile(kDAQ, "DELAYS", str->GetName());
275         if (fileName.Length()>0){
276           AliInfo(Form("Got the file %s, now we can extract some values.", fileName.Data()));
277
278           daqFile = new TFile(fileName.Data(),"READ");
279           if (fh2) delete fh2;
280           fh2 = (TH2S*) daqFile->Get("htoftot");
281           if (!fh2){
282             Log("some problems occurred:: No histo retrieved, TOF exiting from Shuttle");
283             delete daqFile;
284             return 7; //return error code for histograms not existing/junky
285           }
286           else {
287             static const Int_t kSize=fh2->GetNbinsX();
288             static const Int_t kNBins=fh2->GetNbinsY();
289             static const Double_t kXBinmin=fh2->GetYaxis()->GetBinLowEdge(1);
290             if (kSize != fNChannels){
291               Log(" number of bins along x different from number of pads, found only a subset of the histograms, TOF exiting from Shuttle");
292               delete daqFile;
293               return 7; //return error code for histograms not existing/junky
294             }
295             Int_t nNotStatistics = 0; // number of channel with not enough statistics
296             for (Int_t ich=0;ich<kSize;ich++){
297               TH1S *h1 = new TH1S("h1","h1",kNBins,kXBinmin-0.5,kNBins*1.+kXBinmin-0.5);
298               for (Int_t ibin=0;ibin<kNBins;ibin++){
299                 h1->SetBinContent(ibin+1,fh2->GetBinContent(ich+1,ibin+1));
300               }
301               if(h1->Integral()<fgkIntegralThr) {
302                 nNotStatistics++;
303                 if (!fFDRFlag)  Log(Form(" Not enough statistics for bin %i, skipping this channel",ich));  // printing message only if not in FDR runs
304                 delete h1;
305                 h1=0x0;
306                 continue;
307               }
308               if (!fFDRFlag) {  // not computing delays if in FDR runs
309                 Bool_t found=kFALSE; 
310                 Float_t minContent=h1->Integral()*fgkThrPar; 
311                 Int_t nbinsX = h1->GetNbinsX();
312                 Int_t startBin=1;
313                 for (Int_t j=1; j<=nbinsX; j++){
314                   if ((
315                        h1->GetBinContent(j) +     
316                        h1->GetBinContent(j+1)+
317                        h1->GetBinContent(j+2)+ 
318                        h1->GetBinContent(j+3))>minContent){
319                     found=kTRUE;
320                     startBin=j;
321                     break;
322                   }
323                 }
324                 if(!found) AliInfo(Form("WARNING!!! no start of fit found for histo # %i",ich));
325                 // Now calculate the mean over the interval. 
326                 Double_t mean = 0;
327                 Double_t sumw2 = 0;
328                 Double_t nent = 0;
329                 for(Int_t k=0;k<fgkBinRangeAve;k++){
330                   mean=mean+h1->GetBinCenter(startBin+k)*h1->GetBinContent(startBin+k);                 
331                   nent=nent+h1->GetBinContent(startBin+k);                 
332                   sumw2=sumw2+(h1->GetBinCenter(startBin+k))*(h1->GetBinCenter(startBin+k))*(h1->GetBinContent(startBin+k));
333                 }
334                 mean= mean/nent; //<x>
335                 sumw2=sumw2/nent; //<x^2>
336                 Double_t rmsmean= 0;
337                 rmsmean = TMath::Sqrt((sumw2-mean*mean)/nent);
338                 if (ich<fNChannels) {
339                   AliTOFChannelOnline * ch = (AliTOFChannelOnline *)fCal->At(ich);
340                   ch->SetDelay((Double_t)mean*AliTOFGeometry::TdcBinWidth()*1.E-3);  // delay in ns
341                   // ch->SetStatus(1);  // calibrated channel, removed for the time being from AliTOFChannelOnline
342                 }
343               }
344             delete h1;
345             h1=0x0;
346             }
347             if (nNotStatistics!=0) Log(Form("Too little statistics for %d channels!",nNotStatistics)); 
348           }
349           daqFile->Close();
350           delete daqFile;
351         }
352         else{
353           Log("The Cumulative data file from DAQ does not exist, TOF exiting from Shuttle"); 
354           return 6;//return error code for problems in retrieving DAQ data 
355         }
356       }
357       delete listTot;
358     }
359   else{
360     Log("Problem: no list for Cumulative data file from DAQ was found, TOF exiting from Shuttle");
361     return 6; //return error code for problems in retrieving DAQ data 
362   }
363
364   daqFile=0;
365   AliCDBMetaData metaData;
366   metaData.SetBeamPeriod(0);
367   metaData.SetResponsible("Chiara Zampolli");
368   metaData.SetComment("This preprocessor fills a TObjArray object.");
369   AliInfo("Storing Calibration Data");
370   resultTOFPP = Store("Calib","ParOnline",fCal, &metaData,0,kTRUE);
371   if(!resultTOFPP){
372     Log("Some problems occurred while storing online object resulting from DAQ data processing");
373     return 8;//return error code for problems in storing DAQ data 
374   }
375
376   return 0;
377 }
378 //_____________________________________________________________________________
379
380 UInt_t AliTOFPreprocessor::ProcessPulserData()
381 {
382   // Processing Pulser Run data for TOF channel status
383
384   Log("Processing Pulser");
385
386   TH1::AddDirectory(0);
387
388   Bool_t resultPulserRef=kFALSE;
389   Bool_t resultPulser=kFALSE;
390
391   static const Int_t kSize = AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors();
392   TH1S * htofPulser = new TH1S("hTOFpulser","histo with signals on TOF during pulser", kSize,-0.5,kSize-0.5);
393   for (Int_t ibin =1;ibin<=kSize;ibin++){
394     htofPulser->SetBinContent(ibin,-1);
395   }
396
397   // retrieving last stored pulser object, and copying
398
399   AliCDBEntry *cdbEntry = GetFromOCDB("Calib","Pulser");
400   if (cdbEntry!=0x0){
401     TObjArray *currentCalStatus = (TObjArray*)cdbEntry->GetObject();
402     for (Int_t ich = 0; ich<fNChannels; ich ++){
403       AliTOFChannelOnlineStatus * calChOnlineSt = (AliTOFChannelOnlineStatus*)fCalStatus->At(ich);
404       calChOnlineSt->SetStatus(((AliTOFChannelOnlineStatus*)currentCalStatus->At(ich))->GetStatus());
405     }
406   }
407
408   // processing pulser
409   
410   TFile * daqFile=0x0;
411   TH1S *h1=0x0;
412   
413   //retrieving Pulser data 
414   TList* listPulser = GetFileSources(kDAQ, "PULSER");
415   if (listPulser !=0x0 && listPulser->GetEntries()!=0)
416     {
417       AliInfo("The following sources produced files with the id PULSER");
418       listPulser->Print();
419       for (Int_t jj=0;jj<listPulser->GetEntries();jj++){
420         TObjString * str = dynamic_cast<TObjString*> (listPulser->At(jj));
421         AliInfo(Form("found source %s", str->String().Data()));
422         // file to be stored run per run
423         TString fileNamePulser = GetFile(kDAQ, "PULSER", str->GetName());
424         if (fileNamePulser.Length()>0){
425           // storing refernce data
426           AliInfo(Form("Got the file %s, now we can store the Reference Data for pulser for the current Run.", fileNamePulser.Data()));
427           daqFile = new TFile(fileNamePulser.Data(),"READ");
428           h1 = (TH1S*) daqFile->Get("hTOFpulser");
429           for (Int_t ibin=0;ibin<kSize;ibin++){
430             if ((h1->GetBinContent(ibin+1))!=-1){
431               if ((htofPulser->GetBinContent(ibin+1))==-1){
432                 htofPulser->SetBinContent(ibin+1,h1->GetBinContent(ibin+1));
433               }
434               else {
435                 Log(Form("Something strange occurred during Pulser run, channel %i already read by another LDC, please check!",ibin));
436               }
437             }
438           }
439           // elaborating infos
440           Double_t mean =0;
441           Int_t nread=0;
442           Int_t nreadNotEmpty=0;
443           for (Int_t ientry=1;ientry<=h1->GetNbinsX();ientry++){
444             if (h1->GetBinContent(ientry)==-1) continue;
445             else {
446               if (h1->GetBinContent(ientry)>0) {
447                 nreadNotEmpty++;
448                 AliDebug(1,Form(" channel %i is ok with entry = %f; so far %i channels added ",ientry-1,h1->GetBinContent(ientry),nreadNotEmpty));
449               }
450               mean+=h1->GetBinContent(ientry);
451               nread++;
452             }
453           }
454           mean/=nread;
455           AliDebug(1,Form(" nread =  %i , mean = %f",nread,mean));
456           for (Int_t ich =0;ich<fNChannels;ich++){
457             AliTOFChannelOnlineStatus * chSt = (AliTOFChannelOnlineStatus *)fCalStatus->At(ich);
458             if (h1->GetBinContent(ich+1)==-1) continue;
459             AliDebug(1,Form(" channel %i ",ich));
460             AliDebug(1,Form(" channel status before pulser = %i",(Int_t)chSt->GetStatus()));
461             if (h1->GetBinContent(ich+1)<0.05*mean){
462               chSt->SetStatus(chSt->GetStatus()|AliTOFChannelOnlineStatus::kTOFPulserBad);  // bad status for pulser
463               AliDebug(1,Form(" channel status after pulser = %i",(Int_t)chSt->GetStatus()));
464             }
465             else {
466               chSt->SetStatus(chSt->GetStatus()|AliTOFChannelOnlineStatus::kTOFPulserOk);  // bad status for pulser
467               AliDebug(1,Form(" channel status after pulser = %i",(Int_t)chSt->GetStatus()));
468             }
469           }
470           
471           daqFile->Close();
472           delete daqFile;
473           delete h1;
474         }
475         
476         else{
477           Log("The input data file from DAQ (pulser) was not found, TOF exiting from Shuttle "); 
478           return 10;//return error code for failure in retrieving Ref Data 
479         }
480         
481       }
482       delete listPulser;
483     }
484   
485   else{
486     Log("The input data file list from DAQ (pulser) was not found, TOF exiting from Shuttle "); 
487     return 10;//return error code for failure in retrieving Ref Data 
488   }     
489
490   //storing in OCDB  
491
492   AliCDBMetaData metaData;
493   metaData.SetBeamPeriod(0);
494   metaData.SetResponsible("Chiara Zampolli");
495   metaData.SetComment("This preprocessor fills a TObjArray object for Pulser data.");
496   AliInfo("Storing Calibration Data from Pulser Run");
497   resultPulser = Store("Calib","PulserData",fCalStatus, &metaData,0,kTRUE);
498   if(!resultPulser){
499     Log("Some problems occurred while storing online object resulting from Pulser data processing");
500     return 11;//return error code for problems in storing Pulser data 
501   }
502
503   if(fStoreRefData){
504     
505     AliCDBMetaData metaDataHisto;
506     metaDataHisto.SetBeamPeriod(0);
507     metaDataHisto.SetResponsible("Chiara Zampolli");
508     char comment[200];
509     sprintf(comment,"This preprocessor stores the result of the pulser run");
510     metaDataHisto.SetComment(comment);
511     AliInfo("Storing Reference Data");
512     resultPulserRef = StoreReferenceData("Calib","Pulser",htofPulser, &metaDataHisto);
513     if (!resultPulserRef){
514       Log("some problems occurred::No Reference Data for pulser stored, TOF exiting from Shuttle");
515       return 9;//return error code for failure in storing Ref Data 
516     }
517   }
518   
519   daqFile=0;
520
521   return 0;
522 }
523 //_____________________________________________________________________________
524
525 UInt_t AliTOFPreprocessor::ProcessNoiseData()
526 {
527
528   // Processing Noise Run data for TOF channel status
529
530   Log("Processing Noise");
531
532   TH1::AddDirectory(0);
533
534   Bool_t resultNoiseRef=kFALSE;
535   Bool_t resultNoise=kFALSE;
536
537   static const Int_t kSize = AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors();
538   TH1F * htofNoise = new TH1F("hTOFnoise","histo with signals on TOF during pulser", kSize,-0.5,kSize-0.5);
539   for (Int_t ibin =1;ibin<=kSize;ibin++){
540     htofNoise->SetBinContent(ibin,-1);
541   }
542
543   // retrieving last stored noise object, and copying
544
545   AliCDBEntry *cdbEntry = GetFromOCDB("Calib","Noise");
546   if (cdbEntry!=0x0){
547     TObjArray *currentCalStatus = (TObjArray*)cdbEntry->GetObject();
548     for (Int_t ich = 0; ich<fNChannels; ich ++){
549       AliTOFChannelOnlineStatus * calChOnlineSt = (AliTOFChannelOnlineStatus*)fCalStatus->At(ich);
550       calChOnlineSt->SetStatus(((AliTOFChannelOnlineStatus*)currentCalStatus->At(ich))->GetStatus());
551     }
552   }
553
554   // processing noise
555   
556   TFile * daqFile=0x0;
557   TH1F * h1=0x0;
558   
559   //retrieving Noise data 
560   TList* listNoise = GetFileSources(kDAQ, "NOISE");
561   if (listNoise !=0x0 && listNoise->GetEntries()!=0)
562     {
563       AliInfo("The following sources produced files with the id NOISE");
564       listNoise->Print();
565       for (Int_t jj=0;jj<listNoise->GetEntries();jj++){
566         TObjString * str = dynamic_cast<TObjString*> (listNoise->At(jj));
567         AliInfo(Form("found source %s", str->String().Data()));
568         // file to be stored run per run
569         TString fileNameNoise = GetFile(kDAQ, "NOISE", str->GetName());
570         if (fileNameNoise.Length()>0){
571           // storing refernce data
572           AliInfo(Form("Got the file %s, now we can store the Reference Data for noise for the current Run.", fileNameNoise.Data()));
573           daqFile = new TFile(fileNameNoise.Data(),"READ");
574           h1 = (TH1F*) daqFile->Get("hTOFnoise");
575           for (Int_t ibin=0;ibin<kSize;ibin++){
576             if ((h1->GetBinContent(ibin+1))!=-1){
577               if ((htofNoise->GetBinContent(ibin+1))==-1){
578                 htofNoise->SetBinContent(ibin+1,h1->GetBinContent(ibin+1));
579               }
580               else {
581                 Log(Form("Something strange occurred during Noise run, channel %i already read by another LDC, please check!",ibin));
582               }
583             }
584           }
585           // elaborating infos
586           for (Int_t ich =0;ich<fNChannels;ich++){
587             AliTOFChannelOnlineStatus * chSt = (AliTOFChannelOnlineStatus *)fCalStatus->At(ich);
588             if (h1->GetBinContent(ich+1)==-1) continue;
589             AliDebug(1,Form( " channel %i",ich));
590             AliDebug(1,Form( " channel status before noise = %i",(Int_t)chSt->GetStatus()));
591             if (h1->GetBinContent(ich+1)>=1){  // setting limit for noise to 1 kHz
592               chSt->SetStatus(chSt->GetStatus()|AliTOFChannelOnlineStatus::kTOFNoiseBad);  // bad status for noise
593               AliDebug(1,Form( " channel status after noise = %i",(Int_t)chSt->GetStatus()));
594             }
595             else {
596               chSt->SetStatus(chSt->GetStatus()|AliTOFChannelOnlineStatus::kTOFNoiseOk);  // bad status for noise
597               AliDebug(1,Form(" channel status after noise = %i",(Int_t)chSt->GetStatus()));
598             }
599           }
600  
601           daqFile->Close();
602           delete daqFile;
603           delete h1;
604
605         }
606         
607         else{
608           Log("The input data file from DAQ (noise) was not found, TOF exiting from Shuttle "); 
609           return 13;//return error code for failure in retrieving Ref Data 
610         }
611         
612       }
613       delete listNoise;
614     }
615   else{
616     Log("The input data file list from DAQ (noise) was not found, TOF exiting from Shuttle "); 
617     return 13;//return error code for failure in retrieving Ref Data 
618   }     
619   
620   daqFile=0;
621
622   //storing in OCDB
623
624   AliCDBMetaData metaData;
625   metaData.SetBeamPeriod(0);
626   metaData.SetResponsible("Chiara Zampolli");
627   metaData.SetComment("This preprocessor fills a TObjArray object for Noise data.");
628   AliInfo("Storing Calibration Data from Noise Run");
629   resultNoise = Store("Calib","NoiseData",fCalStatus, &metaData,0,kTRUE);
630   if(!resultNoise){
631     Log("Some problems occurred while storing online object resulting from Noise data processing");
632     return 14;//return error code for problems in storing Noise data 
633   }
634
635   if(fStoreRefData){
636     
637     AliCDBMetaData metaDataHisto;
638     metaDataHisto.SetBeamPeriod(0);
639     metaDataHisto.SetResponsible("Chiara Zampolli");
640     char comment[200];
641     sprintf(comment,"This preprocessor stores the result of the noise run, TOF exiting from Shuttle ");
642     metaDataHisto.SetComment(comment);
643     AliInfo("Storing Reference Data");
644     resultNoiseRef = StoreReferenceData("Calib","Noise",htofNoise, &metaDataHisto);
645     if (!resultNoiseRef){
646       Log("some problems occurred::No Reference Data for noise stored");
647       return 12;//return error code for failure in storing Ref Data 
648     }
649   }
650
651   return 0;
652 }
653 //_____________________________________________________________________________
654
655 UInt_t AliTOFPreprocessor::ProcessHWData()
656 {
657   // Processing Pulser Run data for TOF channel status
658   // dummy for the time being
659
660   Log("Processing HW");
661
662   return 0;
663
664 }
665
666 //_____________________________________________________________________________
667
668 UInt_t AliTOFPreprocessor::Process(TMap* dcsAliasMap)
669 {
670   //
671   //
672   //
673
674   TString runType = GetRunType();
675   Log(Form("RunType %s",runType.Data()));
676   
677   //*((TString*) (0x0)) = "bla";
678
679   // processing 
680
681   if (runType == "PULSER") {
682     Int_t iresultPulser = ProcessPulserData();
683     return iresultPulser; 
684   }
685   
686   if (runType == "NOISE") { // for the time being associating noise runs with pedestal runs; proper run type to be defined 
687     Int_t iresultNoise = ProcessNoiseData();
688     return iresultNoise; 
689   }
690   
691   if (runType == "PHYSICS") {
692     Int_t iresultDCS = ProcessDCSDataPoints(dcsAliasMap);
693     if (iresultDCS != 0) {
694       return iresultDCS;
695     }
696     else { 
697       Int_t iresultDAQ = ProcessOnlineDelays();
698       return iresultDAQ;
699     }
700   }
701
702   // storing
703   return 0;
704 }
705
706