]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliInvmass.cxx
Have ProdProcess return const char*
[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 /*
17 $Log$
18 */
19
20 #include "AliInvmass.h"
21  
22 ClassImp(AliInvmass) // Class implementation to enable ROOT I/O
23  
24 AliInvmass::AliInvmass()
25 {
26 // Creation of an AliInvmass object and initialisation of parameters
27  fPi=acos(-1.);
28  fMode=2;
29  fBkg=0;
30  fNewtheta=1;
31  fNewphi=1;
32  fMinv=0;
33  fMbkg=0;
34 }
35 ////////////////////////////////////////////////////////////////////////////////
36 AliInvmass::~AliInvmass()
37 {
38 // Destructor to delete dynamically allocated memory
39  if (fMinv)
40  {
41   fMinv->Delete();
42   delete fMinv;
43   fMinv=0;
44  }
45
46  if (fMbkg)
47  {
48   fMbkg->Delete();
49   delete fMbkg;
50   fMbkg=0;
51  }
52 }
53 ////////////////////////////////////////////////////////////////////////////////
54 void AliInvmass::SetStorageMode(Int_t m)
55 {
56 // Set storage mode for the result arrays for inv. mass and comb. background
57  fMode=2;
58  if (m==1) fMode=1;
59 }
60 ////////////////////////////////////////////////////////////////////////////////
61 void AliInvmass::SetThetaSwitch(Int_t i)
62 {
63 // Enable/Disable (1/0) switching of theta angle in comb. bkg. reconstruction.
64 // Default : Switching of theta is enabled.
65  fNewtheta=1;
66  if (i==0) fNewtheta=0;
67 }
68 ////////////////////////////////////////////////////////////////////////////////
69 void AliInvmass::SetPhiSwitch(Int_t i)
70 {
71 // Enable/Disable (1/0) switching of phi angle in comb. bkg. reconstruction.
72 // Default : Switching of phi is enabled.
73  fNewphi=1;
74  if (i==0) fNewphi=0;
75 }
76 ////////////////////////////////////////////////////////////////////////////////
77 Int_t AliInvmass::GetStorageMode()
78 {
79 // Provide mode of storage for the result arrays for inv. mass and comb. background
80  return fMode;
81 }
82 ////////////////////////////////////////////////////////////////////////////////
83 Int_t AliInvmass::GetThetaSwitch()
84 {
85 // Provide the theta switching flag
86  return fNewtheta;
87 }
88 ////////////////////////////////////////////////////////////////////////////////
89 Int_t AliInvmass::GetPhiSwitch()
90 {
91 // Provide the phi switching flag
92  return fNewphi;
93 }
94 ////////////////////////////////////////////////////////////////////////////////
95 TObjArray* AliInvmass::Invmass(TObjArray* a1,TObjArray* a2)
96 {
97 // Perform two-particle invariant mass reconstruction
98  fBkg=0;
99  Combine(a1,a2);
100  return fMinv;
101 }
102 ////////////////////////////////////////////////////////////////////////////////
103 TObjArray* AliInvmass::CombBkg(TObjArray* a1,TObjArray* a2)
104 {
105 // Perform two-particle combinatorial background reconstruction
106  fBkg=1;
107  Combine(a1,a2);
108  if (fMode != 1)
109  {
110   return fMbkg;
111  }
112  else
113  {
114   return fMinv;
115  }
116 }
117 ////////////////////////////////////////////////////////////////////////////////
118 void AliInvmass::Combine(TObjArray* a1,TObjArray* a2)
119 {
120 // Perform two-particle invariant mass reconstruction
121  
122  if ((!fBkg || fMode==1) && fMinv)
123  {
124   fMinv->Delete();
125   delete fMinv;
126   fMinv=0;
127  }
128  
129  if (fBkg && (fMode !=1) && fMbkg)
130  {
131   fMbkg->Delete();
132   delete fMbkg;
133   fMbkg=0;
134  }
135
136  Int_t isame; // Indicates whether both lists are identical
137  isame=0;
138  if (a1==a2) isame=1;
139
140  // Index i must loop over the shortest of a1 and a2
141  TObjArray* listi=a1;
142  TObjArray* listj=a2;
143  Int_t ni=a1->GetEntries();
144  Int_t nj=a2->GetEntries();
145  if (nj < ni)
146  {
147   ni=a2->GetEntries();
148   nj=a1->GetEntries();
149   listi=a2;
150   listj=a1;
151  }
152
153  AliTrack* p1=0;  
154  AliTrack* p2=0;
155  AliTrack* px=0;
156  Ali4Vector ptot;
157  AliTrack* t=0;
158  Double_t v2[4],vx[4];
159  Float_t q1,q2;
160  
161  Int_t jmin; // Start index for list j
162  Int_t jx;   // Index for randomly picked particle for comb. bkg. reconstruction 
163  
164  for (Int_t i=0; i<ni; i++) // Select first a particle from list i
165  {
166   p1=(AliTrack*)listi->At(i);
167   p2=0;
168
169   if (!p1) continue;
170
171   jmin=0;
172   if (isame) jmin=i+1;
173   for (Int_t j=jmin; j<nj; j++) // Select also a particle from list j
174   {
175    p2=(AliTrack*)listj->At(j);
176    if (p1==p2) p2=0; // Don't combine particle with itself
177
178    if (!p2) continue;
179
180    p2->GetVector(v2,"sph");
181
182    // Take theta and phi from randomly chosen other list j particle for bkg. reconstr.
183    if (fBkg)
184    {
185     px=0; 
186     if ((!isame && nj>1) || (isame && nj>2))
187     {
188      jx=int(fRndm.Uniform(0,float(nj)));
189      px=(AliTrack*)listj->At(jx);
190      
191      while (!px || px==p2 || px==p1)
192      {
193       jx++;
194       if (jx >= nj) jx=0;
195       px=(AliTrack*)listj->At(jx);
196      }
197
198      px->GetVector(vx,"sph");
199      if (fNewtheta) v2[2]=vx[2]; // Replace the theta angle in the v2 vector
200      if (fNewphi) v2[3]=vx[3]; // Replace the phi angle in the v2 vector
201     }
202    }
203  
204    if ((!fBkg && p2) || (fBkg && px))
205    {
206     // Store the data of this two-particle combination
207     ptot.SetVector(v2,"sph");
208     ptot=(Ali4Vector)(ptot+(*p1));
209     q1=p1->GetCharge();
210     q2=p2->GetCharge();
211     t=new AliTrack;
212     t->Set4Momentum(ptot);
213     t->SetCharge(q1+q2);
214     if (!fBkg || fMode==1) 
215     {
216      if (!fMinv) fMinv=new TObjArray();
217      fMinv->Add(t);
218     }
219     else
220     {
221      if (!fMbkg) fMbkg=new TObjArray();
222      fMbkg->Add(t);
223     }
224    }
225
226   } // End of second particle loop
227
228  } // End of first particle loop
229 }
230 ////////////////////////////////////////////////////////////////////////////////