ROOT logo
/**************************************************************************
 * This file is property of and copyright by the ALICE HLT Project        *
 * ALICE Experiment at CERN, All rights reserved.                         *
 *                                                                        *
 * Primary Author: Svein Lindal <slindal@fys.uio.no>                      *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

/// @file   AliAnaConvCorrBase.cxx
/// @author Svein Lindal
/// @brief  Base class for analysation of conversion particle - track correlations


#include "AliAnaConvCorrBase.h"
#include "AliAODTrack.h"

#include "TClonesArray.h"
#include "TH1F.h"
#include "TH3.h"
#include "TList.h"
#include "AliAODConversionParticle.h"

#include <iostream>

// Gamma - jet correlation analysis task
// Authors: Svein Lindal


using namespace std;
ClassImp(AliAnaConvCorrBase)

//________________________________________________________________________
AliAnaConvCorrBase::AliAnaConvCorrBase(TString name, TString title = "title") : TNamed(name, title),
  fHistograms(NULL),
  fAxesList(), 
  fTrigAxisList(), 
  fTrackAxisList(),
  fAxistPt(),
  fAxiscPt(), 
  fAxisdEta(), 
  fAxisdPhi(),
  fAxisIso(), 
  fAxisCent(), 
  fAxisZ(),
  fAxisTrigEta(), 
  fAxisAssEta(), 
  fAxisMEPhi(),
  fCorrSparse(NULL),
  fTrigSparse(NULL),
  fTrackSparse(NULL)
{
  //Constructor
  fAxesList.SetOwner(kFALSE);
  fTrackAxisList.SetOwner(kFALSE);
  fTrigAxisList.SetOwner(kFALSE);
  SetUpDefaultBins();
 
}


//________________________________________________________________________________
AliAnaConvCorrBase::~AliAnaConvCorrBase() {
  ///destructor
}

//________________________________________________________________________________
void AliAnaConvCorrBase::CreateHistograms() {
  CreateBaseHistograms();
}

///________________________________________________________________________________
void AliAnaConvCorrBase::SetUpDefaultBins() {
  //Set up default bins
  fAxisdEta.Set(32, -1.6, 1.6);
  fAxisdEta.SetNameTitle("dEta", "delta eta");

  fAxisTrigEta.SetNameTitle("tEta", "Eta");
  fAxisTrigEta.Set(320, -0.8, 0.8);

  fAxisAssEta.SetNameTitle("aEta", "Eta");
  fAxisAssEta.Set(360, -0.9, 0.9);

  fAxisdPhi.Set(32, -TMath::PiOver2(), 3*TMath::PiOver2());
  fAxisdPhi.SetNameTitle("dPhi", "delta Phi");

  Double_t tptbins[15] = {2.0, 3.0, 3.5, 4.0, 5.0, 6.0, 8.0, 10.0, 12.5, 15, 20, 25, 30, 50, 100};
  fAxistPt.Set(13, tptbins);
  fAxistPt.SetNameTitle("tPt", "trigger Pt");

  Double_t cptbins[19] = {0.5, 1.0, 1.5, 2.0, 2.5, 3.0, 3.5, 4.0, 5.0, 6.0, 8.0, 10.0, 12.5, 15, 20, 25, 30, 50, 100};
  fAxiscPt.Set(18, cptbins);
  fAxiscPt.SetNameTitle("cPt", "track Pt");

  fAxisCent.SetNameTitle("centrality", "centrality");
  fAxisCent.Set(1, -999, 999);
  fAxisZ.SetNameTitle("vtxz", "vtxz");
  Double_t zbins[6] = {-10, -5, -1.5, 1.5, 5, 10 }; 
  fAxisZ.Set(5, zbins);

  fAxisIso.Set(1, -0.5, 2.5);
  fAxisIso.SetNameTitle("iso", "isolation");

  fAxesList.AddAt(&fAxisdEta, 0);
  fAxesList.AddAt(&fAxisdPhi, 1);
  fAxesList.AddAt(&fAxistPt, 2);
  fAxesList.AddAt(&fAxiscPt, 3);
  fAxesList.AddAt(&fAxisCent, 4);
  fAxesList.AddAt(&fAxisZ, 5);

  fTrackAxisList.AddAt(&fAxisAssEta, 0);
  fTrackAxisList.AddAt(&fAxistPt, 1);
  fTrackAxisList.AddAt(&fAxiscPt, 2);
  fTrackAxisList.AddAt(&fAxisCent, 3);
  fTrackAxisList.AddAt(&fAxisZ, 4);

  fTrigAxisList.AddAt(&fAxisTrigEta, 0);
  fTrigAxisList.AddAt(&fAxistPt, 1);
  fTrigAxisList.AddAt(&fAxisCent, 2);
  fTrigAxisList.AddAt(&fAxisZ, 3);


}


//________________________________________________________________________
void AliAnaConvCorrBase::CreateBaseHistograms() {
  //Create histograms add, to outputlis

  //cout << "Creating histograms for "<< GetName() << endl;

  fHistograms = new TList();
  fHistograms->SetOwner(kTRUE);
  fHistograms->SetName(fName);



  fCorrSparse = CreateSparse(GetName(), GetTitle(), &fAxesList);
  fHistograms->Add(fCorrSparse);

  fTrackSparse = CreateSparse(Form("%s_%s", GetName(), "METrack"), Form("%s %s", GetTitle(), "ME Tracks"), &fTrackAxisList);
  fHistograms->Add(fTrackSparse);

  fTrigSparse = CreateSparse(Form("%s_%s", GetName(), "METrig"), Form("%s %s", GetTitle(), "ME Triggers"), &fTrigAxisList);
  fHistograms->Add(fTrigSparse);

}

///________________________________________________________________________
THnSparseF * AliAnaConvCorrBase::CreateSparse(TString nameString, TString titleString, TList * axesList) {
  //Create sparse
  const Int_t dim = axesList->GetSize();
  
  //cout << nameString << " " << titleString << " " <<   "    dimesion: " << dim << endl;

  TAxis * axes[dim];
  Int_t   bins[dim];
  Double_t min[dim];
  Double_t max[dim];

  for(Int_t i = 0; i<dim; i++) {
	TAxis * axis = dynamic_cast<TAxis*>(axesList->At(i));
	if(axis) axes[i] = axis;
	else {
	  cout << "AliAnalysisTaskdPhi::CreateSparse: Error error, all the axes are not present in axis list" << endl;
	  return NULL;
	}
  }

  for(Int_t i = 0; i<dim; i++) {
	//cout << axes[i]->GetTitle() << endl;
 	bins[i] = axes[i]->GetNbins(); 
	min[i] = axes[i]->GetBinLowEdge(1);
	max[i] = axes[i]->GetBinUpEdge(axes[i]->GetNbins());
  }

  THnSparseF * sparse = new THnSparseF(Form("%s", nameString.Data()), 
									   Form("%s", titleString.Data()), 
									   dim, bins, min, max);
  
  for(Int_t i = 0; i<dim; i++) {
	sparse->GetAxis(i)->SetNameTitle(axes[i]->GetName(), axes[i]->GetTitle() );
	if(axes[i]->GetXbins()->GetSize() > 0) {
	  sparse->SetBinEdges(i, axes[i]->GetXbins()->GetArray() );
	}
  }
  return sparse;
}



//_______________________________________________________________________________
void AliAnaConvCorrBase::FillCounters(TObjArray * particles, TObjArray * tracks, Float_t cent, Float_t vtxz) {
  //Fill ME Counters
  const Int_t nbins = fAxistPt.GetNbins();
  Bool_t tmap[nbins];
  for(Int_t ptbin = 0; ptbin < nbins; ptbin++){
    tmap[ptbin] = kFALSE;
  }


  Double_t trackValues[fTrackAxisList.GetSize()];
  trackValues[3] = cent;
  trackValues[4] = vtxz;

  for(Int_t ip = 0; ip < particles->GetEntriesFast(); ip++){
    AliAODConversionParticle * particle = static_cast<AliAODConversionParticle*>(particles->At(ip));

    Int_t tbin = fAxistPt.FindFixBin(particle->Pt());
    if (tbin > 0 && tbin < nbins + 1) {
      if(tmap[tbin - 1] == kTRUE) {
	continue;
      } else {
	tmap[tbin -1 ] = kTRUE;

	if( fTrackAxisList.GetSize() > 5){
	  trackValues[5] = particle->M();
	}

	for(int ij = 0; ij < tracks->GetEntriesFast(); ij++) {
	  AliVTrack * track = static_cast<AliVTrack*>(tracks->UncheckedAt(ij));
	  trackValues[0] = track->Eta();
	  trackValues[1] = particle->Pt();
	  trackValues[2] = track->Pt();
	  fTrackSparse->Fill(trackValues);	
	}
      }
    }
  }
}

//________________________________________________________________
void AliAnaConvCorrBase::CorrelateWithTracks(AliAODConversionParticle * particle, TObjArray * tracks, Int_t const tIDs[4], Float_t cent, Float_t vtxz) {
  //Correlate particle with tracks


  const Int_t nDim = fAxesList.GetSize();
  Double_t dphivalues[nDim];
  dphivalues[4] = cent;
  dphivalues[5] = vtxz;


  Double_t trigValues[fTrigAxisList.GetSize()];
  trigValues[0] = particle->Eta();
  trigValues[1] = particle->Pt();
  trigValues[2] = cent;
  trigValues[3] = vtxz;
  
  if(nDim > 6) {
    dphivalues[6] = particle->M();
    trigValues[4] = particle->M();
  }

  fTrigSparse->Fill(trigValues);

  for(int ij = 0; ij < tracks->GetEntriesFast(); ij++) {
    AliVTrack * track = static_cast<AliVTrack*>(tracks->UncheckedAt(ij));
    
    Int_t tid = track->GetID();
    
    if((tid > 0) && (tid == tIDs[0] || tid == tIDs[1] || tid == tIDs[2] || tid == tIDs[3]) ) {
      continue;
    }
	
    dphivalues[0] = particle->Eta() - track->Eta();
    dphivalues[1] = GetDPhi(particle->Phi() - track->Phi());
    dphivalues[2] = particle->Pt();
    dphivalues[3] = track->Pt();
    fCorrSparse->Fill(dphivalues);
  }
}

 AliAnaConvCorrBase.cxx:1
 AliAnaConvCorrBase.cxx:2
 AliAnaConvCorrBase.cxx:3
 AliAnaConvCorrBase.cxx:4
 AliAnaConvCorrBase.cxx:5
 AliAnaConvCorrBase.cxx:6
 AliAnaConvCorrBase.cxx:7
 AliAnaConvCorrBase.cxx:8
 AliAnaConvCorrBase.cxx:9
 AliAnaConvCorrBase.cxx:10
 AliAnaConvCorrBase.cxx:11
 AliAnaConvCorrBase.cxx:12
 AliAnaConvCorrBase.cxx:13
 AliAnaConvCorrBase.cxx:14
 AliAnaConvCorrBase.cxx:15
 AliAnaConvCorrBase.cxx:16
 AliAnaConvCorrBase.cxx:17
 AliAnaConvCorrBase.cxx:18
 AliAnaConvCorrBase.cxx:19
 AliAnaConvCorrBase.cxx:20
 AliAnaConvCorrBase.cxx:21
 AliAnaConvCorrBase.cxx:22
 AliAnaConvCorrBase.cxx:23
 AliAnaConvCorrBase.cxx:24
 AliAnaConvCorrBase.cxx:25
 AliAnaConvCorrBase.cxx:26
 AliAnaConvCorrBase.cxx:27
 AliAnaConvCorrBase.cxx:28
 AliAnaConvCorrBase.cxx:29
 AliAnaConvCorrBase.cxx:30
 AliAnaConvCorrBase.cxx:31
 AliAnaConvCorrBase.cxx:32
 AliAnaConvCorrBase.cxx:33
 AliAnaConvCorrBase.cxx:34
 AliAnaConvCorrBase.cxx:35
 AliAnaConvCorrBase.cxx:36
 AliAnaConvCorrBase.cxx:37
 AliAnaConvCorrBase.cxx:38
 AliAnaConvCorrBase.cxx:39
 AliAnaConvCorrBase.cxx:40
 AliAnaConvCorrBase.cxx:41
 AliAnaConvCorrBase.cxx:42
 AliAnaConvCorrBase.cxx:43
 AliAnaConvCorrBase.cxx:44
 AliAnaConvCorrBase.cxx:45
 AliAnaConvCorrBase.cxx:46
 AliAnaConvCorrBase.cxx:47
 AliAnaConvCorrBase.cxx:48
 AliAnaConvCorrBase.cxx:49
 AliAnaConvCorrBase.cxx:50
 AliAnaConvCorrBase.cxx:51
 AliAnaConvCorrBase.cxx:52
 AliAnaConvCorrBase.cxx:53
 AliAnaConvCorrBase.cxx:54
 AliAnaConvCorrBase.cxx:55
 AliAnaConvCorrBase.cxx:56
 AliAnaConvCorrBase.cxx:57
 AliAnaConvCorrBase.cxx:58
 AliAnaConvCorrBase.cxx:59
 AliAnaConvCorrBase.cxx:60
 AliAnaConvCorrBase.cxx:61
 AliAnaConvCorrBase.cxx:62
 AliAnaConvCorrBase.cxx:63
 AliAnaConvCorrBase.cxx:64
 AliAnaConvCorrBase.cxx:65
 AliAnaConvCorrBase.cxx:66
 AliAnaConvCorrBase.cxx:67
 AliAnaConvCorrBase.cxx:68
 AliAnaConvCorrBase.cxx:69
 AliAnaConvCorrBase.cxx:70
 AliAnaConvCorrBase.cxx:71
 AliAnaConvCorrBase.cxx:72
 AliAnaConvCorrBase.cxx:73
 AliAnaConvCorrBase.cxx:74
 AliAnaConvCorrBase.cxx:75
 AliAnaConvCorrBase.cxx:76
 AliAnaConvCorrBase.cxx:77
 AliAnaConvCorrBase.cxx:78
 AliAnaConvCorrBase.cxx:79
 AliAnaConvCorrBase.cxx:80
 AliAnaConvCorrBase.cxx:81
 AliAnaConvCorrBase.cxx:82
 AliAnaConvCorrBase.cxx:83
 AliAnaConvCorrBase.cxx:84
 AliAnaConvCorrBase.cxx:85
 AliAnaConvCorrBase.cxx:86
 AliAnaConvCorrBase.cxx:87
 AliAnaConvCorrBase.cxx:88
 AliAnaConvCorrBase.cxx:89
 AliAnaConvCorrBase.cxx:90
 AliAnaConvCorrBase.cxx:91
 AliAnaConvCorrBase.cxx:92
 AliAnaConvCorrBase.cxx:93
 AliAnaConvCorrBase.cxx:94
 AliAnaConvCorrBase.cxx:95
 AliAnaConvCorrBase.cxx:96
 AliAnaConvCorrBase.cxx:97
 AliAnaConvCorrBase.cxx:98
 AliAnaConvCorrBase.cxx:99
 AliAnaConvCorrBase.cxx:100
 AliAnaConvCorrBase.cxx:101
 AliAnaConvCorrBase.cxx:102
 AliAnaConvCorrBase.cxx:103
 AliAnaConvCorrBase.cxx:104
 AliAnaConvCorrBase.cxx:105
 AliAnaConvCorrBase.cxx:106
 AliAnaConvCorrBase.cxx:107
 AliAnaConvCorrBase.cxx:108
 AliAnaConvCorrBase.cxx:109
 AliAnaConvCorrBase.cxx:110
 AliAnaConvCorrBase.cxx:111
 AliAnaConvCorrBase.cxx:112
 AliAnaConvCorrBase.cxx:113
 AliAnaConvCorrBase.cxx:114
 AliAnaConvCorrBase.cxx:115
 AliAnaConvCorrBase.cxx:116
 AliAnaConvCorrBase.cxx:117
 AliAnaConvCorrBase.cxx:118
 AliAnaConvCorrBase.cxx:119
 AliAnaConvCorrBase.cxx:120
 AliAnaConvCorrBase.cxx:121
 AliAnaConvCorrBase.cxx:122
 AliAnaConvCorrBase.cxx:123
 AliAnaConvCorrBase.cxx:124
 AliAnaConvCorrBase.cxx:125
 AliAnaConvCorrBase.cxx:126
 AliAnaConvCorrBase.cxx:127
 AliAnaConvCorrBase.cxx:128
 AliAnaConvCorrBase.cxx:129
 AliAnaConvCorrBase.cxx:130
 AliAnaConvCorrBase.cxx:131
 AliAnaConvCorrBase.cxx:132
 AliAnaConvCorrBase.cxx:133
 AliAnaConvCorrBase.cxx:134
 AliAnaConvCorrBase.cxx:135
 AliAnaConvCorrBase.cxx:136
 AliAnaConvCorrBase.cxx:137
 AliAnaConvCorrBase.cxx:138
 AliAnaConvCorrBase.cxx:139
 AliAnaConvCorrBase.cxx:140
 AliAnaConvCorrBase.cxx:141
 AliAnaConvCorrBase.cxx:142
 AliAnaConvCorrBase.cxx:143
 AliAnaConvCorrBase.cxx:144
 AliAnaConvCorrBase.cxx:145
 AliAnaConvCorrBase.cxx:146
 AliAnaConvCorrBase.cxx:147
 AliAnaConvCorrBase.cxx:148
 AliAnaConvCorrBase.cxx:149
 AliAnaConvCorrBase.cxx:150
 AliAnaConvCorrBase.cxx:151
 AliAnaConvCorrBase.cxx:152
 AliAnaConvCorrBase.cxx:153
 AliAnaConvCorrBase.cxx:154
 AliAnaConvCorrBase.cxx:155
 AliAnaConvCorrBase.cxx:156
 AliAnaConvCorrBase.cxx:157
 AliAnaConvCorrBase.cxx:158
 AliAnaConvCorrBase.cxx:159
 AliAnaConvCorrBase.cxx:160
 AliAnaConvCorrBase.cxx:161
 AliAnaConvCorrBase.cxx:162
 AliAnaConvCorrBase.cxx:163
 AliAnaConvCorrBase.cxx:164
 AliAnaConvCorrBase.cxx:165
 AliAnaConvCorrBase.cxx:166
 AliAnaConvCorrBase.cxx:167
 AliAnaConvCorrBase.cxx:168
 AliAnaConvCorrBase.cxx:169
 AliAnaConvCorrBase.cxx:170
 AliAnaConvCorrBase.cxx:171
 AliAnaConvCorrBase.cxx:172
 AliAnaConvCorrBase.cxx:173
 AliAnaConvCorrBase.cxx:174
 AliAnaConvCorrBase.cxx:175
 AliAnaConvCorrBase.cxx:176
 AliAnaConvCorrBase.cxx:177
 AliAnaConvCorrBase.cxx:178
 AliAnaConvCorrBase.cxx:179
 AliAnaConvCorrBase.cxx:180
 AliAnaConvCorrBase.cxx:181
 AliAnaConvCorrBase.cxx:182
 AliAnaConvCorrBase.cxx:183
 AliAnaConvCorrBase.cxx:184
 AliAnaConvCorrBase.cxx:185
 AliAnaConvCorrBase.cxx:186
 AliAnaConvCorrBase.cxx:187
 AliAnaConvCorrBase.cxx:188
 AliAnaConvCorrBase.cxx:189
 AliAnaConvCorrBase.cxx:190
 AliAnaConvCorrBase.cxx:191
 AliAnaConvCorrBase.cxx:192
 AliAnaConvCorrBase.cxx:193
 AliAnaConvCorrBase.cxx:194
 AliAnaConvCorrBase.cxx:195
 AliAnaConvCorrBase.cxx:196
 AliAnaConvCorrBase.cxx:197
 AliAnaConvCorrBase.cxx:198
 AliAnaConvCorrBase.cxx:199
 AliAnaConvCorrBase.cxx:200
 AliAnaConvCorrBase.cxx:201
 AliAnaConvCorrBase.cxx:202
 AliAnaConvCorrBase.cxx:203
 AliAnaConvCorrBase.cxx:204
 AliAnaConvCorrBase.cxx:205
 AliAnaConvCorrBase.cxx:206
 AliAnaConvCorrBase.cxx:207
 AliAnaConvCorrBase.cxx:208
 AliAnaConvCorrBase.cxx:209
 AliAnaConvCorrBase.cxx:210
 AliAnaConvCorrBase.cxx:211
 AliAnaConvCorrBase.cxx:212
 AliAnaConvCorrBase.cxx:213
 AliAnaConvCorrBase.cxx:214
 AliAnaConvCorrBase.cxx:215
 AliAnaConvCorrBase.cxx:216
 AliAnaConvCorrBase.cxx:217
 AliAnaConvCorrBase.cxx:218
 AliAnaConvCorrBase.cxx:219
 AliAnaConvCorrBase.cxx:220
 AliAnaConvCorrBase.cxx:221
 AliAnaConvCorrBase.cxx:222
 AliAnaConvCorrBase.cxx:223
 AliAnaConvCorrBase.cxx:224
 AliAnaConvCorrBase.cxx:225
 AliAnaConvCorrBase.cxx:226
 AliAnaConvCorrBase.cxx:227
 AliAnaConvCorrBase.cxx:228
 AliAnaConvCorrBase.cxx:229
 AliAnaConvCorrBase.cxx:230
 AliAnaConvCorrBase.cxx:231
 AliAnaConvCorrBase.cxx:232
 AliAnaConvCorrBase.cxx:233
 AliAnaConvCorrBase.cxx:234
 AliAnaConvCorrBase.cxx:235
 AliAnaConvCorrBase.cxx:236
 AliAnaConvCorrBase.cxx:237
 AliAnaConvCorrBase.cxx:238
 AliAnaConvCorrBase.cxx:239
 AliAnaConvCorrBase.cxx:240
 AliAnaConvCorrBase.cxx:241
 AliAnaConvCorrBase.cxx:242
 AliAnaConvCorrBase.cxx:243
 AliAnaConvCorrBase.cxx:244
 AliAnaConvCorrBase.cxx:245
 AliAnaConvCorrBase.cxx:246
 AliAnaConvCorrBase.cxx:247
 AliAnaConvCorrBase.cxx:248
 AliAnaConvCorrBase.cxx:249
 AliAnaConvCorrBase.cxx:250
 AliAnaConvCorrBase.cxx:251
 AliAnaConvCorrBase.cxx:252
 AliAnaConvCorrBase.cxx:253
 AliAnaConvCorrBase.cxx:254
 AliAnaConvCorrBase.cxx:255
 AliAnaConvCorrBase.cxx:256
 AliAnaConvCorrBase.cxx:257
 AliAnaConvCorrBase.cxx:258
 AliAnaConvCorrBase.cxx:259
 AliAnaConvCorrBase.cxx:260
 AliAnaConvCorrBase.cxx:261
 AliAnaConvCorrBase.cxx:262
 AliAnaConvCorrBase.cxx:263
 AliAnaConvCorrBase.cxx:264
 AliAnaConvCorrBase.cxx:265
 AliAnaConvCorrBase.cxx:266
 AliAnaConvCorrBase.cxx:267
 AliAnaConvCorrBase.cxx:268
 AliAnaConvCorrBase.cxx:269
 AliAnaConvCorrBase.cxx:270
 AliAnaConvCorrBase.cxx:271
 AliAnaConvCorrBase.cxx:272
 AliAnaConvCorrBase.cxx:273
 AliAnaConvCorrBase.cxx:274
 AliAnaConvCorrBase.cxx:275
 AliAnaConvCorrBase.cxx:276
 AliAnaConvCorrBase.cxx:277
 AliAnaConvCorrBase.cxx:278