]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ANALYSIS/AliReaderAOD.cxx
Remove memory leak (A. Zinchenko)
[u/mrichter/AliRoot.git] / ANALYSIS / AliReaderAOD.cxx
CommitLineData
dd2b6810 1#include "AliReaderAOD.h"
cd319790 2//______________________________________________________________________________
3////////////////////////////////////////////////////////////////////////////////
4// //
5// class AliReaderAOD //
6// //
7// Reader and Writer for AOD format. //
8// AODs are stored in a tree named by the variable fgkTreeName. //
9// There is stored 1 or 2 branches. Each of them stores AOD objects //
10// First branch is named by the variable fgkReconstructedDataBranchName //
11// ("reconstructed.") and keeps reconstructed data. //
12// Second branch is called by the variable fgkSimulatedDataBranchName //
13// ("simulated.") and stores Monte carlo truth. If both branches are present //
14// AODs are parallel, i.e. nth particle in one branch corresponds to the nth //
15// particle in the other one. //
16// //
17// Since we accept different formats of particles that are stored in AODs //
18// reader must take care of that fact: clean buffer if the next file contains //
19// different particle type. //
20// //
21// Piotr.Skowronski@cern.ch //
22// //
23////////////////////////////////////////////////////////////////////////////////
dd2b6810 24
25ClassImp(AliReaderAOD)
26
27#include <TError.h>
28#include <TFile.h>
29#include <TTree.h>
cd319790 30#include <TH1.h>
dd2b6810 31#include "AliAOD.h"
32
efdb0cc9 33
34const TString AliReaderAOD::fgkTreeName("TAOD");
cd319790 35const TString AliReaderAOD::fgkReconstructedDataBranchName("reconstructed.");
efdb0cc9 36const TString AliReaderAOD::fgkSimulatedDataBranchName("simulated.");
37
38AliReaderAOD::AliReaderAOD(const Char_t* aodfilename):
39 fFileName(aodfilename),
40 fReadSim(kFALSE),
cd319790 41 fReadRec(kTRUE),
efdb0cc9 42 fTree(0x0),
43 fFile(0x0),
44 fSimBuffer(0x0),
45 fRecBuffer(0x0)
46{
47 //ctor
48}
49/********************************************************************/
50
51AliReaderAOD::~AliReaderAOD()
52{
53//dtor
54 if (fEventSim == fSimBuffer )
55 {
56 fEventSim = 0x0;
57 fEventRec = 0x0;
58 }
59 delete fSimBuffer;
60 delete fRecBuffer;
61
62 delete fTree;
63 delete fFile;
64}
65/********************************************************************/
66
67void AliReaderAOD::Rewind()
68{
69//Rewinds reading
70 delete fTree;
71 fTree = 0x0;
72 delete fFile;
73 fFile = 0x0;
74 fCurrentDir = 0;
75 fNEventsRead= 0;
76}
77/********************************************************************/
78Int_t AliReaderAOD::ReadNext()
79{
80//Reads next event
81
d9122a01 82 Info("ReadNext","Entered");
efdb0cc9 83 do //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
84 {
85 if (fFile == 0x0)
86 {
cd319790 87 Int_t openfailed = OpenFile(fCurrentDir);//rl is opened here
88 if (openfailed)
efdb0cc9 89 {
90 //Error("ReadNext","Error Occured while opening directory number %d",fCurrentDir);
91 fCurrentDir++;
92 continue;
93 }
94 fCurrentEvent = 0;
95 }
cd319790 96 //Tree must exist because OpenFile would reuturn error in the other case
97 if ( fCurrentEvent >= fTree->GetEntries() )
98 {
99 delete fTree;
100 fTree = 0x0;
101 delete fFile;
102 fFile = 0x0;
103 fSimBuffer = 0x0;
104 fRecBuffer = 0x0;
105 fCurrentDir++;
106 continue;
107 }
108
d9122a01 109 Info("ReadNext","Getting event %d",fCurrentEvent);
efdb0cc9 110 fTree->GetEvent(fCurrentEvent);
d9122a01 111 Info("ReadNext","Getting event %d Done",fCurrentEvent);
efdb0cc9 112
cd319790 113 Int_t retval = 0;
114 if (fReadRec && fReadSim)
efdb0cc9 115 {
cd319790 116 retval = ReadRecAndSim();
efdb0cc9 117 }
cd319790 118 else
119 {
120 if (fReadRec) retval = ReadRec();
121 if (fReadSim) retval = ReadSim();
122 }
efdb0cc9 123
cd319790 124 fCurrentEvent++;
125 if (retval == 0) fNEventsRead++;
efdb0cc9 126
cd319790 127 return retval;//success -> read one event
efdb0cc9 128
129 }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array
130
131 return 1; //no more directories to read
132
133
134}
135/********************************************************************/
136
cd319790 137Int_t AliReaderAOD::ReadRecAndSim()
138{
139//Reads raconstructed and simulated data
140
141 Info("ReadRecAndSim","Found %d reconstructed tracks and %d simulated particles",
142 fRecBuffer->GetNumberOfParticles(),fSimBuffer->GetNumberOfParticles());
143
144 if (fCuts->GetEntriesFast() == 0x0)
145 {//if there is no cuts we return pointer to the buffer
146 if (fEventRec != fRecBuffer)
147 {
148 delete fEventRec;
149 delete fEventSim;
150 }
151 fEventRec = fRecBuffer;//fEventRec is the pointer that the user gets when he asks about an event
152 fEventSim = fSimBuffer;
153 }
154 else
155 {//if there are cuts specified
156 if ( (fEventRec == 0x0) || (fEventRec == fRecBuffer) )
157 {//we need to create a new event, if it is not existing or it is the same as branch buffer
158 fEventRec = new AliAOD();
159 fEventSim = new AliAOD();
160
161 fEventRec->SetParticleClass( fRecBuffer->GetParticleClass() );
162 fEventSim->SetParticleClass( fSimBuffer->GetParticleClass() );
163 }
164 else
165 {//or simply reset it in case it already exists
166 fEventRec->Reset();
167 fEventSim->Reset();
168 }
169
170 Int_t npart = fRecBuffer->GetNumberOfParticles();
171 for (Int_t i = 0; i < npart; i++)
172 {
173 AliVAODParticle* prec = fRecBuffer->GetParticle(i);
174 if (Rejected(prec)) continue;//we make cuts only on simulated data
175
176 fEventRec->AddParticle(prec);
177 fEventSim->AddParticle( fSimBuffer->GetParticle(i));
178 }
179 }
180
181 Info("ReadRecAndSim","Read %d reconstructed tracks and %d simulated particles",
182 fEventRec->GetNumberOfParticles(),fEventSim->GetNumberOfParticles());
183
184 fTrackCounter->Fill(fEventRec->GetNumberOfParticles());
185
186 return 0;
187}
188/********************************************************************/
189
190Int_t AliReaderAOD::ReadRec()
191{
192//Reads reconstructed data only
193
194 Info("ReadRec","Found %d reconstructed tracks",fRecBuffer->GetNumberOfParticles());
195
196 if (fCuts->GetEntriesFast() == 0x0)
197 {//if there is no cuts we return pointer to the buffer
198 if (fEventRec != fRecBuffer)
199 {
200 delete fEventRec;
201 }
202 fEventRec = fRecBuffer;//fEventRec is the pointer that the user gets when he asks about an event
203 }
204 else
205 {//if there are cuts specified
206 if ( (fEventRec == 0x0) || (fEventRec == fRecBuffer) )
207 {//we need to create a new event, if it is not existing or it is the same as branch buffer
208 fEventRec = new AliAOD();
209
210 fEventRec->SetParticleClass( fRecBuffer->GetParticleClass() );
211 }
212 else
213 {//or simply reset it in case it already exists
214 fEventRec->Reset();
215 }
216
217 Int_t npart = fRecBuffer->GetNumberOfParticles();
218 for (Int_t i = 0; i < npart; i++)
219 {
220 AliVAODParticle* prec = fRecBuffer->GetParticle(i);
221 if (Rejected(prec)) continue;//we make cuts only on simulated data
222
223 fEventRec->AddParticle(prec);
224 }
225 }
226
227 Info("ReadRec","Read %d reconstructed tracks",fEventRec->GetNumberOfParticles());
228 fTrackCounter->Fill(fEventRec->GetNumberOfParticles());
229
230 return 0;
231}
232/********************************************************************/
233
234Int_t AliReaderAOD::ReadSim()
235{
236//Reads simulated data only
237
238 Info("ReadSim","Found %d simulated particles",fSimBuffer->GetNumberOfParticles());
239
240 if (fCuts->GetEntriesFast() == 0x0)
241 {//if there is no cuts we return pointer to the buffer
242 if (fEventSim != fSimBuffer)
243 {
244 delete fEventSim;
245 }
246 fEventSim = fSimBuffer;
247 }
248 else
249 {//if there are cuts specified
250 if ( (fEventSim == 0x0) || (fEventSim == fSimBuffer) )
251 {//we need to create a new event, if it is not existing or it is the same as branch buffer
252 fEventSim = new AliAOD();
253
254 fEventSim->SetParticleClass( fSimBuffer->GetParticleClass() );
255 }
256 else
257 {//or simply reset it in case it already exists
258 fEventSim->Reset();
259 }
260
261 Int_t npart = fSimBuffer->GetNumberOfParticles();
262 for (Int_t i = 0; i < npart; i++)
263 {
264 AliVAODParticle* prec = fSimBuffer->GetParticle(i);
265 if (Rejected(prec)) continue;//we make cuts only on simulated data
266 fEventSim->AddParticle(prec);
267 }
268 }
269
270 Info("ReadSim","Read %d simulated particles",fEventSim->GetNumberOfParticles());
271 fTrackCounter->Fill(fEventSim->GetNumberOfParticles());
272
273
274 return 0;
275}
276/********************************************************************/
277
efdb0cc9 278Int_t AliReaderAOD::OpenFile(Int_t n)
279{
280//opens fFile with tree
281
d9122a01 282// Info("ReadNext","Opening File %d",n);
283 const TString dirname = GetDirName(n);
efdb0cc9 284 if (dirname == "")
285 {
286 if (AliVAODParticle::GetDebug() > 2 )
287 {
288 Info("OpenFile","Got empty string as a directory name.");
289 }
290 return 1;
291 }
292
293 TString filename = dirname +"/"+ fFileName;
294 fFile = TFile::Open(filename.Data());
295 if ( fFile == 0x0)
296 {
297 Error("OpenFile","Can't open fFile %s",filename.Data());
298 return 2;
299 }
300 if (!fFile->IsOpen())
301 {
302 Error("OpenFile","Can't open fFile %s",filename.Data());
303 delete fFile;
304 fFile = 0x0;
305 return 3;
306 }
d9122a01 307
f4df64d6 308 Info("ReadNext","File %s Is Opened, Getting the TREE",filename.Data());
efdb0cc9 309
310 fTree = dynamic_cast<TTree*>(fFile->Get(fgkTreeName));
311 if (fTree == 0x0)
312 {
313 if (AliVAODParticle::GetDebug() > 2 )
314 {
315 Info("ReadNext","Can not find TTree object named %s",fgkTreeName.Data());
316 }
cd319790 317 delete fFile;
efdb0cc9 318 fFile = 0x0;
319 return 4;
320 }
321
d9122a01 322// Info("ReadNext","Got TREE, Setting branch addresses");
323
cd319790 324 if (fReadRec)
325 {
326 TBranch* branch = fTree->GetBranch(fgkReconstructedDataBranchName);
327 if (branch == 0x0)
328 {
329 Error("OpenFile","Can not find branch %s in file %s",
330 fgkReconstructedDataBranchName.Data(),filename.Data());
331
332 delete fTree;
333 fTree = 0x0;
334 delete fFile;
335 fFile = 0x0;
336 return 5;
337 }
338 fTree->SetBranchAddress(fgkReconstructedDataBranchName,&fRecBuffer);
339 }
efdb0cc9 340
cd319790 341
342 if (fReadSim)
343 {
344 TBranch* branch = fTree->GetBranch(fgkSimulatedDataBranchName);
345 if (branch == 0x0)
346 {
347 Error("OpenFile","Can not find branch %s in file %s",
348 fgkSimulatedDataBranchName.Data(),filename.Data());
349
350 delete fTree;
351 fTree = 0x0;
352 delete fFile;
353 fFile = 0x0;
354 return 6;
355 }
356 fTree->SetBranchAddress(fgkSimulatedDataBranchName,&fSimBuffer);
357 }
d9122a01 358// Info("ReadNext","Got TREE, Addresses are set.");
359// Info("ReadNext","Quitting the method.");
360
efdb0cc9 361 return 0;
362
363}
efdb0cc9 364/********************************************************************/
365
beb1c41d 366Int_t AliReaderAOD::WriteAOD(AliReader* reader, const char* outfilename, const char* pclassname, Bool_t /*multcheck*/)
dd2b6810 367{
368//reads tracks from runs and writes them to file
369 ::Info("AliReaderAOD::Write","________________________________________________________");
370 ::Info("AliReaderAOD::Write","________________________________________________________");
371 ::Info("AliReaderAOD::Write","________________________________________________________");
372
373 if (reader == 0x0)
374 {
375 ::Error("AliReaderAOD::Write","Input Reader is NULL");
376 return -1;
377 }
378 TFile *outfile = TFile::Open(outfilename,"recreate");
379 if (outfile == 0x0)
380 {
381 ::Error("AliReaderAOD::Write","Can not open output file %s",outfilename);
382 return -1;
383 }
384
efdb0cc9 385 TTree *tree = new TTree(fgkTreeName,"Tree with tracks");
dd2b6810 386
387 TBranch *recbranch = 0x0, *simbranch = 0x0;
dd2b6810 388
beb1c41d 389 AliAOD* eventrec = new AliAOD();//must be created before Branch is called. Otherwise clones array is not splitted
390 AliAOD* eventsim = new AliAOD();//AOD together with fParticles clones array knowing exact type of particles
391
392 eventrec->SetParticleClassName(pclassname);
393 eventsim->SetParticleClassName(pclassname);
efdb0cc9 394
395 AliAOD* recbuffer = eventrec;
396 AliAOD* simbuffer = eventsim;
beb1c41d 397
cd319790 398 if (reader->ReadsRec()) recbranch = tree->Branch(fgkReconstructedDataBranchName,"AliAOD",&recbuffer,32000,99);
efdb0cc9 399 if (reader->ReadsSim()) simbranch = tree->Branch(fgkSimulatedDataBranchName,"AliAOD",&simbuffer,32000,99);
dd2b6810 400
401 reader->Rewind();
402 while (reader->Next() == kFALSE)
403 {
404
dd2b6810 405 if (reader->ReadsRec())
efdb0cc9 406 {//here we can get AOD that has different particle type
407 AliAOD* event = reader->GetEventRec();
408 if ( eventrec->GetParticleClass() != event->GetParticleClass() )
409 {//if class type is not what what we whant we copy particles
410 eventrec->CopyData(event);
411 recbuffer = eventrec;
412 }
413 else
414 {//else just pointer to event from input reader is passed
415 recbuffer = event;
416 }
beb1c41d 417 }
418
419 if (reader->ReadsSim())
420 {
d9122a01 421 AliAOD* event = reader->GetEventSim();
efdb0cc9 422 if ( eventsim->GetParticleClass() != event->GetParticleClass() )
423 {//if class type is not what what we whant we copy particles
424 eventsim->CopyData(event);
425 simbuffer = eventrec;
426 }
427 else
428 {//else just pointer to event from input reader is passed
429 simbuffer = event;
430 }
dd2b6810 431 }
432 tree->Fill();
dd2b6810 433 }
434
435 ::Info("AliReaderAOD::Write","Written %d events",tree->GetEntries());
436 outfile->cd();
437 tree->Write();
efdb0cc9 438
439 delete eventsim;
440 delete eventrec;
441
dd2b6810 442 delete tree;
443 delete outfile;
444 return 0;
445}
446