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