]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliUA1JetFinderV2.cxx
Initialize TRefArray with correct pid
[u/mrichter/AliRoot.git] / JETAN / AliUA1JetFinderV2.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17  
18 //---------------------------------------------------------------------
19 // UA1 Cone Algorithm Jet finder for charged + neutral jet studies
20 // manages the search for jets using charged particle momentum and 
21 // neutral cell energy information
22 // Based on UA1 V1 (from R. Diaz)
23 // Author: magali.estienne@subatech.in2p3.fr
24 //---------------------------------------------------------------------
25
26 #include <TClonesArray.h>
27 #include <TH1F.h>
28 #include <TH2F.h>
29 #include <TLorentzVector.h>
30 #include <TMath.h>
31 #include <TRefArray.h>
32 #include "TFile.h"
33
34 #include "AliUA1JetFinderV2.h"
35 #include "AliUA1JetHeaderV1.h"
36 #include "AliJetUnitArray.h"
37 #include "AliJetReaderHeader.h"
38 #include "AliJetReader.h"
39 #include "AliJetHeader.h"
40
41 class TArrayF;
42 class TFile;
43 class AliJetReader;
44 class AliAODJet;
45
46 ClassImp(AliUA1JetFinderV2)
47
48
49 ////////////////////////////////////////////////////////////////////////
50 AliUA1JetFinderV2::AliUA1JetFinderV2() :
51   AliJetFinder(),
52   fLego(0),  
53   fDebug(0),
54   fOpt(0)
55 {
56   //
57   // Constructor
58   //
59 }
60
61 ////////////////////////////////////////////////////////////////////////
62 AliUA1JetFinderV2::~AliUA1JetFinderV2()
63 {
64   //
65   // Destructor
66   //
67 }
68
69 ////////////////////////////////////////////////////////////////////////
70 void AliUA1JetFinderV2::FindJetsC()
71
72   // 
73   //  Used to find jets using charged particle momentum information
74   //
75   //  1) Fill cell map array
76   //  2) calculate total energy and fluctuation level
77   //  3) Run algorithm
78   //     3.1) look centroides in cell map
79   //     3.2) calculate total energy in cones
80   //     3.3) flag as a possible jet
81   //     3.4) reorder cones by energy
82   //  4) subtract backg in accepted jets
83   //  5) fill AliJet list
84
85   //  Transform input to pt,eta,phi plus lego
86     
87   AliUA1JetHeaderV1* header  = (AliUA1JetHeaderV1*) fHeader;
88   TClonesArray*      lvArray = fReader->GetMomentumArray();
89   Int_t              nIn     = lvArray->GetEntries();
90   
91   if (nIn == 0) return;
92   
93   // local arrays for input
94   Float_t* ptT    = new Float_t[nIn];
95   Float_t* etaT   = new Float_t[nIn];
96   Float_t* phiT   = new Float_t[nIn];
97   Float_t* cFlagT = new Float_t[nIn]; // Temporarily added
98   Float_t* sFlagT = new Float_t[nIn]; // Temporarily added
99   Int_t*   injet  = new Int_t[nIn];
100   
101   //total energy in array
102   Float_t  etbgTotal = 0.0;
103   TH1F*    hPtTotal  = new TH1F("hPt","Pt distribution of all particles ",100,0.0,15.0);
104   
105   // load input vectors and calculate total energy in array
106   for (Int_t i = 0; i < nIn; i++){
107     TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
108     ptT[i]  = lv->Pt();
109     etaT[i] = lv->Eta();
110     phiT[i] = ((lv->Phi() < 0) ? (lv->Phi()) + 2 * TMath::Pi() : lv->Phi());
111     cFlagT[i] = fReader->GetCutFlag(i); 
112     sFlagT[i] = fReader->GetSignalFlag(i); 
113     
114     if (fReader->GetCutFlag(i) != 1) continue;
115     fLego->Fill(etaT[i], phiT[i], ptT[i]);
116     hPtTotal->Fill(ptT[i]);
117     etbgTotal+= ptT[i];
118   }
119   
120   
121   // calculate total energy and fluctuation in map
122   Double_t meanpt   = hPtTotal->GetMean();
123   Double_t ptRMS    = hPtTotal->GetRMS();
124   Double_t npart    = hPtTotal->GetEntries();
125   Double_t dEtTotal = (TMath::Sqrt(npart))*TMath::Sqrt(meanpt * meanpt + ptRMS*ptRMS);
126   
127   // arrays to hold jets
128   Float_t* etaJet    = new Float_t[30];  // eta jet
129   Float_t* phiJet    = new Float_t[30];  // phi jet
130   Float_t* etJet     = new Float_t[30];  // et jet
131   Float_t* etsigJet  = new Float_t[30];  // signal et in jet
132   Float_t* etallJet  = new Float_t[30];  // total et in jet (tmp variable)
133   Int_t*   ncellsJet = new Int_t[30];
134   Int_t*   multJet   = new Int_t[30];
135   //--- Added for jet reordering at the end of the jet finding procedure
136   Float_t* etaJetOk    = new Float_t[30];
137   Float_t* phiJetOk    = new Float_t[30];
138   Float_t* etJetOk     = new Float_t[30];
139   Float_t* etsigJetOk  = new Float_t[30];  // signal et in jet
140   Float_t* etallJetOk  = new Float_t[30];  // total et in jet (tmp variable)
141   Int_t*   ncellsJetOk = new Int_t[30];
142   Int_t*   multJetOk   = new Int_t[30];
143   //--------------------------
144   Int_t nJets; // to hold number of jets found by algorithm
145   Int_t nj;    // number of jets accepted
146   Float_t prec  = header->GetPrecBg();
147   Float_t bgprec = 1;
148   while(bgprec > prec){
149     //reset jet arrays in memory
150     memset(etaJet,0,sizeof(Float_t)*30);
151     memset(phiJet,0,sizeof(Float_t)*30);
152     memset(etJet,0,sizeof(Float_t)*30);
153     memset(etallJet,0,sizeof(Float_t)*30);
154     memset(etsigJet,0,sizeof(Float_t)*30);
155     memset(ncellsJet,0,sizeof(Int_t)*30);
156     memset(multJet,0,sizeof(Int_t)*30);
157     //--- Added for jet reordering at the end of the jet finding procedure
158     memset(etaJetOk,0,sizeof(Float_t)*30);
159     memset(phiJetOk,0,sizeof(Float_t)*30);
160     memset(etJetOk,0,sizeof(Float_t)*30);
161     memset(etallJetOk,0,sizeof(Float_t)*30);
162     memset(etsigJetOk,0,sizeof(Float_t)*30);
163     memset(ncellsJetOk,0,sizeof(Int_t)*30);
164     memset(multJetOk,0,sizeof(Int_t)*30);
165     //--------------------------
166     nJets = 0;
167     nj = 0;
168     
169     // reset particles-jet array in memory
170     memset(injet,-1,sizeof(Int_t)*nIn);
171     //run cone algorithm finder
172     RunAlgoritmC(etbgTotal,dEtTotal,nJets,etJet,etaJet,phiJet,etallJet,ncellsJet);
173     
174     //run background subtraction
175     if(nJets > header->GetNAcceptJets()) // limited number of accepted jets per event
176       nj = header->GetNAcceptJets();
177     else
178       nj = nJets;
179     //subtract background
180     Float_t etbgTotalN = 0.0; //new background
181     if(header->GetBackgMode() == 1) // standard
182       SubtractBackgC(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
183     if(header->GetBackgMode() == 2) //cone
184       SubtractBackgCone(nIn,nj,etbgTotalN,ptT,etaT,phiT,cFlagT,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
185     if(header->GetBackgMode() == 3) //ratio
186       SubtractBackgRatio(nIn,nj,etbgTotalN,ptT,etaT,phiT,cFlagT,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
187     if(header->GetBackgMode() == 4) //statistic
188       SubtractBackgStat(nIn,nj,etbgTotalN,ptT,etaT,phiT,cFlagT,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
189     //calc precision
190     if(etbgTotalN != 0.0)
191       bgprec = (etbgTotal - etbgTotalN)/etbgTotalN;
192     else
193       bgprec = 0;
194     etbgTotal = etbgTotalN; // update with new background estimation
195   } //end while
196   
197   // add jets to list
198   Int_t* idxjets = new Int_t[nj];
199   Int_t nselectj = 0;
200   printf("Found %d jets \n", nj);
201
202   // Reorder jets by et in cone
203   Int_t * idx  = new Int_t[nJets];
204   TMath::Sort(nJets, etJet, idx);
205   for(Int_t p = 0; p < nJets; p++){
206     etaJetOk[p]    = etaJet[idx[p]];
207     phiJetOk[p]    = phiJet[idx[p]];
208     etJetOk[p]     = etJet[idx[p]];
209     etallJetOk[p]  = etJet[idx[p]];
210     etsigJetOk[p]  = etsigJet[idx[p]];
211     ncellsJetOk[p] = ncellsJet[idx[p]];
212     multJetOk[p]   = multJet[idx[p]];
213   }
214   
215   for(Int_t kj=0; kj<nj; kj++)
216     {
217       if ((etaJetOk[kj] > (header->GetJetEtaMax())) ||
218           (etaJetOk[kj] < (header->GetJetEtaMin())) ||
219           (etJetOk[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
220       Float_t px, py,pz,en; // convert to 4-vector
221       px = etJetOk[kj] * TMath::Cos(phiJetOk[kj]);
222       py = etJetOk[kj] * TMath::Sin(phiJetOk[kj]);
223       pz = etJetOk[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJetOk[kj])));
224       en = TMath::Sqrt(px * px + py * py + pz * pz);
225       
226       AliAODJet jet(px, py, pz, en);
227       jet.Print("");
228       
229       AddJet(jet);
230       
231       idxjets[nselectj] = kj;
232       nselectj++;
233     }
234
235   //add signal percentage and total signal  in AliJets for analysis tool
236   Float_t* percentage  = new Float_t[nselectj];
237   Int_t* ncells      = new Int_t[nselectj];
238   Int_t* mult        = new Int_t[nselectj];
239   for(Int_t i = 0; i< nselectj; i++)
240     {
241       percentage[i] = etsigJetOk[idxjets[i]]/etJetOk[idxjets[i]];
242       ncells[i] = ncellsJetOk[idxjets[i]];
243       mult[i] = multJetOk[idxjets[i]];
244     }
245
246   //add particle-injet relationship ///
247   for(Int_t bj = 0; bj < nIn; bj++)
248     {
249       if(injet[bj] == -1) continue; //background particle
250       Int_t bflag = 0;
251       for(Int_t ci = 0; ci< nselectj; ci++){
252         if(injet[bj] == idxjets[ci]){
253           injet[bj]= ci;
254           bflag++;
255           break;
256         }
257       }
258       if(bflag == 0) injet[bj] = -1; // set as background particle
259     }
260
261  
262   //delete
263   delete[] ptT;
264   delete[] etaT;
265   delete[] phiT;
266   delete[] cFlagT;
267   delete[] sFlagT;
268   delete[] injet;
269   delete[] hPtTotal;
270   delete[] etaJet;
271   delete[] phiJet;
272   delete[] etJet;
273   delete[] etsigJet;
274   delete[] etallJet;
275   delete[] ncellsJet;
276   delete[] multJet;
277   delete[] idxjets;
278   delete[] percentage;
279   delete[] ncells;
280   delete[] mult;
281   //--- Added for jet reordering
282   delete etaJetOk;
283   delete phiJetOk;
284   delete etJetOk;
285   delete etsigJetOk;
286   delete etallJetOk;
287   delete ncellsJetOk;
288   delete multJetOk;
289   //--------------------------
290
291 }
292
293 ////////////////////////////////////////////////////////////////////////
294 void AliUA1JetFinderV2::FindJets()
295 {
296   // 
297   //  Used to find jets using charged particle momentum information 
298   //  & neutral energy from calo cells
299   //
300   //  1) Fill cell map array
301   //  2) calculate total energy and fluctuation level
302   //  3) Run algorithm
303   //     3.1) look centroides in cell map
304   //     3.2) calculate total energy in cones
305   //     3.3) flag as a possible jet
306   //     3.4) reorder cones by energy
307   //  4) subtract backg in accepted jets
308   //  5) fill AliJet list
309
310   // transform input to pt,eta,phi plus lego
311     
312   AliUA1JetHeaderV1* header   = (AliUA1JetHeaderV1*) fHeader;
313   TClonesArray*      fUnit    = fReader->GetUnitArray();
314   Int_t              nCand    = fReader->GetNumCandidate();
315   Int_t              nCandCut = fReader->GetNumCandidateCut();
316   Int_t              nIn      = fUnit->GetEntries();
317   Float_t            ptMin   = fReader->GetReaderHeader()->GetPtCut();
318
319   fDebug   = fReader->GetReaderHeader()->GetDebug();
320   if (nIn == 0) return;
321
322   Int_t nCandidateCut = 0;
323   Int_t nCandidate = 0;
324   
325   nCandidate = nCand;
326   nCandidateCut = nCandCut;
327
328   // local arrays for input No Cuts
329   // Both pt < ptMin and pt > ptMin
330   Float_t*   ptT       = new Float_t[nCandidate];
331   Float_t*   en2T      = new Float_t[nCandidate];
332   Float_t*   pt2T      = new Float_t[nCandidate];
333   Int_t*     detT      = new Int_t[nCandidate]; 
334   Float_t*   etaT      = new Float_t[nCandidate];
335   Float_t*   phiT      = new Float_t[nCandidate];
336   Float_t*   cFlagT    = new Float_t[nCandidate];
337   Float_t*   cFlag2T   = new Float_t[nCandidate];
338   Float_t*   sFlagT    = new Float_t[nCandidate];
339   Float_t*   cClusterT = new Float_t[nCandidate];
340   Int_t*     vectT     = new Int_t[nCandidate];
341   Int_t      loop1     = 0;
342   Int_t*     injet     = new Int_t[nCandidate];
343   Int_t*     sflag     = new Int_t[nCandidate];
344   TRefArray* trackRef  = new TRefArray();
345
346   //total energy in array
347   Float_t  etbgTotal = 0.0;
348   TH1F*    hPtTotal  = new TH1F("hPt","Pt distribution of all particles ",100,0.0,15.0);
349
350   // Input cell info
351   Float_t *etCell    = new Float_t[nIn];   //! Cell Energy - Extracted from UnitArray
352   Float_t *etaCell   = new Float_t[nIn];  //! Cell eta - Extracted from UnitArray
353   Float_t *phiCell   = new Float_t[nIn];  //! Cell phi - Extracted from UnitArray
354   Int_t   *flagCell  = new Int_t[nIn];   //! Cell phi - Extracted from UnitArray
355   Float_t *etCell2   = new Float_t[nIn];  //! Cell Energy - Extracted from UnitArray
356   Float_t *etaCell2  = new Float_t[nIn]; //! Cell eta - Extracted from UnitArray
357   Float_t *phiCell2  = new Float_t[nIn]; //! Cell phi - Extracted from UnitArray
358   Int_t   *flagCell2 = new Int_t[nIn];  //! Cell phi - Extracted from UnitArray
359
360   // Information extracted from fUnitArray
361   // Load input vectors and calculate total energy in array
362   for(Int_t i=0; i<nIn; i++) 
363     {
364       // Recover particle information from UnitArray
365       
366       AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(i);
367       TRefArray* ref = uArray->GetUnitTrackRef();
368       Int_t nRef = ref->GetEntries();
369
370       if(uArray->GetUnitEnergy()>0.){
371
372         for(Int_t jpart=0; jpart<nRef;jpart++)
373           trackRef->Add((AliVTrack*)ref->At(jpart));
374         ptT[loop1]   = uArray->GetUnitEnergy();
375         detT[loop1]  = uArray->GetUnitDetectorFlag();
376         etaT[loop1]  = uArray->GetUnitEta();
377         phiT[loop1]  = uArray->GetUnitPhi();
378         cFlagT[loop1]= uArray->GetUnitCutFlag();   // pt cut tpc
379         cFlag2T[loop1]= uArray->GetUnitCutFlag2(); // pt cut emcal
380         sFlagT[loop1]= uArray->GetUnitSignalFlag();
381         vectT[loop1] = nRef;
382         if(cFlagT[loop1] == 1 || cFlag2T[loop1] == 1) {
383           pt2T[loop1] = 0.;
384           en2T[loop1] = 0.;
385           if(detT[loop1]==1){
386             en2T[loop1] = ptT[loop1] - header->GetMinCellEt();
387             if(en2T[loop1] < 0) en2T[loop1]=0;
388             hPtTotal->Fill(en2T[loop1]);
389             etbgTotal += en2T[loop1];
390           }
391           if(detT[loop1]==0){ // TPC+ITS
392             Float_t pt = 0.;
393             for(Int_t j=0; j<nRef;j++){
394               Float_t x=0.;  Float_t y=0.;  Float_t z=0.;
395               x = ((AliVTrack*)ref->At(j))->Px();
396               y = ((AliVTrack*)ref->At(j))->Py();
397               z = ((AliVTrack*)ref->At(j))->Pz();
398               pt = TMath::Sqrt(x*x+y*y);
399               if(pt>ptMin) {
400                 pt2T[loop1] += pt;
401                 en2T[loop1] += pt;
402                 hPtTotal->Fill(pt);
403                 etbgTotal+= pt;
404               }
405             }
406           }
407           if(detT[loop1]==2) { // EMCal
408             Float_t ptCTot = 0.;
409             Float_t pt = 0.;
410             Float_t enC = 0.;
411             for(Int_t j=0; j<nRef;j++){
412               Float_t x=0.;  Float_t y=0.;  Float_t z=0.;
413               x = ((AliVTrack*)ref->At(j))->Px();
414               y = ((AliVTrack*)ref->At(j))->Py();
415               z = ((AliVTrack*)ref->At(j))->Pz();
416               pt = TMath::Sqrt(x*x+y*y);
417               if(pt>ptMin) {
418                 pt2T[loop1]+=pt;
419                 en2T[loop1]+=pt;
420                 hPtTotal->Fill(pt);
421                 etbgTotal+= pt;
422               } 
423               ptCTot += pt;
424             }
425             enC = ptT[loop1] - ptCTot - header->GetMinCellEt();
426             if(enC < 0.) enC=0.;
427             en2T[loop1] += enC;
428             hPtTotal->Fill(enC);
429             etbgTotal+= enC;
430           }
431         }       
432         loop1++;
433       }
434
435       if(uArray->GetUnitCutFlag()==1) {
436         if(uArray->GetUnitDetectorFlag()==1){ // EMCal case
437           etCell[i] = uArray->GetUnitEnergy() - header->GetMinCellEt();
438           if ((uArray->GetUnitEnergy() - header->GetMinCellEt()) < 0.0) etCell[i]=0.;
439           etaCell[i] = uArray->GetUnitEta();
440           phiCell[i] = uArray->GetUnitPhi();
441           flagCell[i] = 0; // default
442           etCell2[i] = etCell[i];
443           etaCell2[i] = uArray->GetUnitEta();
444           phiCell2[i] = uArray->GetUnitPhi();
445           flagCell2[i] = 0; // default
446         }
447         if(uArray->GetUnitDetectorFlag()==0){ // TPC case
448           Float_t pt = 0.; Float_t et1 = 0.; Float_t et2 = 0.;
449           for(Int_t j=0; j<nRef;j++)
450             {
451               Float_t x=0.;  Float_t y=0.;  Float_t z=0.;
452               x = ((AliVTrack*)ref->At(j))->Px();
453               y = ((AliVTrack*)ref->At(j))->Py();
454               z = ((AliVTrack*)ref->At(j))->Pz();
455               pt = TMath::Sqrt(x*x+y*y);
456               if(pt>ptMin) {
457                 et1 += pt;
458                 et2 += pt;
459               }
460             }
461           etCell[i] = et1;
462           etCell2[i] = et2;
463           if(et1 < 0.) etCell[i] = etCell2[i] = 0.;
464           etaCell[i] = uArray->GetUnitEta();
465           phiCell[i] = uArray->GetUnitPhi();
466           flagCell[i] = 0; // default
467           etaCell2[i] = uArray->GetUnitEta();
468           phiCell2[i] = uArray->GetUnitPhi();
469           flagCell2[i] = 0; // default
470         }
471         if(uArray->GetUnitDetectorFlag()==2){ // TPC + EMCal case
472           Float_t ptCTot = 0.;
473           Float_t pt = 0.; Float_t et1 = 0.; Float_t et2 = 0.;
474           Float_t enC = 0.;
475           for(Int_t j=0; j<nRef;j++)
476             {
477               Float_t x=0.;  Float_t y=0.;  Float_t z=0.;
478               x = ((AliVTrack*)ref->At(j))->Px();
479               y = ((AliVTrack*)ref->At(j))->Py();
480               z = ((AliVTrack*)ref->At(j))->Pz();
481               pt = TMath::Sqrt(x*x+y*y);
482               if(pt>ptMin) {
483                 et1 += pt;
484                 et2 += pt;
485               }
486               ptCTot += pt;
487             }
488           enC = uArray->GetUnitEnergy() - ptCTot;
489           etCell[i] = et1 + enC - header->GetMinCellEt();
490           etCell2[i] = et2 + enC - header->GetMinCellEt();
491           if((enC + et1 - header->GetMinCellEt()) < 0.) etCell[i] = etCell2[i] = 0.;
492           etaCell[i] = uArray->GetUnitEta();
493           phiCell[i] = uArray->GetUnitPhi();
494           flagCell[i] = 0; // default
495           etaCell2[i] = uArray->GetUnitEta();
496           phiCell2[i] = uArray->GetUnitPhi();
497           flagCell2[i] = 0; // default
498         }
499       }
500       else {
501         etCell[i]  = 0.;
502         etaCell[i] = uArray->GetUnitEta();
503         phiCell[i] = uArray->GetUnitPhi();
504         flagCell[i] = 0;
505         etCell2[i]  = 0.;
506         etaCell2[i] = uArray->GetUnitEta();
507         phiCell2[i] = uArray->GetUnitPhi();
508         flagCell2[i] = 0;
509       }
510     } // end loop on nCandidate
511
512
513   // calculate total energy and fluctuation in map
514   Double_t meanpt = hPtTotal->GetMean();
515   Double_t ptRMS = hPtTotal->GetRMS();
516   Double_t npart = hPtTotal->GetEntries();
517   Double_t dEtTotal = (TMath::Sqrt(npart))*TMath::Sqrt(meanpt * meanpt + ptRMS*ptRMS);
518
519   // arrays to hold jets
520   Float_t* etaJet    = new Float_t[30];
521   Float_t* phiJet    = new Float_t[30];
522   Float_t* etJet     = new Float_t[30];
523   Float_t* etsigJet  = new Float_t[30]; //signal et in jet
524   Float_t* etallJet  = new Float_t[30];  // total et in jet (tmp variable)
525   Int_t*   ncellsJet = new Int_t[30];
526   Int_t*   multJet   = new Int_t[30];
527   //--- Added by me for jet reordering at the end of the jet finding procedure
528   Float_t* etaJetOk    = new Float_t[30];
529   Float_t* phiJetOk    = new Float_t[30];
530   Float_t* etJetOk     = new Float_t[30];
531   Float_t* etsigJetOk  = new Float_t[30]; //signal et in jet
532   Float_t* etallJetOk  = new Float_t[30];  // total et in jet (tmp variable)
533   Int_t*   ncellsJetOk = new Int_t[30];
534   Int_t*   multJetOk   = new Int_t[30];
535   //--------------------------
536   Int_t    nJets; // to hold number of jets found by algorithm
537   Int_t    nj;    // number of jets accepted
538   Float_t  prec  = header->GetPrecBg();
539   Float_t  bgprec = 1;
540
541   while(bgprec > prec){
542
543     //reset jet arrays in memory
544     memset(etaJet,0,sizeof(Float_t)*30);
545     memset(phiJet,0,sizeof(Float_t)*30);
546     memset(etJet,0,sizeof(Float_t)*30);
547     memset(etallJet,0,sizeof(Float_t)*30);
548     memset(etsigJet,0,sizeof(Float_t)*30);
549     memset(ncellsJet,0,sizeof(Int_t)*30);
550     memset(multJet,0,sizeof(Int_t)*30);
551     //--- Added by me for jet reordering at the end of the jet finding procedure
552     memset(etaJetOk,0,sizeof(Float_t)*30);
553     memset(phiJetOk,0,sizeof(Float_t)*30);
554     memset(etJetOk,0,sizeof(Float_t)*30);
555     memset(etallJetOk,0,sizeof(Float_t)*30);
556     memset(etsigJetOk,0,sizeof(Float_t)*30);
557     memset(ncellsJetOk,0,sizeof(Int_t)*30);
558     memset(multJetOk,0,sizeof(Int_t)*30);
559
560     nJets = 0;  
561     nj = 0;
562
563     // reset particles-jet array in memory
564     memset(injet,-1,sizeof(Int_t)*nCandidate);
565     //run cone algorithm finder
566     RunAlgoritm(nIn,etCell,etaCell,phiCell,flagCell,etCell2,etaCell2,phiCell2,
567                 flagCell2,etbgTotal,dEtTotal,nJets,etJet,etaJet,phiJet,
568                 etallJet,ncellsJet);
569
570     //run background subtraction
571     if(nJets > header->GetNAcceptJets()) // limited number of accepted jets per event
572       nj = header->GetNAcceptJets();
573     else
574       nj = nJets;
575
576     //subtract background
577     Float_t etbgTotalN = 0.0; //new background
578     if(header->GetBackgMode() == 1) // standard
579       SubtractBackg(nCandidate,nj,etbgTotalN,en2T,vectT,etaT,phiT,cFlagT,cFlag2T,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
580     // To be modified ------------------------
581     if(header->GetBackgMode() == 2) //cone
582       SubtractBackgCone(nCandidate,nj,etbgTotalN,ptT,etaT,phiT,cFlagT,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
583     if(header->GetBackgMode() == 3) //ratio
584       SubtractBackgRatio(nCandidate,nj,etbgTotalN,ptT,etaT,phiT,cFlagT,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
585     if(header->GetBackgMode() == 4) //statistic
586       SubtractBackgStat(nCandidate,nj,etbgTotalN,ptT,etaT,phiT,cFlagT,sFlagT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
587     //----------------------------------------
588     //calc precision
589     if(etbgTotalN != 0.0)
590       bgprec = (etbgTotal - etbgTotalN)/etbgTotalN;
591     else
592       bgprec = 0;
593     etbgTotal = etbgTotalN; // update with new background estimation
594   } //end while
595   
596   // add jets to list
597   Int_t* idxjets = new Int_t[nj];
598   Int_t nselectj = 0;
599   printf("Found %d jets \n", nj);
600   
601   // Reorder jets by et in cone
602   // Sort jets by energy
603   Int_t * idx  = new Int_t[nJets];
604   TMath::Sort(nJets, etJet, idx);
605   for(Int_t p = 0; p < nJets; p++)
606     {
607       etaJetOk[p]    = etaJet[idx[p]];
608       phiJetOk[p]    = phiJet[idx[p]];
609       etJetOk[p]     = etJet[idx[p]];
610       etallJetOk[p]  = etJet[idx[p]];
611       etsigJetOk[p]  = etsigJet[idx[p]];
612       ncellsJetOk[p] = ncellsJet[idx[p]];
613       multJetOk[p]   = multJet[idx[p]];
614     }
615
616   for(Int_t kj=0; kj<nj; kj++)
617     {
618       if ((etaJetOk[kj] > (header->GetJetEtaMax())) ||
619           (etaJetOk[kj] < (header->GetJetEtaMin())) ||
620           (etJetOk[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
621       Float_t px, py,pz,en; // convert to 4-vector
622       px = etJetOk[kj] * TMath::Cos(phiJetOk[kj]);
623       py = etJetOk[kj] * TMath::Sin(phiJetOk[kj]);
624       pz = etJetOk[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJetOk[kj])));
625       en = TMath::Sqrt(px * px + py * py + pz * pz);
626
627       AliAODJet jet(px, py, pz, en);
628       jet.Print("");
629       
630       AddJet(jet);
631       
632       idxjets[nselectj] = kj;
633       nselectj++;
634     }
635
636   //add signal percentage and total signal  in AliJets for analysis tool
637   Float_t* percentage  = new Float_t[nselectj];
638   Int_t* ncells      = new Int_t[nselectj];
639   Int_t* mult        = new Int_t[nselectj];
640   for(Int_t i = 0; i< nselectj; i++)
641     {
642       percentage[i] = etsigJetOk[idxjets[i]]/etJetOk[idxjets[i]];
643       ncells[i]     = ncellsJetOk[idxjets[i]];
644       mult[i]       = multJetOk[idxjets[i]];
645     }
646
647   //add particle-injet relationship ///
648   for(Int_t bj = 0; bj < nCandidate; bj++)
649     {
650       if(injet[bj] == -1) continue; //background particle
651       Int_t bflag = 0;
652       for(Int_t ci = 0; ci< nselectj; ci++){
653         if(injet[bj] == idxjets[ci]){
654           injet[bj]= ci;
655           bflag++;
656           break;
657         }
658       }
659     if(bflag == 0) injet[bj] = -1; // set as background particle
660   }
661
662
663   //delete
664   delete ptT;
665   delete en2T;
666   delete pt2T;
667   delete etaT;
668   delete phiT;
669   trackRef->Delete();
670   delete trackRef;
671   delete detT;
672   delete cFlagT;
673   delete cFlag2T;
674   delete sFlagT;
675   delete cClusterT;
676   delete vectT;
677   delete injet;
678   delete sflag;
679   delete hPtTotal;
680   delete etCell;
681   delete etaCell;
682   delete phiCell;
683   delete flagCell;
684   delete etCell2;
685   delete etaCell2;
686   delete phiCell2;
687   delete flagCell2;
688   delete etaJet;
689   delete phiJet;
690   delete etJet;
691   delete etsigJet;
692   delete etallJet;
693   delete ncellsJet;
694   delete multJet;
695   //--- Added for jet reordering
696   delete etaJetOk;
697   delete phiJetOk;
698   delete etJetOk;
699   delete etsigJetOk;
700   delete etallJetOk;
701   delete ncellsJetOk;
702   delete multJetOk;
703   //--------------------------
704   delete idxjets;
705   delete percentage;
706   delete ncells;
707   delete mult;
708
709 }
710
711 ////////////////////////////////////////////////////////////////////////
712 void AliUA1JetFinderV2::RunAlgoritm(Int_t nIn, Float_t* etCell, Float_t* etaCell, Float_t* phiCell, 
713                                     Int_t* flagCell, const Float_t* etCell2, const Float_t* etaCell2, 
714                                     const Float_t* phiCell2, const Int_t* flagCell2, Float_t etbgTotal, 
715                                     Double_t dEtTotal, Int_t& nJets, Float_t* etJet,Float_t* etaJet, 
716                                     Float_t* phiJet, Float_t* etallJet, Int_t* ncellsJet)
717 {
718   //
719   // Main method for jet finding
720   // UA1 base cone finder
721   //
722
723   Int_t nCell  = nIn; 
724   fDebug = fReader->GetReaderHeader()->GetDebug();
725
726   // Dump lego
727   // Check enough space! *to be done*
728   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
729   for(Int_t i=0; i<nCell; i++){
730     etCell[i]   = etCell2[i];
731     etaCell[i]  = etaCell2[i];
732     phiCell[i]  = phiCell2[i];
733     flagCell[i] = flagCell2[i];
734   }
735
736   // Parameters from header
737   Float_t minmove = header->GetMinMove();
738   Float_t maxmove = header->GetMaxMove();
739   Float_t rc      = header->GetRadius();
740   Float_t etseed  = header->GetEtSeed();
741
742   // Tmp array of jets form algoritm
743   Float_t etaAlgoJet[30];
744   Float_t phiAlgoJet[30];
745   Float_t etAlgoJet[30];
746   Int_t   ncellsAlgoJet[30];
747
748   // Run algorithm//
749
750   // Sort cells by et
751   Int_t * index  = new Int_t[nCell];
752   TMath::Sort(nCell, etCell, index);
753
754   // Variable used in centroide loop
755   Float_t eta   = 0.0;
756   Float_t phi   = 0.0;
757   Float_t eta0  = 0.0;
758   Float_t phi0  = 0.0;
759   Float_t etab  = 0.0;
760   Float_t phib  = 0.0;
761   Float_t etas  = 0.0;
762   Float_t phis  = 0.0;
763   Float_t ets   = 0.0;
764   Float_t deta  = 0.0;
765   Float_t dphi  = 0.0;
766   Float_t dr    = 0.0;
767   Float_t etsb  = 0.0;
768   Float_t etasb = 0.0;
769   Float_t phisb = 0.0;
770   Float_t dphib = 0.0;
771
772   for(Int_t icell = 0; icell < nCell; icell++)
773     {
774       Int_t jcell = index[icell];
775       if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
776       if(flagCell[jcell] != 0) continue; // if cell was used before
777
778       eta  = etaCell[jcell];
779       phi  = phiCell[jcell];
780       eta0 = eta;
781       phi0 = phi;
782       etab = eta;
783       phib = phi;
784       ets  = etCell[jcell];
785       etas = 0.0;
786       phis = 0.0;
787       etsb = ets;
788       etasb = 0.0;
789       phisb = 0.0;
790       for(Int_t kcell =0; kcell < nCell; kcell++)
791         {
792           Int_t lcell = index[kcell];
793           if(lcell == jcell) continue; // cell itself
794           if(flagCell[lcell] != 0) continue; // cell used before
795           if(etCell[lcell] > etCell[jcell]) continue;  // can this happen
796           //calculate dr
797           deta = etaCell[lcell] - eta;
798           dphi = TMath::Abs(phiCell[lcell] - phi);
799           if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
800           dr = TMath::Sqrt(deta * deta + dphi * dphi);
801           if(dr <= rc){
802             // calculate offset from initiate cell
803             deta = etaCell[lcell] - eta0;
804             dphi = phiCell[lcell] - phi0;
805             if (dphi < -TMath::Pi()) dphi= dphi + 2.0 * TMath::Pi();
806             if (dphi > TMath::Pi()) dphi = dphi - 2.0 * TMath::Pi();
807             etas = etas + etCell[lcell]*deta;
808             phis = phis + etCell[lcell]*dphi;
809             ets = ets + etCell[lcell];
810             //new weighted eta and phi including this cell
811             eta = eta0 + etas/ets;
812             phi = phi0 + phis/ets;
813             // if cone does not move much, just go to next step
814             dphib = TMath::Abs(phi - phib);
815             if (dphib > TMath::Pi()) dphib = 2. * TMath::Pi() - dphib;
816             dr = TMath::Sqrt((eta-etab)*(eta-etab) + dphib * dphib);
817             if(dr <= minmove) break;
818             // cone should not move more than max_mov
819             dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
820             if(dr > maxmove){
821               eta = etab;
822               phi = phib;
823               ets = etsb;
824               etas = etasb;
825               phis = phisb;
826             } else { // store this loop information
827               etab  = eta;
828               phib  = phi;
829               etsb  = ets;
830               etasb = etas;
831               phisb = phis;
832             }
833           } // inside cone
834         }//end of cells loop looking centroide
835
836         //avoid cones overloap (to be implemented in the future)
837
838         //flag cells in Rc, estimate total energy in cone
839       Float_t etCone = 0.0;
840       Int_t   nCellIn = 0;
841       Int_t   nCellOut = 0;
842       rc = header->GetRadius();
843
844       for(Int_t ncell =0; ncell < nCell; ncell++)
845         {
846           if(flagCell[ncell] != 0) continue; // cell used before
847           //calculate dr
848           deta = etaCell[ncell] - eta;
849           //        if(deta <= rc){ // Added to improve velocity -> to be tested
850           dphi = phiCell[ncell] - phi;
851           if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
852           if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
853           //          if(dphi <= rc){ // Added to improve velocity -> to be tested
854           dr = TMath::Sqrt(deta * deta + dphi * dphi);
855           if(dr <= rc){  // cell in cone
856             flagCell[ncell] = -1;
857             etCone+=etCell[ncell];
858             nCellIn++;
859           }
860           else nCellOut++;
861           //    } // end deta <= rc
862           //        } // end dphi <= rc
863         }
864
865       // select jets with et > background
866       // estimate max fluctuation of background in cone
867       Double_t ncellin = (Double_t)nCellIn;
868       Double_t ntcell  = (Double_t)nCell;
869       Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/(ntcell));
870       // min cone et
871       Double_t etcmin = etCone ;  // could be used etCone - etmin !!
872       //decisions !! etbmax < etcmin
873
874       for(Int_t mcell =0; mcell < nCell; mcell++)
875         {
876           if(flagCell[mcell] == -1){
877             if(etbmax < etcmin)
878               flagCell[mcell] = 1; //flag cell as used
879             else
880               flagCell[mcell] = 0; // leave it free
881           }
882         }
883       //store tmp jet info !!!
884       if(etbmax < etcmin) 
885         {
886           etaAlgoJet[nJets]    = eta;
887           phiAlgoJet[nJets]    = phi;
888           etAlgoJet[nJets]     = etCone;
889           ncellsAlgoJet[nJets] = nCellIn;
890           nJets++;
891         }
892
893     } // end of cells loop
894
895   for(Int_t p = 0; p < nJets; p++)
896     {
897       etaJet[p]    = etaAlgoJet[p];
898       phiJet[p]    = phiAlgoJet[p];
899       etJet[p]     = etAlgoJet[p];
900       etallJet[p]  = etAlgoJet[p];
901       ncellsJet[p] = ncellsAlgoJet[p];
902     }
903
904   //delete
905   delete index;
906
907 }
908
909 ////////////////////////////////////////////////////////////////////////
910 void AliUA1JetFinderV2::RunAlgoritmC(Float_t etbgTotal, Double_t dEtTotal, Int_t& nJets,
911                                   Float_t* etJet,Float_t* etaJet, Float_t* phiJet,
912                                   Float_t* etallJet, Int_t* ncellsJet)
913 {
914   // Dump lego
915   // Check enough space! *to be done*
916   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
917   Float_t etCell[60000];   //! Cell Energy
918   Float_t etaCell[60000];  //! Cell eta
919   Float_t phiCell[60000];  //! Cell phi
920   Int_t   flagCell[60000]; //! Cell flag
921   
922   Int_t nCell = 0;
923   TAxis* xaxis = fLego->GetXaxis();
924   TAxis* yaxis = fLego->GetYaxis();
925   Float_t e = 0.0;
926   for (Int_t i = 1; i <= header->GetLegoNbinEta(); i++) 
927     {
928       for (Int_t j = 1; j <= header->GetLegoNbinPhi(); j++)
929         {
930           e = fLego->GetBinContent(i,j);
931           if (e < 0.0) continue; // don't include this cells
932           Float_t eta  = xaxis->GetBinCenter(i);
933           Float_t phi  = yaxis->GetBinCenter(j);
934           etCell[nCell]  = e;
935           etaCell[nCell] = eta;
936           phiCell[nCell] = phi;
937           flagCell[nCell] = 0; //default
938           nCell++;
939         }
940     }
941
942   // Parameters from header
943   Float_t minmove = header->GetMinMove();
944   Float_t maxmove = header->GetMaxMove();
945   Float_t rc      = header->GetRadius();
946   Float_t etseed  = header->GetEtSeed();
947
948   // Tmp array of jets form algoritm
949   Float_t etaAlgoJet[30];
950   Float_t phiAlgoJet[30];
951   Float_t etAlgoJet[30];
952   Int_t   ncellsAlgoJet[30];
953
954   // Run algorithm//
955
956   // Sort cells by et
957   Int_t * index  = new Int_t[nCell];
958   TMath::Sort(nCell, etCell, index);
959   // variable used in centroide loop
960   Float_t eta   = 0.0;
961   Float_t phi   = 0.0;
962   Float_t eta0  = 0.0;
963   Float_t phi0  = 0.0;
964   Float_t etab  = 0.0;
965   Float_t phib  = 0.0;
966   Float_t etas  = 0.0;
967   Float_t phis  = 0.0;
968   Float_t ets   = 0.0;
969   Float_t deta  = 0.0;
970   Float_t dphi  = 0.0;
971   Float_t dr    = 0.0;
972   Float_t etsb  = 0.0;
973   Float_t etasb = 0.0;
974   Float_t phisb = 0.0;
975   Float_t dphib = 0.0;
976
977   for(Int_t icell = 0; icell < nCell; icell++)
978     {
979       Int_t jcell = index[icell];
980       if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
981       if(flagCell[jcell] != 0) continue; // if cell was used before
982       
983       eta  = etaCell[jcell];
984       phi  = phiCell[jcell];
985       eta0 = eta;
986       phi0 = phi;
987       etab = eta;
988       phib = phi;
989       ets  = etCell[jcell];
990       etas = 0.0;
991       phis = 0.0;
992       etsb = ets;
993       etasb = 0.0;
994       phisb = 0.0;
995       for(Int_t kcell =0; kcell < nCell; kcell++)
996         {
997           Int_t lcell = index[kcell];
998           if(lcell == jcell) continue; // cell itself
999           if(flagCell[lcell] != 0) continue; // cell used before
1000           if(etCell[lcell] > etCell[jcell]) continue; // can this happen
1001           //calculate dr
1002           deta = etaCell[lcell] - eta;
1003           dphi = TMath::Abs(phiCell[lcell] - phi);
1004           if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1005           dr = TMath::Sqrt(deta * deta + dphi * dphi);
1006           if(dr <= rc)
1007             {
1008               // calculate offset from initiate cell
1009               deta = etaCell[lcell] - eta0;
1010               dphi = phiCell[lcell] - phi0;
1011               if (dphi < -TMath::Pi()) dphi= dphi + 2.0 * TMath::Pi();
1012               if (dphi > TMath::Pi()) dphi = dphi - 2.0 * TMath::Pi();
1013               etas = etas + etCell[lcell]*deta;
1014               phis = phis + etCell[lcell]*dphi;
1015               ets = ets + etCell[lcell];
1016               //new weighted eta and phi including this cell
1017               eta = eta0 + etas/ets;
1018               phi = phi0 + phis/ets;
1019               // if cone does not move much, just go to next step
1020               dphib = TMath::Abs(phi - phib);
1021               if (dphib > TMath::Pi()) dphib = 2. * TMath::Pi() - dphib;
1022               dr = TMath::Sqrt((eta-etab)*(eta-etab) + dphib * dphib);
1023               if(dr <= minmove) break;
1024               // cone should not move more than max_mov
1025               dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
1026               if(dr > maxmove){
1027                 eta = etab;
1028                 phi = phib;
1029                 ets = etsb;
1030                 etas = etasb;
1031                 phis = phisb;
1032               } else { // store this loop information
1033                 etab=eta;
1034                 phib=phi;
1035                 etsb = ets;
1036                 etasb = etas;
1037                 phisb = phis;
1038                 }
1039             } // inside cone
1040         }//end of cells loop looking centroide
1041       
1042       // Avoid cones overloap (to be implemented in the future)
1043       
1044       // Flag cells in Rc, estimate total energy in cone
1045       Float_t etCone   = 0.0;
1046       Int_t   nCellIn  = 0;
1047       Int_t   nCellOut = 0;
1048       rc = header->GetRadius();
1049       for(Int_t ncell =0; ncell < nCell; ncell++)
1050         {
1051           if(flagCell[ncell] != 0) continue; // cell used before
1052           //calculate dr
1053           deta = etaCell[ncell] - eta;
1054           dphi = phiCell[ncell] - phi;
1055           if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
1056           if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1057           dr = TMath::Sqrt(deta * deta + dphi * dphi);
1058           if(dr <= rc){  // cell in cone
1059             flagCell[ncell] = -1;
1060             etCone+=etCell[ncell];
1061             nCellIn++;
1062           }
1063           else nCellOut++;
1064         }
1065       
1066       // Select jets with et > background
1067       // estimate max fluctuation of background in cone
1068       Double_t ncellin = (Double_t)nCellIn;
1069       Double_t ntcell  = (Double_t)nCell;
1070       Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/ntcell);
1071       // min cone et
1072       Double_t etcmin = etCone ;  // could be used etCone - etmin !!
1073       //decisions !! etbmax < etcmin
1074       
1075       for(Int_t mcell =0; mcell < nCell; mcell++){
1076         if(flagCell[mcell] == -1){
1077           if(etbmax < etcmin)
1078             flagCell[mcell] = 1; //flag cell as used
1079           else
1080             flagCell[mcell] = 0; // leave it free
1081         }
1082       }
1083       //store tmp jet info !!!
1084       
1085       if(etbmax < etcmin) {
1086         etaAlgoJet[nJets]    = eta;
1087         phiAlgoJet[nJets]    = phi;
1088         etAlgoJet[nJets]     = etCone;
1089         ncellsAlgoJet[nJets] = nCellIn;
1090         nJets++;
1091       }
1092       
1093     } // end of cells loop
1094
1095   //reorder jets by et in cone
1096   //sort jets by energy
1097   Int_t * idx  = new Int_t[nJets];
1098   TMath::Sort(nJets, etAlgoJet, idx);
1099   for(Int_t p = 0; p < nJets; p++)
1100     {
1101       etaJet[p] = etaAlgoJet[idx[p]];
1102       phiJet[p] = phiAlgoJet[idx[p]];
1103       etJet[p] = etAlgoJet[idx[p]];
1104       etallJet[p] = etAlgoJet[idx[p]];
1105       ncellsJet[p] = ncellsAlgoJet[idx[p]];
1106     }
1107   
1108   //delete
1109   delete index;
1110   delete idx;
1111
1112 }
1113
1114 ////////////////////////////////////////////////////////////////////////
1115 void AliUA1JetFinderV2::SubtractBackg(const Int_t& nIn, const Int_t&nJ, Float_t&etbgTotalN, const Float_t* ptT, 
1116                                       const Int_t*vectT, const Float_t* etaT, const Float_t* phiT, const Float_t* cFlagT, 
1117                                       const Float_t* cFlag2T, const Float_t* sFlagT, Float_t* etJet, const Float_t* etaJet, 
1118                                       const Float_t* phiJet, Float_t* etsigJet, Int_t* multJet, Int_t* injet)
1119 {
1120   //
1121   // Background subtraction using cone method but without correction in dE/deta distribution
1122   // Cases to take into account the EMCal geometry are included
1123   //
1124   
1125   //calculate energy inside and outside cones
1126   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1127   Int_t              fOpt = fReader->GetReaderHeader()->GetDetector();
1128   fDebug = fReader->GetReaderHeader()->GetDebug();
1129   Float_t rc= header->GetRadius();
1130   Float_t etIn[30];
1131   Float_t etOut = 0;
1132   
1133   for(Int_t j=0;j<30;j++){etIn[j]=0.;}
1134   
1135   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
1136
1137     for(Int_t ijet=0; ijet<nJ; ijet++){
1138
1139       Float_t deta = etaT[jpart] - etaJet[ijet];
1140       Float_t dphi = phiT[jpart] - phiJet[ijet];
1141       if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
1142       if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1143       
1144       Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
1145       if(dr <= rc){ // particles inside this cone
1146         multJet[ijet]+=vectT[jpart];
1147         injet[jpart] = ijet;
1148         
1149         if(cFlagT[jpart] == 1 || cFlag2T[jpart] == 1){ // pt cut
1150           etIn[ijet] += ptT[jpart];
1151           if(sFlagT[jpart] == 1) etsigJet[ijet]+= ptT[jpart];
1152         }
1153         break;
1154       }
1155     }// end jets loop
1156
1157     if(injet[jpart] == -1 && (cFlagT[jpart] == 1 || cFlag2T[jpart] == 1)){
1158       etOut += ptT[jpart]; // particle outside cones and pt cut
1159     }
1160   } //end particle loop
1161
1162   //estimate jets and background areas
1163   // TPC case
1164   if(fOpt == 0 || fOpt == 1){
1165     Float_t areaJet[30];
1166     Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
1167
1168     for(Int_t k=0; k<nJ; k++){
1169       Float_t detamax = etaJet[k] + rc;
1170       Float_t detamin = etaJet[k] - rc;
1171       Float_t accmax = 0.0; Float_t accmin = 0.0;
1172       if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
1173         Float_t h = header->GetLegoEtaMax() - etaJet[k];
1174         accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1175       }
1176       if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
1177         Float_t h = header->GetLegoEtaMax() + etaJet[k];
1178         accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1179       }
1180       areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
1181       areaOut = areaOut - areaJet[k];
1182     }
1183     //subtract background using area method
1184     for(Int_t ljet=0; ljet<nJ; ljet++){
1185       Float_t areaRatio = areaJet[ljet]/areaOut;
1186       etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
1187     }
1188
1189     // estimate new total background
1190     Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
1191     etbgTotalN = etOut*areaT/areaOut;
1192   }
1193   else { // If EMCal included
1194     Float_t areaJet[30];
1195     Float_t areaOut = 2*(header->GetLegoEtaMax())*(header->GetLegoPhiMax() - header->GetLegoPhiMin());
1196     for(Int_t k=0; k<nJ; k++){
1197       Float_t detamax = etaJet[k] + rc;
1198       Float_t detamin = etaJet[k] - rc;
1199       Float_t dphimax = phiJet[k] + rc;
1200       Float_t dphimin = phiJet[k] - rc;
1201       Float_t eMax = header->GetLegoEtaMax();
1202       Float_t eMin = header->GetLegoEtaMin();
1203       Float_t pMax = header->GetLegoPhiMax();
1204       Float_t pMin = header->GetLegoPhiMin();
1205       Float_t accetamax = 0.0; Float_t accetamin = 0.0;
1206       Float_t accphimax = 0.0; Float_t accphimin = 0.0;
1207       if((detamax > eMax && dphimax >= (pMin+2*rc) && dphimax <= pMax )||
1208          (detamax > eMax && dphimin <= (pMax-2*rc) && dphimin >= pMin )){
1209         Float_t h = eMax - etaJet[k];
1210         accetamax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1211       }
1212       if((detamin < eMin && dphimax >= (pMin+2*rc) && dphimax <= pMax )||
1213          (detamin < eMin && dphimin <= (pMax-2*rc) && dphimin >= pMin )){
1214         Float_t h = eMax + etaJet[k];
1215         accetamin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1216       }
1217       if((dphimax > pMax && detamax >= (eMin+2*rc) && detamax <= eMax )||
1218          (dphimax > pMax && detamin <= (eMax-2*rc) && detamin >= eMin )){
1219         Float_t h = pMax - phiJet[k];
1220         accphimax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1221       }
1222       if((dphimin < eMin && detamax >= (eMin+2*rc) && detamax <= eMax )||
1223          (dphimin < eMin && detamin <= (eMax-2*rc) && detamin >= eMin )){
1224         Float_t h = phiJet[k] - pMin;
1225         accphimin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1226       }
1227       
1228       if(detamax > eMax && dphimax > pMax ){
1229         Float_t he = eMax - etaJet[k];
1230         Float_t hp = pMax - phiJet[k];
1231         Float_t rlim = TMath::Sqrt(pow(he,2)+pow(hp,2));
1232         Float_t alphae = TMath::ACos(he/rc);
1233         Float_t alphap = TMath::ACos(hp/rc);
1234         Float_t alphad = (alphae+alphap)/2-TMath::Pi()/4;
1235         if(rlim <= rc){
1236           accetamax = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1237           accphimax = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp);
1238         }
1239         if(rlim > rc){
1240           accetamax = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1241           accphimax = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp)-
1242             ((TMath::Sqrt(pow(rc,2)-pow(he,2))-hp)*(TMath::Sqrt(pow(rc,2)-pow(hp,2))-he))/2+
1243             rc*rc*alphad - rc*rc*TMath::Sin(alphad)*TMath::Cos(alphad);
1244         }
1245       }
1246
1247       if(detamax > eMax && dphimin < pMin ){
1248         Float_t he = eMax - etaJet[k];
1249         Float_t hp = phiJet[k] - pMin;
1250         Float_t rlim = TMath::Sqrt(pow(he,2)+pow(hp,2));
1251         Float_t alphae = TMath::ACos(he/rc);
1252         Float_t alphap = TMath::ACos(hp/rc);
1253         Float_t alphad = (alphae+alphap)/2-TMath::Pi()/4;
1254         if(rlim <= rc){
1255           accetamax = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1256           accphimin = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp);
1257         }
1258         if(rlim > rc){
1259           accetamax = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1260           accphimin = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp)-
1261             ((TMath::Sqrt(pow(rc,2)-pow(he,2))-hp)*(TMath::Sqrt(pow(rc,2)-pow(hp,2))-he))/2+
1262             rc*rc*alphad - rc*rc*TMath::Sin(alphad)*TMath::Cos(alphad);
1263         }
1264       }
1265
1266       if(detamin < eMin && dphimax > pMax ){
1267         Float_t he = eMax + etaJet[k];
1268         Float_t hp = pMax - phiJet[k];
1269         Float_t rlim = TMath::Sqrt(pow(he,2)+pow(hp,2));
1270         Float_t alphae = TMath::ACos(he/rc);
1271         Float_t alphap = TMath::ACos(hp/rc);
1272         Float_t alphad = (alphae+alphap)/2-TMath::Pi()/4;
1273         if(rlim <= rc){
1274           accetamin = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1275           accphimax = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp);
1276         }
1277         if(rlim > rc){
1278           accetamin = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1279           accphimax = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp)-
1280             ((TMath::Sqrt(pow(rc,2)-pow(he,2))-hp)*(TMath::Sqrt(pow(rc,2)-pow(hp,2))-he))/2+
1281             rc*rc*alphad - rc*rc*TMath::Sin(alphad)*TMath::Cos(alphad);
1282         }
1283       }
1284
1285       if(detamin < eMin && dphimin < pMin ){
1286         Float_t he = eMax + etaJet[k];
1287         Float_t hp = phiJet[k] - pMin;
1288         Float_t rlim = TMath::Sqrt(pow(he,2)+pow(hp,2));
1289         Float_t alphae = TMath::ACos(he/rc);
1290         Float_t alphap = TMath::ACos(hp/rc);
1291         Float_t alphad = (alphae+alphap)/2-TMath::Pi()/4;
1292         if(rlim <= rc){
1293           accetamin = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1294           accphimin = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp);
1295         }
1296         if(rlim > rc){
1297           accetamin = rc*rc*alphae - he*TMath::Sqrt(rc*rc - he*he);
1298           accphimin = rc*rc*alphap - hp*TMath::Sqrt(rc*rc - hp*hp)-
1299             ((TMath::Sqrt(pow(rc,2)-pow(he,2))-hp)*(TMath::Sqrt(pow(rc,2)-pow(hp,2))-he))/2+
1300             rc*rc*alphad - rc*rc*TMath::Sin(alphad)*TMath::Cos(alphad);
1301         }
1302       }
1303       areaJet[k] = rc*rc*TMath::Pi() - accetamax - accetamin - accphimax - accphimin;
1304       areaOut = areaOut - areaJet[k];
1305     } // end loop on jets
1306    
1307     //subtract background using area method
1308     for(Int_t ljet=0; ljet<nJ; ljet++){
1309       Float_t areaRatio = areaJet[ljet]/areaOut;
1310       etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
1311     }
1312
1313     // estimate new total background
1314     Float_t areaT = 2*(header->GetLegoEtaMax()*header->GetLegoPhiMax());
1315     etbgTotalN = etOut*areaT/areaOut;
1316   }
1317   
1318 }
1319
1320 ////////////////////////////////////////////////////////////////////////
1321 void AliUA1JetFinderV2::SubtractBackgC(const Int_t& nIn, const Int_t&nJ, Float_t&etbgTotalN,
1322                                        const Float_t* ptT, const Float_t* etaT, const Float_t* phiT,
1323                                        Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
1324                                        Int_t* multJet, Int_t* injet)
1325 {
1326   //background subtraction using cone method but without correction in dE/deta distribution
1327   
1328   //calculate energy inside and outside cones
1329   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1330   Float_t rc= header->GetRadius();
1331   Float_t etIn[30];
1332   Float_t etOut = 0;
1333   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
1334     // if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
1335     for(Int_t ijet=0; ijet<nJ; ijet++){
1336       Float_t deta = etaT[jpart] - etaJet[ijet];
1337       Float_t dphi = phiT[jpart] - phiJet[ijet];
1338       if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
1339       if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1340       Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
1341       if(dr <= rc){ // particles inside this cone
1342         multJet[ijet]++;
1343         injet[jpart] = ijet;
1344         if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
1345           etIn[ijet] += ptT[jpart];
1346           if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet]+= ptT[jpart];
1347         }
1348         break;
1349       }
1350     }// end jets loop
1351     if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
1352       etOut += ptT[jpart]; // particle outside cones and pt cut
1353   } //end particle loop
1354   
1355   //estimate jets and background areas
1356   Float_t areaJet[30];
1357   Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
1358   for(Int_t k=0; k<nJ; k++){
1359     Float_t detamax = etaJet[k] + rc;
1360     Float_t detamin = etaJet[k] - rc;
1361     Float_t accmax = 0.0; Float_t accmin = 0.0;
1362     if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
1363       Float_t h = header->GetLegoEtaMax() - etaJet[k];
1364       accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1365     }
1366     if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
1367       Float_t h = header->GetLegoEtaMax() + etaJet[k];
1368       accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1369     }
1370     areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
1371     areaOut = areaOut - areaJet[k];
1372   }
1373   //subtract background using area method
1374   for(Int_t ljet=0; ljet<nJ; ljet++){
1375     Float_t areaRatio = areaJet[ljet]/areaOut;
1376     etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
1377   }
1378   
1379   // estimate new total background
1380   Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
1381   etbgTotalN = etOut*areaT/areaOut;
1382   
1383 }
1384
1385
1386 ////////////////////////////////////////////////////////////////////////
1387 void AliUA1JetFinderV2::SubtractBackgStat(const Int_t& nIn, const Int_t&nJ,Float_t&etbgTotalN,
1388                       const Float_t* ptT, const Float_t* etaT, const Float_t* phiT, const Float_t* cFlagT, 
1389                       const Float_t* sFlagT, Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, 
1390                       Float_t* etsigJet, Int_t* multJet, Int_t* injet)
1391 {
1392
1393   //background subtraction using statistical method
1394   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1395   Float_t etbgStat = header->GetBackgStat(); // pre-calculated background
1396   
1397   //calculate energy inside
1398   Float_t rc= header->GetRadius();
1399   Float_t etIn[30];
1400   
1401   for(Int_t jpart = 0; jpart < nIn; jpart++)
1402     { // loop for all particles in array
1403       
1404       for(Int_t ijet=0; ijet<nJ; ijet++)
1405         {
1406           Float_t deta = etaT[jpart] - etaJet[ijet];
1407           Float_t dphi = phiT[jpart] - phiJet[ijet];
1408           if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
1409           if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1410           Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
1411           if(dr <= rc){ // particles inside this cone
1412             multJet[ijet]++;
1413             injet[jpart] = ijet;
1414             if(cFlagT[jpart] == 1){ // pt cut
1415               etIn[ijet]+= ptT[jpart];
1416               if(sFlagT[jpart] == 1) etsigJet[ijet] += ptT[jpart];
1417             }
1418             break;
1419           }
1420         }// end jets loop
1421     } //end particle loop
1422   
1423   //calc jets areas
1424   Float_t areaJet[30];
1425   Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
1426   for(Int_t k=0; k<nJ; k++)
1427     {
1428       Float_t detamax = etaJet[k] + rc;
1429       Float_t detamin = etaJet[k] - rc;
1430       Float_t accmax = 0.0; Float_t accmin = 0.0;
1431       if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
1432         Float_t h = header->GetLegoEtaMax() - etaJet[k];
1433         accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1434       }
1435       if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
1436       Float_t h = header->GetLegoEtaMax() + etaJet[k];
1437       accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
1438       }
1439       areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
1440     }
1441
1442   //subtract background using area method
1443   for(Int_t ljet=0; ljet<nJ; ljet++){
1444     Float_t areaRatio = areaJet[ljet]/areaOut;
1445     etJet[ljet] = etIn[ljet]-etbgStat*areaRatio; // subtraction
1446   }
1447   
1448   etbgTotalN = etbgStat;
1449 }
1450
1451 ////////////////////////////////////////////////////////////////////////
1452 void AliUA1JetFinderV2::SubtractBackgCone(const Int_t& nIn, const Int_t&nJ,Float_t& etbgTotalN, Float_t* ptT,
1453                                           Float_t* etaT, const Float_t* phiT, const Float_t* cFlagT, const Float_t* sFlagT,
1454                                           Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
1455                                           Int_t* multJet, Int_t* injet)
1456 {
1457   // Cone background subtraction method taking into acount dEt/deta distribution
1458   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1459   //general
1460   Float_t rc= header->GetRadius();
1461   Float_t etamax = header->GetLegoEtaMax();
1462   Float_t etamin = header->GetLegoEtaMin();
1463   Int_t ndiv = 100;
1464   
1465   // jet energy and area arrays
1466   TH1F* hEtJet[30];
1467   TH1F* hAreaJet[30];
1468   for(Int_t mjet=0; mjet<nJ; mjet++){
1469     char hEtname[256]; char hAreaname[256];
1470     sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
1471     hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
1472     hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
1473   }
1474   // background energy and area
1475   TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
1476   TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
1477
1478   //fill energies
1479   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
1480     for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
1481       Float_t deta = etaT[jpart] - etaJet[ijet];
1482       Float_t dphi = phiT[jpart] - phiJet[ijet];
1483       if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
1484       if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1485       Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
1486       if(dr <= rc){ // particles inside this cone
1487         injet[jpart] = ijet;
1488         multJet[ijet]++;
1489         if(cFlagT[jpart] == 1){// pt cut
1490           hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone
1491           if(sFlagT[jpart] == 1) etsigJet[ijet] += ptT[jpart];
1492         }
1493         break;
1494       }
1495     }// end jets loop
1496
1497     if(injet[jpart] == -1  && cFlagT[jpart] == 1)
1498       hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
1499   } //end particle loop
1500
1501   //calc areas
1502   Float_t eta0 = etamin;
1503   Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
1504   Float_t eta1 = eta0 + etaw;
1505   for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
1506     Float_t etac = eta0 + etaw/2.0;
1507     Float_t areabg = etaw*2.0*TMath::Pi();
1508     for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
1509       Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
1510       Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
1511       Float_t acc0 = 0.0; Float_t acc1 = 0.0;
1512       Float_t areaj = 0.0;
1513       if(deta0 > rc && deta1 < rc){
1514         acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
1515         areaj = acc1;
1516       }
1517       if(deta0 < rc && deta1 > rc){
1518         acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
1519         areaj = acc0;
1520       }
1521       if(deta0 < rc && deta1 < rc){
1522         acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
1523         acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
1524         if(eta1<etaJet[ijet]) areaj = acc1-acc0;  // case 1
1525         if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
1526         if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
1527       }
1528       hAreaJet[ijet]->Fill(etac,areaj);
1529       areabg = areabg - areaj;
1530     } // end jets loop
1531     hAreaBackg->Fill(etac,areabg);
1532     eta0 = eta1;
1533     eta1 = eta1 + etaw;
1534   } // end loop for all eta bins
1535
1536   //subtract background
1537   for(Int_t kjet=0; kjet<nJ; kjet++){
1538     etJet[kjet] = 0.0; // first  clear etJet for this jet
1539     for(Int_t bin = 0; bin< ndiv; bin++){
1540       if(hAreaJet[kjet]->GetBinContent(bin)){
1541         Float_t areab = hAreaBackg->GetBinContent(bin);
1542         Float_t etb = hEtBackg->GetBinContent(bin);
1543         Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
1544         etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR); //subtraction
1545       }
1546     }
1547   }
1548
1549   // calc background total
1550   Double_t etOut = hEtBackg->Integral();
1551   Double_t areaOut = hAreaBackg->Integral();
1552   Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
1553   etbgTotalN = etOut*areaT/areaOut;
1554   
1555   //delete
1556   for(Int_t ljet=0; ljet<nJ; ljet++){  // loop for all jets
1557     delete hEtJet[ljet];
1558     delete hAreaJet[ljet];
1559   }
1560
1561   delete hEtBackg;
1562   delete hAreaBackg;
1563 }
1564
1565 ////////////////////////////////////////////////////////////////////////
1566 void AliUA1JetFinderV2::SubtractBackgRatio(const Int_t& nIn, const Int_t&nJ,Float_t& etbgTotalN,
1567                       Float_t* ptT, Float_t* etaT, const Float_t* phiT, const Float_t* cFlagT, const Float_t* sFlagT,
1568                       Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
1569                       Int_t* multJet, Int_t* injet)
1570 {
1571   // Ratio background subtraction method taking into acount dEt/deta distribution
1572   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1573   //factor F calc before
1574   Float_t bgRatioCut = header->GetBackgCutRatio();
1575   
1576   //general
1577   Float_t rc= header->GetRadius();
1578   Float_t etamax = header->GetLegoEtaMax();
1579   Float_t etamin = header->GetLegoEtaMin();
1580   Int_t ndiv = 100;
1581   
1582   // jet energy and area arrays
1583   TH1F* hEtJet[30];
1584   TH1F* hAreaJet[30];
1585   for(Int_t mjet=0; mjet<nJ; mjet++){
1586     char hEtname[256]; char hAreaname[256];
1587     sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
1588     hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);        // change range
1589     hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);  // change range
1590   }
1591   // background energy and area
1592   TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);         // change range
1593   TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);  // change range
1594
1595   //fill energies
1596   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
1597     for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
1598       Float_t deta = etaT[jpart] - etaJet[ijet];
1599       Float_t dphi = phiT[jpart] - phiJet[ijet];
1600       if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
1601       if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
1602       Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
1603       if(dr <= rc){ // particles inside this cone
1604         multJet[ijet]++;
1605         injet[jpart] = ijet;
1606         if(cFlagT[jpart] == 1){ //pt cut
1607           hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone and pt cut
1608           if(sFlagT[jpart] == 1) etsigJet[ijet] += ptT[jpart];
1609         }
1610         break;
1611       }
1612     }// end jets loop
1613     if(injet[jpart] == -1) hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
1614   } //end particle loop
1615
1616   //calc areas
1617   Float_t eta0 = etamin;
1618   Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
1619   Float_t eta1 = eta0 + etaw;
1620   for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
1621     Float_t etac = eta0 + etaw/2.0;
1622     Float_t areabg = etaw*2.0*TMath::Pi();
1623     for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
1624       Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
1625       Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
1626       Float_t acc0 = 0.0; Float_t acc1 = 0.0;
1627       Float_t areaj = 0.0;
1628       if(deta0 > rc && deta1 < rc){
1629         acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
1630         areaj = acc1;
1631       }
1632       if(deta0 < rc && deta1 > rc){
1633         acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
1634         areaj = acc0;
1635       }
1636       if(deta0 < rc && deta1 < rc){
1637         acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
1638         acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
1639         if(eta1<etaJet[ijet]) areaj = acc1-acc0;  // case 1
1640         if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
1641         if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
1642       }
1643       hAreaJet[ijet]->Fill(etac,areaj);
1644       areabg = areabg - areaj;
1645     } // end jets loop
1646     hAreaBackg->Fill(etac,areabg);
1647     eta0 = eta1;
1648     eta1 = eta1 + etaw;
1649   } // end loop for all eta bins
1650
1651   //subtract background
1652   for(Int_t kjet=0; kjet<nJ; kjet++){
1653     etJet[kjet] = 0.0; // first  clear etJet for this jet
1654     for(Int_t bin = 0; bin< ndiv; bin++){
1655       if(hAreaJet[kjet]->GetBinContent(bin)){
1656         Float_t areab = hAreaBackg->GetBinContent(bin);
1657         Float_t etb = hEtBackg->GetBinContent(bin);
1658         Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
1659         etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR*bgRatioCut); //subtraction
1660       }
1661     }
1662   }
1663
1664   // calc background total
1665   Double_t etOut = hEtBackg->Integral();
1666   Double_t areaOut = hAreaBackg->Integral();
1667   Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
1668   etbgTotalN = etOut*areaT/areaOut;
1669   
1670   //delete
1671   for(Int_t ljet=0; ljet<nJ; ljet++){  // loop for all jets
1672     delete hEtJet[ljet];
1673     delete hAreaJet[ljet];
1674   }
1675   
1676   delete hEtBackg;
1677   delete hAreaBackg;
1678 }
1679
1680 ////////////////////////////////////////////////////////////////////////
1681 void AliUA1JetFinderV2::Reset()
1682 {
1683   fLego->Reset();
1684   AliJetFinder::Reset();
1685 }
1686
1687 ////////////////////////////////////////////////////////////////////////
1688 void AliUA1JetFinderV2::WriteJHeaderToFile()
1689 {
1690   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1691   header->Write();
1692 }
1693
1694 ////////////////////////////////////////////////////////////////////////
1695 void AliUA1JetFinderV2::InitTask(TChain* tree)
1696 {
1697
1698   // initializes some variables
1699   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
1700   // book lego
1701   fLego = new TH2F("legoH","eta-phi",
1702                    header->GetLegoNbinEta(), header->GetLegoEtaMin(),
1703                    header->GetLegoEtaMax(),  header->GetLegoNbinPhi(),
1704                    header->GetLegoPhiMin(),  header->GetLegoPhiMax());
1705   
1706   fDebug = fReader->GetReaderHeader()->GetDebug();
1707   fOpt = fReader->GetReaderHeader()->GetDetector();
1708   
1709   // Tasks initialization
1710   if(fOpt>0)
1711     fReader->CreateTasks(tree);
1712
1713 }