]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliPWG0Helper.cxx
added function from Christian to enable branches with new ESD format
[u/mrichter/AliRoot.git] / PWG0 / AliPWG0Helper.cxx
1 /* $Id$ */
2
3 #include <AliPWG0Helper.h>
4
5 #include <TParticle.h>
6 #include <TParticlePDG.h>
7 #include <TH1.h>
8 #include <TH3.h>
9 #include <TList.h>
10 #include <TTree.h>
11 #include <TBranch.h>
12 #include <TLeaf.h>
13
14 #include <AliHeader.h>
15 #include <AliStack.h>
16 #include <AliLog.h>
17
18 #include <AliLog.h>
19 #include <AliESD.h>
20 #include <AliESDVertex.h>
21
22 #include <AliGenEventHeader.h>
23 #include <AliGenPythiaEventHeader.h>
24 #include <AliGenCocktailEventHeader.h>
25
26 //____________________________________________________________________
27 ClassImp(AliPWG0Helper)
28
29 //____________________________________________________________________
30 Bool_t AliPWG0Helper::IsEventTriggered(const AliESD* aEsd, Trigger trigger)
31 {
32   // see function with ULong64_t argument
33
34   ULong64_t triggerMask = aEsd->GetTriggerMask();
35   return IsEventTriggered(triggerMask, trigger);
36 }
37
38 //____________________________________________________________________
39 Bool_t AliPWG0Helper::IsEventTriggered(ULong64_t triggerMask, Trigger trigger)
40 {
41   // check if the event was triggered
42   //
43   // this function needs the branch fTriggerMask
44   //
45   // MB should be
46   // ITS_SPD_GFO_L0  : 32
47   // VZERO_OR_LEFT   : 1
48   // VZERO_OR_RIGHT  : 2
49
50   switch (trigger)
51   {
52     case kMB1:
53     {
54       if (triggerMask&32 || ((triggerMask&1) || (triggerMask&2)))
55         return kTRUE;
56       break;
57     }
58     case kMB2:
59     {
60       if (triggerMask&32 && ((triggerMask&1) || (triggerMask&2)))
61         return kTRUE;
62       break;
63     }
64   }
65
66   return kFALSE;
67 }
68
69 //____________________________________________________________________
70 Bool_t AliPWG0Helper::IsVertexReconstructed(const AliESD* aEsd)
71 {
72   // see function with AliESDVertex argument
73
74   const AliESDVertex* vtxESD = aEsd->GetVertex();
75   return IsVertexReconstructed(vtxESD);
76 }
77
78 //____________________________________________________________________
79 Bool_t AliPWG0Helper::IsVertexReconstructed(const AliESDVertex* vtxESD)
80 {
81   // checks if the vertex is reasonable
82   //
83   // this function needs the branches fSPDVertex*
84
85   if (!vtxESD)
86     return kFALSE;
87
88   // the vertex should be reconstructed
89   if (strcmp(vtxESD->GetName(), "default")==0)
90     return kFALSE;
91
92   Double_t vtx_res[3];
93   vtx_res[0] = vtxESD->GetXRes();
94   vtx_res[1] = vtxESD->GetYRes();
95   vtx_res[2] = vtxESD->GetZRes();
96
97   if (vtx_res[2]==0 || vtx_res[2]>0.1)
98     return kFALSE;
99
100   // check Ncontributors, if <0 it means error *gna*
101
102   return kTRUE;
103 }
104
105 //____________________________________________________________________
106 Bool_t AliPWG0Helper::IsPrimaryCharged(TParticle* aParticle, Int_t aTotalPrimaries, Bool_t adebug)
107 {
108   //
109   // this function checks if a particle from the event generator (i.e. among the nPrim particles in the stack)
110   // shall be counted as a primary particle
111   //
112   // This function or a equivalent should be available in some common place of AliRoot
113   //
114
115   // if the particle has a daughter primary, we do not want to count it
116   if (aParticle->GetFirstDaughter() != -1 && aParticle->GetFirstDaughter() < aTotalPrimaries)
117   {
118     if (adebug)
119       printf("Dropping particle because it has a daughter among the primaries.\n");
120     return kFALSE;
121   }
122
123   Int_t pdgCode = TMath::Abs(aParticle->GetPdgCode());
124
125   // skip quarks and gluon
126   if (pdgCode <= 10 || pdgCode == 21)
127   {
128     if (adebug)
129       printf("Dropping particle because it is a quark or gluon.\n");
130     return kFALSE;
131   }
132
133   if (strcmp(aParticle->GetName(),"XXX") == 0)
134   {
135     if (adebug)
136       printf("WARNING: There is a particle named XXX.\n");
137     return kFALSE;
138   }
139
140   TParticlePDG* pdgPart = aParticle->GetPDG();
141
142   if (strcmp(pdgPart->ParticleClass(),"Unknown") == 0)
143   {
144     if (adebug)
145       printf("WARNING: There is a particle with an unknown particle class (pdg code %d).\n", pdgCode);
146     return kFALSE;
147   }
148
149   if (pdgPart->Charge() == 0)
150   {
151     if (adebug)
152       printf("Dropping particle because it is not charged.\n");
153     return kFALSE;
154   }
155
156   return kTRUE;
157 }
158
159 //____________________________________________________________________
160 void AliPWG0Helper::CreateProjections(TH3* hist, Bool_t save)
161 {
162   // create projections of 3d hists to all 2d combinations
163   // the histograms are not returned, just use them from memory or use this to create them in a file
164
165   TH1* proj = hist->Project3D("yx");
166   proj->SetXTitle(hist->GetXaxis()->GetTitle());
167   proj->SetYTitle(hist->GetYaxis()->GetTitle());
168   if (save)
169     proj->Write();
170
171   proj = hist->Project3D("zx");
172   proj->SetXTitle(hist->GetXaxis()->GetTitle());
173   proj->SetYTitle(hist->GetZaxis()->GetTitle());
174   if (save)
175     proj->Write();
176
177   proj = hist->Project3D("zy");
178   proj->SetXTitle(hist->GetYaxis()->GetTitle());
179   proj->SetYTitle(hist->GetZaxis()->GetTitle());
180   if (save)
181     proj->Write();
182 }
183
184 //____________________________________________________________________
185 void AliPWG0Helper::CreateDividedProjections(TH3* hist, TH3* hist2, const char* axis, Bool_t putErrors, Bool_t save)
186 {
187   // create projections of the 3d hists divides them
188   // axis decides to which plane, if axis is 0 to all planes
189   // the histograms are not returned, just use them from memory or use this to create them in a file
190
191   if (axis == 0)
192   {
193     CreateDividedProjections(hist, hist2, "yx", putErrors, save);
194     CreateDividedProjections(hist, hist2, "zx", putErrors, save);
195     CreateDividedProjections(hist, hist2, "zy", putErrors, save);
196
197     return;
198   }
199
200   TH1* proj = hist->Project3D(axis);
201
202   if (strlen(axis) == 2)
203   {
204     proj->SetYTitle(GetAxisTitle(hist, axis[0]));
205     proj->SetXTitle(GetAxisTitle(hist, axis[1]));
206   }
207   else if (strlen(axis) == 1)
208     proj->SetXTitle(GetAxisTitle(hist, axis[0]));
209
210   TH1* proj2 = hist2->Project3D(axis);
211   if (strlen(axis) == 2)
212   {
213     proj2->SetYTitle(GetAxisTitle(hist2, axis[0]));
214     proj2->SetXTitle(GetAxisTitle(hist2, axis[1]));
215   }
216   else if (strlen(axis) == 1)
217     proj2->SetXTitle(GetAxisTitle(hist2, axis[0]));
218
219   TH1* division = dynamic_cast<TH1*> (proj->Clone(Form("%s_div_%s", proj->GetName(), proj2->GetName())));
220   //printf("doing axis: %s, x axis has %d %d bins, min %f %f max %f %f\n", axis, division->GetNbinsX(), proj2->GetNbinsX(), division->GetXaxis()->GetBinLowEdge(1), proj2->GetXaxis()->GetBinLowEdge(1), division->GetXaxis()->GetBinUpEdge(division->GetNbinsX()), proj2->GetXaxis()->GetBinUpEdge(proj2->GetNbinsX()));
221   //printf("doing axis: %s, y axis has %d %d bins, min %f %f max %f %f\n", axis, division->GetNbinsY(), proj2->GetNbinsY(), division->GetYaxis()->GetBinLowEdge(1), proj2->GetYaxis()->GetBinLowEdge(1), division->GetYaxis()->GetBinUpEdge(division->GetNbinsY()), proj2->GetYaxis()->GetBinUpEdge(proj2->GetNbinsY()));
222   division->Divide(proj, proj2, 1, 1, "B");
223   division->SetTitle(Form("%s divided %s", proj->GetTitle(), proj2->GetTitle()));
224
225   if (putErrors)
226   {
227     division->Sumw2();
228     if (division->GetDimension() == 1)
229     {
230       Int_t nBins = division->GetNbinsX();
231       for (Int_t i = 1; i <= nBins; ++i)
232         if (proj2->GetBinContent(i) != 0)
233           division->SetBinError(i, TMath::Sqrt(proj->GetBinContent(i)) / proj2->GetBinContent(i));
234     }
235     else if (division->GetDimension() == 2)
236     {
237       Int_t nBinsX = division->GetNbinsX();
238       Int_t nBinsY = division->GetNbinsY();
239       for (Int_t i = 1; i <= nBinsX; ++i)
240         for (Int_t j = 1; j <= nBinsY; ++j)
241           if (proj2->GetBinContent(i, j) != 0)
242             division->SetBinError(i, j, TMath::Sqrt(proj->GetBinContent(i, j)) / proj2->GetBinContent(i, j));
243     }
244   }
245
246   if (save)
247   {
248     proj->Write();
249     proj2->Write();
250     division->Write();
251   }
252 }
253
254 //____________________________________________________________________
255 const char* AliPWG0Helper::GetAxisTitle(TH3* hist, const char axis)
256 {
257   // returns the title of the axis given in axis (x, y, z)
258
259   if (axis == 'x')
260     return hist->GetXaxis()->GetTitle();
261   else if (axis == 'y')
262     return hist->GetYaxis()->GetTitle();
263   else if (axis == 'z')
264     return hist->GetZaxis()->GetTitle();
265
266   return 0;
267 }
268
269 //____________________________________________________________________
270 Int_t AliPWG0Helper::GetPythiaEventProcessType(AliHeader* aHeader, Bool_t adebug) {
271   //
272   // get the process type of the event.
273   //
274
275   // can only read pythia headers, either directly or from cocktalil header
276   AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(aHeader->GenEventHeader());
277
278   if (!pythiaGenHeader) {
279
280     AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(aHeader->GenEventHeader());
281     if (!genCocktailHeader) {
282       printf("AliPWG0Helper::GetProcessType : Unknown header type (not Pythia or Cocktail). \n");
283       return -1;
284     }
285
286     TList* headerList = genCocktailHeader->GetHeaders();
287     if (!headerList) {
288       return -1;
289     }
290
291     for (Int_t i=0; i<headerList->GetEntries(); i++) {
292       pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
293       if (pythiaGenHeader)
294         break;
295     }
296
297     if (!pythiaGenHeader) {
298       printf("AliPWG0Helper::GetProcessType : Could not find Pythia header. \n");
299       return -1;
300     }
301   }
302
303   if (adebug) {
304     printf("AliPWG0Helper::GetProcessType : Pythia process type found: %d \n",pythiaGenHeader->ProcessType());
305   }
306
307   return pythiaGenHeader->ProcessType();
308 }
309
310 //____________________________________________________________________
311 TParticle* AliPWG0Helper::FindPrimaryMother(AliStack* stack, Int_t label)
312 {
313   //
314   // Finds the first mother among the primary particles of the particle identified by <label>,
315   // i.e. the primary that "caused" this particle
316   //
317
318   Int_t motherLabel = FindPrimaryMotherLabel(stack, label);
319   if (motherLabel < 0)
320     return 0;
321
322   return stack->Particle(motherLabel);
323 }
324
325 //____________________________________________________________________
326 Int_t AliPWG0Helper::FindPrimaryMotherLabel(AliStack* stack, Int_t label)
327 {
328   //
329   // Finds the first mother among the primary particles of the particle identified by <label>,
330   // i.e. the primary that "caused" this particle
331   //
332   // returns its label
333   //
334
335   Int_t nPrim  = stack->GetNprimary();
336
337   while (label >= nPrim)
338   {
339     //printf("Particle %d (pdg %d) is not a primary. Let's check its mother %d\n", label, mother->GetPdgCode(), mother->GetMother(0));
340
341     TParticle* particle = stack->Particle(label);
342     if (!particle)
343     {
344       AliDebugGeneral("FindPrimaryMother", AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack.", label));
345       return -1;
346     }
347  
348     // find mother
349     if (particle->GetMother(0) < 0)
350     {
351       AliDebugGeneral("FindPrimaryMother", AliLog::kError, Form("UNEXPECTED: Could not find mother of secondary particle %d.", label));
352       return -1;
353     }
354
355     label = particle->GetMother(0);
356   }
357
358   return label;
359 }
360
361 void AliPWG0Helper::SetBranchStatusRecursive(TTree* tree, char *bname, Bool_t status, Bool_t debug)
362 {
363   // Function  to switch on/off all data members of a top level branch
364   // this is needed for branches without a trailing dot ".", for those
365   // the root functionality with regular expressions does not work.
366   // Usage e.g.
367   // chain->SetBranchStatus("*", 0); 
368   // SetBranchStatusRecursive(chain,"SPDVertex",1);
369   // You need to give the full name of the top level branch zou want to access
370   //==========================================================================
371   // Author Christian.Klein-Boesing@cern.ch
372
373   if (!tree)
374     return;
375
376   TBranch *br = tree->GetBranch(bname);
377   if(!br) {
378     Printf("AliPWG0Helper::SetBranchStatusRecursive: Branch %s not found", bname);
379   }
380
381   TObjArray *leaves = tree->GetListOfLeaves();
382   Int_t nleaves = leaves->GetEntries();
383   TLeaf *leaf = 0;
384   TBranch *branch = 0;
385   TBranch *mother = 0;
386   for (Int_t i=0;i<nleaves;i++)  {
387     // the matched entry e.g. SPDVertex is its own Mother
388     leaf = (TLeaf*)leaves->UncheckedAt(i);
389     branch = (TBranch*)leaf->GetBranch();
390     mother = branch->GetMother();
391     if (mother==br){
392       if (debug)
393         Printf(">>>> AliPWG0Helper::SetBranchStatusRecursive: Setting status %s %s to %d", mother->GetName(), leaf->GetName(), status);
394       if (status) branch->ResetBit(kDoNotProcess);
395       else        branch->SetBit(kDoNotProcess);
396     }
397   }
398 }