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