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