]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAnalysisTaskMEVertexingHF.cxx
Add the libraries needed to run the HF QA task
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / AliAnalysisTaskMEVertexingHF.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  * Permission to use, copy, modify and distribute this software and its   *
7  * documentation strictly for non-commercial purposes is hereby granted   *
8  * without fee, provided that the above copyright notice appears in all   *
9  * copies and that both the copyright notice and this permission notice   *
10  * appear in the supporting documentation. The authors make no claims     *
11  * about the suitability of this software for any purpose. It is          *
12  * provided "as is" without express or implied warranty.                  *
13  **************************************************************************/
14
15 /* $Id$ */
16
17 //*************************************************************************
18 // Class AliAnalysisTaskMEVertexingHF
19 // AliAnalysisTaskME for event mixing, building the background for 
20 // heavy-flavour decay candidates
21 // Author: R.Romita, r.romita@gsi.de
22 //*************************************************************************
23
24
25
26 #include "TH1F.h"
27 #include "TObjArray.h"
28 #include "TList.h"
29 #include "TROOT.h"
30 #include "TSystem.h"
31 #include "TCanvas.h"
32
33 #include "AliVEvent.h"
34 #include "AliVVertex.h"
35 #include "AliAODEvent.h"
36 #include "AliESDEvent.h"
37 #include "AliAnalysisVertexingHF.h"
38 #include "AliMixedEvent.h"
39 #include "AliAnalysisTaskMEVertexingHF.h"
40 #include "AliAnalysisManager.h"
41 #include "AliMultiEventInputHandler.h"
42
43 ClassImp(AliAnalysisTaskMEVertexingHF)
44
45 //________________________________________________________________________
46 AliAnalysisTaskMEVertexingHF::AliAnalysisTaskMEVertexingHF(const char *name) : 
47 AliAnalysisTaskME(name), 
48 fvHF(0), 
49 fMixedEvent(),
50 fVerticesHFTClArr(0),
51 fD0toKpiTClArr(0), 
52 fJPSItoEleTClArr(0),
53 fCharm3ProngTClArr(0),
54 fCharm4ProngTClArr(0),
55 fDstarTClArr(0),
56 fCascadesTClArr(0),
57 fLikeSign2ProngTClArr(0),
58 fLikeSign3ProngTClArr(0)
59 {
60   // Constructor
61 }
62 //________________________________________________________________________
63 void AliAnalysisTaskMEVertexingHF::Init()
64 {
65  // Initialization
66  // Instanciates vHF and loads its parameters
67  // Some parameters are changed
68  
69   if(gROOT->LoadMacro("ConfigVertexingHF.C")) {
70     printf("AnalysisTaskMEVertexingHF::Init() \n Using $ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C\n");
71     gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C");
72   }
73   fvHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
74   fvHF->SetMixEventOn();
75   fvHF->SetInputAOD();
76   fvHF->PrintStatus();
77   if(fvHF->GetLikeSign()) {
78     printf("WARNING: fLikeSign will be switched off!");
79     fvHF->SetLikeSignOff();
80   }
81   if(fvHF->GetRecoPrimVtxSkippingTrks() || fvHF->GetRmTrksFromPrimVtx()){
82     fvHF->UnsetRecoPrimVtxSkippingTrks();
83     printf("WARNING: if on, fRecoPrimVtxSkippingTrks and fRmTrksFromPrimVtx  will be switched off!\n");
84   }
85   
86   AliAnalysisManager::GetAnalysisManager()->RegisterExtraFile("AliAOD.VertexingHF.root");
87
88   return;
89 }
90 //________________________________________________________________________
91 void AliAnalysisTaskMEVertexingHF::UserCreateOutputObjects()
92 {  
93 // Create the output container
94
95
96   if (!AODEvent()) {
97     Fatal("UserCreateOutputObjects", "This task needs an AOD handler");
98     return;
99   }
100   
101   if(!fvHF) {
102     printf("AnalysisTaskMEVertexingHF::UserCreateOutPutData() \n ERROR! no fvHF!\n");
103     return;
104   }
105   fVerticesHFTClArr = new TClonesArray("AliAODVertex", 0);
106   fVerticesHFTClArr->SetName("VerticesHF");
107   AddAODBranch("TClonesArray", &fVerticesHFTClArr);
108   if(fvHF->GetD0toKpi()) {
109     fD0toKpiTClArr = new TClonesArray("AliAODRecoDecayHF2Prong", 0);
110     fD0toKpiTClArr->SetName("D0toKpi");
111     AddAODBranch("TClonesArray", &fD0toKpiTClArr);
112   }
113   if(fvHF->GetJPSItoEle()) {
114     fJPSItoEleTClArr = new TClonesArray("AliAODRecoDecayHF2Prong", 0);
115     fJPSItoEleTClArr->SetName("JPSItoEle");
116     AddAODBranch("TClonesArray", &fJPSItoEleTClArr);
117   }
118   if(fvHF->Get3Prong()) {
119     fCharm3ProngTClArr = new TClonesArray("AliAODRecoDecayHF3Prong", 0);
120     fCharm3ProngTClArr->SetName("Charm3Prong");
121     AddAODBranch("TClonesArray", &fCharm3ProngTClArr);
122   }
123   if(fvHF->Get4Prong()) {
124     fCharm4ProngTClArr = new TClonesArray("AliAODRecoDecayHF4Prong", 0);
125     fCharm4ProngTClArr->SetName("Charm4Prong");
126     AddAODBranch("TClonesArray", &fCharm4ProngTClArr);
127   }
128   
129   if(fvHF->GetDstar()) {
130     fDstarTClArr = new TClonesArray("AliAODRecoCascadeHF", 0);
131     fDstarTClArr->SetName("Dstar");
132     AddAODBranch("TClonesArray", &fDstarTClArr);
133   }
134   
135   if(fvHF->GetCascades()){
136     fCascadesTClArr = new TClonesArray("AliAODRecoCascadeHF", 0);
137     fCascadesTClArr->SetName("CascadesHF");
138     AddAODBranch("TClonesArray", &fCascadesTClArr);
139   }
140
141   if(fvHF->GetLikeSign()) {
142     fLikeSign2ProngTClArr = new TClonesArray("AliAODRecoDecayHF2Prong", 0);
143     fLikeSign2ProngTClArr->SetName("LikeSign2Prong");
144     AddAODBranch("TClonesArray", &fLikeSign2ProngTClArr);
145   }
146  
147   if(fvHF->GetLikeSign() && fvHF->Get3Prong()) {
148     fLikeSign3ProngTClArr = new TClonesArray("AliAODRecoDecayHF3Prong", 0);
149     fLikeSign3ProngTClArr->SetName("LikeSign3Prong");
150     AddAODBranch("TClonesArray", &fLikeSign3ProngTClArr);
151   }
152
153   return;
154 }
155
156 //________________________________________________________________________
157 void AliAnalysisTaskMEVertexingHF::UserExec(Option_t *) 
158 {
159   // Execute analysis for current event:
160   // first build the mixed event, compute the new primary vtx 
161   // then heavy flavor vertexing 
162
163   Int_t nev = fInputHandler->GetBufferSize();
164   fMixedEvent = new AliMixedEvent();
165   fMixedEvent->Reset();
166   TString primTitle;
167   TString primTitleFirst;
168   AliAODVertex *vtxCopy=0;
169
170   TObjArray *vertices=new TObjArray(nev);
171   for (Int_t iev = 0; iev < nev; iev++) {
172     AliAODEvent *evt = (AliAODEvent*)GetEvent(iev);
173     if(!evt) {delete vertices;return;}
174     AliAODVertex *evtVtx=(AliAODVertex*)evt->GetPrimaryVertex();
175     if(!evtVtx) {delete vertices;return;}
176     primTitle = evtVtx->GetTitle();
177     Int_t nContrib=evtVtx->GetNContributors();
178     if(!primTitle.Contains("VertexerTracks") || nContrib<=0) {
179       delete vertices;
180       return;
181     }
182
183     vtxCopy=new AliAODVertex(*evtVtx);
184     primTitleFirst=evtVtx->GetTitle();
185         
186
187     fMixedEvent->AddEvent(evt);
188
189     vertices->AddLast(vtxCopy);
190   }
191
192
193   fMixedEvent->Init();
194   Double_t vtxPos[3]={0.,0.,0.},vtxSigma[3]={0.,0.,0.};
195   Int_t nContributors[1]={0};
196   Double_t chi2=0;
197   Bool_t primaryOk=fMixedEvent->ComputeVtx(vertices,vtxPos,vtxSigma,nContributors);
198   if(!primaryOk) {
199     delete vertices;
200     delete vtxCopy;
201     vtxCopy=NULL;
202     return;
203   }
204   Int_t contribCopy=nContributors[0];
205   Double_t vtxCov[6]={vtxSigma[0]*vtxSigma[0],0,vtxSigma[1]*vtxSigma[1],0,0,vtxSigma[2]*vtxSigma[2]};
206   AliVVertex* newVertex=new AliESDVertex(vtxPos,vtxCov,chi2,contribCopy);
207   newVertex->SetTitle(primTitleFirst.Data());
208   fMixedEvent->SetPrimaryVertex(newVertex);
209
210   delete vertices;
211   delete vtxCopy;
212   vtxCopy=NULL;
213
214   fvHF->FindCandidates(fMixedEvent,
215                        fVerticesHFTClArr,
216                        fD0toKpiTClArr,
217                        fJPSItoEleTClArr,
218                        fCharm3ProngTClArr,
219                        fCharm4ProngTClArr,
220                        fDstarTClArr,
221                        fCascadesTClArr, 
222                        fLikeSign2ProngTClArr,
223                        fLikeSign3ProngTClArr);
224
225   delete newVertex;
226   return;
227 }      
228
229 //________________________________________________________________________
230 void AliAnalysisTaskMEVertexingHF::Terminate(Option_t *) 
231 {
232   // Terminate analysis
233 }