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