]> git.uio.no Git - u/mrichter/AliRoot.git/blame - T0/AliT0QADataMakerRec.cxx
Include informations from rstate module
[u/mrichter/AliRoot.git] / T0 / AliT0QADataMakerRec.cxx
CommitLineData
04236e67 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/* $Id$ */
18
19//---
20// Produces the data needed to calculate the quality assurance.
1d7681da 21// Alla.Maevskaya@cern.ch
04236e67 22//---
23
24// --- ROOT system ---
25#include <TClonesArray.h>
26#include <TFile.h>
27#include <TH1F.h>
1d7681da 28#include <TH2F.h>
04236e67 29#include <TDirectory.h>
30// --- Standard library ---
31
32// --- AliRoot header files ---
33#include "AliESDEvent.h"
34#include "AliLog.h"
35#include "AliT0digit.h"
36#include "AliT0hit.h"
37#include "AliT0RecPoint.h"
38#include "AliT0QADataMakerRec.h"
39#include "AliQAChecker.h"
40#include "AliT0RawReader.h"
41
1d7681da 42#include "Riostream.h"
04236e67 43ClassImp(AliT0QADataMakerRec)
44
45//____________________________________________________________________________
46 AliT0QADataMakerRec::AliT0QADataMakerRec() :
47 AliQADataMakerRec(AliQA::GetDetName(AliQA::kT0), "T0 Quality Assurance Data Maker")
48
49{
50 // ctor
1d7681da 51 for (Int_t i=0; i<6; i++) {
52 fNumTriggers[i]=0;
53 fNumTriggersCal[i]=0;
04236e67 54 }
04236e67 55}
56
57//____________________________________________________________________________
58AliT0QADataMakerRec::AliT0QADataMakerRec(const AliT0QADataMakerRec& qadm) :
59 AliQADataMakerRec()
60{
61 //copy ctor
1d7681da 62 SetName((const char*)qadm.GetName()) ;
04236e67 63 SetTitle((const char*)qadm.GetTitle());
64}
65
66//__________________________________________________________________
67AliT0QADataMakerRec& AliT0QADataMakerRec::operator = (const AliT0QADataMakerRec& qadm )
68{
69 // Equal operator.
70 this->~AliT0QADataMakerRec();
71 new(this) AliT0QADataMakerRec(qadm);
72 return *this;
73}
74//____________________________________________________________________________
92a357bf 75void AliT0QADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
04236e67 76{
77 //Detector specific actions at end of cycle
78 // do the QA checking
c4b97145 79 AliQAChecker::Instance()->Run(AliQA::kT0, task, list) ;
1d7681da 80 Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
81 for (Int_t itr=0; itr<6; itr++) {
82 GetRawsData(197)->Fill(triggers[itr],fNumTriggersCal[itr]);
83 GetRawsData(197)->SetBinContent(itr+1, fNumTriggersCal[itr]);
84 }
85 GetRawsData(205)->SetOption("COLZ");
86 GetRawsData(206)->SetOption("COLZ");
87 GetRawsData(207)->SetOption("COLZ");
88
04236e67 89}
90
91//____________________________________________________________________________
92void AliT0QADataMakerRec::StartOfDetectorCycle()
93{
94 //Detector specific actions at start of cycle
95
96}
97
98//____________________________________________________________________________
99void AliT0QADataMakerRec::InitRaws()
100{
101 // create Raw histograms in Raw subdir
1d7681da 102 TString timename, ampname, qtcname, ledname;
103 TString timeCalname, ampCalname, ledCalname, qtcCalname;
04236e67 104
86c020af 105 TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point", 10,1252170, 1252180);
106 Add2RawsList( fhRefPoint,0);
107
1d7681da 108 TH1F *fhRawCFD[24]; TH1F * fhRawLEDamp[24];
109 TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
110 TH1F *fhRawCFDcal[24]; TH1F * fhRawLEDampcal[24];
111 TH1F *fhRawQTCcal[24]; TH1F * fhRawLEDcal[24];
86c020af 112
04236e67 113 for (Int_t i=0; i<24; i++)
114 {
115 timename ="hRawCFD";
1d7681da 116 ledname = "hRawLED";
04236e67 117 qtcname = "hRawQTC";
1d7681da 118 ampname = "hRawLEDminCFD";
04236e67 119 timename += i;
120 ampname += i;
121 qtcname += i;
1d7681da 122 ledname += i;
123 fhRawCFD[i] = new TH1F(timename.Data(), timename.Data(),500,-250,250);
86c020af 124 Add2RawsList( fhRawCFD[i],i+1);
1d7681da 125 fhRawLED[i] = new TH1F(ledname.Data(), ledname.Data(),2000,-1000,1000);
126 Add2RawsList( fhRawLED[i],i+24+1);
446d6ec4 127 fhRawLEDamp[i] = new TH1F(ampname.Data(), ampname.Data(),100,300,600);
1d7681da 128 Add2RawsList( fhRawLEDamp[i],i+48+1);
86c020af 129 fhRawQTC[i] = new TH1F(qtcname.Data(), qtcname.Data(),1500,1000,7000);
1d7681da 130 Add2RawsList( fhRawQTC[i],i+72+1);
04236e67 131 }
cf66e144 132 const Bool_t saveForCorr = kTRUE;
1d7681da 133 TH1F* fhRawTrigger = new TH1F("hRawTrigger"," phys triggers",5,0,5);
134 Add2RawsList(fhRawTrigger ,97);
135
86c020af 136 TH1F* fhRawMean = new TH1F("hRawMean","online mean signal", 100,2400,2500);
1d7681da 137 Add2RawsList( fhRawMean,98);
86c020af 138 TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal", 100,0,600);
1d7681da 139 Add2RawsList( fhRawVertex,99);
86c020af 140 TH1F* fhRawORA = new TH1F("hRawORA","online OR A", 100,2500,2800);
1d7681da 141 Add2RawsList( fhRawORA,100);
86c020af 142 TH1F* fhRawORC = new TH1F("hRawORC","online OR C", 100,2000,2300);
1d7681da 143 Add2RawsList( fhRawORC,101);
144
5ed41460 145 for (Int_t i=0; i<24; i++)
146 {
147 // for events with trigger CALIBRATION_EVENT
148 timeCalname ="hRawCFDcal";
1d7681da 149 ledCalname = "hRawLEDcal";
150 ampCalname = "hRawLEDminCFDcal";
5ed41460 151 qtcCalname = "hRawQTCcal";
152 timeCalname += i;
1d7681da 153 ledCalname += i;
5ed41460 154 ampCalname += i;
155 qtcCalname += i;
1d7681da 156 fhRawCFDcal[i] = new TH1F(timeCalname.Data(), timeCalname.Data(),2000,-1000,1000);
157 Add2RawsList( fhRawCFDcal[i],101+i+1);
158 fhRawLEDcal[i] = new TH1F(ledCalname.Data(), ledCalname.Data(),2000,-1000,1000);
159 Add2RawsList( fhRawLEDcal[i],101+i+24+1);
160 fhRawLEDampcal[i] = new TH1F(ampCalname.Data(), ampCalname.Data(),300,300,600);
161 Add2RawsList( fhRawLEDampcal[i],101+i+48+1);
5ed41460 162 fhRawQTCcal[i] = new TH1F(qtcCalname.Data(), qtcCalname.Data(),1000,0,7000);
1d7681da 163 Add2RawsList( fhRawQTCcal[i],101+i+72+1);
5ed41460 164 }
1d7681da 165
166 TH1F* fhRawTriggerCal = new TH1F("hRawTriggerCal"," laser triggers",6,0,6);
c4b97145 167 Add2RawsList(fhRawTriggerCal ,197 ,saveForCorr);
1d7681da 168
5ed41460 169 TH1F* fhRawMeanCal = new TH1F("hRawMeanCal","online mean signal, calibration event",
54f32a55 170 10000,0,10000);
1d7681da 171 Add2RawsList( fhRawMeanCal,198);
5ed41460 172 TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration even ",
54f32a55 173 10000,0,10000);
1d7681da 174 Add2RawsList( fhRawVertexCal,199);
54f32a55 175 TH1F* fhRawORAcal = new TH1F("hRawORAcal","online OR A", 10000,0,10000);
1d7681da 176 Add2RawsList( fhRawORAcal,200 );
54f32a55 177 TH1F* fhRawORCcal = new TH1F("hRawORCcal","online OR C", 10000,0,10000);
1d7681da 178 Add2RawsList( fhRawORCcal,201);
54f32a55 179 TH1F* fhMultcal = new TH1F("hMultcal","full mulltiplicity", 10000,0,10000);
1d7681da 180 Add2RawsList( fhMultcal,202 );
5ed41460 181 TH1F* fhMultScal = new TH1F("hMultScal","full multiplicity with semi-central trigger",
54f32a55 182 10000,0,10000);
1d7681da 183 Add2RawsList( fhMultScal,203);
5ed41460 184 TH1F* fhMultCcal = new TH1F("hMultCcal","full multiplicity with central trigger",
185 1000,0,10000);
1d7681da 186 Add2RawsList( fhMultCcal,204);
187
188 TH2F* fhEffCFD = new TH2F("hEffCFD"," CFD time",24, 0 ,24,
189 100,-500,500);
c4b97145 190 Add2RawsList( fhEffCFD,205, saveForCorr);
1d7681da 191 TH2F* fhEffLED = new TH2F("hEffLED","LED time",24, 0 ,24,
192 100,-500,500);
c4b97145 193 Add2RawsList( fhEffLED,206, saveForCorr);
1d7681da 194 TH2F* fhEffQTC = new TH2F("hEffQTC","QTC amplitude",24, 0 ,24,
195 100,0,7000);
c4b97145 196 Add2RawsList( fhEffQTC,207, saveForCorr);
04236e67 197}
198
199//____________________________________________________________________________
200
201void AliT0QADataMakerRec::InitRecPoints()
202{
203 // create cluster histograms in RecPoint subdir
1d7681da 204
205
04236e67 206
207 TString timename,ampname, qtcname;
208 TH1F *fhRecCFD[24]; TH1F *fhRecLEDAmp[24]; TH1F * fhRecQTC[24];
209 for (Int_t i=0; i<24; i++)
210 {
211 timename ="hRecCFD";
212 ampname = "hRecLED";
213 qtcname = "hRecQTC";
214 timename += i;
215 ampname += i;
216 qtcname += i;
1d7681da 217 fhRecCFD[i] = new TH1F(timename.Data(), timename.Data(),2000,-1000,1000);
218 // fhRecCFD[i] = new TH1F(timename.Data(), timename.Data(),6000,1244000,1250000);
219 Add2RecPointsList ( fhRecCFD[i],i);
c724dd64 220 fhRecLEDAmp[i] = new TH1F(ampname.Data(), ampname.Data(),100,0, 100);
1d7681da 221 Add2RecPointsList ( fhRecLEDAmp[i],i+24);
c724dd64 222 fhRecQTC[i] = new TH1F(qtcname.Data(), qtcname.Data(),100,0,100);
1d7681da 223 Add2RecPointsList ( fhRecQTC[i],i+48);
04236e67 224 }
1d7681da 225
226
04236e67 227
1d7681da 228 // TH1F *fhOnlineMean = new TH1F("hOnlineMean","online mean",2000,-1000,1000);
229 TH1F *fhOnlineMean = new TH1F("hOnlineMean","online mean",6000,1244000,1250000);
c724dd64 230 Add2RecPointsList ( fhOnlineMean,72);
1d7681da 231 // TH1F * fhRecMean = new TH1F("hRecMean"," reconstructed mean signal",2000,-1000,1000);
232 TH1F * fhRecMean = new TH1F("hRecMean"," reconstructed mean signal",6000,1244000,1250000);
446d6ec4 233 Add2RecPointsList( fhRecMean,73);
1d7681da 234
04236e67 235
236}
237
238//____________________________________________________________________________
239void AliT0QADataMakerRec::InitESDs()
240{
241 //create ESDs histograms in ESDs subdir
f6c2d5e2 242
446d6ec4 243 TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean",100,2400,2500);
04236e67 244 Add2ESDsList(fhESDMean, 0) ;
1d7681da 245 TH1F * fhESDVertex = new TH1F("hESDvertex","ESD vertex",100,-50,50);
04236e67 246 Add2ESDsList(fhESDVertex, 1) ;
c724dd64 247
04236e67 248
249}
250
251//____________________________________________________________________________
252void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
253{
78328afd 254 rawReader->Reset() ;
04236e67 255 //fills QA histos for RAW
5ed41460 256 Int_t shift=0;
86c020af 257 AliT0RawReader *start = new AliT0RawReader(rawReader);
258 // start->Next();
259 if (! start->Next())
260 AliDebug(1,Form(" no raw data found!!"));
261 else
262 {
1d7681da 263
5ed41460 264 UInt_t type =rawReader->GetType();
265 // cout<<" !!!!! new event type = "<<type<<endl;
86c020af 266 Int_t allData[110][5];
267 for (Int_t i0=0; i0<105; i0++)
268 {
269 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
04236e67 270 }
86c020af 271 for (Int_t i=0; i<105; i++)
272 for (Int_t iHit=0; iHit<5; iHit++)
273 allData[i][iHit]= start->GetData(i,iHit);
274
1d7681da 275 if (allData[0][0]>0) GetRawsData(0) -> Fill( allData[0][0]);
276 // allData[0][0] = allData[0][0] - 7000;
277 if (type == 8) shift=101;
f6c2d5e2 278 if (type == 7) shift=0;
1d7681da 279
86c020af 280 for (Int_t ik = 0; ik<12; ik++){
281 for (Int_t iHt=0; iHt<5; iHt++){
1d7681da 282 //cfd
283 if(allData[ik+1][iHt]>0)
284 GetRawsData(shift+ik+1) ->
285 Fill(allData[ik+1][iHt]-allData[1][0]);
286 //led
287 if(allData[ik+13][iHt] > 0 && allData[13][iHt]>0)
288 GetRawsData(shift+ik+24+1)->
289 Fill(allData[ik+13][iHt]-allData[13][iHt]);
290 //led -cfd
291
292 if(allData[ik+13][iHt] > 0 && allData[ik+1][iHt] >0 )
293 GetRawsData(shift+ik+48+1)->
294 Fill(allData[ik+13][iHt]-allData[ik+1][iHt]);
295 //qtc
296 if(allData[2*ik+25][iHt] > 0 && allData[2*ik+26][iHt] > 0)
297 GetRawsData(shift+ik+72+1)->
298 Fill(allData[2*ik+25][iHt]-allData[2*ik+26][iHt]);
299
5ed41460 300
1d7681da 301 if(type == 8 && allData[ik+1][iHt]>0 )
302 GetRawsData(205)->Fill(ik,allData[ik+1][iHt]-allData[1][0]);
303 if(type == 8 && allData[ik+13][iHt]>0 )
304 GetRawsData(206)->Fill(ik,allData[ik+13][iHt]-allData[13][0]);
305 if(type == 8 && (allData[2*ik+25][iHt]>0 && allData[2*ik+26][iHt]>0) )
306 GetRawsData(207)->Fill(ik,allData[2*ik+25][iHt]-allData[2*ik+26][iHt]);
04236e67 307 }
86c020af 308 }
86c020af 309 for (Int_t ik = 12; ik<24; ik++) {
310 for (Int_t iHt=0; iHt<5; iHt++) {
1d7681da 311 if(allData[ik+45][iHt]>0)
312 //cfd
313 GetRawsData(shift+ik+1)->
314 Fill(allData[ik+45][iHt]-allData[57][0]);
315 //qtc
316 if(allData[2*ik+57][iHt]>0 && allData[2*ik+58][iHt]>0)
317 GetRawsData(shift+ik+72+1)->
318 Fill(allData[2*ik+57][iHt]-allData[2*ik+58][iHt]);
319 //led
320 if(allData[ik+57][iHt] > 0 )
321 GetRawsData(shift+ik+24+1)->
322 Fill(allData[ik+57][iHt]-allData[69][iHt]);
323 //led-cfd
324 if(allData[ik+57][iHt] > 0 &&allData[ik+45][iHt]>0)
325 GetRawsData(shift+ik+48+1)->
326 Fill(allData[ik+57][iHt]-allData[ik+45][iHt]);
327
328 if(type == 8 && allData[ik+45][iHt]>0 )
329 GetRawsData(205)->Fill(ik,allData[ik+45][iHt]-allData[57][0]);
330 if(type == 8 && allData[ik+57][iHt]>0 )
331 GetRawsData(206)->Fill(ik,allData[ik+57][iHt]-allData[69][0]);
332 if(type == 8 && (allData[2*ik+57][iHt]>0 && allData[2*ik+58][iHt]>0) )
333 GetRawsData(207)->Fill(ik,allData[2*ik+57][iHt]-allData[2*ik+58][iHt]);
334
335
86c020af 336 }
337 }
1d7681da 338 Int_t trChannel[6] = {49,50,51,52,55,56};
f6c2d5e2 339 if(type == 7)
5ed41460 340 {
1d7681da 341 for (Int_t iHt=0; iHt<6; iHt++) {
342 for (Int_t itr=0; itr<6; itr++) {
343 if(allData[trChannel[itr]][iHt]>0) fNumTriggers[itr]++;
344 }
f6c2d5e2 345 }
5ed41460 346 }
347 if(type == 8)
348 {
f6c2d5e2 349 for (Int_t iHt=0; iHt<5; iHt++) {
1d7681da 350 for (Int_t itr=0; itr<6; itr++) {
351 if(allData[trChannel[itr]][iHt]>0)
352 {
353
354 GetRawsData(198+itr)->Fill(allData[trChannel[itr]][iHt]-allData[1][0]);
355
356 fNumTriggersCal[itr]++;
357 }
358 }
359 if(allData[53][iHt]>0 && allData[54][iHt]>0)
360 GetRawsData(204)->Fill(allData[53][iHt]-allData[54][iHt]);
f6c2d5e2 361 }
1d7681da 362 }
363
364
365
86c020af 366 delete start;
367 }
04236e67 368}
1d7681da 369
370
04236e67 371
372//____________________________________________________________________________
373void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
374{
375 //fills QA histos for clusters
376
446d6ec4 377 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
04236e67 378 if (!frecpoints) {
379 AliError("Reconstruct Fill ESD >> no recpoints found");
380 return;
381 }
382 TBranch *brRec =clustersTree ->GetBranch("T0");
383 if (brRec) {
384 brRec->SetAddress(&frecpoints);
385 }else{
386 AliError(Form("EXEC Branch T0 rec not found "));
387 return;
388 }
389
390 brRec->GetEntry(0);
391
392 for ( Int_t i=0; i<24; i++) {
1d7681da 393 if(i<12)
394 GetRecPointsData(i) -> Fill( frecpoints -> GetTime(i) - frecpoints -> GetTime(0));
395 if(i>11)
396 GetRecPointsData(i) -> Fill( frecpoints -> GetTime(i) - frecpoints -> GetTime(12));
397 GetRecPointsData(i+24) -> Fill( frecpoints -> GetAmp(i));
398 GetRecPointsData(i+48) -> Fill( frecpoints->AmpLED(i));
04236e67 399 }
1d7681da 400 GetRecPointsData(72) ->Fill(frecpoints->GetOnlineMean());
401 GetRecPointsData(73) ->Fill(frecpoints->GetMeanTime());
04236e67 402
403}
404
405//____________________________________________________________________________
406void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
407{
408 //fills QA histos for ESD
409
410 GetESDsData(0) -> Fill(esd->GetT0());
411 GetESDsData(1)-> Fill(esd->GetT0zVertex());
412
413}