]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliUA1JetFinderV1.cxx
Protection added.
[u/mrichter/AliRoot.git] / JETAN / AliUA1JetFinderV1.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
20 // manages the search for jets
21 // Author: Rafael.Diaz.Valdes@cern.ch
22 // (version in c++)
23 //---------------------------------------------------------------------
24
25 #include <TArrayF.h>
26 #include <TClonesArray.h>
27 #include <TFile.h>
28 #include <TH1F.h>
29 #include <TH2F.h>
30 #include <TLorentzVector.h>
31
32 #include "AliUA1JetFinderV1.h"
33 #include "AliUA1JetHeaderV1.h"
34 #include "AliJetReaderHeader.h"
35 #include "AliJetReader.h"
36 #include "AliJet.h"
37 #include "AliAODJet.h"
38
39
40 ClassImp(AliUA1JetFinderV1)
41
42 /////////////////////////////////////////////////////////////////////
43
44 AliUA1JetFinderV1::AliUA1JetFinderV1() :
45     AliJetFinder(),
46     fLego(0)
47 {
48   // Constructor
49 }
50
51 ////////////////////////////////////////////////////////////////////////
52
53 AliUA1JetFinderV1::~AliUA1JetFinderV1()
54
55 {
56   // destructor
57   delete fLego;
58   fLego = 0;
59 }
60
61 ////////////////////////////////////////////////////////////////////////
62
63
64 void AliUA1JetFinderV1::FindJets()
65
66 {
67   //1) Fill cell map array
68   //2) calculate total energy and fluctuation level
69   //3) Run algorithm
70   //   3.1) look centroides in cell map
71   //   3.2) calculate total energy in cones
72   //   3.3) flag as a possible jet
73   //   3.4) reorder cones by energy
74   //4) subtract backg in accepted jets
75   //5) fill AliJet list
76
77   // transform input to pt,eta,phi plus lego
78     
79   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
80   TClonesArray *lvArray = fReader->GetMomentumArray();
81   Int_t nIn =  lvArray->GetEntries();
82   if (nIn == 0) return;
83
84   // local arrays for input
85   Float_t* ptT   = new Float_t[nIn];
86   Float_t* etaT  = new Float_t[nIn];
87   Float_t* phiT  = new Float_t[nIn];
88   Int_t*   injet = new Int_t[nIn];
89
90   //total energy in array
91   Float_t  etbgTotal = 0.0;
92   TH1F* hPtTotal = new TH1F("hPt","Pt distribution of all particles ",100,0.0,15.0);
93
94   // load input vectors and calculate total energy in array
95   for (Int_t i = 0; i < nIn; i++){
96     TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
97     ptT[i]  = lv->Pt();
98     etaT[i] = lv->Eta();
99     phiT[i] = ((lv->Phi() < 0) ? (lv->Phi()) + 2 * TMath::Pi() : lv->Phi());
100     if (fReader->GetCutFlag(i) != 1) continue;
101     fLego ->Fill(etaT[i], phiT[i], ptT[i]);
102     hPtTotal->Fill(ptT[i]);
103     etbgTotal+= ptT[i];
104   }
105
106   fJets->SetNinput(nIn);
107
108   // calculate total energy and fluctuation in map
109   Double_t meanpt = hPtTotal->GetMean();
110   Double_t ptRMS  = hPtTotal->GetRMS();
111   Double_t npart  = hPtTotal->GetEntries();
112   Double_t dEtTotal = (TMath::Sqrt(npart))*TMath::Sqrt(meanpt * meanpt + ptRMS*ptRMS);
113
114   // arrays to hold jets
115   Float_t* etaJet = new Float_t[30];
116   Float_t* phiJet = new Float_t[30];
117   Float_t* etJet  = new Float_t[30];
118   Float_t* etsigJet  = new Float_t[30]; //signal et in jet
119   Float_t* etallJet = new Float_t[30];  // total et in jet (tmp variable)
120   Int_t* ncellsJet = new Int_t[30];
121   Int_t* multJet  = new Int_t[30];
122   Int_t nJets; // to hold number of jets found by algorithm
123   Int_t nj;    // number of jets accepted
124   Float_t prec  = header->GetPrecBg();
125   Float_t bgprec = 1;
126   while(bgprec > prec){
127      //reset jet arrays in memory
128      memset(etaJet,0,sizeof(Float_t)*30);
129      memset(phiJet,0,sizeof(Float_t)*30);
130      memset(etJet,0,sizeof(Float_t)*30);
131      memset(etallJet,0,sizeof(Float_t)*30);
132      memset(etsigJet,0,sizeof(Float_t)*30);
133      memset(ncellsJet,0,sizeof(Int_t)*30);
134      memset(multJet,0,sizeof(Int_t)*30);
135      nJets = 0;
136      nj = 0;
137      // reset particles-jet array in memory
138      memset(injet,-1,sizeof(Int_t)*nIn);
139      //run cone algorithm finder
140      RunAlgoritm(etbgTotal,dEtTotal,nJets,etJet,etaJet,phiJet,etallJet,ncellsJet);
141      //run background subtraction
142      if(nJets > header->GetNAcceptJets()) // limited number of accepted jets per event
143        nj = header->GetNAcceptJets();
144      else
145        nj = nJets;
146      //subtract background
147      Float_t etbgTotalN = 0.0; //new background
148      if(header->GetBackgMode() == 1) // standar
149         SubtractBackg(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
150      if(header->GetBackgMode() == 2) //cone
151         SubtractBackgCone(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
152      if(header->GetBackgMode() == 3) //ratio
153         SubtractBackgRatio(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
154      if(header->GetBackgMode() == 4) //statistic
155         SubtractBackgStat(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
156      //calc precision
157      if(etbgTotalN != 0.0)
158         bgprec = (etbgTotal - etbgTotalN)/etbgTotalN;
159      else
160         bgprec = 0;
161      etbgTotal = etbgTotalN; // update with new background estimation
162   } //end while
163
164   // add jets to list
165   Int_t* idxjets = new Int_t[nj];
166   Int_t nselectj = 0;
167 //  printf("Found %d jets \n", nj);
168   
169   TRefArray *refs = 0;
170   Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
171   if (fromAod) refs = fReader->GetReferences();
172   for(Int_t kj=0; kj<nj; kj++){
173      if ((etaJet[kj] > (header->GetJetEtaMax())) ||
174           (etaJet[kj] < (header->GetJetEtaMin())) ||
175           (etJet[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
176       Float_t px, py,pz,en; // convert to 4-vector
177       px = etJet[kj] * TMath::Cos(phiJet[kj]);
178       py = etJet[kj] * TMath::Sin(phiJet[kj]);
179       pz = etJet[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJet[kj])));
180       en = TMath::Sqrt(px * px + py * py + pz * pz);
181       fJets->AddJet(px, py, pz, en);
182       AliAODJet jet(px, py, pz, en);
183
184       if (fromAod){
185         for(Int_t jpart = 0; jpart < nIn; jpart++) // loop for all particles in array
186           if (injet[jpart] == kj && fReader->GetCutFlag(jpart) == 1)
187                     jet.AddTrack(refs->At(jpart));  // check if the particle belongs to the jet and add the ref
188       }
189       
190       //jet.Print("");
191       
192       AddJet(jet);
193       
194       idxjets[nselectj] = kj;
195       nselectj++;
196   } //end particle loop
197
198   //add signal percentage and total signal  in AliJets for analysis tool
199   Float_t* percentage  = new Float_t[nselectj];
200   Int_t* ncells      = new Int_t[nselectj];
201   Int_t* mult        = new Int_t[nselectj];
202   for(Int_t i = 0; i< nselectj; i++){
203      percentage[i] = etsigJet[idxjets[i]]/etJet[idxjets[i]];
204      ncells[i] = ncellsJet[idxjets[i]];
205      mult[i] = multJet[idxjets[i]];
206   }
207    //add particle-injet relationship ///
208    for(Int_t bj = 0; bj < nIn; bj++){
209        if(injet[bj] == -1) continue; //background particle
210        Int_t bflag = 0;
211        for(Int_t ci = 0; ci< nselectj; ci++){
212            if(injet[bj] == idxjets[ci]){
213               injet[bj]= ci;
214               bflag++;
215               break;
216            }
217        }
218        if(bflag == 0) injet[bj] = -1; // set as background particle
219    }
220   fJets->SetNCells(ncells);
221   fJets->SetPtFromSignal(percentage);
222   fJets->SetMultiplicities(mult);
223   fJets->SetInJet(injet);
224   fJets->SetEtaIn(etaT);
225   fJets->SetPhiIn(phiT);
226   fJets->SetPtIn(ptT);
227   fJets->SetEtAvg(etbgTotal/(4*(header->GetLegoEtaMax())*TMath::Pi()));
228
229
230   //delete
231   delete [] ptT;
232   delete [] etaT;
233   delete [] phiT;
234   delete [] injet;
235   delete hPtTotal;
236   delete [] etaJet;
237   delete [] phiJet;
238   delete [] etJet;
239   delete [] etsigJet;
240   delete [] etallJet;
241   delete [] ncellsJet;
242   delete [] multJet;
243   delete [] idxjets;
244   delete [] percentage;
245   delete [] ncells;
246   delete [] mult;
247
248
249 }
250
251 ////////////////////////////////////////////////////////////////////////
252
253 void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t& nJets,
254                                   Float_t* etJet,Float_t* etaJet, Float_t* phiJet,
255                                   Float_t* etallJet, Int_t* ncellsJet)
256 {
257
258    //dump lego
259   // check enough space! *to be done*
260   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
261   const Int_t nBinsMax = 70000;
262   
263   Float_t etCell[nBinsMax];   //! Cell Energy
264   Float_t etaCell[nBinsMax];  //! Cell eta
265   Float_t phiCell[nBinsMax];  //! Cell phi
266   Int_t   flagCell[nBinsMax]; //! Cell flag
267
268   Int_t nCell = 0;
269   TAxis* xaxis = fLego->GetXaxis();
270   TAxis* yaxis = fLego->GetYaxis();
271   Float_t e = 0.0;
272   for (Int_t i = 1; i <= header->GetLegoNbinEta(); i++) {
273       for (Int_t j = 1; j <= header->GetLegoNbinPhi(); j++) {
274                e = fLego->GetBinContent(i,j);
275                if (e < 0.0) continue; // don't include this cells
276                Float_t eta  = xaxis->GetBinCenter(i);
277                Float_t phi  = yaxis->GetBinCenter(j);
278                etCell[nCell]  = e;
279                etaCell[nCell] = eta;
280                phiCell[nCell] = phi;
281                flagCell[nCell] = 0; //default
282                nCell++;
283       }
284   }
285
286   // Parameters from header
287   Float_t minmove = header->GetMinMove();
288   Float_t maxmove = header->GetMaxMove();
289   Float_t rc      = header->GetRadius();
290   Float_t etseed  = header->GetEtSeed();
291   //Float_t etmin   = header->GetMinJetEt();
292
293
294
295   // tmp array of jets form algoritm
296   Float_t etaAlgoJet[30];
297   Float_t phiAlgoJet[30];
298   Float_t etAlgoJet[30];
299   Int_t   ncellsAlgoJet[30];
300
301   //run algorithm//
302
303   // sort cells by et
304   Int_t * index  = new Int_t[nCell];
305   TMath::Sort(nCell, etCell, index);
306   // variable used in centroide loop
307   Float_t eta   = 0.0;
308   Float_t phi   = 0.0;
309   Float_t eta0  = 0.0;
310   Float_t phi0  = 0.0;
311   Float_t etab  = 0.0;
312   Float_t phib  = 0.0;
313   Float_t etas  = 0.0;
314   Float_t phis  = 0.0;
315   Float_t ets   = 0.0;
316   Float_t deta  = 0.0;
317   Float_t dphi  = 0.0;
318   Float_t dr    = 0.0;
319   Float_t etsb  = 0.0;
320   Float_t etasb = 0.0;
321   Float_t phisb = 0.0;
322   Float_t dphib = 0.0;
323   
324
325   for(Int_t icell = 0; icell < nCell; icell++){
326         Int_t jcell = index[icell];
327         if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
328         if(flagCell[jcell] != 0) continue; // if cell was used before
329         eta  = etaCell[jcell];
330         phi  = phiCell[jcell];
331         eta0 = eta;
332         phi0 = phi;
333         etab = eta;
334         phib = phi;
335         ets  = etCell[jcell];
336         etas = 0.0;
337         phis = 0.0;
338         etsb = ets;
339         etasb = 0.0;
340         phisb = 0.0;
341         for(Int_t kcell =0; kcell < nCell; kcell++){
342             Int_t lcell = index[kcell];
343             if(lcell == jcell)                 continue; // cell itself
344             if(flagCell[lcell] != 0)           continue; // cell used before
345             if(etCell[lcell] > etCell[jcell])  continue; // can this happen
346             //calculate dr
347             deta = etaCell[lcell] - eta;
348             dphi = TMath::Abs(phiCell[lcell] - phi);
349             if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
350             dr = TMath::Sqrt(deta * deta + dphi * dphi);
351             if(dr <= rc){
352                // calculate offset from initiate cell
353                deta = etaCell[lcell] - eta0;
354                dphi = phiCell[lcell] - phi0;
355                if (dphi < - TMath::Pi()) dphi=  dphi + 2.0 * TMath::Pi();
356                if (dphi >   TMath::Pi()) dphi = dphi - 2.0 * TMath::Pi();
357                
358                etas = etas + etCell[lcell]*deta;
359                phis = phis + etCell[lcell]*dphi;
360                ets = ets + etCell[lcell];
361                //new weighted eta and phi including this cell
362                eta = eta0 + etas/ets;
363                phi = phi0 + phis/ets;
364                // if cone does not move much, just go to next step
365                dphib = TMath::Abs(phi - phib);
366                if (dphib > TMath::Pi()) dphib = 2. * TMath::Pi() - dphib;
367                dr = TMath::Sqrt((eta-etab)*(eta-etab) + dphib * dphib);
368                if(dr <= minmove) break;
369                // cone should not move more than max_mov
370                dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
371                if(dr > maxmove){
372                    eta = etab;
373                    phi = phib;
374                    ets = etsb;
375                    etas = etasb;
376                    phis = phisb;
377                } else { // store this loop information
378                    etab  = eta;
379                    phib  = phi;
380                    etsb  = ets;
381                    etasb = etas;
382                    phisb = phis;
383                }
384             } // inside cone
385         }//end of cells loop looking centroide
386
387         //avoid cones overloap (to be implemented in the future)
388
389         //flag cells in Rc, estimate total energy in cone
390         Float_t etCone = 0.0;
391         Int_t   nCellIn = 0;
392         rc = header->GetRadius();
393         for(Int_t ncell =0; ncell < nCell; ncell++){
394             if(flagCell[ncell] != 0) continue; // cell used before
395            //calculate dr
396             deta = etaCell[ncell] - eta;
397                  dphi = phiCell[ncell] - phi;
398                  if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
399                  if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
400                  dr = TMath::Sqrt(deta * deta + dphi * dphi);
401             if(dr <= rc){  // cell in cone
402                flagCell[ncell] = -1;
403                etCone+=etCell[ncell];
404                nCellIn++;
405             }
406         }
407
408         // select jets with et > background
409         // estimate max fluctuation of background in cone
410         Double_t ncellin = (Double_t)nCellIn;
411         Double_t ntcell  = (Double_t)nCell;
412         Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/ntcell);
413         // min cone et
414         Double_t etcmin = etCone ;  // could be used etCone - etmin !!
415         //desicions !! etbmax < etcmin
416         for(Int_t mcell =0; mcell < nCell; mcell++){
417             if(flagCell[mcell] == -1){
418               if(etbmax < etcmin)
419                  flagCell[mcell] = 1; //flag cell as used
420               else
421                  flagCell[mcell] = 0; // leave it free
422             }
423         }
424         //store tmp jet info !!!
425        if(etbmax < etcmin) {
426              etaAlgoJet[nJets] = eta;
427              phiAlgoJet[nJets] = phi;
428              etAlgoJet[nJets] = etCone;
429              ncellsAlgoJet[nJets] = nCellIn;
430              nJets++;
431         }
432
433   } // end of cells loop
434
435   //reorder jets by et in cone
436   //sort jets by energy
437   Int_t * idx  = new Int_t[nJets];
438   TMath::Sort(nJets, etAlgoJet, idx);
439   for(Int_t p = 0; p < nJets; p++){
440      etaJet[p] = etaAlgoJet[idx[p]];
441      phiJet[p] = phiAlgoJet[idx[p]];
442      etJet[p] = etAlgoJet[idx[p]];
443      etallJet[p] = etAlgoJet[idx[p]];
444      ncellsJet[p] = ncellsAlgoJet[idx[p]];
445   }
446
447
448   //delete
449   delete[] index;
450   delete[] idx;
451
452 }
453 ////////////////////////////////////////////////////////////////////////
454
455 void AliUA1JetFinderV1::SubtractBackg(Int_t& nIn, Int_t&nJ, Float_t&etbgTotalN,
456                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
457                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
458                       Int_t* multJet, Int_t* injet)
459 {
460   //background subtraction using cone method but without correction in dE/deta distribution
461
462   //calculate energy inside and outside cones
463   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
464   Float_t rc= header->GetRadius();
465   Float_t etIn[30];
466   Float_t etOut = 0;
467   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
468      // if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
469      for(Int_t ijet=0; ijet<nJ; ijet++){
470          Float_t deta = etaT[jpart] - etaJet[ijet];
471               Float_t dphi = phiT[jpart] - phiJet[ijet];
472          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
473               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
474               Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
475          if(dr <= rc){ // particles inside this cone
476              multJet[ijet]++;
477              injet[jpart] = ijet;
478              if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
479                 etIn[ijet] += ptT[jpart];
480                 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet]+= ptT[jpart];
481              }
482              break;
483          }
484      }// end jets loop
485      if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
486         etOut += ptT[jpart]; // particle outside cones and pt cut
487   } //end particle loop
488
489   //estimate jets and background areas
490   Float_t areaJet[30];
491   Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
492   for(Int_t k=0; k<nJ; k++){
493       Float_t detamax = etaJet[k] + rc;
494       Float_t detamin = etaJet[k] - rc;
495       Float_t accmax = 0.0; Float_t accmin = 0.0;
496       if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
497          Float_t h = header->GetLegoEtaMax() - etaJet[k];
498          accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
499       }
500       if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
501          Float_t h = header->GetLegoEtaMax() + etaJet[k];
502          accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
503       }
504       areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
505       areaOut = areaOut - areaJet[k];
506   }
507   //subtract background using area method
508   for(Int_t ljet=0; ljet<nJ; ljet++){
509      Float_t areaRatio = areaJet[ljet]/areaOut;
510      etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
511   }
512
513   // estimate new total background
514   Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
515   etbgTotalN = etOut*areaT/areaOut;
516
517
518 }
519
520 ////////////////////////////////////////////////////////////////////////
521
522 void AliUA1JetFinderV1::SubtractBackgStat(Int_t& nIn, Int_t&nJ,Float_t&etbgTotalN,
523                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
524                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
525                       Int_t* multJet, Int_t* injet)
526 {
527
528   //background subtraction using statistical method
529   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
530   Float_t etbgStat = header->GetBackgStat(); // pre-calculated background
531
532   //calculate energy inside
533   Float_t rc= header->GetRadius();
534   Float_t etIn[30];
535
536   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
537      //if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
538      for(Int_t ijet=0; ijet<nJ; ijet++){
539          Float_t deta = etaT[jpart] - etaJet[ijet];
540               Float_t dphi = phiT[jpart] - phiJet[ijet];
541          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
542               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
543               Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
544          if(dr <= rc){ // particles inside this cone
545              multJet[ijet]++;
546              injet[jpart] = ijet;
547              if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
548                 etIn[ijet]+= ptT[jpart];
549                 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
550              }
551              break;
552          }
553      }// end jets loop
554   } //end particle loop
555
556   //calc jets areas
557   Float_t areaJet[30];
558   Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
559   for(Int_t k=0; k<nJ; k++){
560       Float_t detamax = etaJet[k] + rc;
561       Float_t detamin = etaJet[k] - rc;
562       Float_t accmax = 0.0; Float_t accmin = 0.0;
563       if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
564          Float_t h = header->GetLegoEtaMax() - etaJet[k];
565          accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
566       }
567       if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
568          Float_t h = header->GetLegoEtaMax() + etaJet[k];
569          accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
570       }
571       areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
572   }
573
574   //subtract background using area method
575   for(Int_t ljet=0; ljet<nJ; ljet++){
576      Float_t areaRatio = areaJet[ljet]/areaOut;
577      etJet[ljet] = etIn[ljet]-etbgStat*areaRatio; // subtraction
578   }
579
580   etbgTotalN = etbgStat;
581
582 }
583
584 ////////////////////////////////////////////////////////////////////////
585
586 void AliUA1JetFinderV1::SubtractBackgCone(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
587                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
588                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
589                       Int_t* multJet, Int_t* injet)
590 {
591    // Cone background subtraction method taking into acount dEt/deta distribution
592     AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
593    //general
594    Float_t rc= header->GetRadius();
595    Float_t etamax = header->GetLegoEtaMax();
596    Float_t etamin = header->GetLegoEtaMin();
597    Int_t ndiv = 100;
598
599    // jet energy and area arrays
600    TH1F* hEtJet[30];
601    TH1F* hAreaJet[30];
602    for(Int_t mjet=0; mjet<nJ; mjet++){
603      char hEtname[256]; char hAreaname[256];
604      sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
605      hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
606      hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
607   }
608    // background energy and area
609    TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
610    TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
611
612    //fill energies
613    for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
614      for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
615          Float_t deta = etaT[jpart] - etaJet[ijet];
616          Float_t dphi = phiT[jpart] - phiJet[ijet];
617          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
618          if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
619          Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
620          if(dr <= rc){ // particles inside this cone
621              injet[jpart] = ijet;
622              multJet[ijet]++;
623              if((fReader->GetCutFlag(jpart)) == 1){// pt cut
624                 hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone
625                 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
626              }
627              break;
628          }
629      }// end jets loop
630      if(injet[jpart] == -1  && fReader->GetCutFlag(jpart) == 1)
631         hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
632   } //end particle loop
633
634    //calc areas
635    Float_t eta0 = etamin;
636    Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
637    Float_t eta1 = eta0 + etaw;
638    for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
639       Float_t etac = eta0 + etaw/2.0;
640       Float_t areabg = etaw*2.0*TMath::Pi();
641       for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
642           Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
643           Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
644           Float_t acc0 = 0.0; Float_t acc1 = 0.0;
645           Float_t areaj = 0.0;
646           if(deta0 > rc && deta1 < rc){
647              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
648              areaj = acc1;
649           }
650           if(deta0 < rc && deta1 > rc){
651              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
652              areaj = acc0;
653           }
654           if(deta0 < rc && deta1 < rc){
655              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
656              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
657              if(eta1<etaJet[ijet]) areaj = acc1-acc0;  // case 1
658              if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
659              if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
660           }
661           hAreaJet[ijet]->Fill(etac,areaj);
662           areabg = areabg - areaj;
663       } // end jets loop
664       hAreaBackg->Fill(etac,areabg);
665       eta0 = eta1;
666       eta1 = eta1 + etaw;
667    } // end loop for all eta bins
668
669    //subtract background
670    for(Int_t kjet=0; kjet<nJ; kjet++){
671        etJet[kjet] = 0.0; // first  clear etJet for this jet
672        for(Int_t bin = 0; bin< ndiv; bin++){
673            if(hAreaJet[kjet]->GetBinContent(bin)){
674               Float_t areab = hAreaBackg->GetBinContent(bin);
675               Float_t etb = hEtBackg->GetBinContent(bin);
676               Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
677               etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR); //subtraction
678            }
679        }
680    }
681
682    // calc background total
683    Double_t etOut = hEtBackg->Integral();
684    Double_t areaOut = hAreaBackg->Integral();
685    Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
686    etbgTotalN = etOut*areaT/areaOut;
687
688    //delete
689    for(Int_t ljet=0; ljet<nJ; ljet++){  // loop for all jets
690        delete hEtJet[ljet];
691        delete hAreaJet[ljet];
692    }
693
694    delete hEtBackg;
695    delete hAreaBackg;
696 }
697
698 ////////////////////////////////////////////////////////////////////////
699
700
701 void AliUA1JetFinderV1::SubtractBackgRatio(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
702                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
703                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
704                        Int_t* multJet, Int_t* injet)
705 {
706    // Ratio background subtraction method taking into acount dEt/deta distribution
707     AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
708    //factor F calc before
709     Float_t bgRatioCut = header->GetBackgCutRatio();
710
711
712    //general
713    Float_t rc= header->GetRadius();
714    Float_t etamax = header->GetLegoEtaMax();
715    Float_t etamin = header->GetLegoEtaMin();
716    Int_t ndiv = 100;
717
718    // jet energy and area arrays
719    TH1F* hEtJet[30];
720    TH1F* hAreaJet[30];
721    for(Int_t mjet=0; mjet<nJ; mjet++){
722      char hEtname[256]; char hAreaname[256];
723      sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
724      hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);        // change range
725      hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);  // change range
726   }
727    // background energy and area
728    TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);         // change range
729    TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);  // change range
730
731    //fill energies
732    for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
733      //if((fReader->GetCutFlag(jpart)) != 1) continue;
734      for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
735          Float_t deta = etaT[jpart] - etaJet[ijet];
736               Float_t dphi = phiT[jpart] - phiJet[ijet];
737          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
738               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
739               Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
740          if(dr <= rc){ // particles inside this cone
741             multJet[ijet]++;
742             injet[jpart] = ijet;
743             if((fReader->GetCutFlag(jpart)) == 1){ //pt cut
744                hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone and pt cut
745                if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
746             }
747             break;
748          }
749      }// end jets loop
750      if(injet[jpart] == -1) hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
751   } //end particle loop
752
753    //calc areas
754    Float_t eta0 = etamin;
755    Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
756    Float_t eta1 = eta0 + etaw;
757    for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
758       Float_t etac = eta0 + etaw/2.0;
759       Float_t areabg = etaw*2.0*TMath::Pi();
760       for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
761           Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
762           Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
763           Float_t acc0 = 0.0; Float_t acc1 = 0.0;
764           Float_t areaj = 0.0;
765           if(deta0 > rc && deta1 < rc){
766              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
767              areaj = acc1;
768           }
769           if(deta0 < rc && deta1 > rc){
770              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
771              areaj = acc0;
772           }
773           if(deta0 < rc && deta1 < rc){
774              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
775              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
776              if(eta1<etaJet[ijet]) areaj = acc1-acc0;  // case 1
777              if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
778              if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
779           }
780           hAreaJet[ijet]->Fill(etac,areaj);
781           areabg = areabg - areaj;
782       } // end jets loop
783       hAreaBackg->Fill(etac,areabg);
784       eta0 = eta1;
785       eta1 = eta1 + etaw;
786    } // end loop for all eta bins
787
788    //subtract background
789    for(Int_t kjet=0; kjet<nJ; kjet++){
790        etJet[kjet] = 0.0; // first  clear etJet for this jet
791        for(Int_t bin = 0; bin< ndiv; bin++){
792            if(hAreaJet[kjet]->GetBinContent(bin)){
793               Float_t areab = hAreaBackg->GetBinContent(bin);
794               Float_t etb = hEtBackg->GetBinContent(bin);
795               Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
796               etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR*bgRatioCut); //subtraction
797            }
798        }
799    }
800
801    // calc background total
802    Double_t etOut = hEtBackg->Integral();
803    Double_t areaOut = hAreaBackg->Integral();
804    Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
805    etbgTotalN = etOut*areaT/areaOut;
806
807    //delete
808    for(Int_t ljet=0; ljet<nJ; ljet++){  // loop for all jets
809        delete hEtJet[ljet];
810        delete hAreaJet[ljet];
811    }
812
813    delete hEtBackg;
814    delete hAreaBackg;
815 }
816
817 ////////////////////////////////////////////////////////////////////////
818
819
820 void AliUA1JetFinderV1::Reset()
821 {
822   fLego->Reset();
823   fJets->ClearJets();
824   AliJetFinder::Reset();
825 }
826
827 ////////////////////////////////////////////////////////////////////////
828
829 void AliUA1JetFinderV1::WriteJHeaderToFile()
830 {
831   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
832   header->Write();
833 }
834
835 ////////////////////////////////////////////////////////////////////////
836
837 void AliUA1JetFinderV1::Init()
838 {
839   // initializes some variables
840   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
841   fLego = new
842     TH2F("legoH","eta-phi",
843          header->GetLegoNbinEta(), header->GetLegoEtaMin(),
844          header->GetLegoEtaMax(),  header->GetLegoNbinPhi(),
845          header->GetLegoPhiMin(),  header->GetLegoPhiMax());
846   // Do not store in current dir
847   fLego->SetDirectory(0);
848
849 }