]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0QADataMakerSim.cxx
Removing some compilation warnings.
[u/mrichter/AliRoot.git] / T0 / AliT0QADataMakerSim.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 /* $Id$ */
18
19 //---
20 //  Produces the data needed to calculate the quality assurance. 
21 //  Alla.Maevskaya@cern.ch
22 //  
23 //---
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH2F.h> 
30 #include <TDirectory.h>
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliESDEvent.h"
35 #include "AliLog.h"
36 #include "AliT0digit.h" 
37 #include "AliT0hit.h"
38 #include "AliT0RecPoint.h"
39 #include "AliT0QADataMakerSim.h"
40 #include "AliQAChecker.h"
41 #include "AliT0RawReader.h"
42
43 #include <Riostream.h>
44
45 ClassImp(AliT0QADataMakerSim)
46            
47 //____________________________________________________________________________ 
48   AliT0QADataMakerSim::AliT0QADataMakerSim() : 
49   AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kT0), "T0 Quality Assurance Data Maker")
50
51 {
52   // ctor
53  //   fDetectorDir = fOutput->GetDirectory(GetName()) ;  
54 //   if (!fDetectorDir) 
55 //     fDetectorDir = fOutput->mkdir(GetName()) ;  
56 }
57
58 //____________________________________________________________________________ 
59 AliT0QADataMakerSim::AliT0QADataMakerSim(const AliT0QADataMakerSim& qadm) :
60   AliQADataMakerSim() 
61 {
62   //copy ctor 
63   SetName((const char*)qadm.GetName()) ; 
64   SetTitle((const char*)qadm.GetTitle()); 
65 }
66
67 //__________________________________________________________________
68 AliT0QADataMakerSim& AliT0QADataMakerSim::operator = (const AliT0QADataMakerSim& qadm )
69 {
70   // Equal operator.
71   this->~AliT0QADataMakerSim();
72   new(this) AliT0QADataMakerSim(qadm);
73   return *this; 
74 }
75 //____________________________________________________________________________
76 void AliT0QADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
77 {
78   //Detector specific actions at end of cycle
79   // do the QA checking
80   AliQAChecker::Instance()->Run(AliQAv1::kT0, task, list) ;
81 }
82
83 //____________________________________________________________________________
84 void AliT0QADataMakerSim::StartOfDetectorCycle()
85 {
86   //Detector specific actions at start of cycle
87
88 }
89  
90 //____________________________________________________________________________ 
91 void AliT0QADataMakerSim::InitHits()
92 {
93   // create Hits histograms in Hits subdir
94   // create Hits histograms in Hits subdir
95   const Bool_t expert   = kTRUE ; 
96   const Bool_t image    = kTRUE ; 
97
98   TString timename;
99   
100   TH2F *fhHitsTimeA = new TH2F("hHitsTimeA", "Hits Efficiency;#PMT; Time [ns];", 13, 12, 25, 100,12,15 );
101   fhHitsTimeA->SetOption("COLZ");
102  Add2HitsList(fhHitsTimeA,0, !expert, image);
103   TH2F *fhHitsTimeC = new TH2F("hHitsTimeC", "Hits Efficiency;#PMT; Time [ns];", 13, 0, 13, 100,2,5 );
104   fhHitsTimeC->SetOption("COLZ");
105   Add2HitsList(fhHitsTimeC,1, !expert, image);
106 }
107
108 //____________________________________________________________________________ 
109 void AliT0QADataMakerSim::InitDigits()
110 {
111   // create Digits histograms in Digits subdir
112   const Bool_t expert   = kTRUE ; 
113   const Bool_t image    = kTRUE ; 
114   
115   TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD time [#channel]",25,-0.5,24.5,100,0,1000);
116   fhDigCFD->SetOption("COLZ");
117   Add2DigitsList( fhDigCFD,0);
118   TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; amplitude  LED-CFD [#channel]",25,-0.5,24.5,100,100,1000);
119   fhDigLEDamp->SetOption("COLZ");
120   Add2DigitsList( fhDigLEDamp,1, !expert, image);
121   TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; amplitude QTC [#channel]",25,-0.5,24.5,200,500,10000);
122   fhDigQTC->SetOption("COLZ");
123   Add2DigitsList( fhDigQTC,2, !expert, image);
124   
125   
126    
127 }
128
129 //____________________________________________________________________________
130
131 void AliT0QADataMakerSim::MakeHits(TTree *hitTree)
132 {
133   //fills QA histos for Hits
134   if (fHitsArray) 
135     fHitsArray->Clear() ; 
136   else 
137     fHitsArray = new TClonesArray("AliT0hit", 1000);
138   
139   TBranch * branch = hitTree->GetBranch("T0") ;
140   if ( ! branch ) {
141     AliWarning("T0 branch in Hit Tree not found") ;
142   } else {
143
144    if (branch) {
145       branch->SetAddress(&fHitsArray);
146     }else{
147       AliError("Branch T0 hit not found");
148       exit(111);
149     } 
150     Int_t ntracks    = (Int_t) hitTree->GetEntries();
151     
152     if (ntracks<=0) return;
153     // Start loop on tracks in the hits containers
154
155     for (Int_t track=0; track<ntracks;track++) {
156       branch->GetEntry(track);
157       Int_t nhits = fHitsArray->GetEntriesFast();
158       for (Int_t ihit=0;ihit<nhits;ihit++) 
159         {
160           AliT0hit  * startHit   = (AliT0hit*) fHitsArray->UncheckedAt(ihit);
161           if (!startHit) {
162             AliError("The unchecked hit doesn't exist");
163             continue;
164           }
165           Int_t pmt=startHit->Pmt();
166           Float_t time = 0.001 * startHit->Time();
167           if(pmt<13)GetHitsData(1)->Fill(pmt,time) ;
168           if(pmt>12)GetHitsData(0)->Fill(pmt,time) ;
169         }
170     }
171   }
172 }
173
174 //____________________________________________________________________________
175 void AliT0QADataMakerSim::MakeDigits( TTree *digitsTree)
176 {
177   //fills QA histos for Digits
178  
179   TArrayI *digCFD = new TArrayI(24);
180   TArrayI *digLED = new TArrayI(24);
181   TArrayI *digQT0 = new TArrayI(24);
182   TArrayI *digQT1 = new TArrayI(24);
183   Int_t refpoint=0;
184
185   TBranch *brDigits=digitsTree->GetBranch("T0");
186   AliT0digit *fDigits = new AliT0digit() ;
187   if (brDigits) {
188     brDigits->SetAddress(&fDigits);
189   }else{
190     AliError(Form("EXEC Branch T0 digits not found"));
191      return;
192   }
193
194   digitsTree->GetEvent(0);
195   digitsTree->GetEntry(0);
196   brDigits->GetEntry(0);
197   fDigits->GetTimeCFD(*digCFD);
198   fDigits->GetTimeLED(*digLED);
199   fDigits->GetQT0(*digQT0);
200   fDigits->GetQT1(*digQT1);
201   refpoint = fDigits->RefPoint();
202
203    for (Int_t i=0; i<24; i++)
204     {
205       if (digCFD->At(i)>0) {
206         Int_t cfd=digCFD->At(i)- refpoint;
207         GetDigitsData(0) ->Fill(i,cfd);
208         GetDigitsData(1) -> Fill(i,(digLED->At(i) - digCFD->At(i)));
209         GetDigitsData(2) -> Fill(i, (digQT1->At(i) - digQT0->At(i)));
210
211       }
212     }  
213       
214   delete digCFD;
215   delete digLED;
216   delete digQT0;
217   delete digQT1;
218
219 }