]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPhysicsSelectionTask.cxx
fix for generating the outputs_valid file in MC generator mode
[u/mrichter/AliRoot.git] / ANALYSIS / AliPhysicsSelectionTask.cxx
1 /* $Id$ */
2
3 #include "AliPhysicsSelectionTask.h"
4
5 #include <TFile.h>
6 #include <TH1F.h>
7 #include <TH2F.h>
8
9 #include <AliLog.h>
10 #include <AliESDEvent.h>
11 #include <AliHeader.h>
12
13 #include "AliPhysicsSelection.h"
14 #include "AliAnalysisManager.h"
15 #include "AliInputEventHandler.h"
16
17 //#include "AliBackgroundSelection.h"
18
19 ClassImp(AliPhysicsSelectionTask)
20
21 AliPhysicsSelectionTask::AliPhysicsSelectionTask() :
22   AliAnalysisTaskSE("AliPhysicsSelectionTask"),
23   fOutput(0),
24   fOption(""),
25   fUseSpecialOutput(kFALSE),
26   fPhysicsSelection(0)
27 {
28   //
29   // Default event handler
30   //
31 }
32
33 AliPhysicsSelectionTask::AliPhysicsSelectionTask(const char* opt) :
34   AliAnalysisTaskSE("AliPhysicsSelectionTask"),
35   fOutput(0),
36   fOption(opt),
37   fUseSpecialOutput(kFALSE),
38   fPhysicsSelection(new AliPhysicsSelection())
39 {
40   //
41   // Constructor. Initialization of pointers
42   //
43   AliInputEventHandler* handler = dynamic_cast<AliInputEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
44   if (handler) {
45     handler->SetEventSelection(fPhysicsSelection);
46     AliInfo("Physics Event Selection enabled.");
47   } else {
48     AliError("No input event handler connected to analysis manager. No Physics Event Selection.");
49   }
50   //
51   TString opts = opt;
52   opts.ToLower();
53   if (opts.Contains("specialoutput")) fUseSpecialOutput = kTRUE;
54
55   // Define input and output slots here
56   DefineOutput(1, TList::Class());
57   fBranchNames = "ESD:AliESDRun.,AliESDHeader.,AliMultiplicity.,AliESDVZERO.,"
58                  "AliESDZDC.,SPDVertex.,PrimaryVertex.,TPCVertex.,Tracks,SPDPileupVertices";
59   
60   AliLog::SetClassDebugLevel("AliPhysicsSelectionTask", AliLog::kWarning);
61 }
62
63 AliPhysicsSelectionTask::~AliPhysicsSelectionTask()
64 {
65   //
66   // Destructor
67   //
68
69   // histograms are in the output list and deleted when the output
70   // list is deleted by the TSelector dtor
71
72   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
73     delete fOutput;
74     fOutput = 0;
75   }
76 }
77
78 void AliPhysicsSelectionTask::UserCreateOutputObjects()
79 {
80   // create result objects and add to output list
81
82   Printf("AliPhysicsSelectionTask::CreateOutputObjects");
83
84   if (fUseSpecialOutput) OpenFile(1);
85
86   fOutput = new TList;
87   fOutput->SetOwner();
88   
89   if (!fPhysicsSelection)
90     fPhysicsSelection = new AliPhysicsSelection;
91   
92   fOutput->Add(fPhysicsSelection);
93   // All tasks must post data once for all outputs (AG)
94   PostData(1, fOutput);
95 }
96
97 void AliPhysicsSelectionTask::UserExec(Option_t*)
98 {
99   // process the event
100
101   // AliPhysicsSelection::IsCollisionCandidate is called from the event handler
102   // post the data here anyway!
103   PostData(1, fOutput);
104 }
105
106 void AliPhysicsSelectionTask::FinishTaskOutput()
107 {
108 // This gets called at the end of the processing on the worker. It allows dumping
109 // statistics printed by the physics selection object to the statistics message
110 // handled by the analysis manager.
111    if (fPhysicsSelection) fPhysicsSelection->Print("STAT");
112 }
113
114 void AliPhysicsSelectionTask::Terminate(Option_t *)
115 {
116   // The Terminate() function is the last function to be called during
117   // a query. It always runs on the client, it can be used to present
118   // the results graphically or save the results to file.
119
120   fOutput = dynamic_cast<TList*> (GetOutputData(1));
121   if (!fOutput)
122     Printf("ERROR: fOutput not available");
123     
124   if (fOutput)
125   {
126     fPhysicsSelection = dynamic_cast<AliPhysicsSelection*> (fOutput->FindObject("AliPhysicsSelection"));
127   }
128
129   TFile* fout = new TFile("event_stat.root", "RECREATE");
130
131   if (fPhysicsSelection)
132   {
133     fPhysicsSelection->Print();
134     fPhysicsSelection->SaveHistograms();
135   }
136     
137   fout->Write();
138   fout->Close();
139   
140   Printf("Writing result to event_stat.root");
141 }