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