]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAOD.cxx
Process codes corrected/completed.
[u/mrichter/AliRoot.git] / ANALYSIS / AliAOD.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 //
20 // base class for AOD containers
21 //
22 /////////////////////////////////////////////////////////////
23
24 #include <TROOT.h>
25 #include <TParticle.h>
26 #include <TClass.h>
27 #include <TString.h>
28
29 #include "AliAOD.h"
30 #include "AliAODParticle.h"
31 #include "AliTrackPoints.h"
32
33 ClassImp(AliAOD)
34
35 AliAOD::AliAOD():
36  fParticles(0x0),
37  fIsRandomized(kFALSE),
38  fPrimaryVertexX(0.0),
39  fPrimaryVertexY(0.0),
40  fPrimaryVertexZ(0.0),
41  fParticleClass(0x0)
42 {
43  //ctor
44 // Info("AliAOD()","Entered");
45 // SetOwner(kTRUE);
46 // Info("AliAOD()","Exited");
47 }
48 /**************************************************************************/
49
50 AliAOD::AliAOD(const AliAOD& in):
51  TObject(in),
52  fParticles((TClonesArray*)in.fParticles->Clone()),
53   fIsRandomized(in.fIsRandomized),
54   fPrimaryVertexX(fPrimaryVertexX),
55   fPrimaryVertexY(in.fPrimaryVertexY),
56   fPrimaryVertexZ(in.fPrimaryVertexZ),
57   fParticleClass(in.fParticleClass)
58 {
59 //copy constructor
60 }
61 /**************************************************************************/
62
63 AliAOD& AliAOD::operator=(const AliAOD& in)
64 {
65 //assigment operator  
66
67   if (this == &in ) return *this;
68   
69   delete fParticles;
70   fParticles = (TClonesArray*)in.fParticles->Clone();
71   fIsRandomized = in.fIsRandomized ;
72   fPrimaryVertexX = in.fPrimaryVertexX ;
73   fPrimaryVertexY = in.fPrimaryVertexY ;
74   fPrimaryVertexZ = in.fPrimaryVertexZ ;
75   fParticleClass = in.fParticleClass ; //althought it is pointer, this points to object in class list of gROOT
76   return *this;
77 }
78 /**************************************************************************/
79
80 AliAOD::~AliAOD()
81 {
82   //Destructor
83   //fParticleClass does not belong to AliAOD -> Do not delete it
84   delete fParticles;
85   
86 }
87 /**************************************************************************/
88
89 void AliAOD::CopyData(AliAOD* aod)
90 {
91  //Copys all data from aod, but leaves local type of particles
92  if (aod == 0x0) return;
93  if (aod == this) return;
94  
95  AliAOD& in = *this;
96  
97  fIsRandomized = in.fIsRandomized ;
98  fPrimaryVertexX = in.fPrimaryVertexX ;
99  fPrimaryVertexY = in.fPrimaryVertexY ;
100  fPrimaryVertexZ = in.fPrimaryVertexZ ;
101  fParticleClass = in.fParticleClass ; //althought it is pointer, this points to object in class list of gROOT
102
103  
104  if (in.fParticles == 0x0)
105   {//if in obj has null fParticles we delete ours
106     delete fParticles;
107     fParticles = 0x0;
108   }
109  else
110   { 
111     if (fParticles)
112      { //if ours particles were already created
113        if (fParticles->GetClass() != in.fParticles->GetClass())
114         {//if in obj has 
115           delete fParticles;
116           fParticles = (TClonesArray*)in.fParticles->Clone();
117         }
118        else
119         {
120          //it should be faster than cloning
121           Int_t inentr = in.fParticles->GetEntriesFast();
122           Int_t curentr = fParticles->GetEntriesFast();
123
124           TClonesArray& arr = *fParticles;
125
126           //we have to take care about different sizes of arrays
127           if ( curentr < inentr )
128            {
129              for (Int_t i = 0; i < curentr; i++)
130               {
131                 TObject& inobj = *(in.fParticles->At(i));
132                 TObject& obj = *(fParticles->At(i));
133                 obj = inobj;
134               }
135              
136              TClass* partclass = GetParticleClass();
137              if (partclass == 0x0)
138               {
139                 Fatal("CopyData","Can not get particle class");
140                 return;//pro forma
141               }
142               
143              for (Int_t i = curentr; i < inentr; i++)
144               {
145                 TObject& inobj = *(in.fParticles->At(i));
146                 TObject& obj =  *((TObject*)(partclass->New(arr[i])));
147                 obj = inobj;
148               }
149            }
150           else 
151            {
152              for (Int_t i = 0; i < inentr; i++)
153               {
154                 TObject& inobj = *(in.fParticles->At(i));
155                 TObject& obj = *(fParticles->At(i));
156                 obj = inobj;
157               }
158              
159              for (Int_t i = curentr ; i >= inentr ; i--)
160               {
161                 fParticles->RemoveAt(i);
162               }
163            }
164         } 
165      }
166     else
167      {
168        fParticles = (TClonesArray*)in.fParticles->Clone();
169      } 
170   } 
171  
172 }
173 /**************************************************************************/
174
175 void AliAOD::SetParticleClassName(const char* classname)
176 {
177 //Sets type of particle that is going to be stored 
178   if (gROOT == 0x0) Fatal("SetParticleClassName","ROOT System not initialized");
179   TClass* pclass = gROOT->GetClass(classname);
180   if ( pclass == 0x0 )
181    {
182      Error("SetParticleClass","Can not get TClass for class named %s",classname);
183      return;
184    }
185   SetParticleClass(pclass);
186 }
187 /**************************************************************************/
188
189 void AliAOD::SetParticleClass(TClass* pclass)
190 {
191 //Sets type of particle that is going to be stored 
192
193   if ( pclass == 0x0 )
194    {
195      Error("SetParticleClass","Parameter is NULL.");
196      return;
197    }
198    
199   if ( pclass->InheritsFrom("AliVAODParticle") == kFALSE )
200    {
201      Error("SetParticleClass","Class named %s does not inherit from AliVAODParticle",pclass->GetName());
202      return;
203    }
204   if (pclass != fParticleClass)
205    {
206      fParticleClass = pclass;
207      if (fParticleClass) delete fParticles;
208      fParticles = new TClonesArray(fParticleClass);
209    }
210 }
211 /**************************************************************************/
212 TClass* AliAOD::GetParticleClass()
213 {
214 //returns TClass of particle class
215   if (fParticleClass) return fParticleClass;
216   
217   if (fParticles == 0x0) return 0x0;
218   
219   fParticleClass = fParticles->GetClass();
220   return fParticleClass;
221 }
222
223 /**************************************************************************/
224
225 void  AliAOD::AddParticle(TParticle* part, Int_t idx)
226 {
227   //Adds TParticle to event
228   if (part == 0x0) 
229    {
230      Error("AddParticle(TParticle*,Int_t)","pointer to particle is NULL");
231      return;
232    }
233
234   if (fParticles == 0x0) SetParticleClassName("AliAODParticle");
235   AddParticle( new AliAODParticle(*part,idx) );
236 }
237 /**************************************************************************/
238
239 void  AliAOD::AddParticle(AliVAODParticle* particle)
240 {
241  //add particle to AOD
242  //MAKES ITS OWN COPY OF THE PARTICLE!!! (AOD is not going to keep and delete input pointer)
243  
244   if (fParticles == 0x0) SetParticleClassName("AliAODParticle");
245
246   Int_t idx = fParticles->GetLast() + 1;
247   TClonesArray& arr = *fParticles;
248
249   TClass* partclass = GetParticleClass();
250   if (partclass == 0x0)
251    {
252      Error("AddParticle(AliVAODParticle*)","Can not get particle class");
253      return;
254    }
255   
256   AliVAODParticle* pp = (AliVAODParticle*)partclass->New(arr[idx]);
257   pp->operator=(*particle);
258   
259 }
260 /**************************************************************************/
261
262 void  AliAOD::AddParticle(Int_t pdg, Int_t idx,
263                           Double_t px, Double_t py, Double_t pz, Double_t etot,
264                           Double_t vx, Double_t vy, Double_t vz, Double_t time)
265 {
266   //adds particle to event (standard AOD class)
267
268   if (fParticles == 0x0) SetParticleClassName("AliAODParticle");
269
270   Int_t newpartidx = fParticles->GetLast() + 1;
271   TClonesArray& arr = *fParticles;
272
273   TClass* partclass = GetParticleClass();
274   if (partclass == 0x0)
275    {
276      Error("AddParticle(Int_t,...)","Can not get particle class");
277      return;
278    }
279
280   AliVAODParticle* p =  (AliVAODParticle*)partclass->New(arr[newpartidx]);
281   
282   p->SetPdgCode(pdg);
283   p->SetUID(idx);
284   p->SetMomentum(px,py,pz,etot);
285   p->SetProductionVertex(vx,vy,vz,time);
286   
287 }
288 /**************************************************************************/
289
290 void AliAOD::SwapParticles(Int_t i, Int_t j)
291 {
292 //swaps particles positions; used by AliHBTEvent::Blend
293   if ( (i<0) || (i>=GetNumberOfParticles()) ) return;
294   if ( (j<0) || (j>=GetNumberOfParticles()) ) return;
295   
296
297   TClass* partclass = GetParticleClass();
298   if (partclass == 0x0)
299    {
300      Error("SwapParticles","Can not get particle class");
301      return;
302    }
303     
304   AliVAODParticle* tmpobj = (AliVAODParticle*)partclass->New();
305   AliVAODParticle& tmp = *tmpobj;
306   AliVAODParticle& first = *(GetParticle(i));
307   AliVAODParticle& second = *(GetParticle(j));
308   
309   tmp = first;
310   first = second;
311   second = tmp;
312   delete tmpobj;
313 }
314 /**************************************************************************/
315
316 void  AliAOD::Reset()
317 {
318   //deletes all particles from the event
319    if (fParticles) fParticles->Clear("C");
320    
321    fIsRandomized = kFALSE;
322 }
323 /**************************************************************************/
324
325 void AliAOD::GetPrimaryVertex(Double_t&x, Double_t&y, Double_t&z)
326 {
327 //returns positions of the primary vertex
328   x = fPrimaryVertexX;
329   y = fPrimaryVertexY;
330   z = fPrimaryVertexZ;
331 }
332 /**************************************************************************/
333
334 void AliAOD::SetPrimaryVertex(Double_t x, Double_t y, Double_t z)
335 {
336 //Sets positions of the primary vertex 
337   fPrimaryVertexX = x;
338   fPrimaryVertexY = y;
339   fPrimaryVertexZ = z;
340 }
341 /**************************************************************************/
342
343 Int_t AliAOD::GetNumberOfCharged(Double_t etamin, Double_t etamax) const
344 {
345   //reurns number of charged particles within given pseudorapidity range
346   Int_t n = 0;
347   Int_t npart = GetNumberOfParticles();
348   for (Int_t i = 0; i < npart; i++)
349    {
350      AliVAODParticle* p = GetParticle(i);
351      Double_t eta = p->Eta();
352      if ( (eta < etamin) || (eta > etamax) ) continue;
353      if (p->Charge() != 0.0) n++;
354    }
355   return n;
356 }
357 /**************************************************************************/
358
359 void AliAOD::Move(Double_t x, Double_t y, Double_t z)
360 {
361  //moves all spacial coordinates about this vector
362  // vertex
363  // track points
364  // and whatever will be added to AOD and AOD particles that is a space coordinate
365
366   fPrimaryVertexX += x;
367   fPrimaryVertexY += y;
368   fPrimaryVertexZ += z;
369
370   Int_t npart = GetNumberOfParticles();
371   for (Int_t i = 0; i < npart; i++)
372    {
373      AliVAODParticle* p = GetParticle(i);
374      AliTrackPoints* tp  = p->GetTPCTrackPoints();
375      if (tp) tp->Move(x,y,z);
376      tp  = p->GetITSTrackPoints();
377      if (tp) tp->Move(x,y,z);
378    }
379 }
380 /**************************************************************************/
381
382 void AliAOD::Print(const Option_t* /*option*/) const
383 {
384   //Prints AOD
385   TString ts;
386   TString msg("\n");
387   msg+="Particle Class: ";
388   if (fParticleClass)
389    {
390      msg+=fParticleClass->GetName();
391    }
392   else
393    {
394      msg+="Not specified yet";
395    } 
396   msg += "\n";
397   msg += "Vertex position X: ";
398   msg += fPrimaryVertexX;
399   msg += " Y:" ;
400   msg += fPrimaryVertexY;
401   msg += " Z:";
402   msg += fPrimaryVertexZ;
403   msg += "\n";
404   
405   msg += "Randomized: ";
406   msg += fIsRandomized;
407   msg += "\n";
408   
409   Info("Print","%s",msg.Data());
410   
411   Int_t npart = GetNumberOfParticles();
412   Info("Print","Npart: %d",npart);
413   for (Int_t i = 0; i < npart; i++)
414    {
415      Info("Print","Getting particle %d",i);
416      AliVAODParticle* p = GetParticle(i);
417      Info("Print","Printing particle %d, address %#x",i,p);
418      p->Dump();
419      p->Print();
420      Info("Print","particle %d printed",i);
421    }
422 }
423
424 void AliAOD::SetOwner(Bool_t /*owner*/)
425 {
426 //Sets the ownership of particles: if particles should be also deleted if AOD is deleted/reseted
427 //Since fParticles is Clones and not Object Array, it is always the owner and this method does not have sense
428  
429  MayNotUse("SetOwner");
430  //if fParticles->SetOwner(owner);
431  
432 }