]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALJetMicroDst.cxx
Moved the EMCAL back to 60 < phi < 180
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALJetMicroDst.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2002, 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
17 /* $Id$ */
18
19 //*-- Authors: Aleksei Pavlinov (WSU) 
20
21 //*
22
23 #include "AliEMCALJetMicroDst.h"
24 #include "AliRun.h"
25 #include "AliHeader.h"
26 #include "AliGenEventHeader.h"
27 #include "AliGenHijingEventHeader.h"
28 #include "AliEMCALJetFinder.h"
29 #include <TVector3.h>
30 #include <TROOT.h>
31 #include <TBrowser.h>
32 #include <TString.h>
33 #include <TParticle.h>
34 #include "AliMC.h"
35
36 ClassImp(AliEMCALJetMicroDst)
37
38 TString nameTree("jetMDST"); // 7-feb-2002
39
40 TH1F*  hPtPart, *hNJet, *hPtJet;
41 TH2F*  hEtaPhiPart, *hEtaPhiJet;
42 TH1F*  hNcell, *hCellId, *hCellEt, *hSumEt;
43 TH1F*  hNgrid, *hGridId, *hGridEt, *hSumEtGrForJF;
44
45 extern "C" void sgpdge_(Int_t &i, Int_t &pdggea);
46
47 AliEMCALJetMicroDst::AliEMCALJetMicroDst(char *name, char *tit) : TNamed(name,tit)
48 {
49   fFile  = 0;
50   fTree  = 0;
51   fDebug = 0;
52     
53 //  Don't add histos to the current directory
54 //  TH1::AddDirectory(0);
55   gROOT->cd();
56   hPtPart     = new TH1F("hPtPart","P_{T} for partons", 300, 0., 300.);
57   // 16-jan-2002 - new limit fo phi 
58   hEtaPhiPart = new TH2F("hEtaPhiPart","#eta #phi distr.for partons after HSc", 
59                        28, -0.7, 0.7, 21, 0.0, (2.0/3.0)*TMath::Pi());
60
61   hNJet  = new TH1F("hNJet","number of jets", 11, -0.5, 10.5);
62   hPtJet = new TH1F("hPtJet","P_{T} for jets", 500, 0., 500.);
63   hEtaPhiJet = new TH2F("hEtaPhiJet","#eta #phi distr.for jets (W)", 
64                        28, -0.7, 0.7, 21, 0.0, (2.0/3.0)*TMath::Pi());
65
66   hNcell  = new TH1F("hNcell","#cell with de>0.0 for EMCAL", 1400, 0.0, 14000.);
67   hCellId = new TH1F("hCellId","cell ID with de>0.0 for EMCAL", 1400, 0.0, 14000.);
68   hCellEt = new TH1F("hCellEt","cell Et for EMCAL", 1000, 0.0, 10.);
69   hSumEt  = new TH1F("hSumEt","sum Et for EMCAL", 1000, 0.0, 1000.);
70
71   hNgrid  = new TH1F("hNgrid","#cell with de>0.0 in EMCAL grid for JF", 1400, 0.0, 14000.);
72   hGridId = new TH1F("hGridId","cell ID with de>0.0 in EMCAL grid for JF", 1400, 0.0, 14000.);
73   hGridEt = new TH1F("hGridEt","cell Et in EMCAL grid for JF", 1000, 0.0, 10.);
74   hSumEtGrForJF  = new TH1F("hSumEtGrForJF","sum Et in EMCAL grid for JF", 1000, 0.0, 1000.);
75
76   fListHist = MoveHistsToList("Hist For AliEMCALJetMicroDst", kFALSE); 
77 }
78
79 AliEMCALJetMicroDst::~AliEMCALJetMicroDst()
80 {
81   if(fFile) fFile->Close();
82 }  
83
84 Bool_t AliEMCALJetMicroDst::Create(TFile *file)
85 {
86   if(!file) {
87     Error("Create", "define TFile for output\n");
88     return kFALSE;
89   }
90   fFile     = file;
91   fFileName = fFile->GetName();
92   fFile->cd();
93   fTree = new TTree(nameTree.Data(),"Temporary micro DST for jets analysis");
94   // for jet calibration - 4-mar-2003
95   fTree->Branch("decone", &decone, "decone/F");
96   fTree->Branch("ptcone", &ptcone, "ptcone/F");
97   // partons 
98   fTree->Branch("npart", &npart, "npart/I");
99   fTree->Branch("xpt",  xpt,  "xpt[npart]/F");
100   fTree->Branch("xeta", xeta, "xeta[npart]/F");
101   fTree->Branch("xphi", xphi, "xphi[npart]/F");
102   // jets 
103   fTree->Branch("njet", &njet, "njet/I");
104   fTree->Branch("jet",  jet,  "jet[njet]/F");
105   fTree->Branch("jetaw", jetaw, "jetaw[njet]/F");
106   fTree->Branch("jphiw", jphiw, "jphiw[njet]/F");
107   fTree->Branch("jetal", jetal, "jetal[njet]/F");
108   fTree->Branch("jphil", jphil, "jphil[njet]/F");
109
110   // Et in EMCAL itself 
111   fTree->Branch("ncell", &ncell, "ncell/I");
112   fTree->Branch("idcell", idcell, "idcell[ncell]/I");
113   fTree->Branch("etcell", etcell, "etcell[ncell]/F");
114
115   // Et in EMCAL grid for JF 
116   fTree->Branch("ngrid", &ngrid,  "ngrid/I");
117   fTree->Branch("idgrid", idgrid, "idgrid[ngrid]/I");
118   fTree->Branch("etgrid", etgrid, "etgrid[ngrid]/F");
119
120   // charge particle which hit to EMCAL
121   fTree->Branch("nchp", &nchp, "nchp/I");
122   fTree->Branch("pid", pid, "pid[nchp]/I");
123   fTree->Branch("ppt", ppt, "ppt[nchp]/F");
124   fTree->Branch("peta", peta, "peta[nchp]/F");
125   fTree->Branch("pphi", pphi, "pphi[nchp]/F");
126
127   return kTRUE;
128 }
129
130 Bool_t AliEMCALJetMicroDst::Create(const char *fname)
131 {
132   TFile *file = new TFile(fname, "RECREATE");
133   if(file) {
134     //    fNameFile = fname;
135     return Create(file);
136   } else   return kFALSE;
137 }
138
139 Bool_t AliEMCALJetMicroDst::Open(const char *fname)
140 {
141   if(fFile && fFile->IsOpen()) fFile->Close();
142   if(strlen(fname)) fName = fname;
143   TFile *file = new TFile(fName.Data(), "READ");
144   if(file) {
145     Bool_t ini =  Initialize(file);
146     Info("Open", "open file %s : initialize TTree %i",fName.Data(), Int_t(ini)); 
147     return ini;
148   } else {
149     Error("Open", "can not open file %s",fName.Data()); 
150     return kFALSE;
151   }
152 }
153
154 const Char_t* AliEMCALJetMicroDst::DefineName(const Int_t mode)
155 {
156   static TString dir, name;
157   //  dir = "jetDST/"; // 24-jan-2003
158   dir = "/auto/alice/pavlinov/jet/microDST/"; // 24-jan-2003
159   switch (mode) {
160   case 1: // for characteristic of BG
161     name = dir + "mDst1_1.root"; // Bg 2000 - first version 
162     SetTitle("Bg2000");
163     break;
164   case 2: // for characteristic of BG
165     name = dir + "mDst2_1.root"; // Bg 4000 - first version 
166     SetTitle("Bg4000");
167     break;
168   case 3: // for characteristic of BG
169     name = dir + "mDst3_1.root"; // Bg 8000 - first version 
170     SetTitle("Bg8000");
171     break;
172   case 4: // Central Hijing - 18-mar-2003
173     name  = dir  + "march/";
174     name += "jF_R0.50MinCell1.0PtCut0.0EtSeed4.0MinEt40.0BGSubtr0SF11.6Smear0Eff0HijingCentral.root";
175     SetTitle("HijingCentral");
176     break;
177   case 5: // Para Hijing (Dn/Dy=8000) - 21-mar-2003
178     name  = dir  + "march/";
179     name += "jF_R0.50MinCell1.0PtCut0.0EtSeed4.0MinEt40.0BGSubtr0SF11.6Smear0Eff0ParaHijing8000.root";
180     SetTitle("HIJINGparaDnDy8000");
181     break;
182   case 6: // Para Hijing (Dn/Dy=4000) - 21-mar-2003
183     name  = dir  + "march/";
184     name += "jF_R0.50MinCell1.0PtCut0.0EtSeed4.0MinEt40.0BGSubtr0SF11.6Smear0Eff0ParaHijing4000.root";
185     SetTitle("HIJINGparaDnDy4000");
186     break;
187   case 7: // Para Hijing (Dn/Dy=2000) - 21-mar-2003
188     name  = dir  + "march/";
189     name += "jF_R0.50MinCell1.0PtCut0.0EtSeed4.0MinEt40.0BGSubtr0SF11.6Smear0Eff0ParaHijing2000.root";
190     SetTitle("HIJINGparaDnDy2000");
191     break;
192   case 11: // pure PYTHIA with default value of parameters
193     name = dir + "jF_R0.50MinCell0.0PtCut0.0EtSeed8.0MinEt40.0BGSubtr0SF11.6.root";
194     break;
195   case 12: // 0 + background 
196     name = dir + "jF_R0.50MinCell0.0PtCut0.0EtSeed8.0MinEt40.0BGSubtr0SF11.6kBackground2000.root";
197     break;
198     // calibration case 
199   case 101:
200     name  = dir  + "march/";
201     name += "Pythia100_1.root";
202     SetTitle("Pythia100_1");
203     break;
204   case 102: // 2-apr-2003
205     name  = dir  + "march/";
206     name += "Pythia50_1.root";
207     SetTitle("Pythia50_1");
208     //    name = "microDst3th.root"; // 101 + (smearing and eff) - 14-mar-2003
209     break;
210   case 103:// 4-apr-2003
211     name  = dir  + "march/";
212     name += "Pythia200_1.root";
213     SetTitle("Pythia200_1");
214     //    name = "microDst4th.root"; // 102 + MinCell= 1.0
215     break;
216   default:
217     Fatal("DefineName", "NO  D E F A U L T : mode %i\n", mode);
218   }
219   Info("DefineName", "mode %5i file : %s : Title %s\n", mode, name.Data(), GetTitle());
220   return name.Data();
221 }
222
223 Bool_t AliEMCALJetMicroDst::Initialize(TFile *file)
224 {
225   if(file) fFile = file; 
226   fFile->cd();
227   fTree = (TTree*)fFile->Get(nameTree.Data());
228   if(!fTree) return kFALSE;
229   // for jet calibration - 4-mar-2003
230   fTree->SetBranchAddress("decone",&decone);
231   fTree->SetBranchAddress("ptcone",&ptcone);
232   // partons
233   fTree->SetBranchAddress("npart",&npart);
234   fTree->SetBranchAddress("xpt",   xpt);
235   fTree->SetBranchAddress("xeta",  xeta);
236   fTree->SetBranchAddress("xphi",  xphi);
237   // jets 
238   fTree->SetBranchAddress("njet", &njet);
239   fTree->SetBranchAddress("jet",   jet);
240   fTree->SetBranchAddress("jetaw", jetaw);
241   fTree->SetBranchAddress("jphiw", jphiw);
242   fTree->SetBranchAddress("jetal", jetal);
243   fTree->SetBranchAddress("jphil", jphil);
244   // eT in EMCAL
245   fTree->SetBranchAddress("ncell", &ncell);
246   fTree->SetBranchAddress("idcell", idcell);
247   fTree->SetBranchAddress("etcell", etcell);
248   // eT in EMCAL grid for JF
249   fTree->SetBranchAddress("ngrid", &ngrid);
250   fTree->SetBranchAddress("idgrid", idgrid);
251   fTree->SetBranchAddress("etgrid", etgrid);
252   // 28-jan-2003
253   fTree->SetBranchAddress("nchp", &nchp);
254   fTree->SetBranchAddress("pid", pid);
255   fTree->SetBranchAddress("ppt", ppt);
256   fTree->SetBranchAddress("peta", peta);
257   fTree->SetBranchAddress("pphi", pphi);
258
259   return kTRUE;
260 }
261
262 void AliEMCALJetMicroDst::Print(Option_t* option) const
263 {
264   if(option);
265   if(fFile) {
266     fFile->Print();
267     if(fTree) fTree->Print();
268     else Info("Print", "TRee is zero\n");
269   } else {
270      Info("Print", "File with TRee is closed \n Name of file %s(from fFileName", fFileName.Data());
271   }
272
273   Info("Print", "******* Current(last) event *****");
274   printf("#partons %2i \n", npart);
275   for(Int_t i=0; i<npart; i++){
276     printf("     %1i Pt %7.1f eta  %7.4f phi  %7.4f \n",
277     i, xpt[i],xeta[i],xphi[i]);  
278   }
279   printf("#jets    %2i \n", njet);
280   for(Int_t i=0; i<njet; i++){
281     printf("     %1i Et %7.1f etaw %7.4f phiw %7.4f \n",
282     i,jet[i],jetaw[i],jphiw[i]);  
283   }
284   printf(" Title %s \n", GetTitle());
285 }
286
287 void AliEMCALJetMicroDst::Fill(AliRun *run, AliEMCALJetFinder* jetFinder, Int_t modeFilling)
288 {//  modeFilling >=1 - fill info for EMCAL grid
289   if(!run) run = gAlice;
290   AliGenEventHeader* evHeader = run->GetHeader()->GenEventHeader();
291   TString tmp(evHeader->ClassName());
292   if(tmp.Contains("Hijing")) {
293      AliGenHijingEventHeader *hijEvHeader = (AliGenHijingEventHeader*)evHeader;
294      FillPartons(hijEvHeader);
295   } else if(tmp.Contains("Pythia")) {
296      FillPartons();     
297   } else {
298     Error("Fill", "Wrong type of generator -> %s \n Info about partons will be absent",tmp.Data()); 
299   }
300
301   FillJets(jetFinder);
302
303   if(modeFilling >= 1) {
304      FillEtForEMCAL(jetFinder);
305      FillEtForGrid(jetFinder);
306      FillChargeParticles(jetFinder);
307   } else {
308      ncell = 0; // 27-jan-2003
309      ngrid = 0; // 27-jan-2003
310      nchp  = 0;
311      // negative - no signal 
312      decone = -1.;
313      ptcone = -1.;
314   }
315
316   FillJetsControl();  //24-jan-2003
317
318   fTree->Fill();
319 }
320
321 void AliEMCALJetMicroDst::FillPartons(AliGenHijingEventHeader *header)
322 {
323   TLorentzVector parton[4];
324   header->GetJets(parton[0], parton[1], parton[2], parton[3]);
325
326   npart = 4; // 
327   for(Int_t i=0; i<4; i++){
328      xpt[i]  = parton[i].Pt();
329      xeta[i] = parton[i].Eta();
330      xphi[i] = parton[i].Phi();
331   }
332 }
333
334 void AliEMCALJetMicroDst::FillPartons()
335 {// for case of Pythia -> get info from full event record
336
337   npart = 2;
338   TParticle *MPart;
339   Int_t ind;
340   for(Int_t i=6; i<8; i++){
341      MPart = gAlice->GetMCApp()->Particle(i);
342      ind   = i-6;
343      xpt[ind]  = MPart->Pt();
344      xeta[ind] = MPart->Eta();
345      xphi[ind] = MPart->Phi();
346   }
347 }
348
349 void AliEMCALJetMicroDst::FillJets(AliEMCALJetFinder* jetFinder)
350 {
351   njet = 0;
352   if(fDebug>1) Info("FillJets", "Debug"); 
353   if(!jetFinder) {
354     if(fDebug>1) Info("FillJets", "jetFinder is zero"); 
355     return;
356   }
357   njet = jetFinder->Njets();
358   if(njet>10) {
359     if(fDebug>1) Warning("FillJets", "wrong value of jetFinder->Njets() %i ", njet); 
360     njet = 10;
361   }
362   //  hNJet->Fill(njet);
363   if(fDebug>1) Info("FillJets", "njet %i", njet); 
364   if(njet){
365     for(Int_t i=0; i<njet; i++){
366       jet[i]   = jetFinder->JetEnergy(i);
367       jetaw[i] = jetFinder->JetEtaW(i);
368       jphiw[i] = jetFinder->JetPhiW(i);
369       jetal[i] = jetFinder->JetEtaL(i);
370       jphil[i] = jetFinder->JetPhiL(i);
371     }
372   }
373 }
374
375 void AliEMCALJetMicroDst::FillEtForEMCAL(AliEMCALJetFinder* jetFinder)
376 {
377    ncell = 0;
378    TH2F *hid = jetFinder->GetLegoEMCAL();
379    if(!hid) return;
380
381    Double_t de = 0.;
382    Int_t neta = hid->GetNbinsX(), nphi = hid->GetNbinsY();
383    for(Int_t ieta=1; ieta<=neta; ieta++) {
384       for(Int_t iphi=1; iphi<=nphi; iphi++) {
385          de = hid->GetBinContent(ieta,iphi);
386          if(de > 0.0) {
387            etcell[ncell] = Float_t(de);
388            idcell[ncell] = nphi*(ieta-1) + iphi;
389            ncell++;
390            if(ncell >= 13824) break; 
391            // Info("FillEtForEMCAL", " ncell %i6 id %i6 de %f \n", ncell, idcell[ncell], etcell[ncell]); 
392          }
393       }
394    }
395    if(njet == 1) {
396      // jet energy calculate around LP direction !!! - 10-mar-2003 
397       decone = jetFinder->EMCALConeEnergy(jetal[0],jphil[0]); 
398       ptcone = jetFinder->TrackConeEnergy(jetal[0],jphil[0]); // get from lego plot fo ch.part
399       Info("FillEtForEMCAL", " njet %i Emcal in cone %f pt ch.part in cone %f\n", njet, decone, ptcone); 
400       Info("FillEtForEMCAL", " jet - decone - ptcone : %9.2f\n", jet[0]-decone-ptcone);
401    } else {
402       decone = -1.;
403       ptcone = -1.;
404    }
405
406    Info("FillEtForEMCAL", "neta %3i nphi %3i # array size %i Sum.Et %f\n", 
407    neta,nphi, ncell, hid->Integral());
408 }
409
410 void AliEMCALJetMicroDst::FillEtForGrid(AliEMCALJetFinder* jetFinder)
411 {
412    TH2F *hid = jetFinder->GetLego();
413    if(!hid) return;
414
415    FillArrays(hid, ngrid, idgrid, etgrid);
416 }
417
418 void AliEMCALJetMicroDst::FillArrays(TH2* hid, Int_t &n, Int_t *id, Float_t *et)
419 {
420    n = 0;
421    Double_t de = 0.;
422    Int_t neta = hid->GetNbinsX(), nphi = hid->GetNbinsY();
423    for(Int_t ieta=1; ieta<=neta; ieta++) {
424       for(Int_t iphi=1; iphi<=nphi; iphi++) {
425          de = hid->GetBinContent(ieta,iphi);
426          if(de > 0.0) {
427            et[n] = Float_t(de);
428            id[n] = nphi*(ieta-1) + iphi;
429            n++;
430            if(n >= 13824) break; 
431          }
432       }
433    }
434    Info("FillArrays", "neta %3i nphi %3i # array size %i Sum.Et %f\n", 
435    neta, nphi, n, hid->Integral());
436 }
437
438 void AliEMCALJetMicroDst::FillChargeParticles(AliEMCALJetFinder* jetFinder)
439 {// 28-jan-2003 for fullness ; 18-mar - sometimes 
440   nchp = 0;
441   Int_t gid=0;
442   for(Int_t i=0; i<jetFinder->fNt; i++) {
443     //     fPdgT[i];
444      if(jetFinder->fTrackList[i] >= 1) {
445         sgpdge_(jetFinder->fPdgT[i], gid);
446         pid[nchp]  = gid; 
447         ppt[nchp]  = jetFinder->fPtT[i];
448         peta[nchp] = jetFinder->fEtaT[i];
449         pphi[nchp] = jetFinder->fPhiT[i];
450         nchp++;
451      }
452   }
453   Info("FillChargedParticles", "fNtS %i : nchp %i -> %i\n", jetFinder->fNtS, nchp, jetFinder->fNtS - nchp);
454 }
455
456 void AliEMCALJetMicroDst::FillJetsControl()
457 {  // see FillJets(AliEMCALJetFinder* jetFinder) and FillPartons
458   hNJet->Fill(njet);
459   for(Int_t i=0; i<njet; i++){
460      hPtJet->Fill(jet[i]);
461      hEtaPhiJet->Fill(jetaw[i],jphiw[i]); 
462   }
463
464   for(Int_t i=0; i < npart; i++){
465      hEtaPhiPart->Fill(xeta[i], xphi[i]);
466      hPtPart->Fill(xpt[i]);
467   }
468
469   Double_t sum = 0.0;
470   hNcell->Fill(ncell);
471   for(Int_t i=0; i < ncell; i++){
472     hCellId->Fill(idcell[i]);
473     hCellEt->Fill(etcell[i]);
474     sum += Double_t(etcell[i]);
475   }
476   hSumEt->Fill(sum);
477
478   sum = 0.0;
479   hNgrid->Fill(ngrid);
480   for(Int_t i=0; i < ngrid; i++){
481     hGridId->Fill(idgrid[i]);
482     hGridEt->Fill(etgrid[i]);
483     sum += Double_t(etgrid[i]);
484   }
485   hSumEtGrForJF->Fill(sum);
486 }
487
488 Int_t AliEMCALJetMicroDst::GetEntry(Int_t entry)
489 { // Read contents of entry.
490    if (!fTree) {
491       Error("GetEntry", "define TTree");
492       return -1;
493    }
494    return fTree->GetEntry(entry);
495 }
496
497 Bool_t AliEMCALJetMicroDst::GetParton(Int_t i, Float_t& pt, Float_t& eta, Float_t& phi)
498 {
499   if(i>=0 && i<npart) {
500     pt  = xpt[i];
501     eta = xeta[i];
502     phi = xphi[i];
503     return kTRUE;
504   } else return kFALSE; 
505 }
506
507 Bool_t AliEMCALJetMicroDst::GetParton(Int_t i, TVector3& vec)
508 {
509   static Float_t pt, eta, phi;
510
511   if(!GetParton(i, pt, eta, phi)) return kFALSE;
512
513   FillVector(pt, eta, phi, vec);
514   return kTRUE;
515 }
516
517 Bool_t AliEMCALJetMicroDst::GetJet(Int_t i,Int_t mode,Float_t& pt, Float_t& eta, Float_t& phi)
518 {// mode=1(W) mode=any(L)
519   if(i>=0 && i<njet) {
520     pt  = jet[i];
521     if(mode==1) {
522       eta = jetaw[i];
523       phi = jphiw[i];
524     } else {
525       eta = jetal[i];
526       phi = jphil[i];
527     }
528     return kTRUE;
529   } else return kFALSE; 
530 }
531
532 Bool_t AliEMCALJetMicroDst::GetJet(Int_t i, Int_t mode, TVector3& vec)
533 {
534   static Float_t pt, eta, phi;
535
536   if(!GetJet(i, mode, pt, eta, phi)) return kFALSE;
537   FillVector(pt, eta, phi, vec);
538   return kTRUE;
539 }
540
541 void AliEMCALJetMicroDst::Test()
542 {
543   if(!fFile || !fTree ) {
544     Info("Test", "define file with proper TTree !");
545     return;
546   }
547   Int_t nbytes=0, nb=0, nentries=Int_t(fTree->GetEntries());
548   for(Int_t i=0; i<nentries; i++){
549     nb = fTree->GetEntry(i);  
550     nbytes += nb;
551     for(Int_t j=0; j<npart; j++){
552       hEtaPhiPart->Fill(xeta[j], xphi[j]);
553       hPtPart->Fill(xpt[j]);
554     }
555
556     hNJet->Fill(njet);
557     if(njet){
558       for(Int_t j=0; j<njet; j++) {
559         hPtJet->Fill(jet[j]);
560         hEtaPhiJet->Fill(jetaw[j],jphiw[j]); 
561       }
562     }
563   }
564   Info("Test", "Entries %5i Bytes %10i\n", nentries, nbytes);
565 }
566
567 void AliEMCALJetMicroDst::FillVector(Float_t pt, Float_t eta, Float_t phi, TVector3& vec)
568 { // service function 
569   static Float_t px, py, pz;
570
571   px = pt*TMath::Cos(phi);
572   py = pt*TMath::Sin(phi);
573   pz = pt*TMath::SinH(eta);  // sinh(eta) = cot(theta)
574
575   vec.SetXYZ(px, py, pz);
576 }
577
578 void AliEMCALJetMicroDst::GetEtaPhi(Int_t id, Double_t &eta, Double_t &phi)
579 { // see AliEMCALGeometry 
580   static Int_t ieta, iphi, nphi=144, neta=96;
581   static Double_t phiMax=(2.0/3.0)*TMath::Pi(), phiMin=0.0;
582   static Double_t phiStep=(phiMax-phiMin)/nphi, phiBeg = phiMin + phiStep/2.; 
583   static Double_t etaMax=0.7, etaMin=-etaMax;
584   static Double_t etaStep=(etaMax-etaMin)/neta, etaBeg = etaMin + etaStep/2.;
585
586   ieta = (id-1)/nphi + 1; // id = nphi*(ieta-1) + iphi
587   iphi = id - nphi*(ieta-1);
588   if(ieta<=0 || ieta>neta) {
589     Fatal("GetEtaPhi", "wrong id %i(ieta %i,iphi %i) : nphi %i neta %i", id,iphi,ieta, nphi,neta);
590
591   }
592
593   eta  = etaBeg + etaStep*(ieta-1);
594   phi  = phiBeg + phiStep*(iphi-1);
595 }
596
597 TVector3& AliEMCALJetMicroDst::GetCellVector(Int_t i)
598 {
599   static Double_t eta,phi;
600   static TVector3 vec;
601   vec.SetXYZ(0.,0.,0.);
602   if(i>=0 && i<ncell) {
603      GetEtaPhi(idcell[i], eta, phi);
604      vec.SetPtEtaPhi(Double_t(etcell[i]),eta,phi);
605   }
606   return vec;
607 }
608
609 TVector3& AliEMCALJetMicroDst::GetGridVector(Int_t i)
610 {
611   static Double_t eta,phi;
612   static TVector3 vec;
613   vec.SetXYZ(0.,0.,0.);
614   if(i>=0 && i<ngrid) {
615      GetEtaPhi(idgrid[i], eta, phi);
616      vec.SetPtEtaPhi(Double_t(etgrid[i]),eta,phi);
617   }
618   return vec;
619 }
620
621 Double_t AliEMCALJetMicroDst::GetSumInCone(TVector3 &jet,Int_t nc, Float_t *et,Float_t *eta,Float_t *phi, Double_t cellEtCut, Double_t rJet)
622 {
623   static Double_t sum=0.;
624   static TVector3 cell(0., 0., 0.);
625   if(nc<=0 || et==0 || eta==0 || phi==0) {
626     Error("GetSumInCone", "nc %d %f %f %f ", nc, et, eta, phi);
627     return -1.;
628   }
629
630   sum=0.;
631   //  jet.SetPtEtaPhi(jet[0],jetaw[0],jphiw[0]); // must be one jet !!
632   Info("GetSumInCone", "jet.Mag() %f : njet %i\n", jet.Mag(), njet);
633   for(Int_t i=0; i<nc; i++){
634     if(et[i] < cellEtCut)    continue;
635     cell.SetPtEtaPhi(et[i], eta[i], phi[i]);
636     if(jet.DeltaR(cell) > rJet) continue;
637     sum += et[i];
638   }
639   Info("GetSumCone", "Sum %f \n", sum);
640   return sum;
641 }
642
643 Double_t AliEMCALJetMicroDst::GetEmcalEtInCone(TVector3 &jet, Double_t cellEtCut, Double_t rJet)
644 {
645   Int_t nc = ncell;
646   if(nc<=0) return 0.;
647
648   Float_t *et=etcell, *eta=new Float_t[nc], *phi=new Float_t[nc];
649   Double_t etaCell=0., phiCell=0., ET=0;
650
651   for(Int_t i=0; i<nc; i++) {
652     GetEtaPhi(idcell[i], etaCell, phiCell);
653     eta[i] = etaCell;
654     phi[i] = phiCell;
655   }
656
657   ET = GetSumInCone(jet, nc, et,eta,phi, cellEtCut,rJet);
658   delete [] eta;
659   delete [] phi;
660
661   return ET;
662 }
663
664 Double_t AliEMCALJetMicroDst::GetTpcPtInCone(TVector3 &jet,Double_t cellEtCut, Double_t rJet)
665 {
666   if(nchp<=0) return 0.;
667   return GetSumInCone(jet, nchp, ppt,peta,pphi, cellEtCut,rJet);
668 }
669
670 Double_t AliEMCALJetMicroDst::GetSum(Int_t n, Float_t *ar, Double_t cut)
671 { // 25-apr-2003
672   Double_t sum=0.0;
673   if(n<=0 || ar==0) return sum;
674   for(Int_t i=0; i<n; i++) {if(ar[i]>=cut) sum += Double_t(ar[i]);}
675   return sum;
676 }
677
678 void AliEMCALJetMicroDst::Close()
679 {
680   fFile->Write(); 
681   fTree->Print(); 
682   fFile->Close();
683   fFile = 0;
684   fTree = 0;
685 }
686
687 void AliEMCALJetMicroDst::Browse(TBrowser* b)
688 {
689    if(fTree)      b->Add((TObject*)fTree);
690    if(fListHist)  b->Add((TObject*)fListHist);
691    //   TObject::Browse(b);
692 }
693
694 Bool_t  AliEMCALJetMicroDst::IsPythiaDst()
695 {
696   TString st(GetTitle());
697   if(st.Contains("py",TString::kIgnoreCase)||!st.Contains("hijing", TString::kIgnoreCase)) return kTRUE;
698   else return kFALSE;
699 }
700
701 Bool_t  AliEMCALJetMicroDst::IsFolder() const
702 {
703   if(fTree || fListHist) return kTRUE;
704   else                   return kFALSE;
705 }
706
707 TList* AliEMCALJetMicroDst::MoveHistsToList(char* name, Bool_t putToBrowser)
708 {
709   gROOT->cd();
710   TIter nextHist(gDirectory->GetList());
711   TList *list = new TList;
712   list->SetName(name);
713   TObject *objHist;
714   while((objHist=nextHist())){
715     if (!objHist->InheritsFrom("TH1")) continue;
716     ((TH1*)objHist)->SetDirectory(0); // Remove from gROOT
717     list->Add(objHist);
718   }
719   if(putToBrowser) gROOT->GetListOfBrowsables()->Add((TObject*)list);
720   return list;
721 }