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