]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/UserTasks/CaloCellQA/AliAnalysisTaskCaloCellsQA.cxx
1) static variable were replaced by local variables and data members
[u/mrichter/AliRoot.git] / PWG4 / UserTasks / CaloCellQA / AliAnalysisTaskCaloCellsQA.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 // Container class for bad channels & bad runs identification
16 //
17 // By default, pileup events are not processed, use SetAvoidPileup()
18 // to change this.
19 //
20 // Clusters containing a bad cell can be rejected, use SetBadCells().
21 //
22 // See AddTaskCaloCellsQA.C for usage example.
23 //
24 //----
25 //  Author: Olga Driga (SUBATECH)
26
27 // --- ROOT system ---
28 #include <TFile.h>
29 #include <TObjArray.h>
30
31 // --- AliRoot header files ---
32 #include <AliAnalysisTaskCaloCellsQA.h>
33 #include <AliVEvent.h>
34 #include <AliVCaloCells.h>
35 #include <AliVCluster.h>
36 #include <AliVVertex.h>
37
38 ClassImp(AliAnalysisTaskCaloCellsQA)
39
40 //________________________________________________________________
41 AliAnalysisTaskCaloCellsQA::AliAnalysisTaskCaloCellsQA(const char *name) : AliAnalysisTaskSE(name),
42   fkAvoidPileup(kTRUE),
43   fCellsQA(0),
44   fOutfile(new TString),
45   fBadCells(0),
46   fNBad(0)
47 {
48 }
49
50 //________________________________________________________________
51 AliAnalysisTaskCaloCellsQA::~AliAnalysisTaskCaloCellsQA()
52 {
53   if (fCellsQA) delete fCellsQA;
54   delete fOutfile;
55   if (fBadCells) delete [] fBadCells;
56 }
57
58 //________________________________________________________________
59 void AliAnalysisTaskCaloCellsQA::InitCaloCellsQA(char* fname, Int_t nmods, Int_t det)
60 {
61   // Must be called at the very beginning.
62   //
63   // fname -- output file name;
64   // nmods -- number of supermodules + 1;
65   // det -- detector;
66
67   if (det == kEMCAL)
68     fCellsQA = new AliCaloCellsQA(nmods, AliCaloCellsQA::kEMCAL);
69   else if (det == kPHOS)
70     fCellsQA = new AliCaloCellsQA(nmods, AliCaloCellsQA::kPHOS);
71   else
72     Fatal("AliAnalysisTaskCaloCellsQA::InitCellsQA", "Wrong detector provided");
73
74   *fOutfile = fname;
75 }
76
77 //________________________________________________________________
78 void AliAnalysisTaskCaloCellsQA::UserCreateOutputObjects()
79 {
80   // Per run histograms cannot be initialized here
81
82   fCellsQA->InitSummaryHistograms();
83 }
84
85 //________________________________________________________________
86 void AliAnalysisTaskCaloCellsQA::UserExec(Option_t *)
87 {
88   // Does the job for one event
89
90   // event
91   AliVEvent *event = InputEvent();
92   if (!event) {
93     Warning("AliAnalysisTaskCaloCellsQA::UserExec", "Could not get event");
94     return;
95   }
96
97   // pileup;  FIXME: why AliVEvent does not allow a version without arguments?
98   if (fkAvoidPileup && event->IsPileupFromSPD(3,0.8,3.,2.,5.))
99     return;
100
101   // cells
102   AliVCaloCells *cells;
103   if (fCellsQA->GetDetector() == AliCaloCellsQA::kEMCAL)
104     cells = event->GetEMCALCells();
105   else
106     cells = event->GetPHOSCells();
107
108   if (!cells) {
109     Warning("AliAnalysisTaskCaloCellsQA::UserExec", "Could not get cells");
110     return;
111   }
112
113   // primary vertex
114   AliVVertex *vertex = (AliVVertex*) event->GetPrimaryVertex();
115   if (!vertex) {
116     Warning("AliAnalysisTaskCaloCellsQA::UserExec", "Could not get primary vertex");
117     return;
118   }
119
120   // collect clusters
121   TObjArray clusArray;
122   for (Int_t i = 0; i < event->GetNumberOfCaloClusters(); i++) {
123     AliVCluster *clus = event->GetCaloCluster(i);
124     if (!clus) {
125       Warning("AliAnalysisTaskCaloCellsQA::UserExec", "Could not get cluster");
126       return;
127     }
128
129     // basic filtering
130     if (fCellsQA->GetDetector() == AliCaloCellsQA::kEMCAL && !clus->IsEMCAL()) continue;
131     if (fCellsQA->GetDetector() == AliCaloCellsQA::kPHOS && !clus->IsPHOS()) continue;
132     if (IsClusterBad(clus)) continue;
133
134     clusArray.Add(clus);
135   }
136
137   Double_t vertexXYZ[3];
138   vertex->GetXYZ(vertexXYZ);
139   fCellsQA->Fill(event->GetRunNumber(), &clusArray, cells, vertexXYZ);
140 }
141
142 //________________________________________________________________
143 void AliAnalysisTaskCaloCellsQA::Terminate(Option_t*)
144 {
145   // The AliCaloCellsQA analysis output should not be saved through
146   // AliAnalysisManager, because it will write with TObject::kSingleKey
147   // option. Such an output cannot be merged later: we have histograms
148   // which are run-dependent, i.e. not the same between every analysis
149   // instance.
150
151   TFile f(fOutfile->Data(), "RECREATE");
152   fCellsQA->GetListOfHistos()->Write();
153 }
154
155 //____________________________________________________________
156 void AliAnalysisTaskCaloCellsQA::SetBadCells(Int_t badcells[], Int_t nbad)
157 {
158   // Set absId numbers for bad cells;
159   // clusters which contain a bad cell will be rejected.
160
161   // switch off bad cells, if asked
162   if (nbad <= 0) {
163     if (fBadCells) delete [] fBadCells;
164     fNBad = 0;
165     return;
166   }
167
168   fNBad = nbad;
169   fBadCells = new Int_t[nbad];
170
171   for (Int_t i = 0; i < nbad; i++)
172     fBadCells[i] = badcells[i];
173 }
174
175 //________________________________________________________________
176 Bool_t AliAnalysisTaskCaloCellsQA::IsClusterBad(AliVCluster *clus)
177 {
178   // Returns true if cluster contains a bad cell
179
180   for (Int_t b = 0; b < fNBad; b++)
181     for (Int_t c = 0; c < clus->GetNCells(); c++)
182       if (clus->GetCellAbsId(c) == fBadCells[b])
183         return kTRUE;
184
185   return kFALSE;
186 }