]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - FMD/flow/AliFMDFlowResolution.cxx
- Adapt default setting for ESDpid object
[u/mrichter/AliRoot.git] / FMD / flow / AliFMDFlowResolution.cxx
index b37d2d95e6bbf812ba8008aa3d02b71b70a52015..934bf46d5892d929408842751e483c2df4e2a229 100644 (file)
 #include <TH1.h>
 #include <TH2.h>
 #include <iostream>
+#include <TBrowser.h>
 //#include <cmath>
 
 //====================================================================
+AliFMDFlowResolution::AliFMDFlowResolution(UShort_t n) 
+  : fOrder(n), 
+    fContrib("contrib", Form("cos(%d(#Psi_{A} - #Psi_{B}))", fOrder), 100,-1,1)
+{
+  fContrib.SetDirectory(0);
+  fContrib.SetXTitle(Form("cos(%d(#Psi_{A} - #Psi_{B}))", fOrder));
+}
+
+//____________________________________________________________________
 AliFMDFlowResolution::AliFMDFlowResolution(const AliFMDFlowResolution& o)
   : AliFMDFlowStat(o), 
-    fOrder(o.fOrder)
+    fOrder(o.fOrder),
+    fContrib(o.fContrib)
 {
   // Copy constructor 
   // Parameters: 
   //   o   Object to copy from. 
+  fContrib.SetDirectory(0);
+  fContrib.SetXTitle(Form("cos(%d(#Psi_{A} - #Psi_{B}", fOrder));
 }
 //____________________________________________________________________
 AliFMDFlowResolution&
@@ -53,6 +66,10 @@ AliFMDFlowResolution::operator=(const AliFMDFlowResolution& o)
   // 
   AliFMDFlowStat::operator=(o);
   fOrder = o.fOrder;
+
+  fContrib.Reset();
+  fContrib.Add(&o.fContrib);
+
   return *this;
 }
 //____________________________________________________________________
@@ -66,6 +83,7 @@ AliFMDFlowResolution::Add(Double_t psiA, Double_t psiB)
   Double_t diff    = NormalizeAngle(fOrder * (psiA - psiB));
   Double_t contrib = cos(diff);
   AliFMDFlowStat::Add(contrib);
+  fContrib.Fill(contrib);
 }
 
 //____________________________________________________________________
@@ -89,7 +107,7 @@ AliFMDFlowResolution::Correction(UShort_t, Double_t& e2) const
   //  e2 The square error on the correction 
   // Returns <cos(n(psi - psi_R))>
   e2 = fSqVar / fN;
-  return sqrt(2) * sqrt(fabs(fAverage));
+  return sqrt(Double_t(2)) * sqrt(fabs(fAverage));
 }
 
 //____________________________________________________________________
@@ -102,7 +120,7 @@ AliFMDFlowResolution::Draw(Option_t* option)
   TGraph* g = new TGraph(100);
   for (UShort_t i = 0; i < g->GetN(); i++) { 
     Double_t x = -1. + 2. / 100 * i;
-    Double_t y = sqrt(2) * sqrt(fabs(x));
+    Double_t y = sqrt(Double_t(2)) * sqrt(fabs(x));
     g->SetPoint(i, x, y);
   }
   g->SetName("naive_res");
@@ -113,6 +131,12 @@ AliFMDFlowResolution::Draw(Option_t* option)
   g->Draw(Form("lh %s", option));
 }
 
+//____________________________________________________________________
+void
+AliFMDFlowResolution::Browse(TBrowser* b)
+{
+  b->Add(&fContrib);
+}
 
 //====================================================================
 Double_t 
@@ -127,7 +151,7 @@ AliFMDFlowResolutionStar::Correction(UShort_t k, Double_t& e2) const
   Double_t delta = 0;
   Double_t chi   = Chi(fAverage, k, delta);
   Double_t dr    = 0;
-  Double_t res   = Res(sqrt(2) * chi, k, dr);
+  Double_t res   = Res(sqrt(Double_t(2)) * chi, k, dr);
   e2             = pow(dr * delta,2);
   return res;
 }
@@ -202,7 +226,7 @@ AliFMDFlowResolutionStar::Res(Double_t chi, UShort_t k, Double_t& dr) const
   //            2   
   // 
   Double_t y  = chi * chi / 4;
