]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetBkg.cxx
Add data members for jet background, enable filling (Leticia Mendez)
[u/mrichter/AliRoot.git] / JETAN / AliJetBkg.cxx
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 #include <Riostream.h> 
17 #include <TChain.h>
18 #include <TFile.h>
19 #include <TList.h>
20 #include <TArrayF.h>
21 #include <TClonesArray.h>
22 #include <TF1.h>
23 #include <TString.h>
24
25 #include "AliJetHeader.h"
26 #include "AliJetReader.h"
27 #include "AliJetReaderHeader.h"
28 #include "AliFastJetFinder.h"
29 #include "AliFastJetHeaderV1.h"
30 #include "AliJetReaderHeader.h"
31 #include "AliJetReader.h"
32 #include "AliJetUnitArray.h"
33 #include "AliFastJetInput.h"
34 #include "AliESDEvent.h"
35
36
37 #include "fastjet/PseudoJet.hh"
38 #include "fastjet/ClusterSequenceArea.hh"
39 #include "fastjet/AreaDefinition.hh"
40 #include "fastjet/JetDefinition.hh"
41 // get info on how fastjet was configured
42 #include "fastjet/config.h"
43
44 #ifdef ENABLE_PLUGIN_SISCONE
45 #include "fastjet/SISConePlugin.hh"
46 #endif
47
48 #include<sstream>  // needed for internal io
49 #include<vector> 
50 #include <cmath> 
51
52 using namespace std;
53 #include "AliJetBkg.h"
54
55 ClassImp(AliJetBkg)
56
57 ////////////////////////////////////////////////////////////////////////
58
59 AliJetBkg::AliJetBkg():
60   TObject(),
61   fReader(0),
62   fHeader(0),
63   fInputFJ(0)
64 {
65   // Default constructor
66
67 }
68 //______________________________________________________________________
69 AliJetBkg::AliJetBkg(const AliJetBkg& input):
70   TObject(input),
71     fReader(input.fReader),
72     fHeader(input.fHeader),
73     fInputFJ(input.fInputFJ)
74 {
75   // copy constructor
76
77 }
78 //______________________________________________________________________
79 AliJetBkg& AliJetBkg::operator=(const AliJetBkg& source){
80     // Assignment operator. 
81     this->~AliJetBkg();
82     new(this) AliJetBkg(source);
83     return *this;
84
85 }
86
87
88 //___________________________________________________________________
89   void AliJetBkg::BkgFastJetb(Double_t& rho,Double_t& sigma, 
90 Double_t& meanarea){
91   
92   cout<<"===============  AliJetBkg::BkgFastJetb()  =========== "<<endl;
93   vector<fastjet::PseudoJet> inputParticles=fInputFJ->GetInputParticles();
94   
95   //cout<<"printing inputParticles for BKG "<<inputParticles.size()<<endl;
96   
97
98    
99   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
100   double rParamBkg = header->GetRparamBkg(); //Radius for background calculation
101
102   Double_t medianb,sigmab,meanareab;
103   CalcRhob(medianb,sigmab,meanareab,inputParticles,rParamBkg,"All");
104   rho=medianb;
105   sigma=sigmab;
106   meanarea=meanareab;
107  
108 }
109 //_________________________________________________________________
110
111   void AliJetBkg::BkgFastJetWoHardest(Double_t& rho,Double_t& sigma, 
112 Double_t& meanarea){
113   
114   cout<<"===============  AliJetBkg::BkgWoHardest()  =========== "<<endl;
115   vector<fastjet::PseudoJet> inputParticles=fInputFJ->GetInputParticles();
116   
117   //cout<<"printing inputParticles for BKG "<<inputParticles.size()<<endl;
118   
119
120    
121   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
122   double rParamBkg = header->GetRparamBkg(); //Radius for background calculation  
123   Double_t medianb,sigmab,meanareab;
124   CalcRhoWoHardest(medianb,sigmab,meanareab,inputParticles,rParamBkg,"All");
125   rho=medianb;
126   sigma=sigmab;
127   meanarea=meanareab;
128
129   }
130
131 //____________________________________________________________________
132
133  void AliJetBkg::CalcRhob(Double_t& median,Double_t& 
134 sigma,Double_t& 
135 meanarea,vector<fastjet::PseudoJet> inputParticles,Double_t 
136 rParamBkg,TString method){
137   //calculate rho using the fastjet method
138
139   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
140   Bool_t debug  = header->GetDebug();     // debug option
141
142   fastjet::Strategy strategy = header->GetStrategy();
143   fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
144   fastjet::JetAlgorithm algorithm = header->GetAlgorithm(); 
145   fastjet::JetDefinition jetDef(algorithm, rParamBkg, recombScheme, strategy);
146
147   // create an object that specifies how we to define the area
148   fastjet::AreaDefinition areaDef;
149   double ghostEtamax = header->GetGhostEtaMax(); 
150   double ghostArea   = header->GetGhostArea(); 
151   int    activeAreaRepeats = header->GetActiveAreaRepeats(); 
152
153   // now create the object that holds info about ghosts
154
155   
156
157   fastjet::GhostedAreaSpec ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea);
158   // and from that get an area definition
159   fastjet::AreaType areaType = header->GetAreaType();
160   areaDef = fastjet::AreaDefinition(areaType,ghost_spec);
161   
162   //fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef);
163   fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef,areaDef);
164   TString comment = "Running FastJet algorithm for BKG calculation with the following setup. ";
165   comment+= "Jet definition: ";
166   comment+= TString(jetDef.description());
167   //  comment+= ". Area definition: ";
168   //comment+= TString(areaDef.description());
169   comment+= ". Strategy adopted by FastJet: ";
170   comment+= TString(clust_seq.strategy_string());
171   header->SetComment(comment);
172   if(debug){
173     cout << "--------------------------------------------------------" << endl;
174     cout << comment << endl;
175     cout << "--------------------------------------------------------" << endl;
176   }
177
178   vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(0.);
179   vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets); 
180
181  
182   //cout<<"# of BKG jets = "<<jets.size()<<endl;
183       //  for (size_t j = 0; j < jets.size(); j++) { // loop for jets   
184       //printf("BKG Jet found %5d %9.5f %8.5f %10.3f %4.4f 
185       //\n",(Int_t)j,jets[j].rap(),jets[j].phi(),jets[j].perp(),
186       //clust_seq.area(jets[j]));
187       // }
188   
189   
190   double phiMin = 0, phiMax = 0, rapMin = 0, rapMax = 0;
191
192  
193     phiMin = 0;
194     phiMax = 2*TMath::Pi();
195  
196   rapMax = ghostEtamax - rParamBkg;
197   rapMin = - ghostEtamax + rParamBkg;
198
199
200   fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
201
202
203   double medianb, sigmab, meanareab;
204   clust_seq.get_median_rho_and_sigma(inclusiveJets, range, false, medianb, sigmab, meanareab, false);
205   median=medianb;
206   sigma=sigmab;
207   meanarea=meanareab; 
208
209  
210  }
211
212 //____________________________________________________________________
213
214
215  void AliJetBkg::CalcRhoWoHardest(Double_t& median,Double_t& 
216 sigma,Double_t& 
217 meanarea,vector<fastjet::PseudoJet> inputParticles,Double_t 
218 rParamBkg,TString method){
219   //calculate rho using the fastjet method
220
221   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
222   Bool_t debug  = header->GetDebug();     // debug option
223
224   fastjet::Strategy strategy = header->GetStrategy();
225   fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
226   fastjet::JetAlgorithm algorithm = header->GetAlgorithm(); 
227   fastjet::JetDefinition jetDef(algorithm, rParamBkg, recombScheme, strategy);
228
229   // create an object that specifies how we to define the area
230   fastjet::AreaDefinition areaDef;
231   double ghostEtamax = header->GetGhostEtaMax(); 
232   double ghostArea   = header->GetGhostArea(); 
233   int    activeAreaRepeats = header->GetActiveAreaRepeats(); 
234
235   // now create the object that holds info about ghosts
236
237  
238
239   fastjet::GhostedAreaSpec ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea);
240   // and from that get an area definition
241   fastjet::AreaType areaType = header->GetAreaType();
242   areaDef = fastjet::AreaDefinition(areaType,ghost_spec);
243   //cout<<"rParamBkg="<<rParamBkg<<"  ghostEtamax="<<ghostEtamax<<"  ghostArea="<<ghostArea<<" areadef="<<TString(areaDef.description())<<endl;
244   //fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef);
245   fastjet::ClusterSequenceArea clust_seq(inputParticles,jetDef,areaDef);
246   TString comment = "Running FastJet algorithm for BKG calculation with the following setup. ";
247   comment+= "Jet definition: ";
248   comment+= TString(jetDef.description());
249   //  comment+= ". Area definition: ";
250   //comment+= TString(areaDef.description());
251   comment+= ". Strategy adopted by FastJet: ";
252   comment+= TString(clust_seq.strategy_string());
253   header->SetComment(comment);
254   if(debug){
255     cout << "--------------------------------------------------------" << endl;
256     cout << comment << endl;
257     cout << "--------------------------------------------------------" << endl;
258   }
259
260   vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(0.);
261   vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets); 
262   vector<fastjet::PseudoJet> jets2=sorted_by_pt(inclusiveJets);
263   if(jets2.size()>=2) jets2.erase(jets2.begin(),jets2.begin()+1);
264   
265    
266   
267   double phiMin = 0, phiMax = 0, rapMin = 0, rapMax = 0;
268
269  
270     phiMin = 0;
271     phiMax = 2*TMath::Pi();
272  
273   rapMax = ghostEtamax - rParamBkg;
274   rapMin = - ghostEtamax + rParamBkg;
275
276
277   fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
278
279
280   double medianb, sigmab, meanareab;
281   clust_seq.get_median_rho_and_sigma(jets2, range, false, medianb, sigmab, 
282 meanareab, false);
283   median=medianb;
284   sigma=sigmab;
285   meanarea=meanareab; 
286
287   
288  }
289
290
291 //___________________________________________________________________
292 Float_t AliJetBkg::BkgFastJet(){
293   
294   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
295   Bool_t debug  = header->GetDebug();     // debug option
296
297   if(debug)cout<<"===============  AliJetBkg::BkgFastJet()  =========== "<<endl;
298   vector<fastjet::PseudoJet> inputParticles=fInputFJ->GetInputParticles();
299   
300   if(debug)cout<<"printing inputParticles for BKG "<<inputParticles.size()<<endl;
301   
302   for(UInt_t i=0;i<inputParticles.size();i++){
303     //  cout<<"                "<<inputParticles[i].px()<<" "<<inputParticles[i].py()<<" "<<inputParticles[i].pz()<<endl;
304     
305   }
306    
307   double rParamBkg = header->GetRparamBkg(); //Radius for background calculation
308   Double_t rho=CalcRho(inputParticles,rParamBkg,"All");
309   if(debug)cout<<"-------- rho (from all part)="<<rho<<endl; 
310   return rho;
311  
312 }
313 //___________________________________________________________________
314 Float_t  AliJetBkg::BkgChargedFastJet(){
315
316   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
317   Bool_t debug  = header->GetDebug();     // debug option
318
319   if(debug)cout<<"===============  AliJetBkg::BkgChargedFastJet()  =========== "<<endl;
320
321   vector<fastjet::PseudoJet> inputParticlesCharged=fInputFJ->GetInputParticlesCh();
322   
323   if(debug)cout<<"printing CHARGED inputParticles for BKG "<<inputParticlesCharged.size()<<endl;
324
325   for(UInt_t i=0;i<inputParticlesCharged.size();i++){
326     // cout<<"                "<<inputParticlesCharged[i].px()<<" "<<inputParticlesCharged[i].py()<<" "<<inputParticlesCharged[i].pz()<<endl;
327
328   } 
329
330   double rParam = header->GetRparam();
331
332   Double_t rho=CalcRho(inputParticlesCharged,rParam,"Charg");
333
334   cout<<"-------- rho (from CHARGED part)="<<rho<<endl; 
335   return rho;
336 }
337
338
339
340 Float_t AliJetBkg::BkgStat()
341 {
342   //background subtraction using statistical method
343  
344   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
345   Bool_t debug  = header->GetDebug();     // debug option
346
347   if(debug)cout<<"==============AliJetBkg::BkgStat()============="<<endl;
348   //TO BE IMPLEMENTED 
349   // Int_t nTracks= fReader->GetESD()->GetNumberOfTracks();
350   Int_t nTracks= 0;
351   TF1 fun("fun",BkgFunction,0,800,1);
352   Double_t enTot=fun.Eval(nTracks);
353   Double_t accEMCal=2*0.7*110./180*TMath::Pi();//2.68 area of EMCal
354   return enTot/accEMCal;
355
356 }
357
358 ////////////////////////////////////////////////////////////////////////
359 Float_t AliJetBkg::BkgFastJetCone(TClonesArray* fAODJets)
360 {
361
362   // Cone background subtraction method applied on the fastjet: REmove the particles of the
363   // two largest jets with the given R from the estimation of new rho. 
364
365   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
366   Bool_t debug  = header->GetDebug();     // debug option
367
368   if(debug)cout<<"==============AliJetBkg::SubtractFastJetBackgCone()============="<<endl;
369
370   Float_t rc= header->GetRparam();
371
372   //Hard wired Calorimeter area (get it later from the AliJetReaderHeader.h)
373   Double_t accEMCal=2*0.7*110./180*TMath::Pi();//2.68 area of EMCal
374
375   Int_t nJ=fAODJets->GetEntries(); //this must be the # of jets... 
376   if(debug)cout<<"nJets:  "<<nJ<<endl;
377   
378
379   //begin unit array 
380   TClonesArray* fUnit = fReader->GetUnitArray();
381   if(fUnit == 0) { cout << "Could not get the momentum array" << endl; return -99; }
382
383   Int_t nIn = fUnit->GetEntries();
384   if(nIn == 0) { cout << "entries = 0 ; Event empty !!!" << endl ; return -99; }
385   
386   // Information extracted from fUnitArray
387   // load input vectors and calculate total energy in array
388   Float_t pt,eta,phi;
389   Float_t jeteta = 0,jetphi = 0,jeteta1 = 0, jetphi1 = 0;
390   Float_t rhoback=0.0;
391
392   Float_t ptallback=0.0; //particles without the jet
393   Float_t restarea=accEMCal; //initial area set 
394   Bool_t acc=0;
395   Bool_t acc1=0;
396   Float_t rCone=0.4;
397   
398   if(nJ==1) { 
399     AliAODJet *jettmp = dynamic_cast<AliAODJet*>(fAODJets->At(0));
400     jeteta=jettmp->Eta();
401     jetphi=jettmp->Phi();
402     acc=EmcalAcceptance(jeteta,jetphi,rCone);
403     if(acc==1)restarea= accEMCal-TMath::Pi()*rc*rc;
404     if(debug)cout<<" acc  "<<acc<<endl;
405   }
406
407   
408   if(nJ>=2) { 
409     AliAODJet *jettmp = dynamic_cast<AliAODJet*>(fAODJets->At(0));
410     AliAODJet *jettmp1 = dynamic_cast<AliAODJet*>(fAODJets->At(1));
411     jeteta=jettmp->Eta();
412     jetphi=jettmp->Phi();
413     jeteta1=jettmp1->Eta();
414     jetphi1=jettmp1->Phi(); 
415     acc=EmcalAcceptance(jeteta,jetphi,rCone);
416     acc1=EmcalAcceptance(jeteta1,jetphi1,rCone);
417     if(acc1==1 && acc==1)restarea= accEMCal-2*TMath::Pi()*rc*rc;
418     if(acc1==1 && acc==0)restarea= accEMCal-TMath::Pi()*rc*rc;
419     if(acc1==0 && acc==1)restarea= accEMCal-TMath::Pi()*rc*rc;
420
421     if(debug)cout<<" acc1="<<acc<<"  acc2="<<acc1<<"  restarea="<<restarea<<endl;
422
423   }
424   
425   // cout<<" nIn = "<<nIn<<endl;
426   Float_t sumpt=0;
427   for(Int_t i=0; i<nIn; i++) 
428     { //Unit Array Loop
429       AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(i);
430       if(uArray->GetUnitEnergy()>0.){
431         
432         pt    = uArray->GetUnitEnergy();
433         eta   = uArray->GetUnitEta();
434         phi   = uArray->GetUnitPhi();
435
436         //cout<<"test emcal acceptance for particles "<<EmcalAcceptance(eta,phi,0.)<<endl;
437         
438         Float_t deta=0.0, dphi=0.0, dr=100.0;
439         Float_t deta1=0.0, dphi1=0.0, dr1=100.0;
440         
441         //cout<<i<<"  pt="<<pt<<"  eta="<<eta<<"  phi="<<phi<<endl;
442         if(phi>1.396 && phi<3.316 && eta>-0.7 && eta<0.7){
443           sumpt+=pt;
444           //if(i<30)cout<<i<<" pt unit = "<<pt<<endl;
445
446           if(nJ==1 && acc==1) { 
447             deta = eta - jeteta;
448             dphi = phi - jetphi;
449             if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
450             if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
451             dr = TMath::Sqrt(deta * deta + dphi * dphi);
452             if(dr<=rc)sumpt-=pt;
453           }
454         
455           if(nJ>=2) { 
456             if(acc==1){
457               deta = eta - jeteta;
458               dphi = phi - jetphi;
459               if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
460               if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
461               dr = TMath::Sqrt(deta * deta + dphi * dphi);
462               if(dr<=rc)sumpt-=pt;
463             }
464             if(acc1==1){
465               deta1 = eta - jeteta1;
466               dphi1 = phi - jetphi1;
467               if (dphi1 < -TMath::Pi()) dphi1= -dphi1 - 2.0 * TMath::Pi();
468               if (dphi1 > TMath::Pi()) dphi1 = 2.0 * TMath::Pi() - dphi1;
469               dr1 = TMath::Sqrt(deta1 * deta1 + dphi1 * dphi1);
470               if(dr1<=rc)sumpt-=pt;
471             }
472
473           }
474
475           if(dr >= rc && dr1 >=rc) { 
476             // particles outside both cones
477           
478             //cout<<" out of the cone  "<<dr<<"    "<<deta<<"  deltaeta  "<<dphi<<"  dphi "<<i<<"  particle  "<<endl;
479             //cout<<" out of the cone  "<<dr1<<"    "<<deta1<<"  deltaeta1  "<<dphi1<<"  dphi1 "<<i<<"  particle  "<<endl;
480             ptallback+=pt;
481           }
482         }
483         //cout<<"  ipart  "<<ipart<<"  rhointegral  "<<rhoback <<endl;
484       }         
485     } // End loop on UnitArray 
486   
487   if(debug)cout<<"total area left "<<restarea<<endl;
488   if(debug)cout<<"sumpt="<<sumpt<<endl;
489   // if(acc==1 || acc1==1) rhoback= ptallback/restarea;
490   //else rhoback=ptallback;
491
492   rhoback= ptallback/restarea;
493   if(debug)cout<<"rhoback    "<<rhoback<<"     "<<nJ<<"   "<<endl;
494
495   return rhoback;
496    
497 }
498
499
500 Double_t AliJetBkg::CalcRho(vector<fastjet::PseudoJet> inputParticles,Double_t rParamBkg,TString method){
501   //calculate rho using the fastjet method
502
503   AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
504   Bool_t debug  = header->GetDebug();     // debug option
505
506   fastjet::Strategy strategy = header->GetStrategy();
507   fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
508   fastjet::JetAlgorithm algorithm = header->GetAlgorithm(); 
509   fastjet::JetDefinition jetDef(algorithm, rParamBkg, recombScheme, strategy);
510
511   // create an object that specifies how we to define the area
512   fastjet::AreaDefinition areaDef;
513   double ghostEtamax = header->GetGhostEtaMax(); 
514   double ghostArea   = header->GetGhostArea(); 
515   int    activeAreaRepeats = header->GetActiveAreaRepeats(); 
516   
517   // now create the object that holds info about ghosts
518
519   if (method.Contains("Charg"))ghostEtamax=0.9;
520
521   fastjet::GhostedAreaSpec ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea);
522   // and from that get an area definition
523   fastjet::AreaType areaType = header->GetAreaType();
524   areaDef = fastjet::AreaDefinition(areaType,ghost_spec);
525   if(debug)cout<<"rParamBkg="<<rParamBkg<<"  ghostEtamax="<<ghostEtamax<<"  ghostArea="<<ghostArea<<" areadef="<<TString(areaDef.description())<<endl;
526   //fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef);
527   fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef,areaDef);
528   TString comment = "Running FastJet algorithm for BKG calculation with the following setup. ";
529   comment+= "Jet definition: ";
530   comment+= TString(jetDef.description());
531   //  comment+= ". Area definition: ";
532   //comment+= TString(areaDef.description());
533   comment+= ". Strategy adopted by FastJet: ";
534   comment+= TString(clust_seq.strategy_string());
535   header->SetComment(comment);
536   if(debug){
537     cout << "--------------------------------------------------------" << endl;
538     cout << comment << endl;
539     cout << "--------------------------------------------------------" << endl;
540   }
541
542   double ptmin = header->GetPtMin(); 
543   vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
544   vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets); 
545
546   if (debug) {
547       cout<<"# of BKG jets = "<<jets.size()<<endl;
548       for (size_t j = 0; j < jets.size(); j++) { // loop for jets   
549           printf("BKG Jet found %5d %9.5f %8.5f %10.3f %4.4f \n",(Int_t)j,jets[j].rap(),jets[j].phi(),jets[j].perp(),clust_seq.area(jets[j]));
550       }
551   }
552   
553   double phiMin = 0, phiMax = 0, rapMin = 0, rapMax = 0;
554
555   if (method.Contains("All")){
556     phiMin = 80.*TMath::Pi()/180+rParamBkg;
557     phiMax = 190.*TMath::Pi()/180-rParamBkg;
558   }
559   if (method.Contains("Charg")){
560     phiMin = 0;
561     phiMax = 2*TMath::Pi();
562   }
563   rapMax = ghostEtamax - rParamBkg;
564   rapMin = - ghostEtamax + rParamBkg;
565
566
567   fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
568
569
570   Double_t rho=clust_seq.median_pt_per_unit_area(range);
571   //    double median, sigma, meanArea;
572   //clust_seq.get_median_rho_and_sigma(inclusiveJets, range, false, median, sigma, meanArea, true);
573   //fastjet::ActiveAreaSpec area_spec(ghostEtamax,activeAreaRepeats,ghostArea);
574
575   // fastjet::ClusterSequenceActiveArea clust_seq_bkg(inputParticles, jetDef,area_spec);
576
577
578   if(debug)cout<<"bkg in R="<<rParamBkg<<"  : "<<rho<<" range: Rap="<<rapMin<<","<<rapMax<<" --  phi="<<phiMin<<","<<phiMax<<endl;
579   return rho;
580 }
581
582 Float_t  AliJetBkg::EtaToTheta(Float_t arg)
583 {
584   //  return (180./TMath::Pi())*2.*atan(exp(-arg));
585   return 2.*atan(exp(-arg));
586 }
587
588
589 Double_t  AliJetBkg::BkgFunction(Double_t */*x*/,Double_t */*par*/){
590   //to be implemented--- (pT + Energy in EMCal Acceptance vs Multiplicity)
591   return 1;
592 }
593
594
595 Bool_t AliJetBkg::EmcalAcceptance(const Float_t eta, const Float_t phi, const Float_t radius){
596  
597   Float_t meanPhi=190./180.*TMath::Pi()-110./180.*TMath::Pi()/2;
598   Float_t deltaphi=110./180.*TMath::Pi();
599   Float_t phicut=deltaphi/2.-radius;
600   Float_t etacut=0.7-radius;
601   //cout<<"  eta    "<<eta<<"  phi    "<<phi<<endl;
602   //cout<<etacut<<"    "<<phicut<<"    "<<meanPhi<<"    "<<endl;
603   if(TMath::Abs(eta)<etacut && TMath::Abs(phi-meanPhi)<phicut) return 1;
604   else return 0; 
605 }