]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Updated version of the flattening task, adding A, C and A&C event planes.
authorcvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 19 Mar 2012 12:55:43 +0000 (12:55 +0000)
committercvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 19 Mar 2012 12:55:43 +0000 (12:55 +0000)
PWGPP/VZERO/AliAnaVZEROEPFlatenning.cxx
PWGPP/VZERO/AliAnaVZEROEPFlatenning.h

index 26da564c4e1f25a9d560bd5fa4e431af83492734..e9bb9780a61d55cd000ab888734ec0092843cec8 100644 (file)
@@ -27,15 +27,20 @@ AliAnaVZEROEPFlatenning::AliAnaVZEROEPFlatenning()
   : AliAnalysisTaskSE("AliAnaVZEROEPFlatenning"), fEvent(0), fOutputList(0),
   fMBTrigName("CPBI"),
   fUsePhysSel(kFALSE),
-  fPsiAC(NULL),
-  fPsiACOrg(NULL)
+  fPsiARawCentr(NULL),
+  fPsiAFlatCentr(NULL),
+  fPsiCRawCentr(NULL),
+  fPsiCFlatCentr(NULL),
+  fPsiACRawCentr(NULL),
+  fPsiACFlatCentr(NULL)
 {
   // Default constructor
   // Init pointers
-  for(Int_t i = 0; i < 8; ++i) fX2[i] = fY2[i] = fX2Y2[i] = fCos8Psi[i] = NULL;
+  for(Int_t i = 0; i < 11; ++i) fX2[i] = fY2[i] = fX2Y2[i] = fCos8Psi[i] = NULL;
   for(Int_t i = 0; i < 8; ++i) fX2In[i] = fY2In[i] = fX2Y2In[i] = fCos8PsiIn[i] = NULL;
   for(Int_t i = 0; i < 8; ++i) fPsiRingRawCentr[i] = fPsiRingFlatCentr[i] = fPsiRingFlatFinalCentr[i] = NULL;
   for(Int_t i = 0; i < 8; ++i) fC2[i] = fS2[i] = fC4[i] = fS4[i] = NULL;
+  for(Int_t i = 0; i < 11; ++i) fX2Corr[i] = fY2Corr[i] = fX2Y2Corr[i] = NULL;
   // Define input and output slots here
   // Input slot #0 works with a TChain
   DefineInput(0, TChain::Class());
@@ -49,15 +54,20 @@ AliAnaVZEROEPFlatenning::AliAnaVZEROEPFlatenning(const char *name)
   : AliAnalysisTaskSE(name), fEvent(0), fOutputList(0),
   fMBTrigName("CPBI"),
   fUsePhysSel(kFALSE),
-  fPsiAC(NULL),
-  fPsiACOrg(NULL)
+  fPsiARawCentr(NULL),
+  fPsiAFlatCentr(NULL),
+  fPsiCRawCentr(NULL),
+  fPsiCFlatCentr(NULL),
+  fPsiACRawCentr(NULL),
+  fPsiACFlatCentr(NULL)
 {
   // Constructor
   // Init pointers
-  for(Int_t i = 0; i < 8; ++i) fX2[i] = fY2[i] = fX2Y2[i] = fCos8Psi[i] = NULL;
+  for(Int_t i = 0; i < 11; ++i) fX2[i] = fY2[i] = fX2Y2[i] = fCos8Psi[i] = NULL;
   for(Int_t i = 0; i < 8; ++i) fX2In[i] = fY2In[i] = fX2Y2In[i] = fCos8PsiIn[i] = NULL;
   for(Int_t i = 0; i < 8; ++i) fPsiRingRawCentr[i] = fPsiRingFlatCentr[i] = fPsiRingFlatFinalCentr[i] = NULL;
   for(Int_t i = 0; i < 8; ++i) fC2[i] = fS2[i] = fC4[i] = fS4[i] = NULL;
+  for(Int_t i = 0; i < 11; ++i) fX2Corr[i] = fY2Corr[i] = fX2Y2Corr[i] = NULL;
   // Define input and output slots here
   // Input slot #0 works with a TChain
   DefineInput(0, TChain::Class());
@@ -98,7 +108,7 @@ void AliAnaVZEROEPFlatenning::UserCreateOutputObjects()
   fOutputList = new TList();
   fOutputList->SetOwner(kTRUE);
 
-  for(Int_t i = 0; i < 8; ++i) {
+  for(Int_t i = 0; i < 11; ++i) {
     fX2[i] = new TProfile(Form("fX2_%d",i),"",21,0,105,"s");
     fOutputList->Add(fX2[i]);
     fY2[i] = new TProfile(Form("fY2_%d",i),"",21,0,105,"s");
@@ -108,6 +118,14 @@ void AliAnaVZEROEPFlatenning::UserCreateOutputObjects()
     fCos8Psi[i] = new TProfile(Form("fCos8Psi_%d",i),"",21,0,105,"s");
     fOutputList->Add(fCos8Psi[i]);
   }
+  for(Int_t i = 0; i < 11; ++i) {
+    fX2Corr[i] = new TProfile(Form("fX2Corr_%d",i),"",21,0,105,"s");
+    fOutputList->Add(fX2Corr[i]);
+    fY2Corr[i] = new TProfile(Form("fY2Corr_%d",i),"",21,0,105,"s");
+    fOutputList->Add(fY2Corr[i]);
+    fX2Y2Corr[i] = new TProfile(Form("fX2Y2Corr_%d",i),"",21,0,105,"s");
+    fOutputList->Add(fX2Y2Corr[i]);
+  }
 
   for(Int_t i = 0; i < 8; ++i) {
     fPsiRingRawCentr[i] = new TH2F(Form("fPsiRingRawCentr_%d",i),"",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
@@ -117,6 +135,18 @@ void AliAnaVZEROEPFlatenning::UserCreateOutputObjects()
     fPsiRingFlatFinalCentr[i] = new TH2F(Form("fPsiRingFlatFinalCentr_%d",i),"",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
     fOutputList->Add(fPsiRingFlatFinalCentr[i]);
   }
+  fPsiARawCentr = new TH2F("fPsiARawCentr","",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
+  fOutputList->Add(fPsiARawCentr);
+  fPsiAFlatCentr = new TH2F("fPsiAFlatCentr","",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
+  fOutputList->Add(fPsiAFlatCentr);
+  fPsiCRawCentr = new TH2F("fPsiCRawCentr","",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
+  fOutputList->Add(fPsiCRawCentr);
+  fPsiCFlatCentr = new TH2F("fPsiCFlatCentr","",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
+  fOutputList->Add(fPsiCFlatCentr);
+  fPsiACRawCentr = new TH2F("fPsiACRawCentr","",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
+  fOutputList->Add(fPsiACRawCentr);
+  fPsiACFlatCentr = new TH2F("fPsiACFlatCentr","",105,0,105,100,-TMath::Pi()/2,TMath::Pi()/2);
+  fOutputList->Add(fPsiACFlatCentr);
 
   for(Int_t i = 0; i < 8; ++i) {
     fC2[i] = new TProfile(Form("fC2_%d",i),"",21,0,105,"s");
@@ -129,11 +159,6 @@ void AliAnaVZEROEPFlatenning::UserCreateOutputObjects()
     fOutputList->Add(fS4[i]);
   }
 
-  fPsiAC = new TH2F("fPsiAC","",100,-TMath::Pi()/2,TMath::Pi()/2,100,-TMath::Pi()/2,TMath::Pi()/2);
-  fOutputList->Add(fPsiAC);
-  fPsiACOrg = new TH2F("fPsiACOrg","",100,-TMath::Pi()/2,TMath::Pi()/2,100,-TMath::Pi()/2,TMath::Pi()/2);
-  fOutputList->Add(fPsiACOrg);
-
   PostData(1, fOutputList);
 }
 
@@ -212,12 +237,13 @@ void AliAnaVZEROEPFlatenning::UserExec(Option_t *)
   }
 
   if (goodEvent) {
-    Double_t qxTierce[8],qyTierce[8];
+    Double_t totMult = 0, totMultA = 0, totMultC = 0;
+    Double_t c2A = 0, c2C = 0, s2A = 0, s2C = 0;
+    Double_t qxA = 0, qyA = 0;
+    Double_t qxC = 0, qyC = 0;
     for(Int_t iring = 0; iring < 8; ++iring) {
-      qxTierce[iring] = qyTierce[iring] = 0.;
-      Double_t c2 = 0;
-      Double_t s2 = 0;
-      Double_t totMult = 0;
+      Double_t ringMult = 0;
+      Double_t c2 = 0, s2 = 0;
       for(Int_t iCh = iring*8; iCh < (iring+1)*8; ++iCh) {
        Double_t phi = TMath::Pi()/8. + TMath::Pi()/4.*(iCh%8);
        c2 += fEvent->GetVZEROEqMultiplicity(iCh)*TMath::Cos(2.*phi);
@@ -228,10 +254,20 @@ void AliAnaVZEROEPFlatenning::UserExec(Option_t *)
          fC4[iring]->Fill(spdPercentile,TMath::Cos(4.*phi),fEvent->GetVZEROEqMultiplicity(iCh));
          fS4[iring]->Fill(spdPercentile,TMath::Sin(4.*phi),fEvent->GetVZEROEqMultiplicity(iCh));
        }
-       totMult += fEvent->GetVZEROEqMultiplicity(iCh);
+       ringMult += fEvent->GetVZEROEqMultiplicity(iCh);
+      }
+      if (ringMult < 1e-6) continue;
+      totMult += ringMult;
+      if (iring >= 4) {
+       totMultA += ringMult;
+       c2A += c2;
+       s2A += s2;
+      }
+      else {
+       totMultC += ringMult;
+       c2C += c2;
+       s2C += s2;
       }
-      if (totMult < 1e-6) continue;
-
       fX2[iring]->Fill(spdPercentile,c2);
       fY2[iring]->Fill(spdPercentile,s2);
       fX2Y2[iring]->Fill(spdPercentile,c2*s2);
@@ -239,26 +275,67 @@ void AliAnaVZEROEPFlatenning::UserExec(Option_t *)
       Double_t qxOut = 0, qyOut = 0;
       Double_t psiRingRaw = fEvent->GetEventplane()->CalculateVZEROEventPlane(fEvent,iring,iring,2,qxOut,qyOut);
       fPsiRingRawCentr[iring]->Fill(spdPercentile,psiRingRaw);
-      Double_t psiRingFlat2 = CalculateVZEROEventPlane(fEvent,iring,spdPercentile,qxTierce[iring],qyTierce[iring]);
-      fPsiRingFlatCentr[iring]->Fill(spdPercentile,psiRingFlat2);
-      fCos8Psi[iring]->Fill(spdPercentile, 2./4.*TMath::Cos(2.*4.*psiRingFlat2));
+      fCos8Psi[iring]->Fill(spdPercentile, 2./4.*TMath::Cos(2.*4.*psiRingRaw));
+      Double_t qxTierce = 0, qyTierce = 0;
+      Double_t psiRingFlat = CalculateVZEROEventPlane(fEvent,iring,spdPercentile,qxTierce,qyTierce);
+      fPsiRingFlatCentr[iring]->Fill(spdPercentile,psiRingFlat);
       Int_t ibin = fCos8PsiIn[iring]->FindBin(spdPercentile);
-      Double_t psiRingFlatFinal = psiRingFlat2 + (fCos8PsiIn[iring]->GetBinContent(ibin)*TMath::Sin(2.*4.*psiRingFlat2))/2.;
+      Double_t psiRingFlatFinal = psiRingFlat + (fCos8PsiIn[iring]->GetBinContent(ibin)*TMath::Sin(2.*4.*psiRingFlat))/2.;
       if (psiRingFlatFinal > TMath::Pi()/2) psiRingFlatFinal -= TMath::Pi();
       if (psiRingFlatFinal <-TMath::Pi()/2) psiRingFlatFinal += TMath::Pi();
       fPsiRingFlatFinalCentr[iring]->Fill(spdPercentile,psiRingFlatFinal);
+      fX2Corr[iring]->Fill(spdPercentile,qxTierce);
+      fY2Corr[iring]->Fill(spdPercentile,qyTierce);
+      fX2Y2Corr[iring]->Fill(spdPercentile,qxTierce*qyTierce);
+      if (iring >= 4) {
+       qxA += qxTierce;
+       qyA += qyTierce;
+      }
+      else {
+       qxC += qxTierce;
+       qyC += qyTierce;
+     }
+    }
+
+    if (totMultA > 1e-6) {
+      fX2[8]->Fill(spdPercentile,c2A);
+      fY2[8]->Fill(spdPercentile,s2A);
+      fX2Y2[8]->Fill(spdPercentile,c2A*s2A);
+      fX2Corr[8]->Fill(spdPercentile,qxA);
+      fY2Corr[8]->Fill(spdPercentile,qyA);
+      fX2Y2Corr[8]->Fill(spdPercentile,qxA*qyA);
+      Double_t psiA = TMath::ATan2(qyA,qxA)/2.;
+      fPsiAFlatCentr->Fill(spdPercentile,psiA);
+      Double_t psiAOrg = fEvent->GetEventplane()->GetEventplane("V0A",fEvent,2);
+      fPsiARawCentr->Fill(spdPercentile,psiAOrg);
+      fCos8Psi[8]->Fill(spdPercentile, 2./4.*TMath::Cos(2.*4.*psiAOrg));
+    }
+    if (totMultC > 1e-6) {
+      fX2[9]->Fill(spdPercentile,c2C);
+      fY2[9]->Fill(spdPercentile,s2C);
+      fX2Y2[9]->Fill(spdPercentile,c2C*s2C);
+      fX2Corr[9]->Fill(spdPercentile,qxC);
+      fY2Corr[9]->Fill(spdPercentile,qyC);
+      fX2Y2Corr[9]->Fill(spdPercentile,qxC*qyC);
+      Double_t psiC = TMath::ATan2(qyC,qxC)/2.;
+      fPsiCFlatCentr->Fill(spdPercentile,psiC);
+      Double_t psiCOrg = fEvent->GetEventplane()->GetEventplane("V0C",fEvent,2);
+      fPsiCRawCentr->Fill(spdPercentile,psiCOrg);
+      fCos8Psi[9]->Fill(spdPercentile, 2./4.*TMath::Cos(2.*4.*psiCOrg));
+   }
+    if (totMult > 1e-6) {
+      fX2[10]->Fill(spdPercentile,c2A+c2C);
+      fY2[10]->Fill(spdPercentile,s2A+s2C);
+      fX2Y2[10]->Fill(spdPercentile,(c2A+c2C)*(s2A+s2C));
+      fX2Corr[10]->Fill(spdPercentile,qxA+qxC);
+      fY2Corr[10]->Fill(spdPercentile,qyA+qyC);
+      fX2Y2Corr[10]->Fill(spdPercentile,(qxA+qxC)*(qyA+qyC));
+      Double_t psiAC = TMath::ATan2(qyA+qyC,qxA+qxC)/2.;
+      fPsiACFlatCentr->Fill(spdPercentile,psiAC);
+      Double_t psiACOrg = fEvent->GetEventplane()->GetEventplane("V0",fEvent,2);
+      fPsiACRawCentr->Fill(spdPercentile,psiACOrg);
+      fCos8Psi[10]->Fill(spdPercentile, 2./4.*TMath::Cos(2.*4.*psiACOrg));
     }
-    Double_t qxA = qxTierce[4] + qxTierce[5] + qxTierce[6] + qxTierce[7];
-    Double_t qyA = qyTierce[4] + qyTierce[5] + qyTierce[6] + qyTierce[7];
-    Double_t qxC = qxTierce[0] + qxTierce[1] + qxTierce[2] + qxTierce[3];
-    Double_t qyC = qyTierce[0] + qyTierce[1] + qyTierce[2] + qyTierce[3];
-    Double_t psiA = TMath::ATan2(qyA,qxA)/2.;
-    Double_t psiC = TMath::ATan2(qyC,qxC)/2.;
-
-    fPsiAC->Fill(psiA,psiC);
-    Double_t psiAOrg = fEvent->GetEventplane()->GetEventplane("V0A",fEvent,2);
-    Double_t psiCOrg = fEvent->GetEventplane()->GetEventplane("V0C",fEvent,2);
-    fPsiACOrg->Fill(psiAOrg,psiCOrg);
   }
 
   PostData(1, fOutputList);
@@ -310,12 +387,6 @@ Double_t AliAnaVZEROEPFlatenning::CalculateVZEROEventPlane(const AliVEvent *  ev
   Double_t lambdaPlus = b/aPlus;
   Double_t lambdaMinus = b/aMinus;
 
-  Double_t trans[2][2];
-  trans[0][0] = 1./(1.-lambdaPlus*lambdaMinus);
-  trans[0][1] = -lambdaMinus/(1.-lambdaPlus*lambdaMinus);
-  trans[1][0] = -lambdaPlus/(1.-lambdaPlus*lambdaMinus);
-  trans[1][1] = 1./(1.-lambdaPlus*lambdaMinus);
-
   Double_t qx=0., qy=0.;
   for(Int_t iCh = ring*8; iCh < (ring+1)*8; ++iCh) {
     Double_t phi = TMath::Pi()/8. + TMath::Pi()/4.*(iCh%8);
@@ -327,6 +398,11 @@ Double_t AliAnaVZEROEPFlatenning::CalculateVZEROEventPlane(const AliVEvent *  ev
   Double_t qxPrime = qx - meanX2;
   Double_t qyPrime = qy - meanY2;
   // Twist
+  Double_t trans[2][2];
+  trans[0][0] = 1./(1.-lambdaPlus*lambdaMinus);
+  trans[0][1] = -lambdaMinus/(1.-lambdaPlus*lambdaMinus);
+  trans[1][0] = -lambdaPlus/(1.-lambdaPlus*lambdaMinus);
+  trans[1][1] = 1./(1.-lambdaPlus*lambdaMinus);
   Double_t qxSeconde = trans[0][0]*qxPrime + trans[0][1]*qyPrime;
   Double_t qySeconde = trans[1][0]*qxPrime + trans[1][1]*qyPrime;
   // Rescaling
index 7c0f8ae5999c557799b2bed0adb82d643d2166f0..134c83d84fa80850560828a0a0074873f088262e 100644 (file)
@@ -30,14 +30,17 @@ class AliAnaVZEROEPFlatenning : public AliAnalysisTaskSE {
   TString fMBTrigName; // MB trigger name (for evt sel)
   Bool_t  fUsePhysSel; // Use or not phys sel
 
-  TProfile *fX2[8]; //! Profile histogram for Q^2_x
-  TProfile *fY2[8]; //! Profile histogram for Q^2_y
-  TProfile *fX2Y2[8]; //! Profile histogram for Q^2_x*Q^2_y
-  TProfile *fCos8Psi[8]; //! Profile histogram for Cos(8*Psi)
+  TProfile *fX2[11]; //! Profile histogram for Q^2_x
+  TProfile *fY2[11]; //! Profile histogram for Q^2_y
+  TProfile *fX2Y2[11]; //! Profile histogram for Q^2_x*Q^2_y
+  TProfile *fCos8Psi[11]; //! Profile histogram for Cos(8*Psi)
   TProfile *fC2[8]; //! Profile histogram for Cos(2*phi)
   TProfile *fS2[8]; //! Profile histogram for Sin(2*phi)
   TProfile *fC4[8]; //! Profile histogram for Cos(4*phi)
   TProfile *fS4[8]; //! Profile histogram for Sin(4*phi)
+  TProfile *fX2Corr[11]; //! Profile histogram for Q^2_x
+  TProfile *fY2Corr[11]; //! Profile histogram for Q^2_y
+  TProfile *fX2Y2Corr[11]; //! Profile histogram for Q^2_x*Q^2_y
 
   TProfile *fX2In[8]; // Profile histogram for Q^2_x (read from input file)
   TProfile *fY2In[8]; // Profile histogram for Q^2_y (read from input file)
@@ -47,9 +50,12 @@ class AliAnaVZEROEPFlatenning : public AliAnalysisTaskSE {
   TH2F *fPsiRingRawCentr[8]; //! Raw VZERO event plane vs centrality (ring-by-ring)
   TH2F *fPsiRingFlatCentr[8]; //! Flatenned with corrections on cumulants VZERO event plane vs centrality (ring-by-ring)
   TH2F *fPsiRingFlatFinalCentr[8]; //! Flatenned with corrections on cumulants and fourier VZERO event plane vs centrality (ring-by-ring)
-
-  TH2F *fPsiAC; //! Correlation between V0A and V0C event-plane angles
-  TH2F *fPsiACOrg; //! Correlation between V0A and V0C event-plane angles with all weight equal 1
+  TH2F *fPsiARawCentr; //! Raw VZEROA event plane vs centrality
+  TH2F *fPsiAFlatCentr; //! Flatenned with corrections on cumulants VZEROA event plane vs centrality 
+  TH2F *fPsiCRawCentr; //! Raw VZEROC event plane vs centrality
+  TH2F *fPsiCFlatCentr; //! Flatenned with corrections on cumulants VZEROC event plane vs centrality 
+  TH2F *fPsiACRawCentr; //! Raw VZERO event plane vs centrality
+  TH2F *fPsiACFlatCentr; //! Flatenned with corrections on cumulants VZERO event plane vs centrality 
 
   AliAnaVZEROEPFlatenning(const AliAnaVZEROEPFlatenning&); // not implemented
   AliAnaVZEROEPFlatenning& operator=(const AliAnaVZEROEPFlatenning&); // not implemented