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