1 /**************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /////////////////////////////////////////////////////////////
18 // AliAnalysisTaskSE for the reconstruction of heavy flavor
19 // decays, using the class AliAnalysisVertexingHF.
21 // Author: A.Dainese, andrea.dainese@lnl.infn.it
22 /////////////////////////////////////////////////////////////
26 #include <TClonesArray.h>
28 #include "AliVEvent.h"
29 #include "AliAODEvent.h"
30 #include "AliESDEvent.h"
31 #include "AliAnalysisVertexingHF.h"
32 #include "AliAnalysisTaskSE.h"
33 #include "AliAnalysisManager.h"
34 #include "AliAnalysisTaskSEVertexingHF.h"
36 ClassImp(AliAnalysisTaskSEVertexingHF)
39 //________________________________________________________________________
40 AliAnalysisTaskSEVertexingHF::AliAnalysisTaskSEVertexingHF():
46 fCharm3ProngTClArr(0),
47 fCharm4ProngTClArr(0),
49 fLikeSign2ProngTClArr(0),
50 fLikeSign3ProngTClArr(0)
52 // Default constructor
55 //________________________________________________________________________
56 AliAnalysisTaskSEVertexingHF::AliAnalysisTaskSEVertexingHF(const char *name):
57 AliAnalysisTaskSE(name),
62 fCharm3ProngTClArr(0),
63 fCharm4ProngTClArr(0),
65 fLikeSign2ProngTClArr(0),
66 fLikeSign3ProngTClArr(0)
68 // Default constructor
71 //________________________________________________________________________
72 AliAnalysisTaskSEVertexingHF::~AliAnalysisTaskSEVertexingHF()
77 //________________________________________________________________________
78 void AliAnalysisTaskSEVertexingHF::Init()
81 // Instanciates vHF and loads its parameters
83 if(fDebug > 1) printf("AnalysisTaskSEVertexingHF::Init() \n");
85 if(gROOT->LoadMacro("ConfigVertexingHF.C")) {
86 printf("AnalysisTaskSEVertexingHF::Init() \n Using $ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C\n");
87 gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C");
90 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
92 AliAnalysisManager::GetAnalysisManager()->RegisterExtraFile("AliAOD.VertexingHF.root");
97 //________________________________________________________________________
98 void AliAnalysisTaskSEVertexingHF::UserCreateOutputObjects()
100 // Create the output container
102 if(fDebug > 1) printf("AnalysisTaskSEVertexingHF::UserCreateOutPutData() \n");
103 // Support both the case when the AOD + deltaAOD are produced in an ESD
104 // analysis or if the deltaAOD is produced on an analysis on AOD's. (A.G. 27/04/09)
106 Fatal("UserCreateOutputObjects", "This task needs an AOD handler");
109 TString filename = "AliAOD.VertexingHF.root";
110 if (!IsStandardAOD()) filename = "";
112 printf("AnalysisTaskSEVertexingHF::UserCreateOutPutData() \n ERROR! no fvHF!\n");
116 fVerticesHFTClArr = new TClonesArray("AliAODVertex", 0);
117 fVerticesHFTClArr->SetName("VerticesHF");
118 AddAODBranch("TClonesArray", &fVerticesHFTClArr, filename);
120 if(fVHF->GetD0toKpi()) {
121 fD0toKpiTClArr = new TClonesArray("AliAODRecoDecayHF2Prong", 0);
122 fD0toKpiTClArr->SetName("D0toKpi");
123 AddAODBranch("TClonesArray", &fD0toKpiTClArr, filename);
126 if(fVHF->GetJPSItoEle()) {
127 fJPSItoEleTClArr = new TClonesArray("AliAODRecoDecayHF2Prong", 0);
128 fJPSItoEleTClArr->SetName("JPSItoEle");
129 AddAODBranch("TClonesArray", &fJPSItoEleTClArr, filename);
132 if(fVHF->Get3Prong()) {
133 fCharm3ProngTClArr = new TClonesArray("AliAODRecoDecayHF3Prong", 0);
134 fCharm3ProngTClArr->SetName("Charm3Prong");
135 AddAODBranch("TClonesArray", &fCharm3ProngTClArr, filename);
138 if(fVHF->Get4Prong()) {
139 fCharm4ProngTClArr = new TClonesArray("AliAODRecoDecayHF4Prong", 0);
140 fCharm4ProngTClArr->SetName("Charm4Prong");
141 AddAODBranch("TClonesArray", &fCharm4ProngTClArr, filename);
144 if(fVHF->GetDstar()) {
145 fDstarTClArr = new TClonesArray("AliAODRecoCascadeHF", 0);
146 fDstarTClArr->SetName("Dstar");
147 AddAODBranch("TClonesArray", &fDstarTClArr, filename);
150 if(fVHF->GetLikeSign()) {
151 fLikeSign2ProngTClArr = new TClonesArray("AliAODRecoDecayHF2Prong", 0);
152 fLikeSign2ProngTClArr->SetName("LikeSign2Prong");
153 AddAODBranch("TClonesArray", &fLikeSign2ProngTClArr, filename);
156 if(fVHF->GetLikeSign() && fVHF->Get3Prong()) {
157 fLikeSign3ProngTClArr = new TClonesArray("AliAODRecoDecayHF3Prong", 0);
158 fLikeSign3ProngTClArr->SetName("LikeSign3Prong");
159 AddAODBranch("TClonesArray", &fLikeSign3ProngTClArr, filename);
165 //________________________________________________________________________
166 void AliAnalysisTaskSEVertexingHF::UserExec(Option_t */*option*/)
168 // Execute analysis for current event:
169 // heavy flavor vertexing
171 AliVEvent *event = dynamic_cast<AliVEvent*> (InputEvent());
172 // In case there is an AOD handler writing a standard AOD, use the AOD
173 // event in memory rather than the input (ESD) event. (A.G. 27/04/09)
174 if (AODEvent() && IsStandardAOD()) event = dynamic_cast<AliVEvent*> (AODEvent());
176 // heavy flavor vertexing
177 fVHF->FindCandidates(event,
184 fLikeSign2ProngTClArr,
185 fLikeSign3ProngTClArr);
190 //________________________________________________________________________
191 void AliAnalysisTaskSEVertexingHF::Terminate(Option_t */*option*/)
193 // Terminate analysis
195 if(fDebug > 1) printf("AnalysisTaskSEVertexingHF: Terminate() \n");