]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFPreprocessor.cxx
Adding RAWDatarec before STEER (Laurent)
[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   delete fData;
76   fData = 0;
77   delete fh2;
78   fh2 = 0;
79   fCal->Clear();
80   delete fCal;
81   fCal = 0;
82 }
83
84 //______________________________________________________________________________
85 void AliTOFPreprocessor::Initialize(Int_t run, UInt_t startTime,
86         UInt_t endTime)
87 {
88   // Creates AliTOFDataDCS object
89
90   AliPreprocessor::Initialize(run, startTime, endTime);
91
92         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
93                 TTimeStamp(startTime).AsString(),
94                 TTimeStamp(endTime).AsString()));
95
96         fData = new AliTOFDataDCS(fRun, fStartTime, fEndTime);
97         fh2 = 0x0;
98         fNChannels = AliTOFGeometry::NSectors()*(2*(AliTOFGeometry::NStripC()+AliTOFGeometry::NStripB())+AliTOFGeometry::NStripA())*AliTOFGeometry::NpadZ()*AliTOFGeometry::NpadX();
99         fCal = new TObjArray(fNChannels);
100         fCal->SetOwner();
101         for (Int_t ich = 0; ich<fNChannels; ich ++){
102           AliTOFChannelOnline * calChOnline = new AliTOFChannelOnline();
103           fCal->AddAt(calChOnline,ich);
104         }
105 }
106
107 //_____________________________________________________________________________
108
109 UInt_t AliTOFPreprocessor::Process(TMap* dcsAliasMap)
110 {
111   // Fills data into a AliTOFDataDCS object
112   // return codes:
113   // return=0 : all ok
114   // return=1 : no DCS input data Map
115   // return=2 : no DCS input data processing
116   // return=3 : no DCS processed data was stored
117   // return=4 : no DAQ input for Ref Data
118   // return=5 : failed to store Ref data
119   // return=6 : failed to retrieve DAQ data for calibration 
120   // return=7 : problems in histos in the input DAQ file 
121   // return=8 : failed to store Online Delays
122
123   TH1::AddDirectory(0);
124
125   Bool_t resultDCSMap=kFALSE;
126   Bool_t resultDCSStore=kFALSE;
127   Bool_t resultDAQ=kFALSE;
128   Bool_t resultDAQRef=kFALSE;
129
130   // processing DCS
131
132   if (!dcsAliasMap){
133     Log("No DCS map found: TOF exiting from Shuttle");
134     return 1;// return error Code for DCS input data not found 
135   }
136   else {
137   // The processing of the DCS input data is forwarded to AliTOFDataDCS
138     resultDCSMap=fData->ProcessData(*dcsAliasMap);
139     if(!resultDCSMap){
140       Log("Some problems occurred while processing DCS data, TOF exiting from Shuttle");
141       return 2;// return error Code for processed DCS data not stored 
142     }
143     else{
144       AliCDBMetaData metaDataDCS;
145       metaDataDCS.SetBeamPeriod(0);
146       metaDataDCS.SetResponsible("Chiara Zampolli");
147       metaDataDCS.SetComment("This preprocessor fills an AliTOFDataDCS object.");
148       AliInfo("Storing DCS Data");
149       resultDCSStore = Store("Calib","DCSData",fData, &metaDataDCS);
150       if (!resultDCSStore){
151         Log("Some problems occurred while storing DCS data results in OCDB, TOF exiting from Shuttle");
152         return 3;// return error Code for processed DCS data not stored 
153         
154       }
155     }
156   }
157   
158   // processing DAQ
159   
160   TFile * daqFile=0x0;
161   
162   if(fStoreRefData){
163     //retrieving data at Run level
164     TList* list = GetFileSources(kDAQ, "RUNLevel");
165     if (list)
166       {
167         AliInfo("The following sources produced files with the id RUNLevel");
168         list->Print();
169         for (Int_t jj=0;jj<list->GetEntries();jj++){
170           TObjString * str = dynamic_cast<TObjString*> (list->At(jj));
171           AliInfo(Form("found source %s", str->String().Data()));
172           // file to be stored run per run
173           TString fileNameRun = GetFile(kDAQ, "RUNLevel", str->GetName());
174           if (fileNameRun.Length()>0){
175             AliInfo(Form("Got the file %s, now we can store the Reference Data for the current Run.", fileNameRun.Data()));
176             daqFile = new TFile(fileNameRun.Data(),"READ");
177             fh2 = (TH2S*) daqFile->Get("htof");
178             AliCDBMetaData metaDataHisto;
179             metaDataHisto.SetBeamPeriod(0);
180             metaDataHisto.SetResponsible("Chiara Zampolli");
181             metaDataHisto.SetComment("This preprocessor stores the array of histos object as Reference Data.");
182             AliInfo("Storing Reference Data");
183             resultDAQRef = StoreReferenceData("Calib","DAQData",fh2, &metaDataHisto);
184             if (!resultDAQRef){
185               Log("some problems occurred::No Reference Data stored, TOF exiting from Shuttle");
186               return 5;//return error code for failure in storing Ref Data 
187             }
188             daqFile->Close();
189             delete daqFile;
190           }
191           
192           else{
193             Log("The input data file from DAQ (run-level) was not found, TOF exiting from Shuttle "); 
194             return 4;//return error code for failure in retrieving Ref Data 
195           }
196         }
197       }
198     else{
199       Log("The input data file list from DAQ (run-level) was not found, TOF exiting from Shuttle "); 
200       return 4;//return error code for failure in retrieving Ref Data 
201     }   
202   }
203
204
205 //Total files, with cumulative histos
206   
207   TList* listTot = GetFileSources(kDAQ, "DELAYS");
208   if (listTot)
209     {
210       AliInfo("The following sources produced files with the id DELAYS");
211       listTot->Print();
212       for (Int_t jj=0;jj<listTot->GetEntries();jj++){
213         TObjString * str = dynamic_cast<TObjString*> (listTot->At(jj));
214         AliInfo(Form("found source %s", str->String().Data()));
215
216         // file with summed histos, to extract calib params
217         TString fileName = GetFile(kDAQ, "DELAYS", str->GetName());
218         if (fileName.Length()>0){
219           AliInfo(Form("Got the file %s, now we can extract some values.", fileName.Data()));
220
221           daqFile = new TFile(fileName.Data(),"READ");
222           fh2 = (TH2S*) daqFile->Get("htoftot");
223           if (!fh2){
224             Log("some problems occurred:: No histo retrieved, TOF exiting from Shuttle");
225             delete daqFile;
226             return 7; //return error code for histograms not existing/junky
227           }
228           else {
229             static const Int_t kSize=fh2->GetNbinsX();
230             static const Int_t kNBins=fh2->GetNbinsY();
231             static const Double_t kXBinmin=fh2->GetYaxis()->GetBinLowEdge(1);
232             if (kSize != fNChannels){
233               Log(" number of bins along x different from number of pads, found only a subset of the histograms, TOF exiting from Shuttle");
234               delete daqFile;
235               return 7; //return error code for histograms not existing/junky
236             }
237             for (Int_t ich=0;ich<kSize;ich++){
238               TH1S *h1 = new TH1S("h1","h1",kNBins,kXBinmin-0.5,kNBins*1.+kXBinmin-0.5);
239               for (Int_t ibin=0;ibin<kNBins;ibin++){
240                 h1->SetBinContent(ibin+1,fh2->GetBinContent(ich+1,ibin+1));
241               }
242               Bool_t found=kFALSE; 
243               Float_t minContent=h1->Integral()*fgkThrPar; 
244               Int_t nbinsX = h1->GetNbinsX();
245               Int_t startBin=1;
246               for (Int_t j=1; j<=nbinsX; j++){
247                 if ((
248                      h1->GetBinContent(j) +     
249                      h1->GetBinContent(j+1)+
250                      h1->GetBinContent(j+2)+ 
251                      h1->GetBinContent(j+3))>minContent){
252                   found=kTRUE;
253                   startBin=j;
254                   break;
255                 }
256               }
257               if(!found) AliInfo(Form("WARNING!!! no start of fit found for histo # %i",ich));
258               // Now calculate the mean over the interval. 
259               Double_t mean = 0;
260               Double_t sumw2 = 0;
261               Double_t nent = 0;
262               for(Int_t k=0;k<fgkBinRangeAve;k++){
263                 mean=mean+h1->GetBinCenter(startBin+k)*h1->GetBinContent(startBin+k);                 
264                 nent=nent+h1->GetBinContent(startBin+k);                 
265                 sumw2=sumw2+(h1->GetBinCenter(startBin+k))*(h1->GetBinCenter(startBin+k))*(h1->GetBinContent(startBin+k));
266               }
267               mean= mean/nent; //<x>
268               sumw2=sumw2/nent; //<x^2>
269               Double_t rmsmean= 0;
270               rmsmean = TMath::Sqrt((sumw2-mean*mean)/nent);
271               if (ich<fNChannels) {
272                 AliTOFChannelOnline * ch = (AliTOFChannelOnline *)fCal->At(ich);
273                 ch->SetDelay((Double_t)mean*AliTOFGeometry::TdcBinWidth()*1.E-3);  // delay in ns
274               }
275             delete h1;
276             h1=0x0;
277             }
278           }
279           daqFile->Close();
280           delete daqFile;
281           AliCDBMetaData metaData;
282           metaData.SetBeamPeriod(0);
283           metaData.SetResponsible("Chiara Zampolli");
284           metaData.SetComment("This preprocessor fills an AliTOFCal object.");
285           AliInfo("Storing Calibration Data");
286           resultDAQ = Store("Calib","OnlineDelay",fCal, &metaData);
287           if(!resultDAQ){
288             Log("Some problems occurred while storing DAQ data processing results");
289             return 8;//return error code for problems in storing DAQ data 
290           }
291         }
292         else{
293           Log("The Cumulative data file from DAQ does not exist, TOF exiting from Shuttle"); 
294           return 6;//return error code for problems in retrieving DAQ data 
295         }
296       }
297     }
298   else{
299     Log("Problem: no list for Cumulative data file from DAQ was found, TOF exiting from Shuttle");
300     return 6; //return error code for problems in retrieving DAQ data 
301   }
302
303   daqFile=0;
304
305   return 0;
306 }
307
308