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 //*-- Authors: Andreas Morsch (CERN)
20 //* Aleksei Pavlinov (WSU)
29 #include <TBranchElement.h>
31 #include <TClonesArray.h>
38 #include <TParticle.h>
39 #include <TParticlePDG.h>
40 #include <TPaveText.h>
41 #include <TPythia6Calls.h>
49 #include "AliEMCALDigit.h"
50 #include "AliEMCALDigitizer.h"
51 #include "AliEMCALFast.h"
52 #include "AliEMCALGeometry.h"
53 #include "AliEMCALHadronCorrection.h"
54 #include "AliEMCALHit.h"
55 #include "AliEMCALJetFinder.h"
56 #include "AliEMCALJetMicroDst.h"
57 #include "AliHeader.h"
59 #include "AliMagFCM.h"
61 #include "AliGenerator.h"
62 #include "AliEMCALGetter.h"
63 // Interface to FORTRAN
68 ClassImp(AliEMCALJetFinder)
70 //____________________________________________________________________________
71 AliEMCALJetFinder::AliEMCALJetFinder()
73 // Default constructor
102 SetParametersForBgSubtraction();
105 AliEMCALJetFinder::AliEMCALJetFinder(const char* name, const char *title)
109 // Title is used in method GetFileNameForParameters();
111 fJets = new TClonesArray("AliEMCALJet",10000);
113 for (Int_t i = 0; i < 30000; i++)
135 fHadronCorrector = 0;
144 SetMomentumSmearing();
147 SetHadronCorrection();
151 SetParametersForBgSubtraction();
154 void AliEMCALJetFinder::SetParametersForBgSubtraction
155 (Int_t mode, Float_t minMove, Float_t maxMove, Float_t precBg)
157 // see file /home/pavlinov/cosmic/pythia/jetFinderParamData.inc
158 // at WSU Linux cluster - 11-feb-2002
159 // These parameters must be tuned more carefull !!!
166 //____________________________________________________________________________
167 AliEMCALJetFinder::~AliEMCALJetFinder()
179 delete fhLegoHadrCorr;
182 delete fhCellEMCALEt;
184 delete fhTrackPtBcut;
185 delete fhChPartMultInTpc;
193 delete[] fTrackListB;
201 # define jet_finder_ua1 jet_finder_ua1_
203 # define type_of_call
206 # define jet_finder_ua1 JET_FINDER_UA1
208 # define type_of_call _stdcall
211 extern "C" void type_of_call
212 jet_finder_ua1(Int_t& ncell, Int_t& ncell_tot,
213 Float_t etc[30000], Float_t etac[30000],
215 Float_t& min_move, Float_t& max_move, Int_t& mode,
216 Float_t& prec_bg, Int_t& ierror);
218 extern "C" void type_of_call hf1(Int_t& id, Float_t& x, Float_t& wgt);
221 void AliEMCALJetFinder::Init()
224 // Geometry and I/O initialization
228 // Get geometry parameters from EMCAL
232 //AliEMCAL* pEMCAL = (AliEMCAL*) gAlice->GetModule("EMCAL");
233 // AliEMCALGeometry* geom =
234 // AliEMCALGeometry::GetInstance(pEMCAL->GetTitle(), "");
235 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
236 AliEMCALGeometry* geom = gime->EMCALGeometry() ;
238 // SetSamplingFraction(geom->GetSampling());
240 fNbinEta = geom->GetNZ();
241 fNbinPhi = geom->GetNPhi();
242 fPhiMin = geom->GetArm1PhiMin()*TMath::Pi()/180.;
243 fPhiMax = geom->GetArm1PhiMax()*TMath::Pi()/180.;
244 fEtaMin = geom->GetArm1EtaMin();
245 fEtaMax = geom->GetArm1EtaMax();
246 fDphi = (fPhiMax-fPhiMin)/fNbinPhi;
247 fDeta = (fEtaMax-fEtaMin)/fNbinEta;
248 fNtot = fNbinPhi*fNbinEta;
249 fWeightingMethod = kFALSE;
252 SetCellSize(fDeta, fDphi);
255 if (fOutFileName) fOutFile = new TFile(fOutFileName, "recreate");
260 void AliEMCALJetFinder::Find(Int_t ncell, Int_t ncelltot, Float_t etc[30000],
261 Float_t etac[30000], Float_t phic[30000],
262 Float_t minmove, Float_t maxmove, Int_t mode,
263 Float_t precbg, Int_t ierror)
265 // Wrapper for fortran coded jet finder
266 // Full argument list
267 jet_finder_ua1(ncell, ncelltot, etc, etac, phic,
268 minmove, maxmove, mode, precbg, ierror);
269 // Write result to output
270 if(fWrite) WriteJets();
274 void AliEMCALJetFinder::Find()
276 // Wrapper for fortran coded jet finder using member data for
279 Float_t minmove = fMinMove;
280 Float_t maxmove = fMaxMove;
282 Float_t precbg = fPrecBg;
284 ResetJets(); // 4-feb-2002 by PAI
286 jet_finder_ua1(fNcell, fNtot, fEtCell, fEtaCell, fPhiCell,
287 minmove, maxmove, mode, precbg, ierror);
289 // Write result to output
290 Int_t njet = Njets();
292 for (Int_t nj=0; nj<njet; nj++)
294 if (fWeightingMethod)
296 fJetT[nj] = new AliEMCALJet(WeightedJetEnergy( JetEtaW(nj),JetPhiW(nj) ),
302 fJetT[nj] = new AliEMCALJet(JetEnergy(nj),
306 fJetT[nj]->SetIsWeightedEnergy(fWeightingMethod);
307 fJetT[nj]->SetEMCALEnergy( EMCALConeEnergy(JetEtaW(nj),JetPhiW(nj)) );
308 fJetT[nj]->SetTrackEnergy(TrackConeEnergy( JetEtaW(nj),JetPhiW(nj) ));
312 FindTracksInJetCone();
313 if(fWrite) WriteJets();
318 Float_t AliEMCALJetFinder::EMCALConeEnergy(Float_t eta, Float_t phi)
320 // Calculate the energy in the cone
321 Float_t newenergy = 0.0;
322 Float_t bineta,binphi;
323 TAxis *x = fhLegoEMCAL->GetXaxis();
324 TAxis *y = fhLegoEMCAL->GetYaxis();
325 for (Int_t i = 0 ; i < fNbinEta ; i++) // x coord
327 for (Int_t j = 0 ; j < fNbinPhi ; j++) // y coord
329 bineta = x->GetBinCenter(i);
330 binphi = y->GetBinCenter(j);
331 if ( (bineta-eta)*(bineta-eta) + (binphi-phi)*(binphi-phi) < fConeRadius*fConeRadius)
333 newenergy += fhLegoEMCAL->GetBinContent(i,j);
340 Float_t AliEMCALJetFinder::TrackConeEnergy(Float_t eta, Float_t phi)
342 // Calculate the track energy in the cone
343 Float_t newenergy = 0.0;
344 Float_t bineta,binphi;
345 TAxis *x = fhLegoTracks->GetXaxis();
346 TAxis *y = fhLegoTracks->GetYaxis();
347 for (Int_t i = 0 ; i < fNbinEta ; i++) // x coord
349 for (Int_t j = 0 ; j < fNbinPhi ; j++) // y coord
351 bineta = x->GetBinCenter(i);
352 binphi = y->GetBinCenter(j);
353 if ( (bineta-eta)*(bineta-eta) + (binphi-phi)*(binphi-phi) < fConeRadius*fConeRadius)
355 newenergy += fhLegoTracks->GetBinContent(i,j);
363 Float_t AliEMCALJetFinder::WeightedJetEnergy(Float_t eta, Float_t phi)
365 // Calculate the weighted jet energy
367 Float_t newenergy = 0.0;
368 Float_t bineta,binphi;
369 TAxis *x = fhLegoEMCAL->GetXaxis();
370 TAxis *y = fhLegoEMCAL->GetYaxis();
373 for (Int_t i = 0 ; i < fNbinEta ; i++) // x coord
375 for (Int_t j = 0 ; j < fNbinPhi ; j++) // y coord
377 bineta = x->GetBinCenter(i);
378 binphi = y->GetBinCenter(j);
379 if ( (bineta-eta)*(bineta-eta) + (binphi-phi)*(binphi-phi) < fConeRadius*fConeRadius)
381 newenergy += (fEMCALWeight)* fhLegoEMCAL->GetBinContent(i,j) + (fTrackWeight)* fhLegoTracks->GetBinContent(i,j);
391 Int_t AliEMCALJetFinder::Njets() const
393 // Get number of reconstructed jets
394 return EMCALJETS.njet;
397 Float_t AliEMCALJetFinder::JetEnergy(Int_t i) const
399 // Get reconstructed jet energy
400 return EMCALJETS.etj[i];
403 Float_t AliEMCALJetFinder::JetPhiL(Int_t i) const
405 // Get reconstructed jet phi from leading particle
406 return EMCALJETS.phij[0][i];
409 Float_t AliEMCALJetFinder::JetPhiW(Int_t i) const
411 // Get reconstructed jet phi from weighting
412 return EMCALJETS.phij[1][i];
415 Float_t AliEMCALJetFinder::JetEtaL(Int_t i) const
417 // Get reconstructed jet eta from leading particles
418 return EMCALJETS.etaj[0][i];
422 Float_t AliEMCALJetFinder::JetEtaW(Int_t i) const
424 // Get reconstructed jet eta from weighting
425 return EMCALJETS.etaj[1][i];
428 void AliEMCALJetFinder::SetCellSize(Float_t eta, Float_t phi)
430 // Set grid cell size
431 EMCALCELLGEO.etaCellSize = eta;
432 EMCALCELLGEO.phiCellSize = phi;
435 void AliEMCALJetFinder::SetConeRadius(Float_t par)
437 // Set jet cone radius
438 EMCALJETPARAM.coneRad = par;
442 void AliEMCALJetFinder::SetEtSeed(Float_t par)
444 // Set et cut for seeds
445 EMCALJETPARAM.etSeed = par;
449 void AliEMCALJetFinder::SetMinJetEt(Float_t par)
451 // Set minimum jet et
452 EMCALJETPARAM.ejMin = par;
456 void AliEMCALJetFinder::SetMinCellEt(Float_t par)
458 // Set et cut per cell
459 EMCALJETPARAM.etMin = par;
463 void AliEMCALJetFinder::SetPtCut(Float_t par)
465 // Set pT cut on charged tracks
470 void AliEMCALJetFinder::Test()
473 // Test the finder call
475 const Int_t knmax = 30000;
477 Int_t ncelltot = 100;
489 Find(ncell, ncelltot, etc, etac, phic,
490 minmove, maxmove, mode, precbg, ierror);
498 void AliEMCALJetFinder::AddJet(const AliEMCALJet& jet)
503 TClonesArray &lrawcl = *fJets;
504 new(lrawcl[fNjets++]) AliEMCALJet(jet);
507 void AliEMCALJetFinder::ResetJets()
516 void AliEMCALJetFinder::WriteJets()
519 // Add all jets to the list
521 const Int_t kBufferSize = 4000;
522 const char* file = 0;
524 Int_t njet = Njets();
526 for (Int_t nj = 0; nj < njet; nj++)
533 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
536 // output written to input file
538 AliEMCAL* pEMCAL = (AliEMCAL* )gAlice->GetModule("EMCAL");
539 TTree* pK = gAlice->TreeK();
540 file = (pK->GetCurrentFile())->GetName();
541 TBranch * jetBranch ;
543 printf("Make Branch - TreeR address %p %p\n",gAlice->TreeR(), pEMCAL);
544 //if (fJets && gAlice->TreeR()) {
545 if (fJets && gime->TreeR()) {
546 // pEMCAL->MakeBranchInTree(gAlice->TreeR(),
547 jetBranch = gime->TreeR()->Branch("EMCALJets", &fJets, kBufferSize, 0) ;
548 //pEMCAL->MakeBranchInTree(gime->TreeR(),
554 //Int_t nev = gAlice->GetHeader()->GetEvent();
555 //gAlice->TreeR()->Fill();
558 //sprintf(hname,"TreeR%d", nev);
559 //gAlice->TreeR()->Write(hname);
560 //gAlice->TreeR()->Reset();
561 gime->WriteRecPoints("OVERWRITE");
565 // Output written to user specified output file
567 //TTree* pK = gAlice->TreeK();
568 TTree* pK = gAlice->TreeK();
569 fInFile = pK->GetCurrentFile();
573 sprintf(hname,"TreeR%d", fEvent);
574 TTree* treeJ = new TTree(hname, "EMCALJets");
575 treeJ->Branch("EMCALJets", &fJets, kBufferSize);
583 void AliEMCALJetFinder::BookLego()
586 // Book histo for discretization
590 // Don't add histos to the current directory
591 if(fDebug) printf("\n AliEMCALJetFinder::BookLego() \n");
593 // TH2::AddDirectory(0); // hists wil be put to the list from gROOT
594 // TH1::AddDirectory(0);
598 fLego = new TH2F("legoH","eta-phi",
599 fNbinEta, fEtaMin, fEtaMax,
600 fNbinPhi, fPhiMin, fPhiMax);
603 fLegoB = new TH2F("legoB","eta-phi for BG event",
604 fNbinEta, fEtaMin, fEtaMax,
605 fNbinPhi, fPhiMin, fPhiMax);
608 fhLegoTracks = new TH2F("hLegoTracks","eta-phi for Tracks",
609 fNbinEta, fEtaMin, fEtaMax, fNbinPhi, fPhiMin, fPhiMax);
611 fhLegoEMCAL = new TH2F("hLegoEMCAL","eta-phi for EMCAL",
612 fNbinEta, fEtaMin, fEtaMax, fNbinPhi, fPhiMin, fPhiMax);
613 // Hadron correction map
614 fhLegoHadrCorr = new TH2F("hLegoHadrCorr","eta-phi for Hadron. Correction",
615 fNbinEta, fEtaMin, fEtaMax, fNbinPhi, fPhiMin, fPhiMax);
616 // Hists. for tuning jet finder
617 fhEff = new TH2F("hEff","#epsilon vs momentum ", 100,0.,20., 50,0.5,1.);
621 for(Int_t i=1; i<=1000; i++) eTmp[i] = 0.1*i; // step 100 mev
622 for(Int_t i=1001; i<=1100; i++) eTmp[i] = eTmp[1000] + 1.0*(i-1000); // step 1GeV
624 fhCellEt = new TH1F("hCellEt","Cell E_{T} from fLego",
625 eTmp.GetSize()-1, eTmp.GetArray());
626 fhCellEMCALEt = new TH1F("hCellEMCALEt","Cell E_{T} for EMCAL itself",
627 eTmp.GetSize()-1, eTmp.GetArray());
628 fhTrackPt = new TH1F("hTrackPt","Ch.particles P_{T} ",
629 eTmp.GetSize()-1, eTmp.GetArray());
630 fhTrackPtBcut = new TH1F("hTrackPtBcut","Ch.particles P_{T} + magnetic field cut",
631 eTmp.GetSize()-1, eTmp.GetArray());
633 fhChPartMultInTpc = new TH1F("hChPartMultInTpc",
634 "Charge partilcle multiplicity in |%eta|<0.9", 2000, 0, 20000);
636 fhSinTheta = new TH1F("fhSinTheta","sin(theta)", fNbinEta, fEtaMin, fEtaMax);
637 TAxis *xax = fhSinTheta->GetXaxis();
638 for(Int_t i=1; i<=fNbinEta; i++) {
639 Double_t eta = xax->GetBinCenter(i);
640 fhSinTheta->Fill(eta, 1./TMath::CosH(eta)); // cosh(eta) = 1./sin(theta)
643 //! first canvas for drawing
644 fHistsList=AliEMCALJetMicroDst::MoveHistsToList("Hists from AliEMCALJetFinder", kFALSE);
647 void AliEMCALJetFinder::DumpLego()
650 // Dump lego histo into array
653 TAxis* xaxis = fLego->GetXaxis();
654 TAxis* yaxis = fLego->GetYaxis();
655 // fhCellEt->Clear();
657 for (Int_t i = 1; i <= fNbinEta; i++) {
658 for (Int_t j = 1; j <= fNbinPhi; j++) {
660 e = fLego->GetBinContent(i,j);
662 if (gRandom->Rndm() < 0.5) {
663 Float_t ebg = 0.28 * TMath::Abs(gRandom->Gaus(0.,1.));
667 if (e > 0.0) e -= fMinCellEt;
669 Float_t eta = xaxis->GetBinCenter(i);
670 Float_t phi = yaxis->GetBinCenter(j);
672 fEtaCell[fNcell] = eta;
673 fPhiCell[fNcell] = phi;
677 eH = fhLegoEMCAL->GetBinContent(i,j);
678 if(eH > 0.0) fhCellEMCALEt->Fill(eH);
686 void AliEMCALJetFinder::ResetMap()
689 // Reset eta-phi array
691 for (Int_t i=0; i<30000; i++)
700 void AliEMCALJetFinder::FillFromTracks(Int_t flag, Int_t ich)
704 const char* name = gAlice->Generator()->GetName();
705 enum {kPythia, kHijing, kHijingPara};
708 if (!strcmp(name, "Hijing")){
710 } else if (!strcmp(name, "Pythia")) {
712 } else if (!strcmp(name, "HIJINGpara") ||!strcmp(name, "HIGINGpara")) {
713 genType = kHijingPara;
716 printf("Fill tracks generated by %s %d\n", name, genType);
720 // Fill Cells with track information
723 printf("\n AliEMCALJetFinder::FillFromTracks(%i,%i) ",flag,ich);
728 if (!fLego) BookLego();
730 if (flag == 0) fLego->Reset();
732 // Access particle information
733 Int_t npart = (gAlice->GetHeader())->GetNprimary();
734 Int_t ntr = (gAlice->GetHeader())->GetNtrack();
735 printf(" : #primary particles %i # tracks %i : (before) Sum.Et %f\n",
736 npart, ntr, fLego->Integral());
741 // 1: selected for jet finding
744 if (fTrackList) delete[] fTrackList;
745 if (fPtT) delete[] fPtT;
746 if (fEtaT) delete[] fEtaT;
747 if (fPhiT) delete[] fPhiT;
748 if (fPdgT) delete[] fPdgT;
750 fTrackList = new Int_t [npart];
751 fPtT = new Float_t[npart];
752 fEtaT = new Float_t[npart];
753 fPhiT = new Float_t[npart];
754 fPdgT = new Int_t[npart];
758 Float_t chTmp=0.0; // charge of current particle
761 // this is for Pythia ??
762 for (Int_t part = 0; part < npart; part++) {
763 TParticle *mPart = gAlice->GetMCApp()->Particle(part);
764 Int_t mpart = mPart->GetPdgCode();
765 Int_t child1 = mPart->GetFirstDaughter();
766 Float_t pT = mPart->Pt();
767 Float_t p = mPart->P();
768 Float_t phi = mPart->Phi();
770 if(pT > 0.001) eta = mPart->Eta();
771 Float_t theta = mPart->Theta();
772 if (fDebug>=2 && mPart->GetStatusCode()==1) {
773 printf("ind %7i part %7i -> child1 %5i child2 %5i Status %5i\n",
774 part, mpart, child1, mPart->GetLastDaughter(), mPart->GetStatusCode());
777 if (fDebug >= 2 && genType == kPythia) {
778 if (part == 6 || part == 7)
780 printf("\n Simulated Jet (pt, eta, phi): %d %f %f %f\n",
781 part-5, pT, eta, phi);
785 fTrackList[part] = 0;
786 fPtT[part] = pT; // must be change after correction for resolution !!!
791 // final state particles only
793 if (genType == kPythia) {
794 if (mPart->GetStatusCode() != 1) continue;
795 } else if (genType == kHijing) {
796 if (child1 >= 0 && child1 < npart) continue;
800 TParticlePDG* pdgP = 0;
801 // charged or neutral
802 pdgP = mPart->GetPDG();
803 chTmp = pdgP->Charge() / 3.; // 13-feb-2001!!
810 if (mpart != kNeutron &&
811 mpart != kNeutronBar &&
812 mpart != kK0Long) continue;
815 } else if (ich == 2) {
816 if (mpart == kNeutron ||
817 mpart == kNeutronBar ||
818 mpart == kK0Long) continue;
821 if (TMath::Abs(eta)<=0.9) fNChTpc++;
823 if (child1 >= 0 && child1 < npart) continue;
825 if (eta > fEtaMax || eta < fEtaMin) continue;
826 if (phi > fPhiMax || phi < fPhiMin ) continue;
829 printf("\n=>nsel:%5d mpart %5d child1 %5d eta %6.2f phi %6.2f pT %6.2f ch %3.0f ",
830 part, mpart, child1, eta, phi, pT, chTmp);
833 // Momentum smearing goes here ...
837 if (fSmear && TMath::Abs(chTmp)) {
838 pw = AliEMCALFast::SmearMomentum(1,p);
839 // p changed - take into account when calculate pT,
842 if(fDebug >= 4) printf("\n Smearing : p %8.4f change to %8.4f ", p, pw);
846 // Tracking Efficiency and TPC acceptance goes here ...
848 if (fEffic && TMath::Abs(chTmp)) {
849 eff = 0.9; // AliEMCALFast::Efficiency(2,p);
850 if(fhEff) fhEff->Fill(p, eff);
851 if (AliEMCALFast::RandomReject(eff)) {
852 if(fDebug >= 5) printf(" reject due to unefficiency ");
857 // Correction of Hadronic Energy goes here ...
860 // phi propagation for hadronic correction
862 Bool_t curls = kFALSE; // hit two the EMCAL (no curl)
863 Float_t phiHC=0.0, dpH=0.0, dphi=0.0, eTdpH=0;
864 if(TMath::Abs(chTmp)) {
865 // hadr. correction only for charge particle
866 dphi = PropagatePhi(pT, chTmp, curls);
869 printf("\n Delta phi %f pT %f ", dphi, pT);
870 if (curls) printf("\n !! Track is curling");
872 if(!curls) fhTrackPtBcut->Fill(pT);
874 if (fHCorrection && !curls) {
875 if (!fHadronCorrector)
876 Fatal("AliEMCALJetFinder",
877 "Hadronic energy correction required but not defined !");
879 dpH = fHadronCorrector->GetEnergy(p, eta, 7);
880 eTdpH = dpH*TMath::Sin(theta);
882 if (fDebug >= 7) printf(" phi %f phiHC %f eTcorr %f\n",
883 phi, phiHC, -eTdpH); // correction is negative
885 xbin = fLego->GetXaxis()->FindBin(eta);
886 ybin = fLego->GetYaxis()->FindBin(phiHC);
887 cout <<"Hadron Correction affected bin - contents before correction : "<<fLego->GetBinContent(xbin,ybin)<<endl;
888 fLego->Fill(eta, phi, -fSamplingF*eTdpH );
889 cout <<"Hadron Correction affected bin - contents after correction : "<<fLego->GetBinContent(xbin,ybin)<<endl;
890 fhLegoHadrCorr->Fill(eta, phi, fSamplingF*eTdpH);
894 // More to do ? Just think about it !
896 if (phi > fPhiMax || phi < fPhiMin ) continue;
898 if(TMath::Abs(chTmp) ) { // charge particle
899 if (pT > fPtCut && !curls) {
900 if (fDebug >= 8) printf("Charge : fLego->Fill(%5.2f, %5.2f, %6.2f, %d)\n",
901 eta , phi, pT, fNtS);
902 fLego->Fill(eta, phi, pT);
903 fhLegoTracks->Fill(eta, phi, pT); // 20-feb for checking
904 fTrackList[part] = 1;
907 } else if(ich > 0 || fK0N) {
908 // case of n, nbar and K0L
909 if (fDebug >= 9) printf("Neutral : fLego->Fill(%5.2f, %5.2f, %6.2f, %d)\n",
910 eta , phi, pT, fNtS);
911 fLego->Fill(eta, phi, pT);
912 fTrackList[part] = 1;
917 for(Int_t i=0; i<fLego->GetSize(); i++) {
918 Float_t etc = (*fLego)[i];
919 if (etc > fMinCellEt) etsum += etc;
922 printf("FillFromTracks: Sum above threshold %f -> %f (%f)\n", fMinCellEt, etsum, fLego->Integral());
923 printf(" Track selected(fNtS) %i \n", fNtS);
928 void AliEMCALJetFinder::FillFromHits(Int_t flag)
931 // Fill Cells with hit information
935 printf("\n AliEMCALJetFinder::FillFromHits(%i)\n",flag);
939 if (!fLego) BookLego();
940 // Reset eta-phi maps if needed
941 if (flag == 0) { // default behavior
943 fhLegoTracks->Reset();
944 fhLegoEMCAL->Reset();
945 fhLegoHadrCorr->Reset();
947 // Initialize from background event if available
949 // Access hit information
950 AliEMCAL* pEMCAL = (AliEMCAL*) gAlice->GetModule("EMCAL");
951 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
952 TTree *treeH = gime->TreeH();
953 Int_t ntracks = (Int_t) treeH->GetEntries();
958 // Double_t etH = 0.0;
960 for (Int_t track=0; track<ntracks;track++) {
962 nbytes += treeH->GetEvent(track);
966 for(AliEMCALHit* mHit=(AliEMCALHit*) pEMCAL->FirstHit(-1);
968 mHit=(AliEMCALHit*) pEMCAL->NextHit())
970 Float_t x = mHit->X(); // x-pos of hit
971 Float_t y = mHit->Y(); // y-pos
972 Float_t z = mHit->Z(); // z-pos
973 Float_t eloss = mHit->GetEnergy(); // deposited energy
975 Float_t r = TMath::Sqrt(x*x+y*y);
976 Float_t theta = TMath::ATan2(r,z);
977 Float_t eta = -TMath::Log(TMath::Tan(theta/2.));
978 Float_t phi = TMath::ATan2(y,x);
980 if (fDebug >= 21) printf("\n Hit %f %f %f %f %f %f %f %f", x, y, z, eloss, r, eta, phi, fSamplingF);
982 // etH = fSamplingF*eloss*TMath::Sin(theta);
983 fLego->Fill(eta, phi, eloss);
987 // Transition from deposit energy to eT (eT = de*SF*sin(theta))
989 for(Int_t i=1; i<=fLego->GetNbinsX(); i++){ // eta
990 Double_t sinTheta = fhSinTheta->GetBinContent(i), eT=0;
991 for(Int_t j=1; j<=fLego->GetNbinsY(); j++){ // phi
992 eT = fLego->GetBinContent(i,j)*fSamplingF*sinTheta;
993 fLego->SetBinContent(i,j,eT);
994 // copy content of fLego to fhLegoEMCAL (fLego and fhLegoEMCAL are identical)
995 fhLegoEMCAL->SetBinContent(i,j,eT);
996 if (eT > fMinCellEt) etsum += eT;
1000 // for(Int_t i=0; i<fLego->GetSize(); i++) {
1001 // (*fhLegoEMCAL)[i] = (*fLego)[i];
1002 // Float_t etc = (*fLego)[i];
1003 // if (etc > fMinCellEt) etsum += etc;
1006 printf("FillFromHits: Sum above threshold %f -> %f \n ", fMinCellEt, etsum);
1010 void AliEMCALJetFinder::FillFromDigits(Int_t flag)
1013 // Fill Cells with digit information
1018 if (!fLego) BookLego();
1019 if (flag == 0) fLego->Reset();
1026 TClonesArray* digs = new TClonesArray("AliEMCALDigit", 10000);
1027 TTree *treeD = gAlice->TreeD();
1028 TBranchElement* branchDg = (TBranchElement*)
1029 treeD->GetBranch("EMCAL");
1031 if (!branchDg) Fatal("AliEMCALJetFinder",
1032 "Reading digits requested but no digits in file !");
1034 branchDg->SetAddress(&digs);
1035 Int_t nent = (Int_t) branchDg->GetEntries();
1037 // Connect digitizer
1039 AliEMCALDigitizer* digr = new AliEMCALDigitizer();
1040 TBranchElement* branchDr = (TBranchElement*)
1041 treeD->GetBranch("AliEMCALDigitizer");
1042 branchDr->SetAddress(&digr);
1045 nbytes += branchDg->GetEntry(0);
1046 nbytes += branchDr->GetEntry(0);
1048 // Get digitizer parameters
1049 Float_t ecADCped = digr->GetECApedestal();
1050 Float_t ecADCcha = digr->GetECAchannel();
1052 AliEMCAL* pEMCAL = (AliEMCAL*) gAlice->GetModule("EMCAL");
1053 AliEMCALGeometry* geom =
1054 AliEMCALGeometry::GetInstance(pEMCAL->GetTitle(), "");
1057 Int_t ndig = digs->GetEntries();
1058 Info("FillFromDigits","Number of Digits: %d %d\n Parameters: EC: %f %f\n Geometry: %d %d",
1059 ndig, nent, ecADCped, ecADCcha, geom->GetNEta(), geom->GetNPhi());
1066 while ((sdg = (AliEMCALDigit*)(next())))
1068 Double_t pedestal = 0.;
1069 Double_t channel = 0.;
1070 if (geom->IsInECA(sdg->GetId())) {
1071 pedestal = ecADCped;
1075 Fatal("FillFromDigits", "unexpected digit is number!") ;
1077 Float_t eta = sdg->GetEta();
1078 Float_t phi = sdg->GetPhi() * TMath::Pi()/180.;
1079 Float_t amp = (Float_t) (channel*(sdg->GetAmp())-pedestal);
1082 Info("FillFromDigits", "Digit: eta %8.3f phi %8.3f amp %8.3f %8d",
1083 eta, phi, amp, sdg->GetAmp());
1085 fLego->Fill(eta, phi, fSamplingF*amp);
1093 void AliEMCALJetFinder::FillFromHitFlaggedTracks(Int_t flag)
1096 // Fill Cells with hit information
1101 if (!fLego) BookLego();
1103 if (flag == 0) fLego->Reset();
1105 // Flag charged tracks passing through TPC which
1106 // are associated to EMCAL Hits
1107 BuildTrackFlagTable();
1110 // Access particle information
1111 TTree *treeK = gAlice->TreeK();
1112 Int_t ntracks = (Int_t) treeK->GetEntries();
1114 if (fPtT) delete[] fPtT;
1115 if (fEtaT) delete[] fEtaT;
1116 if (fPhiT) delete[] fPhiT;
1117 if (fPdgT) delete[] fPdgT;
1119 fPtT = new Float_t[ntracks];
1120 fEtaT = new Float_t[ntracks];
1121 fPhiT = new Float_t[ntracks];
1122 fPdgT = new Int_t[ntracks];
1127 for (Int_t track = 0; track < ntracks; track++) {
1128 TParticle *mPart = gAlice->GetMCApp()->Particle(track);
1129 Float_t pT = mPart->Pt();
1130 Float_t phi = mPart->Phi();
1131 Float_t eta = mPart->Eta();
1133 if(fTrackList[track]) {
1137 fPdgT[track] = mPart->GetPdgCode();
1139 if (track < 2) continue; //Colliding particles?
1140 if (pT == 0 || pT < fPtCut) continue;
1142 fLego->Fill(eta, phi, pT);
1148 void AliEMCALJetFinder::FillFromParticles()
1150 // 26-feb-2002 PAI - for checking all chain
1151 // Work on particles level; accept all particle (not neutrino )
1153 Double_t pX=0, pY=0, pZ=0, energy=0; // checking conservation law
1157 if (!fLego) BookLego();
1160 // Access particles information
1161 Int_t npart = (gAlice->GetHeader())->GetNprimary();
1162 if (fDebug >= 2 || npart<=0) {
1163 printf(" AliEMCALJetFinder::FillFromParticles : npart %i\n", npart);
1164 if(npart<=0) return;
1168 RearrangeParticlesMemory(npart);
1170 // Go through the particles
1171 Int_t mpart, child1, child2, geantPdg;
1172 Float_t pT, phi, eta, e=0, px=0, py=0, pz=0;
1174 for (Int_t part = 0; part < npart; part++) {
1176 fTrackList[part] = 0;
1178 mPart = gAlice->GetMCApp()->Particle(part);
1179 mpart = mPart->GetPdgCode();
1180 child1 = mPart->GetFirstDaughter();
1181 child2 = mPart->GetLastDaughter();
1189 e = mPart->Energy();
1191 // see pyedit in Pythia's text
1193 if (IsThisPartonsOrDiQuark(mpart)) continue;
1194 printf("%5i: %5i(%2i) px %5.1f py %5.1f pz %6.1f e %6.1f childs %5i,%5i \n",
1195 part, mpart, geantPdg, px, py, pz, e, child1, child2);
1197 // exclude partons (21 - gluon, 92 - string)
1200 // exclude neutrinous also ??
1201 if (fDebug >= 11 && pT>0.01)
1202 printf("\n part:%5d mpart %5d eta %9.2f phi %9.2f pT %9.2f ",
1203 part, mpart, eta, phi, pT);
1208 fPdgT[part] = mpart;
1212 if (child1 >= 0 && child1 < npart) continue;
1214 // printf("%4i -> %5i(%3i) px %6.1f py %6.1f pz %7.1f e %8.2f child1 %5i %s\n",
1215 // part, mpart, geantPdg, px, py, pz, e, child1, name.Data());
1222 if (TMath::Abs(eta) <= 0.9) fNChTpc++;
1224 if (eta > fEtaMax || eta < fEtaMin) continue;
1225 if (phi > fPhiMax || phi < fPhiMin ) continue;
1227 if(fK0N==0 ) { // exclude neutral hadrons
1228 if (mpart == kNeutron || mpart == kNeutronBar || mpart == kK0Long) continue;
1230 fTrackList[part] = 1;
1231 fLego->Fill(eta, phi, pT);
1234 printf("\n PX %8.2f PY %8.2f PZ %8.2f E %8.2f \n",
1235 pX, pY, pZ, energy);
1237 if(fhChPartMultInTpc) fhChPartMultInTpc->Fill(fNChTpc);
1240 void AliEMCALJetFinder::FillFromPartons()
1242 // function under construction - 13-feb-2002 PAI
1245 printf("\n AliEMCALJetFinder::FillFromPartons()\n");
1249 if (!fLego) BookLego();
1252 // Access particle information
1253 Int_t npart = (gAlice->GetHeader())->GetNprimary();
1254 if (fDebug >= 2 || npart<=0)
1255 printf("\n AliEMCALJetFinder::FillFromPartons : npart %i\n", npart);
1256 fNt = 0; // for FindTracksInJetCone
1259 // Go through the partons
1261 for (Int_t part = 8; part < npart; part++) {
1262 TParticle *mPart = gAlice->GetMCApp()->Particle(part);
1263 Int_t mpart = mPart->GetPdgCode();
1264 // Int_t child1 = MPart->GetFirstDaughter();
1265 Float_t pT = mPart->Pt();
1266 // Float_t p = MPart->P();
1267 Float_t phi = mPart->Phi();
1268 Float_t eta = mPart->Eta();
1269 // Float_t theta = MPart->Theta();
1270 statusCode = mPart->GetStatusCode();
1272 // accept partons (21 - gluon, 92 - string)
1273 if (!(TMath::Abs(mpart) <= 6 || mpart == 21 ||mpart == 92)) continue;
1274 if (fDebug > 1 && pT>0.01)
1275 printf("\n part:%5d mpart %5d status %5d eta %8.2f phi %8.2f pT %8.2f ",
1276 part, mpart, statusCode, eta, phi, pT);
1277 // if (fDebug >= 3) MPart->Print();
1278 // accept partons before fragmentation - p.57 in Pythia manual
1279 // if(statusCode != 1) continue;
1281 if (eta > fEtaMax || eta < fEtaMin) continue;
1282 if (phi > fPhiMax || phi < fPhiMin ) continue;
1284 // if (child1 >= 0 && child1 < npart) continue;
1287 fLego->Fill(eta, phi, pT);
1293 void AliEMCALJetFinder::BuildTrackFlagTable() {
1295 // Method to generate a lookup table for TreeK particles
1296 // which are linked to hits in the EMCAL
1298 // --Author: J.L. Klay
1300 // Access hit information
1301 AliEMCAL* pEMCAL = (AliEMCAL*) gAlice->GetModule("EMCAL");
1303 TTree *treeK = gAlice->TreeK(); // Get the number of entries in the kine tree
1304 Int_t nKTrks = (Int_t) treeK->GetEntries(); // (Number of particles created somewhere)
1306 if(fTrackList) delete[] fTrackList; //Make sure we get rid of the old one
1307 fTrackList = new Int_t[nKTrks]; //before generating a new one
1309 for (Int_t i = 0; i < nKTrks; i++) { //Initialize members to 0
1313 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
1314 // TTree *treeH = gAlice->TreeH();
1315 TTree *treeH = gime->TreeH();
1316 Int_t ntracks = (Int_t) treeH->GetEntries();
1322 for (Int_t track=0; track<ntracks;track++) {
1323 gAlice->ResetHits();
1324 nbytes += treeH->GetEvent(track);
1330 for(AliEMCALHit* mHit=(AliEMCALHit*) pEMCAL->FirstHit(-1);
1332 mHit=(AliEMCALHit*) pEMCAL->NextHit())
1334 Int_t iTrk = mHit->Track(); // track number
1335 Int_t idprim = mHit->GetPrimary(); // primary particle
1337 //Determine the origin point of this particle - it made a hit in the EMCAL
1338 TParticle *trkPart = gAlice->GetMCApp()->Particle(iTrk);
1339 TParticlePDG *trkPdg = trkPart->GetPDG();
1340 Int_t trkCode = trkPart->GetPdgCode();
1342 if (trkCode < 10000) { //Big Ions cause problems for
1343 trkChg = trkPdg->Charge(); //this function. Since they aren't
1344 } else { //likely to make it very far, set
1345 trkChg = 0.0; //their charge to 0 for the Flag test
1347 Float_t vxTrk = trkPart->Vx();
1348 Float_t vyTrk = trkPart->Vy();
1349 Float_t vrTrk = TMath::Sqrt(vxTrk*vxTrk+vyTrk*vyTrk);
1350 fTrackList[iTrk] = SetTrackFlag(vrTrk,trkCode,trkChg);
1352 //Loop through the ancestry of the EMCAL entrance particles
1353 Int_t ancestor = trkPart->GetFirstMother(); //Get track's Mother
1354 while (ancestor != -1) {
1355 TParticle *ancPart = gAlice->GetMCApp()->Particle(ancestor); //get particle info on ancestor
1356 TParticlePDG *ancPdg = ancPart->GetPDG();
1357 Int_t ancCode = ancPart->GetPdgCode();
1359 if (ancCode < 10000) {
1360 ancChg = ancPdg->Charge();
1364 Float_t vxAnc = ancPart->Vx();
1365 Float_t vyAnc = ancPart->Vy();
1366 Float_t vrAnc = TMath::Sqrt(vxAnc*vxAnc+vyAnc*vyAnc);
1367 fTrackList[ancestor] = SetTrackFlag(vrAnc,ancCode,ancChg);
1368 ancestor = ancPart->GetFirstMother(); //Get the next ancestor
1371 //Determine the origin point of the primary particle associated with the hit
1372 TParticle *primPart = gAlice->GetMCApp()->Particle(idprim);
1373 TParticlePDG *primPdg = primPart->GetPDG();
1374 Int_t primCode = primPart->GetPdgCode();
1376 if (primCode < 10000) {
1377 primChg = primPdg->Charge();
1381 Float_t vxPrim = primPart->Vx();
1382 Float_t vyPrim = primPart->Vy();
1383 Float_t vrPrim = TMath::Sqrt(vxPrim*vxPrim+vyPrim*vyPrim);
1384 fTrackList[idprim] = SetTrackFlag(vrPrim,primCode,primChg);
1390 Int_t AliEMCALJetFinder
1391 ::SetTrackFlag(Float_t radius, Int_t code, Double_t charge) {
1392 // Set the flag for the track
1397 if (charge == 0) neutral = 1;
1399 if (TMath::Abs(code) <= 6 ||
1401 code == 92) parton = 1;
1403 //It's not a parton, it's charged and it went through the TPC
1404 if (!parton && !neutral && radius <= 84.0) flag = 1;
1411 void AliEMCALJetFinder::SaveBackgroundEvent(Char_t *name)
1413 // Saves the eta-phi lego and the tracklist and name of file with BG events
1417 (*fLegoB) = (*fLegoB) + (*fLego);
1419 printf("\n AliEMCALJetFinder::SaveBackgroundEvent() (fLegoB) %f = %f(fLego) \n",
1420 fLegoB->Integral(), fLego->Integral());
1423 if (fPtB) delete[] fPtB;
1424 if (fEtaB) delete[] fEtaB;
1425 if (fPhiB) delete[] fPhiB;
1426 if (fPdgB) delete[] fPdgB;
1427 if (fTrackListB) delete[] fTrackListB;
1429 fPtB = new Float_t[fNtS];
1430 fEtaB = new Float_t[fNtS];
1431 fPhiB = new Float_t[fNtS];
1432 fPdgB = new Int_t [fNtS];
1433 fTrackListB = new Int_t [fNtS];
1437 for (Int_t i = 0; i < fNt; i++) {
1438 if (!fTrackList[i]) continue;
1439 fPtB [fNtB] = fPtT [i];
1440 fEtaB[fNtB] = fEtaT[i];
1441 fPhiB[fNtB] = fPhiT[i];
1442 fPdgB[fNtB] = fPdgT[i];
1443 fTrackListB[fNtB] = 1;
1447 printf(" fNtB %i => fNtS %i #particles %i \n", fNtB, fNtS, fNt);
1449 if(strlen(name) == 0) {
1450 TSeqCollection *li = gROOT->GetListOfFiles();
1452 for(Int_t i=0; i<li->GetSize(); i++) {
1453 nf = ((TFile*)li->At(i))->GetName();
1454 if(nf.Contains("backgorund")) break;
1460 printf("BG file name is \n %s\n", fBGFileName.Data());
1463 void AliEMCALJetFinder::InitFromBackground()
1467 if (fDebug) printf("\n AliEMCALJetFinder::InitFromBackground() ");
1471 (*fLego) = (*fLego) + (*fLegoB);
1473 printf("\n AliEMCALJetFinder::InitBackgroundEvent() (fLego) %f = %f(fLegoB) \n",
1474 fLego->Integral(), fLegoB->Integral());
1476 printf(" => fLego undefined \n");
1481 void AliEMCALJetFinder::FindTracksInJetCone()
1484 // Build list of tracks inside jet cone
1487 Int_t njet = Njets();
1488 for (Int_t nj = 0; nj < njet; nj++)
1490 Float_t etaj = JetEtaW(nj);
1491 Float_t phij = JetPhiW(nj);
1492 Int_t nT = 0; // number of associated tracks
1494 // Loop over particles in current event
1495 for (Int_t part = 0; part < fNt; part++) {
1496 if (!fTrackList[part]) continue;
1497 Float_t phi = fPhiT[part];
1498 Float_t eta = fEtaT[part];
1499 Float_t dr = TMath::Sqrt((etaj-eta)*(etaj-eta) +
1500 (phij-phi)*(phij-phi));
1501 if (dr < fConeRadius) {
1502 fTrackList[part] = nj+2;
1507 // Same for background event if available
1510 for (Int_t part = 0; part < fNtB; part++) {
1511 Float_t phi = fPhiB[part];
1512 Float_t eta = fEtaB[part];
1513 Float_t dr = TMath::Sqrt((etaj-eta)*(etaj-eta) +
1514 (phij-phi)*(phij-phi));
1515 fTrackListB[part] = 1;
1517 if (dr < fConeRadius) {
1518 fTrackListB[part] = nj+2;
1522 } // Background available ?
1524 Int_t nT0 = nT + nTB;
1525 printf("Total number of tracks %d\n", nT0);
1527 if (nT0 > 1000) nT0 = 1000;
1529 Float_t* ptT = new Float_t[nT0];
1530 Float_t* etaT = new Float_t[nT0];
1531 Float_t* phiT = new Float_t[nT0];
1532 Int_t* pdgT = new Int_t[nT0];
1537 for (Int_t part = 0; part < fNt; part++) {
1538 if (fTrackList[part] == nj+2) {
1540 for (j=iT-1; j>=0; j--) {
1541 if (fPtT[part] > ptT[j]) {
1546 for (j=iT-1; j>=index; j--) {
1547 ptT [j+1] = ptT [j];
1548 etaT[j+1] = etaT[j];
1549 phiT[j+1] = phiT[j];
1550 pdgT[j+1] = pdgT[j];
1552 ptT [index] = fPtT [part];
1553 etaT[index] = fEtaT[part];
1554 phiT[index] = fPhiT[part];
1555 pdgT[index] = fPdgT[part];
1557 } // particle associated
1558 if (iT > nT0) break;
1562 for (Int_t part = 0; part < fNtB; part++) {
1563 if (fTrackListB[part] == nj+2) {
1565 for (j=iT-1; j>=0; j--) {
1566 if (fPtB[part] > ptT[j]) {
1572 for (j=iT-1; j>=index; j--) {
1573 ptT [j+1] = ptT [j];
1574 etaT[j+1] = etaT[j];
1575 phiT[j+1] = phiT[j];
1576 pdgT[j+1] = pdgT[j];
1578 ptT [index] = fPtB [part];
1579 etaT[index] = fEtaB[part];
1580 phiT[index] = fPhiB[part];
1581 pdgT[index] = fPdgB[part];
1583 } // particle associated
1584 if (iT > nT0) break;
1586 } // Background available ?
1588 fJetT[nj]->SetTrackList(nT0, ptT, etaT, phiT, pdgT);
1597 Float_t AliEMCALJetFinder::PropagatePhi(Float_t pt, Float_t charge, Bool_t& curls)
1599 // Propagates phi angle to EMCAL radius
1601 static Float_t b = 0.0, rEMCAL = -1.0;
1603 b = gAlice->Field()->SolenoidField();
1604 // Get EMCAL radius in cm
1605 rEMCAL = AliEMCALGeometry::GetInstance()->GetIPDistance();
1613 Float_t rB = 3335.6 * pt / b; // [cm] (case of |charge|=1)
1615 // check if particle is curling below EMCAL
1616 if (2.*rB < rEMCAL) {
1621 // if not calculate delta phi
1622 Float_t phi = TMath::ACos(1.-rEMCAL*rEMCAL/(2.*rB*rB));
1623 dPhi = TMath::ATan2(1.-TMath::Cos(phi), TMath::Sin(phi));
1624 dPhi = -TMath::Sign(dPhi, charge);
1629 void hf1(Int_t& , Float_t& , Float_t& )
1631 // dummy for hbook calls
1635 void AliEMCALJetFinder::DrawLego(Char_t *opt)
1638 if(fLego) fLego->Draw(opt);
1641 void AliEMCALJetFinder::DrawLegoBackground(Char_t *opt)
1643 // Draw background lego
1644 if(fLegoB) fLegoB->Draw(opt);
1647 void AliEMCALJetFinder::DrawLegoEMCAL(Char_t *opt)
1650 if(fhLegoEMCAL) fhLegoEMCAL->Draw(opt);
1653 void AliEMCALJetFinder::DrawHistsForTuning(Int_t mode)
1656 static TPaveText *varLabel=0;
1658 fC1 = new TCanvas("C1","Hists. for tunning", 0,25,600,800);
1668 fhCellEMCALEt->Draw();
1673 fhTrackPtBcut->SetLineColor(2);
1674 fhTrackPtBcut->Draw("same");
1679 varLabel = new TPaveText(0.05,0.5,0.95,0.95,"NDC");
1680 varLabel->SetTextAlign(12);
1681 varLabel->SetFillColor(19); // see TAttFill
1682 TString tmp(GetTitle());
1683 varLabel->ReadFile(GetFileNameForParameters());
1687 if(mode) { // for saving picture to the file
1688 TString stmp(GetFileNameForParameters());
1689 stmp.ReplaceAll("_Par.txt",".ps");
1690 fC1->Print(stmp.Data());
1694 void AliEMCALJetFinder::PrintParameters(Int_t mode)
1696 // Print all parameters out
1699 if(mode==0) file = stdout; // output to terminal
1701 file = fopen(GetFileNameForParameters(),"w");
1702 if(file==0) file = stdout;
1704 fprintf(file,"==== Filling lego ==== \n");
1705 fprintf(file,"Sampling fraction %6.3f ", fSamplingF);
1706 fprintf(file,"Smearing %6i ", fSmear);
1707 fprintf(file,"Efficiency %6i\n", fEffic);
1708 fprintf(file,"Hadr.Correct. %6i ", fHCorrection);
1709 fprintf(file,"P_{T} Cut of ch.par. %6.2f\n", fPtCut);
1710 fprintf(file,"==== Jet finding ==== \n");
1711 fprintf(file,"Cone radius %6.2f ", fConeRadius);
1712 fprintf(file,"Seed E_{T} %6.1f\n", fEtSeed);
1713 fprintf(file,"Min E_{T} of cell %6.1f ", fMinCellEt);
1714 fprintf(file,"Min E_{T} of jet %6.1f\n", fMinJetEt);
1716 fprintf(file,"==== Bg subtraction ==== \n");
1717 fprintf(file,"BG subtraction %6i ", fMode);
1718 fprintf(file,"Min cone move %6.2f\n", fMinMove);
1719 fprintf(file,"Max cone move %6.2f ", fMaxMove);
1720 fprintf(file,"%% change for BG %6.4f\n", fPrecBg);
1722 fprintf(file,"==== No Bg subtraction ==== \n");
1724 printf("BG file name is %s \n", fBGFileName.Data());
1725 if(file != stdout) fclose(file);
1728 void AliEMCALJetFinder::DrawLegos()
1732 fC1 = new TCanvas("C1","Hists. for tunning", 0,25,600,800);
1736 gStyle->SetOptStat(111111);
1738 Int_t nent1, nent2, nent3, nent4;
1739 Double_t int1, int2, int3, int4;
1740 nent1 = (Int_t)fLego->GetEntries();
1741 int1 = fLego->Integral();
1743 if(int1) fLego->Draw("lego");
1745 nent2 = (Int_t)fhLegoTracks->GetEntries();
1746 int2 = fhLegoTracks->Integral();
1748 if(int2) fhLegoTracks->Draw("lego");
1750 nent3 = (Int_t)fhLegoEMCAL->GetEntries();
1751 int3 = fhLegoEMCAL->Integral();
1753 if(int3) fhLegoEMCAL->Draw("lego");
1755 nent4 = (Int_t)fhLegoHadrCorr->GetEntries();
1756 int4 = fhLegoHadrCorr->Integral();
1758 if(int4) fhLegoHadrCorr->Draw("lego");
1760 // just for checking
1761 printf(" Integrals \n");
1762 printf("lego %10.3f \ntrack %10.3f \nhits %10.3f \nHCorr %10.3f\n-- %10.3f(must be 0)\n",
1763 int1, int2, int3, int4, int1 - (int2 + int3 - int4));
1766 const Char_t* AliEMCALJetFinder::GetFileNameForParameters(Char_t* dir)
1768 // Get paramters from a file
1770 if(strlen(dir)) tmp = dir;
1776 void AliEMCALJetFinder::RearrangeParticlesMemory(Int_t npart)
1778 // See FillFromTracks() - npart must be positive
1779 if (fTrackList) delete[] fTrackList;
1780 if (fPtT) delete[] fPtT;
1781 if (fEtaT) delete[] fEtaT;
1782 if (fPhiT) delete[] fPhiT;
1783 if (fPdgT) delete[] fPdgT;
1786 fTrackList = new Int_t [npart];
1787 fPtT = new Float_t[npart];
1788 fEtaT = new Float_t[npart];
1789 fPhiT = new Float_t[npart];
1790 fPdgT = new Int_t[npart];
1792 printf("AliEMCALJetFinder::RearrangeParticlesMemory : npart = %d\n", npart);
1796 Bool_t AliEMCALJetFinder::IsThisPartonsOrDiQuark(Int_t pdg)
1798 // Return quark info
1799 Int_t absPdg = TMath::Abs(pdg);
1800 if(absPdg<=6) return kTRUE; // quarks
1801 if(pdg == 21) return kTRUE; // gluon
1802 if(pdg == 92) return kTRUE; // string
1804 // see p.51 of Pythia Manual
1805 // Not include diquarks with c and b quark - 4-mar-2002
1806 // ud_0,sd_0,su_0; dd_1,ud_1,uu_1; sd_1,su_1,ss_1
1807 static Int_t diquark[9]={2101,3101,3201, 1103,2103,2203, 3103,3203,3303};
1808 for(Int_t i=0; i<9; i++) if(absPdg == diquark[i]) return kTRUE; // diquarks
1813 void AliEMCALJetFinder::FindChargedJet()
1816 // Find jet from charged particle information only
1831 for (part = 0; part < fNt; part++) {
1832 if (!fTrackList[part]) continue;
1833 if (fPtT[part] > fEtSeed) nseed++;
1835 printf("\nThere are %d seeds (%d)\n", nseed, fNtS);
1836 Int_t* iSeeds = new Int_t[nseed];
1839 for (part = 0; part < fNt; part++) {
1840 if (!fTrackList[part]) continue;
1841 if (fPtT[part] > fEtSeed) iSeeds[nseed++] = part;
1852 // Find seed with highest pt
1854 Float_t ptmax = -1.;
1857 for (seed = 0; seed < nseed; seed++) {
1858 if ((pt = fPtT[iSeeds[seed]]) > ptmax && iSeeds[seed] != -1) {
1863 if (ptmax < 0.) break;
1864 jndex = iSeeds[index];
1866 // Remove from the list
1868 printf("\n Next Seed %d %f", jndex, ptmax);
1870 // Find tracks in cone around seed
1872 Float_t phiSeed = fPhiT[jndex];
1873 Float_t etaSeed = fEtaT[jndex];
1879 for (part = 0; part < fNt; part++) {
1880 if (!fTrackList[part]) continue;
1881 Float_t deta = fEtaT[part] - etaSeed;
1882 Float_t dphi = fPhiT[part] - phiSeed;
1883 Float_t dR = TMath::Sqrt(deta * deta + dphi * dphi);
1884 if (dR < fConeRadius) {
1886 Float_t theta = 2. * TMath::ATan(TMath::Exp(-fEtaT[part]));
1887 Float_t px = fPtT[part] * TMath::Cos(fPhiT[part]);
1888 Float_t py = fPtT[part] * TMath::Sin(fPhiT[part]);
1889 Float_t pz = fPtT[part] / TMath::Tan(theta);
1894 // if seed, remove it
1896 for (seed = 0; seed < nseed; seed++) {
1897 if (part == iSeeds[seed]) iSeeds[seed] = -1;
1903 // Estimate of jet direction
1904 Float_t phiJ = TMath::ATan2(pyJ, pxJ);
1905 Float_t thetaJ = TMath::ATan2(TMath::Sqrt(pxJ * pxJ + pyJ * pyJ), pzJ);
1906 Float_t etaJ = TMath::Log(TMath::Tan(thetaJ / 2.));
1907 // Float_t ptJ = TMath::Sqrt(pxJ * pxJ + pyJ * pyJ);
1910 // Sum up all energy
1912 Int_t iPhi0 = Int_t((phiJ - fPhiMin) / fDphi);
1913 Int_t iEta0 = Int_t((etaJ - fEtaMin) / fDeta);
1914 Int_t dIphi = Int_t(fConeRadius / fDphi);
1915 Int_t dIeta = Int_t(fConeRadius / fDeta);
1918 for (iPhi = iPhi0 -dIphi; iPhi < iPhi0 + dIphi; iPhi++) {
1919 for (iEta = iEta0 - dIeta; iEta < iEta0 + dIeta; iEta++) {
1920 if (iPhi < 0 || iEta < 0) continue;
1921 Float_t dPhi = fPhiMin + iPhi * fDphi;
1922 Float_t dEta = fEtaMin + iEta * fDeta;
1923 if (TMath::Sqrt(dPhi * dPhi + dEta * dEta) < fConeRadius) continue;
1924 sumE += fLego->GetBinContent(iEta, iPhi);
1930 fJetT[njets++] = new AliEMCALJet(sumE, phiJ, etaJ);
1931 FindTracksInJetCone();
1932 printf("\n Jet Energy %f %f %f %f %d\n", eT, sumE, fPtT[6], fPtT[7], njets);
1933 printf("\n Jet Phi %f %f %f \n", phiJ, fPhiT[6], fPhiT[7]);
1934 printf("\n Jet Eta %f %f %f \n", etaJ, fEtaT[6], fEtaT[7]);
1936 EMCALJETS.njet = njets;
1937 if (fWrite) WriteJets();
1940 // 16-jan-2003 - just for convenience
1941 void AliEMCALJetFinder::Browse(TBrowser* b)
1944 if(fHistsList) b->Add((TObject*)fHistsList);
1947 Bool_t AliEMCALJetFinder::IsFolder() const
1949 // Return folder status
1950 if(fHistsList) return kTRUE;
1954 const Char_t* AliEMCALJetFinder::GetNameOfVariant()
1956 // generate the literal string with info about jet finder
1958 sprintf(name, "jF_R%3.2fMinCell%4.1fPtCut%4.1fEtSeed%4.1fMinEt%4.1fBGSubtr%iSF%4.1f",
1959 fConeRadius,fMinCellEt,fPtCut,fEtSeed,fMinJetEt, fMode, fSamplingF);
1961 nt.ReplaceAll(" ","");
1962 if(fBGFileName.Length()) {
1963 Int_t i1 = fBGFileName.Index("kBackground");
1964 Int_t i2 = fBGFileName.Index("/0000") - 1;
1965 if(i1>=0 && i2>=0) {
1966 TString bg(fBGFileName(i1,i2-i1+1));
1970 printf("<I> Name of variant %s \n", nt.Data());