]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - JETAN/AliFastJetHeader.cxx
Treatment of bad SDD half-modules (F. Prino)
[u/mrichter/AliRoot.git] / JETAN / AliFastJetHeader.cxx
index ecfe7f8983dc07def7ff7d559d44c42dd62caa33..dcc888edad94a6b3e5208eb1e94af08aa0e5ed32 100644 (file)
  * provided "as is" without express or implied warranty.                  *
  **************************************************************************/
 
+
 //---------------------------------------------------------------------
-//  FastJet header class
-// Stores parameters of particle algoritm
+// FastJet v2.3.4 finder algorithm interface
+// Finder Header Class 
 // Author: Rafael.Diaz.Valdes@cern.ch
 //---------------------------------------------------------------------
 
 #include <Riostream.h>
 #include <TMath.h>
+
+#include "fastjet/ClusterSequenceArea.hh"
+#include "fastjet/AreaDefinition.hh"
+#include "fastjet/JetDefinition.hh"
+
 #include "AliFastJetHeader.h"
 
 ClassImp(AliFastJetHeader)
@@ -29,18 +35,26 @@ ClassImp(AliFastJetHeader)
 
 AliFastJetHeader::AliFastJetHeader():
     AliJetHeader("AliFastJetHeader"),
-    fLegoNbinEta(60),
-    fLegoNbinPhi(210),
-    fLegoEtaMin(-0.9),
-    fLegoEtaMax(0.9),
-    fLegoPhiMin(0.),
-    fLegoPhiMax(2. * TMath::Pi()),
-    fRadius(1.0),
-    fMinJetEt(10.0),
-    fSoftBackg(kTRUE),
-    fPrecBg(0.035) 
+    fRparam(1.0), 
+    fAlgorithm(fastjet::kt_algorithm),
+    fStrategy(fastjet::Best),
+    fRecombScheme(fastjet::BIpt_scheme),
+    fGhostEtaMax(2.0),
+    fGhostArea(0.05),
+    fActiveAreaRepeats(1),
+    fAreaType(fastjet::active_area), 
+    fPtMin(5.0),
+    fRapMax(0),
+    fRapMin(0),
+    fPhiMax(TMath::TwoPi()),
+    fPhiMin(0)
 {
   // Constructor
+  
+  Double_t rapmax = fGhostEtaMax - fRparam;
+  Double_t rapmin = -fGhostEtaMax + fRparam;
+  SetRapRange(rapmin, rapmax);
+  
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -49,6 +63,27 @@ void AliFastJetHeader::PrintParameters() const
 {
   // prints out parameters of jet algorithm
 
-  cout << " FastJet algorithm " << endl;
+  cout << "FastJet algorithm  parameters:" << endl;
+  
+  cout << "-- Jet Definition --- " << endl;
+  cout << "R " << fRparam << endl;
+  cout << "Jet Algorithm " << fAlgorithm << endl; 
+  cout << "Strategy " << fStrategy << endl;  
+  cout << "Recombination Scheme " << fRecombScheme << endl; 
+  
+  cout << "-- Ghosted Area Spec parameters --- " << endl;
+  cout << "Ghost Eta Max " << fGhostEtaMax << endl;
+  cout << "Ghost Area " << fGhostArea << endl;
+  cout << "Active Area Repeats " << fActiveAreaRepeats << endl;
+  
+  cout << "-- Area Definition parameters --- " << endl;
+  cout << "Area Type " << fAreaType << endl; 
+  
+  cout << "-- Cluster Sequence Area parameters --- " << endl;
+  cout << "pt min " << fPtMin << endl; 
+  
+  cout << "-- Range Definition parameters --- " << endl;
+  cout << " bkg rapidity range from  " << fRapMin << " to " << fRapMax << endl;
+  cout << " bkg phi range from " << fPhiMin << " to " << fPhiMax << endl;
 
 }