]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/JCORRAN/AliJJetJtAnalysis.cxx
Install macros and scripts needed for QA
[u/mrichter/AliRoot.git] / PWGCF / Correlations / JCORRAN / AliJJetJtAnalysis.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2014, 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 // Comment describing what this class does needed!
17
18 //===========================================================
19 // Dummy comment, should be replaced by a real one
20 // comment
21 // comment
22 // Simple class for the jt anlyais by Beomkyu Kim and Dongjo Kim
23 //===========================================================
24
25 #include <TRandom.h>
26 #include <TMath.h>
27 #include <TRegexp.h>
28 #include <TVector.h>
29 #include "AliJJet.h"
30 #include "AliJEfficiency.h"
31 #include "AliJJetJtAnalysis.h"
32 #include "AliJHistManager.h"
33 #include "TClonesArray.h"
34
35 AliJJetJtAnalysis::AliJJetJtAnalysis():
36         fInputList(NULL)
37         , fJetList(NULL)
38         , fJetListOfList()
39         //, fJetBgList(NULL)
40         , fJetBgListOfList()
41     , fJetTriggPtBorders(NULL)
42     , fJetConstPtLowLimits(NULL)
43     , fJetAssocPtBorders(NULL)
44     , fDeltaRBorders(NULL)
45     , nJetContainer(0)
46         , fCard(NULL)
47         , fJJetAnalysis(NULL)
48     , fJetFinderName(0)
49     , fConeSizes(0)
50         , fEfficiency(0x0)
51         , cBin(-1)
52         , fcent(-999)
53         , zBin(-1)
54         , zVert(-999)
55     , fTracks(NULL)
56     , fHMG(NULL)
57     , fJetFinderBin()
58     , fJetTriggerBin()
59     , fTrkPtBin()
60     , fTrkLimPtBin()
61     , fdRBin()
62     , fhNumber()
63     , fhKNumber()
64     , fhJetPt()
65     , fhJetPtBin()
66     , fhZ()
67     , fhZBin()
68     , fhJt()
69     , fhJtBin()
70     , fhJtWeightBin()
71     , fhLogJtWeightBin()
72     , fhJtWithPtCutWeightBinBin()
73     , fhLogJtWithPtCutWeightBinBin()
74     , fhJtBinLimBin()
75     , fhJtWeightBinLimBin()
76     , fhLogJtWeightBinLimBin()
77     , fhJetBgPt()
78     , fhJetBgPtBin()
79     , fhBgZ()
80     , fhBgZBin()
81     , fhBgJt()
82     , fhBgJtBin()
83     , fhBgJtWeightBin()
84     , fhBgLogJtWeightBin()
85     , fhBgJtWithPtCutWeightBinBin()
86     , fhBgLogJtWithPtCutWeightBinBin()
87     , fhdeltaE()
88     , fhdeltaN()
89     , fhFullJetEChJetBin()
90     , fhFullChdRChJetBin()
91     , fh2DFullEvsChEdN0()
92     , fh2DFullEvsChEdNnot0()
93 {
94
95 }
96
97 AliJJetJtAnalysis::AliJJetJtAnalysis( AliJCard * card ):
98         fInputList(NULL)
99         , fJetList(NULL)
100         , fJetListOfList()
101         //, fJetBgList(NULL)
102         , fJetBgListOfList()
103     , fJetTriggPtBorders(NULL)
104     , fJetConstPtLowLimits(NULL)
105     , fJetAssocPtBorders(NULL)
106     , fDeltaRBorders(NULL)
107     , nJetContainer(0)
108         , fCard(card)
109         , fJJetAnalysis(NULL)
110     , fJetFinderName(0)
111     , fConeSizes(0)
112         , fEfficiency(0x0)
113         , cBin(-1)
114         , fcent(-999)
115         , zBin(-1)
116         , zVert(-999)
117     , fTracks(NULL)
118     , fHMG(NULL)
119     , fJetFinderBin()
120     , fJetTriggerBin()
121     , fTrkPtBin()
122     , fTrkLimPtBin()
123     , fdRBin()
124     , fhNumber()
125     , fhKNumber()
126     , fhJetPt()
127     , fhJetPtBin()
128     , fhZ()
129     , fhZBin()
130     , fhJt()
131     , fhJtBin()
132     , fhJtWeightBin()
133     , fhLogJtWeightBin()
134     , fhJtWithPtCutWeightBinBin()
135     , fhLogJtWithPtCutWeightBinBin()
136     , fhJtBinLimBin()
137     , fhJtWeightBinLimBin()
138     , fhLogJtWeightBinLimBin()
139     , fhJetBgPt()
140     , fhJetBgPtBin()
141     , fhBgZ()
142     , fhBgZBin()
143     , fhBgJt()
144     , fhBgJtBin()
145     , fhBgJtWeightBin()
146     , fhBgLogJtWeightBin()
147     , fhBgJtWithPtCutWeightBinBin()
148     , fhBgLogJtWithPtCutWeightBinBin()
149     , fhdeltaE()
150     , fhdeltaN()
151     , fhFullJetEChJetBin()
152     , fhFullChdRChJetBin()
153     , fh2DFullEvsChEdN0()
154     , fh2DFullEvsChEdNnot0()
155 {
156
157 }
158
159 AliJJetJtAnalysis::AliJJetJtAnalysis(const AliJJetJtAnalysis& ap) :
160         fInputList(ap.fInputList)
161         , fJetList(ap.fJetList)
162         , fJetListOfList(ap.fJetListOfList)
163         //, fJetBgList(ap.fJetBgList)
164         , fJetBgListOfList(ap.fJetBgListOfList)
165     , fJetTriggPtBorders(ap.fJetTriggPtBorders)
166     , fJetConstPtLowLimits(ap.fJetConstPtLowLimits)
167     , fJetAssocPtBorders(ap.fJetAssocPtBorders)
168     , fDeltaRBorders(ap.fDeltaRBorders)
169     , nJetContainer(ap.nJetContainer)
170         , fCard(ap.fCard)
171         , fJJetAnalysis(ap.fJJetAnalysis)
172     , fJetFinderName(ap.fJetFinderName)
173     , fConeSizes(ap.fConeSizes)
174         , fEfficiency(ap.fEfficiency)
175         , cBin(-1)
176         , fcent(-999)
177         , zBin(-1)
178         , zVert(-999)
179     , fTracks(ap.fTracks)
180     , fHMG(ap.fHMG)
181     , fJetFinderBin(ap.fJetFinderBin)
182     , fJetTriggerBin(ap.fJetTriggerBin)
183     , fTrkPtBin(ap.fTrkPtBin)
184     , fTrkLimPtBin(ap.fTrkLimPtBin)
185     , fdRBin(ap.fdRBin)
186     , fhNumber(ap.fhNumber)
187     , fhKNumber(ap.fhKNumber)
188     , fhJetPt(ap.fhJetPt)
189     , fhJetPtBin(ap.fhJetPtBin)
190     , fhZ(ap.fhZ)
191     , fhZBin(ap.fhZBin)
192     , fhJt(ap.fhJt)
193     , fhJtBin(ap.fhJtBin)
194     , fhJtWeightBin(ap.fhJtWeightBin)
195     , fhLogJtWeightBin(ap.fhLogJtWeightBin)
196     , fhJtWithPtCutWeightBinBin(ap.fhJtWithPtCutWeightBinBin)
197     , fhLogJtWithPtCutWeightBinBin(ap.fhLogJtWithPtCutWeightBinBin)
198     , fhJtBinLimBin(ap.fhJtBinLimBin)
199     , fhJtWeightBinLimBin(ap.fhJtWeightBinLimBin)
200     , fhLogJtWeightBinLimBin(ap.fhLogJtWeightBinLimBin)
201     , fhJetBgPt(ap.fhJetBgPt)
202     , fhJetBgPtBin(ap.fhJetBgPtBin)
203     , fhBgZ(ap.fhBgZ)
204     , fhBgZBin(ap.fhBgZBin)
205     , fhBgJt(ap.fhBgJt)
206     , fhBgJtBin(ap.fhBgJtBin)
207     , fhBgJtWeightBin(ap.fhBgJtWeightBin)
208     , fhBgLogJtWeightBin(ap.fhBgLogJtWeightBin)
209     , fhBgJtWithPtCutWeightBinBin(ap.fhBgJtWithPtCutWeightBinBin)
210     , fhBgLogJtWithPtCutWeightBinBin(ap.fhBgLogJtWithPtCutWeightBinBin)
211     , fhdeltaE(ap.fhdeltaE)
212     , fhdeltaN(ap.fhdeltaN)
213     , fhFullJetEChJetBin(ap.fhFullJetEChJetBin)
214     , fhFullChdRChJetBin(ap.fhFullChdRChJetBin)
215     , fh2DFullEvsChEdN0(ap.fh2DFullEvsChEdN0)
216     , fh2DFullEvsChEdNnot0(ap.fh2DFullEvsChEdNnot0)
217 {
218
219 }
220
221 AliJJetJtAnalysis& AliJJetJtAnalysis::operator = (const AliJJetJtAnalysis& ap)
222 {
223         // assignment operator
224
225         this->~AliJJetJtAnalysis();
226         new(this) AliJJetJtAnalysis(ap);
227         return *this;
228 }
229
230
231 AliJJetJtAnalysis::~AliJJetJtAnalysis(){
232
233
234     delete fJJetAnalysis;
235     fJetFinderName.clear();
236     fConeSizes.clear();
237     delete fEfficiency;
238     delete fHMG;
239    
240
241
242 }
243
244
245
246 void AliJJetJtAnalysis::UserCreateOutputObjects(){
247         //fJetListOfList always point one address in the whole time of this analysis.
248     //Thus mustn't be cleared in it's life.     
249     //fJetListOfList.Clear();
250
251
252     fJJetAnalysis = new AliJJetAnalysis();
253
254     fJetTriggPtBorders = fCard->GetVector("JetTriggPtBorders");
255     fJetConstPtLowLimits = fCard->GetVector("JetConstPtLowLimits");
256     fJetAssocPtBorders = fCard->GetVector("JetAssocPtBorders");
257     fDeltaRBorders = fCard->GetVector("DeltaRBorders");
258
259         fEfficiency = new AliJEfficiency();
260         fEfficiency->SetMode( fCard->Get("EfficiencyMode") ); // 0:NoEff, 1:Period 2:RunNum 3:Auto
261         fEfficiency->SetDataPath("alien:///alice/cern.ch/user/d/djkim/legotrain/efficieny/data"); // Efficiency root file location local or alien
262
263     TRegexp reg("R[0-9][0-9][0-9]");
264     TRegexp reg2("[0-9][0-9][0-9]");
265
266     //container name has information of cone size like **R040**
267     //this cone size information will be pulled to a numerical variable
268     nJetContainer = fJetFinderName.size();
269     fJetBgListOfList.resize(nJetContainer, TClonesArray("AliJJet",100));
270         for (int i=0; i<nJetContainer; i++){
271                 //AliJJetJtHistos *histo = new AliJJetJtHistos(fCard);
272                 //histo->CreateJetJtHistos();
273         //      fHistos.push_back( histo  );
274         TString fullNameOfiJetContainer(fJetFinderName[i]);
275         TString coneSizeName (fullNameOfiJetContainer(reg));
276         TString coneSizeValue (coneSizeName(reg2));
277         fConeSizes.push_back( (double) coneSizeValue.Atoi()/100.);
278         }
279
280
281     int NBINS=150;
282     double LogBinsX[NBINS+1], LimL=0.1, LimH=150;
283     double logBW = (log(LimH)-log(LimL))/NBINS;
284     for(int ij=0;ij<=NBINS;ij++) LogBinsX[ij]=LimL*exp(ij*logBW);
285
286
287
288     fHMG = new AliJHistManager( "AliJJetJtHistManager");
289     fJetFinderBin .Set("JetFinderOrder","NFin","NFin:%d", AliJBin::kSingle).SetBin(nJetContainer);
290     fJetTriggerBin .Set("JetTriggerBin","JetPt","p_{T,jet} : %.1f - %.1f").SetBin(fCard->GetVector("JetTriggPtBorders"));
291     fTrkPtBin .Set("TrkPtBin","TrkPt","p_{T,constituent}:%.1f-%.1f").SetBin(fCard->GetVector("JetAssocPtBorders"));
292     fTrkLimPtBin .Set("TrkLimitPtBin","TrkLimitPt","p_{T,Limit}<%.1f", AliJBin::kSingle).SetBin(fJetConstPtLowLimits->GetNoElements());
293     fdRBin.Set("dRBin","dR","dR : %.1f - %.1f ").SetBin(fCard->GetVector("DeltaRBorders"));
294
295     fhNumber
296         << TH1D("hNumber","Number",6,0,6) << fJetFinderBin
297         <<"END";
298
299     fhKNumber
300         << TH1D("hKNumber","KNumber",17,0,17) << fJetFinderBin
301         <<"END";
302
303     fhJetPt 
304         << TH1D("JetPt","",NBINS, LogBinsX ) << fJetFinderBin
305         <<"END";
306     fhJetPtBin 
307         << TH1D("JetPtBin","",NBINS, LogBinsX ) << fJetFinderBin << fJetTriggerBin
308         <<"END";
309
310     int NBINSZ=150;
311     double LogBinsZ[NBINSZ+1], LimLZ=0.001, LimHZ=1.1;
312     double logBWZ = (TMath::Log(LimHZ)-TMath::Log(LimLZ))/NBINSZ;
313     for(int ij=0;ij<=NBINSZ;ij++) LogBinsZ[ij]=LimLZ*exp(ij*logBWZ);//
314
315     fhZ 
316         << TH1D("Z","",NBINSZ, LogBinsZ ) << fJetFinderBin
317         <<"END";
318     fhZBin 
319         << TH1D("ZBin","",NBINSZ, LogBinsZ ) << fJetFinderBin << fJetTriggerBin
320         <<"END";
321
322     int NBINSJt=150;
323     double LogBinsJt[NBINSJt+1], LimLJt=0.01, LimHJt=10;
324     double logBWJt = (TMath::Log(LimHJt)-TMath::Log(LimLJt))/NBINSJt;
325     for(int ij=0;ij<=NBINSJt;ij++) LogBinsJt[ij]=LimLJt*exp(ij*logBWJt);
326     int NBINSJtW=150;
327     double LimLJtW=TMath::Log(0.01), LimHJtW=TMath::Log(10);
328
329     fhJt 
330         << TH1D("Jt","",NBINSJt, LogBinsJt ) << fJetFinderBin
331         <<"END";
332     fhJtBin 
333         << TH1D("JtBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin
334         <<"END";
335     fhJtWeightBin 
336         << TH1D("JtWeightBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin
337         <<"END";
338     fhLogJtWeightBin 
339         << TH1D("LogJtWeightBin","",NBINSJtW, LimLJtW, LimHJtW ) << fJetFinderBin << fJetTriggerBin
340         <<"END";
341
342     fhJtWithPtCutWeightBinBin 
343         << TH1D("JtWithPtCutWeightBinBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin << fTrkPtBin
344         <<"END";
345     fhLogJtWithPtCutWeightBinBin 
346         << TH1D("LogJtWeightBinBin","",NBINSJtW, LimLJtW, LimHJtW ) << fJetFinderBin << fJetTriggerBin << fTrkPtBin
347         <<"END";
348
349     fhJtBinLimBin 
350         << TH1D("JtBinLimBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin << fTrkLimPtBin
351         <<"END";
352     fhJtWeightBinLimBin 
353         << TH1D("JtWeightBinLimBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin << fTrkLimPtBin
354         <<"END";
355     fhLogJtWeightBinLimBin 
356         << TH1D("LogJtWeightBinLimBin","",NBINSJtW, LimLJtW, LimHJtW ) << fJetFinderBin << fJetTriggerBin << fTrkLimPtBin
357         <<"END";
358
359     fhJetBgPt 
360         << TH1D("JetBgPt","",NBINS, LogBinsX ) << fJetFinderBin
361         <<"END";
362     fhJetBgPtBin 
363         << TH1D("JetBgPtBin","",NBINS, LogBinsX ) << fJetFinderBin << fJetTriggerBin
364         <<"END";
365     fhBgZ 
366         << TH1D("BgZ","",NBINSZ, LogBinsZ ) << fJetFinderBin
367         <<"END";
368     fhBgZBin 
369         << TH1D("BgZBin","",NBINSZ, LogBinsZ ) << fJetFinderBin << fJetTriggerBin
370         <<"END";
371
372
373     fhBgJt 
374         << TH1D("BgJt","",NBINSJt, LogBinsJt ) << fJetFinderBin
375         <<"END";
376     fhBgJtBin 
377         << TH1D("BgJtBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin
378         <<"END";
379     fhBgJtWeightBin 
380         << TH1D("BgJtWeightBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin
381         <<"END";
382     fhBgLogJtWeightBin 
383         << TH1D("BgLogJtWeightBin","",NBINSJtW, LimLJtW, LimHJtW ) << fJetFinderBin << fJetTriggerBin
384         <<"END";
385
386     fhBgJtWithPtCutWeightBinBin 
387         << TH1D("BgJtWithPtCutWeightBinBin","",NBINSJt, LogBinsJt ) << fJetFinderBin << fJetTriggerBin << fTrkPtBin
388         <<"END";
389     fhBgLogJtWithPtCutWeightBinBin 
390         << TH1D("BgLogJtWeightBin","",NBINSJtW, LimLJtW, LimHJtW ) << fJetFinderBin << fJetTriggerBin << fTrkPtBin
391         <<"END";
392
393     
394     int NBINSdeltaN=40;
395     double LimLdeltaN=-19.5, LimHdeltaN=19.5;
396     fhdeltaN
397     << TH1D("hdeltaN","",NBINSdeltaN,LimLdeltaN,LimHdeltaN )
398     << fJetTriggerBin << fdRBin <<"END";
399
400     int NBINSdeltaE=400;
401     double LimLdeltaE=-20, LimHdeltaE=20;
402     fhdeltaE
403     << TH1D("hdeltaE","",NBINSdeltaE,LimLdeltaE,LimHdeltaE )
404     << fJetTriggerBin << fdRBin <<"END";
405
406
407     fhFullJetEChJetBin 
408         << TH1D("hFullJetEChJetBin","",NBINS, LogBinsX )  << fJetTriggerBin
409         <<"END";
410
411     int nDR = 1000;double xDR0= -10; double xDR1 = 10;
412     fhFullChdRChJetBin 
413         << TH1D("hFullChdRChJetBin","",nDR,xDR0,xDR1)  << fJetTriggerBin
414         <<"END";
415
416     fh2DFullEvsChEdN0
417         << TH2D("h2DFullEvsChEdN0","",NBINS, LogBinsX, NBINS, LogBinsX )  
418         <<"END";
419
420     fh2DFullEvsChEdNnot0
421         << TH2D("h2DFullEvsChEdNnot0","",NBINS, LogBinsX, NBINS, LogBinsX )  
422         <<"END";
423     fHMG->Print();
424     fHMG->WriteConfig();
425
426
427
428
429 }
430
431 void AliJJetJtAnalysis::ClearBeforeEvent(){
432     //fJetListOfList.Clear();
433
434
435 }
436
437 void AliJJetJtAnalysis::UserExec(){
438     for( int i=0;i<fJetListOfList.GetEntries();i++ ){
439         TObjArray * Jets = (TObjArray*) fJetListOfList[i];
440         if(!Jets) {
441             continue;
442         }
443         this->FillJtHistogram(Jets,i);
444     }
445
446
447     int iS1 = 0; //full 0.4
448     int iS2 = 3; //Ch   0.4
449     TObjArray * jetfinder1 = (TObjArray*) fJetListOfList[iS1];
450     TObjArray * jetfinder2 = (TObjArray*) fJetListOfList[iS2];
451     AliJJet *jet1 = NULL;
452     AliJJet *jet2 = NULL;
453     double deltaeta; 
454     int chEbin=-1, rbin=-1;
455     int dN=-1000;
456     double dE=-1000.;
457     for (int ijet = 0; ijet<jetfinder1->GetEntriesFast(); ijet++){
458         jet1 = dynamic_cast<AliJJet*>( jetfinder1->At(ijet) );
459         if (!jet1) continue;
460         for (int jjet = 0; jjet<jetfinder2->GetEntriesFast(); jjet++){
461             cout<<"Check new system is working"<<endl;
462             cout<<"Check new system is working2"<<endl;
463             jet2 = dynamic_cast<AliJJet*>( jetfinder2->At(jjet) );
464             if (!jet2) continue;
465             chEbin = GetBin(fJetTriggPtBorders,jet2->E());
466             deltaeta = TMath::Abs(jet1->Eta()-jet2->Eta());
467             rbin   = GetBin(fDeltaRBorders,deltaeta); 
468             fJJetAnalysis->CompareTwoJets(jet1, jet2, dE, dN);
469             if (chEbin < 0 || rbin < 0 ) continue;
470             fhdeltaE[chEbin][rbin]->Fill(dE);
471             fhdeltaN[chEbin][rbin]->Fill(dN);
472             if (dN ==0) { 
473                 fhFullJetEChJetBin[chEbin]->Fill(jet1->E());
474                 fhFullChdRChJetBin[chEbin]->Fill(jet1->DeltaR(*jet2));
475                 fh2DFullEvsChEdN0->Fill(jet1->E(), jet2->E());
476             } else { 
477                 fh2DFullEvsChEdNnot0->Fill(jet1->E(), jet2->E());
478             }
479
480         }
481     }
482
483
484 }
485
486 void AliJJetJtAnalysis::WriteHistograms(){
487
488
489     TDirectory * cwd = gDirectory;
490     //const int nJetContainer = fJetListOfList.GetEntries();
491
492
493     for (int i=0; i<nJetContainer; i++){
494         TDirectory *nwd = gDirectory->mkdir(fJetFinderName[i]);
495         nwd->cd();
496         //fHistos[i]->WriteHistograms();
497         cwd->cd();
498     }
499
500
501 }
502
503
504
505 void AliJJetJtAnalysis::FillJtHistogram( TObjArray *Jets , int iContainer)
506 {       
507
508
509
510
511
512     int iBin, iptaBin=0;
513     int jBin=0;
514     double pT = 0;
515     double conPtMax =0;
516
517     double z; double jt;
518     double pta;
519     //double Y , deltaY = 0;
520     //double Phi, deltaPhi;
521     //double deltaR= 0;
522     //cout<<"histogram filling number of jets : "<<Jets->GetEntriesFast()<<endl;
523
524     TLorentzVector  vOrtho;
525     fJetBgListOfList[iContainer].Clear();
526     TClonesArray & bgjets = fJetBgListOfList[iContainer];
527
528
529
530     int k = 0;
531     double deltaR = -1;
532     double deltaEta = -999;
533     double deltaPhi = -999;
534     double effCorrection = -1;
535     double thisConeSize = fConeSizes[iContainer] ;
536     int iBgJet = 0;
537
538     // iJet loop for an event
539     for (int i = 0; i<Jets->GetEntries(); i++){
540         AliJJet *jet = dynamic_cast<AliJJet*>( Jets->At(i) );
541         pT = jet->Pt();
542         if (pT<(*fJetTriggPtBorders)[1]) continue;
543         iBin = GetBin(fJetTriggPtBorders,pT); // fill jetPt histos
544         if( iBin < 0 ) continue;
545         fhJetPt[iContainer]->Fill( pT );
546         fhJetPtBin[iContainer][iBin]->Fill( pT );
547
548         for (int icon = 0; icon<jet->GetConstituents()->GetEntries(); icon++){
549             AliJBaseTrack *con = jet->GetConstituent(icon);
550             if (con->Pt()>conPtMax) conPtMax = con->Pt();
551         }
552
553         for (int ii = fJetConstPtLowLimits->GetNoElements(); ii >= 1 ; ii--)   {       // could also be done using GetBin( ... )
554             if (conPtMax > (*fJetConstPtLowLimits)[ii]) {               // if JetConstPtLowLimits={a,...,b} -> ConPtBinBorders={a,...,b,c}
555                 jBin = ii-1;                                                              // where c(>>b) is ''sufficiently'' high
556                 break;
557             }
558         }
559
560
561         //iConstituent loop for the iJet
562         //jt, z are calcualted and filled  
563         for (int icon = 0; icon<jet->GetConstituents()->GetEntries(); icon++){
564             AliJBaseTrack *constituent = jet->GetConstituent(icon);
565             z = (constituent->Vect()*jet->Vect().Unit())/jet->P();
566             pta = constituent->Pt();
567             constituent->SetTrackEff( fEfficiency->GetCorrection( pta, 5, fcent) );
568             effCorrection = 1.0/constituent->GetTrackEff();
569             iptaBin = GetBin(fJetAssocPtBorders, pta);
570             if( iptaBin < 0 ) continue;
571
572
573             fhZ[iContainer]->Fill( z , effCorrection);
574             fhZBin[iContainer][iBin]->Fill( z , effCorrection);
575             jt = (constituent->Vect()-z*jet->Vect()).Mag();
576             fhJt[iContainer]->Fill( jt , effCorrection);
577             fhJtBin[iContainer][iBin]->Fill( jt , effCorrection);
578             fhJtWeightBin[iContainer][iBin]->Fill( jt, 1.0/jt * effCorrection );
579             fhLogJtWeightBin[iContainer][iBin]->Fill( TMath::Log(jt), 1.0/jt * effCorrection );
580
581
582
583             if (iptaBin < 0) continue;
584             fhJtWithPtCutWeightBinBin[iContainer][iBin][iptaBin]->Fill( jt, 1.0/jt * effCorrection );
585             fhLogJtWithPtCutWeightBinBin[iContainer][iBin][iptaBin]->Fill( TMath::Log(jt), 1.0/jt * effCorrection);
586
587             for (int jj = 0; jj <= jBin ; jj++) {
588                 fhJtBinLimBin[iContainer][iBin][jj]->Fill( jt, effCorrection );
589                 fhJtWeightBinLimBin[iContainer][iBin][jj]->Fill( jt, 1.0/jt * effCorrection );
590                 fhLogJtWeightBinLimBin[iContainer][iBin][jj]->Fill( TMath::Log(jt), 1.0/jt * effCorrection );
591                 //histos->fHistosJT[0][0][iBin][jj]->Fill( TMath::Log(jt), 1.0/jt );
592             }
593
594         }
595
596
597
598         vOrtho.SetVect(jet->Vect().Orthogonal());
599         vOrtho.SetE(jet->E());
600         //if (Log) cout<<"Before R caluation, R = "<<TMath::Sqrt(it->area()/TMath::Pi())<<endl;
601         //R_area = TMath::Sqrt(it->area()/TMath::Pi())*Rs[order]/R;
602         //if (Log )cout<<"Rs[order] = "<<Rs[order]<<" R = "<<R<<" Bg R area : "<<R_area<<endl;
603
604         //Background jet (iBgJet) will be produced. This background jet is orthogonal to the iJet.  
605         //If there is another jJet, then iBgJet will be consecutevely moved not to have jJet in the cone size. 
606         if (Jets->GetEntries()>1){
607             fhNumber[iContainer]->Fill(3.5);
608             for (int j = 0; j<Jets->GetEntries(); j++){
609                 if (i == j) continue;
610                 AliJJet *jet2 = dynamic_cast<AliJJet*>( Jets->At(j) );
611
612                 if (k>15) {
613                     fhNumber[iContainer]->Fill(5.5);
614                     break;
615                 }
616
617
618                 deltaEta = vOrtho.Eta() - jet2->Eta();
619                 deltaPhi = vOrtho.Phi() - jet2->Phi();
620                 deltaR   = TMath::Sqrt(deltaEta*deltaEta + deltaPhi*deltaPhi);
621                 if ( deltaR < thisConeSize) {
622
623                     vOrtho.Rotate(TMath::Pi()/8, jet->Vect());
624                     j=0;
625                     k++;
626                     fhNumber[iContainer]->Fill(4.5);
627                 }
628
629             }
630             fhKNumber[iContainer]->Fill(k);
631         }
632
633         // Filling iBgJet,  Bgjt and Bgz
634         // "k<16" means that we will select a iBgJet which hasn't moved 
635         // more than 16 times by the process above
636         double maxconpt = 0;
637         if ( k<16 ){
638             new (bgjets[iBgJet]) AliJJet(vOrtho.Px(),vOrtho.Py(), vOrtho.Pz(), vOrtho.E(), i ,0,0);
639             AliJJet * jbg = (AliJJet*) fJetBgListOfList[iContainer][iBgJet];
640             iBgJet++;
641
642             pT = vOrtho.Pt(); 
643             if (pT<(*fJetTriggPtBorders)[1]) continue;
644
645             fhJetBgPt[iContainer]->Fill( pT );
646             //bbfHistos[iContainer]->fhJetBgPtWeight->Fill( pT, 1./pT);
647             iBin = GetBin(fJetTriggPtBorders, pT);
648             if( iBin < 0 ) continue;
649             fhJetBgPtBin[iContainer][iBin]->Fill( pT );
650
651
652             for (int icon = 0; icon<fTracks->GetEntries(); icon++){
653                 AliJBaseTrack *track = dynamic_cast<AliJBaseTrack*>(fTracks->At(icon));
654                 if (!track) continue;
655                 deltaEta = vOrtho.Eta() - track->Eta();
656                 deltaPhi = vOrtho.Phi() - track->Phi();
657                 deltaR   = TMath::Sqrt(deltaEta*deltaEta + deltaPhi*deltaPhi);
658                 if ( deltaR > thisConeSize) continue;
659         
660                 jbg->AddConstituent(track);
661
662                 pta = track->Pt();
663                 if (pta > maxconpt) maxconpt = pta;
664                 track->SetTrackEff( fEfficiency->GetCorrection( pta, 5, fcent) );
665                 effCorrection = 1.0/track->GetTrackEff();
666                 iptaBin = GetBin(fJetAssocPtBorders, pta);
667                 if( iptaBin < 0 ) continue;
668
669
670                 z = (track->Vect()*vOrtho.Vect().Unit())/vOrtho.P();
671                 fhBgZ[iContainer]->Fill( z , effCorrection);
672                 fhBgZBin[iContainer][iBin]->Fill( z , effCorrection);
673
674                 jt = (track->Vect()-z*vOrtho.Vect()).Mag();
675                 fhBgJt[iContainer]->Fill( jt , effCorrection);
676                 fhBgJtBin[iContainer][iBin]->Fill( jt , effCorrection);
677                 fhBgJtWeightBin[iContainer][iBin]->Fill( jt, 1.0/jt * effCorrection );
678                 fhBgLogJtWeightBin[iContainer][iBin]->Fill( TMath::Log(jt), 1.0/jt * effCorrection );
679
680                 if (iptaBin < 0) continue;
681                 fhBgJtWithPtCutWeightBinBin[iContainer][iBin][iptaBin]->Fill( jt, 1.0/jt * effCorrection );
682                 fhBgLogJtWithPtCutWeightBinBin[iContainer][iBin][iptaBin]->Fill( TMath::Log(jt), 1.0/jt * effCorrection );
683             } 
684         }
685
686
687     }
688 }
689
690
691