]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
add QA histograms for HBT/conversion cuts
authormiweber <miweber@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 21 Jun 2012 07:20:16 +0000 (07:20 +0000)
committermiweber <miweber@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 21 Jun 2012 07:20:16 +0000 (07:20 +0000)
PWGCF/EBYE/BalanceFunctions/AliAnalysisTaskBF.cxx
PWGCF/EBYE/BalanceFunctions/AliBalance.cxx
PWGCF/EBYE/BalanceFunctions/AliBalance.h
PWGCF/EBYE/macros/readBalanceFunction.C

index 66c58c9a6810686d83d7f4cb5933cbd9ffeef453..a3ce34fd198c9ea9debded96c7cba23e91f078a4 100755 (executable)
@@ -265,6 +265,12 @@ void AliAnalysisTaskBF::UserCreateOutputObjects() {
     fHistRefTracks->GetXaxis()->SetBinLabel(i,gRefTrackName[i-1].Data());
   fList->Add(fHistRefTracks);
 
+  // QA histograms for HBTinspired and Conversion cuts
+  fList->Add(fBalance->GetQAHistHBTbefore());
+  fList->Add(fBalance->GetQAHistHBTafter());
+  fList->Add(fBalance->GetQAHistConversionbefore());
+  fList->Add(fBalance->GetQAHistConversionafter());
+
   // Balance function histograms
   // Initialize histograms if not done yet
   if(!fBalance->GetHistNp(0)){
@@ -297,7 +303,7 @@ void AliAnalysisTaskBF::UserCreateOutputObjects() {
       fListBFS->Add(fShuffledBalance->GetHistNpp(a));
       fListBFS->Add(fShuffledBalance->GetHistNnp(a));
     }  
-  }
+  }  
 
   if(fESDtrackCuts) fList->Add(fESDtrackCuts);
 
index f0b4f9a8144714e4c617030db382909fbf797418..67ccc90d3b6dddf9c293d807149004a8f2189e27 100644 (file)
@@ -93,6 +93,13 @@ AliBalance::AliBalance() :
     fHistNN[i] = NULL;
 
   }
+
+  //QA histograms
+  fHistHBTbefore = NULL;
+  fHistHBTafter  = NULL;
+  fHistConversionbefore = NULL;
+  fHistConversionafter  = NULL;
+
 }
 
 
@@ -222,6 +229,13 @@ void AliBalance::InitHistograms() {
     if(fCentralityId) histName += fCentralityId.Data();
     fHistNN[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,fNumberOfBins[iAnalysisType],fP2Start[iAnalysisType],fP2Stop[iAnalysisType]);
   }
+
+  // QA histograms
+  fHistHBTbefore        = new TH2D("fHistHBTbefore","before HBT cut",200,0,2,200,0,200);
+  fHistHBTafter         = new TH2D("fHistHBTafter","after HBT cut",200,0,2,200,0,200);
+  fHistConversionbefore = new TH2D("fHistConversionbefore","before Conversion cut",200,0,2,200,0,200);
+  fHistConversionafter  = new TH2D("fHistConversionafter","after Conversion cut",200,0,2,200,0,200);
+
 }
 
 //____________________________________________________________________//
@@ -401,6 +415,8 @@ void AliBalance::CalculateBalance(Float_t fCentrality,vector<Double_t> **chargeV
          // VERSION 2 (Taken from DPhiCorrelations)
          // the variables & cuthave been developed by the HBT group 
          // see e.g. https://indico.cern.ch/materialDisplay.py?contribId=36&sessionId=6&materialId=slides&confId=142700
+
+         fHistHBTbefore->Fill(deta,dphi);
          
          // optimization
          if (TMath::Abs(deta) < 0.02 * 2.5 * 3) //twoTrackEfficiencyCutValue = 0.02 [default for dphicorrelations]
@@ -440,12 +456,16 @@ void AliBalance::CalculateBalance(Float_t fCentrality,vector<Double_t> **chargeV
                    }
                }
            }
