]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliMuonsHFHeader.cxx
993cfd4781c0c6cc54d65f8210dd07d50b0bb412
[u/mrichter/AliRoot.git] / PWG3 / muon / AliMuonsHFHeader.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2006, 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 /* $Id$ */
17
18 /////////////////////////////////////////////////////////////
19 //
20 // class used to extract and store info at event level
21 //
22 // Author: X-M. Zhang, zhang@clermont.in2p3.fr
23 //                     zhangxm@iopp.ccnu.edu.cn
24 /////////////////////////////////////////////////////////////
25
26 #include <TMath.h>
27 #include <TH1.h>
28 #include <TH1D.h>
29 #include <TList.h>
30
31 #include "AliVEvent.h"
32 #include "AliMuonInfoStoreRD.h"
33 #include "AliMuonInfoStoreMC.h"
34 #include "AliDimuInfoStoreRD.h"
35 #include "AliDimuInfoStoreMC.h"
36 #include "AliMuonsHFHeader.h"
37
38 class TNamed;
39 class AliVVertex;
40
41 ClassImp(AliMuonsHFHeader)
42
43 const TString AliMuonsHFHeader::fgkStdBranchName("MuEvsH");
44 Int_t         AliMuonsHFHeader::fgAnaMode = 0;
45 Bool_t        AliMuonsHFHeader::fgIsMC    = kFALSE;
46 Double_t      AliMuonsHFHeader::fgCuts[5] = { -999999., 999999., 999999., -999999., 999999. };
47
48 //_____________________________________________________________________________
49 AliMuonsHFHeader::AliMuonsHFHeader() :
50 TNamed(),
51 fSelMask(AliVEvent::kAny),
52 fIsMB(kFALSE),
53 fIsMU(kFALSE),
54 fVtxContrsN(0),
55 fFiredTriggerClass(),
56 fCentrality(0.)
57 {
58   //
59   // default constructor
60   //
61   for (Int_t i=3; i--;) fVtx[i] = 0.;
62 }
63
64 //_____________________________________________________________________________
65 AliMuonsHFHeader::AliMuonsHFHeader(const AliMuonsHFHeader &src) :
66 TNamed(),
67 fSelMask(src.fSelMask),
68 fIsMB(src.fIsMB),
69 fIsMU(src.fIsMU),
70 fVtxContrsN(src.fVtxContrsN),
71 fFiredTriggerClass(src.fFiredTriggerClass),
72 fCentrality(src.fCentrality)
73 {
74   //
75   // copy constructor
76   //
77   for (Int_t i=3; i--;) fVtx[i] = src.fVtx[i];
78 }
79
80 //_____________________________________________________________________________
81 AliMuonsHFHeader& AliMuonsHFHeader::operator=(const AliMuonsHFHeader &src)
82 {
83   //
84   // assignment constructor
85   //
86
87   if(&src==this) return *this;
88
89   fSelMask           = src.fSelMask;
90   fIsMB              = src.fIsMB;
91   fIsMU              = src.fIsMU;
92   fVtxContrsN        = src.fVtxContrsN;
93   fFiredTriggerClass = src.fFiredTriggerClass;
94   fCentrality        = src.fCentrality;
95   for (Int_t i=3; i--;) fVtx[i] = src.fVtx[i];
96
97   return *this;
98 }
99
100 //_____________________________________________________________________________
101 AliMuonsHFHeader::~AliMuonsHFHeader()
102 {
103   //
104   // default destructor
105   //
106 }
107
108 //_____________________________________________________________________________
109 void AliMuonsHFHeader::SetVertex(AliVVertex *vertex)
110 {
111   // extract event info from AOD event
112
113   vertex->GetXYZ(fVtx);
114   fVtxContrsN = vertex->GetNContributors();
115   this->SetTitle(vertex->GetTitle());
116   return;
117 }
118
119 void AliMuonsHFHeader::SetFiredTriggerClass(TString trigger)
120 {
121   fFiredTriggerClass = trigger;
122   fIsMB=(fFiredTriggerClass.Contains("CINT1B-ABCE-NOPF-ALL")  ||   // p-p   min. bias trigger before Aug/2010
123          fFiredTriggerClass.Contains("CINT1-B-NOPF-ALLNOTRD") ||   // p-p   min. bias trigger from   Aug/2010
124          fFiredTriggerClass.Contains("CMBS2A-B-NOPF-ALL")     ||   // Pb-Pb min. bias trigger 2-out-of-3
125          fFiredTriggerClass.Contains("CMBS2C-B-NOPF-ALL")     ||   // Pb-Pb min. bias trigger 2-out-of-3
126          fFiredTriggerClass.Contains("CMBAC-B-NOPF-ALL")      ||   // Pb-Pb min. bias trigger 2-out-of-3
127          fFiredTriggerClass.Contains("CMBACS2-B-NOPF-ALL")    ||   // Pb-Pb min. bias trigger 3-out-of-3 (early)
128          fFiredTriggerClass.Contains("CMBACS2-B-NOPF-ALLNOTRD"));  // Pb-Pb min. bias trigger 3-out-of-3 (late 2010)
129   fIsMU=(fFiredTriggerClass.Contains("CMUS1B-ABCE-NOPF-MUON") ||   // p-p MUON trigger before Aug/2010
130          fFiredTriggerClass.Contains("CMUS1-B-NOPF-ALLNOTRD"));    // p-p MUON trigger from   Aug/2010 
131   return;
132 }
133
134 //_____________________________________________________________________________
135 Bool_t AliMuonsHFHeader::IsSelected()
136 {
137   // select event according to the event selection cuts
138   if (this->VtxContrsN()<fgCuts[0])       return kFALSE;
139   if (TMath::Abs(this->Vz())>fgCuts[1])   return kFALSE;
140   if (this->Vt()>fgCuts[2])               return kFALSE;
141
142   // centrality selection
143   Float_t centr = this->Centrality();
144   if (centr<fgCuts[3] || centr>fgCuts[4]) return kFALSE;
145   return kTRUE;
146 }
147
148 //_____________________________________________________________________________
149 void AliMuonsHFHeader::CreateHistograms(TList *list)
150 {
151   // create output histos of muon analysis according to the analysis mode & MC flag
152
153   if (fgIsMC) {
154     this->CreateHistosEvnH(list);
155     if (fgAnaMode!=2) {
156       TString sName[7] = { "Unidentified", "Hadron", "SecondaryMu", "PrimaryMu", "CharmMu", "BottomMu", "" };
157       for (Int_t i=7; i--;) this->CreateHistosMuon(list, sName[i]);
158     }
159     if (fgAnaMode!=1) {
160       TString sName[7] = { "Uncorr", "Resonance", "DDsame", "DDdiff", "BBsame", "BBdiff", "" };
161       for (Int_t i=7; i--;) this->CreateHistosDimu(list, sName[i]);
162     }
163     return;
164   }
165
166   this->CreateHistosEvnH(list,"MB"); this->CreateHistosEvnH(list,"MU");
167   if (fgAnaMode!=2) { this->CreateHistosMuon(list,"MB"); this->CreateHistosMuon(list,"MU"); }
168   if (fgAnaMode!=1) { this->CreateHistosDimu(list,"MB"); this->CreateHistosDimu(list,"MU"); }
169   return;
170 }
171
172 //_____________________________________________________________________________
173 void AliMuonsHFHeader::CreateHistosEvnH(TList *list, TString sName)
174 {
175   // create histograms at event level
176
177   if (!list) list = new TList();
178   list->SetOwner();
179   Bool_t oldStatus = TH1::AddDirectoryStatus();
180   TH1::AddDirectory(kFALSE);
181
182   const Int_t nhs    = 3;
183   TString tName[nhs] = {  "Vz",  "Vt",  "VtxNcontr" };
184   Int_t   nbins[nhs] = {  800 ,   40 ,        202   };
185   Double_t xlow[nhs] = {  -40.,    0.,         -2.5 };
186   Double_t  xup[nhs] = {   40.,    4.,        199.5 };
187
188   TH1D *histo = 0;
189   for (Int_t i=0; i<nhs; i++) {
190     char *hName = Form("h%s_%s", sName.Data(), tName[i].Data());
191     histo = new TH1D(hName, hName, nbins[i], xlow[i], xup[i]);
192     histo->Sumw2(); list->Add(histo); histo = 0;
193   }
194
195   TH1::AddDirectory(oldStatus);
196   return;
197 }
198
199 //_____________________________________________________________________________
200 void AliMuonsHFHeader::CreateHistosMuon(TList *list, TString sName)
201 {
202   // create histograms for single muon
203
204   if (!list) list = new TList();
205   list->SetOwner();
206   Bool_t oldStatus = TH1::AddDirectoryStatus();
207   TH1::AddDirectory(kFALSE);
208
209   const Int_t nhs    = 7;
210   TString tName[nhs] = {   "P",  "Pt",  "Eta",  "DCA",  "TrM",  "Charge", "Rabs" };
211   Int_t   nbins[nhs] = { 1500 ,  300 ,    15 ,  1000 ,    4  ,       3  ,    48  };
212   Double_t xlow[nhs] = {    0.,    0.,   -4.0,     0.,   -0.5,      -1.5,   17.6 };
213   Double_t  xup[nhs] = {  150.,   30.,   -2.5,   500.,    3.5,       1.5,   80.0 };
214
215   TH1D *histo = 0;
216   for (Int_t i=0; i<nhs; i++) {
217     char *hName = Form("h%s_%s", sName.Data(), tName[i].Data());
218     histo = new TH1D(hName, hName, nbins[i], xlow[i], xup[i]);
219     histo->Sumw2(); list->Add(histo); histo = 0;
220   }
221
222   TH1::AddDirectory(oldStatus);
223   return;
224 }
225
226 //_____________________________________________________________________________
227 void AliMuonsHFHeader::CreateHistosDimu(TList *list, TString sName)
228 {
229   // create histograms for dimuon
230
231   if (!list) list = new TList();
232   list->SetOwner();
233   Bool_t oldStatus = TH1::AddDirectoryStatus();
234   TH1::AddDirectory(kFALSE);
235
236   TH1D *histo = 0;
237   const Int_t nhs    = 3;
238   TString tName[nhs] = {   "P",  "Pt",  "InvM"   };
239   Int_t   nbins[nhs] = { 1500 ,  300 ,    300    };
240   Double_t xlow[nhs] = {    0.,    0.,      0.   };
241   Double_t  xup[nhs] = {  150.,   30.,     30.   };
242   TString dimuName[3] = { "DimuNN", "DimuNP", "DimuPP" };
243   for (Int_t i=0; i<3; i++) {
244     for (Int_t j=0; j<nhs; j++) {
245       char *hName = Form("h%s_%s_%s", sName.Data(), dimuName[i].Data(), tName[j].Data());
246       histo = new TH1D(hName, hName, nbins[j], xlow[j], xup[j]);
247       histo->Sumw2(); list->Add(histo); histo = 0;
248     }
249   }
250
251   TH1::AddDirectory(oldStatus);
252   return;
253 }
254
255 //_____________________________________________________________________________
256 void AliMuonsHFHeader::FillHistosEvnH(TList *list)
257 {
258   // fill histograms at event level according to event selection cuts
259
260   if (!list)                   return;
261   if (!this->IsSelected()) return;
262
263   const Int_t nhs    = 3;
264   TString tName[nhs] = {       "Vz",       "Vt",        "VtxNcontr" };
265   Double_t dist[nhs] = { this->Vz(), this->Vt(), this->VtxContrsN() };
266   if (fgIsMC && (fSelMask & AliVEvent::kAny)) {
267     for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h_%s",tName[i].Data())))->Fill(dist[i]);
268   } else {
269     if (fIsMB && (fSelMask & AliVEvent::kMB))   { for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s","MB",tName[i].Data())))->Fill(dist[i]); }
270     if (fIsMU && (fSelMask & AliVEvent::kMUON)) { for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s","MU",tName[i].Data())))->Fill(dist[i]); }
271   }
272   return;
273 }
274
275 //_____________________________________________________________________________
276 void AliMuonsHFHeader::FillHistosMuon(TList *list, AliMuonInfoStoreRD* const infoStore, Int_t s)
277 {
278   // fill histograms for single muon according to event & muon track selection cuts
279
280   if (!list)                    return;
281   if (!this->IsSelected())      return;
282   if (!infoStore->IsSelected()) return;
283
284   const Int_t nhs    = 7;
285   TString tName[nhs] = { "P", "Pt", "Eta", "DCA", "TrM", "Charge", "Rabs" };
286   Double_t dist[nhs] = { infoStore->Momentum().Mag(),
287                              infoStore->Momentum().Pt(),
288                              infoStore->Momentum().Eta(),
289                              infoStore->DCA(),
290                              infoStore->MatchTrigger(),
291                              infoStore->Charge(),
292                              infoStore->RabsEnd() };
293
294   if (fgIsMC && (fSelMask & AliVEvent::kAny)) {
295     TString sName[7] = { "BottomMu", "CharmMu", "PrimaryMu", "SecondaryMu", "Hadron", "Unidentified", "" };
296     for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s",sName[6].Data(),tName[i].Data())))->Fill(dist[i]);
297     for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s",sName[s].Data(),tName[i].Data())))->Fill(dist[i]);
298   } else {
299     if (fIsMB && (fSelMask & AliVEvent::kMB))   { for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s","MB",tName[i].Data())))->Fill(dist[i]); }
300     if (fIsMU && (fSelMask & AliVEvent::kMUON)) { for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s","MU",tName[i].Data())))->Fill(dist[i]); }
301   }
302
303   return; 
304 }
305
306 //_____________________________________________________________________________
307 void AliMuonsHFHeader::FillHistosDimu(TList *list, AliDimuInfoStoreRD* const infoStore, Int_t s)
308 {
309   // fill histograms for dimuon according to evnet & dimuon candidates selection cuts
310
311   if (!list)                    return;
312   if (!this->IsSelected())      return;
313   if (!infoStore->IsSelected()) return;
314
315   TString dimuName = "DimuNN";
316   if (infoStore->Charge()==0)     dimuName = "DimuNP";
317   else if (infoStore->Charge()>0) dimuName = "DimuPP";
318
319   const Int_t nhs    = 3;
320   TString tName[nhs] = { "P", "Pt", "InvM" };
321   Double_t dist[nhs] = { infoStore->Momentum().Mag(),
322                              infoStore->Momentum().Pt(),
323                              infoStore->InvM() };
324
325   if (fgIsMC && (fSelMask & AliVEvent::kAny)) {
326     TString sName[7] = { "BBdiff", "BBsame", "DDdiff", "DDsame", "Resonance", "Uncorr", "" };
327     for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s_%s",sName[6].Data(),dimuName.Data(),tName[i].Data())))->Fill(dist[i]);
328     for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s_%s",sName[s].Data(),dimuName.Data(),tName[i].Data())))->Fill(dist[i]);
329   } else {
330     if (fIsMB && (fSelMask & AliVEvent::kMB)) {
331       for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s_%s","MB",dimuName.Data(),tName[i].Data())))->Fill(dist[i]);
332     }
333     if (fIsMU && (fSelMask & AliVEvent::kMUON)) {
334       for (Int_t i=nhs; i--;) ((TH1D*)list->FindObject(Form("h%s_%s_%s","MU",dimuName.Data(),tName[i].Data())))->Fill(dist[i]);
335     }
336   }
337
338   return;
339 }