]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemto/AliFemtoSimpleAnalysis.cxx
Reaction Plane analysis updates
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemto / AliFemtoSimpleAnalysis.cxx
1 ///////////////////////////////////////////////////////////////////////////
2 //                                                                       //
3 // AliFemtoSimpleAnalysis - the most basic analysis there is. All other        //
4 // inherit from this one. Provides basic functionality for the analysis. //
5 // To properly set up the analysis the following steps should be taken:  //
6 //                                                                       //
7 // - create particle cuts and add them via SetFirstParticleCut and       //
8 //   SetSecondParticleCut. If one analyzes identical particle            //
9 //   correlations, the first particle cut must be also the second        //
10 //   particle cut.                                                       //
11 //                                                                       //
12 // - create pair cuts and add them via SetPairCut                        //
13 //                                                                       //
14 // - create one or many correlation functions and add them via           //
15 //   AddCorrFctn method.                                                 //
16 //                                                                       //
17 // - specify how many events are to be strored in the mixing buffer for  //
18 //   background construction                                             //
19 //                                                                       //
20 // Then, when the analysis is run, for each event, the EventBegin is     //
21 // called before any processing is done, then the ProcessEvent is called //
22 // which takes care of creating real and mixed pairs and sending them    //
23 // to all the registered correlation functions. At the end of each event,//
24 // after all pairs are processed, EventEnd is called. After the whole    //
25 // analysis finishes (there is no more events to process) Finish() is    //
26 // called.                                                               //
27 //                                                                       //
28 ///////////////////////////////////////////////////////////////////////////
29 #include "AliFemtoSimpleAnalysis.h"
30 #include "AliFemtoTrackCut.h"
31 #include "AliFemtoV0Cut.h"
32 #include "AliFemtoKinkCut.h"
33 #include <string>
34 #include <iostream>
35
36 // blah blah
37
38 #ifdef __ROOT__ 
39 ClassImp(AliFemtoSimpleAnalysis)
40 #endif
41
42 AliFemtoEventCut*    copyTheCut(AliFemtoEventCut*);
43 AliFemtoParticleCut* copyTheCut(AliFemtoParticleCut*);
44 AliFemtoPairCut*     copyTheCut(AliFemtoPairCut*);
45 AliFemtoCorrFctn*    copyTheCorrFctn(AliFemtoCorrFctn*);
46
47 // this little function used to apply ParticleCuts (TrackCuts or V0Cuts) and fill ParticleCollections of picoEvent
48 //  it is called from AliFemtoSimpleAnalysis::ProcessEvent()
49 void FillHbtParticleCollection(AliFemtoParticleCut*         partCut,
50                                AliFemtoEvent*               hbtEvent,
51                                AliFemtoParticleCollection*  partCollection)
52 {
53   // Fill particle collections from the event
54   // by the particles that pass all the cuts
55   switch (partCut->Type()) {
56   case hbtTrack:       // cut is cutting on Tracks
57     {
58       AliFemtoTrackCut* pCut = (AliFemtoTrackCut*) partCut;
59       AliFemtoTrack* pParticle;
60       AliFemtoTrackIterator pIter;
61       AliFemtoTrackIterator startLoop = hbtEvent->TrackCollection()->begin();
62       AliFemtoTrackIterator endLoop   = hbtEvent->TrackCollection()->end();
63       for (pIter=startLoop;pIter!=endLoop;pIter++){
64         pParticle = *pIter;
65         bool tmpPassParticle = pCut->Pass(pParticle);
66         pCut->FillCutMonitor(pParticle, tmpPassParticle);
67         if (tmpPassParticle){   
68           AliFemtoParticle* particle = new AliFemtoParticle(pParticle,pCut->Mass());
69           partCollection->push_back(particle);
70         }
71       }
72       break;
73     }
74   case hbtV0:          // cut is cutting on V0s
75     {
76       AliFemtoV0Cut* pCut = (AliFemtoV0Cut*) partCut;
77       AliFemtoV0* pParticle;
78       AliFemtoV0Iterator pIter;
79       AliFemtoV0Iterator startLoop = hbtEvent->V0Collection()->begin();
80       AliFemtoV0Iterator endLoop   = hbtEvent->V0Collection()->end();
81       // this following "for" loop is identical to the one above, but because of scoping, I can's see how to avoid repitition...
82       for (pIter=startLoop;pIter!=endLoop;pIter++){
83         pParticle = *pIter; 
84         bool tmpPassV0 = pCut->Pass(pParticle);
85         pCut->FillCutMonitor(pParticle,tmpPassV0);
86         if (tmpPassV0){
87           AliFemtoParticle* particle = new AliFemtoParticle(pParticle,partCut->Mass());
88           partCollection->push_back(particle);
89         }
90       }
91       pCut->FillCutMonitor(hbtEvent,partCollection);// Gael 19/06/02
92       
93       break;
94     }
95   case hbtKink:          // cut is cutting on Kinks  -- mal 25May2001
96     {
97       AliFemtoKinkCut* pCut = (AliFemtoKinkCut*) partCut;
98       AliFemtoKink* pParticle;
99       AliFemtoKinkIterator pIter;
100       AliFemtoKinkIterator startLoop = hbtEvent->KinkCollection()->begin();
101       AliFemtoKinkIterator endLoop   = hbtEvent->KinkCollection()->end();
102       // this following "for" loop is identical to the one above, but because of scoping, I can's see how to avoid repitition...
103       for (pIter=startLoop;pIter!=endLoop;pIter++){
104         pParticle = *pIter; 
105         bool tmpPass = pCut->Pass(pParticle);
106         pCut->FillCutMonitor(pParticle,tmpPass);
107         if (tmpPass){
108           AliFemtoParticle* particle = new AliFemtoParticle(pParticle,partCut->Mass());
109           partCollection->push_back(particle);
110         }
111       }
112       break;
113     }
114   default:
115     cout << "FillHbtParticleCollection function (in AliFemtoSimpleAnalysis.cxx) - undefined Particle Cut type!!! \n";
116   }
117 }
118 //____________________________
119 AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis() :
120   fPicoEventCollectionVectorHideAway(0), 
121   fPairCut(0),            
122   fCorrFctnCollection(0), 
123   fEventCut(0),           
124   fFirstParticleCut(0),   
125   fSecondParticleCut(0),  
126   fMixingBuffer(0),       
127   fPicoEvent(0),          
128   fNumEventsToMix(0),                     
129   fNeventsProcessed(0),                   
130   fMinSizePartCollection(0)
131 {
132   // Default constructor
133   //  mControlSwitch     = 0;
134   fCorrFctnCollection = new AliFemtoCorrFctnCollection;
135   fMixingBuffer = new AliFemtoPicoEventCollection;
136 }
137 //____________________________
138
139 AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) : 
140   AliFemtoAnalysis(),
141   fPicoEventCollectionVectorHideAway(0), 
142   fPairCut(0),            
143   fCorrFctnCollection(0), 
144   fEventCut(0),           
145   fFirstParticleCut(0),   
146   fSecondParticleCut(0),  
147   fMixingBuffer(0),       
148   fPicoEvent(0),          
149   fNumEventsToMix(0),                     
150   fNeventsProcessed(0),                   
151   fMinSizePartCollection(0)
152 {
153   // Copy constructor
154   //AliFemtoSimpleAnalysis();
155   fCorrFctnCollection = new AliFemtoCorrFctnCollection;
156   fMixingBuffer = new AliFemtoPicoEventCollection;
157
158   // find the right event cut
159   fEventCut = a.fEventCut->Clone();
160   // find the right first particle cut
161   fFirstParticleCut = a.fFirstParticleCut->Clone();
162   // find the right second particle cut
163   if (a.fFirstParticleCut==a.fSecondParticleCut) 
164     SetSecondParticleCut(fFirstParticleCut); // identical particle hbt
165   else
166   fSecondParticleCut = a.fSecondParticleCut->Clone();
167
168   fPairCut = a.fPairCut->Clone();
169   
170   if ( fEventCut ) {
171       SetEventCut(fEventCut); // this will set the myAnalysis pointer inside the cut
172       cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - event cut set " << endl;
173   }
174   if ( fFirstParticleCut ) {
175       SetFirstParticleCut(fFirstParticleCut); // this will set the myAnalysis pointer inside the cut
176       cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - first particle cut set " << endl;
177   }
178   if ( fSecondParticleCut ) {
179       SetSecondParticleCut(fSecondParticleCut); // this will set the myAnalysis pointer inside the cut
180       cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - second particle cut set " << endl;
181   }  if ( fPairCut ) {
182       SetPairCut(fPairCut); // this will set the myAnalysis pointer inside the cut
183       cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - pair cut set " << endl;
184   }
185
186   AliFemtoCorrFctnIterator iter;
187   for (iter=a.fCorrFctnCollection->begin(); iter!=a.fCorrFctnCollection->end();iter++){
188     cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - looking for correlation functions " << endl;
189     AliFemtoCorrFctn* fctn = (*iter)->Clone();
190     if (fctn) AddCorrFctn(fctn);
191     else cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - correlation function not found " << endl;
192   }
193
194   fNumEventsToMix = a.fNumEventsToMix;
195
196   fMinSizePartCollection = a.fMinSizePartCollection;  // minimum # particles in ParticleCollection
197
198   cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - analysis copied " << endl;
199
200 }
201 //____________________________
202 AliFemtoSimpleAnalysis::~AliFemtoSimpleAnalysis(){
203   // destructor
204   cout << " AliFemtoSimpleAnalysis::~AliFemtoSimpleAnalysis()" << endl;
205   if (fEventCut) delete fEventCut; fEventCut=0;
206   if (fFirstParticleCut == fSecondParticleCut) fSecondParticleCut=0;
207   if (fFirstParticleCut)  delete fFirstParticleCut; fFirstParticleCut=0;
208   if (fSecondParticleCut) delete fSecondParticleCut; fSecondParticleCut=0;
209   if (fPairCut) delete fPairCut; fPairCut=0;
210   // now delete every CorrFunction in the Collection, and then the Collection itself
211   AliFemtoCorrFctnIterator iter;
212   for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
213     delete *iter;
214   }
215   delete fCorrFctnCollection;
216   // now delete every PicoEvent in the EventMixingBuffer and then the Buffer itself
217   if (fMixingBuffer) {
218     AliFemtoPicoEventIterator piter;
219     for (piter=fMixingBuffer->begin();piter!=fMixingBuffer->end();piter++){
220       delete *piter;
221     }
222     delete fMixingBuffer;
223   }
224 }
225 //______________________
226 AliFemtoSimpleAnalysis& AliFemtoSimpleAnalysis::operator=(const AliFemtoSimpleAnalysis& aAna) 
227 {
228   // Assignment operator
229   if (this == &aAna)
230     return *this;
231
232   if (fCorrFctnCollection) delete fCorrFctnCollection;
233   fCorrFctnCollection = new AliFemtoCorrFctnCollection;
234   if (fMixingBuffer) delete fMixingBuffer;
235   fMixingBuffer = new AliFemtoPicoEventCollection;
236
237   // find the right event cut
238   if (fEventCut) delete fEventCut;
239   fEventCut = aAna.fEventCut->Clone();
240   // find the right first particle cut
241   if (fFirstParticleCut) delete fFirstParticleCut;
242   fFirstParticleCut = aAna.fFirstParticleCut->Clone();
243   // find the right second particle cut
244   if (fSecondParticleCut) delete fSecondParticleCut;
245   if (aAna.fFirstParticleCut==aAna.fSecondParticleCut) 
246     SetSecondParticleCut(fFirstParticleCut); // identical particle hbt
247   else
248     fSecondParticleCut = aAna.fSecondParticleCut->Clone();
249
250   if (fPairCut) delete fPairCut;
251   fPairCut = aAna.fPairCut->Clone();
252   
253   if ( fEventCut ) {
254     SetEventCut(fEventCut); // this will set the myAnalysis pointer inside the cut
255   }
256   if ( fFirstParticleCut ) {
257     SetFirstParticleCut(fFirstParticleCut); // this will set the myAnalysis pointer inside the cut
258   }
259   if ( fSecondParticleCut ) {
260     SetSecondParticleCut(fSecondParticleCut); // this will set the myAnalysis pointer inside the cut
261   }  
262   if ( fPairCut ) {
263     SetPairCut(fPairCut); // this will set the myAnalysis pointer inside the cut
264   }
265
266   AliFemtoCorrFctnIterator iter;
267   for (iter=aAna.fCorrFctnCollection->begin(); iter!=aAna.fCorrFctnCollection->end();iter++){
268     AliFemtoCorrFctn* fctn = (*iter)->Clone();
269     if (fctn) AddCorrFctn(fctn);
270   }
271
272   fNumEventsToMix = aAna.fNumEventsToMix;
273
274   fMinSizePartCollection = aAna.fMinSizePartCollection;  // minimum # particles in ParticleCollection
275
276   return *this;
277 }
278 //______________________
279 AliFemtoCorrFctn* AliFemtoSimpleAnalysis::CorrFctn(int n){  
280   // return pointer to n-th correlation function
281   if ( n<0 || n > (int)fCorrFctnCollection->size() )
282     return NULL;
283   AliFemtoCorrFctnIterator iter=fCorrFctnCollection->begin();
284   for (int i=0; i<n ;i++){
285     iter++;
286   }
287   return *iter;
288 }
289 //____________________________
290 AliFemtoString AliFemtoSimpleAnalysis::Report()
291 {
292   // Create a simple report from the analysis execution
293   cout << "AliFemtoSimpleAnalysis - constructing Report..."<<endl;
294   string temp = "-----------\nHbt Analysis Report:\n";
295   temp += "\nEvent Cuts:\n";
296   temp += fEventCut->Report();
297   temp += "\nParticle Cuts - First Particle:\n";
298   temp += fFirstParticleCut->Report();
299   temp += "\nParticle Cuts - Second Particle:\n";
300   temp += fSecondParticleCut->Report();
301   temp += "\nPair Cuts:\n";
302   temp += fPairCut->Report();
303   temp += "\nCorrelation Functions:\n";
304   AliFemtoCorrFctnIterator iter;
305   if ( fCorrFctnCollection->size()==0 ) {
306     cout << "AliFemtoSimpleAnalysis-Warning : no correlations functions in this analysis " << endl;
307   }
308   for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
309     temp += (*iter)->Report();
310     temp += "\n";
311   }
312   temp += "-------------\n";
313   AliFemtoString returnThis=temp;
314   return returnThis;
315 }
316 //_________________________
317 void AliFemtoSimpleAnalysis::ProcessEvent(const AliFemtoEvent* hbtEvent) {
318   // Add event to processed events
319
320   fPicoEvent=0; // we will get a new pico event, if not prevent corr. fctn to access old pico event
321   AddEventProcessed();
322   // startup for EbyE 
323   EventBegin(hbtEvent);  
324   // event cut and event cut monitor
325   bool tmpPassEvent = fEventCut->Pass(hbtEvent);
326   if (!tmpPassEvent) 
327     fEventCut->FillCutMonitor(hbtEvent, tmpPassEvent);
328   if (tmpPassEvent) {
329     //     cout << "AliFemtoSimpleAnalysis::ProcessEvent() - Event has passed cut - build picoEvent from " <<
330 //       hbtEvent->TrackCollection()->size() << " tracks in TrackCollection" << endl;
331     //  cout << "Event has passed cut with " << hbtEvent->TrackCollection()->size() << " tracks" << endl;
332     // OK, analysis likes the event-- build a pico event from it, using tracks the analysis likes...
333     fPicoEvent = new AliFemtoPicoEvent; // this is what we will make pairs from and put in Mixing Buffer
334     // no memory leak. we will delete picoevents when they come out of the mixing buffer
335     FillHbtParticleCollection(fFirstParticleCut,(AliFemtoEvent*)hbtEvent,fPicoEvent->FirstParticleCollection());
336     if ( !(AnalyzeIdenticalParticles()) )
337       FillHbtParticleCollection(fSecondParticleCut,(AliFemtoEvent*)hbtEvent,fPicoEvent->SecondParticleCollection());
338     //cout <<"AliFemtoSimpleAnalysis::ProcessEvent - #particles in First, Second Collections: " <<
339 //       fPicoEvent->FirstParticleCollection()->size() << " " <<
340 //       fPicoEvent->SecondParticleCollection()->size() << endl;
341     
342     // cout << "#particles in Collection 1, 2: " <<
343     //   fPicoEvent->FirstParticleCollection()->size() << " " <<
344     //   fPicoEvent->SecondParticleCollection()->size() << endl;
345
346     fEventCut->FillCutMonitor(fPicoEvent->FirstParticleCollection(),fPicoEvent->SecondParticleCollection()); //MJ!
347     
348     
349     // mal - implement a switch which allows only using events with ParticleCollections containing a minimum
350     // number of entries (jun2002)
351     if ((fPicoEvent->FirstParticleCollection()->size() >= fMinSizePartCollection )
352         && ( AnalyzeIdenticalParticles() || (fPicoEvent->SecondParticleCollection()->size() >= fMinSizePartCollection ))) {
353       fEventCut->FillCutMonitor(hbtEvent, tmpPassEvent);
354       
355
356 //------------------------------------------------------------------------------
357 //   Temporary comment:
358 //      This whole section rewritten so that all pairs are built using the
359 //      same code... easier to read and manage, and MakePairs() can be called by
360 //      derived classes.  Also, the requirement of a full mixing buffer before
361 //      mixing is removed.
362 //                          Dan Magestro, 11/2002
363
364       //------ Make real pairs. If identical, make pairs for one collection ------//
365
366       if (AnalyzeIdenticalParticles()) {
367         MakePairs("real", fPicoEvent->FirstParticleCollection() );
368       }
369       else {
370         MakePairs("real", fPicoEvent->FirstParticleCollection(),
371                           fPicoEvent->SecondParticleCollection() );
372       }
373       // cout << "AliFemtoSimpleAnalysis::ProcessEvent() - reals done ";
374
375       //---- Make pairs for mixed events, looping over events in mixingBuffer ----//
376
377       AliFemtoPicoEvent* storedEvent;
378       AliFemtoPicoEventIterator fPicoEventIter;
379       for (fPicoEventIter=MixingBuffer()->begin();fPicoEventIter!=MixingBuffer()->end();fPicoEventIter++){
380         storedEvent = *fPicoEventIter;
381         if (AnalyzeIdenticalParticles()) {
382           MakePairs("mixed",fPicoEvent->FirstParticleCollection(),
383                             storedEvent->FirstParticleCollection() );
384         }
385         else {
386           MakePairs("mixed",fPicoEvent->FirstParticleCollection(),
387                             storedEvent->SecondParticleCollection() );
388
389           MakePairs("mixed",storedEvent->FirstParticleCollection(),
390                             fPicoEvent->SecondParticleCollection() );
391         }
392       }
393       // cout << " - mixed done   " << endl;
394
395       //--------- If mixing buffer is full, delete oldest event ---------//
396
397       if ( MixingBufferFull() ) {
398         delete MixingBuffer()->back();
399         MixingBuffer()->pop_back();
400       }
401
402       //-------- Add current event (fPicoEvent) to mixing buffer --------//
403
404       MixingBuffer()->push_front(fPicoEvent);
405
406
407 // Temporary comment: End of rewritten section... Dan Magestro, 11/2002
408 //------------------------------------------------------------------------------
409
410
411     }  // if ParticleCollections are big enough (mal jun2002)
412     else{
413       fEventCut->FillCutMonitor(hbtEvent, !tmpPassEvent);
414       delete fPicoEvent;
415     }
416   }   // if currentEvent is accepted by currentAnalysis
417   EventEnd(hbtEvent);  // cleanup for EbyE 
418   //cout << "AliFemtoSimpleAnalysis::ProcessEvent() - return to caller ... " << endl;
419 }
420 //_________________________
421 void AliFemtoSimpleAnalysis::MakePairs(const char* typeIn, AliFemtoParticleCollection *partCollection1,
422                                        AliFemtoParticleCollection *partCollection2){
423 // Build pairs, check pair cuts, and call CFs' AddRealPair() or
424 // AddMixedPair() methods. If no second particle collection is
425 // specfied, make pairs within first particle collection.
426
427   string type = typeIn;
428
429   //  int swpart = ((long int) partCollection1) % 2;
430   int swpart = fNeventsProcessed % 2;
431
432   AliFemtoPair* tPair = new AliFemtoPair;
433
434   AliFemtoCorrFctnIterator tCorrFctnIter;
435
436   AliFemtoParticleIterator tPartIter1, tPartIter2;
437
438   AliFemtoParticleIterator tStartOuterLoop = partCollection1->begin();  // always
439   AliFemtoParticleIterator tEndOuterLoop   = partCollection1->end();    // will be one less if identical
440   AliFemtoParticleIterator tStartInnerLoop;
441   AliFemtoParticleIterator tEndInnerLoop;
442   if (partCollection2) {                        // Two collections:
443     tStartInnerLoop = partCollection2->begin();  //   Full inner & outer loops
444     tEndInnerLoop   = partCollection2->end();    //
445   }
446   else {                                        // One collection:
447     tEndOuterLoop--;                             //   Outer loop goes to next-to-last particle
448     tEndInnerLoop = partCollection1->end() ;     //   Inner loop goes to last particle
449   }
450   for (tPartIter1=tStartOuterLoop;tPartIter1!=tEndOuterLoop;tPartIter1++) {
451     if (!partCollection2){
452       tStartInnerLoop = tPartIter1;
453       tStartInnerLoop++;
454     }
455     tPair->SetTrack1(*tPartIter1);
456     for (tPartIter2 = tStartInnerLoop; tPartIter2!=tEndInnerLoop;tPartIter2++) {
457       tPair->SetTrack2(*tPartIter2);
458
459       // The following lines have to be uncommented if you want pairCutMonitors
460       // they are not in for speed reasons
461       // bool tmpPassPair = fPairCut->Pass(tPair);
462       // fPairCut->FillCutMonitor(tPair, tmpPassPair);
463       // if ( tmpPassPair )
464
465       //---- If pair passes cut, loop over CF's and add pair to real/mixed ----//
466
467       if (!partCollection2) {
468         if (swpart) {
469           tPair->SetTrack1(*tPartIter2);
470           tPair->SetTrack2(*tPartIter1);
471           swpart = 0;
472         }
473         else {
474           tPair->SetTrack1(*tPartIter1);
475           tPair->SetTrack2(*tPartIter2);
476           swpart = 1;
477         }
478       }
479
480       if (fPairCut->Pass(tPair)){
481         for (tCorrFctnIter=fCorrFctnCollection->begin();
482              tCorrFctnIter!=fCorrFctnCollection->end();tCorrFctnIter++){
483           AliFemtoCorrFctn* tCorrFctn = *tCorrFctnIter;
484           if(type == "real")
485             tCorrFctn->AddRealPair(tPair);
486           else if(type == "mixed")
487             tCorrFctn->AddMixedPair(tPair);
488           else
489             cout << "Problem with pair type, type = " << type.c_str() << endl;
490         }
491       }
492
493     }    // loop over second particle
494
495   }      // loop over first particle
496
497   delete tPair;
498
499 }
500 //_________________________
501 void AliFemtoSimpleAnalysis::EventBegin(const AliFemtoEvent* ev){
502   // Perform initialization operations at the beginning of the event processing
503   //cout << " AliFemtoSimpleAnalysis::EventBegin(const AliFemtoEvent* ev) " << endl;
504   fFirstParticleCut->EventBegin(ev);
505   fSecondParticleCut->EventBegin(ev);
506   fPairCut->EventBegin(ev);
507   for (AliFemtoCorrFctnIterator iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
508     (*iter)->EventBegin(ev);
509   }
510 }
511 //_________________________
512 void AliFemtoSimpleAnalysis::EventEnd(const AliFemtoEvent* ev){
513   // Fiinsh operations at the end of event processing
514   fFirstParticleCut->EventEnd(ev);
515   fSecondParticleCut->EventEnd(ev);
516   fPairCut->EventEnd(ev);
517   for (AliFemtoCorrFctnIterator iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
518     (*iter)->EventEnd(ev);
519   }
520 }
521 //_________________________
522 void AliFemtoSimpleAnalysis::Finish(){
523   // Perform finishing operations after all events are processed
524   AliFemtoCorrFctnIterator iter;
525   for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
526     (*iter)->Finish();
527   }
528 }
529 //_________________________
530 void AliFemtoSimpleAnalysis::AddEventProcessed() {
531   // Increase count of processed events
532   fNeventsProcessed++;
533 }
534 //_________________________
535 TList* AliFemtoSimpleAnalysis::ListSettings()
536 {
537   // Collect settings list
538   TList *tListSettings = new TList();
539
540   TList *p1Cut = fFirstParticleCut->ListSettings();
541
542   TListIter nextp1(p1Cut);
543   while (TObject *obj = nextp1.Next()) {
544     TString cuts(obj->GetName());
545     cuts.Prepend("AliFemtoSimpleAnalysis.");
546     tListSettings->Add(new TObjString(cuts.Data()));
547   }
548
549   if (fSecondParticleCut != fFirstParticleCut) {
550     TList *p2Cut = fSecondParticleCut->ListSettings();
551     
552     TIter nextp2(p2Cut);
553     while (TObject *obj = nextp2()) {
554       TString cuts(obj->GetName());
555       cuts.Prepend("AliFemtoSimpleAnalysis.");
556       tListSettings->Add(new TObjString(cuts.Data()));
557     }
558   }
559
560   TList *pairCut = fPairCut->ListSettings();
561
562   TIter nextpair(pairCut);
563   while (TObject *obj = nextpair()) {
564     TString cuts(obj->GetName());
565     cuts.Prepend("AliFemtoSimpleAnalysis.");
566     tListSettings->Add(new TObjString(cuts.Data()));
567   }
568
569   return tListSettings;
570   
571 }
572
573 //_________________________
574 TList* AliFemtoSimpleAnalysis::GetOutputList()
575 {
576   // Collect the list of output objects
577   // to be written 
578   TList *tOutputList = new TList();
579
580   TList *p1Cut = fFirstParticleCut->GetOutputList();
581
582   TListIter nextp1(p1Cut);
583   while (TObject *obj = nextp1.Next()) {
584     tOutputList->Add(obj);
585   }
586
587   if (fSecondParticleCut != fFirstParticleCut) {
588     TList *p2Cut = fSecondParticleCut->GetOutputList();
589     
590     TIter nextp2(p2Cut);
591     while (TObject *obj = nextp2()) {
592       tOutputList->Add(obj);
593     }
594   }
595
596   TList *pairCut = fPairCut->GetOutputList();
597
598   TIter nextpair(pairCut);
599   while (TObject *obj = nextpair()) {
600     tOutputList->Add(obj);
601   }
602
603   TList *eventCut = fEventCut->GetOutputList();
604
605   TIter nextevent(eventCut);
606   while (TObject *obj = nextevent()) {
607     tOutputList->Add(obj);
608   }
609
610   AliFemtoCorrFctnIterator iter;
611   for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
612     TList *tListCf = (*iter)->GetOutputList();
613     
614     TIter nextListCf(tListCf);
615     while (TObject *obj = nextListCf()) {
616       tOutputList->Add(obj);
617     }
618   }
619
620   return tOutputList;
621   
622 }