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