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