]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCQADataMakerRec.cxx
Compilation on Windows/Cygwin
[u/mrichter/AliRoot.git] / ZDC / AliZDCQADataMakerRec.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 // --- ROOT system ---
16 #include <TClonesArray.h>
17 #include <TFile.h> 
18 #include <TH1F.h> 
19 #include <TH2F.h>
20 #include <TProfile.h>
21 #include <Riostream.h>
22 // --- Standard library ---
23
24 // --- AliRoot header files ---
25 #include "AliLog.h"
26 #include "AliQAChecker.h"
27 #include "AliZDCQADataMakerRec.h"
28 #include "AliZDCRawStream.h"
29 #include "AliESDZDC.h"
30 #include "AliESDEvent.h"
31
32 ClassImp(AliZDCQADataMakerRec)
33            
34 //____________________________________________________________________________ 
35   AliZDCQADataMakerRec::AliZDCQADataMakerRec() : 
36   AliQADataMakerRec(AliQA::GetDetName(AliQA::kZDC), "ZDC Quality Assurance Data Maker")
37 {
38   // ctor
39 }
40
41 //____________________________________________________________________________ 
42 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
43   AliQADataMakerRec() 
44 {
45   //copy ctor 
46   SetName((const char*)qadm.GetName()); 
47   SetTitle((const char*)qadm.GetTitle()); 
48 }
49
50 //__________________________________________________________________
51 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
52 {
53   // Equal operator.
54   this->~AliZDCQADataMakerRec();
55   new(this) AliZDCQADataMakerRec(qadm);
56   return *this;
57 }
58
59 //____________________________________________________________________________
60
61 void AliZDCQADataMakerRec::InitRaws()
62 {
63   // create Digits histograms in Digits subdir
64   //
65   TH1F * hRawZNCTot = new TH1F("hRawZNCTot", "Raw signal in ZNC", 100, 0., 6000.);
66   TH1F * hRawZNATot = new TH1F("hRawZNATot", "Raw signal in ZNA", 100, 0., 6000.);
67   TH1F * hRawZPCTot = new TH1F("hRawZPCTot", "Raw signal in ZPC", 100, 0., 10000.);
68   TH1F * hRawZPATot = new TH1F("hRawZPATot", "Raw signal in ZPA", 100, 0., 10000.);
69   Add2RawsList(hRawZNCTot, 0);
70   Add2RawsList(hRawZPCTot, 1);
71   Add2RawsList(hRawZNATot, 2);
72   Add2RawsList(hRawZPATot, 3);
73   //
74   TH1F * hRawSumQZNC = new TH1F("hRawSumQZNC", "Signal in 4 ZNC PMQ[i]",100, 0., 4000.);
75   TH1F * hRawSumQZPC = new TH1F("hRawSumQZPC", "Signal in 4 ZPC PMQ[i]",100, 0., 4000.);
76   TH1F * hRawSumQZNA = new TH1F("hRawSumQZNA", "Signal in 4 ZNA PMQ[i]",100, 0., 4000.);
77   TH1F * hRawSumQZPA = new TH1F("hRawSumQZPA", "Signal in 4 ZPA PMQ[i]",100, 0., 4000.);
78   Add2RawsList(hRawSumQZNC, 4);
79   Add2RawsList(hRawSumQZPC, 5);
80   Add2RawsList(hRawSumQZNA, 6);
81   Add2RawsList(hRawSumQZPA, 7);
82   //
83   TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Signal in 4 ZNC PMQ[i]",100, 0., 4000.);
84   TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Signal in 4 ZPC PMQ[i]",100, 0., 4000.);
85   TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Signal in 4 ZNA PMQ[i]",100, 0., 4000.);
86   TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Signal in 4 ZPA PMQ[i]",100, 0., 4000.);
87   Add2RawsList(hRawPMCZNC, 8);
88   Add2RawsList(hRawPMCZPC, 9);
89   Add2RawsList(hRawPMCZNA, 10);
90   Add2RawsList(hRawPMCZPA, 11);
91   // 
92   // ------------------- LOW GAIN CHAIN ---------------------------
93   TH1F * hRawZNCTotlg = new TH1F("hRawZNCTotlg", "Rawit lg signal in ZNC", 100, 0., 6000.);
94   TH1F * hRawZNATotlg = new TH1F("hRawZNATotlg", "Rawit lg signal in ZNA", 100, 0., 6000.);
95   TH1F * hRawZPCTotlg = new TH1F("hRawZPCTotlg", "Rawit lg signal in ZPC", 100, 0., 10000.);
96   TH1F * hRawZPATotlg = new TH1F("hRawZPATotlg", "Rawit lg signal in ZPA", 100, 0., 10000.);
97   Add2RawsList(hRawZNCTotlg, 12);
98   Add2RawsList(hRawZPCTotlg, 13);
99   Add2RawsList(hRawZNATotlg, 14);
100   Add2RawsList(hRawZPATotlg, 15);
101   //
102   TH1F * hRawSumQZNClg = new TH1F("hRawSumQZNClg", "Signal in 4 ZNC PMQlg[i]",100, 0., 4000.);
103   TH1F * hRawSumQZPClg = new TH1F("hRawSumQZPClg", "Signal in 4 ZPC PMQlg[i]",100, 0., 4000.);
104   TH1F * hRawSumQZNAlg = new TH1F("hRawSumQZNAlg", "Signal in 4 ZNA PMQlg[i]",100, 0., 4000.);
105   TH1F * hRawSumQZPAlg = new TH1F("hRawSumQZPAlg", "Signal in 4 ZPA PMQlg[i]",100, 0., 4000.);
106   Add2RawsList(hRawSumQZNClg, 16);
107   Add2RawsList(hRawSumQZPClg, 17);
108   Add2RawsList(hRawSumQZNAlg, 18);
109   Add2RawsList(hRawSumQZPAlg, 19);
110   //
111   TH1F * hRawPMCZNClg = new TH1F("hRawPMCZNClg", "Signal in 4 ZNC PMQlg[i]",100, 0., 4000.);
112   TH1F * hRawPMCZPClg = new TH1F("hRawPMCZPClg", "Signal in 4 ZPC PMQlg[i]",100, 0., 4000.);
113   TH1F * hRawPMCZNAlg = new TH1F("hRawPMCZNAlg", "Signal in 4 ZNA PMQlg[i]",100, 0., 4000.);
114   TH1F * hRawPMCZPAlg = new TH1F("hRawPMCZPAlg", "Signal in 4 ZPA PMQlg[i]",100, 0., 4000.);
115   Add2RawsList(hRawPMCZNClg, 20);
116   Add2RawsList(hRawPMCZPClg, 21);
117   Add2RawsList(hRawPMCZNAlg, 22);
118   Add2RawsList(hRawPMCZPAlg, 23);
119 }
120
121 //____________________________________________________________________________
122 void AliZDCQADataMakerRec::InitESDs()
123 {
124   //Booking ESDs histograms
125   //
126   TH2F * hZNC  = new TH2F("hZNC", "Centroid in ZNC", 100, -5.,5.,100,-5.,5.);
127   TH2F * hZNA  = new TH2F("hZNA", "Centroid in ZNA", 100, -5.,5.,100,-5.,5.);
128   Add2DigitsList(hZNC, 0);
129   Add2DigitsList(hZNA, 1);
130   //
131   TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "ESD signal in ZNC", 100, 0., 6000.);
132   TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "ESD signal in ZPC", 100, 0., 10000.);
133   TH1F * hESDZNATot = new TH1F("hESDZNATot", "ESD signal in ZNA", 100, 0., 6000.);
134   TH1F * hESDZPATot = new TH1F("hESDZPATot", "ESD signal in ZPA", 100, 0., 10000.);
135   Add2ESDsList(hESDZNCTot, 2);
136   Add2ESDsList(hESDZPCTot, 3);
137   Add2ESDsList(hESDZNATot, 4);
138   Add2ESDsList(hESDZPATot, 5);
139   //
140   TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Signal in 4 ZNC PMQ[i]",100, 0., 4000.);
141   TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Signal in 4 ZPC PMQ[i]",100, 0., 4000.);
142   TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Signal in 4 ZNA PMQ[i]",100, 0., 4000.);
143   TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Signal in 4 ZPA PMQ[i]",100, 0., 4000.);
144   Add2ESDsList(hESDSumQZNC, 6);
145   Add2ESDsList(hESDSumQZPC, 7);
146   Add2ESDsList(hESDSumQZNA, 8);
147   Add2ESDsList(hESDSumQZPA, 9);
148   //
149   TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Signal in 4 ZNC PMQ[i]",100, 0., 4000.);
150   TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Signal in 4 ZPC PMQ[i]",100, 0., 4000.);
151   TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Signal in 4 ZNA PMQ[i]",100, 0., 4000.);
152   TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Signal in 4 ZPA PMQ[i]",100, 0., 4000.);
153   Add2ESDsList(hESDPMCZNC, 10);
154   Add2ESDsList(hESDPMCZPC, 11);
155   Add2ESDsList(hESDPMCZNA, 12);
156   Add2ESDsList(hESDPMCZPA, 13);
157   // 
158   // ------------------- LOW GAIN CHAIN ---------------------------
159   TH1F * hESDZNCTotlg = new TH1F("hESDZNCTotlg", "ESDit lg signal in ZNC", 100, 0., 6000.);
160   TH1F * hESDZNATotlg = new TH1F("hESDZNATotlg", "ESDit lg signal in ZNA", 100, 0., 6000.);
161   TH1F * hESDZPCTotlg = new TH1F("hESDZPCTotlg", "ESDit lg signal in ZPC", 100, 0., 10000.);
162   TH1F * hESDZPATotlg = new TH1F("hESDZPATotlg", "ESDit lg signal in ZPA", 100, 0., 10000.);
163   Add2ESDsList(hESDZNCTotlg, 14);
164   Add2ESDsList(hESDZPCTotlg, 15);
165   Add2ESDsList(hESDZNATotlg, 16);
166   Add2ESDsList(hESDZPATotlg, 17);
167   //
168   TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Signal in 4 ZNC PMQlg[i]",100, 0., 4000.);
169   TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Signal in 4 ZPC PMQlg[i]",100, 0., 4000.);
170   TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Signal in 4 ZNA PMQlg[i]",100, 0., 4000.);
171   TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Signal in 4 ZPA PMQlg[i]",100, 0., 4000.);
172   Add2ESDsList(hESDSumQZNClg, 18);
173   Add2ESDsList(hESDSumQZPClg, 19);
174   Add2ESDsList(hESDSumQZNAlg, 20);
175   Add2ESDsList(hESDSumQZPAlg, 21);
176   //
177   TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in 4 ZNC PMQlg[i]",100, 0., 4000.);
178   TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in 4 ZPC PMQlg[i]",100, 0., 4000.);
179   TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in 4 ZNA PMQlg[i]",100, 0., 4000.);
180   TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in 4 ZPA PMQlg[i]",100, 0., 4000.);
181   Add2ESDsList(hESDPMCZNClg, 22);
182   Add2ESDsList(hESDPMCZPClg, 23);
183   Add2ESDsList(hESDPMCZNAlg, 24);
184   Add2ESDsList(hESDPMCZPAlg, 25);
185 }
186   
187 //____________________________________________________________________________
188
189 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
190 {
191   // Filling Raws QA histos
192   //
193   Float_t Sum_ZNC=0., Sum_ZNA=0., Sum_ZPC=0., Sum_ZPA=0.;
194   Float_t SumQ_ZNC=0., SumQ_ZNA=0., SumQ_ZPC=0., SumQ_ZPA=0.;
195   Float_t Sum_ZNC_lg=0., Sum_ZNA_lg=0., Sum_ZPC_lg=0., Sum_ZPA_lg=0.;
196   Float_t SumQ_ZNC_lg=0., SumQ_ZNA_lg=0., SumQ_ZPC_lg=0., SumQ_ZPA_lg=0.;
197   //
198   AliZDCRawStream stream(rawReader);
199   while(stream.Next()){
200     if(stream.IsADCDataWord() && 
201      (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
202        if(stream.GetSector(0)==1){
203          if(stream.GetADCGain()==0){
204            Sum_ZNC += stream.GetADCValue();
205            if(stream.GetSector(1)!=0) SumQ_ZNC += stream.GetADCValue();
206            else GetRawsData(8)->Fill(stream.GetADCValue());
207          }
208          else{
209            Sum_ZNC_lg += stream.GetADCValue();
210            if(stream.GetSector(1)!=0) SumQ_ZNC_lg += stream.GetADCValue();
211            else GetRawsData(20)->Fill(stream.GetADCValue());
212          }
213        }
214        else if(stream.GetSector(0)==2){
215          if(stream.GetADCGain()==0){
216            Sum_ZPC += stream.GetADCValue();
217            if(stream.GetSector(1)!=0) SumQ_ZPC += stream.GetADCValue();
218            else GetRawsData(9)->Fill(stream.GetADCValue());
219          }
220          else{
221            Sum_ZPC_lg += stream.GetADCValue();
222            if(stream.GetSector(1)!=0) SumQ_ZPC_lg += stream.GetADCValue();
223            else GetRawsData(21)->Fill(stream.GetADCValue());
224          }
225        }
226        else if(stream.GetSector(0)==4){
227          if(stream.GetADCGain()==0){
228            Sum_ZNA += stream.GetADCValue();
229            if(stream.GetSector(1)!=0) SumQ_ZNA += stream.GetADCValue();
230            else GetRawsData(10)->Fill(stream.GetADCValue());
231          }
232          else{
233            Sum_ZNA_lg += stream.GetADCValue();
234            if(stream.GetSector(1)!=0) SumQ_ZNA_lg += stream.GetADCValue();
235            else GetRawsData(22)->Fill(stream.GetADCValue());
236          }
237        }
238        else if(stream.GetSector(0)==5){
239          if(stream.GetADCGain()==0){
240            Sum_ZPA += stream.GetADCValue();
241            if(stream.GetSector(1)!=0) SumQ_ZPA += stream.GetADCValue();
242            else GetRawsData(11)->Fill(stream.GetADCValue());
243          }
244          else{
245            Sum_ZPA_lg += stream.GetADCValue();
246            if(stream.GetSector(1)!=0) SumQ_ZPA_lg += stream.GetADCValue();
247            else GetRawsData(23)->Fill(stream.GetADCValue());
248          }
249        }
250     }
251   }
252   //
253   GetRawsData(0)->Fill(Sum_ZNC);
254   GetRawsData(1)->Fill(Sum_ZPC);
255   GetRawsData(2)->Fill(Sum_ZNA);
256   GetRawsData(3)->Fill(Sum_ZPA);
257   //
258   GetRawsData(4)->Fill(SumQ_ZNC);
259   GetRawsData(5)->Fill(SumQ_ZPC);
260   GetRawsData(6)->Fill(SumQ_ZNA);
261   GetRawsData(7)->Fill(SumQ_ZPA);
262   //
263   GetRawsData(12)->Fill(Sum_ZNC_lg);
264   GetRawsData(13)->Fill(Sum_ZPC_lg);
265   GetRawsData(14)->Fill(Sum_ZNA_lg);
266   GetRawsData(15)->Fill(Sum_ZPA_lg);
267   //
268   GetRawsData(16)->Fill(SumQ_ZNC_lg);
269   GetRawsData(17)->Fill(SumQ_ZPC_lg);
270   GetRawsData(18)->Fill(SumQ_ZNA_lg);
271   GetRawsData(19)->Fill(SumQ_ZPA_lg);
272   //
273   stream.Delete();
274 }
275
276 //____________________________________________________________________________
277 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
278 {
279   // make QA data from ESDs
280   //
281   AliESDZDC * zdcESD =  esd->GetESDZDC();
282   //
283   const Float_t * Centr_ZNC, * Centr_ZNA;
284   Int_t NSpecnC = (Int_t) (esd->GetZDCN1Energy()/2.7);
285   Int_t NSpecnA = (Int_t) (esd->GetZDCN2Energy()/2.7);
286   Centr_ZNC = zdcESD->GetZNCCentroid(NSpecnC);
287   Centr_ZNA = zdcESD->GetZNACentroid(NSpecnA);
288   GetESDsData(0)->Fill(Centr_ZNC[0], Centr_ZNC[1]);
289   GetESDsData(1)->Fill(Centr_ZNA[0], Centr_ZNA[1]);
290   //
291   GetESDsData(2)->Fill(esd->GetZDCN1Energy());
292   GetESDsData(3)->Fill(esd->GetZDCP1Energy());
293   GetESDsData(4)->Fill(esd->GetZDCN2Energy());
294   GetESDsData(5)->Fill(esd->GetZDCP2Energy());
295   //
296   Double_t ZNCSumQ=0., ZPCSumQ=0., ZNASumQ=0., ZPASumQ=0.;
297   Double_t ZNCSumQ_lg=0., ZPCSumQ_lg=0., ZNASumQ_lg=0., ZPASumQ_lg=0.;
298   //
299   const Double_t *ZNCTow, *ZPCTow, *ZNATow, *ZPATow;
300   const Double_t *ZNCTow_lg, *ZPCTow_lg, *ZNATow_lg, *ZPATow_lg;
301   //
302   ZNCTow = zdcESD->GetZN1TowerEnergy();
303   ZPCTow = zdcESD->GetZP1TowerEnergy();
304   ZNATow = zdcESD->GetZN2TowerEnergy();
305   ZPATow = zdcESD->GetZP2TowerEnergy();
306   //
307   ZNCTow_lg = zdcESD->GetZN1TowerEnergyLR();
308   ZPCTow_lg = zdcESD->GetZP1TowerEnergyLR();
309   ZNATow_lg = zdcESD->GetZN2TowerEnergyLR();
310   ZPATow_lg = zdcESD->GetZP2TowerEnergyLR();
311   //
312   for(Int_t i=0; i<5; i++){
313      if(i==0){
314        GetESDsData(10)->Fill(ZNCTow[i]);
315        GetESDsData(11)->Fill(ZPCTow[i]);
316        GetESDsData(12)->Fill(ZNATow[i]);
317        GetESDsData(13)->Fill(ZPATow[i]);
318        //
319        GetESDsData(22)->Fill(ZNCTow_lg[i]);
320        GetESDsData(23)->Fill(ZPCTow_lg[i]);
321        GetESDsData(24)->Fill(ZNATow_lg[i]);
322        GetESDsData(25)->Fill(ZPATow_lg[i]);
323      }
324      else{
325        ZNCSumQ += ZNCTow[i];
326        ZPCSumQ += ZPCTow[i];
327        ZNASumQ += ZNATow[i];
328        ZPASumQ += ZPATow[i];
329        //
330        ZNCSumQ_lg += ZNCTow_lg[i];
331        ZPCSumQ_lg += ZPCTow_lg[i];
332        ZNASumQ_lg += ZNATow_lg[i];
333        ZPASumQ_lg += ZPATow_lg[i];
334      }
335   }
336   GetESDsData(6)->Fill(ZNCSumQ);
337   GetESDsData(7)->Fill(ZPCSumQ);
338   GetESDsData(8)->Fill(ZNASumQ);
339   GetESDsData(9)->Fill(ZPASumQ);
340   //
341   GetESDsData(18)->Fill(ZNCSumQ_lg);
342   GetESDsData(19)->Fill(ZPCSumQ_lg);
343   GetESDsData(20)->Fill(ZNASumQ_lg);
344   GetESDsData(21)->Fill(ZPASumQ_lg);
345 }
346
347 //____________________________________________________________________________
348 void AliZDCQADataMakerRec::StartOfDetectorCycle()
349 {
350   //Detector specific actions at start of cycle
351   
352 }
353
354 //____________________________________________________________________________ 
355 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
356 {
357   //Detector specific actions at end of cycle
358   // do the QA checking
359   AliQAChecker::Instance()->Run(AliQA::kZDC, task, list) ;  
360 }
361