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