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