]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnAnalysisSimple.cxx
Introducing a new AOD class: AliAODcascade (A.Maire)
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnAnalysisSimple.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //-------------------------------------------------------------------------
17 //                     Class AliRsnAnalysisSimple
18 //-------------------------------------------------------------------------
19 // This class is a manager to process one or more pair analysis with a
20 // given event sample, specified as a TTree passed to this object.
21 // Each single pair analysis must be defined with its specifications, like
22 // histogram binning, cuts and everything else.
23 // This class utiliy consists in defining a unique event sample which is
24 // used for all pairs analysis, and all kinds of event mixing (if required).
25 // All histograms computed in a single execution, are then saved into a file.
26 // This object contains a two TObjArray's:
27 //  - one to contain single event pair analysis objects (signal, like-sign)
28 //  - one to contain all event-mixing pair analysis objects
29 // When a new pair is added, it must be then specified in what container it
30 // must be placed, in order to avoid meaningless results.
31 //
32 // author: A. Pulvirenti
33 // email : alberto.pulvirenti@ct.infn.it
34 //-------------------------------------------------------------------------
35
36 #include <TH1.h>
37 #include <TTree.h>
38 #include <TFile.h>
39 #include <TArray.h>
40 #include <TClonesArray.h>
41
42 #include "AliLog.h"
43 #include "AliRsnPID.h"
44 #include "AliRsnDaughter.h"
45 #include "AliRsnEvent.h"
46 #include "AliRsnEventBuffer.h"
47 #include "AliRsnPairSimple.h"
48 #include "AliRsnAnalyzerSimple.h"
49 #include "AliRsnAnalysisSimple.h"
50
51 ClassImp(AliRsnAnalysisSimple)
52
53 //_____________________________________________________________________________
54 AliRsnAnalysisSimple::AliRsnAnalysisSimple(AliRsnAnalyzerSimple *ana, AliRsnPID *pid) :
55   TObject(),
56   fInitialized(kFALSE),
57   fStep(1000),
58   fTree(0x0),
59   fPID(pid),
60   fAnalyzer(ana)
61 {
62 //
63 // Constructor
64 // Initializes all pointers and collections to NULL.
65 //
66
67     strcpy(fFileName, "default.root");
68 }
69
70
71 //_____________________________________________________________________________
72 void AliRsnAnalysisSimple::Clear(Option_t* /*option*/)
73 {
74 //
75 // Clear heap
76 //
77 }
78
79 //_____________________________________________________________________________
80 void AliRsnAnalysisSimple::SetEventsTree(TTree *tree)
81 {
82 //
83 // Set the tree containing the events to be processed.
84 // Counts also the number of events and stores it in a private datamember.
85 // This can avoid the time-wasting entries count in a long TChain.
86 //
87     fTree = tree;
88     AliInfo(Form("Total number of events to be processed: %d", fTree->GetEntries()));
89 }
90
91 //_____________________________________________________________________________
92 Bool_t AliRsnAnalysisSimple::Initialize()
93 {
94 //
95 // Various initialization processes
96 //
97     // check process objects
98     if (!fPID) {
99         AliError("PID not initialized");
100         return kFALSE;
101     }
102     if (!fAnalyzer) {
103         AliError("Analyzer not initialized");
104         return kFALSE;
105     }
106
107     // set reference to PID method to all pairs
108     AliRsnPairSimple *pair = 0;
109     TObjArrayIter pairIterator(fAnalyzer->GetPairs());
110     while ( (pair = (AliRsnPairSimple*)pairIterator.Next()) ) {
111         pair->SetPIDMethod(fPID->GetMethod());
112     }
113     if (fAnalyzer->GetMixPairs() && !fAnalyzer->GetMixPairs()->IsEmpty()) {
114         TObjArrayIter mixIterator(fAnalyzer->GetMixPairs());
115         while ( (pair = (AliRsnPairSimple*)mixIterator.Next()) ) {
116             pair->SetPIDMethod(fPID->GetMethod());
117         }
118     }
119
120     // initialize analyzer
121     fAnalyzer->Init();
122
123     // at the end, update flag for initialization
124     fInitialized = kTRUE;
125
126     return kTRUE;
127 }
128
129 //_____________________________________________________________________________
130 Stat_t AliRsnAnalysisSimple::Process()
131 {
132 //
133 // Computes all invariant mass distributions defined in the AliRsnPair objects collected.
134 // Depending on the kind of background evaluation method, computes also this one.
135 //
136
137     // check initialization
138     if (!fInitialized) {
139         AliError("Analysis not initialized. Use method 'Initialize()'");
140         return 0.0;
141     }
142
143     // set cursor object
144     AliRsnEvent *event = 0x0;
145     fTree->SetBranchAddress("rsnEvents", &event);
146
147     // output counter
148     Stat_t nPairs = 0.0;
149
150         // loop on events
151         Int_t i, nEvents = (Int_t)fTree->GetEntries();
152     for (i = 0; i < nEvents; i++) {
153         // message
154         if ((i % fStep) == 0) AliInfo(Form("Processing event %d", i));
155         // get entry
156         fTree->GetEntry(i);
157         if (!event) continue;
158         fPID->Identify(event);
159         nPairs += fAnalyzer->Process(event);
160     }
161
162     return nPairs;
163 }
164
165 //_____________________________________________________________________________
166 void AliRsnAnalysisSimple::SaveOutput() const
167 {
168 //
169 // Writes histograms in current directory
170 //
171     TFile *file = TFile::Open(fFileName, "RECREATE");
172     AliRsnPairSimple *pair = 0;
173     TH1D *hist = 0;
174     TObjArrayIter pairIterator(fAnalyzer->GetPairs());
175     while ( (pair = (AliRsnPairSimple*)pairIterator.Next()) ) {
176         hist = pair->GetHistogram();
177         if (hist) hist->Write();
178         hist = pair->GetHistogramMC();
179         if (hist) hist->Write();
180     }
181     if (fAnalyzer->GetMixPairs() && !fAnalyzer->GetMixPairs()->IsEmpty()) {
182         TObjArrayIter mixIterator(fAnalyzer->GetMixPairs());
183         while ( (pair = (AliRsnPairSimple*)mixIterator.Next()) ) {
184             hist = pair->GetHistogram();
185             if (hist) hist->Write();
186             hist = pair->GetHistogramMC();
187             if (hist) hist->Write();
188         }
189         }
190         file->Close();
191 }