ROOT logo
///////////////////////////////////////////////////////////////////////////
//                                                                       //
// AliFemtoSimpleAnalysis - the most basic analysis there is. All other  //
// inherit from this one. Provides basic functionality for the analysis. //
// To properly set up the analysis the following steps should be taken:  //
//                                                                       //
// - create particle cuts and add them via SetFirstParticleCut and       //
//   SetSecondParticleCut. If one analyzes identical particle            //
//   correlations, the first particle cut must be also the second        //
//   particle cut.                                                       //
//                                                                       //
// - create pair cuts and add them via SetPairCut                        //
//                                                                       //
// - create one or many correlation functions and add them via           //
//   AddCorrFctn method.                                                 //
//                                                                       //
// - specify how many events are to be strored in the mixing buffer for  //
//   background construction                                             //
//                                                                       //
// Then, when the analysis is run, for each event, the EventBegin is     //
// called before any processing is done, then the ProcessEvent is called //
// which takes care of creating real and mixed pairs and sending them    //
// to all the registered correlation functions. At the end of each event,//
// after all pairs are processed, EventEnd is called. After the whole    //
// analysis finishes (there is no more events to process) Finish() is    //
// called.                                                               //
//                                                                       //
///////////////////////////////////////////////////////////////////////////

#include "AliFemtoSimpleAnalysis.h"
#include "AliFemtoTrackCut.h"
#include "AliFemtoV0Cut.h"
#include "AliFemtoKinkCut.h"
#include <string>
#include <iostream>

// blah blah

#ifdef __ROOT__
ClassImp(AliFemtoSimpleAnalysis)
#endif

AliFemtoEventCut*    copyTheCut(AliFemtoEventCut*);
AliFemtoParticleCut* copyTheCut(AliFemtoParticleCut*);
AliFemtoPairCut*     copyTheCut(AliFemtoPairCut*);
AliFemtoCorrFctn*    copyTheCorrFctn(AliFemtoCorrFctn*);

// this little function used to apply ParticleCuts (TrackCuts or V0Cuts) and fill ParticleCollections of picoEvent
//  it is called from AliFemtoSimpleAnalysis::ProcessEvent()
void FillHbtParticleCollection(AliFemtoParticleCut*         partCut,
			       AliFemtoEvent*               hbtEvent,
			       AliFemtoParticleCollection*  partCollection,
			       bool performSharedDaughterCut=kFALSE)
{
  // Fill particle collections from the event
  // by the particles that pass all the cuts
  switch (partCut->Type()) {
  case hbtTrack:       // cut is cutting on Tracks
  {
    AliFemtoTrackCut* pCut = (AliFemtoTrackCut*) partCut;
    AliFemtoTrack* pParticle;
    AliFemtoTrackIterator pIter;
    AliFemtoTrackIterator startLoop = hbtEvent->TrackCollection()->begin();
    AliFemtoTrackIterator endLoop   = hbtEvent->TrackCollection()->end();
    for (pIter=startLoop;pIter!=endLoop;pIter++){
      pParticle = *pIter;
      bool tmpPassParticle = pCut->Pass(pParticle);
      pCut->FillCutMonitor(pParticle, tmpPassParticle);
      if (tmpPassParticle){
        AliFemtoParticle* particle = new AliFemtoParticle(pParticle,pCut->Mass());
        partCollection->push_back(particle);
      }
    }
    break;
  }
  case hbtV0:          // cut is cutting on V0s
    {
      AliFemtoV0Cut* pCut = (AliFemtoV0Cut*) partCut;
      AliFemtoV0* pParticle;
      AliFemtoV0Iterator pIter;

      if(performSharedDaughterCut) {
        AliFemtoV0Collection V0CorrectedCollection;
        AliFemtoV0SharedDaughterCut sharedDaughterCut;
        V0CorrectedCollection = sharedDaughterCut.AliFemtoV0SharedDaughterCutCollection(hbtEvent->V0Collection(), pCut);

        AliFemtoV0Iterator startLoop = V0CorrectedCollection.begin();
        AliFemtoV0Iterator endLoop   = V0CorrectedCollection.end();

        for (pIter=startLoop;pIter!=endLoop;pIter++){
          pParticle = *pIter;
          AliFemtoParticle* particle = new AliFemtoParticle(pParticle,partCut->Mass());
          partCollection->push_back(particle);
        }
      }
      else { //previous, untouched loop:
        AliFemtoV0Iterator startLoop = hbtEvent->V0Collection()->begin();
        AliFemtoV0Iterator endLoop   = hbtEvent->V0Collection()->end();
        for (pIter=startLoop;pIter!=endLoop;pIter++){
          pParticle = *pIter;
          bool tmpPassV0 = pCut->Pass(pParticle);
          pCut->FillCutMonitor(pParticle,tmpPassV0);
          if(tmpPassV0) {
            AliFemtoParticle* particle = new AliFemtoParticle(pParticle,partCut->Mass());
            partCollection->push_back(particle);
          }
        }
      }

      pCut->FillCutMonitor(hbtEvent,partCollection);// Gael 19/06/02

      break;
    }
  case hbtKink:          // cut is cutting on Kinks  -- mal 25May2001
    {
      AliFemtoKinkCut* pCut = (AliFemtoKinkCut*) partCut;
      AliFemtoKink* pParticle;
      AliFemtoKinkIterator pIter;
      AliFemtoKinkIterator startLoop = hbtEvent->KinkCollection()->begin();
      AliFemtoKinkIterator endLoop   = hbtEvent->KinkCollection()->end();
      // this following "for" loop is identical to the one above, but because of scoping, I can's see how to avoid repitition...
      for (pIter=startLoop;pIter!=endLoop;pIter++){
	pParticle = *pIter;
	bool tmpPass = pCut->Pass(pParticle);
	pCut->FillCutMonitor(pParticle,tmpPass);
	if (tmpPass){
	  AliFemtoParticle* particle = new AliFemtoParticle(pParticle,partCut->Mass());
	  partCollection->push_back(particle);
	}
      }
      break;
    }
  default:
    cout << "FillHbtParticleCollection function (in AliFemtoSimpleAnalysis.cxx) - undefined Particle Cut type!!! \n";
  }
}
//____________________________
AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis() :
  fPicoEventCollectionVectorHideAway(0),
  fPairCut(0),
  fCorrFctnCollection(0),
  fEventCut(0),
  fFirstParticleCut(0),
  fSecondParticleCut(0),
  fMixingBuffer(0),
  fPicoEvent(0),
  fNumEventsToMix(0),
  fNeventsProcessed(0),
  fMinSizePartCollection(0),
  fVerbose(kTRUE),
  fPerformSharedDaughterCut(kFALSE)
{
  // Default constructor
  //  mControlSwitch     = 0;
  fCorrFctnCollection = new AliFemtoCorrFctnCollection;
  fMixingBuffer = new AliFemtoPicoEventCollection;
}
//____________________________

AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) :
  AliFemtoAnalysis(),
  fPicoEventCollectionVectorHideAway(0),
  fPairCut(0),
  fCorrFctnCollection(0),
  fEventCut(0),
  fFirstParticleCut(0),
  fSecondParticleCut(0),
  fMixingBuffer(0),
  fPicoEvent(0),
  fNumEventsToMix(0),
  fNeventsProcessed(0),
  fMinSizePartCollection(0),
  fVerbose(kTRUE),
  fPerformSharedDaughterCut(kFALSE)
{
  // Copy constructor
  //AliFemtoSimpleAnalysis();
  fCorrFctnCollection = new AliFemtoCorrFctnCollection;
  fMixingBuffer = new AliFemtoPicoEventCollection;

  // find the right event cut
  fEventCut = a.fEventCut->Clone();
  // find the right first particle cut
  fFirstParticleCut = a.fFirstParticleCut->Clone();
  // find the right second particle cut
  if (a.fFirstParticleCut==a.fSecondParticleCut)
    SetSecondParticleCut(fFirstParticleCut); // identical particle hbt
  else
  fSecondParticleCut = a.fSecondParticleCut->Clone();

  fPairCut = a.fPairCut->Clone();

  if ( fEventCut ) {
      SetEventCut(fEventCut); // this will set the myAnalysis pointer inside the cut
      cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - event cut set " << endl;
  }
  if ( fFirstParticleCut ) {
      SetFirstParticleCut(fFirstParticleCut); // this will set the myAnalysis pointer inside the cut
      cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - first particle cut set " << endl;
  }
  if ( fSecondParticleCut ) {
      SetSecondParticleCut(fSecondParticleCut); // this will set the myAnalysis pointer inside the cut
      cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - second particle cut set " << endl;
  }  if ( fPairCut ) {
      SetPairCut(fPairCut); // this will set the myAnalysis pointer inside the cut
      cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - pair cut set " << endl;
  }

  AliFemtoCorrFctnIterator iter;
  for (iter=a.fCorrFctnCollection->begin(); iter!=a.fCorrFctnCollection->end();iter++){
    cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - looking for correlation functions " << endl;
    AliFemtoCorrFctn* fctn = (*iter)->Clone();
    if (fctn) AddCorrFctn(fctn);
    else cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - correlation function not found " << endl;
  }

  fNumEventsToMix = a.fNumEventsToMix;

  fMinSizePartCollection = a.fMinSizePartCollection;  // minimum # particles in ParticleCollection

  cout << " AliFemtoSimpleAnalysis::AliFemtoSimpleAnalysis(const AliFemtoSimpleAnalysis& a) - analysis copied " << endl;

}
//____________________________
AliFemtoSimpleAnalysis::~AliFemtoSimpleAnalysis(){
  // destructor
  cout << " AliFemtoSimpleAnalysis::~AliFemtoSimpleAnalysis()" << endl;
  if (fEventCut) delete fEventCut; fEventCut=0;
  if (fFirstParticleCut == fSecondParticleCut) fSecondParticleCut=0;
  if (fFirstParticleCut)  delete fFirstParticleCut; fFirstParticleCut=0;
  if (fSecondParticleCut) delete fSecondParticleCut; fSecondParticleCut=0;
  if (fPairCut) delete fPairCut; fPairCut=0;
  // now delete every CorrFunction in the Collection, and then the Collection itself
  AliFemtoCorrFctnIterator iter;
  for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
    delete *iter;
  }
  delete fCorrFctnCollection;
  // now delete every PicoEvent in the EventMixingBuffer and then the Buffer itself
  if (fMixingBuffer) {
    AliFemtoPicoEventIterator piter;
    for (piter=fMixingBuffer->begin();piter!=fMixingBuffer->end();piter++){
      delete *piter;
    }
    delete fMixingBuffer;
  }
}
//______________________
AliFemtoSimpleAnalysis& AliFemtoSimpleAnalysis::operator=(const AliFemtoSimpleAnalysis& aAna)
{
  // Assignment operator
  if (this == &aAna)
    return *this;

  if (fCorrFctnCollection) delete fCorrFctnCollection;
  fCorrFctnCollection = new AliFemtoCorrFctnCollection;
  if (fMixingBuffer) delete fMixingBuffer;
  fMixingBuffer = new AliFemtoPicoEventCollection;

  // find the right event cut
  if (fEventCut) delete fEventCut;
  fEventCut = aAna.fEventCut->Clone();
  // find the right first particle cut
  if (fFirstParticleCut) delete fFirstParticleCut;
  fFirstParticleCut = aAna.fFirstParticleCut->Clone();
  // find the right second particle cut
  if (fSecondParticleCut) delete fSecondParticleCut;
  if (aAna.fFirstParticleCut==aAna.fSecondParticleCut)
    SetSecondParticleCut(fFirstParticleCut); // identical particle hbt
  else
    fSecondParticleCut = aAna.fSecondParticleCut->Clone();

  if (fPairCut) delete fPairCut;
  fPairCut = aAna.fPairCut->Clone();

  if ( fEventCut ) {
    SetEventCut(fEventCut); // this will set the myAnalysis pointer inside the cut
  }
  if ( fFirstParticleCut ) {
    SetFirstParticleCut(fFirstParticleCut); // this will set the myAnalysis pointer inside the cut
  }
  if ( fSecondParticleCut ) {
    SetSecondParticleCut(fSecondParticleCut); // this will set the myAnalysis pointer inside the cut
  }
  if ( fPairCut ) {
    SetPairCut(fPairCut); // this will set the myAnalysis pointer inside the cut
  }

  AliFemtoCorrFctnIterator iter;
  for (iter=aAna.fCorrFctnCollection->begin(); iter!=aAna.fCorrFctnCollection->end();iter++){
    AliFemtoCorrFctn* fctn = (*iter)->Clone();
    if (fctn) AddCorrFctn(fctn);
  }

  fNumEventsToMix = aAna.fNumEventsToMix;

  fMinSizePartCollection = aAna.fMinSizePartCollection;  // minimum # particles in ParticleCollection

  fVerbose = aAna.fVerbose;

  fPerformSharedDaughterCut = aAna.fPerformSharedDaughterCut;

  return *this;
}
//______________________
AliFemtoCorrFctn* AliFemtoSimpleAnalysis::CorrFctn(int n){
  // return pointer to n-th correlation function
  if ( n<0 || n > (int)fCorrFctnCollection->size() )
    return NULL;
  AliFemtoCorrFctnIterator iter=fCorrFctnCollection->begin();
  for (int i=0; i<n ;i++){
    iter++;
  }
  return *iter;
}
//____________________________
AliFemtoString AliFemtoSimpleAnalysis::Report()
{
  // Create a simple report from the analysis execution
  cout << "AliFemtoSimpleAnalysis - constructing Report..."<<endl;
  string temp = "-----------\nHbt Analysis Report:\n";
  temp += "\nEvent Cuts:\n";
  temp += fEventCut->Report();
  temp += "\nParticle Cuts - First Particle:\n";
  temp += fFirstParticleCut->Report();
  temp += "\nParticle Cuts - Second Particle:\n";
  temp += fSecondParticleCut->Report();
  temp += "\nPair Cuts:\n";
  temp += fPairCut->Report();
  temp += "\nCorrelation Functions:\n";
  AliFemtoCorrFctnIterator iter;
  if ( fCorrFctnCollection->size()==0 ) {
    cout << "AliFemtoSimpleAnalysis-Warning : no correlations functions in this analysis " << endl;
  }
  for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
    temp += (*iter)->Report();
    temp += "\n";
  }
  temp += "-------------\n";
  AliFemtoString returnThis=temp;
  return returnThis;
}
//_________________________
void AliFemtoSimpleAnalysis::ProcessEvent(const AliFemtoEvent* hbtEvent) {
  // Add event to processed events

  fPicoEvent=0; // we will get a new pico event, if not prevent corr. fctn to access old pico event
  AddEventProcessed();
  // startup for EbyE
  EventBegin(hbtEvent);
  // event cut and event cut monitor
  bool tmpPassEvent = fEventCut->Pass(hbtEvent);
  if (!tmpPassEvent)
    fEventCut->FillCutMonitor(hbtEvent, tmpPassEvent);
  if (tmpPassEvent) {
    //  cout << "AliFemtoSimpleAnalysis::ProcessEvent() - Event has passed cut - build picoEvent from " <<
    //   hbtEvent->TrackCollection()->size() << " tracks in TrackCollection" << endl;
    //  cout << "Event has passed cut with " << hbtEvent->TrackCollection()->size() << " tracks" << endl;
    // OK, analysis likes the event-- build a pico event from it, using tracks the analysis likes...
    fPicoEvent = new AliFemtoPicoEvent; // this is what we will make pairs from and put in Mixing Buffer
    // no memory leak. we will delete picoevents when they come out of the mixing buffer
    FillHbtParticleCollection(fFirstParticleCut,(AliFemtoEvent*)hbtEvent,fPicoEvent->FirstParticleCollection(), fPerformSharedDaughterCut);
    if ( !(AnalyzeIdenticalParticles()) )
      FillHbtParticleCollection(fSecondParticleCut,(AliFemtoEvent*)hbtEvent,fPicoEvent->SecondParticleCollection(),
				fPerformSharedDaughterCut);
    //cout <<"AliFemtoSimpleAnalysis::ProcessEvent - #particles in First, Second Collections: " <<
//       fPicoEvent->FirstParticleCollection()->size() << " " <<
//       fPicoEvent->SecondParticleCollection()->size() << endl;

    if (fVerbose)
    cout << "#particles in Collection 1, 2: " <<
       fPicoEvent->FirstParticleCollection()->size() << " " <<
       fPicoEvent->SecondParticleCollection()->size() << endl;
    fEventCut->FillCutMonitor(fPicoEvent->FirstParticleCollection(),fPicoEvent->SecondParticleCollection()); //MJ!


    // mal - implement a switch which allows only using events with ParticleCollections containing a minimum
    // number of entries (jun2002)
    if ((fPicoEvent->FirstParticleCollection()->size() >= fMinSizePartCollection )
	&& ( AnalyzeIdenticalParticles() || (fPicoEvent->SecondParticleCollection()->size() >= fMinSizePartCollection ))) {
      fEventCut->FillCutMonitor(hbtEvent, tmpPassEvent);


//------------------------------------------------------------------------------
//   Temporary comment:
//      This whole section rewritten so that all pairs are built using the
//      same code... easier to read and manage, and MakePairs() can be called by
//      derived classes.  Also, the requirement of a full mixing buffer before
//      mixing is removed.
//                          Dan Magestro, 11/2002

      //------ Make real pairs. If identical, make pairs for one collection ------//

      if (AnalyzeIdenticalParticles()) {
        MakePairs("real", fPicoEvent->FirstParticleCollection() );
      }
      else {
        MakePairs("real", fPicoEvent->FirstParticleCollection(),
                          fPicoEvent->SecondParticleCollection() );
      }

      if (fVerbose)
       cout << "AliFemtoSimpleAnalysis::ProcessEvent() - reals done ";

      //---- Make pairs for mixed events, looping over events in mixingBuffer ----//

      AliFemtoPicoEvent* storedEvent;
      AliFemtoPicoEventIterator fPicoEventIter;
      for (fPicoEventIter=MixingBuffer()->begin();fPicoEventIter!=MixingBuffer()->end();fPicoEventIter++){
        storedEvent = *fPicoEventIter;
        if (AnalyzeIdenticalParticles()) {
          MakePairs("mixed",fPicoEvent->FirstParticleCollection(),
                            storedEvent->FirstParticleCollection() );
        }
        else {
          MakePairs("mixed",fPicoEvent->FirstParticleCollection(),
                            storedEvent->SecondParticleCollection() );

          MakePairs("mixed",storedEvent->FirstParticleCollection(),
                            fPicoEvent->SecondParticleCollection() );
        }
      }

      if (fVerbose)
       cout << " - mixed done   " << endl;

      //--------- If mixing buffer is full, delete oldest event ---------//

      if ( MixingBufferFull() ) {
        delete MixingBuffer()->back();
        MixingBuffer()->pop_back();
      }

      //-------- Add current event (fPicoEvent) to mixing buffer --------//

      MixingBuffer()->push_front(fPicoEvent);


// Temporary comment: End of rewritten section... Dan Magestro, 11/2002
//------------------------------------------------------------------------------


    }  // if ParticleCollections are big enough (mal jun2002)
    else{
      fEventCut->FillCutMonitor(hbtEvent, !tmpPassEvent);
      delete fPicoEvent;
    }
  }   // if currentEvent is accepted by currentAnalysis
  EventEnd(hbtEvent);  // cleanup for EbyE
  //cout << "AliFemtoSimpleAnalysis::ProcessEvent() - return to caller ... " << endl;
}
//_________________________
void AliFemtoSimpleAnalysis::MakePairs(const char* typeIn, AliFemtoParticleCollection *partCollection1,
				       AliFemtoParticleCollection *partCollection2){
// Build pairs, check pair cuts, and call CFs' AddRealPair() or
// AddMixedPair() methods. If no second particle collection is
// specfied, make pairs within first particle collection.

  string type = typeIn;

  //  int swpart = ((long int) partCollection1) % 2;
  int swpart = fNeventsProcessed % 2;

  AliFemtoPair* tPair = new AliFemtoPair;

  AliFemtoCorrFctnIterator tCorrFctnIter;

  AliFemtoParticleIterator tPartIter1, tPartIter2;

  AliFemtoParticleIterator tStartOuterLoop = partCollection1->begin();  // always
  AliFemtoParticleIterator tEndOuterLoop   = partCollection1->end();    // will be one less if identical
  AliFemtoParticleIterator tStartInnerLoop;
  AliFemtoParticleIterator tEndInnerLoop;
  if (partCollection2) {                        // Two collections:
    tStartInnerLoop = partCollection2->begin();  //   Full inner & outer loops
    tEndInnerLoop   = partCollection2->end();    //
  }
  else {                                        // One collection:
    tEndOuterLoop--;                             //   Outer loop goes to next-to-last particle
    tEndInnerLoop = partCollection1->end() ;     //   Inner loop goes to last particle
  }
  for (tPartIter1=tStartOuterLoop;tPartIter1!=tEndOuterLoop;tPartIter1++) {
    if (!partCollection2){
      tStartInnerLoop = tPartIter1;
      tStartInnerLoop++;
    }
    tPair->SetTrack1(*tPartIter1);
    for (tPartIter2 = tStartInnerLoop; tPartIter2!=tEndInnerLoop;tPartIter2++) {
      tPair->SetTrack2(*tPartIter2);

      // The following lines have to be uncommented if you want pairCutMonitors
      // they are not in // for speed reasons
      //bool tmpPassPair = fPairCut->Pass(tPair);
      //fPairCut->FillCutMonitor(tPair, tmpPassPair);
      // // if ( tmpPassPair )

      //---- If pair passes cut, loop over CF's and add pair to real/mixed ----//

      if (!partCollection2) {
	if (swpart) {
 	  tPair->SetTrack1(*tPartIter2);
 	  tPair->SetTrack2(*tPartIter1);
 	  swpart = 0;
 	}
 	else {
 	  tPair->SetTrack1(*tPartIter1);
 	  tPair->SetTrack2(*tPartIter2);
 	  swpart = 1;
	}
      }

      if (fPairCut->Pass(tPair)){
        for (tCorrFctnIter=fCorrFctnCollection->begin();
             tCorrFctnIter!=fCorrFctnCollection->end();tCorrFctnIter++){
          AliFemtoCorrFctn* tCorrFctn = *tCorrFctnIter;
          if(type == "real")
            tCorrFctn->AddRealPair(tPair);
	  else if(type == "mixed")
            tCorrFctn->AddMixedPair(tPair);
          else
            cout << "Problem with pair type, type = " << type.c_str() << endl;
        }
      }

    }    // loop over second particle

  }      // loop over first particle

  delete tPair;

}
//_________________________
void AliFemtoSimpleAnalysis::EventBegin(const AliFemtoEvent* ev){
  // Perform initialization operations at the beginning of the event processing
  //cout << " AliFemtoSimpleAnalysis::EventBegin(const AliFemtoEvent* ev) " << endl;
  fFirstParticleCut->EventBegin(ev);
  fSecondParticleCut->EventBegin(ev);
  fPairCut->EventBegin(ev);
  for (AliFemtoCorrFctnIterator iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
    (*iter)->EventBegin(ev);
  }
}
//_________________________
void AliFemtoSimpleAnalysis::EventEnd(const AliFemtoEvent* ev){
  // Fiinsh operations at the end of event processing
  fFirstParticleCut->EventEnd(ev);
  fSecondParticleCut->EventEnd(ev);
  fPairCut->EventEnd(ev);
  for (AliFemtoCorrFctnIterator iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
    (*iter)->EventEnd(ev);
  }
}
//_________________________
void AliFemtoSimpleAnalysis::Finish(){
  // Perform finishing operations after all events are processed
  AliFemtoCorrFctnIterator iter;
  for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
    (*iter)->Finish();
  }
}
//_________________________
void AliFemtoSimpleAnalysis::AddEventProcessed() {
  // Increase count of processed events
  fNeventsProcessed++;
}
//_________________________
TList* AliFemtoSimpleAnalysis::ListSettings()
{
  // Collect settings list
  TList *tListSettings = new TList();

  TList *p1Cut = fFirstParticleCut->ListSettings();

  TListIter nextp1(p1Cut);
  while (TObject *obj = nextp1.Next()) {
    TString cuts(obj->GetName());
    cuts.Prepend("AliFemtoSimpleAnalysis.");
    tListSettings->Add(new TObjString(cuts.Data()));
  }

  if (fSecondParticleCut != fFirstParticleCut) {
    TList *p2Cut = fSecondParticleCut->ListSettings();

    TIter nextp2(p2Cut);
    while (TObject *obj = nextp2()) {
      TString cuts(obj->GetName());
      cuts.Prepend("AliFemtoSimpleAnalysis.");
      tListSettings->Add(new TObjString(cuts.Data()));
    }
  }

  TList *pairCut = fPairCut->ListSettings();

  TIter nextpair(pairCut);
  while (TObject *obj = nextpair()) {
    TString cuts(obj->GetName());
    cuts.Prepend("AliFemtoSimpleAnalysis.");
    tListSettings->Add(new TObjString(cuts.Data()));
  }

  return tListSettings;

}

