]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliInvmass.cxx
Protection against wrong name added
[u/mrichter/AliRoot.git] / RALICE / AliInvmass.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 // $Id$
17
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
102 //- Modified: NvE $Date$ UU-SAP Utrecht
103 ////////////////////////////////////////////////////////////////////////////////
104
105 #include "AliInvmass.h"
106  
107 ClassImp(AliInvmass) // Class implementation to enable ROOT I/O
108  
109 AliInvmass::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 ////////////////////////////////////////////////////////////////////////////////
121 AliInvmass::~AliInvmass()
122 {
123 // Destructor to delete dynamically allocated memory
124  if (fMinv)
125  {
126   delete fMinv;
127   fMinv=0;
128  }
129
130  if (fMbkg)
131  {
132   delete fMbkg;
133   fMbkg=0;
134  }
135 }
136 ////////////////////////////////////////////////////////////////////////////////
137 void AliInvmass::SetStorageMode(Int_t m)
138 {
139 // Set storage mode for the result arrays for inv. mass and comb. background
140  fMode=2;
141  if (m==1) fMode=1;
142 }
143 ////////////////////////////////////////////////////////////////////////////////
144 void AliInvmass::SetThetaSwitch(Int_t i)
145 {
146 // Enable/Disable (1/0) switching of theta angle in comb. bkg. reconstruction.
147 // Default : Switching of theta is enabled.
148  fNewtheta=1;
149  if (i==0) fNewtheta=0;
150 }
151 ////////////////////////////////////////////////////////////////////////////////
152 void AliInvmass::SetPhiSwitch(Int_t i)
153 {
154 // Enable/Disable (1/0) switching of phi angle in comb. bkg. reconstruction.
155 // Default : Switching of phi is enabled.
156  fNewphi=1;
157  if (i==0) fNewphi=0;
158 }
159 ////////////////////////////////////////////////////////////////////////////////
160 Int_t AliInvmass::GetStorageMode()
161 {
162 // Provide mode of storage for the result arrays for inv. mass and comb. background
163  return fMode;
164 }
165 ////////////////////////////////////////////////////////////////////////////////
166 Int_t AliInvmass::GetThetaSwitch()
167 {
168 // Provide the theta switching flag
169  return fNewtheta;
170 }
171 ////////////////////////////////////////////////////////////////////////////////
172 Int_t AliInvmass::GetPhiSwitch()
173 {
174 // Provide the phi switching flag
175  return fNewphi;
176 }
177 ////////////////////////////////////////////////////////////////////////////////
178 TObjArray* AliInvmass::Invmass(TObjArray* a1,TObjArray* a2)
179 {
180 // Perform two-particle invariant mass reconstruction
181  fBkg=0;
182  Combine(a1,a2);
183  return fMinv;
184 }
185 ////////////////////////////////////////////////////////////////////////////////
186 TObjArray* AliInvmass::CombBkg(TObjArray* a1,TObjArray* a2)
187 {
188 // Perform two-particle combinatorial background reconstruction
189  fBkg=1;
190  Combine(a1,a2);
191  if (fMode != 1)
192  {
193   return fMbkg;
194  }
195  else
196  {
197   return fMinv;
198  }
199 }
200 ////////////////////////////////////////////////////////////////////////////////
201 void AliInvmass::Combine(TObjArray* a1,TObjArray* a2)
202 {
203 // Perform two-particle invariant mass reconstruction
204  
205  if ((!fBkg || fMode==1) && fMinv)
206  {
207   delete fMinv;
208   fMinv=0;
209  }
210  
211  if (fBkg && (fMode !=1) && fMbkg)
212  {
213   delete fMbkg;
214   fMbkg=0;
215  }
216
217  Int_t isame; // Indicates whether both lists are identical
218  isame=0;
219  if (a1==a2) isame=1;
220
221  // Index i must loop over the shortest of a1 and a2
222  TObjArray* listi=a1;
223  TObjArray* listj=a2;
224  Int_t ni=a1->GetEntries();
225  Int_t nj=a2->GetEntries();
226  if (nj < ni)
227  {
228   ni=a2->GetEntries();
229   nj=a1->GetEntries();
230   listi=a2;
231   listj=a1;
232  }
233
234  AliTrack* p1=0;  
235  AliTrack* p2=0;
236  AliTrack* px=0;
237  Ali4Vector ptot;
238  AliTrack* t=0;
239  Double_t v2[4],vx[4];
240  Float_t q1,q2;
241  
242  Int_t jmin; // Start index for list j
243  Int_t jx;   // Index for randomly picked particle for comb. bkg. reconstruction 
244  
245  for (Int_t i=0; i<ni; i++) // Select first a particle from list i
246  {
247   p1=(AliTrack*)listi->At(i);
248   p2=0;
249
250   if (!p1) continue;
251
252   jmin=0;
253   if (isame) jmin=i+1;
254   for (Int_t j=jmin; j<nj; j++) // Select also a particle from list j
255   {
256    p2=(AliTrack*)listj->At(j);
257    if (p1==p2) p2=0; // Don't combine particle with itself
258
259    if (!p2) continue;
260
261    p2->GetVector(v2,"sph");
262
263    // Take theta and phi from randomly chosen other list j particle for bkg. reconstr.
264    if (fBkg)
265    {
266     px=0; 
267     if ((!isame && nj>1) || (isame && nj>2))
268     {
269      jx=int(fRndm.Uniform(0,float(nj)));
270      px=(AliTrack*)listj->At(jx);
271      
272      while (!px || px==p2 || px==p1)
273      {
274       jx++;
275       if (jx >= nj) jx=0;
276       px=(AliTrack*)listj->At(jx);
277      }
278
279      px->GetVector(vx,"sph");
280      if (fNewtheta) v2[2]=vx[2]; // Replace the theta angle in the v2 vector
281      if (fNewphi) v2[3]=vx[3]; // Replace the phi angle in the v2 vector
282     }
283    }
284  
285    if ((!fBkg && p2) || (fBkg && px))
286    {
287     // Store the data of this two-particle combination
288     ptot.SetVector(v2,"sph");
289     ptot=(Ali4Vector)(ptot+(*p1));
290     q1=p1->GetCharge();
291     q2=p2->GetCharge();
292     t=new AliTrack;
293     t->Set4Momentum(ptot);
294     t->SetCharge(q1+q2);
295     if (!fBkg || fMode==1) 
296     {
297      if (!fMinv)
298      {
299       fMinv=new TObjArray();
300       fMinv->SetOwner();
301      }
302      fMinv->Add(t);
303     }
304     else
305     {
306      if (!fMbkg)
307      {
308       fMbkg=new TObjArray();
309       fMbkg->SetOwner();
310      }
311      fMbkg->Add(t);
312     }
313    }
314
315   } // End of second particle loop
316
317  } // End of first particle loop
318 }
319 ////////////////////////////////////////////////////////////////////////////////