]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ACORDE/AliACORDEQADataMakerSim.cxx
Update for ACORDE-QA (Mario & Luciano)
[u/mrichter/AliRoot.git] / ACORDE / AliACORDEQADataMakerSim.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
18 //---
19 //  Produces the data needed to calculate the quality assurance. 
20 //  All data must be mergeable objects.
21
22 //  Authors:
23 //
24 //  Luciano Diaz Gonzalez <luciano.diaz@nucleares.unam.mx> (ICN-UNAM)
25 //  Mario Rodriguez Cahuantzi <mrodrigu@mail.cern.ch> (FCFM-BUAP)
26 //  Arturo Fernandez Tellez <afernan@mail.cern.ch (FCFM-BUAP)
27 //
28 //  Created: June 13th 2008
29 //---
30
31 // --- ROOT system ---
32 #include <TClonesArray.h>
33 #include <TFile.h> 
34 #include <TH1F.h> 
35 #include <TH2F.h>
36 #include <TDirectory.h>
37 // --- Standard library ---
38
39 // --- AliRoot header files ---
40 #include "AliESDEvent.h"
41 #include "AliLog.h"
42 #include "AliACORDEdigit.h" 
43 #include "AliACORDEhit.h"
44 #include "AliACORDEQADataMakerSim.h"
45 #include "AliQAChecker.h"
46 #include "AliACORDERawReader.h"
47 ClassImp(AliACORDEQADataMakerSim)
48            
49 //____________________________________________________________________________ 
50 AliACORDEQADataMakerSim::AliACORDEQADataMakerSim():AliQADataMakerSim(AliQA::GetDetName(AliQA::kACORDE), "ACORDE Quality Assurance Data Maker")
51 {
52 }
53 //____________________________________________________________________________ 
54 AliACORDEQADataMakerSim::AliACORDEQADataMakerSim(const AliACORDEQADataMakerSim& qadm) :
55   AliQADataMakerSim() 
56 {
57   SetName((const char*)qadm.GetName()) ; 
58   SetTitle((const char*)qadm.GetTitle()); 
59 }
60 //__________________________________________________________________
61 AliACORDEQADataMakerSim& AliACORDEQADataMakerSim::operator = (const AliACORDEQADataMakerSim& qadm )
62 {
63   // Equal operator.
64   this->~AliACORDEQADataMakerSim();
65   new(this) AliACORDEQADataMakerSim(qadm);
66   return *this;
67 }
68 //____________________________________________________________________________
69 void AliACORDEQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
70 {
71   //Detector specific actions at end of cycle
72   // do the QA checking
73    printf("ACORDE---->Detector specific actions at END of cycle\n................\n");
74
75   AliQAChecker::Instance()->Run(AliQA::kACORDE, task, list) ;
76 }
77 //____________________________________________________________________________
78 void AliACORDEQADataMakerSim::StartOfDetectorCycle()
79 {
80   //Detector specific actions at start of cycle
81   printf("ACORDE---->Detector specific actions at START of cycle\n................\n");
82 }
83 //____________________________________________________________________________ 
84 void AliACORDEQADataMakerSim::InitHits()
85 {
86   // create Hits histograms in Hits subdir
87         TH1F *fAHitsACORDE[8];
88
89         fAHitsACORDE[0] = new TH1F("hACORDEloss" ,"Energy Loss ",1000,0.,1500.);
90         fAHitsACORDE[1] = new TH1F("hACORDEPolar" ," Polar Angle ",90,0.,90.);
91         fAHitsACORDE[2] = new TH1F("hACORDEAzimuth" ,"Azimuth Angle  ",360,-180.,180.);
92         fAHitsACORDE[3] = new TH1F("hACORDEPx" ,"Px Distribution ",60,-30.,30.);
93         fAHitsACORDE[4] = new TH1F("hACORDEPy" ,"Py Distribution ",60,-30.,30.);
94         fAHitsACORDE[5] = new TH1F("hACORDEPz" ,"Pz Distribution ",60,-30.,30.);
95         fAHitsACORDE[6] = new TH1F("hACORDEPt" ,"Pt Distribution ",60,0.,50.);
96         fAHitsACORDE[7] = new TH1F("hACORDEpxpz" ,"Pt Distribution ",100,-50.,50.);
97
98         TH2F *hACORDExy = new TH2F("hACORDExy" ,"Dist. xy",2800,-2400.,1400.,200,-4805.,4825.);
99         TH2F *hACORDExz = new TH2F("hACORDExz" ,"Dist.xz ",900,-1500.,2850.,1200,-1000.,4000.);
100         TH2F *hACORDEyz = new TH2F("hACORDEyz" ,"Dist.yz ",5,817.,819.,1200,-600.,600.);
101         TH2F *hACORDEAzimPol =  new TH2F("hACORDEAzimPol" ,"Azimuth vs Polar ",360,-180.,180.,180,0.,180.);
102
103         for(Int_t i=0; i<8; i++)
104            Add2HitsList(fAHitsACORDE[i],i);
105  
106          Add2HitsList(hACORDExy,8);
107          Add2HitsList(hACORDExz,9);
108          Add2HitsList(hACORDEyz,10);
109          Add2HitsList(hACORDEAzimPol,11);
110
111
112
113 }
114 //____________________________________________________________________________ 
115 void AliACORDEQADataMakerSim::InitDigits()
116 {
117   // create Digits histograms in Digits subdir
118
119    TH1F *    fhDigitsModule;
120    TString   modulename;
121    modulename = "hDigitsModule";
122    fhDigitsModule = new TH1F(modulename.Data(),"hDigitsModuleSingle",60,0,60);
123    Add2DigitsList( fhDigitsModule,0);
124
125 }
126 //____________________________________________________________________________
127
128 void AliACORDEQADataMakerSim::MakeHits(TTree *hitTree)
129 {
130   // Here we fill the QA histos for Hits declared above
131
132         printf("Estamos en make Hits");
133         TClonesArray * hits = new TClonesArray("AliACORDEhit",1000);
134         TBranch * branch = hitTree->GetBranch("ACORDE");
135         if (!branch)
136         {
137                 AliWarning("ACORDE branch in Hit Tree not found");
138         }else
139         {
140                 if (branch)
141                 {
142                         branch->SetAddress(&hits);
143                 }else
144                 {
145                         AliError("Branch ACORDE hit not found");
146                         exit(111);
147                 }
148                 Int_t ntracks = (Int_t)hitTree->GetEntries();
149                 if (ntracks<=0) return;
150                 for(Int_t track=0;track<ntracks;track++)
151                 {
152                         branch->GetEntry(track);
153                         Int_t nhits = hits->GetEntriesFast();
154                         for(Int_t ihit=0;ihit<nhits;ihit++)
155                         {
156                                 AliACORDEhit *AcoHit = (AliACORDEhit*) hits->UncheckedAt(ihit);
157                                 if (!AcoHit)
158                                 {
159                                         AliError("The unchecked hit doesn't exist");
160                                         break;
161                                 }
162                                 GetHitsData(0)->Fill(AcoHit->Eloss());
163                                 GetHitsData(1)->Fill(AcoHit->PolarAngle());
164                                 GetHitsData(2)->Fill(AcoHit->AzimuthAngle());
165                                 GetHitsData(3)->Fill(AcoHit->Px());
166                                 GetHitsData(4)->Fill(AcoHit->Py());
167                                 GetHitsData(5)->Fill(AcoHit->Pz());
168                                 GetHitsData(6)->Fill(TMath::Sqrt( (AcoHit->Px())*(AcoHit->Px())+
169                                              (AcoHit->Py())*(AcoHit->Py())));
170                                 if((AcoHit->Py()) != 0.0 ) GetHitsData(7)->Fill(TMath::ATan(AcoHit->Px()/AcoHit->Py()));
171                                 GetHitsData(8)->Fill( (Float_t)(AcoHit->X()),(Float_t)(AcoHit->Y()) );
172                                 GetHitsData(9)->Fill( (Float_t)(AcoHit->X()),(Float_t)(AcoHit->Z()) );
173                                 GetHitsData(10)->Fill( (Float_t)(AcoHit->Y()),(Float_t)(AcoHit->Z()) );
174                                 GetHitsData(11)->Fill( (Float_t)(AcoHit->AzimuthAngle()),
175                                 (Float_t)(AcoHit->PolarAngle()));
176                         }
177                 }
178         }
179
180 }
181 //____________________________________________________________________________
182 void AliACORDEQADataMakerSim::MakeDigits( TTree *digitsTree)
183 {
184   //fills QA histos for Digits
185
186
187         TClonesArray * digits = new TClonesArray("AliACORDEdigit",1000);
188         TBranch * branch = digitsTree->GetBranch("ACORDEdigit");
189         if (!branch)
190         {
191                 AliWarning("ACORDE branch in Digits Tree not found");
192         }else
193         {
194                 if (branch)
195                 {
196                         branch->SetAddress(&digits);
197                 }else
198                 {
199                         AliError("Branch ACORDE digit not found");
200                         exit(111);
201                 }
202                 Int_t ntracks = (Int_t)digitsTree->GetEntries();
203                 if (ntracks<=0) return;
204                 printf("Entries in DigitsTree:%d\n",ntracks);
205                 for(Int_t track=0;track<ntracks;track++)
206                 {
207                         branch->GetEntry(track);
208                         Int_t ndigits = digits->GetEntriesFast();
209                         for(Int_t idigit=0;idigit<ndigits;idigit++)
210                         {
211                                 AliACORDEdigit *AcoDigit = (AliACORDEdigit*) digits->UncheckedAt(idigit);
212                                 if (!AcoDigit)
213                                 {
214                                         AliError("The unchecked digit doesn't exist");
215                                         break;
216                                 }
217                                 GetDigitsData(0)->Fill(AcoDigit->GetModule()-1);
218                         }
219                 }
220
221
222         }
223
224 }