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"
36 #include "AliJetHeader.h"
39 #include "AliAODJet.h"
43 ClassImp(AliUA1JetFinderV1)
45 /////////////////////////////////////////////////////////////////////
47 AliUA1JetFinderV1::AliUA1JetFinderV1() :
54 for(int i = 0;i < kMaxJets;i++){
55 fhAreaJet[i] = fhEtJet[i] = 0;
59 ////////////////////////////////////////////////////////////////////////
61 AliUA1JetFinderV1::~AliUA1JetFinderV1()
67 if(fhEtBackg)delete fhEtBackg;
69 if( fhAreaBackg) delete fhAreaBackg;
71 for(int i = 0;i < kMaxJets;i++){
72 if(fhAreaJet[i])delete fhAreaJet[i];
73 if(fhEtJet[i]) delete fhEtJet[i];
74 fhAreaJet[i] = fhEtJet[i] = 0;
79 ////////////////////////////////////////////////////////////////////////
82 void AliUA1JetFinderV1::FindJets()
85 //1) Fill cell map array
86 //2) calculate total energy and fluctuation level
88 // 3.1) look centroides in cell map
89 // 3.2) calculate total energy in cones
90 // 3.3) flag as a possible jet
91 // 3.4) reorder cones by energy
92 //4) subtract backg in accepted jets
95 // transform input to pt,eta,phi plus lego
97 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
98 TClonesArray *lvArray = fReader->GetMomentumArray();
99 Int_t nIn = lvArray->GetEntries();
100 if (nIn <= 0) return;
102 // local arrays for input
103 // ToDo: check memory fragmentation, maybe better to
104 // define them globally and resize as needed
105 // Fragementation should be worse for low mult...
106 Float_t* ptT = new Float_t[nIn];
107 Float_t* etaT = new Float_t[nIn];
108 Float_t* phiT = new Float_t[nIn];
109 Int_t* injet = new Int_t[nIn];
111 memset(ptT,0,sizeof(Float_t)*nIn);
112 memset(etaT,0,sizeof(Float_t)*nIn);
113 memset(phiT,0,sizeof(Float_t)*nIn);
116 // load input vectors and calculate total energy in array
118 //total energy in array
119 Float_t etbgTotal = 0.0;
123 for (Int_t i = 0; i < nIn; i++){
124 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
127 phiT[i] = ((lv->Phi() < 0) ? (lv->Phi()) + 2 * TMath::Pi() : lv->Phi());
128 if (fReader->GetCutFlag(i) != 1) continue;
129 fLego ->Fill(etaT[i], phiT[i], ptT[i]);
132 etbg2 += ptT[i]*ptT[i];
135 // calculate total energy and fluctuation in map
139 meanpt = etbgTotal/npart;
141 if(etbg2>(meanpt*meanpt)){// prenent NAN, should only happen due to numerical instabilities
142 ptRMS = TMath::Sqrt(etbg2-meanpt*meanpt);
145 Double_t dEtTotal = (TMath::Sqrt(npart))*TMath::Sqrt(meanpt * meanpt + ptRMS*ptRMS);
147 // arrays to hold jets
148 Float_t etaJet[kMaxJets];
149 Float_t phiJet[kMaxJets];
150 Float_t etJet[kMaxJets];
151 Float_t etsigJet[kMaxJets]; //signal et in jet
152 Float_t etallJet[kMaxJets]; // total et in jet (tmp variable)
153 Int_t ncellsJet[kMaxJets];
154 Int_t multJet[kMaxJets];
155 Int_t nJets; // to hold number of jets found by algorithm
156 Int_t nj; // number of jets accepted
157 Float_t prec = header->GetPrecBg();
159 while(bgprec > prec){
160 //reset jet arrays in memory
161 memset(etaJet,0,sizeof(Float_t)*kMaxJets);
162 memset(phiJet,0,sizeof(Float_t)*kMaxJets);
163 memset(etJet,0,sizeof(Float_t)*kMaxJets);
164 memset(etallJet,0,sizeof(Float_t)*kMaxJets);
165 memset(etsigJet,0,sizeof(Float_t)*kMaxJets);
166 memset(ncellsJet,0,sizeof(Int_t)*kMaxJets);
167 memset(multJet,0,sizeof(Int_t)*kMaxJets);
170 // reset particles-jet array in memory
171 memset(injet,-1,sizeof(Int_t)*nIn);
172 //run cone algorithm finder
173 RunAlgoritm(etbgTotal,dEtTotal,nJets,etJet,etaJet,phiJet,etallJet,ncellsJet);
174 //run background subtraction
175 if(nJets > header->GetNAcceptJets()) // limited number of accepted jets per event
176 nj = header->GetNAcceptJets();
179 //subtract background
180 Float_t etbgTotalN = 0.0; //new background
181 if(header->GetBackgMode() == 1) // standar
182 SubtractBackg(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
183 if(header->GetBackgMode() == 2) //cone
184 SubtractBackgCone(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
185 if(header->GetBackgMode() == 3) //ratio
186 SubtractBackgRatio(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
187 if(header->GetBackgMode() == 4) //statistic
188 SubtractBackgStat(nIn,nj,etbgTotalN,ptT,etaT,phiT,etJet,etaJet,phiJet,etsigJet,multJet,injet);
190 if(etbgTotalN != 0.0)
191 bgprec = (etbgTotal - etbgTotalN)/etbgTotalN;
194 etbgTotal = etbgTotalN; // update with new background estimation
197 // add tracks to the jet if it wasn't yet done
198 if (header->GetBackgMode() == 0){
199 Float_t rc= header->GetRadius();
200 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
201 for(Int_t ijet=0; ijet<nj; ijet++){
202 Float_t deta = etaT[jpart] - etaJet[ijet];
203 Float_t dphi = phiT[jpart] - phiJet[ijet];
204 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
205 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
206 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
207 if(dr <= rc){ // particles inside this cone
212 } //end particle loop
216 Int_t idxjets[kMaxJets];
220 Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
221 if (fromAod) refs = fReader->GetReferences();
222 for(Int_t kj=0; kj<nj; kj++){
223 if ((etaJet[kj] > (header->GetJetEtaMax())) ||
224 (etaJet[kj] < (header->GetJetEtaMin())) ||
225 (etJet[kj] < header->GetMinJetEt())) continue; // acceptance eta range and etmin
226 Float_t px, py,pz,en; // convert to 4-vector
227 px = etJet[kj] * TMath::Cos(phiJet[kj]);
228 py = etJet[kj] * TMath::Sin(phiJet[kj]);
229 pz = etJet[kj] / TMath::Tan(2.0 * TMath::ATan(TMath::Exp(-etaJet[kj])));
230 en = TMath::Sqrt(px * px + py * py + pz * pz);
232 AliAODJet jet(px, py, pz, en);
235 for(Int_t jpart = 0; jpart < nIn; jpart++) // loop for all particles in array
236 if (injet[jpart] == kj && fReader->GetCutFlag(jpart) == 1)
237 jet.AddTrack(refs->At(jpart)); // check if the particle belongs to the jet and add the ref
244 idxjets[nselectj] = kj;
246 } //end particle loop
248 //add signal percentage and total signal in AliJets for analysis tool
249 Float_t percentage[kMaxJets];
250 Int_t ncells[kMaxJets];
251 Int_t mult[kMaxJets];
252 for(Int_t i = 0; i< nselectj; i++){
253 percentage[i] = etsigJet[idxjets[i]]/etJet[idxjets[i]];
254 ncells[i] = ncellsJet[idxjets[i]];
255 mult[i] = multJet[idxjets[i]];
257 //add particle-injet relationship ///
258 for(Int_t bj = 0; bj < nIn; bj++){
259 if(injet[bj] == -1) continue; //background particle
261 for(Int_t ci = 0; ci< nselectj; ci++){
262 if(injet[bj] == idxjets[ci]){
268 if(bflag == 0) injet[bj] = -1; // set as background particle
278 ////////////////////////////////////////////////////////////////////////
280 void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t& nJets,
281 Float_t* etJet,Float_t* etaJet, Float_t* phiJet,
282 Float_t* etallJet, Int_t* ncellsJet)
286 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
287 const Int_t nBinsMax = 120000; // we use a fixed array not to fragment memory
289 const Int_t nBinEta = header->GetLegoNbinEta();
290 const Int_t nBinPhi = header->GetLegoNbinPhi();
291 if((nBinPhi*nBinEta)>nBinsMax){
292 AliError("Too many bins of the ETA-PHI histogram");
295 Float_t etCell[nBinsMax]; //! Cell Energy
296 Float_t etaCell[nBinsMax]; //! Cell eta
297 Float_t phiCell[nBinsMax]; //! Cell phi
298 Short_t flagCell[nBinsMax]; //! Cell flag
301 TAxis* xaxis = fLego->GetXaxis();
302 TAxis* yaxis = fLego->GetYaxis();
304 for (Int_t i = 1; i <= nBinEta; i++) {
305 for (Int_t j = 1; j <= nBinPhi; j++) {
306 e = fLego->GetBinContent(i,j);
307 if (e < 0.0) continue; // don't include this cells
308 Float_t eta = xaxis->GetBinCenter(i);
309 Float_t phi = yaxis->GetBinCenter(j);
311 etaCell[nCell] = eta;
312 phiCell[nCell] = phi;
313 flagCell[nCell] = 0; //default
318 // Parameters from header
319 Float_t minmove = header->GetMinMove();
320 Float_t maxmove = header->GetMaxMove();
321 Float_t rc = header->GetRadius();
322 Float_t etseed = header->GetEtSeed();
323 //Float_t etmin = header->GetMinJetEt();
327 // tmp array of jets form algoritm
328 Float_t etaAlgoJet[kMaxJets];
329 Float_t phiAlgoJet[kMaxJets];
330 Float_t etAlgoJet[kMaxJets];
331 Int_t ncellsAlgoJet[kMaxJets];
336 Int_t index[nBinsMax];
337 TMath::Sort(nCell, etCell, index);
338 // variable used in centroide loop
357 for(Int_t icell = 0; icell < nCell; icell++){
358 Int_t jcell = index[icell];
359 if(etCell[jcell] <= etseed) continue; // if cell energy is low et seed
360 if(flagCell[jcell] != 0) continue; // if cell was used before
361 eta = etaCell[jcell];
362 phi = phiCell[jcell];
373 for(Int_t kcell =0; kcell < nCell; kcell++){
374 Int_t lcell = index[kcell];
375 if(lcell == jcell) continue; // cell itself
376 if(flagCell[lcell] != 0) continue; // cell used before
377 if(etCell[lcell] > etCell[jcell]) continue; // can this happen
379 deta = etaCell[lcell] - eta;
380 dphi = TMath::Abs(phiCell[lcell] - phi);
381 if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
382 dr = TMath::Sqrt(deta * deta + dphi * dphi);
384 // calculate offset from initiate cell
385 deta = etaCell[lcell] - eta0;
386 dphi = phiCell[lcell] - phi0;
387 if (dphi < - TMath::Pi()) dphi= dphi + 2.0 * TMath::Pi();
388 if (dphi > TMath::Pi()) dphi = dphi - 2.0 * TMath::Pi();
390 etas = etas + etCell[lcell]*deta;
391 phis = phis + etCell[lcell]*dphi;
392 ets = ets + etCell[lcell];
393 //new weighted eta and phi including this cell
394 eta = eta0 + etas/ets;
395 phi = phi0 + phis/ets;
396 // if cone does not move much, just go to next step
397 dphib = TMath::Abs(phi - phib);
398 if (dphib > TMath::Pi()) dphib = 2. * TMath::Pi() - dphib;
399 dr = TMath::Sqrt((eta-etab)*(eta-etab) + dphib * dphib);
400 if(dr <= minmove) break;
401 // cone should not move more than max_mov
402 dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
409 } else { // store this loop information
417 }//end of cells loop looking centroide
419 //avoid cones overloap (to be implemented in the future)
421 //flag cells in Rc, estimate total energy in cone
422 Float_t etCone = 0.0;
424 rc = header->GetRadius();
425 for(Int_t ncell =0; ncell < nCell; ncell++){
426 if(flagCell[ncell] != 0) continue; // cell used before
428 deta = etaCell[ncell] - eta;
429 dphi = phiCell[ncell] - phi;
430 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
431 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
432 dr = TMath::Sqrt(deta * deta + dphi * dphi);
433 if(dr <= rc){ // cell in cone
434 flagCell[ncell] = -1;
435 etCone+=etCell[ncell];
440 // select jets with et > background
441 // estimate max fluctuation of background in cone
442 Double_t ncellin = (Double_t)nCellIn;
443 Double_t ntcell = (Double_t)nCell;
444 Double_t etbmax = (etbgTotal + dEtTotal )*(ncellin/ntcell);
446 Double_t etcmin = etCone ; // could be used etCone - etmin !!
447 //desicions !! etbmax < etcmin
448 for(Int_t mcell =0; mcell < nCell; mcell++){
449 if(flagCell[mcell] == -1){
451 flagCell[mcell] = 1; //flag cell as used
453 flagCell[mcell] = 0; // leave it free
456 //store tmp jet info !!!
457 if(etbmax < etcmin) {
459 etaAlgoJet[nJets] = eta;
460 phiAlgoJet[nJets] = phi;
461 etAlgoJet[nJets] = etCone;
462 ncellsAlgoJet[nJets] = nCellIn;
466 AliError(Form("Too many jets (> %d) found by UA1JetFinder, adapt your cuts",kMaxJets));
470 } // end of cells loop
472 //reorder jets by et in cone
473 //sort jets by energy
475 TMath::Sort(nJets, etAlgoJet, idx); // sort only the found jets
476 for(Int_t p = 0; p < nJets; p++){
477 etaJet[p] = etaAlgoJet[idx[p]];
478 phiJet[p] = phiAlgoJet[idx[p]];
479 etJet[p] = etAlgoJet[idx[p]];
480 etallJet[p] = etAlgoJet[idx[p]];
481 ncellsJet[p] = ncellsAlgoJet[idx[p]];
485 ////////////////////////////////////////////////////////////////////////
487 void AliUA1JetFinderV1::SubtractBackg(const Int_t& nIn, const Int_t&nJ, Float_t&etbgTotalN,
488 const Float_t* ptT, const Float_t* etaT, const Float_t* phiT,
489 Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
490 Int_t* multJet, Int_t* injet)
492 //background subtraction using cone method but without correction in dE/deta distribution
494 //calculate energy inside and outside cones
495 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
496 Float_t rc= header->GetRadius();
497 Float_t etIn[kMaxJets];
499 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
500 // if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
501 for(Int_t ijet=0; ijet<nJ; ijet++){
502 Float_t deta = etaT[jpart] - etaJet[ijet];
503 Float_t dphi = phiT[jpart] - phiJet[ijet];
504 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
505 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
506 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
507 if(dr <= rc){ // particles inside this cone
510 if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
511 etIn[ijet] += ptT[jpart];
512 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet]+= ptT[jpart];
517 if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
518 etOut += ptT[jpart]; // particle outside cones and pt cut
519 } //end particle loop
521 //estimate jets and background areas
522 Float_t areaJet[kMaxJets];
523 Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
524 for(Int_t k=0; k<nJ; k++){
525 Float_t detamax = etaJet[k] + rc;
526 Float_t detamin = etaJet[k] - rc;
527 Float_t accmax = 0.0; Float_t accmin = 0.0;
528 if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
529 Float_t h = header->GetLegoEtaMax() - etaJet[k];
530 accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
532 if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
533 Float_t h = header->GetLegoEtaMax() + etaJet[k];
534 accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
536 areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
537 areaOut = areaOut - areaJet[k];
539 //subtract background using area method
540 for(Int_t ljet=0; ljet<nJ; ljet++){
541 Float_t areaRatio = areaJet[ljet]/areaOut;
542 etJet[ljet] = etIn[ljet]-etOut*areaRatio; // subtraction
545 // estimate new total background
546 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
547 etbgTotalN = etOut*areaT/areaOut;
552 ////////////////////////////////////////////////////////////////////////
554 void AliUA1JetFinderV1::SubtractBackgStat(const Int_t& nIn, const Int_t&nJ,Float_t&etbgTotalN,
555 const Float_t* ptT, const Float_t* etaT, const Float_t* phiT,
556 Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
557 Int_t* multJet, Int_t* injet)
560 //background subtraction using statistical method
561 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
562 Float_t etbgStat = header->GetBackgStat(); // pre-calculated background
564 //calculate energy inside
565 Float_t rc= header->GetRadius();
566 Float_t etIn[kMaxJets];
568 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
569 //if((fReader->GetCutFlag(jpart)) != 1) continue; // pt cut
570 for(Int_t ijet=0; ijet<nJ; ijet++){
571 Float_t deta = etaT[jpart] - etaJet[ijet];
572 Float_t dphi = phiT[jpart] - phiJet[ijet];
573 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
574 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
575 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
576 if(dr <= rc){ // particles inside this cone
579 if((fReader->GetCutFlag(jpart)) == 1){ // pt cut
580 etIn[ijet]+= ptT[jpart];
581 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
586 } //end particle loop
589 Float_t areaJet[kMaxJets];
590 Float_t areaOut = 4*(header->GetLegoEtaMax())*TMath::Pi();
591 for(Int_t k=0; k<nJ; k++){
592 Float_t detamax = etaJet[k] + rc;
593 Float_t detamin = etaJet[k] - rc;
594 Float_t accmax = 0.0; Float_t accmin = 0.0;
595 if(detamax > header->GetLegoEtaMax()){ // sector outside etamax
596 Float_t h = header->GetLegoEtaMax() - etaJet[k];
597 accmax = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
599 if(detamin < header->GetLegoEtaMin()){ // sector outside etamin
600 Float_t h = header->GetLegoEtaMax() + etaJet[k];
601 accmin = rc*rc*TMath::ACos(h/rc) - h*TMath::Sqrt(rc*rc - h*h);
603 areaJet[k] = rc*rc*TMath::Pi() - accmax - accmin;
606 //subtract background using area method
607 for(Int_t ljet=0; ljet<nJ; ljet++){
608 Float_t areaRatio = areaJet[ljet]/areaOut;
609 etJet[ljet] = etIn[ljet]-etbgStat*areaRatio; // subtraction
612 etbgTotalN = etbgStat;
616 ////////////////////////////////////////////////////////////////////////
618 void AliUA1JetFinderV1::SubtractBackgCone(const Int_t& nIn, const Int_t&nJ,Float_t& etbgTotalN,
619 const Float_t* ptT, const Float_t* etaT, const Float_t* phiT,
620 Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
621 Int_t* multJet, Int_t* injet)
623 // Cone background subtraction method taking into acount dEt/deta distribution
624 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
626 Float_t rc= header->GetRadius();
627 Float_t etamax = header->GetLegoEtaMax();
628 Float_t etamin = header->GetLegoEtaMin();
631 // jet energy and area arrays
632 char hEtname[256];char hAreaname[256];
634 for(Int_t mjet=0; mjet<nJ; mjet++){
636 sprintf(hEtname, "hEtJet%d", mjet);
637 fhEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
639 if(!fhAreaJet[mjet]){
640 sprintf(hAreaname, "hAreaJet%d", mjet);
641 fhAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
643 fhEtJet[mjet]->Reset();
644 fhAreaJet[mjet]->Reset();
646 // background energy and area
647 if(!fhEtBackg)fhEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
649 if(!fhAreaBackg) fhAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
650 fhAreaBackg->Reset();
653 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
654 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
655 Float_t deta = etaT[jpart] - etaJet[ijet];
656 Float_t dphi = phiT[jpart] - phiJet[ijet];
657 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
658 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
659 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
660 if(dr <= rc){ // particles inside this cone
663 if((fReader->GetCutFlag(jpart)) == 1){// pt cut
664 fhEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone
665 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
670 if(injet[jpart] == -1 && fReader->GetCutFlag(jpart) == 1)
671 fhEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
672 } //end particle loop
675 Float_t eta0 = etamin;
676 Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
677 Float_t eta1 = eta0 + etaw;
678 for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
679 Float_t etac = eta0 + etaw/2.0;
680 Float_t areabg = etaw*2.0*TMath::Pi();
681 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
682 Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
683 Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
684 Float_t acc0 = 0.0; Float_t acc1 = 0.0;
686 if(deta0 > rc && deta1 < rc){
687 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
690 if(deta0 < rc && deta1 > rc){
691 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
694 if(deta0 < rc && deta1 < rc){
695 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
696 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
697 if(eta1<etaJet[ijet]) areaj = acc1-acc0; // case 1
698 if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
699 if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
701 fhAreaJet[ijet]->Fill(etac,areaj);
702 areabg = areabg - areaj;
704 fhAreaBackg->Fill(etac,areabg);
707 } // end loop for all eta bins
709 //subtract background
710 for(Int_t kjet=0; kjet<nJ; kjet++){
711 etJet[kjet] = 0.0; // first clear etJet for this jet
712 for(Int_t bin = 0; bin< ndiv; bin++){
713 if(fhAreaJet[kjet]->GetBinContent(bin)){
714 Float_t areab = fhAreaBackg->GetBinContent(bin);
715 Float_t etb = fhEtBackg->GetBinContent(bin);
716 Float_t areaR = (fhAreaJet[kjet]->GetBinContent(bin))/areab;
717 etJet[kjet] = etJet[kjet] + ((fhEtJet[kjet]->GetBinContent(bin)) - etb*areaR); //subtraction
722 // calc background total
723 Double_t etOut = fhEtBackg->Integral();
724 Double_t areaOut = fhAreaBackg->Integral();
725 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
726 etbgTotalN = etOut*areaT/areaOut;
729 ////////////////////////////////////////////////////////////////////////
732 void AliUA1JetFinderV1::SubtractBackgRatio(const Int_t& nIn, const Int_t&nJ, Float_t& etbgTotalN,
733 const Float_t* ptT, const Float_t* etaT, const Float_t* phiT,
734 Float_t* etJet, const Float_t* etaJet, const Float_t* phiJet, Float_t* etsigJet,
735 Int_t* multJet, Int_t* injet)
737 // Ratio background subtraction method taking into acount dEt/deta distribution
738 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
739 //factor F calc before
740 Float_t bgRatioCut = header->GetBackgCutRatio();
744 Float_t rc= header->GetRadius();
745 Float_t etamax = header->GetLegoEtaMax();
746 Float_t etamin = header->GetLegoEtaMin();
749 // jet energy and area arrays
750 // jet energy and area arrays
751 char hEtname[256];char hAreaname[256];
753 for(Int_t mjet=0; mjet<nJ; mjet++){
755 sprintf(hEtname, "hEtJet%d", mjet);
756 fhEtJet[mjet] = new TH1F(hEtname,"et dist in eta ",ndiv,etamin,etamax);
758 if(!fhAreaJet[mjet]){
759 sprintf(hAreaname, "hAreaJet%d", mjet);
760 fhAreaJet[mjet] = new TH1F(hAreaname,"area dist in eta ",ndiv,etamin,etamax);
762 fhEtJet[mjet]->Reset();
763 fhAreaJet[mjet]->Reset();
765 // background energy and area
766 if(!fhEtBackg)fhEtBackg = new TH1F("hEtBackg"," backg et dist in eta ",ndiv,etamin,etamax);
768 if(!fhAreaBackg) fhAreaBackg = new TH1F("hAreaBackg","backg area dist in eta ",ndiv,etamin,etamax);
769 fhAreaBackg->Reset();
772 for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array
773 //if((fReader->GetCutFlag(jpart)) != 1) continue;
774 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
775 Float_t deta = etaT[jpart] - etaJet[ijet];
776 Float_t dphi = phiT[jpart] - phiJet[ijet];
777 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
778 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
779 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
780 if(dr <= rc){ // particles inside this cone
783 if((fReader->GetCutFlag(jpart)) == 1){ //pt cut
784 fhEtJet[ijet]->Fill(etaT[jpart],ptT[jpart]); //particle inside cone and pt cut
785 if(fReader->GetSignalFlag(jpart) == 1) etsigJet[ijet] += ptT[jpart];
790 if(injet[jpart] == -1) fhEtBackg->Fill(etaT[jpart],ptT[jpart]); // particle outside cones
791 } //end particle loop
794 Float_t eta0 = etamin;
795 Float_t etaw = (etamax - etamin)/((Float_t)ndiv);
796 Float_t eta1 = eta0 + etaw;
797 for(Int_t etabin = 0; etabin< ndiv; etabin++){ // loop for all eta bins
798 Float_t etac = eta0 + etaw/2.0;
799 Float_t areabg = etaw*2.0*TMath::Pi();
800 for(Int_t ijet=0; ijet<nJ; ijet++){ // loop for all jets
801 Float_t deta0 = TMath::Abs(eta0 - etaJet[ijet]);
802 Float_t deta1 = TMath::Abs(eta1 - etaJet[ijet]);
803 Float_t acc0 = 0.0; Float_t acc1 = 0.0;
805 if(deta0 > rc && deta1 < rc){
806 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
809 if(deta0 < rc && deta1 > rc){
810 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
813 if(deta0 < rc && deta1 < rc){
814 acc0 = rc*rc*TMath::ACos(deta0/rc) - deta0*TMath::Sqrt(rc*rc - deta0*deta0);
815 acc1 = rc*rc*TMath::ACos(deta1/rc) - deta1*TMath::Sqrt(rc*rc - deta1*deta1);
816 if(eta1<etaJet[ijet]) areaj = acc1-acc0; // case 1
817 if((eta0 < etaJet[ijet]) && (etaJet[ijet]<eta1)) areaj = rc*rc*TMath::Pi() - acc1 -acc0; // case 2
818 if(etaJet[ijet] < eta0) areaj = acc0 -acc1; // case 3
820 fhAreaJet[ijet]->Fill(etac,areaj);
821 areabg = areabg - areaj;
823 fhAreaBackg->Fill(etac,areabg);
826 } // end loop for all eta bins
828 //subtract background
829 for(Int_t kjet=0; kjet<nJ; kjet++){
830 etJet[kjet] = 0.0; // first clear etJet for this jet
831 for(Int_t bin = 0; bin< ndiv; bin++){
832 if(fhAreaJet[kjet]->GetBinContent(bin)){
833 Float_t areab = fhAreaBackg->GetBinContent(bin);
834 Float_t etb = fhEtBackg->GetBinContent(bin);
835 Float_t areaR = (fhAreaJet[kjet]->GetBinContent(bin))/areab;
836 etJet[kjet] = etJet[kjet] + ((fhEtJet[kjet]->GetBinContent(bin)) - etb*areaR*bgRatioCut); //subtraction
841 // calc background total
842 Double_t etOut = fhEtBackg->Integral();
843 Double_t areaOut = fhAreaBackg->Integral();
844 Float_t areaT = 4*(header->GetLegoEtaMax())*TMath::Pi();
845 etbgTotalN = etOut*areaT/areaOut;
848 ////////////////////////////////////////////////////////////////////////
851 void AliUA1JetFinderV1::Reset()
854 AliJetFinder::Reset();
857 ////////////////////////////////////////////////////////////////////////
859 void AliUA1JetFinderV1::WriteJHeaderToFile() const
861 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
865 ////////////////////////////////////////////////////////////////////////
867 void AliUA1JetFinderV1::Init()
869 // initializes some variables
870 AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
872 TH2F("legoH","eta-phi",
873 header->GetLegoNbinEta(), header->GetLegoEtaMin(),
874 header->GetLegoEtaMax(), header->GetLegoNbinPhi(),
875 header->GetLegoPhiMin(), header->GetLegoPhiMax());
876 // Do not store in current dir
877 fLego->SetDirectory(0);