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 for(Int_t kj=0; kj<nj; kj++){
165 if ((etaJet[kj] > (header->GetJetEtaMax())) ||
166 (etaJet[kj] < (header->GetJetEtaMin())) ||
167 (etJet[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
168 Float_t px, py,pz,en; // convert to 4-vector
169 px = etJet[kj] * TMath::Cos(phiJet[kj]);
170 py = etJet[kj] * TMath::Sin(phiJet[kj]);
171 pz = etJet[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJet[kj])));
172 en = TMath::Sqrt(px * px + py * py + pz * pz);
173 fJets->AddJet(px, py, pz, en);
174 idxjets[nselectj] = kj;
177 //add signal percentage and total signal in AliJets for analysis tool
178 Float_t* percentage = new Float_t[nselectj];
179 Int_t* ncells = new Int_t[nselectj];
180 Int_t* mult = new Int_t[nselectj];
181 for(Int_t i = 0; i< nselectj; i++){
182 percentage[i] = etsigJet[idxjets[i]]/etJet[idxjets[i]];
183 ncells[i] = ncellsJet[idxjets[i]];
184 mult[i] = multJet[idxjets[i]];
186 //add particle-injet relationship ///
187 for(Int_t bj = 0; bj < nIn; bj++){
188 if(injet[bj] == -1) continue; //background particle
190 for(Int_t ci = 0; ci< nselectj; ci++){
191 if(injet[bj] == idxjets[ci]){
197 if(bflag == 0) injet[bj] = -1; // set as background particle
199 fJets->SetNCells(ncells);
200 fJets->SetPtFromSignal(percentage);
201 fJets->SetMultiplicities(mult);
202 fJets->SetInJet(injet);
203 fJets->SetEtaIn(etaT);
204 fJets->SetPhiIn(phiT);
206 fJets->SetEtAvg(etbgTotal/(4*(header->GetLegoEtaMax())*TMath::Pi()));
230 ////////////////////////////////////////////////////////////////////////
232 void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t& nJets,
233 Float_t* etJet,Float_t* etaJet, Float_t* phiJet,
234 Float_t* etallJet, Int_t* ncellsJet)
238 // check enough space! *to be done*
239 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
240 Float_t etCell[60000]; //! Cell Energy
241 Float_t etaCell[60000]; //! Cell eta
242 Float_t phiCell[60000]; //! Cell phi
243 Int_t flagCell[60000]; //! Cell flag
246 TAxis* xaxis = fLego->GetXaxis();
247 TAxis* yaxis = fLego->GetYaxis();
249 for (Int_t i = 1; i <= header->GetLegoNbinEta(); i++) {
250 for (Int_t j = 1; j <= header->GetLegoNbinPhi(); j++) {
251 e = fLego->GetBinContent(i,j);
252 if (e < 0.0) continue; // don't include this cells
253 Float_t eta = xaxis->GetBinCenter(i);
254 Float_t phi = yaxis->GetBinCenter(j);
256 etaCell[nCell] = eta;
257 phiCell[nCell] = phi;
258 flagCell[nCell] = 0; //default
263 // Parameters from header
264 Float_t minmove = header->GetMinMove();
265 Float_t maxmove = header->GetMaxMove();
266 Float_t rc = header->GetRadius();
267 Float_t etseed = header->GetEtSeed();
268 //Float_t etmin = header->GetMinJetEt();
272 // tmp array of jets form algoritm
273 Float_t etaAlgoJet[30];
274 Float_t phiAlgoJet[30];
275 Float_t etAlgoJet[30];
276 Int_t ncellsAlgoJet[30];
281 Int_t * index = new Int_t[nCell];
282 TMath::Sort(nCell, etCell, index);
283 // variable used in centroide loop
301 for(Int_t icell = 0; icell < nCell; icell++){
302 Int_t jcell = index[icell];
303 if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
304 if(flagCell[jcell] != 0) continue; // if cell was used before
305 eta = etaCell[jcell];
306 phi = phiCell[jcell];
317 for(Int_t kcell =0; kcell < nCell; kcell++){
318 Int_t lcell = index[kcell];
319 if(lcell == jcell) continue; // cell itself
320 if(flagCell[lcell] != 0) continue; // cell used before
321 if(etCell[lcell] > etCell[jcell]) continue;
323 deta = etaCell[lcell] - eta;
324 dphi = phiCell[lcell] - phi;
325 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
326 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
327 dr = TMath::Sqrt(deta * deta + dphi * dphi);
329 // calculate offset from initiate cell
330 deta = etaCell[lcell] - eta0;
331 dphi = phiCell[lcell] - phi0;
332 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
333 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
334 etas = etas + etCell[lcell]*deta;
335 phis = phis + etCell[lcell]*dphi;
336 ets = ets + etCell[lcell];
337 //new weighted eta and phi including this cell
338 eta = eta0 + etas/ets;
339 phi = phi0 + phis/ets;
340 // if cone does not move much, just go to next step
341 dr = TMath::Sqrt((eta-etab)*(eta-etab) + (phi-phib)*(phi-phib));
342 if(dr <= minmove) break;
343 // cone should not move more than max_mov
344 dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
351 }else{ // store this loop information
359 }//end of cells loop looking centroide
361 //avoid cones overloap (to be implemented in the future)
363 //flag cells in Rc, estimate total energy in cone
364 Float_t etCone = 0.0;
366 rc = header->GetRadius();
367 for(Int_t ncell =0; ncell < nCell; ncell++){
368 if(flagCell[ncell] != 0) continue; // cell used before
370 deta = etaCell[ncell] - eta;
371 dphi = phiCell[ncell] - phi;
372 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
373 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
374 dr = TMath::Sqrt(deta * deta + dphi * dphi);
375 if(dr <= rc){ // cell in cone
376 flagCell[ncell] = -1;
377 etCone+=etCell[ncell];
382 // select jets with et > background
383 // estimate max fluctuation of background in cone
384 Double_t ncellin = (Double_t)nCellIn;
385 Double_t ntcell = (Double_t)nCell;
386 Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/ntcell);
388 Double_t etcmin = etCone ; // could be used etCone - etmin !!
389 //desicions !! etbmax < etcmin
390 for(Int_t mcell =0; mcell < nCell; mcell++){
391 if(flagCell[mcell] == -1){
393 flagCell[mcell] = 1; //flag cell as used
395 flagCell[mcell] = 0; // leave it free
398 //store tmp jet info !!!
399 if(etbmax < etcmin) {
400 etaAlgoJet[nJets] = eta;
401 phiAlgoJet[nJets] = phi;
402 etAlgoJet[nJets] = etCone;
403 ncellsAlgoJet[nJets] = nCellIn;
407 } // end of cells loop
409 //reorder jets by et in cone
410 //sort jets by energy
411 Int_t * idx = new Int_t[nJets];
412 TMath::Sort(nJets, etAlgoJet, idx);
413 for(Int_t p = 0; p < nJets; p++){
414 etaJet[p] = etaAlgoJet[idx[p]];
415 phiJet[p] = phiAlgoJet[idx[p]];
416 etJet[p] = etAlgoJet[idx[p]];
417 etallJet[p] = etAlgoJet[idx[p]];
418 ncellsJet[p] = ncellsAlgoJet[idx[p]];
427 ////////////////////////////////////////////////////////////////////////
429 void AliUA1JetFinderV1::SubtractBackg(Int_t& nIn, Int_t&nJ, Float_t&etbgTotalN,
430 Float_t* ptT, Float_t* etaT, Float_t* phiT,
431 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
432 Int_t* multJet, Int_t* injet)
434 //background subtraction using cone method but without correction in dE/deta distribution
436 //calculate energy inside and outside cones
437 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
438 Float_t rc= header->GetRadius();
441 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
442 // if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
443 for(Int_t ijet=0; ijet<nJ; ijet++){
444 Float_t deta = etaT[jpart] - etaJet[ijet];
445 Float_t dphi = phiT[jpart] - phiJet[ijet];
446 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
447 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
448 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
449 if(dr <= rc){ // particles inside this cone
452 if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
453 etIn[ijet] += ptT[jpart];
454 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet]+= ptT[jpart];
459 if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
460 etOut += ptT[jpart]; // particle outside cones and pt cut
461 } //end particle loop
463 //estimate jets and background areas
465 Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
466 for(Int_t k=0; k<nJ; k++){
467 Float_t detamax = etaJet[k] + rc;
468 Float_t detamin = etaJet[k] - rc;
469 Float_t accmax = 0.0; Float_t accmin = 0.0;
470 if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
471 Float_t h = header->GetLegoEtaMax() - etaJet[k];
472 accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
474 if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
475 Float_t h = header->GetLegoEtaMax() + etaJet[k];
476 accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
478 areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
479 areaOut = areaOut - areaJet[k];
481 //subtract background using area method
482 for(Int_t ljet=0; ljet<nJ; ljet++){
483 Float_t areaRatio = areaJet[ljet]/areaOut;
484 etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
487 // estimate new total background
488 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
489 etbgTotalN = etOut*areaT/areaOut;
494 ////////////////////////////////////////////////////////////////////////
496 void AliUA1JetFinderV1::SubtractBackgStat(Int_t& nIn, Int_t&nJ,Float_t&etbgTotalN,
497 Float_t* ptT, Float_t* etaT, Float_t* phiT,
498 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
499 Int_t* multJet, Int_t* injet)
502 //background subtraction using statistical method
503 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
504 Float_t etbgStat = header->GetBackgStat(); // pre-calculated background
506 //calculate energy inside
507 Float_t rc= header->GetRadius();
510 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
511 //if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
512 for(Int_t ijet=0; ijet<nJ; ijet++){
513 Float_t deta = etaT[jpart] - etaJet[ijet];
514 Float_t dphi = phiT[jpart] - phiJet[ijet];
515 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
516 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
517 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
518 if(dr <= rc){ // particles inside this cone
521 if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
522 etIn[ijet]+= ptT[jpart];
523 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
528 } //end particle loop
532 Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
533 for(Int_t k=0; k<nJ; k++){
534 Float_t detamax = etaJet[k] + rc;
535 Float_t detamin = etaJet[k] - rc;
536 Float_t accmax = 0.0; Float_t accmin = 0.0;
537 if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
538 Float_t h = header->GetLegoEtaMax() - etaJet[k];
539 accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
541 if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
542 Float_t h = header->GetLegoEtaMax() + etaJet[k];
543 accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
545 areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
548 //subtract background using area method
549 for(Int_t ljet=0; ljet<nJ; ljet++){
550 Float_t areaRatio = areaJet[ljet]/areaOut;
551 etJet[ljet] = etIn[ljet]-etbgStat*areaRatio; // subtraction
554 etbgTotalN = etbgStat;
558 ////////////////////////////////////////////////////////////////////////
560 void AliUA1JetFinderV1::SubtractBackgCone(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
561 Float_t* ptT, Float_t* etaT, Float_t* phiT,
562 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
563 Int_t* multJet, Int_t* injet)
565 // Cone background subtraction method taking into acount dEt/deta distribution
566 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
568 Float_t rc= header->GetRadius();
569 Float_t etamax = header->GetLegoEtaMax();
570 Float_t etamin = header->GetLegoEtaMin();
573 // jet energy and area arrays
576 for(Int_t mjet=0; mjet<nJ; mjet++){
577 char hEtname[256]; char hAreaname[256];
578 sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
579 hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
580 hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
582 // background energy and area
583 TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
584 TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
587 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
588 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
589 Float_t deta = etaT[jpart] - etaJet[ijet];
590 Float_t dphi = phiT[jpart] - phiJet[ijet];
591 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
592 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
593 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
594 if(dr <= rc){ // particles inside this cone
597 if((fReader->GetCutFlag(jpart)) == 1){// pt cut
598 hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone
599 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
604 if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
605 hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
606 } //end particle loop
609 Float_t eta0 = etamin;
610 Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
611 Float_t eta1 = eta0 + etaw;
612 for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
613 Float_t etac = eta0 + etaw/2.0;
614 Float_t areabg = etaw*2.0*TMath::Pi();
615 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
616 Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
617 Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
618 Float_t acc0 = 0.0; Float_t acc1 = 0.0;
620 if(deta0 > rc && deta1 < rc){
621 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
624 if(deta0 < rc && deta1 > rc){
625 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
628 if(deta0 < rc && deta1 < rc){
629 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
630 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
631 if(eta1<etaJet[ijet]) areaj = acc1-acc0; // case 1
632 if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
633 if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
635 hAreaJet[ijet]->Fill(etac,areaj);
636 areabg = areabg - areaj;
638 hAreaBackg->Fill(etac,areabg);
641 } // end loop for all eta bins
643 //subtract background
644 for(Int_t kjet=0; kjet<nJ; kjet++){
645 etJet[kjet] = 0.0; // first clear etJet for this jet
646 for(Int_t bin = 0; bin< ndiv; bin++){
647 if(hAreaJet[kjet]->GetBinContent(bin)){
648 Float_t areab = hAreaBackg->GetBinContent(bin);
649 Float_t etb = hEtBackg->GetBinContent(bin);
650 Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
651 etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR); //subtraction
656 // calc background total
657 Double_t etOut = hEtBackg->Integral();
658 Double_t areaOut = hAreaBackg->Integral();
659 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
660 etbgTotalN = etOut*areaT/areaOut;
663 for(Int_t ljet=0; ljet<nJ; ljet++){ // loop for all jets
665 delete hAreaJet[ljet];
672 ////////////////////////////////////////////////////////////////////////
675 void AliUA1JetFinderV1::SubtractBackgRatio(Int_t& nIn, Int_t&nJ,Float_t& etbgTotalN,
676 Float_t* ptT, Float_t* etaT, Float_t* phiT,
677 Float_t* etJet,Float_t* etaJet, Float_t* phiJet, Float_t* etsigJet,
678 Int_t* multJet, Int_t* injet)
680 // Ratio background subtraction method taking into acount dEt/deta distribution
681 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
682 //factor F calc before
683 Float_t bgRatioCut = header->GetBackgCutRatio();
687 Float_t rc= header->GetRadius();
688 Float_t etamax = header->GetLegoEtaMax();
689 Float_t etamin = header->GetLegoEtaMin();
692 // jet energy and area arrays
695 for(Int_t mjet=0; mjet<nJ; mjet++){
696 char hEtname[256]; char hAreaname[256];
697 sprintf(hEtname, "hEtJet%d", mjet); sprintf(hAreaname, "hAreaJet%d", mjet);
698 hEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax); // change range
699 hAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax); // change range
701 // background energy and area
702 TH1F* hEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax); // change range
703 TH1F* hAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax); // change range
706 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
707 //if((fReader->GetCutFlag(jpart)) != 1) continue;
708 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
709 Float_t deta = etaT[jpart] - etaJet[ijet];
710 Float_t dphi = phiT[jpart] - phiJet[ijet];
711 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
712 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
713 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
714 if(dr <= rc){ // particles inside this cone
717 if((fReader->GetCutFlag(jpart)) == 1){ //pt cut
718 hEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone and pt cut
719 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
724 if(injet[jpart] == -1) hEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
725 } //end particle loop
728 Float_t eta0 = etamin;
729 Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
730 Float_t eta1 = eta0 + etaw;
731 for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
732 Float_t etac = eta0 + etaw/2.0;
733 Float_t areabg = etaw*2.0*TMath::Pi();
734 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
735 Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
736 Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
737 Float_t acc0 = 0.0; Float_t acc1 = 0.0;
739 if(deta0 > rc && deta1 < rc){
740 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
743 if(deta0 < rc && deta1 > rc){
744 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
747 if(deta0 < rc && deta1 < rc){
748 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
749 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
750 if(eta1<etaJet[ijet]) areaj = acc1-acc0; // case 1
751 if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
752 if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
754 hAreaJet[ijet]->Fill(etac,areaj);
755 areabg = areabg - areaj;
757 hAreaBackg->Fill(etac,areabg);
760 } // end loop for all eta bins
762 //subtract background
763 for(Int_t kjet=0; kjet<nJ; kjet++){
764 etJet[kjet] = 0.0; // first clear etJet for this jet
765 for(Int_t bin = 0; bin< ndiv; bin++){
766 if(hAreaJet[kjet]->GetBinContent(bin)){
767 Float_t areab = hAreaBackg->GetBinContent(bin);
768 Float_t etb = hEtBackg->GetBinContent(bin);
769 Float_t areaR = (hAreaJet[kjet]->GetBinContent(bin))/areab;
770 etJet[kjet] = etJet[kjet] + ((hEtJet[kjet]->GetBinContent(bin)) - etb*areaR*bgRatioCut); //subtraction
775 // calc background total
776 Double_t etOut = hEtBackg->Integral();
777 Double_t areaOut = hAreaBackg->Integral();
778 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
779 etbgTotalN = etOut*areaT/areaOut;
782 for(Int_t ljet=0; ljet<nJ; ljet++){ // loop for all jets
784 delete hAreaJet[ljet];
791 ////////////////////////////////////////////////////////////////////////
794 void AliUA1JetFinderV1::Reset()
800 ////////////////////////////////////////////////////////////////////////
802 void AliUA1JetFinderV1::WriteJHeaderToFile()
804 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
809 ////////////////////////////////////////////////////////////////////////
811 void AliUA1JetFinderV1::Init()
813 // initializes some variables
814 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
817 TH2F("legoH","eta-phi",
818 header->GetLegoNbinEta(), header->GetLegoEtaMin(),
819 header->GetLegoEtaMax(), header->GetLegoNbinPhi(),
820 header->GetLegoPhiMin(), header->GetLegoPhiMax());