]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliUA1JetFinderV1.cxx
Possibility to instantiate CDB manager from prefetched cache. It is used on the proof...
[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   for(Int_t kj=0; kj<nj; kj++){
170      if ((etaJet[kj] > (header->GetJetEtaMax())) ||
171           (etaJet[kj] < (header->GetJetEtaMin())) ||
172           (etJet[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
173       Float_t px, py,pz,en; // convert to 4-vector
174       px = etJet[kj] * TMath::Cos(phiJet[kj]);
175       py = etJet[kj] * TMath::Sin(phiJet[kj]);
176       pz = etJet[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJet[kj])));
177       en = TMath::Sqrt(px * px + py * py + pz * pz);
178       fJets->AddJet(px, py, pz, en);
179       AliAODJet jet(px, py, pz, en);
180       //jet.Print("");
181       
182       AddJet(jet);
183       
184       idxjets[nselectj] = kj;
185       nselectj++;
186   }
187   //add signal percentage and total signal  in AliJets for analysis tool
188   Float_t* percentage  = new Float_t[nselectj];
189   Int_t* ncells      = new Int_t[nselectj];
190   Int_t* mult        = new Int_t[nselectj];
191   for(Int_t i = 0; i< nselectj; i++){
192      percentage[i] = etsigJet[idxjets[i]]/etJet[idxjets[i]];
193      ncells[i] = ncellsJet[idxjets[i]];
194      mult[i] = multJet[idxjets[i]];
195   }
196    //add particle-injet relationship ///
197    for(Int_t bj = 0; bj < nIn; bj++){
198        if(injet[bj] == -1) continue; //background particle
199        Int_t bflag = 0;
200        for(Int_t ci = 0; ci< nselectj; ci++){
201            if(injet[bj] == idxjets[ci]){
202               injet[bj]= ci;
203               bflag++;
204               break;
205            }
206        }
207        if(bflag == 0) injet[bj] = -1; // set as background particle
208    }
209   fJets->SetNCells(ncells);
210   fJets->SetPtFromSignal(percentage);
211   fJets->SetMultiplicities(mult);
212   fJets->SetInJet(injet);
213   fJets->SetEtaIn(etaT);
214   fJets->SetPhiIn(phiT);
215   fJets->SetPtIn(ptT);
216   fJets->SetEtAvg(etbgTotal/(4*(header->GetLegoEtaMax())*TMath::Pi()));
217
218
219   //delete
220   delete [] ptT;
221   delete [] etaT;
222   delete [] phiT;
223   delete [] injet;
224   delete hPtTotal;
225   delete [] etaJet;
226   delete [] phiJet;
227   delete [] etJet;
228   delete [] etsigJet;
229   delete [] etallJet;
230   delete [] ncellsJet;
231   delete [] multJet;
232   delete [] idxjets;
233   delete [] percentage;
234   delete [] ncells;
235   delete [] mult;
236
237
238 }
239
240 ////////////////////////////////////////////////////////////////////////
241
242 void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t& nJets,
243                                   Float_t* etJet,Float_t* etaJet, Float_t* phiJet,
244                                   Float_t* etallJet, Int_t* ncellsJet)
245 {
246
247    //dump lego
248   // check enough space! *to be done*
249   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
250   const Int_t nBinsMax = 70000;
251   
252   Float_t etCell[nBinsMax];   //! Cell Energy
253   Float_t etaCell[nBinsMax];  //! Cell eta
254   Float_t phiCell[nBinsMax];  //! Cell phi
255   Int_t   flagCell[nBinsMax]; //! Cell flag
256
257   Int_t nCell = 0;
258   TAxis* xaxis = fLego->GetXaxis();
259   TAxis* yaxis = fLego->GetYaxis();
260   Float_t e = 0.0;
261   for (Int_t i = 1; i <= header->GetLegoNbinEta(); i++) {
262       for (Int_t j = 1; j <= header->GetLegoNbinPhi(); j++) {
263                e = fLego->GetBinContent(i,j);
264                if (e < 0.0) continue; // don't include this cells
265                Float_t eta  = xaxis->GetBinCenter(i);
266                Float_t phi  = yaxis->GetBinCenter(j);
267                etCell[nCell]  = e;
268                etaCell[nCell] = eta;
269                phiCell[nCell] = phi;
270                flagCell[nCell] = 0; //default
271                nCell++;
272       }
273   }
274
275   // Parameters from header
276   Float_t minmove = header->GetMinMove();
277   Float_t maxmove = header->GetMaxMove();
278   Float_t rc      = header->GetRadius();
279   Float_t etseed  = header->GetEtSeed();
280   //Float_t etmin   = header->GetMinJetEt();
281
282
283
284   // tmp array of jets form algoritm
285   Float_t etaAlgoJet[30];
286   Float_t phiAlgoJet[30];
287   Float_t etAlgoJet[30];
288   Int_t   ncellsAlgoJet[30];
289
290   //run algorithm//
291
292   // sort cells by et
293   Int_t * index  = new Int_t[nCell];
294   TMath::Sort(nCell, etCell, index);
295   // variable used in centroide loop
296   Float_t eta   = 0.0;
297   Float_t phi   = 0.0;
298   Float_t eta0  = 0.0;
299   Float_t phi0  = 0.0;
300   Float_t etab  = 0.0;
301   Float_t phib  = 0.0;
302   Float_t etas  = 0.0;
303   Float_t phis  = 0.0;
304   Float_t ets   = 0.0;
305   Float_t deta  = 0.0;
306   Float_t dphi  = 0.0;
307   Float_t dr    = 0.0;
308   Float_t etsb  = 0.0;
309   Float_t etasb = 0.0;
310   Float_t phisb = 0.0;
311   Float_t dphib = 0.0;
312   
313
314   for(Int_t icell = 0; icell < nCell; icell++){
315         Int_t jcell = index[icell];
316         if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
317         if(flagCell[jcell] != 0) continue; // if cell was used before
318         eta  = etaCell[jcell];
319         phi  = phiCell[jcell];
320         eta0 = eta;
321         phi0 = phi;
322         etab = eta;
323         phib = phi;
324         ets  = etCell[jcell];
325         etas = 0.0;
326         phis = 0.0;
327         etsb = ets;
328         etasb = 0.0;
329         phisb = 0.0;
330         for(Int_t kcell =0; kcell < nCell; kcell++){
331             Int_t lcell = index[kcell];
332             if(lcell == jcell)                 continue; // cell itself
333             if(flagCell[lcell] != 0)           continue; // cell used before
334             if(etCell[lcell] > etCell[jcell])  continue; // can this happen
335             //calculate dr
336             deta = etaCell[lcell] - eta;
337             dphi = TMath::Abs(phiCell[lcell] - phi);
338             if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
339             dr = TMath::Sqrt(deta * deta + dphi * dphi);
340             if(dr <= rc){
341                // calculate offset from initiate cell
342                deta = etaCell[lcell] - eta0;
343                dphi = phiCell[lcell] - phi0;
344                if (dphi < - TMath::Pi()) dphi=  dphi + 2.0 * TMath::Pi();
345                if (dphi >   TMath::Pi()) dphi = dphi - 2.0 * TMath::Pi();
346                
347                etas = etas + etCell[lcell]*deta;
348                phis = phis + etCell[lcell]*dphi;
349                ets = ets + etCell[lcell];
350                //new weighted eta and phi including this cell
351                eta = eta0 + etas/ets;
352                phi = phi0 + phis/ets;
353                // if cone does not move much, just go to next step
354                dphib = TMath::Abs(phi - phib);
355                if (dphib > TMath::Pi()) dphib = 2. * TMath::Pi() - dphib;
356                dr = TMath::Sqrt((eta-etab)*(eta-etab) + dphib * dphib);
357                if(dr <= minmove) break;
358                // cone should not move more than max_mov
359                dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
360                if(dr > maxmove){
361                    eta = etab;
362                    phi = phib;
363                    ets = etsb;
364                    etas = etasb;
365                    phis = phisb;
366                } else { // store this loop information
367                    etab  = eta;
368                    phib  = phi;
369                    etsb  = ets;
370                    etasb = etas;
371                    phisb = phis;
372                }
373             } // inside cone
374         }//end of cells loop looking centroide
375
376         //avoid cones overloap (to be implemented in the future)
377
378         //flag cells in Rc, estimate total energy in cone
379         Float_t etCone = 0.0;
380         Int_t   nCellIn = 0;
381         rc = header->GetRadius();
382         for(Int_t ncell =0; ncell < nCell; ncell++){
383             if(flagCell[ncell] != 0) continue; // cell used before
384            //calculate dr
385             deta = etaCell[ncell] - eta;
386                  dphi = phiCell[ncell] - phi;
387                  if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
388                  if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
389                  dr = TMath::Sqrt(deta * deta + dphi * dphi);
390             if(dr <= rc){  // cell in cone
391                flagCell[ncell] = -1;
392                etCone+=etCell[ncell];
393                nCellIn++;
394             }
395         }
396
397         // select jets with et > background
398         // estimate max fluctuation of background in cone
399         Double_t ncellin = (Double_t)nCellIn;
400         Double_t ntcell  = (Double_t)nCell;
401         Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/ntcell);
402         // min cone et
403         Double_t etcmin = etCone ;  // could be used etCone - etmin !!
404         //desicions !! etbmax < etcmin
405         for(Int_t mcell =0; mcell < nCell; mcell++){
406             if(flagCell[mcell] == -1){
407               if(etbmax < etcmin)
408                  flagCell[mcell] = 1; //flag cell as used
409               else
410                  flagCell[mcell] = 0; // leave it free
411             }
412         }
413         //store tmp jet info !!!
414        if(etbmax < etcmin) {
415              etaAlgoJet[nJets] = eta;
416              phiAlgoJet[nJets] = phi;
417              etAlgoJet[nJets] = etCone;
418              ncellsAlgoJet[nJets] = nCellIn;
419              nJets++;
420         }
421
422   } // end of cells loop
423
424   //reorder jets by et in cone
425   //sort jets by energy
426   Int_t * idx  = new Int_t[nJets];
427   TMath::Sort(nJets, etAlgoJet, idx);
428   for(Int_t p = 0; p < nJets; p++){
429      etaJet[p] = etaAlgoJet[idx[p]];
430      phiJet[p] = phiAlgoJet[idx[p]];
431      etJet[p] = etAlgoJet[idx[p]];
432      etallJet[p] = etAlgoJet[idx[p]];
433      ncellsJet[p] = ncellsAlgoJet[idx[p]];
434   }
435
436
437   //delete
438   delete[] index;
439   delete[] idx;
440
441 }
442 ////////////////////////////////////////////////////////////////////////
443
444 void AliUA1JetFinderV1::SubtractBackg(Int_t& nIn, Int_t&nJ, Float_t&etbgTotalN,
445                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
446                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
447                       Int_t* multJet, Int_t* injet)
448 {
449   //background subtraction using cone method but without correction in dE/deta distribution
450
451   //calculate energy inside and outside cones
452   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
453   Float_t rc= header->GetRadius();
454   Float_t etIn[30];
455   Float_t etOut = 0;
456   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
457      // if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
458      for(Int_t ijet=0; ijet<nJ; ijet++){
459          Float_t deta = etaT[jpart] - etaJet[ijet];
460               Float_t dphi = phiT[jpart] - phiJet[ijet];
461          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
462               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
463               Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
464          if(dr <= rc){ // particles inside this cone
465              multJet[ijet]++;
466              injet[jpart] = ijet;
467              if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
468                 etIn[ijet] += ptT[jpart];
469                 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet]+= ptT[jpart];
470              }
471              break;
472          }
473      }// end jets loop
474      if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
475         etOut += ptT[jpart]; // particle outside cones and pt cut
476   } //end particle loop
477
478   //estimate jets and background areas
479   Float_t areaJet[30];
480   Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
481   for(Int_t k=0; k<nJ; k++){
482       Float_t detamax = etaJet[k] + rc;
483       Float_t detamin = etaJet[k] - rc;
484       Float_t accmax = 0.0; Float_t accmin = 0.0;
485       if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
486          Float_t h = header->GetLegoEtaMax() - etaJet[k];
487          accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
488       }
489       if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
490          Float_t h = header->GetLegoEtaMax() + etaJet[k];
491          accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
492       }
493       areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
494       areaOut = areaOut - areaJet[k];
495   }
496   //subtract background using area method
497   for(Int_t ljet=0; ljet<nJ; ljet++){
498      Float_t areaRatio = areaJet[ljet]/areaOut;
499      etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
500   }
501
502   // estimate new total background
503   Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
504   etbgTotalN = etOut*areaT/areaOut;
505
506
507 }
508
509 ////////////////////////////////////////////////////////////////////////
510
511 void AliUA1JetFinderV1::SubtractBackgStat(Int_t& nIn, Int_t&nJ,Float_t&etbgTotalN,
512                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
513                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
514                       Int_t* multJet, Int_t* injet)
515 {
516
517   //background subtraction using statistical method
518   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
519   Float_t etbgStat = header->GetBackgStat(); // pre-calculated background
520
521   //calculate energy inside
522   Float_t rc= header->GetRadius();
523   Float_t etIn[30];
524
525   for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
526      //if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
527      for(Int_t ijet=0; ijet<nJ; ijet++){
528          Float_t deta = etaT[jpart] - etaJet[ijet];
529               Float_t dphi = phiT[jpart] - phiJet[ijet];
530          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
531               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
532               Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
533          if(dr <= rc){ // particles inside this cone
534              multJet[ijet]++;
535              injet[jpart] = ijet;
536              if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
537                 etIn[ijet]+= ptT[jpart];
538                 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
539              }
540              break;
541          }
542      }// end jets loop
543   } //end particle loop
544
545   //calc jets areas
546   Float_t areaJet[30];
547   Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
548   for(Int_t k=0; k<nJ; k++){
549       Float_t detamax = etaJet[k] + rc;
550       Float_t detamin = etaJet[k] - rc;
551       Float_t accmax = 0.0; Float_t accmin = 0.0;
552       if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
553          Float_t h = header->GetLegoEtaMax() - etaJet[k];
554          accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
555       }
556       if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
557          Float_t h = header->GetLegoEtaMax() + etaJet[k];
558          accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
559       }
560       areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
561   }
562
563   //subtract background using area method
564   for(Int_t ljet=0; ljet<nJ; ljet++){
565      Float_t areaRatio = areaJet[ljet]/areaOut;
566      etJet[ljet] = etIn[ljet]-etbgStat*areaRatio; // subtraction
567   }
568
569   etbgTotalN = etbgStat;
570
571 }
572
573 ////////////////////////////////////////////////////////////////////////
574
575 void AliUA1JetFinderV1::SubtractBackgCone(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
576                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
577                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
578                       Int_t* multJet, Int_t* injet)
579 {
580    // Cone background subtraction method taking into acount dEt/deta distribution
581     AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
582    //general
583    Float_t rc= header->GetRadius();
584    Float_t etamax = header->GetLegoEtaMax();
585    Float_t etamin = header->GetLegoEtaMin();
586    Int_t ndiv = 100;
587
588    // jet energy and area arrays
589    TH1F* hEtJet[30];
590    TH1F* hAreaJet[30];
591    for(Int_t mjet=0; mjet<nJ; mjet++){
592      char hEtname[256]; char hAreaname[256];
593      sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
594      hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
595      hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
596   }
597    // background energy and area
598    TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
599    TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
600
601    //fill energies
602    for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
603      for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
604          Float_t deta = etaT[jpart] - etaJet[ijet];
605          Float_t dphi = phiT[jpart] - phiJet[ijet];
606          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
607          if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
608          Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
609          if(dr <= rc){ // particles inside this cone
610              injet[jpart] = ijet;
611              multJet[ijet]++;
612              if((fReader->GetCutFlag(jpart)) == 1){// pt cut
613                 hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone
614                 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
615              }
616              break;
617          }
618      }// end jets loop
619      if(injet[jpart] == -1  && fReader->GetCutFlag(jpart) == 1)
620         hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
621   } //end particle loop
622
623    //calc areas
624    Float_t eta0 = etamin;
625    Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
626    Float_t eta1 = eta0 + etaw;
627    for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
628       Float_t etac = eta0 + etaw/2.0;
629       Float_t areabg = etaw*2.0*TMath::Pi();
630       for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
631           Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
632           Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
633           Float_t acc0 = 0.0; Float_t acc1 = 0.0;
634           Float_t areaj = 0.0;
635           if(deta0 > rc && deta1 < rc){
636              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
637              areaj = acc1;
638           }
639           if(deta0 < rc && deta1 > rc){
640              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
641              areaj = acc0;
642           }
643           if(deta0 < rc && deta1 < rc){
644              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
645              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
646              if(eta1<etaJet[ijet]) areaj = acc1-acc0;  // case 1
647              if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
648              if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
649           }
650           hAreaJet[ijet]->Fill(etac,areaj);
651           areabg = areabg - areaj;
652       } // end jets loop
653       hAreaBackg->Fill(etac,areabg);
654       eta0 = eta1;
655       eta1 = eta1 + etaw;
656    } // end loop for all eta bins
657
658    //subtract background
659    for(Int_t kjet=0; kjet<nJ; kjet++){
660        etJet[kjet] = 0.0; // first  clear etJet for this jet
661        for(Int_t bin = 0; bin< ndiv; bin++){
662            if(hAreaJet[kjet]->GetBinContent(bin)){
663               Float_t areab = hAreaBackg->GetBinContent(bin);
664               Float_t etb = hEtBackg->GetBinContent(bin);
665               Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
666               etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR); //subtraction
667            }
668        }
669    }
670
671    // calc background total
672    Double_t etOut = hEtBackg->Integral();
673    Double_t areaOut = hAreaBackg->Integral();
674    Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
675    etbgTotalN = etOut*areaT/areaOut;
676
677    //delete
678    for(Int_t ljet=0; ljet<nJ; ljet++){  // loop for all jets
679        delete hEtJet[ljet];
680        delete hAreaJet[ljet];
681    }
682
683    delete hEtBackg;
684    delete hAreaBackg;
685 }
686
687 ////////////////////////////////////////////////////////////////////////
688
689
690 void AliUA1JetFinderV1::SubtractBackgRatio(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
691                       Float_t* ptT, Float_t* etaT, Float_t* phiT,
692                       Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
693                        Int_t* multJet, Int_t* injet)
694 {
695    // Ratio background subtraction method taking into acount dEt/deta distribution
696     AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
697    //factor F calc before
698     Float_t bgRatioCut = header->GetBackgCutRatio();
699
700
701    //general
702    Float_t rc= header->GetRadius();
703    Float_t etamax = header->GetLegoEtaMax();
704    Float_t etamin = header->GetLegoEtaMin();
705    Int_t ndiv = 100;
706
707    // jet energy and area arrays
708    TH1F* hEtJet[30];
709    TH1F* hAreaJet[30];
710    for(Int_t mjet=0; mjet<nJ; mjet++){
711      char hEtname[256]; char hAreaname[256];
712      sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
713      hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);        // change range
714      hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);  // change range
715   }
716    // background energy and area
717    TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);         // change range
718    TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);  // change range
719
720    //fill energies
721    for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
722      //if((fReader->GetCutFlag(jpart)) != 1) continue;
723      for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
724          Float_t deta = etaT[jpart] - etaJet[ijet];
725               Float_t dphi = phiT[jpart] - phiJet[ijet];
726          if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
727               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
728               Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
729          if(dr <= rc){ // particles inside this cone
730             multJet[ijet]++;
731             injet[jpart] = ijet;
732             if((fReader->GetCutFlag(jpart)) == 1){ //pt cut
733                hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone and pt cut
734                if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
735             }
736             break;
737          }
738      }// end jets loop
739      if(injet[jpart] == -1) hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
740   } //end particle loop
741
742    //calc areas
743    Float_t eta0 = etamin;
744    Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
745    Float_t eta1 = eta0 + etaw;
746    for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
747       Float_t etac = eta0 + etaw/2.0;
748       Float_t areabg = etaw*2.0*TMath::Pi();
749       for(Int_t ijet=0; ijet<nJ; ijet++){  // loop for all jets
750           Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
751           Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
752           Float_t acc0 = 0.0; Float_t acc1 = 0.0;
753           Float_t areaj = 0.0;
754           if(deta0 > rc && deta1 < rc){
755              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
756              areaj = acc1;
757           }
758           if(deta0 < rc && deta1 > rc){
759              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
760              areaj = acc0;
761           }
762           if(deta0 < rc && deta1 < rc){
763              acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
764              acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
765              if(eta1<etaJet[ijet]) areaj = acc1-acc0;  // case 1
766              if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
767              if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
768           }
769           hAreaJet[ijet]->Fill(etac,areaj);
770           areabg = areabg - areaj;
771       } // end jets loop
772       hAreaBackg->Fill(etac,areabg);
773       eta0 = eta1;
774       eta1 = eta1 + etaw;
775    } // end loop for all eta bins
776
777    //subtract background
778    for(Int_t kjet=0; kjet<nJ; kjet++){
779        etJet[kjet] = 0.0; // first  clear etJet for this jet
780        for(Int_t bin = 0; bin< ndiv; bin++){
781            if(hAreaJet[kjet]->GetBinContent(bin)){
782               Float_t areab = hAreaBackg->GetBinContent(bin);
783               Float_t etb = hEtBackg->GetBinContent(bin);
784               Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
785               etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR*bgRatioCut); //subtraction
786            }
787        }
788    }
789
790    // calc background total
791    Double_t etOut = hEtBackg->Integral();
792    Double_t areaOut = hAreaBackg->Integral();
793    Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
794    etbgTotalN = etOut*areaT/areaOut;
795
796    //delete
797    for(Int_t ljet=0; ljet<nJ; ljet++){  // loop for all jets
798        delete hEtJet[ljet];
799        delete hAreaJet[ljet];
800    }
801
802    delete hEtBackg;
803    delete hAreaBackg;
804 }
805
806 ////////////////////////////////////////////////////////////////////////
807
808
809 void AliUA1JetFinderV1::Reset()
810 {
811   fLego->Reset();
812   fJets->ClearJets();
813   AliJetFinder::Reset();
814 }
815
816 ////////////////////////////////////////////////////////////////////////
817
818 void AliUA1JetFinderV1::WriteJHeaderToFile()
819 {
820   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
821   header->Write();
822 }
823
824 ////////////////////////////////////////////////////////////////////////
825
826 void AliUA1JetFinderV1::Init()
827 {
828   // initializes some variables
829   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
830   fLego = new
831     TH2F("legoH","eta-phi",
832          header->GetLegoNbinEta(), header->GetLegoEtaMin(),
833          header->GetLegoEtaMax(),  header->GetLegoNbinPhi(),
834          header->GetLegoPhiMin(),  header->GetLegoPhiMax());
835   // Do not store in current dir
836   fLego->SetDirectory(0);
837
838 }