]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliCdfJetFinder.cxx
Treatment of bad SDD half-modules (F. Prino)
[u/mrichter/AliRoot.git] / JETAN / AliCdfJetFinder.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 //---------------------------------------------------------------------
17 // Jet Finder based on CDF algortihm
18 // Charged jet evolution and the underlying event in proton-antiproton collisions at 1.8 TeV
19 // Physical Review D, vol. 65, Issue 9, id. 092002
20 // http://www.phys.ufl.edu/~rfield/cdf/chgjet/chgjet_intro.html
21 // Authors : Adrian.Sevcenco@cern.ch (adriansev@spacescience.ro )
22 //           Daniel.Felea@cern.ch    (dfelea@spacescience.ro)
23 //           Ciprian.Mihai.Mitu@cern.ch  (mcm@spacescience.ro)
24 //---------------------------------------------------------------------
25
26 /*
27 Changelog
28
29
30
31 */
32
33 #include <Riostream.h>
34 #include <TROOT.h>
35 #include <TMath.h>
36 #include <TBits.h>
37 #include <TFile.h>
38 #include <TCanvas.h>
39 #include <TClonesArray.h>
40 #include <TLorentzVector.h>
41 #include <TH1F.h>
42 #include <TH2F.h>
43 #include <TProfile.h>
44 #include <TArrayF.h>
45 #include <TVector2.h>
46
47 #include "AliJetReader.h"
48 #include "AliJetReaderHeader.h"
49 #include "AliAODJet.h"
50 #include "AliAODEvent.h"
51 #include "AliJetFinder.h"
52
53 #include "AliCdfJetFinder.h"
54 #include "AliCdfJetHeader.h"
55
56 ClassImp(AliCdfJetFinder)
57
58 //______________________________________________________________________________
59 AliCdfJetFinder::AliCdfJetFinder():
60     AliJetFinder(),
61     fHistos(0),
62     fDebug(0),
63     fFromAod(0),
64     fAODwrite(0),
65     fAODtracksWrite(0),
66     fAnalyseJets(0),
67     fRefArr (NULL),
68     fNJets(-9999),
69     fNPart(-9999),
70     fRadius(0.7),
71     fMinJetParticles(1),
72     fJetPtCut(0.),
73     fVectParticle(NULL),
74     fVectJet(NULL),
75     fPtArray(NULL),
76     fIdxArray(NULL)
77   {  /* Constructor */  }
78
79 //______________________________________________________________________________
80 AliCdfJetFinder::~AliCdfJetFinder()
81   {
82   // destructor
83   }
84
85 //______________________________________________________________________________
86 void AliCdfJetFinder::CreateOutputObjects(TList * const histos)
87 {
88   // Create the list of histograms. Only the list is owned.
89   fHistos = histos;
90
91 //  gStyle->SetOptStat(11111111);
92
93   TH1F *h1 = new TH1F ("histo1", "Pt distribution of jets", 200, 0,200);  // 1GeV/bin
94   h1->SetStats(kTRUE);
95   h1->GetXaxis()->SetTitle("P_{T} of jets");
96   h1->GetYaxis()->SetTitle("Number of jets");
97   h1->GetXaxis()->SetTitleColor(1);
98   h1->SetMarkerStyle(kFullCircle);
99   fHistos->Add(h1);
100
101   TH1F *h2 = new TH1F ("histo2", "Eta distribution of jets", 240, -1.2,1.2); // 1 unit of rapidity / 100 bin
102   h2->SetStats(kTRUE);
103   h2->GetXaxis()->SetTitle("Eta of jets");
104   h2->GetYaxis()->SetTitle("Number of jets");
105   h2->GetXaxis()->SetTitleColor(1);
106   h2->SetMarkerStyle(kFullCircle);
107   fHistos->Add(h2);
108
109   TH1F *h3 = new TH1F ("histo3", "Phi distribution of jets", 400, -4,4);
110   h3->SetStats(kTRUE);
111   h3->GetXaxis()->SetTitle("Phi of jets");
112   h3->GetYaxis()->SetTitle("Number of jets");
113   h3->GetXaxis()->SetTitleColor(1);
114   h3->SetMarkerStyle(kFullCircle);
115   fHistos->Add(h3);
116
117   TH1F *h4 = new TH1F ("histo4", "Multiplicity of jets", 40, 0,40);  // 1 unit of multiplicity /bin
118   h4->SetStats(kTRUE);
119   h4->GetXaxis()->SetTitle("Particles in jets");
120   h4->GetYaxis()->SetTitle("Number of jets");
121   h4->GetXaxis()->SetTitleColor(1);
122   h4->SetMarkerStyle(kFullCircle);
123   fHistos->Add(h4);
124
125   TH1F *h5 = new TH1F ("histo5", "Distribution of jets in events", 100, 0,100);
126   h5->SetStats(kTRUE);
127   h5->GetXaxis()->SetTitle("Number of jets");
128   h5->GetYaxis()->SetTitle("Number of events");
129   h5->GetXaxis()->SetTitleColor(1);
130   h5->SetMarkerStyle(kFullCircle);
131   fHistos->Add(h5);
132
133   TH1F *h6 = new TH1F ("histo6", "Jet1 Charged Multiplicity Distribution", 30, 0,30);
134   h6->SetStats(kTRUE);
135   h6->GetXaxis()->SetTitle("N_{chg}");
136   h6->GetYaxis()->SetTitle("Number of jets");
137   h6->GetXaxis()->SetTitleColor(1);
138   h6->SetMarkerStyle(kFullCircle);
139   fHistos->Add(h6);
140
141   TProfile * h7 = new TProfile ("histo7","N_{chg}(jet1) vs P_{T}(charged jet1)", 200, 0. ,200. , 0.,200. ) ;
142   h7->SetStats(kTRUE);
143   h7->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
144   h7->GetYaxis()->SetTitle("<N_{chg}(jet1)> in 1 GeV/c bin");
145   h7->GetXaxis()->SetTitleColor(1);
146   h7->SetMarkerStyle(kFullCircle);
147   fHistos->Add(h7);
148
149   TH1F *h8 = new TH1F ("histo8", "Charge momentum distribution for leading jet", 120, 0 , 1.2);
150   h8->SetStats(kTRUE);
151   h8->GetXaxis()->SetTitle("Jets");
152   h8->GetYaxis()->SetTitle("Particle distribution");
153   h8->GetXaxis()->SetTitleColor(1);
154   h8->SetMarkerStyle(kFullCircle);
155   fHistos->Add(h8);
156
157   TProfile *h9 = new TProfile ("histo9", "N_{chg} vs the Azimuthal Angle from Charged Jet1", 50 , 0. , 180. , 0 , 20 );
158   h9->SetStats(kTRUE);
159   h9->GetXaxis()->SetTitle("#Delta#phi (degrees)");
160   h9->GetYaxis()->SetTitle("<N_{chg}> in 3.6 degree bin");
161   h9->GetXaxis()->SetTitleColor(1);
162   h9->SetMarkerStyle(kFullCircle);
163   fHistos->Add(h9);
164
165   TProfile *h10 = new TProfile ("histo10", "P_{T} sum vs the Azimuthal Angle from Charged Jet1", 50 , 0. , 180. , 0 , 100 );
166   h10->SetStats(kTRUE);
167   h10->GetXaxis()->SetTitle("#Delta#phi (degrees)");
168   h10->GetYaxis()->SetTitle("<P_{T} sum> in 3.6 degree bin");
169   h10->GetXaxis()->SetTitleColor(1);
170   h10->SetMarkerStyle(kFullCircle);
171   fHistos->Add(h10);
172
173   TH1F *h11 = new TH1F ("histo11", " \"Transverse\" Pt Distribution ", 70, 0 , 14);
174   h11->SetStats(kTRUE);
175   h11->GetXaxis()->SetTitle("P_{T} (GeV/c)");
176   h11->GetYaxis()->SetTitle("dN_{chg}/dP_{T} (1/GeV/c)");
177   h11->GetXaxis()->SetTitleColor(1);
178   h11->SetMarkerStyle(kFullCircle);
179   fHistos->Add(h11);
180
181   TH1F *h20 = new TH1F ("histo20", "Distribution of R in leading jet", 400, 0.,4.);
182   h20->SetStats(kTRUE);
183   h20->GetXaxis()->SetTitle("R [formula]");
184   h20->GetYaxis()->SetTitle("dN/dR");
185   h20->GetXaxis()->SetTitleColor(1);
186   h20->SetMarkerStyle(kFullCircle);
187   fHistos->Add(h20);
188
189   TProfile * h21 = new TProfile ("histo21","N_{chg}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0., 50. , 0., 30. ) ;
190   h21->SetStats(kTRUE);
191   h21->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
192   h21->GetYaxis()->SetTitle("<N_{chg}(in the event - including jet1)> in 1 GeV/c bin");
193   h21->GetXaxis()->SetTitleColor(1);
194   h21->SetMarkerStyle(kFullCircle);
195   fHistos->Add(h21);
196
197   TProfile * h22 = new TProfile ("histo22","PT_{sum}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0. , 50. , 0., 50. ) ;
198   h22->SetStats(kTRUE);
199   h22->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
200   h22->GetYaxis()->SetTitle("<PT_{sum}(in the event - including jet1)> in 1 GeV/c bin");
201   h22->GetXaxis()->SetTitleColor(1);
202   h22->SetMarkerStyle(kFullCircle);
203   fHistos->Add(h22);
204
205   TProfile * h21Toward = new TProfile ("histo21_toward","N_{chg}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0., 50. , 0., 12. ) ;
206   h21Toward->SetStats(kTRUE);
207   h21Toward->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
208   h21Toward->GetYaxis()->SetTitle("<N_{chg}(in the event - including jet1)> in 1 GeV/c bin");
209   h21Toward->GetXaxis()->SetTitleColor(1);
210   h21Toward->SetMarkerStyle(kFullCircle);
211   fHistos->Add(h21Toward);
212
213   TProfile * h21Transverse = new TProfile ("histo21_transverse","N_{chg}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0., 50. , 0., 12. ) ;
214   h21Transverse->SetStats(kTRUE);
215   h21Transverse->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
216   h21Transverse->GetYaxis()->SetTitle("<N_{chg}(in the event - including jet1)> in 1 GeV/c bin");
217   h21Transverse->GetXaxis()->SetTitleColor(1);
218   h21Transverse->SetMarkerStyle(kFullCircle);
219   fHistos->Add(h21Transverse);
220
221   TProfile * h21Away = new TProfile ("histo21_away","N_{chg}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0., 50. , 0., 12. ) ;
222   h21Away->SetStats(kTRUE);
223   h21Away->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
224   h21Away->GetYaxis()->SetTitle("<N_{chg}(in the event - including jet1)> in 1 GeV/c bin");
225   h21Away->GetXaxis()->SetTitleColor(1);
226   h21Away->SetMarkerStyle(kFullCircle);
227   fHistos->Add(h21Away);
228
229   TProfile * h22Toward = new TProfile ("histo22_toward","PT_{sum}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0. , 50. , 0., 50. ) ;
230   h22Toward->SetStats(kTRUE);
231   h22Toward->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
232   h22Toward->GetYaxis()->SetTitle("<PT_{sum}(in the event - including jet1)> in 1 GeV/c bin");
233   h22Toward->GetXaxis()->SetTitleColor(1);
234   h22Toward->SetMarkerStyle(kFullCircle);
235   fHistos->Add(h22Toward);
236
237   TProfile * h22Transverse = new TProfile ("histo22_transverse","PT_{sum}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0. , 50. , 0., 50. ) ;
238   h22Transverse->SetStats(kTRUE);
239   h22Transverse->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
240   h22Transverse->GetYaxis()->SetTitle("<PT_{sum}(in the event - including jet1)> in 1 GeV/c bin");
241   h22Transverse->GetXaxis()->SetTitleColor(1);
242   h22Transverse->SetMarkerStyle(kFullCircle);
243   fHistos->Add(h22Transverse);
244
245   TProfile * h22Away = new TProfile ("histo22_away","PT_{sum}(in the event - including jet1) vs P_{T}(charged jet1)", 200, 0. , 50. , 0., 50. ) ;
246   h22Away->SetStats(kTRUE);
247   h22Away->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
248   h22Away->GetYaxis()->SetTitle("<PT_{sum}(in the event - including jet1)> in 1 GeV/c bin");
249   h22Away->GetXaxis()->SetTitleColor(1);
250   h22Away->SetMarkerStyle(kFullCircle);
251   fHistos->Add(h22Away);
252
253   TH1F *h23Toward = new TH1F ("histo23_toward","'Toward' Pt Distribution of charged particles", 200, 0., 14.);
254   h23Toward->SetStats(kTRUE);
255   h23Toward->GetXaxis()->SetTitle("P_{T} (charged) (GeV/c)");
256   h23Toward->GetYaxis()->SetTitle("dN_{chg}/dP_{T} (1/GeV/c)");
257   h23Toward->GetXaxis()->SetTitleColor(1);
258   h23Toward->SetMarkerStyle(kFullCircle);
259   fHistos->Add(h23Toward);
260
261   TH1F *h23Transverse = new TH1F ("histo23_transverse","'Transverse' Pt Distribution of charged particles", 200, 0., 14.);
262   h23Transverse->SetStats(kTRUE);
263   h23Transverse->GetXaxis()->SetTitle("P_{T} (charged) (GeV/c)");
264   h23Transverse->GetYaxis()->SetTitle("dN_{chg}/dP_{T} (1/GeV/c)");
265   h23Transverse->GetXaxis()->SetTitleColor(1);
266   h23Transverse->SetMarkerStyle(kFullCircle);
267   fHistos->Add(h23Transverse);
268
269   TH1F *h23Away = new TH1F ("histo23_away","'Away' Pt Distribution of charged particles", 200, 0., 14.);
270   h23Away->SetStats(kTRUE);
271   h23Away->GetXaxis()->SetTitle("P_{T} (charged) (GeV/c)");
272   h23Away->GetYaxis()->SetTitle("dN_{chg}/dP_{T} (1/GeV/c)");
273   h23Away->GetXaxis()->SetTitleColor(1);
274   h23Away->SetMarkerStyle(kFullCircle);
275   fHistos->Add(h23Away);
276
277   TProfile * h24 = new TProfile ("histo24","Jet1 Size vs P_{T}(charged jet1)", 200, 0., 50. , 0., 0.5) ;
278   h24->SetStats(kTRUE);
279   h24->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
280   h24->GetYaxis()->SetTitle("<R(chgjet1)> in 1 GeV/c bin");
281   h24->GetXaxis()->SetTitleColor(1);
282   h24->SetMarkerStyle(kFullCircle);
283   fHistos->Add(h24);
284
285   TProfile * h25 = new TProfile ("histo25","Jet1 Size vs P_{T}(charged jet1)", 200, 0., 50. , 0., 0.5) ;
286   h25->SetStats(kTRUE);
287   h25->GetXaxis()->SetTitle("P_{T} (charged jet1) (GeV/c)");
288   h25->GetYaxis()->SetTitle("<R(chgjet1)> in 1 GeV/c bin");
289   h25->GetXaxis()->SetTitleColor(1);
290   h25->SetMarkerStyle(kFullCircle);
291   fHistos->Add(h25);
292
293   TProfile *h26 = new TProfile ("histo26", "N_{chg} vs the Distance R from Charged Jet1", 30, 0., 0.6, 0., 0.8);
294   h26->SetStats(kTRUE);
295   h26->GetXaxis()->SetTitle("Distance R");
296   h26->GetYaxis()->SetTitle("<N_{chg}> in 0.02 bin");
297   h26->GetXaxis()->SetTitleColor(1);
298   h26->SetMarkerStyle(kFullCircle);
299   fHistos->Add(h26);
300
301   TProfile *h27 = new TProfile ("histo27", "N_{chg} vs the Distance R from Charged Jet1", 30, 0., 0.6, 0., 0.8);
302   h27->SetStats(kTRUE);
303   h27->GetXaxis()->SetTitle("Distance R");
304   h27->GetYaxis()->SetTitle("<N_{chg}> in 0.02 bin");
305   h27->GetXaxis()->SetTitleColor(1);
306   h27->SetMarkerStyle(kFullCircle);
307   fHistos->Add(h27);
308
309   TProfile *h28 = new TProfile ("histo28", "PT_{sum} vs the Distance R from Charged Jet1", 30, 0., 0.6, 0.01, 10.);
310   h28->SetStats(kTRUE);
311   h28->GetXaxis()->SetTitle("Distance R");
312   h28->GetYaxis()->SetTitle("<PT_{sum} (GeV/c)> in 0.02 bin");
313   h28->GetXaxis()->SetTitleColor(1);
314   h28->SetMarkerStyle(kFullCircle);
315   fHistos->Add(h28);
316
317   TProfile *h29 = new TProfile ("histo29", "PT_{sum} vs the Distance R from Charged Jet1", 30, 0., 0.6, 0.01, 10.);
318   h29->SetStats(kTRUE);
319   h29->GetXaxis()->SetTitle("Distance R");
320   h29->GetYaxis()->SetTitle("<PT_{sum} (GeV/c)> in 0.02 bin");
321   h29->GetXaxis()->SetTitleColor(1);
322   h29->SetMarkerStyle(kFullCircle);
323   fHistos->Add(h29);
324
325 }
326
327 //______________________________________________________________________________
328 void AliCdfJetFinder::FindJets()
329 {
330 // Jet Algorithm:
331 //  * Order all charged particles according to their PT.
332 //  * Start with the highest PT particle and include in the "jet" all particles within the "radius" R = 0.7
333 //     (considering each particle in the order of decreasing PT and recalculating the centroid of the jet after
334 //     each new particle is added to the jet).
335 //  * Go to the next highest PT particle (not already included in a jet) and include in the "jet" all particles
336 //     (not already included in a jet) within the radius R =0.7.
337 //  * Continue until all particles are in a "jet".
338
339 AliCdfJetHeader *header = (AliCdfJetHeader*)fHeader;
340
341   if (header)
342     {
343     fDebug            = header->IsDebugCDF();
344     fAODwrite         = header->IsAODwrite() ;       // write jets to AOD
345     fAODtracksWrite   = header->IsAODtracksWrite() ; // write jet tracks to AOD
346     fRadius           = header->GetRadius();      // get Radius from jet finder header
347     fMinJetParticles  = header->GetMinPartJet (); // get minimum multiplicity of an jet
348     fJetPtCut         = header->GetJetPtCut ();   // get minimum of jet pt
349     }
350   else
351     { cout << "Header not found" << endl; return; }
352
353 if (fAODwrite)
354   {
355   fFromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
356   if (fFromAod) { fRefArr = fReader->GetReferences(); }
357   }
358
359 InitData();
360
361 if (!fNPart) { if (fDebug) {cout << "entries = 0 ; Event empty !!!" << endl ;} return; } // if event empty then exit
362
363 FindCones();
364
365 ComputeConesWeight();
366
367 if (fAODwrite) { cout << "Writing AOD" << endl ; WriteJets(); }
368
369 if (fAnalyseJets) AnalizeJets();
370
371 Clean();
372
373 }
374
375 //______________________________________________________________________________
376 void AliCdfJetFinder::InitData()
377 {
378 // initialisation of variables and data members
379
380   TClonesArray * vectArray = fReader->GetMomentumArray() ;
381   if ( vectArray == 0 ) { cout << "Could not get the momentum array" << endl; return; }
382
383   fNPart = vectArray->GetEntries()  ; // n particles in this event;
384
385   if ( fNPart == 0 ) { return ; } // if event empty then exit
386
387   fVectParticle = new varContainer* [fNPart]; // container for Particles
388
389   fPtArray  = new Double_t [fNPart] ; // momentum array
390   fIdxArray = new Int_t    [fNPart] ; // index array of sorted pts
391
392   // initialisation of momentum and index arrays
393   for (  Int_t i = 0 ; i < fNPart ; i++ )
394     {// SORTING STEP :: fPtArray with data from TClonesArray of TLorentzVector
395     TLorentzVector * lv = (TLorentzVector*) vectArray->At(i);
396
397    // INITIALISATION of local arrays for temporary storage
398     varContainer *aParticle = new varContainer;
399     aParticle->pt   = lv->Pt();
400     aParticle->eta  = lv->Eta();
401     aParticle->phi  = TVector2::Phi_mpi_pi ( lv->Phi() ); // normalize to -pi,pi
402     aParticle->njet = -999;
403
404     fVectParticle[i] = aParticle;  // vector of Particles
405
406     // initializing arrays
407     fIdxArray [i] = -999 ;
408      fPtArray [i] = aParticle->pt ;
409     }
410
411   TMath::Sort ( fNPart, fPtArray, fIdxArray ) ; // get a sorted array of indexes with TClonesArray.Size()
412
413 }
414
415
416 //______________________________________________________________________________
417 void AliCdfJetFinder::FindCones()
418 {
419 // parsing of particles in event and estlabish jets (label them with jet index)
420
421   Double_t  ptSeed = 0. , etaSeed = 0. , phiSeed = 0. ; // leading particle params
422   Double_t pttmp = 0. , etatmp = 0. , phitmp = 0. ; // temporary variables to be used in various calculations
423   Double_t deta = 0. , dphi = 0. , dcomputed = 0. ;
424   Bool_t injet = 0 ;
425
426   fNJets = -1 ; // n jets in this event
427   Int_t idxPtSort = -1 ;  // index of array of sorted pt indexes
428
429   if (fDebug) { cout << "\n\n\n\n\n\n------------------\nBegin Event Analysis\n------------------\n\n" << endl ;}
430
431   cout << "fNPart = " << fNPart << endl;
432
433   TBits lkupTable ( fNPart ) ;  // bit container ; 1-to-1 corespondence with fIdxArray
434
435   while ( lkupTable.CountBits() != (UInt_t)fNPart )
436     { // loop over particles in event until all flags are set
437     UInt_t firstnonflagged = lkupTable.FirstNullBit() ; // set the index to the first NON flagged bit ; less conditions
438
439     cout << "\n\nfirst_non_flagged : " << firstnonflagged << endl;
440
441     ++fNJets; // incrementing the jet counter
442     if (fDebug) { printf("JET %d \n", fNJets); }
443
444     ptSeed = 0. ; etaSeed = 0. ; phiSeed = 0. ;  // reseting leading particle params
445
446     for (  UInt_t ipart = firstnonflagged ; ipart < (UInt_t)fNPart ; ipart++ )
447       {// iteration over particles in event
448       // the loop is done over sorted array of pt
449       idxPtSort = fIdxArray[ipart] ;  // index of particle ! fIdxArray is an index list pt sorted
450
451       if ( lkupTable.TestBitNumber(ipart) ) { continue; } // if 4vector is already flagged skip it
452
453       //init computed and used vars
454       pttmp = 0. ; etatmp = 0. ; phitmp = 0. ;
455       deta = 0. ; dphi = 0. ; dcomputed = 0. ; injet = 0 ;
456
457       //taking info from fVectParticle ;
458        pttmp = fVectParticle[idxPtSort]->pt ;
459       etatmp = fVectParticle[idxPtSort]->eta ;
460       phitmp = fVectParticle[idxPtSort]->phi ;
461
462       if ( ipart == firstnonflagged )
463         {// this is first particle in event; leading particle
464         // begin the search around this particle in a fRadius
465
466         // CENTRE OF THE JET
467         ptSeed = pttmp ; etaSeed = etatmp ; phiSeed = phitmp ; // seeding the jet with first particle idxPtSort
468
469         lkupTable.SetBitNumber ( ipart ) ; // flag the index of particle in lkup_table
470         fVectParticle[idxPtSort]->njet = fNJets ; // associate particle with current jet number
471
472         if (fDebug) { printf("\nLeading particle :: particle index = %d ;  at sorted index %d ; in jet %d \n", idxPtSort, ipart, fNJets); }
473         if (fDebug) { printf("pt= %g ; eta= %g ; phi = %g \n", pttmp, etatmp, phitmp) ; }
474         if (fDebug) { lkupTable.Print() ;}
475
476         continue ; // skip to next particle
477         }
478
479       // condition to be in jet
480       deta = etatmp - etaSeed ;
481       dphi = TVector2::Phi_mpi_pi (phitmp - phiSeed) ; // computing dphi and normalizing to (0,2pi) interval in one step
482
483       dcomputed = TMath::Hypot(deta, dphi) ; // Distance(fRadius) to (eta,phi) seed
484
485       injet = ( ( fRadius - dcomputed ) >= 0.000000001 ) ? 1 : 0 ; // if r_computed is within jet_r in_jet == 1 else 0
486
487       if ( injet )
488         { // calculus of jet variables
489         lkupTable.SetBitNumber ( ipart ) ;  // flag the index of particle in lkup_table
490         fVectParticle[idxPtSort]->njet = fNJets ; // setting in particle list the associated jet
491
492         if (fDebug) { printf("\njet particle :: particle index = %d ; at sorted index %d ; in jet %d ; found at radius %g ;  \n", idxPtSort, ipart, fNJets, dcomputed); }
493         if (fDebug) { printf("pt= %g ; eta= %g ; phi = %g \n", pttmp, etatmp, phitmp) ; }
494         if (fDebug) { lkupTable.Print() ;}
495
496         continue ; // skip to next particle
497         }
498
499       }
500       // end of iteration over event; one jet definition of content ; jet parameters to be computed later
501     }
502 }
503
504
505 //______________________________________________________________________________
506 void AliCdfJetFinder::ComputeConesWeight()
507 {
508 // computing of jets Pt, Eta and Phi (centre of weight in (eta,phi) plane)
509 // rescan the vector of particles by identify them by asociate jet number for computing of weight centre
510
511 // JET CONTAINER
512 fVectJet      = new varContainer* [fNJets]; // container for Jets
513
514 Double_t ptJet, ptJet2 , etaJet , phiJet ; Int_t npartJet ;
515 Double_t pttmp = 0. , etatmp = 0. , phitmp = 0. ; // temporary variables to be used in various calculations
516 Int_t idxPtSort = -999 ;  // index of array of sorted pt indexes
517
518 for(  Int_t jet = 0 ; jet < fNJets ; jet++ )
519   {
520   if (fDebug) { printf("\n\n--- Computing weight of Jet %d \n", jet ); }
521   npartJet = 0 ; ptJet = 0. ; etaJet = 0. ; phiJet = 0. ; // reset variables for a new computation
522
523   for (  Int_t ipart = 0 ; ipart < fNPart ; ipart++ )
524     {// iteration over particles in event
525     // the loop is done over sorted array of pt
526     idxPtSort = fIdxArray[ipart] ;  // index of particle ! fIdxArray is an index list pt sorted
527
528     if ( fVectParticle[idxPtSort]->njet == jet )
529       {
530       ++npartJet; // incrementing the counter of jet particles
531
532       //taking info from fVectParticle ;
533        pttmp = fVectParticle[idxPtSort]->pt ;
534       etatmp = fVectParticle[idxPtSort]->eta ;
535       phitmp = TVector2::Phi_mpi_pi (fVectParticle[idxPtSort]->phi) ;
536
537 //      jet_new_angular_coordinate = jet_old_angular_coordinate * jet_old_pt / jet_new_pt +
538 //                                    part[i]_angular_coordinate * part[i]_pt/jet_new_pt
539
540       ptJet2 = ptJet + pttmp ;
541
542       etaJet = etaJet * ptJet / ptJet2 +  etatmp * pttmp / ptJet2 ;
543       phiJet = phiJet * ptJet / ptJet2 +  phitmp * pttmp / ptJet2 ;
544
545       ptJet = ptJet2 ;
546
547       }
548       // add a particle and recalculation of centroid
549     }
550     // end of 1 jet computation
551
552     varContainer *aJet = new varContainer;  // Jet container
553     aJet->pt = ptJet; aJet->eta = etaJet; aJet->phi = phiJet; aJet->njet = npartJet; // setting jet vars in container
554     fVectJet[jet] = aJet;   // store the number of the jet(fNJets) and increment afterwards
555
556     if (fDebug) { printf ("=== current jet %d : npartjet= %d ; pt_jet= %g ; eta_jet = %g ; phi_jet = %g \n\n\n",
557                                        jet,     npartJet,      ptJet,      etaJet,       phiJet ) ; }
558
559   }
560   //end loop over jets
561
562 }
563
564
565 //______________________________________________________________________________
566 void AliCdfJetFinder::WriteJets()  
567
568 // Writing AOD jets and AOD tracks
569
570 for(  Int_t jetnr = 0 ; jetnr < fNJets ; jetnr++ )
571   {
572   Double_t pt = 0., eta = 0., phi = 0., // jet variables
573             px = 0., py = 0., pz = 0., en = 0.; // convert to 4-vector
574   pt  = fVectJet[ jetnr ]->pt   ; // pt  of jet
575   eta = fVectJet[ jetnr ]->eta  ; // eta of jet
576   phi = fVectJet[ jetnr ]->phi  ; // phi of jet
577
578   px = pt * TMath::Cos ( phi ) ;
579   py = pt * TMath::Sin ( phi ) ;
580   pz = pt / TMath::Tan ( 2.0 * TMath::ATan ( TMath::Exp ( -eta ) ) ) ;
581   en = TMath::Sqrt ( px * px + py * py + pz * pz );
582
583   AliAODJet jet (px, py, pz, en);
584
585
586   if (fDebug) jet.Print("");
587
588   if (fFromAod && fAODtracksWrite)
589   {
590       for (  Int_t jetTrack = 0; jetTrack < fNPart; jetTrack++ )
591       { if ( fVectParticle[jetTrack]->njet == jetnr ) { jet.AddTrack(fRefArr->At(jetTrack)) ; } }
592   }
593   // tracks REFs written in AOD
594   AddJet(jet);
595   }
596 //jets vector parsed and written to AOD
597 }
598
599
600 //______________________________________________________________________________
601 void AliCdfJetFinder::AnalizeJets()
602 {
603 // analyzing of jets and filling of histograms
604
605     const Double_t kPI = TMath::Pi();
606     
607   //persistent pointer to histo20
608   TH1F *hR = (TH1F*)fHistos->FindObject("histo20");
609
610   Int_t   *jetsptidx = 0;     // sorted array of jets pt
611   Double_t    *jetspt = 0;     // array of jets pts
612   Int_t leadingjetindex = -1 ;   // index of leading jet from fVectJet
613   Int_t partleadjet = 0 ; // number of particles in leading jet
614   Double_t   ptleadjet = 0. ; // pt  of leading jet
615   Double_t  etaleadjet = 0. ; // eta of leading jet
616   Double_t  phileadjet = 0. ; // phi of leading jet
617
618   jetsptidx = new Int_t    [fNJets] ;
619   jetspt    = new Double_t [fNJets] ;
620
621 //________________________________________________________________________________________
622 //  Jet sorting and finding the leading jet that coresponds to cuts in pt and multiplicity
623 //________________________________________________________________________________________
624
625   // filing the idx_ptjets array
626   if (fDebug) printf("List of unsorted jets:\n");
627   for(  Int_t i = 0 ; i < fNJets ; i++ )
628     {
629     jetsptidx [i] = 0 ;
630     jetspt    [i] = fVectJet[i]->pt ;
631     if (fDebug) { cout << "   jet found: " << i << " npartjet=" << fVectJet[i]->njet << " ; jets_pt = " << jetspt[i] << endl; }
632     }
633
634   TMath::Sort ( fNJets, jetspt , jetsptidx ) ; // sorting pt of jets
635
636   // selection of leading jet
637   // looping over jets searching for __first__ one that coresponds to cuts
638   for(  Int_t i = 0 ; i < fNJets ; i++ )
639     {
640     if ( ( fVectJet[ jetsptidx[i] ]->njet >= fMinJetParticles ) && ( fVectJet[ jetsptidx[i] ]->pt >= fJetPtCut ) )
641       {
642       leadingjetindex = jetsptidx[i] ;
643       partleadjet = fVectJet[ leadingjetindex ]->njet ; // number of particles in leading jet
644         ptleadjet = fVectJet[ leadingjetindex ]->pt   ; // pt  of leading jet
645        etaleadjet = fVectJet[ leadingjetindex ]->eta  ; // eta of leading jet
646        phileadjet = fVectJet[ leadingjetindex ]->phi  ; // phi of leading jet
647
648       if (fDebug)
649       { printf("Leading jet %d : npart= %d ; pt= %g ; eta = %g ; phi = %g \n", leadingjetindex, partleadjet, ptleadjet, etaleadjet, phileadjet ); }
650
651       break ;
652       }
653     }
654     // end of selection of leading jet
655
656
657
658 //////////////////////////////////////////////////
659 ////  Computing of values used in histograms
660 //////////////////////////////////////////////////
661
662 //___________________________________________________________________________
663 // pt_sum of all particles in event
664 //___________________________________________________________________________
665 cout << "Computing sum of pt in event" << endl ;
666 Double_t ptsumevent = 0.;
667 for (  Int_t i = 0 ; i< fNPart ; i++ ) { ptsumevent += fVectParticle[i]->pt ; }
668 printf ("Sum of all Pt in event : pt_sum_event = %g", ptsumevent) ;
669
670 //___________________________________________________________________________
671 // Filling an array with indexes of leading jet particles
672 //___________________________________________________________________________
673 Int_t * idxpartLJ = new Int_t [partleadjet] ;
674 Int_t counterpartleadjet = 0;
675
676 cout << "Filling an array with indexes of leading jet particles" << endl;
677
678 for( Int_t i = 0 ; i < fNPart ; i++ )
679   {
680   if ( fVectParticle[i]->njet == leadingjetindex )
681     {  idxpartLJ[counterpartleadjet++] = i ; }
682   }
683
684 if ( (counterpartleadjet-1) > partleadjet ) { cout << " Counter_part_lead_jet > part_leadjet !!!!" << endl;}
685
686
687 //___________________________________________________________________________
688 // Calculus of part distribution in leading jet
689 //___________________________________________________________________________
690 Double_t z = 0. ;
691 Double_t *zpartljet = new Double_t [ partleadjet ] ; // array of z of particles in leading jet
692
693 cout << "Entering loop of calculus of part distribution in leading jet" << endl ;
694
695 for( Int_t j = 0 ; j < partleadjet ; j++ )
696   {
697   Double_t zj = fVectParticle[idxpartLJ[j]]->pt ;
698   z =  zj / ptleadjet ;
699   zpartljet [j] = z ;
700   cout << "idx_leadjet_part[j] = " << idxpartLJ[j]
701       << " p of particle = " << zj
702       << " pt lead jet = " << ptleadjet
703       << " Z = " << z << endl;
704   }
705
706
707 //___________________________________________________________________________
708 // array of delta phi's between phi of particles and leading jet phi
709 //___________________________________________________________________________
710 cout << "array of delta phi's between phi of particles and leading jet phi" << endl;
711 Double_t dphipartLJ = 0. ;
712 Double_t *dphipartljet = new Double_t [fNPart];
713 for(  Int_t part = 0 ; part < fNPart ; part++ )
714   {
715   dphipartLJ = fVectParticle[part]->phi - phileadjet ;
716   dphipartLJ = TVector2::Phi_mpi_pi (dphipartLJ) ; // restrict the delta phi to (-pi,pi) interval
717   dphipartljet [part] = dphipartLJ ;
718   printf("part= %d ; dphi_partLJ = %g  \n", part, dphipartLJ );
719   }
720
721
722 //______________________________________________________________________________
723 //  Pt distribution for all particles
724 //______________________________________________________________________________
725 TH1F * hpt = (TH1F*)fHistos->FindObject("histo11");
726 if ( hpt ) { for (  Int_t i = 0 ; i < fNPart ; i++ ) { hpt->Fill( fVectParticle[i]->pt ); } }
727
728 //___________________________________________________________________________
729 // Recomputing of radius of particles in leading jet
730 //___________________________________________________________________________
731 if (fDebug) { printf("   Searching particles with jet index %d\n", leadingjetindex); }
732
733 Double_t ddeta = 0. , ddphi = 0. , rpart = 0. ;
734
735 for( Int_t j = 0 ; j < partleadjet ; j++ )
736   {
737   ddeta = etaleadjet - fVectParticle[idxpartLJ[j]]->eta;
738
739   Double_t phitmp = fVectParticle[idxpartLJ[j]]->phi ;
740
741   ddphi = TVector2::Phi_mpi_pi ( phileadjet - phitmp ) ; // restrict the delta phi to (-pi,pi) interval
742
743   rpart = TMath::Hypot (ddeta, ddphi) ;
744
745   printf ("Particle %d with Re-Computed radius = %f ", idxpartLJ[j], rpart) ;
746   if ( (rpart - fRadius) >= 0.00000001 )
747     { printf ("    bigger than selected radius of %f\n", fRadius ); }
748   else
749     { printf ("\n") ; }
750
751   if (hR) hR->Fill(rpart);
752
753   }
754
755
756
757 //_______________________________________________________________________
758 // Computing of radius that contain 80% of Leading Jet ( PT and multiplicity )
759 //_______________________________________________________________________
760 Double_t corepartleadjet = 0.8 * partleadjet ;
761 Double_t coreptleadjet = 0.8 * ptleadjet ;
762 Int_t countercorepart = 0 ;
763 Double_t countercorept = 0. ;
764 Int_t sortedindex = -1 ;
765
766 TProfile * hprof24 = (TProfile*)fHistos->FindObject("histo24");
767 TProfile * hprof25 = (TProfile*)fHistos->FindObject("histo25");
768
769 TProfile * hprof26 = (TProfile*)fHistos->FindObject("histo26");
770 TProfile * hprof27 = (TProfile*)fHistos->FindObject("histo27");
771 TProfile * hprof28 = (TProfile*)fHistos->FindObject("histo28");
772 TProfile * hprof29 = (TProfile*)fHistos->FindObject("histo29");
773
774
775 if ((hprof24) && (hprof25) && (hprof26) && (hprof27) && (hprof28) && (hprof29) )
776 {
777 for(  Int_t part = 0 ; part < fNPart ; part++ )
778   {
779   Double_t pttmp = 0. ; Double_t etatmp = 0. ; Double_t phitmp = 0. ; // temporary variables
780   Double_t dpart = 0. ;
781   sortedindex = fIdxArray[part] ;
782
783   if ( fVectParticle [ sortedindex ]->njet == leadingjetindex )
784     {
785      pttmp = fVectParticle[sortedindex]->pt ;
786     etatmp = fVectParticle[sortedindex]->eta ;
787     phitmp = fVectParticle[sortedindex]->phi ;
788
789     ++countercorepart ;
790     countercorept += pttmp ;
791
792     dpart = TMath::Hypot ( etaleadjet - etatmp, TVector2::Phi_mpi_pi (phileadjet - phitmp) ) ;
793
794     if ( countercorepart <=  corepartleadjet ) { hprof24->Fill(ptleadjet, dpart); }
795     if ( countercorept <= coreptleadjet ) { hprof25->Fill(ptleadjet, dpart); }
796
797     if (ptleadjet >  5.) { hprof26->Fill(dpart, countercorepart); hprof28->Fill(dpart, countercorept); }
798     if (ptleadjet > 30.) { hprof27->Fill(dpart, countercorepart); hprof29->Fill(dpart, countercorept); }
799
800     }
801   }
802 }
803
804   TH1F *hjetpt = (TH1F*)fHistos->FindObject("histo1");
805   TH1F *hjeteta = (TH1F*)fHistos->FindObject("histo2");
806   TH1F *hjetphi = (TH1F*)fHistos->FindObject("histo3");
807   TH1F *hjetnjet = (TH1F*)fHistos->FindObject("histo4");
808
809   for(  Int_t jet = 0 ; jet < fNJets ; jet++ )
810     {
811     if (hjetpt)   hjetpt   ->Fill ( fVectJet[jet]->pt   ) ;
812     if (hjeteta)  hjeteta  ->Fill ( fVectJet[jet]->eta  ) ;
813     if (hjetphi)  hjetphi  ->Fill ( fVectJet[jet]->phi  ) ;
814     if (hjetnjet) hjetnjet ->Fill ( fVectJet[jet]->njet ) ;
815     }
816
817   TH1F *hjets = (TH1F*)fHistos->FindObject("histo5");
818   if (hjets) hjets->Fill(fNJets);
819
820   TH1F *hleadpart = (TH1F*)fHistos->FindObject("histo6");
821   if (hleadpart) hleadpart->Fill(partleadjet);
822
823   TProfile * hprof = (TProfile*)fHistos->FindObject("histo7");
824   if (hprof) hprof->Fill(ptleadjet,partleadjet);
825
826   TH1F *hMD = (TH1F*)fHistos->FindObject("histo8");
827    for(  Int_t k = 0  ; k < partleadjet ; k++)
828      { hMD->Fill( zpartljet[k] ); }
829
830   TProfile * hphi = (TProfile*)fHistos->FindObject("histo9");
831     for(  Int_t k = 0  ; k < partleadjet ; k++)
832         { hphi->Fill( TMath::RadToDeg() * dphipartljet [k] , fNPart ) ; }
833
834   TProfile * htpd = (TProfile*)fHistos->FindObject("histo10");
835     for(  Int_t k = 0  ; k < fNPart ; k++)
836         { htpd->Fill( TMath::RadToDeg() * dphipartljet [k] , ptsumevent ) ; }
837
838
839   TProfile * hprof1 = (TProfile*)fHistos->FindObject("histo21");
840   if (hprof1) hprof1->Fill(ptleadjet, fNPart);
841
842   TProfile * hprof2 = (TProfile*)fHistos->FindObject("histo22");
843   if (hprof2) hprof2->Fill(ptleadjet, ptsumevent);
844
845   TProfile * hprof1toward = (TProfile*)fHistos->FindObject("histo21_toward");
846   TProfile * hprof1transverse = (TProfile*)fHistos->FindObject("histo21_transverse");
847   TProfile * hprof1away = (TProfile*)fHistos->FindObject("histo21_away");
848   TProfile * hprof2toward = (TProfile*)fHistos->FindObject("histo22_toward");
849   TProfile * hprof2transverse = (TProfile*)fHistos->FindObject("histo22_transverse");
850   TProfile * hprof2away = (TProfile*)fHistos->FindObject("histo22_away");
851   TH1F * hpttoward = (TH1F*)fHistos->FindObject("histo23_toward");
852   TH1F * hpttransverse = (TH1F*)fHistos->FindObject("histo23_transverse");
853   TH1F * hptaway = (TH1F*)fHistos->FindObject("histo23_away");
854
855   if ( (hprof1toward) && (hprof1transverse) && (hprof1away) && (hprof2toward) && (hprof2transverse) && (hprof2away) )
856     {
857     for( Int_t part = 0  ; part < fNPart ; part++)
858       {
859       Double_t ptpart = fVectParticle[part]->pt ; // pt of particle
860       if ( ( dphipartljet[part] >=0.) && ( dphipartljet[part] < kPI/3. ) )
861         {
862         hprof1toward->Fill( ptleadjet, fNPart );
863         hprof2toward->Fill( ptleadjet, ptsumevent);
864         hpttoward->Fill( ptpart );
865         }
866       else
867       if ( ( dphipartljet[part] >= (kPI/3.)) && ( dphipartljet[part] < (2.*kPI/3.)) )
868         {
869         hprof1transverse->Fill( ptleadjet, fNPart );
870         hprof2transverse->Fill( ptleadjet, ptsumevent);
871         hpttransverse->Fill( ptpart );
872         }
873       else
874       if ( ( dphipartljet[part] >= ( 2.*kPI/3.)) && ( dphipartljet[part] < kPI ) )
875         {
876         hprof1away->Fill( ptleadjet, fNPart );
877         hprof2away->Fill( ptleadjet, ptsumevent);
878         hptaway->Fill( ptpart );
879         }
880       }
881     }
882
883 }
884
885
886 //______________________________________________________________________________
887 void AliCdfJetFinder::Clean()
888 {
889 // CLEANING SECTION
890   delete [] fVectParticle;
891   delete [] fVectJet;
892   delete [] fPtArray;
893   delete [] fIdxArray;
894
895   Reset();
896 }
897
898 //______________________________________________________________________________
899 void AliCdfJetFinder::FinishRun()
900 {
901 // do i need this?
902 }
903