+         fHistHBTafter->Fill(deta,dphi);
        }
        
        // conversions
        if(bConversionCut){
          if (charge1 * charge2 < 0)
            {
+
+             fHistConversionbefore->Fill(deta,dphi);
+
              Float_t m0 = 0.510e-3;
              Float_t tantheta1 = 1e10;
 
@@ -469,6 +489,7 @@ void AliBalance::CalculateBalance(Float_t fCentrality,vector<Double_t> **chargeV
                //AliInfo(Form("Conversion: Removed track pair %d %d with [[%f %f] %f %f] %d %d <- %f %f  %f %f   %f %f ", i, j, deta, dphi, masssqu, charge1, charge2,eta1,eta2,phi1,phi2,pt1,pt2));
                continue;
              }
+             fHistConversionafter->Fill(deta,dphi);
            }
        }
        
index 36f1fdc176ae48f0a91a49f102966fe376bf0759..3fb901cf0c9029122f4b0ec3c50bbb4a4d738fcd 100644 (file)
@@ -92,6 +92,11 @@ class AliBalance : public TObject {
   TH2D *GetHistNpp(Int_t iAnalysisType) { return fHistPP[iAnalysisType];}
   TH2D *GetHistNnn(Int_t iAnalysisType) { return fHistNN[iAnalysisType];}
 
+  TH2D *GetQAHistHBTbefore()         {return fHistHBTbefore;};
+  TH2D *GetQAHistHBTafter()          {return fHistHBTafter;};
+  TH2D *GetQAHistConversionbefore()  {return fHistConversionbefore;};
+  TH2D *GetQAHistConversionafter()   {return fHistConversionafter;};
+
   void PrintAnalysisSettings();
   TGraphErrors *DrawBalance(Int_t fAnalysisType);
 
@@ -149,6 +154,12 @@ class AliBalance : public TObject {
   TH2D *fHistPP[ANALYSIS_TYPES]; //N++
   TH2D *fHistNN[ANALYSIS_TYPES]; //N--
 
+  //QA histograms
+  TH2D *fHistHBTbefore; // Delta Eta vs. Delta Phi before HBT inspired cuts
+  TH2D *fHistHBTafter; // Delta Eta vs. Delta Phi after HBT inspired cuts
+  TH2D *fHistConversionbefore; // Delta Eta vs. Delta Phi before Conversion cuts
+  TH2D *fHistConversionafter; // Delta Eta vs. Delta Phi before Conversion cuts
+
 
 
   AliBalance & operator=(const AliBalance & ) {return *this;}
index a69991958dd0fc96cf3dbcd084cb27d1a4c3c867..82ccfdb82779b72b3b0ec11497dbcde2eab68e91 100644 (file)
@@ -34,6 +34,7 @@ void drawBF(Bool_t bHistos = kFALSE, TString inFile = "AnalysisResults.root", In
   Int_t maximumCanvases = 13;
   Int_t iCanvas         = 0;
   Int_t iList           = -1;
+  TCanvas *cQACuts[13][10];
   TCanvas *cQA[13][10];
   TCanvas *cQAV0M = new TCanvas("cQAV0M","V0M multiplicities");
   cQAV0M->Divide(4,3);
@@ -229,6 +230,41 @@ void drawBF(Bool_t bHistos = kFALSE, TString inFile = "AnalysisResults.root", In
        histRef->Draw("colz");
       }
     }
+
+
+    cQACuts[iList][iCanvas] = new TCanvas(Form("%sCuts",listName.Data()),Form("%sCuts",listName.Data()));
+    cQACuts[iList][iCanvas]->Divide(2,2);
+
+    cQACuts[iList][iCanvas]->cd(1);
+    TH2D* histHBTcuts = (TH2D*)list->FindObject("fHistHBTafter");
+    if(histHBTcuts){
+      histHBTcuts->DrawCopy("colz");
+    }
+
+    cQACuts[iList][iCanvas]->cd(2);
+    TH2D* histHBTall = (TH2D*)list->FindObject("fHistHBTbefore");
+    if(histHBTall){
+      histHBTcuts->Divide(histHBTall);
+      histHBTcuts->SetMaximum(1.5);
+      histHBTcuts->SetMinimum(0.5);
+      histHBTcuts->DrawCopy("colz");
+    }
+
+    cQACuts[iList][iCanvas]->cd(3);
+    TH2D* histConversioncuts = (TH2D*)list->FindObject("fHistConversionafter");
+    if(histConversioncuts){
+      histConversioncuts->DrawCopy("colz");
+    }
+
+    cQACuts[iList][iCanvas]->cd(4);
+    TH2D* histConversionall = (TH2D*)list->FindObject("fHistConversionbefore");
+    if(histConversionall){
+      histConversioncuts->Divide(histConversionall);
+      histConversioncuts->SetMaximum(1.5);
+      histConversioncuts->SetMinimum(0.5);
+      histConversioncuts->DrawCopy("colz");
+    }
+
     // ----------------------------------------------------
 
     // ----------------------------------------------------
@@ -311,16 +347,29 @@ void drawBF(Bool_t bHistos = kFALSE, TString inFile = "AnalysisResults.root", In
          }
          else{
            fHistPN[iList][a]->SetLineColor(2);
-           fHistPN[iList][a]->ProjectionY(Form("pn%d",a))->DrawCopy();
+           TH1D *fHistTmp = (TH1D*)fHistPN[iList][a]->ProjectionY(Form("pn%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1]);
+           fHistTmp->Scale(0.5);
+           fHistTmp->DrawCopy();
            fHistPP[iList][a]->SetLineColor(1);
-           fHistPP[iList][a]->ProjectionY(Form("pp%d",a))->DrawCopy("same");
+           fHistPP[iList][a]->ProjectionY(Form("pp%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
            fHistNP[iList][a]->SetLineColor(4);
-           fHistNP[iList][a]->ProjectionY(Form("np%d",a))->DrawCopy("same");
+           fHistNP[iList][a]->ProjectionY(Form("np%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
            fHistNN[iList][a]->SetLineColor(8);
-           fHistNN[iList][a]->ProjectionY(Form("nn%d",a))->DrawCopy("same");
+           fHistNN[iList][a]->ProjectionY(Form("nn%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
          }
        }
       }
+
+      if(bHistos){
+       for(iCanvas = 0; iCanvas < nrOfCentralities; iCanvas++){
+         cBF[iList][iCanvas]->cd(8);
+         fHistP[iList][1]->ProjectionY(Form("p%d",1),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy();
+         fHistN[iList][1]->ProjectionY(Form("n%d",1),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
+         cBF[iList][iCanvas]->cd(9);
+         fHistP[iList][6]->ProjectionY(Form("p%d",6),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy();
+         fHistN[iList][6]->ProjectionY(Form("n%d",6),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
+       }
+      }
     }
     // ----------------------------------------------------
 
@@ -400,16 +449,29 @@ void drawBF(Bool_t bHistos = kFALSE, TString inFile = "AnalysisResults.root", In
          }
          else{
            fHistPNS[iList][a]->SetLineColor(2);
-           fHistPNS[iList][a]->ProjectionY(Form("pns%d",a))->DrawCopy();
+           TH1D *fHistTmp = (TH1D*)fHistPN[iList][a]->ProjectionY(Form("pns%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1]);
+           fHistTmp->Scale(0.5);
+           fHistTmp->DrawCopy();
            fHistPPS[iList][a]->SetLineColor(1);
-           fHistPPS[iList][a]->ProjectionY(Form("pps%d",a))->DrawCopy("same");
+           fHistPPS[iList][a]->ProjectionY(Form("pps%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
            fHistNPS[iList][a]->SetLineColor(4);
-           fHistNPS[iList][a]->ProjectionY(Form("nps%d",a))->DrawCopy("same");
+           fHistNPS[iList][a]->ProjectionY(Form("nps%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
            fHistNNS[iList][a]->SetLineColor(8);
-           fHistNNS[iList][a]->ProjectionY(Form("nns%d",a))->DrawCopy("same");
+           fHistNNS[iList][a]->ProjectionY(Form("nns%d",a),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
          }
        }
       }
+
+      if(bHistos){
+       for(iCanvas = 0; iCanvas < nrOfCentralities; iCanvas++){
+         cBFS[iList][iCanvas]->cd(8);
+         fHistPS[iList][1]->ProjectionY(Form("pS%d",1),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy();
+         fHistNS[iList][1]->ProjectionY(Form("nS%d",1),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
+         cBFS[iList][iCanvas]->cd(9);
+         fHistPS[iList][6]->ProjectionY(Form("pS%d",6),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy();
+         fHistNS[iList][6]->ProjectionY(Form("nS%d",6),centralityArray[iCanvas]+1,centralityArray[iCanvas+1])->DrawCopy("same");
+       }
+      }
     }
     // ----------------------------------------------------
   }