]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RALICE/AliInvmass.cxx
27-may-2001 NvE New class AliEvent introduced and RALICELinkDef.h & co. updated accor...
[u/mrichter/AliRoot.git] / RALICE / AliInvmass.cxx
CommitLineData
4c039060 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
f531a546 16// $Id$
4c039060 17
959fbac5 18////////////////////////////////////////////////////////////////////////////////
19// Class AliInvmass
20// Construction of invariant mass and combinatorial background.
21//
22// Example :
23// ---------
24//
25// TObjArray* photons=new TObjArray(); // Array with photon tracks for pi0 rec.
26//
27// // Code to create some photon tracks from pi0 decays
28// Int_t ntracks=200;
29// for (Int_t i=0; i<ntracks; i++)
30// {
31// photons->Add(new Alitrack);
32// ...
33// ...
34// ...
35// }
36//
37// // Perform the invariant mass and comb. bkg. reconstruction
38//
39// TObjArray* allm=q.Invmass(photons,photons); // All reconstructed invariant masses
40//
41// TH1F* hall=new TH1F("hall","hall",200,0,2); // Histo with M_inv of all combinations
42//
43// Int_t nall=0;
44// if (allm) nall=allm->GetEntries();
45//
46// AliTrack* t;
47// Float_t minv;
48// for (Int_t j=0; j<nall; j++)
49// {
50// t=(AliTrack*)allm->At(j);
51// if (t)
52// {
53// minv=t->GetMass();
54// hall->Fill(minv);
55// }
56// }
57//
58// TObjArray* bkgm=q.CombBkg(photons,photons); // Reconstructed comb. background
59//
60// TH1F* hbkg=new TH1F("hbkg","hbkg",200,0,2); // Histo with M_inv. of comb. background
61//
62// Int_t nbkg=0;
63// if (bkgm) nbkg=bkgm->GetEntries();
64//
65// for (Int_t j=0; j<nbkg; j++)
66// {
67// t=(AliTrack*)bkgm->At(j);
68// if (t)
69// {
70// minv=t->GetMass();
71// hbkg->Fill(minv);
72// }
73// }
74//
75// TH1F* hsig=new TH1F("sig","sig",200,0,2); // Histo with the bkg. subtracted signal
76// hsig->Sumw2();
77// hsig->Add(hall,hbkg,1,-1);
78//
79//
80// Note : By default the storage of the reconstructed information is performed
81// in separate TObjArrays for the signal and comb. background resp.
82// In order to limit the memory usage, AliInvmass::SetStorageMode(1) may be
83// used to activate only a single TObjArray to store the reconstructed information.
84// Consequently, the following statements
85//
86// TObjArray* allm=q.Invmass(photons,photons);
87// TObjArray* bkgm=q.CombBkg(photons,photons);
88//
89// will result in the fact that after he invokation of CombBkg
90// the information of "allm" is lost due to the fact that the storage is
91// is re-used for "bkgm" in case the "single storage" option has been selected.
92// Usage of the, in that case invalid, pointer "allm" may cause your
93// program to crash.
94//
95// * Thus : In case of single storage usage, all invokations of the returned
96// array pointer have to be completed before invoking any memberfunction
97// of the same AliInvmass object again.
98//
99//
100//
101//--- Author: Nick van Eijndhoven 12-apr-1999 UU-SAP Utrecht
f531a546 102//- Modified: NvE $Date$ UU-SAP Utrecht
959fbac5 103////////////////////////////////////////////////////////////////////////////////
104
d88f97cc 105#include "AliInvmass.h"
106
107ClassImp(AliInvmass) // Class implementation to enable ROOT I/O
108
109AliInvmass::AliInvmass()
110{
111// Creation of an AliInvmass object and initialisation of parameters
112 fPi=acos(-1.);
113 fMode=2;
114 fBkg=0;
115 fNewtheta=1;
116 fNewphi=1;
117 fMinv=0;
118 fMbkg=0;
119}
120////////////////////////////////////////////////////////////////////////////////
121AliInvmass::~AliInvmass()
122{
123// Destructor to delete dynamically allocated memory
124 if (fMinv)
125 {
126 fMinv->Delete();
127 delete fMinv;
128 fMinv=0;
129 }
130
131 if (fMbkg)
132 {
133 fMbkg->Delete();
134 delete fMbkg;
135 fMbkg=0;
136 }
137}
138////////////////////////////////////////////////////////////////////////////////
139void AliInvmass::SetStorageMode(Int_t m)
140{
141// Set storage mode for the result arrays for inv. mass and comb. background
142 fMode=2;
143 if (m==1) fMode=1;
144}
145////////////////////////////////////////////////////////////////////////////////
146void AliInvmass::SetThetaSwitch(Int_t i)
147{
148// Enable/Disable (1/0) switching of theta angle in comb. bkg. reconstruction.
149// Default : Switching of theta is enabled.
150 fNewtheta=1;
151 if (i==0) fNewtheta=0;
152}
153////////////////////////////////////////////////////////////////////////////////
154void AliInvmass::SetPhiSwitch(Int_t i)
155{
156// Enable/Disable (1/0) switching of phi angle in comb. bkg. reconstruction.
157// Default : Switching of phi is enabled.
158 fNewphi=1;
159 if (i==0) fNewphi=0;
160}
161////////////////////////////////////////////////////////////////////////////////
162Int_t AliInvmass::GetStorageMode()
163{
164// Provide mode of storage for the result arrays for inv. mass and comb. background
165 return fMode;
166}
167////////////////////////////////////////////////////////////////////////////////
168Int_t AliInvmass::GetThetaSwitch()
169{
170// Provide the theta switching flag
171 return fNewtheta;
172}
173////////////////////////////////////////////////////////////////////////////////
174Int_t AliInvmass::GetPhiSwitch()
175{
176// Provide the phi switching flag
177 return fNewphi;
178}
179////////////////////////////////////////////////////////////////////////////////
180TObjArray* AliInvmass::Invmass(TObjArray* a1,TObjArray* a2)
181{
182// Perform two-particle invariant mass reconstruction
183 fBkg=0;
184 Combine(a1,a2);
185 return fMinv;
186}
187////////////////////////////////////////////////////////////////////////////////
188TObjArray* AliInvmass::CombBkg(TObjArray* a1,TObjArray* a2)
189{
190// Perform two-particle combinatorial background reconstruction
191 fBkg=1;
192 Combine(a1,a2);
193 if (fMode != 1)
194 {
195 return fMbkg;
196 }
197 else
198 {
199 return fMinv;
200 }
201}
202////////////////////////////////////////////////////////////////////////////////
203void AliInvmass::Combine(TObjArray* a1,TObjArray* a2)
204{
205// Perform two-particle invariant mass reconstruction
206
207 if ((!fBkg || fMode==1) && fMinv)
208 {
209 fMinv->Delete();
210 delete fMinv;
211 fMinv=0;
212 }
213
214 if (fBkg && (fMode !=1) && fMbkg)
215 {
216 fMbkg->Delete();
217 delete fMbkg;
218 fMbkg=0;
219 }
220
221 Int_t isame; // Indicates whether both lists are identical
222 isame=0;
223 if (a1==a2) isame=1;
224
225 // Index i must loop over the shortest of a1 and a2
226 TObjArray* listi=a1;
227 TObjArray* listj=a2;
228 Int_t ni=a1->GetEntries();
229 Int_t nj=a2->GetEntries();
230 if (nj < ni)
231 {
232 ni=a2->GetEntries();
233 nj=a1->GetEntries();
234 listi=a2;
235 listj=a1;
236 }
237
238 AliTrack* p1=0;
239 AliTrack* p2=0;
240 AliTrack* px=0;
241 Ali4Vector ptot;
242 AliTrack* t=0;
243 Double_t v2[4],vx[4];
244 Float_t q1,q2;
245
246 Int_t jmin; // Start index for list j
247 Int_t jx; // Index for randomly picked particle for comb. bkg. reconstruction
248
249 for (Int_t i=0; i<ni; i++) // Select first a particle from list i
250 {
251 p1=(AliTrack*)listi->At(i);
252 p2=0;
253
254 if (!p1) continue;
255
256 jmin=0;
257 if (isame) jmin=i+1;
258 for (Int_t j=jmin; j<nj; j++) // Select also a particle from list j
259 {
260 p2=(AliTrack*)listj->At(j);
261 if (p1==p2) p2=0; // Don't combine particle with itself
262
263 if (!p2) continue;
264
265 p2->GetVector(v2,"sph");
266
267 // Take theta and phi from randomly chosen other list j particle for bkg. reconstr.
268 if (fBkg)
269 {
270 px=0;
271 if ((!isame && nj>1) || (isame && nj>2))
272 {
273 jx=int(fRndm.Uniform(0,float(nj)));
274 px=(AliTrack*)listj->At(jx);
275
276 while (!px || px==p2 || px==p1)
277 {
278 jx++;
279 if (jx >= nj) jx=0;
280 px=(AliTrack*)listj->At(jx);
281 }
282
283 px->GetVector(vx,"sph");
284 if (fNewtheta) v2[2]=vx[2]; // Replace the theta angle in the v2 vector
285 if (fNewphi) v2[3]=vx[3]; // Replace the phi angle in the v2 vector
286 }
287 }
288
289 if ((!fBkg && p2) || (fBkg && px))
290 {
291 // Store the data of this two-particle combination
292 ptot.SetVector(v2,"sph");
293 ptot=(Ali4Vector)(ptot+(*p1));
294 q1=p1->GetCharge();
295 q2=p2->GetCharge();
296 t=new AliTrack;
297 t->Set4Momentum(ptot);
298 t->SetCharge(q1+q2);
299 if (!fBkg || fMode==1)
300 {
301 if (!fMinv) fMinv=new TObjArray();
302 fMinv->Add(t);
303 }
304 else
305 {
306 if (!fMbkg) fMbkg=new TObjArray();
307 fMbkg->Add(t);
308 }
309 }
310
311 } // End of second particle loop
312
313 } // End of first particle loop
314}
315////////////////////////////////////////////////////////////////////////////////