]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliFastJetInput.cxx
New histograms for centrality and multiplcity checks (Gian Michele)
[u/mrichter/AliRoot.git] / JETAN / AliFastJetInput.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 #include <Riostream.h> 
17 #include <TChain.h>
18 #include <TFile.h>
19 #include <TF1.h>
20 #include <TRandom.h>
21 #include <TList.h>
22 #include <TLorentzVector.h>
23 #include <TArrayF.h>
24 #include <TClonesArray.h>
25
26 #include "AliJetHeader.h"
27 #include "AliJetReader.h"
28 #include "AliJetReaderHeader.h"
29 #include "AliJetHistos.h"
30
31 #include "AliFastJetFinder.h"
32 #include "AliFastJetHeaderV1.h"
33 #include "AliJetReaderHeader.h"
34 #include "AliJetReader.h"
35 #include "AliJetESDReader.h"
36 #include "AliJetUnitArray.h"
37 #include "AliFastJetInput.h"
38
39 #include "fastjet/PseudoJet.hh"
40 #include "fastjet/ClusterSequenceArea.hh"
41 #include "fastjet/AreaDefinition.hh"
42 #include "fastjet/JetDefinition.hh"
43 // get info on how fastjet was configured
44 #include "fastjet/config.h"
45
46 #ifdef ENABLE_PLUGIN_SISCONE
47 #include "fastjet/SISConePlugin.hh"
48 #endif
49
50 #include<sstream>  // needed for internal io
51 #include<vector> 
52 #include <cmath> 
53
54 using namespace std;
55
56
57 ClassImp(AliFastJetInput)
58
59 ////////////////////////////////////////////////////////////////////////
60
61 AliFastJetInput::AliFastJetInput():
62     fReader(0),
63     fHeader(0),
64     fInputParticles(0),
65     fInputParticlesCh(0)
66 {
67   // Default constructor
68 }
69 AliFastJetInput::AliFastJetInput(const AliFastJetInput &input):
70     TObject(),
71     fReader(input.fReader),
72     fHeader(input.fHeader),
73     fInputParticles(input.fInputParticles),
74     fInputParticlesCh(input.fInputParticlesCh)
75 {
76   // copy constructor
77 }
78 //______________________________________________________________________
79 AliFastJetInput& AliFastJetInput::operator=(const AliFastJetInput& source){
80     // Assignment operator. 
81     this->~AliFastJetInput();
82     new(this) AliFastJetInput(source);
83     return *this;
84
85 }
86 //___________________________________________________________
87 void AliFastJetInput::FillInput(){
88   //cout<<"-------- AliFastJetInput::FillInput()  ----------------"<<endl;
89
90   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
91   fInputParticles.clear();
92   fInputParticlesCh.clear();
93
94   Int_t debug  = header->GetDebug();     // debug option
95   Int_t fOpt   = fReader->GetReaderHeader()->GetDetector();
96
97   // RUN ALGORITHM  
98   // read input particles -----------------------------
99   vector<fastjet::PseudoJet> inputParticles;
100   //  cout<<"=============== AliFastJetInput::FillInput()  =========== fOpt="<<fOpt<<endl;
101   //cout<<"pointers --> fReader="<<fReader<<" header="<<header<<" fHeader="<<fHeader<<endl;
102   //cout<<"Rparam="<<Rparam<<"  ghost_etamax="<<ghost_etamax<<"  ghost_area="<<ghost_area<<endl;
103   //cout<<fReader->ClassName()<<endl;
104
105   if(fOpt==0)
106     {
107       TClonesArray *lvArray = fReader->GetMomentumArray();
108       if(lvArray == 0) { cout << "Could not get the momentum array" << endl; return; }
109       Int_t nIn =  lvArray->GetEntries();
110       if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
111       Float_t px,py,pz,en;
112       // load input vectors
113       for(Int_t i = 0; i < nIn; i++){ // loop for all input particles
114         TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
115         px = lv->Px();
116         py = lv->Py();
117         pz = lv->Pz();
118         en = lv->Energy();
119         //      cout<<"in FillInput...... "<<i<<" "<<px<<" "<<py<<"  "<<pz<<" "<<en<<endl;
120         fastjet::PseudoJet inputPart(px,py,pz,en); // create PseudoJet object
121         inputPart.set_user_index(i); //label the particle into Fastjet algortihm
122         fInputParticles.push_back(inputPart);  // back of the inputParticles vector  
123       } // end loop 
124     }
125   else {
126     TClonesArray* fUnit = fReader->GetUnitArray();
127     if(fUnit == 0) { cout << "Could not get the momentum array" << endl; return; }
128     Int_t         nIn = fUnit->GetEntries();
129     if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
130    
131     // Information extracted from fUnitArray
132     // load input vectors and calculate total energy in array
133     Float_t pt,eta,phi,theta,px,py,pz,en;
134     Int_t ipart = 0;
135     Int_t countUnit=0,countUnitNonZero=0;
136     //cout<<" nIn = "<<nIn<<endl;
137
138     for(Int_t i=0; i<nIn; i++) 
139       {
140         AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(i);
141         if(uArray->GetUnitEnergy()>0.){
142           countUnit++;
143           // It is not necessary anymore to cut on particle pt
144           pt    = uArray->GetUnitEnergy();
145           eta   = uArray->GetUnitEta();
146           phi   = uArray->GetUnitPhi();
147           theta = EtaToTheta(eta);
148           en    = (TMath::Abs(TMath::Sin(theta)) == 0) ? pt : pt/TMath::Abs(TMath::Sin(theta));
149           px    = TMath::Cos(phi)*pt;
150           py    = TMath::Sin(phi)*pt;
151           pz    = en*TMath::TanH(eta);
152           if(debug) cout << "pt: " << pt << ", eta: " << eta << ", phi: " << phi << ", en: " << en << ", px: " << px << ", py: " << py << ", pz: " << pz << endl;
153           //cout << i<<" "<<"pt: " << pt << ", eta: " << eta << ", phi: " << phi << ", en: " << en << ", px: " << px << ", py: " << py << ", pz: " << pz << endl;
154           //cout<<"in FillInput...... "<<i<<" "<<px<<" "<<py<<"  "<<pz<<" "<<en<<endl;
155
156           fastjet::PseudoJet inputPart(px,py,pz,en); // create PseudoJet object
157           inputPart.set_user_index(ipart); //label the particle into Fastjet algortihm
158           fInputParticles.push_back(inputPart);  // back of the inputParticles vector 
159           ipart++;
160
161         
162         
163           //only for charged particles (TPC+ITS)
164           TRefArray* ref = uArray->GetUnitTrackRef();
165           Int_t nRef = ref->GetEntries();
166           for(Int_t j=0; j<nRef;j++){
167             Float_t pxj=0.;  Float_t pyj=0.;  Float_t pzj=0.;Float_t enj=0.;
168             pxj = ((AliVTrack*)ref->At(j))->Px();
169             pyj = ((AliVTrack*)ref->At(j))->Py();
170             pzj = ((AliVTrack*)ref->At(j))->Pz();
171             enj=TMath::Sqrt(pxj*pxj+pyj*pyj+pzj*pzj);
172             fastjet::PseudoJet inputPartCh(pxj,pyj,pzj,enj); // create PseudoJet object
173             inputPartCh.set_user_index(((AliVTrack*)ref->At(j))->GetID()); //label the particle into Fastjet algortihm
174             fInputParticlesCh.push_back(inputPartCh);  // back of the inputParticles vector 
175
176           }
177         }
178       } // End loop on UnitArray 
179     if (debug) cout<<"countUnit(En>0) = "<<countUnit<<"  countUnit with Non ZeroSize = "<<countUnitNonZero<<endl;
180   }
181
182 }
183
184 //_____________________________________________________________________
185 Float_t  AliFastJetInput::EtaToTheta(Float_t arg)
186 {
187   //  return (180./TMath::Pi())*2.*atan(exp(-arg));
188   return 2.*atan(exp(-arg));
189
190
191 }
192 Double_t AliFastJetInput::Thermalspectrum(const Double_t *x, const Double_t *par){
193
194   return x[0]*TMath::Exp(-x[0]/par[0]);
195
196 }