]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFPreprocessor.cxx
Corrected a bug in kalman tracking (final parameters and covariances
[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::Process(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   Bool_t resultDAQ=kFALSE;
134   Bool_t resultDAQRef=kFALSE;
135
136   // processing DCS
137
138   if (!dcsAliasMap){
139     Log("No DCS map found: TOF exiting from Shuttle");
140     return 1;// return error Code for DCS input data not found 
141   }
142   else {
143   // The processing of the DCS input data is forwarded to AliTOFDataDCS
144     resultDCSMap=fData->ProcessData(*dcsAliasMap);
145     if(!resultDCSMap){
146       Log("Some problems occurred while processing DCS data, TOF exiting from Shuttle");
147       return 2;// return error Code for processed DCS data not stored 
148     }
149     else{
150       AliCDBMetaData metaDataDCS;
151       metaDataDCS.SetBeamPeriod(0);
152       metaDataDCS.SetResponsible("Chiara Zampolli");
153       metaDataDCS.SetComment("This preprocessor fills an AliTOFDataDCS object.");
154       AliInfo("Storing DCS Data");
155       resultDCSStore = StoreReferenceData("Calib","DCSData",fData, &metaDataDCS);
156       if (!resultDCSStore){
157         Log("Some problems occurred while storing DCS data results in Reference Data, TOF exiting from Shuttle");
158         return 3;// return error Code for processed DCS data not stored 
159                  // in reference data
160         
161       }
162     }
163   }
164   
165   // processing DAQ
166   
167   TFile * daqFile=0x0;
168   
169   if(fStoreRefData){
170     //retrieving data at Run level
171     TList* list = GetFileSources(kDAQ, "RUNLevel");
172     if (list)
173       {
174         AliInfo("The following sources produced files with the id RUNLevel");
175         list->Print();
176         for (Int_t jj=0;jj<list->GetEntries();jj++){
177           TObjString * str = dynamic_cast<TObjString*> (list->At(jj));
178           AliInfo(Form("found source %s", str->String().Data()));
179           // file to be stored run per run
180           TString fileNameRun = GetFile(kDAQ, "RUNLevel", str->GetName());
181           if (fileNameRun.Length()>0){
182             AliInfo(Form("Got the file %s, now we can store the Reference Data for the current Run.", fileNameRun.Data()));
183             daqFile = new TFile(fileNameRun.Data(),"READ");
184             fh2 = (TH2S*) daqFile->Get("htof");
185             AliCDBMetaData metaDataHisto;
186             metaDataHisto.SetBeamPeriod(0);
187             metaDataHisto.SetResponsible("Chiara Zampolli");
188             metaDataHisto.SetComment("This preprocessor stores the array of histos object as Reference Data.");
189             AliInfo("Storing Reference Data");
190             resultDAQRef = StoreReferenceData("Calib","DAQData",fh2, &metaDataHisto);
191             if (!resultDAQRef){
192               Log("some problems occurred::No Reference Data stored, TOF exiting from Shuttle");
193               return 5;//return error code for failure in storing Ref Data 
194             }
195             daqFile->Close();
196             delete daqFile;
197           }
198           
199           else{
200             Log("The input data file from DAQ (run-level) was not found, TOF exiting from Shuttle "); 
201             return 4;//return error code for failure in retrieving Ref Data 
202           }
203         }
204       }
205     else{
206       Log("The input data file list from DAQ (run-level) was not found, TOF exiting from Shuttle "); 
207       return 4;//return error code for failure in retrieving Ref Data 
208     }   
209   }
210
211
212 //Total files, with cumulative histos
213   
214   TList* listTot = GetFileSources(kDAQ, "DELAYS");
215   if (listTot)
216     {
217       AliInfo("The following sources produced files with the id DELAYS");
218       listTot->Print();
219       for (Int_t jj=0;jj<listTot->GetEntries();jj++){
220         TObjString * str = dynamic_cast<TObjString*> (listTot->At(jj));
221         AliInfo(Form("found source %s", str->String().Data()));
222
223         // file with summed histos, to extract calib params
224         TString fileName = GetFile(kDAQ, "DELAYS", str->GetName());
225         if (fileName.Length()>0){
226           AliInfo(Form("Got the file %s, now we can extract some values.", fileName.Data()));
227
228           daqFile = new TFile(fileName.Data(),"READ");
229           if (fh2) delete fh2;
230           fh2 = (TH2S*) daqFile->Get("htoftot");
231           if (!fh2){
232             Log("some problems occurred:: No histo retrieved, TOF exiting from Shuttle");
233             delete daqFile;
234             return 7; //return error code for histograms not existing/junky
235           }
236           else {
237             static const Int_t kSize=fh2->GetNbinsX();
238             static const Int_t kNBins=fh2->GetNbinsY();
239             static const Double_t kXBinmin=fh2->GetYaxis()->GetBinLowEdge(1);
240             if (kSize != fNChannels){
241               Log(" number of bins along x different from number of pads, found only a subset of the histograms, TOF exiting from Shuttle");
242               delete daqFile;
243               return 7; //return error code for histograms not existing/junky
244             }
245             for (Int_t ich=0;ich<kSize;ich++){
246               TH1S *h1 = new TH1S("h1","h1",kNBins,kXBinmin-0.5,kNBins*1.+kXBinmin-0.5);
247               for (Int_t ibin=0;ibin<kNBins;ibin++){
248                 h1->SetBinContent(ibin+1,fh2->GetBinContent(ich+1,ibin+1));
249               }
250               Bool_t found=kFALSE; 
251               Float_t minContent=h1->Integral()*fgkThrPar; 
252               Int_t nbinsX = h1->GetNbinsX();
253               Int_t startBin=1;
254               for (Int_t j=1; j<=nbinsX; j++){
255                 if ((
256                      h1->GetBinContent(j) +     
257                      h1->GetBinContent(j+1)+
258                      h1->GetBinContent(j+2)+ 
259                      h1->GetBinContent(j+3))>minContent){
260                   found=kTRUE;
261                   startBin=j;
262                   break;
263                 }
264               }
265               if(!found) AliInfo(Form("WARNING!!! no start of fit found for histo # %i",ich));
266               // Now calculate the mean over the interval. 
267               Double_t mean = 0;
268               Double_t sumw2 = 0;
269               Double_t nent = 0;
270               for(Int_t k=0;k<fgkBinRangeAve;k++){
271                 mean=mean+h1->GetBinCenter(startBin+k)*h1->GetBinContent(startBin+k);                 
272                 nent=nent+h1->GetBinContent(startBin+k);                 
273                 sumw2=sumw2+(h1->GetBinCenter(startBin+k))*(h1->GetBinCenter(startBin+k))*(h1->GetBinContent(startBin+k));
274               }
275               mean= mean/nent; //<x>
276               sumw2=sumw2/nent; //<x^2>
277               Double_t rmsmean= 0;
278               rmsmean = TMath::Sqrt((sumw2-mean*mean)/nent);
279               if (ich<fNChannels) {
280                 AliTOFChannelOnline * ch = (AliTOFChannelOnline *)fCal->At(ich);
281                 ch->SetDelay((Double_t)mean*AliTOFGeometry::TdcBinWidth()*1.E-3);  // delay in ns
282               }
283             delete h1;
284             h1=0x0;
285             }
286           }
287           daqFile->Close();
288           delete daqFile;
289           AliCDBMetaData metaData;
290           metaData.SetBeamPeriod(0);
291           metaData.SetResponsible("Chiara Zampolli");
292           metaData.SetComment("This preprocessor fills an AliTOFCal object.");
293           AliInfo("Storing Calibration Data");
294           resultDAQ = Store("Calib","ParOnline",fCal, &metaData);
295           if(!resultDAQ){
296             Log("Some problems occurred while storing DAQ data processing results");
297             return 8;//return error code for problems in storing DAQ data 
298           }
299         }
300         else{
301           Log("The Cumulative data file from DAQ does not exist, TOF exiting from Shuttle"); 
302           return 6;//return error code for problems in retrieving DAQ data 
303         }
304       }
305     }
306   else{
307     Log("Problem: no list for Cumulative data file from DAQ was found, TOF exiting from Shuttle");
308     return 6; //return error code for problems in retrieving DAQ data 
309   }
310
311   daqFile=0;
312
313   return 0;
314 }
315
316