1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //---------------------------------------------------------------------
19 // UA1 Cone Algorithm Jet finder
20 // manages the search for jets
21 // Author: Rafael.Diaz.Valdes@cern.ch
23 //---------------------------------------------------------------------
26 #include <TClonesArray.h>
30 #include <TLorentzVector.h>
32 #include "AliUA1JetFinderV1.h"
33 #include "AliUA1JetHeaderV1.h"
34 #include "AliJetReaderHeader.h"
35 #include "AliJetReader.h"
39 ClassImp(AliUA1JetFinderV1)
41 ////////////////////////////////////////////////////////////////////////
43 AliUA1JetFinderV1::AliUA1JetFinderV1()
51 ////////////////////////////////////////////////////////////////////////
53 AliUA1JetFinderV1::~AliUA1JetFinderV1()
59 ////////////////////////////////////////////////////////////////////////
62 void AliUA1JetFinderV1::FindJets()
65 //1) Fill cell map array
66 //2) calculate total energy and fluctuation level
68 // 3.1) look centroides in cell map
69 // 3.2) calculate total energy in cones
70 // 3.3) flag as a possible jet
71 // 3.4) reorder cones by energy
72 //4) subtract backg in accepted jets
75 // transform input to pt,eta,phi plus lego
77 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
78 TClonesArray *lvArray = fReader->GetMomentumArray();
79 Int_t nIn = lvArray->GetEntries();
82 // local arrays for input
83 Float_t* ptT = new Float_t[nIn];
84 Float_t* etaT = new Float_t[nIn];
85 Float_t* phiT = new Float_t[nIn];
86 Int_t* injet = new Int_t[nIn];
88 //total energy in array
89 Float_t etbgTotal = 0.0;
90 TH1F* hPtTotal = new TH1F("hPt","Pt distribution of all particles ",100,0.0,15.0);
92 // load input vectors and calculate total energy in array
93 for (Int_t i = 0; i < nIn; i++){
94 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
97 phiT[i] = ((lv->Phi() < 0) ? (lv->Phi()) + 2 * TMath::Pi() : lv->Phi());
98 if (fReader->GetCutFlag(i) != 1) continue;
99 fLego->Fill(etaT[i], phiT[i], ptT[i]);
100 hPtTotal->Fill(ptT[i]);
103 fJets->SetNinput(nIn);
105 // calculate total energy and fluctuation in map
106 Double_t meanpt = hPtTotal->GetMean();
107 Double_t ptRMS = hPtTotal->GetRMS();
108 Double_t npart = hPtTotal->GetEntries();
109 Double_t dEtTotal = (TMath::Sqrt(npart))*TMath::Sqrt(meanpt * meanpt + ptRMS*ptRMS);
111 // arrays to hold jets
112 Float_t* etaJet = new Float_t[30];
113 Float_t* phiJet = new Float_t[30];
114 Float_t* etJet = new Float_t[30];
115 Float_t* etsigJet = new Float_t[30]; //signal et in jet
116 Float_t* etallJet = new Float_t[30]; // total et in jet (tmp variable)
117 Int_t* ncellsJet = new Int_t[30];
118 Int_t* multJet = new Int_t[30];
119 Int_t nJets; // to hold number of jets found by algorithm
120 Int_t nj; // number of jets accepted
121 Float_t prec = header->GetPrecBg();
123 while(bgprec > prec){
124 //reset jet arrays in memory
125 memset(etaJet,0,sizeof(Float_t)*30);
126 memset(phiJet,0,sizeof(Float_t)*30);
127 memset(etJet,0,sizeof(Float_t)*30);
128 memset(etallJet,0,sizeof(Float_t)*30);
129 memset(etsigJet,0,sizeof(Float_t)*30);
130 memset(ncellsJet,0,sizeof(Int_t)*30);
131 memset(multJet,0,sizeof(Int_t)*30);
134 // reset particles-jet array in memory
135 memset(injet,-1,sizeof(Int_t)*nIn);
136 //run cone algorithm finder
137 RunAlgoritm(etbgTotal,dEtTotal,nJets,etJet,etaJet,phiJet,etallJet,ncellsJet);
138 //run background subtraction
139 if(nJets > header->GetNAcceptJets()) // limited number of accepted jets per event
140 nj = header->GetNAcceptJets();
143 //subtract background
144 Float_t etbgTotalN = 0.0; //new background
145 if(header->GetBackgMode() == 1) // standar
146 SubtractBackg(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
147 if(header->GetBackgMode() == 2) //cone
148 SubtractBackgCone(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
149 if(header->GetBackgMode() == 3) //ratio
150 SubtractBackgRatio(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
151 if(header->GetBackgMode() == 4) //statistic
152 SubtractBackgStat(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
154 if(etbgTotalN != 0.0)
155 bgprec = (etbgTotal - etbgTotalN)/etbgTotalN;
158 etbgTotal = etbgTotalN; // update with new background estimation
162 Int_t* idxjets = new Int_t[nj];
164 printf("Found %d jets \n", nj);
166 for(Int_t kj=0; kj<nj; kj++){
167 if ((etaJet[kj] > (header->GetJetEtaMax())) ||
168 (etaJet[kj] < (header->GetJetEtaMin())) ||
169 (etJet[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
170 Float_t px, py,pz,en; // convert to 4-vector
171 px = etJet[kj] * TMath::Cos(phiJet[kj]);
172 py = etJet[kj] * TMath::Sin(phiJet[kj]);
173 pz = etJet[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJet[kj])));
174 en = TMath::Sqrt(px * px + py * py + pz * pz);
175 fJets->AddJet(px, py, pz, en);
176 idxjets[nselectj] = kj;
179 //add signal percentage and total signal in AliJets for analysis tool
180 Float_t* percentage = new Float_t[nselectj];
181 Int_t* ncells = new Int_t[nselectj];
182 Int_t* mult = new Int_t[nselectj];
183 for(Int_t i = 0; i< nselectj; i++){
184 percentage[i] = etsigJet[idxjets[i]]/etJet[idxjets[i]];
185 ncells[i] = ncellsJet[idxjets[i]];
186 mult[i] = multJet[idxjets[i]];
188 //add particle-injet relationship ///
189 for(Int_t bj = 0; bj < nIn; bj++){
190 if(injet[bj] == -1) continue; //background particle
192 for(Int_t ci = 0; ci< nselectj; ci++){
193 if(injet[bj] == idxjets[ci]){
199 if(bflag == 0) injet[bj] = -1; // set as background particle
201 fJets->SetNCells(ncells);
202 fJets->SetPtFromSignal(percentage);
203 fJets->SetMultiplicities(mult);
204 fJets->SetInJet(injet);
205 fJets->SetEtaIn(etaT);
206 fJets->SetPhiIn(phiT);
208 fJets->SetEtAvg(etbgTotal/(4*(header->GetLegoEtaMax())*TMath::Pi()));
232 ////////////////////////////////////////////////////////////////////////
234 void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t& nJets,
235 Float_t* etJet,Float_t* etaJet, Float_t* phiJet,
236 Float_t* etallJet, Int_t* ncellsJet)
240 // check enough space! *to be done*
241 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
242 Float_t etCell[60000]; //! Cell Energy
243 Float_t etaCell[60000]; //! Cell eta
244 Float_t phiCell[60000]; //! Cell phi
245 Int_t flagCell[60000]; //! Cell flag
248 TAxis* xaxis = fLego->GetXaxis();
249 TAxis* yaxis = fLego->GetYaxis();
251 for (Int_t i = 1; i <= header->GetLegoNbinEta(); i++) {
252 for (Int_t j = 1; j <= header->GetLegoNbinPhi(); j++) {
253 e = fLego->GetBinContent(i,j);
254 if (e < 0.0) continue; // don't include this cells
255 Float_t eta = xaxis->GetBinCenter(i);
256 Float_t phi = yaxis->GetBinCenter(j);
258 etaCell[nCell] = eta;
259 phiCell[nCell] = phi;
260 flagCell[nCell] = 0; //default
265 // Parameters from header
266 Float_t minmove = header->GetMinMove();
267 Float_t maxmove = header->GetMaxMove();
268 Float_t rc = header->GetRadius();
269 Float_t etseed = header->GetEtSeed();
270 //Float_t etmin = header->GetMinJetEt();
274 // tmp array of jets form algoritm
275 Float_t etaAlgoJet[30];
276 Float_t phiAlgoJet[30];
277 Float_t etAlgoJet[30];
278 Int_t ncellsAlgoJet[30];
283 Int_t * index = new Int_t[nCell];
284 TMath::Sort(nCell, etCell, index);
285 // variable used in centroide loop
303 for(Int_t icell = 0; icell < nCell; icell++){
304 Int_t jcell = index[icell];
305 if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
306 if(flagCell[jcell] != 0) continue; // if cell was used before
307 eta = etaCell[jcell];
308 phi = phiCell[jcell];
319 for(Int_t kcell =0; kcell < nCell; kcell++){
320 Int_t lcell = index[kcell];
321 if(lcell == jcell) continue; // cell itself
322 if(flagCell[lcell] != 0) continue; // cell used before
323 if(etCell[lcell] > etCell[jcell]) continue;
325 deta = etaCell[lcell] - eta;
326 dphi = phiCell[lcell] - phi;
327 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
328 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
329 dr = TMath::Sqrt(deta * deta + dphi * dphi);
331 // calculate offset from initiate cell
332 deta = etaCell[lcell] - eta0;
333 dphi = phiCell[lcell] - phi0;
334 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
335 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
336 etas = etas + etCell[lcell]*deta;
337 phis = phis + etCell[lcell]*dphi;
338 ets = ets + etCell[lcell];
339 //new weighted eta and phi including this cell
340 eta = eta0 + etas/ets;
341 phi = phi0 + phis/ets;
342 // if cone does not move much, just go to next step
343 dr = TMath::Sqrt((eta-etab)*(eta-etab) + (phi-phib)*(phi-phib));
344 if(dr <= minmove) break;
345 // cone should not move more than max_mov
346 dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
353 }else{ // store this loop information
361 }//end of cells loop looking centroide
363 //avoid cones overloap (to be implemented in the future)
365 //flag cells in Rc, estimate total energy in cone
366 Float_t etCone = 0.0;
368 rc = header->GetRadius();
369 for(Int_t ncell =0; ncell < nCell; ncell++){
370 if(flagCell[ncell] != 0) continue; // cell used before
372 deta = etaCell[ncell] - eta;
373 dphi = phiCell[ncell] - phi;
374 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
375 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
376 dr = TMath::Sqrt(deta * deta + dphi * dphi);
377 if(dr <= rc){ // cell in cone
378 flagCell[ncell] = -1;
379 etCone+=etCell[ncell];
384 // select jets with et > background
385 // estimate max fluctuation of background in cone
386 Double_t ncellin = (Double_t)nCellIn;
387 Double_t ntcell = (Double_t)nCell;
388 Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/ntcell);
390 Double_t etcmin = etCone ; // could be used etCone - etmin !!
391 //desicions !! etbmax < etcmin
392 for(Int_t mcell =0; mcell < nCell; mcell++){
393 if(flagCell[mcell] == -1){
395 flagCell[mcell] = 1; //flag cell as used
397 flagCell[mcell] = 0; // leave it free
400 //store tmp jet info !!!
401 if(etbmax < etcmin) {
402 etaAlgoJet[nJets] = eta;
403 phiAlgoJet[nJets] = phi;
404 etAlgoJet[nJets] = etCone;
405 ncellsAlgoJet[nJets] = nCellIn;
409 } // end of cells loop
411 //reorder jets by et in cone
412 //sort jets by energy
413 Int_t * idx = new Int_t[nJets];
414 TMath::Sort(nJets, etAlgoJet, idx);
415 for(Int_t p = 0; p < nJets; p++){
416 etaJet[p] = etaAlgoJet[idx[p]];
417 phiJet[p] = phiAlgoJet[idx[p]];
418 etJet[p] = etAlgoJet[idx[p]];
419 etallJet[p] = etAlgoJet[idx[p]];
420 ncellsJet[p] = ncellsAlgoJet[idx[p]];
429 ////////////////////////////////////////////////////////////////////////
431 void AliUA1JetFinderV1::SubtractBackg(Int_t& nIn, Int_t&nJ, Float_t&etbgTotalN,
432 Float_t* ptT, Float_t* etaT, Float_t* phiT,
433 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
434 Int_t* multJet, Int_t* injet)
436 //background subtraction using cone method but without correction in dE/deta distribution
438 //calculate energy inside and outside cones
439 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
440 Float_t rc= header->GetRadius();
443 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
444 // if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
445 for(Int_t ijet=0; ijet<nJ; ijet++){
446 Float_t deta = etaT[jpart] - etaJet[ijet];
447 Float_t dphi = phiT[jpart] - phiJet[ijet];
448 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
449 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
450 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
451 if(dr <= rc){ // particles inside this cone
454 if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
455 etIn[ijet] += ptT[jpart];
456 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet]+= ptT[jpart];
461 if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
462 etOut += ptT[jpart]; // particle outside cones and pt cut
463 } //end particle loop
465 //estimate jets and background areas
467 Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
468 for(Int_t k=0; k<nJ; k++){
469 Float_t detamax = etaJet[k] + rc;
470 Float_t detamin = etaJet[k] - rc;
471 Float_t accmax = 0.0; Float_t accmin = 0.0;
472 if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
473 Float_t h = header->GetLegoEtaMax() - etaJet[k];
474 accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
476 if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
477 Float_t h = header->GetLegoEtaMax() + etaJet[k];
478 accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
480 areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
481 areaOut = areaOut - areaJet[k];
483 //subtract background using area method
484 for(Int_t ljet=0; ljet<nJ; ljet++){
485 Float_t areaRatio = areaJet[ljet]/areaOut;
486 etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
489 // estimate new total background
490 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
491 etbgTotalN = etOut*areaT/areaOut;
496 ////////////////////////////////////////////////////////////////////////
498 void AliUA1JetFinderV1::SubtractBackgStat(Int_t& nIn, Int_t&nJ,Float_t&etbgTotalN,
499 Float_t* ptT, Float_t* etaT, Float_t* phiT,
500 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
501 Int_t* multJet, Int_t* injet)
504 //background subtraction using statistical method
505 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
506 Float_t etbgStat = header->GetBackgStat(); // pre-calculated background
508 //calculate energy inside
509 Float_t rc= header->GetRadius();
512 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
513 //if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
514 for(Int_t ijet=0; ijet<nJ; ijet++){
515 Float_t deta = etaT[jpart] - etaJet[ijet];
516 Float_t dphi = phiT[jpart] - phiJet[ijet];
517 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
518 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
519 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
520 if(dr <= rc){ // particles inside this cone
523 if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
524 etIn[ijet]+= ptT[jpart];
525 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
530 } //end particle loop
534 Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
535 for(Int_t k=0; k<nJ; k++){
536 Float_t detamax = etaJet[k] + rc;
537 Float_t detamin = etaJet[k] - rc;
538 Float_t accmax = 0.0; Float_t accmin = 0.0;
539 if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
540 Float_t h = header->GetLegoEtaMax() - etaJet[k];
541 accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
543 if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
544 Float_t h = header->GetLegoEtaMax() + etaJet[k];
545 accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
547 areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
550 //subtract background using area method
551 for(Int_t ljet=0; ljet<nJ; ljet++){
552 Float_t areaRatio = areaJet[ljet]/areaOut;
553 etJet[ljet] = etIn[ljet]-etbgStat*areaRatio; // subtraction
556 etbgTotalN = etbgStat;
560 ////////////////////////////////////////////////////////////////////////
562 void AliUA1JetFinderV1::SubtractBackgCone(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
563 Float_t* ptT, Float_t* etaT, Float_t* phiT,
564 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
565 Int_t* multJet, Int_t* injet)
567 // Cone background subtraction method taking into acount dEt/deta distribution
568 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
570 Float_t rc= header->GetRadius();
571 Float_t etamax = header->GetLegoEtaMax();
572 Float_t etamin = header->GetLegoEtaMin();
575 // jet energy and area arrays
578 for(Int_t mjet=0; mjet<nJ; mjet++){
579 char hEtname[256]; char hAreaname[256];
580 sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
581 hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
582 hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
584 // background energy and area
585 TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
586 TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
589 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
590 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
591 Float_t deta = etaT[jpart] - etaJet[ijet];
592 Float_t dphi = phiT[jpart] - phiJet[ijet];
593 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
594 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
595 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
596 if(dr <= rc){ // particles inside this cone
599 if((fReader->GetCutFlag(jpart)) == 1){// pt cut
600 hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone
601 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
606 if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
607 hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
608 } //end particle loop
611 Float_t eta0 = etamin;
612 Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
613 Float_t eta1 = eta0 + etaw;
614 for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
615 Float_t etac = eta0 + etaw/2.0;
616 Float_t areabg = etaw*2.0*TMath::Pi();
617 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
618 Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
619 Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
620 Float_t acc0 = 0.0; Float_t acc1 = 0.0;
622 if(deta0 > rc && deta1 < rc){
623 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
626 if(deta0 < rc && deta1 > rc){
627 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
630 if(deta0 < rc && deta1 < rc){
631 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
632 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
633 if(eta1<etaJet[ijet]) areaj = acc1-acc0; // case 1
634 if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
635 if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
637 hAreaJet[ijet]->Fill(etac,areaj);
638 areabg = areabg - areaj;
640 hAreaBackg->Fill(etac,areabg);
643 } // end loop for all eta bins
645 //subtract background
646 for(Int_t kjet=0; kjet<nJ; kjet++){
647 etJet[kjet] = 0.0; // first clear etJet for this jet
648 for(Int_t bin = 0; bin< ndiv; bin++){
649 if(hAreaJet[kjet]->GetBinContent(bin)){
650 Float_t areab = hAreaBackg->GetBinContent(bin);
651 Float_t etb = hEtBackg->GetBinContent(bin);
652 Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
653 etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR); //subtraction
658 // calc background total
659 Double_t etOut = hEtBackg->Integral();
660 Double_t areaOut = hAreaBackg->Integral();
661 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
662 etbgTotalN = etOut*areaT/areaOut;
665 for(Int_t ljet=0; ljet<nJ; ljet++){ // loop for all jets
667 delete hAreaJet[ljet];
674 ////////////////////////////////////////////////////////////////////////
677 void AliUA1JetFinderV1::SubtractBackgRatio(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
678 Float_t* ptT, Float_t* etaT, Float_t* phiT,
679 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
680 Int_t* multJet, Int_t* injet)
682 // Ratio background subtraction method taking into acount dEt/deta distribution
683 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
684 //factor F calc before
685 Float_t bgRatioCut = header->GetBackgCutRatio();
689 Float_t rc= header->GetRadius();
690 Float_t etamax = header->GetLegoEtaMax();
691 Float_t etamin = header->GetLegoEtaMin();
694 // jet energy and area arrays
697 for(Int_t mjet=0; mjet<nJ; mjet++){
698 char hEtname[256]; char hAreaname[256];
699 sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
700 hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax); // change range
701 hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax); // change range
703 // background energy and area
704 TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax); // change range
705 TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax); // change range
708 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
709 //if((fReader->GetCutFlag(jpart)) != 1) continue;
710 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
711 Float_t deta = etaT[jpart] - etaJet[ijet];
712 Float_t dphi = phiT[jpart] - phiJet[ijet];
713 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
714 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
715 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
716 if(dr <= rc){ // particles inside this cone
719 if((fReader->GetCutFlag(jpart)) == 1){ //pt cut
720 hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone and pt cut
721 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
726 if(injet[jpart] == -1) hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
727 } //end particle loop
730 Float_t eta0 = etamin;
731 Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
732 Float_t eta1 = eta0 + etaw;
733 for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
734 Float_t etac = eta0 + etaw/2.0;
735 Float_t areabg = etaw*2.0*TMath::Pi();
736 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
737 Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
738 Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
739 Float_t acc0 = 0.0; Float_t acc1 = 0.0;
741 if(deta0 > rc && deta1 < rc){
742 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
745 if(deta0 < rc && deta1 > rc){
746 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
749 if(deta0 < rc && deta1 < rc){
750 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
751 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
752 if(eta1<etaJet[ijet]) areaj = acc1-acc0; // case 1
753 if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
754 if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
756 hAreaJet[ijet]->Fill(etac,areaj);
757 areabg = areabg - areaj;
759 hAreaBackg->Fill(etac,areabg);
762 } // end loop for all eta bins
764 //subtract background
765 for(Int_t kjet=0; kjet<nJ; kjet++){
766 etJet[kjet] = 0.0; // first clear etJet for this jet
767 for(Int_t bin = 0; bin< ndiv; bin++){
768 if(hAreaJet[kjet]->GetBinContent(bin)){
769 Float_t areab = hAreaBackg->GetBinContent(bin);
770 Float_t etb = hEtBackg->GetBinContent(bin);
771 Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
772 etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR*bgRatioCut); //subtraction
777 // calc background total
778 Double_t etOut = hEtBackg->Integral();
779 Double_t areaOut = hAreaBackg->Integral();
780 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
781 etbgTotalN = etOut*areaT/areaOut;
784 for(Int_t ljet=0; ljet<nJ; ljet++){ // loop for all jets
786 delete hAreaJet[ljet];
793 ////////////////////////////////////////////////////////////////////////
796 void AliUA1JetFinderV1::Reset()
802 ////////////////////////////////////////////////////////////////////////
804 void AliUA1JetFinderV1::WriteJHeaderToFile()
806 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
810 ////////////////////////////////////////////////////////////////////////
812 void AliUA1JetFinderV1::Init()
814 // initializes some variables
815 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
818 TH2F("legoH","eta-phi",
819 header->GetLegoNbinEta(), header->GetLegoEtaMin(),
820 header->GetLegoEtaMax(), header->GetLegoNbinPhi(),
821 header->GetLegoPhiMin(), header->GetLegoPhiMax());