]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFPreprocessor.cxx
Bug fix by Alex
[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
9edefa04 16/* $Id$ */
9d883ed9 17
9edefa04 18#include <Riostream.h>
19#include <stdio.h>
20#include <stdlib.h>
c9fe8530 21
9edefa04 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>
c9fe8530 31
32#include "AliCDBMetaData.h"
33#include "AliLog.h"
c9fe8530 34#include "AliTOFChannelOnline.h"
9edefa04 35#include "AliTOFDataDCS.h"
10056aa6 36#include "AliTOFGeometry.h"
9edefa04 37#include "AliTOFPreprocessor.h"
c9fe8530 38
39class TF1;
40class AliDCSValue;
41class 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
51ClassImp(AliTOFPreprocessor)
52
53const Int_t AliTOFPreprocessor::fgkBinRangeAve = 13; // number of bins where to calculate the mean
9d883ed9 54const Double_t AliTOFPreprocessor::fgkThrPar = 0.013; // parameter used to trigger the calculation of the delay
c9fe8530 55
56//_____________________________________________________________________________
57
708db10b 58AliTOFPreprocessor::AliTOFPreprocessor(AliShuttleInterface* shuttle) :
59 AliPreprocessor("TOF", shuttle),
c9fe8530 60 fData(0),
708db10b 61 fh2(0),
c9fe8530 62 fCal(0),
10056aa6 63 fNChannels(0),
5936ab02 64 fStoreRefData(kTRUE)
c9fe8530 65{
66 // constructor
10056aa6 67
c9fe8530 68}
69
70//_____________________________________________________________________________
71
72AliTOFPreprocessor::~AliTOFPreprocessor()
73{
74 // destructor
9d883ed9 75 delete fData;
76 fData = 0;
77 delete fh2;
78 fh2 = 0;
10056aa6 79 fCal->Clear();
9d883ed9 80 delete fCal;
81 fCal = 0;
c9fe8530 82}
83
84//______________________________________________________________________________
85void 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);
708db10b 97 fh2 = 0x0;
10056aa6 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 }
c9fe8530 105}
106
107//_____________________________________________________________________________
108
109UInt_t AliTOFPreprocessor::Process(TMap* dcsAliasMap)
110{
111 // Fills data into a AliTOFDataDCS object
5936ab02 112 // return codes:
113 // return=0 : all ok
7885c291 114 // return=1 : no DCS input data Map
115 // return=2 : no DCS input data processing
b0d8f3fc 116 // return=3 : no DCS processed data was stored in Ref Data
7885c291 117 // return=4 : no DAQ input for Ref Data
b0d8f3fc 118 // return=5 : failed to store Ref Data
7885c291 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
c9fe8530 122
708db10b 123 TH1::AddDirectory(0);
5936ab02 124
7885c291 125 Bool_t resultDCSMap=kFALSE;
126 Bool_t resultDCSStore=kFALSE;
127 Bool_t resultDAQ=kFALSE;
128 Bool_t resultDAQRef=kFALSE;
c9fe8530 129
130 // processing DCS
131
132 if (!dcsAliasMap){
5936ab02 133 Log("No DCS map found: TOF exiting from Shuttle");
134 return 1;// return error Code for DCS input data not found
c9fe8530 135 }
136 else {
137 // The processing of the DCS input data is forwarded to AliTOFDataDCS
7885c291 138 resultDCSMap=fData->ProcessData(*dcsAliasMap);
139 if(!resultDCSMap){
140 Log("Some problems occurred while processing DCS data, TOF exiting from Shuttle");
5936ab02 141 return 2;// return error Code for processed DCS data not stored
c9fe8530 142 }
7885c291 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");
b0d8f3fc 149 resultDCSStore = StoreReferenceData("Calib","DCSData",fData, &metaDataDCS);
7885c291 150 if (!resultDCSStore){
b0d8f3fc 151 Log("Some problems occurred while storing DCS data results in Reference Data, TOF exiting from Shuttle");
7885c291 152 return 3;// return error Code for processed DCS data not stored
b0d8f3fc 153 // in reference data
7885c291 154
155 }
156 }
c9fe8530 157 }
7885c291 158
c9fe8530 159 // processing DAQ
7885c291 160
c9fe8530 161 TFile * daqFile=0x0;
7885c291 162
5936ab02 163 if(fStoreRefData){
164 //retrieving data at Run level
165 TList* list = GetFileSources(kDAQ, "RUNLevel");
166 if (list)
167 {
168 AliInfo("The following sources produced files with the id RUNLevel");
169 list->Print();
170 for (Int_t jj=0;jj<list->GetEntries();jj++){
171 TObjString * str = dynamic_cast<TObjString*> (list->At(jj));
172 AliInfo(Form("found source %s", str->String().Data()));
173 // file to be stored run per run
174 TString fileNameRun = GetFile(kDAQ, "RUNLevel", str->GetName());
175 if (fileNameRun.Length()>0){
176 AliInfo(Form("Got the file %s, now we can store the Reference Data for the current Run.", fileNameRun.Data()));
177 daqFile = new TFile(fileNameRun.Data(),"READ");
178 fh2 = (TH2S*) daqFile->Get("htof");
179 AliCDBMetaData metaDataHisto;
180 metaDataHisto.SetBeamPeriod(0);
181 metaDataHisto.SetResponsible("Chiara Zampolli");
182 metaDataHisto.SetComment("This preprocessor stores the array of histos object as Reference Data.");
183 AliInfo("Storing Reference Data");
184 resultDAQRef = StoreReferenceData("Calib","DAQData",fh2, &metaDataHisto);
185 if (!resultDAQRef){
186 Log("some problems occurred::No Reference Data stored, TOF exiting from Shuttle");
7885c291 187 return 5;//return error code for failure in storing Ref Data
5936ab02 188 }
189 daqFile->Close();
190 delete daqFile;
191 }
192
193 else{
194 Log("The input data file from DAQ (run-level) was not found, TOF exiting from Shuttle ");
7885c291 195 return 4;//return error code for failure in retrieving Ref Data
c9fe8530 196 }
c9fe8530 197 }
9bc469d1 198 }
5936ab02 199 else{
200 Log("The input data file list from DAQ (run-level) was not found, TOF exiting from Shuttle ");
7885c291 201 return 4;//return error code for failure in retrieving Ref Data
5936ab02 202 }
9bc469d1 203 }
204
9bc469d1 205
5936ab02 206//Total files, with cumulative histos
207
9bc469d1 208 TList* listTot = GetFileSources(kDAQ, "DELAYS");
209 if (listTot)
210 {
211 AliInfo("The following sources produced files with the id DELAYS");
212 listTot->Print();
213 for (Int_t jj=0;jj<listTot->GetEntries();jj++){
214 TObjString * str = dynamic_cast<TObjString*> (listTot->At(jj));
215 AliInfo(Form("found source %s", str->String().Data()));
216
c9fe8530 217 // file with summed histos, to extract calib params
9bc469d1 218 TString fileName = GetFile(kDAQ, "DELAYS", str->GetName());
219 if (fileName.Length()>0){
220 AliInfo(Form("Got the file %s, now we can extract some values.", fileName.Data()));
c9fe8530 221
9bc469d1 222 daqFile = new TFile(fileName.Data(),"READ");
708db10b 223 fh2 = (TH2S*) daqFile->Get("htoftot");
224 if (!fh2){
9d883ed9 225 Log("some problems occurred:: No histo retrieved, TOF exiting from Shuttle");
226 delete daqFile;
7885c291 227 return 7; //return error code for histograms not existing/junky
c9fe8530 228 }
c9fe8530 229 else {
3a3ece53 230 static const Int_t kSize=fh2->GetNbinsX();
231 static const Int_t kNBins=fh2->GetNbinsY();
232 static const Double_t kXBinmin=fh2->GetYaxis()->GetBinLowEdge(1);
10056aa6 233 if (kSize != fNChannels){
9d883ed9 234 Log(" number of bins along x different from number of pads, found only a subset of the histograms, TOF exiting from Shuttle");
235 delete daqFile;
7885c291 236 return 7; //return error code for histograms not existing/junky
9bc469d1 237 }
3a3ece53 238 for (Int_t ich=0;ich<kSize;ich++){
239 TH1S *h1 = new TH1S("h1","h1",kNBins,kXBinmin-0.5,kNBins*1.+kXBinmin-0.5);
240 for (Int_t ibin=0;ibin<kNBins;ibin++){
708db10b 241 h1->SetBinContent(ibin+1,fh2->GetBinContent(ich+1,ibin+1));
242 }
c9fe8530 243 Bool_t found=kFALSE;
9bc469d1 244 Float_t minContent=h1->Integral()*fgkThrPar;
c9fe8530 245 Int_t nbinsX = h1->GetNbinsX();
246 Int_t startBin=1;
247 for (Int_t j=1; j<=nbinsX; j++){
248 if ((
249 h1->GetBinContent(j) +
250 h1->GetBinContent(j+1)+
251 h1->GetBinContent(j+2)+
252 h1->GetBinContent(j+3))>minContent){
253 found=kTRUE;
254 startBin=j;
255 break;
256 }
257 }
708db10b 258 if(!found) AliInfo(Form("WARNING!!! no start of fit found for histo # %i",ich));
c9fe8530 259 // Now calculate the mean over the interval.
260 Double_t mean = 0;
261 Double_t sumw2 = 0;
262 Double_t nent = 0;
263 for(Int_t k=0;k<fgkBinRangeAve;k++){
264 mean=mean+h1->GetBinCenter(startBin+k)*h1->GetBinContent(startBin+k);
265 nent=nent+h1->GetBinContent(startBin+k);
266 sumw2=sumw2+(h1->GetBinCenter(startBin+k))*(h1->GetBinCenter(startBin+k))*(h1->GetBinContent(startBin+k));
267 }
c9fe8530 268 mean= mean/nent; //<x>
269 sumw2=sumw2/nent; //<x^2>
270 Double_t rmsmean= 0;
271 rmsmean = TMath::Sqrt((sumw2-mean*mean)/nent);
10056aa6 272 if (ich<fNChannels) {
273 AliTOFChannelOnline * ch = (AliTOFChannelOnline *)fCal->At(ich);
274 ch->SetDelay((Double_t)mean*AliTOFGeometry::TdcBinWidth()*1.E-3); // delay in ns
c9fe8530 275 }
708db10b 276 delete h1;
277 h1=0x0;
c9fe8530 278 }
c9fe8530 279 }
280 daqFile->Close();
281 delete daqFile;
282 AliCDBMetaData metaData;
283 metaData.SetBeamPeriod(0);
284 metaData.SetResponsible("Chiara Zampolli");
285 metaData.SetComment("This preprocessor fills an AliTOFCal object.");
5936ab02 286 AliInfo("Storing Calibration Data");
a7e9db2a 287 resultDAQ = Store("Calib","ParOnline",fCal, &metaData);
9d883ed9 288 if(!resultDAQ){
289 Log("Some problems occurred while storing DAQ data processing results");
7885c291 290 return 8;//return error code for problems in storing DAQ data
9d883ed9 291 }
c9fe8530 292 }
293 else{
9d883ed9 294 Log("The Cumulative data file from DAQ does not exist, TOF exiting from Shuttle");
7885c291 295 return 6;//return error code for problems in retrieving DAQ data
c9fe8530 296 }
297 }
298 }
299 else{
9d883ed9 300 Log("Problem: no list for Cumulative data file from DAQ was found, TOF exiting from Shuttle");
7885c291 301 return 6; //return error code for problems in retrieving DAQ data
c9fe8530 302 }
303
c9fe8530 304 daqFile=0;
9d883ed9 305
5936ab02 306 return 0;
c9fe8530 307}
308
309