-  Double_t c  = sqrt(M_PI) * exp(-y) / (2 * sqrt(2));   
+  Double_t c  = sqrt(M_PI) * exp(-y) / (2 * sqrt(Double_t(2)));   
   
   // i[0] = I[(k-1)/2-1], i[1] = I[(k-1)/2], i[2] = I[(k-1)/2+1]
   Double_t i[3], di[3];
@@ -249,7 +273,7 @@ AliFMDFlowResolutionStar::Draw(Option_t* option)
       if (chi) y  = c;
       else { 
        Double_t dr = 0;
-       y           = Res(sqrt(2) * c, k, dr);
+       y           = Res(sqrt(Double_t(2)) * c, k, dr);
        e2          = pow(dr * e2,2);
       }
       g->SetPoint(i, x, y);
@@ -267,11 +291,29 @@ AliFMDFlowResolutionStar::Draw(Option_t* option)
 }
 
 //====================================================================
+AliFMDFlowResolutionTDR::AliFMDFlowResolutionTDR(UShort_t n) 
+    : AliFMDFlowResolution(n), 
+      fLarge(0),
+      fNK("nk", "Number of events with (#Psi_{}A-#Psi_{B})>90^{o}", 2, 0, 2)
+{
+  fNK.SetDirectory(0);
+  fNK.SetYTitle("# events");
+  fNK.GetXaxis()->SetBinLabel(1, "k");
+  fNK.GetXaxis()->SetBinLabel(2, "N");
+}
+
+//____________________________________________________________________
 AliFMDFlowResolutionTDR::AliFMDFlowResolutionTDR(const 
                                                 AliFMDFlowResolutionTDR& o)
   : AliFMDFlowResolution(o), 
-    fLarge(o.fLarge)
-{}
+    fLarge(o.fLarge),
+    fNK(o.fNK)
+{
+  fNK.SetDirectory(0);
+  fNK.SetYTitle("# events");
+  fNK.GetXaxis()->SetBinLabel(1, "k");
+  fNK.GetXaxis()->SetBinLabel(2, "N");
+}
 //____________________________________________________________________
 AliFMDFlowResolutionTDR&
 AliFMDFlowResolutionTDR::operator=(const AliFMDFlowResolutionTDR& o)
@@ -282,6 +324,8 @@ AliFMDFlowResolutionTDR::operator=(const AliFMDFlowResolutionTDR& o)
   // 
   AliFMDFlowResolution::operator=(o);
   fLarge = o.fLarge;
+  fNK.Reset();
+  fNK.Add(&fNK);
   return *this;
 }
 //____________________________________________________________________
@@ -303,8 +347,13 @@ AliFMDFlowResolutionTDR::Add(Double_t psiA, Double_t psiB)
   // Parameters: 
   //  psiA   A sub-event plane angle Psi_A in [0,2pi]
   //  psiB   B sub-event plane angle Psi_B in [0,2pi]
+  AliFMDFlowResolution::Add(psiA, psiB);
   Double_t a = fabs(psiA - psiB);
-  if (a >= .5 * M_PI) fLarge++;
+  if (a >= .5 * M_PI) { 
+    fNK.Fill(.5);
+    fLarge++;
+  }
+  fNK.Fill(1.5);
   fN++;
 }
 //____________________________________________________________________
@@ -323,11 +372,11 @@ AliFMDFlowResolutionTDR::Correction(UShort_t k, Double_t& e2) const
   // where z = chi^2 / 2
   //
   if (fLarge == 0) { 
-    std::cerr << "TDR: K = 0" << std::endl;
+    // std::cerr << "TDR: K = 0" << std::endl;
     return -1;
   }
   if (fN == 0) { 
-    std::cerr << "TDR: N = 0" << std::endl;
+    // std::cerr << "TDR: N = 0" << std::endl;
     return -1;
   }
   Double_t r     = Double_t(fLarge) / fN;
@@ -546,6 +595,14 @@ AliFMDFlowResolutionTDR::Draw(Option_t* option)
   }
 }
 
+//____________________________________________________________________
+void
+AliFMDFlowResolutionTDR::Browse(TBrowser* b)
+{
+  AliFMDFlowResolution::Browse(b);
+  b->Add(&fNK);
+}
+
 //____________________________________________________________________
 //
 // EOF