]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/DEV/AliFastJetInput.cxx
adding JETAN and FASTJETAN development libs for new i/o of tracks/particles for the...
[u/mrichter/AliRoot.git] / JETAN / DEV / 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 /* $Id$ */
17
18 //---------------------------------------------------------------------
19 // Class for input particles
20 // manages the search for jets
21 // Authors: Elena Bruna elena.bruna@yale.edu
22 //
23 // ** 2011 magali.estienne@subatech.in2p3.fr &  alexandre.shabetai@cern.ch
24 // Modified accordingly to reader/finder splitting and new handling of neutral information
25 //---------------------------------------------------------------------
26
27 #include <Riostream.h> 
28
29 #include "AliJetHeader.h"
30 #include "AliFastJetHeaderV1.h"
31 #include "AliFastJetInput.h"
32 #include "AliJetCalTrk.h"
33
34 #include "fastjet/PseudoJet.hh"
35
36 #include<vector> 
37
38 using namespace std;
39
40 ClassImp(AliFastJetInput)
41
42 ////////////////////////////////////////////////////////////////////////
43
44 AliFastJetInput::AliFastJetInput():
45   fHeader(0x0),
46   fCalTrkEvent(0x0),
47   fInputParticles(0),
48   fInputParticlesCh(0)
49 {
50   // Default constructor
51 }
52
53 //______________________________________________________________________
54 AliFastJetInput::AliFastJetInput(const AliFastJetInput &input):
55   TObject(input),
56   fHeader(input.fHeader),
57   fCalTrkEvent(input.fCalTrkEvent),
58   fInputParticles(input.fInputParticles),
59   fInputParticlesCh(input.fInputParticlesCh)
60 {
61   // copy constructor
62 }
63
64 //______________________________________________________________________
65 AliFastJetInput& AliFastJetInput::operator=(const AliFastJetInput& source)
66 {
67   // Assignment operator. 
68   if(this!=&source){
69    TObject::operator=(source);
70    fHeader = source.fHeader;
71    fCalTrkEvent = source.fCalTrkEvent;
72    fInputParticles = source.fInputParticles;
73    fInputParticlesCh = source.fInputParticlesCh;
74   }
75
76   return *this;
77
78 }
79
80 //___________________________________________________________
81 void AliFastJetInput::FillInput()
82 {
83   // fills input particles for FASTJET based analysis
84   
85   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
86   Int_t debug  = header->GetDebug();     // debug option
87
88   if(debug>0) cout<<"-------- AliFastJetInput::FillInput()  ----------------"<<endl;
89
90   fInputParticles.clear();
91   fInputParticlesCh.clear();
92
93   // RUN ALGORITHM  
94   // read input particles -----------------------------
95   vector<fastjet::PseudoJet> inputParticles;
96
97   if(fCalTrkEvent == 0) { cout << "Could not get the CalTrk Event" << endl; return; }
98   Int_t nIn =  fCalTrkEvent->GetNCalTrkTracks() ;
99   if(nIn == 0) { if (debug>0) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
100
101   // Information extracted from fCalTrkEvent
102   // load input vectors and calculate total energy in array
103   Float_t px = -999., py = -999., pz = -999., en = -999.; 
104  
105   // Fill charged tracks
106   for(Int_t i = 0; i < fCalTrkEvent->GetNCalTrkTracks(); i++)
107     { // loop for all input particles
108       if (fCalTrkEvent->GetCalTrkTrack(i)->GetCutFlag() != 1) continue;
109       px =  fCalTrkEvent->GetCalTrkTrack(i)->GetPx();
110       py =  fCalTrkEvent->GetCalTrkTrack(i)->GetPy();
111       pz =  fCalTrkEvent->GetCalTrkTrack(i)->GetPz();
112       en =  fCalTrkEvent->GetCalTrkTrack(i)->GetP();
113
114       fastjet::PseudoJet inputPart(px,py,pz,en);  // create PseudoJet object
115       inputPart.set_user_index(i);      //label the particle into Fastjet algortihm
116       fInputParticles.push_back(inputPart);       // back of the inputParticles vector
117  
118       // only for charged particles (TPC+ITS)
119       fastjet::PseudoJet inputPartCh(px,py,pz,en); // create PseudoJet object
120       inputPartCh.set_user_index(i);               //label the particle into Fastjet algortihm
121       fInputParticlesCh.push_back(inputPartCh);    // back of the inputParticles vector
122     } // End loop on CalTrk
123
124 }
125
126 //_____________________________________________________________________
127 Double_t AliFastJetInput::Thermalspectrum(const Double_t *x, const Double_t *par)
128 {
129   // compute an exponential function
130   return x[0]*TMath::Exp(-x[0]/par[0]);
131
132 }