]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliDielectronCutQA.cxx
update from pr task : sjena
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliDielectronCutQA.cxx
1 /*************************************************************************
2  * Copyright(c) 1998-2009, 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 //                           CutQA                                      //
18 //                                                                      //
19 /*
20    Allow to monitor how many tracks,pair,events pass the selection criterion 
21    in any of the cuts added to the corresponding filters. All you need to 
22    add to your config is the following:
23
24    dielectron->SetCutQA();
25
26
27 */
28 //                                                                      //
29 //////////////////////////////////////////////////////////////////////////
30
31 #include "AliDielectronCutQA.h"
32
33 #include <TList.h>
34 #include <TCollection.h>
35
36 #include "AliDielectronCutGroup.h"
37 #include "AliAnalysisCuts.h"
38
39 #include "AliVEvent.h"
40 #include "AliVParticle.h"
41 #include "AliVTrack.h"
42 #include "AliDielectronPair.h"
43
44
45 ClassImp(AliDielectronCutQA)
46
47
48 AliDielectronCutQA::AliDielectronCutQA() :
49   TNamed(),
50   fQAHistArray()
51 {
52   //
53   // Default constructor
54   //
55   for(Int_t itype=0; itype<kNtypes; itype++) {
56     fCutQA[itype]=0x0;
57     fNCuts[itype]=1;
58     for(Int_t i=0; i<20; i++) {
59       fCutNames[i][itype]="";
60     }
61   }
62   fTypeKeys[kTrack] = "Track";
63   fTypeKeys[kPair]  = "Pair";
64   fTypeKeys[kEvent] = "Event";
65   fQAHistArray.SetOwner();
66 }
67
68 //_____________________________________________________________________
69 AliDielectronCutQA::AliDielectronCutQA(const char* name, const char* title) :
70   TNamed(name, title),
71   fQAHistArray()
72 {
73   //
74   // Named Constructor
75   //
76   for(Int_t itype=0; itype<kNtypes; itype++) {
77     fCutQA[itype]=0x0;
78     fNCuts[itype]=1;
79     for(Int_t i=0; i<20; i++) {
80       fCutNames[i][itype]="";
81     }
82   }
83   fTypeKeys[kTrack] = "Track";
84   fTypeKeys[kPair]  = "Pair";
85   fTypeKeys[kEvent] = "Event";
86   fQAHistArray.SetOwner();
87 }
88
89 //_____________________________________________________________________
90 AliDielectronCutQA::~AliDielectronCutQA() 
91 {
92   //
93   //Default Destructor
94   //
95   fQAHistArray.Delete();
96 }
97
98 //_____________________________________________________________________
99 void AliDielectronCutQA::Init()
100 {
101
102   fQAHistArray.SetName(Form("%s",GetName()));
103
104   // loop over all types
105   for(Int_t itype=0; itype<kNtypes; itype++) {
106     //    printf("\n type: %d\n",itype);
107     fCutNames[0][itype]="no cuts";
108
109     // create histogram based on added cuts
110     fCutQA[itype] = new TH1F(fTypeKeys[itype],
111                              Form("%sQA;cuts;# passed %ss",fTypeKeys[itype],fTypeKeys[itype]),
112                              fNCuts[itype],0,fNCuts[itype]);
113     // loop over all cuts
114     for(Int_t i=0; i<fNCuts[itype]; i++) {
115       fCutQA[itype]->GetXaxis()->SetBinLabel(i+1,fCutNames[i][itype]);
116       //      printf(" %s \n",fCutNames[i][itype]);
117     }
118     fQAHistArray.AddLast(fCutQA[itype]);
119   }
120
121 }
122
123 //_____________________________________________________________________
124 void AliDielectronCutQA::AddTrackFilter(AliAnalysisFilter *trackFilter)
125 {
126   //
127   // add track filter cuts to the qa histogram
128   //
129
130
131   TIter listIterator(trackFilter->GetCuts());
132   while (AliAnalysisCuts *thisCut = (AliAnalysisCuts*) listIterator()) {
133     Bool_t addCut=kTRUE;
134
135     // add new cut class to the array
136     if(addCut) {
137       fCutNames[fNCuts[kTrack]][kTrack]=thisCut->GetTitle();
138       //      printf("add cut %s to %d \n",thisCut->GetTitle(),fNCuts[kTrack]);
139       fNCuts[kTrack]++;
140     }
141
142   } // pair filter loop
143
144 }
145
146
147 //_____________________________________________________________________
148 void AliDielectronCutQA::AddPairFilter(AliAnalysisFilter *pairFilter)
149 {
150   //
151   // add track filter cuts to the qa histogram
152   //
153
154
155   TIter listIterator(pairFilter->GetCuts());
156   while (AliAnalysisCuts *thisCut = (AliAnalysisCuts*) listIterator()) {
157     Bool_t addCut=kTRUE;
158
159     // add new cut class to the array
160     if(addCut) {
161       fCutNames[fNCuts[kPair]][kPair]=thisCut->GetTitle();
162       //  printf("add cut %s to %d \n",thisCut->GetTitle(),fNCuts[kPair]);
163       fNCuts[kPair]++;
164     }
165
166   } // trk filter loop
167
168 }
169
170 //_____________________________________________________________________
171 void AliDielectronCutQA::AddEventFilter(AliAnalysisFilter *eventFilter)
172 {
173   //
174   // add track filter cuts to the qa histogram
175   //
176
177
178   TIter listIterator(eventFilter->GetCuts());
179   while (AliAnalysisCuts *thisCut = (AliAnalysisCuts*) listIterator()) {
180     Bool_t addCut=kTRUE;
181
182     // add new cut class to the array
183     if(addCut) {
184       fCutNames[fNCuts[kEvent]][kEvent]=thisCut->GetTitle();
185       //      printf("add cut %s to %d \n",thisCut->GetTitle(),fNCuts[kEvent]);
186       fNCuts[kEvent]++;
187     }
188
189   } // trk filter loop
190
191 }
192
193 //_____________________________________________________________________
194 void AliDielectronCutQA::Fill(UInt_t mask, TObject *obj)
195 {
196   //
197   // fill the corresponding step in the qa histogram
198   //
199
200   UInt_t idx = GetObjIndex(obj);
201
202   Int_t cutstep=1;
203   for (Int_t iCut=0; iCut<fNCuts[idx]-1;++iCut) {
204     //    UInt_t cutMask=1<<iCut;         // for each cut
205     UInt_t cutMask=(1<<(iCut+1))-1; // increasing cut match
206
207     if ((mask&cutMask)==cutMask) {
208       fCutQA[idx]->Fill(cutstep);
209       ++cutstep;
210     }
211
212   }
213
214 }
215
216 //_____________________________________________________________________
217 void AliDielectronCutQA::FillAll(TObject *obj)
218 {
219   //
220   // fill the corresponding step in the qa histogram
221   //
222
223   UInt_t idx = GetObjIndex(obj);
224   fCutQA[idx]->Fill(0);
225
226 }
227
228 //______________________________________________________________________
229 UInt_t AliDielectronCutQA::GetObjIndex(TObject *obj)
230 {
231   //
232   // return the corresponding idex
233   //
234   //  printf("INFO: object type is a %s \n", obj->IsA()->GetName());
235   if(obj->InheritsFrom(AliDielectronPair::Class())) return kPair;
236   if(obj->InheritsFrom(AliVTrack::Class())        ) return kTrack;
237   if(obj->InheritsFrom(AliVParticle::Class())     ) return kTrack;
238   if(obj->InheritsFrom(AliVEvent::Class())        ) return kEvent;
239   printf("FATAL: object type %s not yet supported, please let the author know\n", obj->IsA()->GetName());
240   return -1;
241   //TODO complete
242
243 }
244
245
246
247