1 /**************************************************************************
2 * Copyright(c) 1998-1999, 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 //=========================================================================
17 // Class AliRsnSimpleAnalyzer
19 // Implementation of the event processing which returns histograms of
20 // invariant mass for resonances and backgrounds evaluation.
22 // author: A. Pulvirenti (email: alberto.pulvirenti@ct.infn.it)
23 //=========================================================================
28 #include <TClonesArray.h>
29 #include <TDirectory.h>
33 #include "AliRsnDaughter.h"
34 #include "AliRsnEvent.h"
35 #include "AliRsnEventBuffer.h"
36 #include "AliRsnSimpleFunction.h"
38 #include "AliRsnSimpleAnalyzer.h"
40 ClassImp(AliRsnSimpleAnalyzer)
42 //_____________________________________________________________________________
43 AliRsnSimpleAnalyzer::AliRsnSimpleAnalyzer(Int_t bufferSize) :
44 TNamed("RsnSimpleAnalyzer", ""),
45 fBufferSize(bufferSize),
55 // Initializes all pointers and collections to NULL.
56 // Adds this object to the global Directory.
58 gDirectory->Append(this, kTRUE);
62 //_____________________________________________________________________________
63 void AliRsnSimpleAnalyzer::Clear(Option_t *option)
69 fSingle->Clear(option);
75 //_____________________________________________________________________________
76 void AliRsnSimpleAnalyzer::Add(AliRsnSimpleFunction *fcn)
79 // Add a pair of particle types to be analyzed.
80 // Second argument tells if the added object is for event mixing.
83 Bool_t mixing = fcn->MixFlag();
84 TObjArray* &target = mixing ? fMix : fSingle;
85 if (!target) target = new TObjArray(0);
89 //_____________________________________________________________________________
90 void AliRsnSimpleAnalyzer::Init()
93 // Initialize what needs to.
97 fBuffer = new AliRsnEventBuffer(fBufferSize, kFALSE);
100 AliRsnSimpleFunction *fcn = 0;
102 TObjArrayIter iter(fSingle);
103 while ( (fcn = (AliRsnSimpleFunction*)iter.Next()) ) {
108 TObjArrayIter iter(fMix);
109 while ( (fcn = (AliRsnSimpleFunction*)iter.Next()) ) {
115 //_____________________________________________________________________________
116 void AliRsnSimpleAnalyzer::Process(AliRsnEvent *event)
119 // Computes all invariant mass distributions defined in the AliRsnPair objects collected.
120 // Depending on the kind of background evaluation method, computes also this one.
123 // loop over the collection of pair defs
124 ProcessEvents(fSingle, event, 0x0);
127 // variables for mixing
128 Int_t mult1, mult2, i, j, count = 0;
130 // add this event to buffer
131 fBuffer->AddEvent(event);
133 vz1 = event->GetPrimaryVertexZ();
134 mult1 = event->GetMultiplicity();
135 if (!fBuffer->NEmptySlots()) {
136 i = fBuffer->GetEventsBufferIndex();
139 if (count >= fNMix) break;
140 if (j == fBufferSize) j = 0;
142 AliRsnEvent *evmix = fBuffer->GetEvent(j);
143 if (!evmix) continue;
144 vz2 = evmix->GetPrimaryVertexZ();
145 mult2 = evmix->GetMultiplicity();
146 if (TMath::Abs(vz1 - vz2) <= fMixVzCut && TMath::Abs(mult1- mult2) <= fMixMultCut) {
147 // loop over the collection of pair defs
148 ProcessEvents(fMix, event, evmix);
149 ProcessEvents(fMix, evmix, event);
157 //_____________________________________________________________________________
158 void AliRsnSimpleAnalyzer::ProcessEvents
159 (TObjArray *array, AliRsnEvent *event1, AliRsnEvent *event2)
162 // Takes all AliRsnSimpleFunction objects in the passed array
163 // and processes the given pair of AliRsnEvents with all of them
164 // If the array is NULL nothing is done
169 AliRsnSimpleFunction *fcn = 0;
170 TObjArrayIter iter(array);
171 while ( (fcn = (AliRsnSimpleFunction*)iter.Next()) ) {
172 if (event2) fcn->ProcessTwo(event1, event2);
173 else fcn->ProcessOne(event1);