2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 /////////////////////////////////////////////////////////////
21 // base class for AOD containers
23 /////////////////////////////////////////////////////////////
26 #include <TParticle.h>
29 #include "AliAODParticle.h"
30 #include "AliTrackPoints.h"
36 fIsRandomized(kFALSE),
43 // Info("AliAOD()","Entered");
45 // Info("AliAOD()","Exited");
47 /**************************************************************************/
52 //fParticleClass does not belong to AliAOD -> Do not delete it
56 /**************************************************************************/
58 void AliAOD::SetParticleClassName(const char* classname)
60 //Sets type of particle that is going to be stored
61 if (gROOT == 0x0) Fatal("SetParticleClassName","ROOT System not initialized");
62 TClass* pclass = gROOT->GetClass(classname);
65 Error("SetParticleClass","Can not get TClass for class named %s",classname);
68 SetParticleClass(pclass);
70 /**************************************************************************/
72 void AliAOD::SetParticleClass(TClass* pclass)
74 //Sets type of particle that is going to be stored
78 Error("SetParticleClass","Parameter is NULL.");
82 if ( pclass->InheritsFrom("AliVAODParticle") == kFALSE )
84 Error("SetParticleClass","Class named %s does not inherit from AliVAODParticle",pclass->GetName());
87 if (pclass != fParticleClass)
89 fParticleClass = pclass;
90 if (fParticleClass) delete fParticles;
91 fParticles = new TClonesArray(fParticleClass);
95 /**************************************************************************/
97 void AliAOD::AddParticle(TParticle* part, Int_t idx)
99 //Adds TParticle to event
102 Error("AddParticle(TParticle*,Int_t)","pointer to particle is NULL");
106 if (fParticles == 0x0) SetParticleClassName("AliAODParticle");
107 AddParticle( new AliAODParticle(*part,idx) );
109 /**************************************************************************/
111 void AliAOD::AddParticle(AliVAODParticle* particle)
113 //add particle to AOD
114 //MAKES ITS OWN COPY OF THE PARTICLE!!! (AOD is not going to keep and delete input pointer)
116 if (fParticles == 0x0) SetParticleClassName("AliAODParticle");
118 Int_t idx = fParticles->GetLast() + 1;
119 TClonesArray& arr = *fParticles;
121 AliVAODParticle* pp = (AliVAODParticle*)fParticleClass->New(arr[idx]);
122 pp->operator=(*particle);
125 /**************************************************************************/
127 void AliAOD::AddParticle(Int_t pdg, Int_t idx,
128 Double_t px, Double_t py, Double_t pz, Double_t etot,
129 Double_t vx, Double_t vy, Double_t vz, Double_t time)
131 //adds particle to event (standard AOD class)
133 if (fParticles == 0x0) SetParticleClassName("AliAODParticle");
135 Int_t newpartidx = fParticles->GetLast() + 1;
136 TClonesArray& arr = *fParticles;
138 AliVAODParticle* p = (AliVAODParticle*)fParticleClass->New(arr[newpartidx]);
142 p->SetMomentum(px,py,pz,etot);
143 p->SetProductionVertex(vx,vy,vz,time);
146 /**************************************************************************/
148 void AliAOD::SwapParticles(Int_t i, Int_t j)
150 //swaps particles positions; used by AliHBTEvent::Blend
151 if ( (i<0) || (i>=GetNumberOfParticles()) ) return;
152 if ( (j<0) || (j>=GetNumberOfParticles()) ) return;
154 AliVAODParticle* tmpobj = (AliVAODParticle*)fParticleClass->New();
155 AliVAODParticle& tmp = *tmpobj;
156 AliVAODParticle& first = *(GetParticle(i));
157 AliVAODParticle& second = *(GetParticle(j));
164 /**************************************************************************/
168 //deletes all particles from the event
169 if (fParticles) fParticles->Clear("C");
171 fIsRandomized = kFALSE;
173 /**************************************************************************/
175 void AliAOD::GetPrimaryVertex(Double_t&x, Double_t&y, Double_t&z)
177 //returns positions of the primary vertex
182 /**************************************************************************/
184 void AliAOD::SetPrimaryVertex(Double_t x, Double_t y, Double_t z)
186 //Sets positions of the primary vertex
191 /**************************************************************************/
193 Int_t AliAOD::GetNumberOfCharged(Double_t etamin, Double_t etamax) const
195 //reurns number of charged particles within given pseudorapidity range
197 Int_t npart = GetNumberOfParticles();
198 for (Int_t i = 0; i < npart; i++)
200 AliVAODParticle* p = GetParticle(i);
201 Double_t eta = p->Eta();
202 if ( (eta < etamin) || (eta > etamax) ) continue;
203 if (p->Charge() != 0.0) n++;
207 /**************************************************************************/
209 void AliAOD::Move(Double_t x, Double_t y, Double_t z)
211 //moves all spacial coordinates about this vector
214 // and whatever will be added to AOD and AOD particles that is a space coordinate
216 fPrimaryVertexX += x;
217 fPrimaryVertexY += y;
218 fPrimaryVertexZ += z;
220 Int_t npart = GetNumberOfParticles();
221 for (Int_t i = 0; i < npart; i++)
223 AliVAODParticle* p = GetParticle(i);
224 AliTrackPoints* tp = p->GetTPCTrackPoints();
225 if (tp) tp->Move(x,y,z);
226 tp = p->GetITSTrackPoints();
227 if (tp) tp->Move(x,y,z);
231 void AliAOD::Print(Option_t* /*option*/)
236 msg+="Particle Class: ";
239 msg+=fParticleClass->GetName();
243 msg+="Not specified yet";
246 msg += "Vertex position X: ";
247 msg += fPrimaryVertexX;
249 msg += fPrimaryVertexY;
251 msg += fPrimaryVertexZ;
254 msg += "Randomized: ";
255 msg += fIsRandomized;
258 Info("Print","%s",msg.Data());
260 Int_t npart = GetNumberOfParticles();
261 Info("Print","Npart: %d",npart);
262 for (Int_t i = 0; i < npart; i++)
264 Info("Print","Getting particle %d",i);
265 AliVAODParticle* p = GetParticle(i);
266 Info("Print","Printing particle %d, address %#x",i,p);
269 Info("Print","particle %d printed",i);
273 void AliAOD::SetOwner(Bool_t /*owner*/)
275 //Sets the ownership of particles: if particles should be also deleted if AOD is deleted/reseted
276 //Since fParticles is Clones and not Object Array, it is always the owner and this method does not have sense
278 MayNotUse("SetOwner");
279 //if fParticles->SetOwner(owner);