//_________________________
TList* AliFemtoSimpleAnalysis::GetOutputList()
{
  // Collect the list of output objects
  // to be written
  TList *tOutputList = new TList();

  TList *p1Cut = fFirstParticleCut->GetOutputList();

  TListIter nextp1(p1Cut);
  while (TObject *obj = nextp1.Next()) {
    tOutputList->Add(obj);
  }
  delete p1Cut;

  if (fSecondParticleCut != fFirstParticleCut) {
    TList *p2Cut = fSecondParticleCut->GetOutputList();

    TIter nextp2(p2Cut);
    while (TObject *obj = nextp2()) {
      tOutputList->Add(obj);
    }
    delete p2Cut;
  }

  TList *pairCut = fPairCut->GetOutputList();

  TIter nextpair(pairCut);
  while (TObject *obj = nextpair()) {
    tOutputList->Add(obj);
  }
  delete pairCut;

  TList *eventCut = fEventCut->GetOutputList();

  TIter nextevent(eventCut);
  while (TObject *obj = nextevent()) {
    tOutputList->Add(obj);
  }
  delete eventCut;

  AliFemtoCorrFctnIterator iter;
  for (iter=fCorrFctnCollection->begin(); iter!=fCorrFctnCollection->end();iter++){
    TList *tListCf = (*iter)->GetOutputList();

    TIter nextListCf(tListCf);
    while (TObject *obj = nextListCf()) {
      tOutputList->Add(obj);
    }
    delete tListCf;
  }

  return tOutputList;

}
 AliFemtoSimpleAnalysis.cxx:1
 AliFemtoSimpleAnalysis.cxx:2
 AliFemtoSimpleAnalysis.cxx:3
 AliFemtoSimpleAnalysis.cxx:4
 AliFemtoSimpleAnalysis.cxx:5
 AliFemtoSimpleAnalysis.cxx:6
 AliFemtoSimpleAnalysis.cxx:7
 AliFemtoSimpleAnalysis.cxx:8
 AliFemtoSimpleAnalysis.cxx:9
 AliFemtoSimpleAnalysis.cxx:10
 AliFemtoSimpleAnalysis.cxx:11
 AliFemtoSimpleAnalysis.cxx:12
 AliFemtoSimpleAnalysis.cxx:13
 AliFemtoSimpleAnalysis.cxx:14
 AliFemtoSimpleAnalysis.cxx:15
 AliFemtoSimpleAnalysis.cxx:16
 AliFemtoSimpleAnalysis.cxx:17
 AliFemtoSimpleAnalysis.cxx:18
 AliFemtoSimpleAnalysis.cxx:19
 AliFemtoSimpleAnalysis.cxx:20
 AliFemtoSimpleAnalysis.cxx:21
 AliFemtoSimpleAnalysis.cxx:22
 AliFemtoSimpleAnalysis.cxx:23
 AliFemtoSimpleAnalysis.cxx:24
 AliFemtoSimpleAnalysis.cxx:25
 AliFemtoSimpleAnalysis.cxx:26
 AliFemtoSimpleAnalysis.cxx:27
 AliFemtoSimpleAnalysis.cxx:28
 AliFemtoSimpleAnalysis.cxx:29
 AliFemtoSimpleAnalysis.cxx:30
 AliFemtoSimpleAnalysis.cxx:31
 AliFemtoSimpleAnalysis.cxx:32
 AliFemtoSimpleAnalysis.cxx:33
 AliFemtoSimpleAnalysis.cxx:34
 AliFemtoSimpleAnalysis.cxx:35
 AliFemtoSimpleAnalysis.cxx:36
 AliFemtoSimpleAnalysis.cxx:37
 AliFemtoSimpleAnalysis.cxx:38
 AliFemtoSimpleAnalysis.cxx:39
 AliFemtoSimpleAnalysis.cxx:40
 AliFemtoSimpleAnalysis.cxx:41
 AliFemtoSimpleAnalysis.cxx:42
 AliFemtoSimpleAnalysis.cxx:43
 AliFemtoSimpleAnalysis.cxx:44
 AliFemtoSimpleAnalysis.cxx:45
 AliFemtoSimpleAnalysis.cxx:46
 AliFemtoSimpleAnalysis.cxx:47
 AliFemtoSimpleAnalysis.cxx:48
 AliFemtoSimpleAnalysis.cxx:49
 AliFemtoSimpleAnalysis.cxx:50
 AliFemtoSimpleAnalysis.cxx:51
 AliFemtoSimpleAnalysis.cxx:52
 AliFemtoSimpleAnalysis.cxx:53
 AliFemtoSimpleAnalysis.cxx:54
 AliFemtoSimpleAnalysis.cxx:55
 AliFemtoSimpleAnalysis.cxx:56
 AliFemtoSimpleAnalysis.cxx:57
 AliFemtoSimpleAnalysis.cxx:58
 AliFemtoSimpleAnalysis.cxx:59
 AliFemtoSimpleAnalysis.cxx:60
 AliFemtoSimpleAnalysis.cxx:61
 AliFemtoSimpleAnalysis.cxx:62
 AliFemtoSimpleAnalysis.cxx:63
 AliFemtoSimpleAnalysis.cxx:64
 AliFemtoSimpleAnalysis.cxx:65
 AliFemtoSimpleAnalysis.cxx:66
 AliFemtoSimpleAnalysis.cxx:67
 AliFemtoSimpleAnalysis.cxx:68
 AliFemtoSimpleAnalysis.cxx:69
 AliFemtoSimpleAnalysis.cxx:70
 AliFemtoSimpleAnalysis.cxx:71
 AliFemtoSimpleAnalysis.cxx:72
 AliFemtoSimpleAnalysis.cxx:73
 AliFemtoSimpleAnalysis.cxx:74
 AliFemtoSimpleAnalysis.cxx:75
 AliFemtoSimpleAnalysis.cxx:76
 AliFemtoSimpleAnalysis.cxx:77
 AliFemtoSimpleAnalysis.cxx:78
 AliFemtoSimpleAnalysis.cxx:79
 AliFemtoSimpleAnalysis.cxx:80
 AliFemtoSimpleAnalysis.cxx:81
 AliFemtoSimpleAnalysis.cxx:82
 AliFemtoSimpleAnalysis.cxx:83
 AliFemtoSimpleAnalysis.cxx:84
 AliFemtoSimpleAnalysis.cxx:85
 AliFemtoSimpleAnalysis.cxx:86
 AliFemtoSimpleAnalysis.cxx:87
 AliFemtoSimpleAnalysis.cxx:88
 AliFemtoSimpleAnalysis.cxx:89
 AliFemtoSimpleAnalysis.cxx:90
 AliFemtoSimpleAnalysis.cxx:91
 AliFemtoSimpleAnalysis.cxx:92
 AliFemtoSimpleAnalysis.cxx:93
 AliFemtoSimpleAnalysis.cxx:94
 AliFemtoSimpleAnalysis.cxx:95
 AliFemtoSimpleAnalysis.cxx:96
 AliFemtoSimpleAnalysis.cxx:97
 AliFemtoSimpleAnalysis.cxx:98
 AliFemtoSimpleAnalysis.cxx:99
 AliFemtoSimpleAnalysis.cxx:100
 AliFemtoSimpleAnalysis.cxx:101
 AliFemtoSimpleAnalysis.cxx:102
 AliFemtoSimpleAnalysis.cxx:103
 AliFemtoSimpleAnalysis.cxx:104
 AliFemtoSimpleAnalysis.cxx:105
 AliFemtoSimpleAnalysis.cxx:106
 AliFemtoSimpleAnalysis.cxx:107
 AliFemtoSimpleAnalysis.cxx:108
 AliFemtoSimpleAnalysis.cxx:109
 AliFemtoSimpleAnalysis.cxx:110
 AliFemtoSimpleAnalysis.cxx:111
 AliFemtoSimpleAnalysis.cxx:112
 AliFemtoSimpleAnalysis.cxx:113
 AliFemtoSimpleAnalysis.cxx:114
 AliFemtoSimpleAnalysis.cxx:115
 AliFemtoSimpleAnalysis.cxx:116
 AliFemtoSimpleAnalysis.cxx:117
 AliFemtoSimpleAnalysis.cxx:118
 AliFemtoSimpleAnalysis.cxx:119
 AliFemtoSimpleAnalysis.cxx:120
 AliFemtoSimpleAnalysis.cxx:121
 AliFemtoSimpleAnalysis.cxx:122
 AliFemtoSimpleAnalysis.cxx:123
 AliFemtoSimpleAnalysis.cxx:124
 AliFemtoSimpleAnalysis.cxx:125
 AliFemtoSimpleAnalysis.cxx:126
 AliFemtoSimpleAnalysis.cxx:127
 AliFemtoSimpleAnalysis.cxx:128
 AliFemtoSimpleAnalysis.cxx:129
 AliFemtoSimpleAnalysis.cxx:130
 AliFemtoSimpleAnalysis.cxx:131
 AliFemtoSimpleAnalysis.cxx:132
 AliFemtoSimpleAnalysis.cxx:133
 AliFemtoSimpleAnalysis.cxx:134
 AliFemtoSimpleAnalysis.cxx:135
 AliFemtoSimpleAnalysis.cxx:136
 AliFemtoSimpleAnalysis.cxx:137
 AliFemtoSimpleAnalysis.cxx:138
 AliFemtoSimpleAnalysis.cxx:139
 AliFemtoSimpleAnalysis.cxx:140
 AliFemtoSimpleAnalysis.cxx:141
 AliFemtoSimpleAnalysis.cxx:142
 AliFemtoSimpleAnalysis.cxx:143
 AliFemtoSimpleAnalysis.cxx:144
 AliFemtoSimpleAnalysis.cxx:145
 AliFemtoSimpleAnalysis.cxx:146
 AliFemtoSimpleAnalysis.cxx:147
 AliFemtoSimpleAnalysis.cxx:148
 AliFemtoSimpleAnalysis.cxx:149
 AliFemtoSimpleAnalysis.cxx:150
 AliFemtoSimpleAnalysis.cxx:151
 AliFemtoSimpleAnalysis.cxx:152
 AliFemtoSimpleAnalysis.cxx:153
 AliFemtoSimpleAnalysis.cxx:154
 AliFemtoSimpleAnalysis.cxx:155
 AliFemtoSimpleAnalysis.cxx:156
 AliFemtoSimpleAnalysis.cxx:157
 AliFemtoSimpleAnalysis.cxx:158
 AliFemtoSimpleAnalysis.cxx:159
 AliFemtoSimpleAnalysis.cxx:160
 AliFemtoSimpleAnalysis.cxx:161
 AliFemtoSimpleAnalysis.cxx:162
 AliFemtoSimpleAnalysis.cxx:163
 AliFemtoSimpleAnalysis.cxx:164
 AliFemtoSimpleAnalysis.cxx:165
 AliFemtoSimpleAnalysis.cxx:166
 AliFemtoSimpleAnalysis.cxx:167
 AliFemtoSimpleAnalysis.cxx:168
 AliFemtoSimpleAnalysis.cxx:169
 AliFemtoSimpleAnalysis.cxx:170
 AliFemtoSimpleAnalysis.cxx:171
 AliFemtoSimpleAnalysis.cxx:172
 AliFemtoSimpleAnalysis.cxx:173
 AliFemtoSimpleAnalysis.cxx:174
 AliFemtoSimpleAnalysis.cxx:175
 AliFemtoSimpleAnalysis.cxx:176
 AliFemtoSimpleAnalysis.cxx:177
 AliFemtoSimpleAnalysis.cxx:178
 AliFemtoSimpleAnalysis.cxx:179
 AliFemtoSimpleAnalysis.cxx:180
 AliFemtoSimpleAnalysis.cxx:181
 AliFemtoSimpleAnalysis.cxx:182
 AliFemtoSimpleAnalysis.cxx:183
 AliFemtoSimpleAnalysis.cxx:184
 AliFemtoSimpleAnalysis.cxx:185
 AliFemtoSimpleAnalysis.cxx:186
 AliFemtoSimpleAnalysis.cxx:187
 AliFemtoSimpleAnalysis.cxx:188
 AliFemtoSimpleAnalysis.cxx:189
 AliFemtoSimpleAnalysis.cxx:190
 AliFemtoSimpleAnalysis.cxx:191
 AliFemtoSimpleAnalysis.cxx:192
 AliFemtoSimpleAnalysis.cxx:193
 AliFemtoSimpleAnalysis.cxx:194
 AliFemtoSimpleAnalysis.cxx:195
 AliFemtoSimpleAnalysis.cxx:196
 AliFemtoSimpleAnalysis.cxx:197
 AliFemtoSimpleAnalysis.cxx:198
 AliFemtoSimpleAnalysis.cxx:199
 AliFemtoSimpleAnalysis.cxx:200
 AliFemtoSimpleAnalysis.cxx:201
 AliFemtoSimpleAnalysis.cxx:202
 AliFemtoSimpleAnalysis.cxx:203
 AliFemtoSimpleAnalysis.cxx:204
 AliFemtoSimpleAnalysis.cxx:205
 AliFemtoSimpleAnalysis.cxx:206
 AliFemtoSimpleAnalysis.cxx:207
 AliFemtoSimpleAnalysis.cxx:208
 AliFemtoSimpleAnalysis.cxx:209
 AliFemtoSimpleAnalysis.cxx:210
 AliFemtoSimpleAnalysis.cxx:211
 AliFemtoSimpleAnalysis.cxx:212
 AliFemtoSimpleAnalysis.cxx:213
 AliFemtoSimpleAnalysis.cxx:214
 AliFemtoSimpleAnalysis.cxx:215
 AliFemtoSimpleAnalysis.cxx:216
 AliFemtoSimpleAnalysis.cxx:217
 AliFemtoSimpleAnalysis.cxx:218
 AliFemtoSimpleAnalysis.cxx:219
 AliFemtoSimpleAnalysis.cxx:220
 AliFemtoSimpleAnalysis.cxx:221
 AliFemtoSimpleAnalysis.cxx:222
 AliFemtoSimpleAnalysis.cxx:223
 AliFemtoSimpleAnalysis.cxx:224
 AliFemtoSimpleAnalysis.cxx:225
 AliFemtoSimpleAnalysis.cxx:226
 AliFemtoSimpleAnalysis.cxx:227
 AliFemtoSimpleAnalysis.cxx:228
 AliFemtoSimpleAnalysis.cxx:229
 AliFemtoSimpleAnalysis.cxx:230
 AliFemtoSimpleAnalysis.cxx:231
 AliFemtoSimpleAnalysis.cxx:232
 AliFemtoSimpleAnalysis.cxx:233
 AliFemtoSimpleAnalysis.cxx:234
 AliFemtoSimpleAnalysis.cxx:235
 AliFemtoSimpleAnalysis.cxx:236
 AliFemtoSimpleAnalysis.cxx:237
 AliFemtoSimpleAnalysis.cxx:238
 AliFemtoSimpleAnalysis.cxx:239
 AliFemtoSimpleAnalysis.cxx:240
 AliFemtoSimpleAnalysis.cxx:241
 AliFemtoSimpleAnalysis.cxx:242
 AliFemtoSimpleAnalysis.cxx:243
 AliFemtoSimpleAnalysis.cxx:244
 AliFemtoSimpleAnalysis.cxx:245
 AliFemtoSimpleAnalysis.cxx:246
 AliFemtoSimpleAnalysis.cxx:247
 AliFemtoSimpleAnalysis.cxx:248
 AliFemtoSimpleAnalysis.cxx:249
 AliFemtoSimpleAnalysis.cxx:250
 AliFemtoSimpleAnalysis.cxx:251
 AliFemtoSimpleAnalysis.cxx:252
 AliFemtoSimpleAnalysis.cxx:253
 AliFemtoSimpleAnalysis.cxx:254
 AliFemtoSimpleAnalysis.cxx:255
 AliFemtoSimpleAnalysis.cxx:256
 AliFemtoSimpleAnalysis.cxx:257
 AliFemtoSimpleAnalysis.cxx:258
 AliFemtoSimpleAnalysis.cxx:259
 AliFemtoSimpleAnalysis.cxx:260
 AliFemtoSimpleAnalysis.cxx:261
 AliFemtoSimpleAnalysis.cxx:262
 AliFemtoSimpleAnalysis.cxx:263
 AliFemtoSimpleAnalysis.cxx:264
 AliFemtoSimpleAnalysis.cxx:265
 AliFemtoSimpleAnalysis.cxx:266
 AliFemtoSimpleAnalysis.cxx:267
 AliFemtoSimpleAnalysis.cxx:268
 AliFemtoSimpleAnalysis.cxx:269
 AliFemtoSimpleAnalysis.cxx:270
 AliFemtoSimpleAnalysis.cxx:271
 AliFemtoSimpleAnalysis.cxx:272
 AliFemtoSimpleAnalysis.cxx:273
 AliFemtoSimpleAnalysis.cxx:274
 AliFemtoSimpleAnalysis.cxx:275
 AliFemtoSimpleAnalysis.cxx:276
 AliFemtoSimpleAnalysis.cxx:277
 AliFemtoSimpleAnalysis.cxx:278
 AliFemtoSimpleAnalysis.cxx:279
 AliFemtoSimpleAnalysis.cxx:280
 AliFemtoSimpleAnalysis.cxx:281
 AliFemtoSimpleAnalysis.cxx:282
 AliFemtoSimpleAnalysis.cxx:283
 AliFemtoSimpleAnalysis.cxx:284
 AliFemtoSimpleAnalysis.cxx:285
 AliFemtoSimpleAnalysis.cxx:286
 AliFemtoSimpleAnalysis.cxx:287
 AliFemtoSimpleAnalysis.cxx:288
 AliFemtoSimpleAnalysis.cxx:289
 AliFemtoSimpleAnalysis.cxx:290
 AliFemtoSimpleAnalysis.cxx:291
 AliFemtoSimpleAnalysis.cxx:292
 AliFemtoSimpleAnalysis.cxx:293
 AliFemtoSimpleAnalysis.cxx:294
 AliFemtoSimpleAnalysis.cxx:295
 AliFemtoSimpleAnalysis.cxx:296
 AliFemtoSimpleAnalysis.cxx:297
 AliFemtoSimpleAnalysis.cxx:298
 AliFemtoSimpleAnalysis.cxx:299
 AliFemtoSimpleAnalysis.cxx:300
 AliFemtoSimpleAnalysis.cxx:301
 AliFemtoSimpleAnalysis.cxx:302
 AliFemtoSimpleAnalysis.cxx:303
 AliFemtoSimpleAnalysis.cxx:304
 AliFemtoSimpleAnalysis.cxx:305
 AliFemtoSimpleAnalysis.cxx:306
 AliFemtoSimpleAnalysis.cxx:307
 AliFemtoSimpleAnalysis.cxx:308
 AliFemtoSimpleAnalysis.cxx:309
 AliFemtoSimpleAnalysis.cxx:310
 AliFemtoSimpleAnalysis.cxx:311
 AliFemtoSimpleAnalysis.cxx:312
 AliFemtoSimpleAnalysis.cxx:313
 AliFemtoSimpleAnalysis.cxx:314
 AliFemtoSimpleAnalysis.cxx:315
 AliFemtoSimpleAnalysis.cxx:316
 AliFemtoSimpleAnalysis.cxx:317
 AliFemtoSimpleAnalysis.cxx:318
 AliFemtoSimpleAnalysis.cxx:319
 AliFemtoSimpleAnalysis.cxx:320
 AliFemtoSimpleAnalysis.cxx:321
 AliFemtoSimpleAnalysis.cxx:322
 AliFemtoSimpleAnalysis.cxx:323
 AliFemtoSimpleAnalysis.cxx:324
 AliFemtoSimpleAnalysis.cxx:325
 AliFemtoSimpleAnalysis.cxx:326
 AliFemtoSimpleAnalysis.cxx:327
 AliFemtoSimpleAnalysis.cxx:328
 AliFemtoSimpleAnalysis.cxx:329
 AliFemtoSimpleAnalysis.cxx:330
 AliFemtoSimpleAnalysis.cxx:331
 AliFemtoSimpleAnalysis.cxx:332
 AliFemtoSimpleAnalysis.cxx:333
 AliFemtoSimpleAnalysis.cxx:334
 AliFemtoSimpleAnalysis.cxx:335
 AliFemtoSimpleAnalysis.cxx:336
 AliFemtoSimpleAnalysis.cxx:337
 AliFemtoSimpleAnalysis.cxx:338
 AliFemtoSimpleAnalysis.cxx:339
 AliFemtoSimpleAnalysis.cxx:340
 AliFemtoSimpleAnalysis.cxx:341
 AliFemtoSimpleAnalysis.cxx:342
 AliFemtoSimpleAnalysis.cxx:343
 AliFemtoSimpleAnalysis.cxx:344
 AliFemtoSimpleAnalysis.cxx:345
 AliFemtoSimpleAnalysis.cxx:346
 AliFemtoSimpleAnalysis.cxx:347
 AliFemtoSimpleAnalysis.cxx:348
 AliFemtoSimpleAnalysis.cxx:349
 AliFemtoSimpleAnalysis.cxx:350
 AliFemtoSimpleAnalysis.cxx:351
 AliFemtoSimpleAnalysis.cxx:352
 AliFemtoSimpleAnalysis.cxx:353
 AliFemtoSimpleAnalysis.cxx:354
 AliFemtoSimpleAnalysis.cxx:355
 AliFemtoSimpleAnalysis.cxx:356
 AliFemtoSimpleAnalysis.cxx:357
 AliFemtoSimpleAnalysis.cxx:358
 AliFemtoSimpleAnalysis.cxx:359
 AliFemtoSimpleAnalysis.cxx:360
 AliFemtoSimpleAnalysis.cxx:361
 AliFemtoSimpleAnalysis.cxx:362
 AliFemtoSimpleAnalysis.cxx:363
 AliFemtoSimpleAnalysis.cxx:364
 AliFemtoSimpleAnalysis.cxx:365
 AliFemtoSimpleAnalysis.cxx:366
 AliFemtoSimpleAnalysis.cxx:367
 AliFemtoSimpleAnalysis.cxx:368
 AliFemtoSimpleAnalysis.cxx:369
 AliFemtoSimpleAnalysis.cxx:370
 AliFemtoSimpleAnalysis.cxx:371
 AliFemtoSimpleAnalysis.cxx:372
 AliFemtoSimpleAnalysis.cxx:373
 AliFemtoSimpleAnalysis.cxx:374
 AliFemtoSimpleAnalysis.cxx:375
 AliFemtoSimpleAnalysis.cxx:376
 AliFemtoSimpleAnalysis.cxx:377
 AliFemtoSimpleAnalysis.cxx:378
 AliFemtoSimpleAnalysis.cxx:379
 AliFemtoSimpleAnalysis.cxx:380
 AliFemtoSimpleAnalysis.cxx:381
 AliFemtoSimpleAnalysis.cxx:382
 AliFemtoSimpleAnalysis.cxx:383
 AliFemtoSimpleAnalysis.cxx:384
 AliFemtoSimpleAnalysis.cxx:385
 AliFemtoSimpleAnalysis.cxx:386
 AliFemtoSimpleAnalysis.cxx:387
 AliFemtoSimpleAnalysis.cxx:388
 AliFemtoSimpleAnalysis.cxx:389
 AliFemtoSimpleAnalysis.cxx:390
 AliFemtoSimpleAnalysis.cxx:391
 AliFemtoSimpleAnalysis.cxx:392
 AliFemtoSimpleAnalysis.cxx:393
 AliFemtoSimpleAnalysis.cxx:394
 AliFemtoSimpleAnalysis.cxx:395
 AliFemtoSimpleAnalysis.cxx:396
 AliFemtoSimpleAnalysis.cxx:397
 AliFemtoSimpleAnalysis.cxx:398
 AliFemtoSimpleAnalysis.cxx:399
 AliFemtoSimpleAnalysis.cxx:400
 AliFemtoSimpleAnalysis.cxx:401
 AliFemtoSimpleAnalysis.cxx:402
 AliFemtoSimpleAnalysis.cxx:403
 AliFemtoSimpleAnalysis.cxx:404
 AliFemtoSimpleAnalysis.cxx:405
 AliFemtoSimpleAnalysis.cxx:406
 AliFemtoSimpleAnalysis.cxx:407
 AliFemtoSimpleAnalysis.cxx:408
 AliFemtoSimpleAnalysis.cxx:409
 AliFemtoSimpleAnalysis.cxx:410
 AliFemtoSimpleAnalysis.cxx:411
 AliFemtoSimpleAnalysis.cxx:412
 AliFemtoSimpleAnalysis.cxx:413
 AliFemtoSimpleAnalysis.cxx:414
 AliFemtoSimpleAnalysis.cxx:415
 AliFemtoSimpleAnalysis.cxx:416
 AliFemtoSimpleAnalysis.cxx:417
 AliFemtoSimpleAnalysis.cxx:418
 AliFemtoSimpleAnalysis.cxx:419
 AliFemtoSimpleAnalysis.cxx:420
 AliFemtoSimpleAnalysis.cxx:421
 AliFemtoSimpleAnalysis.cxx:422
 AliFemtoSimpleAnalysis.cxx:423
 AliFemtoSimpleAnalysis.cxx:424
 AliFemtoSimpleAnalysis.cxx:425
 AliFemtoSimpleAnalysis.cxx:426
 AliFemtoSimpleAnalysis.cxx:427
 AliFemtoSimpleAnalysis.cxx:428
 AliFemtoSimpleAnalysis.cxx:429
 AliFemtoSimpleAnalysis.cxx:430
 AliFemtoSimpleAnalysis.cxx:431
 AliFemtoSimpleAnalysis.cxx:432
 AliFemtoSimpleAnalysis.cxx:433
 AliFemtoSimpleAnalysis.cxx:434
 AliFemtoSimpleAnalysis.cxx:435
 AliFemtoSimpleAnalysis.cxx:436
 AliFemtoSimpleAnalysis.cxx:437
 AliFemtoSimpleAnalysis.cxx:438
 AliFemtoSimpleAnalysis.cxx:439
 AliFemtoSimpleAnalysis.cxx:440
 AliFemtoSimpleAnalysis.cxx:441
 AliFemtoSimpleAnalysis.cxx:442
 AliFemtoSimpleAnalysis.cxx:443
 AliFemtoSimpleAnalysis.cxx:444
 AliFemtoSimpleAnalysis.cxx:445
 AliFemtoSimpleAnalysis.cxx:446
 AliFemtoSimpleAnalysis.cxx:447
 AliFemtoSimpleAnalysis.cxx:448
 AliFemtoSimpleAnalysis.cxx:449
 AliFemtoSimpleAnalysis.cxx:450
 AliFemtoSimpleAnalysis.cxx:451
 AliFemtoSimpleAnalysis.cxx:452
 AliFemtoSimpleAnalysis.cxx:453
 AliFemtoSimpleAnalysis.cxx:454
 AliFemtoSimpleAnalysis.cxx:455
 AliFemtoSimpleAnalysis.cxx:456
 AliFemtoSimpleAnalysis.cxx:457
 AliFemtoSimpleAnalysis.cxx:458
 AliFemtoSimpleAnalysis.cxx:459
 AliFemtoSimpleAnalysis.cxx:460
 AliFemtoSimpleAnalysis.cxx:461
 AliFemtoSimpleAnalysis.cxx:462
 AliFemtoSimpleAnalysis.cxx:463
 AliFemtoSimpleAnalysis.cxx:464
 AliFemtoSimpleAnalysis.cxx:465
 AliFemtoSimpleAnalysis.cxx:466
 AliFemtoSimpleAnalysis.cxx:467
 AliFemtoSimpleAnalysis.cxx:468
 AliFemtoSimpleAnalysis.cxx:469
 AliFemtoSimpleAnalysis.cxx:470
 AliFemtoSimpleAnalysis.cxx:471
 AliFemtoSimpleAnalysis.cxx:472
 AliFemtoSimpleAnalysis.cxx:473
 AliFemtoSimpleAnalysis.cxx:474
 AliFemtoSimpleAnalysis.cxx:475
 AliFemtoSimpleAnalysis.cxx:476
 AliFemtoSimpleAnalysis.cxx:477
 AliFemtoSimpleAnalysis.cxx:478
 AliFemtoSimpleAnalysis.cxx:479
 AliFemtoSimpleAnalysis.cxx:480
 AliFemtoSimpleAnalysis.cxx:481
 AliFemtoSimpleAnalysis.cxx:482
 AliFemtoSimpleAnalysis.cxx:483
 AliFemtoSimpleAnalysis.cxx:484
 AliFemtoSimpleAnalysis.cxx:485
 AliFemtoSimpleAnalysis.cxx:486
 AliFemtoSimpleAnalysis.cxx:487
 AliFemtoSimpleAnalysis.cxx:488
 AliFemtoSimpleAnalysis.cxx:489
 AliFemtoSimpleAnalysis.cxx:490
 AliFemtoSimpleAnalysis.cxx:491
 AliFemtoSimpleAnalysis.cxx:492
 AliFemtoSimpleAnalysis.cxx:493
 AliFemtoSimpleAnalysis.cxx:494
 AliFemtoSimpleAnalysis.cxx:495
 AliFemtoSimpleAnalysis.cxx:496
 AliFemtoSimpleAnalysis.cxx:497
 AliFemtoSimpleAnalysis.cxx:498
 AliFemtoSimpleAnalysis.cxx:499
 AliFemtoSimpleAnalysis.cxx:500
 AliFemtoSimpleAnalysis.cxx:501
 AliFemtoSimpleAnalysis.cxx:502
 AliFemtoSimpleAnalysis.cxx:503
 AliFemtoSimpleAnalysis.cxx:504
 AliFemtoSimpleAnalysis.cxx:505
 AliFemtoSimpleAnalysis.cxx:506
 AliFemtoSimpleAnalysis.cxx:507
 AliFemtoSimpleAnalysis.cxx:508
 AliFemtoSimpleAnalysis.cxx:509
 AliFemtoSimpleAnalysis.cxx:510
 AliFemtoSimpleAnalysis.cxx:511
 AliFemtoSimpleAnalysis.cxx:512
 AliFemtoSimpleAnalysis.cxx:513
 AliFemtoSimpleAnalysis.cxx:514
 AliFemtoSimpleAnalysis.cxx:515
 AliFemtoSimpleAnalysis.cxx:516
 AliFemtoSimpleAnalysis.cxx:517
 AliFemtoSimpleAnalysis.cxx:518
 AliFemtoSimpleAnalysis.cxx:519
 AliFemtoSimpleAnalysis.cxx:520
 AliFemtoSimpleAnalysis.cxx:521
 AliFemtoSimpleAnalysis.cxx:522
 AliFemtoSimpleAnalysis.cxx:523
 AliFemtoSimpleAnalysis.cxx:524
 AliFemtoSimpleAnalysis.cxx:525
 AliFemtoSimpleAnalysis.cxx:526
 AliFemtoSimpleAnalysis.cxx:527
 AliFemtoSimpleAnalysis.cxx:528
 AliFemtoSimpleAnalysis.cxx:529
 AliFemtoSimpleAnalysis.cxx:530
 AliFemtoSimpleAnalysis.cxx:531
 AliFemtoSimpleAnalysis.cxx:532
 AliFemtoSimpleAnalysis.cxx:533
 AliFemtoSimpleAnalysis.cxx:534
 AliFemtoSimpleAnalysis.cxx:535
 AliFemtoSimpleAnalysis.cxx:536
 AliFemtoSimpleAnalysis.cxx:537
 AliFemtoSimpleAnalysis.cxx:538
 AliFemtoSimpleAnalysis.cxx:539
 AliFemtoSimpleAnalysis.cxx:540
 AliFemtoSimpleAnalysis.cxx:541
 AliFemtoSimpleAnalysis.cxx:542
 AliFemtoSimpleAnalysis.cxx:543
 AliFemtoSimpleAnalysis.cxx:544
 AliFemtoSimpleAnalysis.cxx:545
 AliFemtoSimpleAnalysis.cxx:546
 AliFemtoSimpleAnalysis.cxx:547
 AliFemtoSimpleAnalysis.cxx:548
 AliFemtoSimpleAnalysis.cxx:549
 AliFemtoSimpleAnalysis.cxx:550
 AliFemtoSimpleAnalysis.cxx:551
 AliFemtoSimpleAnalysis.cxx:552
 AliFemtoSimpleAnalysis.cxx:553
 AliFemtoSimpleAnalysis.cxx:554
 AliFemtoSimpleAnalysis.cxx:555
 AliFemtoSimpleAnalysis.cxx:556
 AliFemtoSimpleAnalysis.cxx:557
 AliFemtoSimpleAnalysis.cxx:558
 AliFemtoSimpleAnalysis.cxx:559
 AliFemtoSimpleAnalysis.cxx:560
 AliFemtoSimpleAnalysis.cxx:561
 AliFemtoSimpleAnalysis.cxx:562
 AliFemtoSimpleAnalysis.cxx:563
 AliFemtoSimpleAnalysis.cxx:564
 AliFemtoSimpleAnalysis.cxx:565
 AliFemtoSimpleAnalysis.cxx:566
 AliFemtoSimpleAnalysis.cxx:567
 AliFemtoSimpleAnalysis.cxx:568
 AliFemtoSimpleAnalysis.cxx:569
 AliFemtoSimpleAnalysis.cxx:570
 AliFemtoSimpleAnalysis.cxx:571
 AliFemtoSimpleAnalysis.cxx:572
 AliFemtoSimpleAnalysis.cxx:573
 AliFemtoSimpleAnalysis.cxx:574
 AliFemtoSimpleAnalysis.cxx:575
 AliFemtoSimpleAnalysis.cxx:576
 AliFemtoSimpleAnalysis.cxx:577
 AliFemtoSimpleAnalysis.cxx:578
 AliFemtoSimpleAnalysis.cxx:579
 AliFemtoSimpleAnalysis.cxx:580
 AliFemtoSimpleAnalysis.cxx:581
 AliFemtoSimpleAnalysis.cxx:582
 AliFemtoSimpleAnalysis.cxx:583
 AliFemtoSimpleAnalysis.cxx:584
 AliFemtoSimpleAnalysis.cxx:585
 AliFemtoSimpleAnalysis.cxx:586
 AliFemtoSimpleAnalysis.cxx:587
 AliFemtoSimpleAnalysis.cxx:588
 AliFemtoSimpleAnalysis.cxx:589
 AliFemtoSimpleAnalysis.cxx:590
 AliFemtoSimpleAnalysis.cxx:591
 AliFemtoSimpleAnalysis.cxx:592
 AliFemtoSimpleAnalysis.cxx:593
 AliFemtoSimpleAnalysis.cxx:594
 AliFemtoSimpleAnalysis.cxx:595
 AliFemtoSimpleAnalysis.cxx:596
 AliFemtoSimpleAnalysis.cxx:597
 AliFemtoSimpleAnalysis.cxx:598
 AliFemtoSimpleAnalysis.cxx:599
 AliFemtoSimpleAnalysis.cxx:600
 AliFemtoSimpleAnalysis.cxx:601
 AliFemtoSimpleAnalysis.cxx:602
 AliFemtoSimpleAnalysis.cxx:603
 AliFemtoSimpleAnalysis.cxx:604
 AliFemtoSimpleAnalysis.cxx:605
 AliFemtoSimpleAnalysis.cxx:606
 AliFemtoSimpleAnalysis.cxx:607
 AliFemtoSimpleAnalysis.cxx:608
 AliFemtoSimpleAnalysis.cxx:609
 AliFemtoSimpleAnalysis.cxx:610
 AliFemtoSimpleAnalysis.cxx:611
 AliFemtoSimpleAnalysis.cxx:612
 AliFemtoSimpleAnalysis.cxx:613
 AliFemtoSimpleAnalysis.cxx:614
 AliFemtoSimpleAnalysis.cxx:615
 AliFemtoSimpleAnalysis.cxx:616
 AliFemtoSimpleAnalysis.cxx:617
 AliFemtoSimpleAnalysis.cxx:618
 AliFemtoSimpleAnalysis.cxx:619
 AliFemtoSimpleAnalysis.cxx:620
 AliFemtoSimpleAnalysis.cxx:621
 AliFemtoSimpleAnalysis.cxx:622
 AliFemtoSimpleAnalysis.cxx:623
 AliFemtoSimpleAnalysis.cxx:624
 AliFemtoSimpleAnalysis.cxx:625
 AliFemtoSimpleAnalysis.cxx:626
 AliFemtoSimpleAnalysis.cxx:627
 AliFemtoSimpleAnalysis.cxx:628
 AliFemtoSimpleAnalysis.cxx:629
 AliFemtoSimpleAnalysis.cxx:630
 AliFemtoSimpleAnalysis.cxx:631
 AliFemtoSimpleAnalysis.cxx:632
 AliFemtoSimpleAnalysis.cxx:633
 AliFemtoSimpleAnalysis.cxx:634
 AliFemtoSimpleAnalysis.cxx:635
 AliFemtoSimpleAnalysis.cxx:636
 AliFemtoSimpleAnalysis.cxx:637
 AliFemtoSimpleAnalysis.cxx:638
 AliFemtoSimpleAnalysis.cxx:639
 AliFemtoSimpleAnalysis.cxx:640
 AliFemtoSimpleAnalysis.cxx:641
 AliFemtoSimpleAnalysis.cxx:642
 AliFemtoSimpleAnalysis.cxx:643
 AliFemtoSimpleAnalysis.cxx:644
 AliFemtoSimpleAnalysis.cxx:645
 AliFemtoSimpleAnalysis.cxx:646
 AliFemtoSimpleAnalysis.cxx:647
 AliFemtoSimpleAnalysis.cxx:648
 AliFemtoSimpleAnalysis.cxx:649
 AliFemtoSimpleAnalysis.cxx:650
 AliFemtoSimpleAnalysis.cxx:651
 AliFemtoSimpleAnalysis.cxx:652
 AliFemtoSimpleAnalysis.cxx:653
 AliFemtoSimpleAnalysis.cxx:654
 AliFemtoSimpleAnalysis.cxx:655
 AliFemtoSimpleAnalysis.cxx:656
 AliFemtoSimpleAnalysis.cxx:657
 AliFemtoSimpleAnalysis.cxx:658
 AliFemtoSimpleAnalysis.cxx:659