]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Fixing conding violation rules. Add task macro that was forgotten last commit (Livio)
authormartinez <martinez@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 21 Apr 2010 16:22:44 +0000 (16:22 +0000)
committermartinez <martinez@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 21 Apr 2010 16:22:44 +0000 (16:22 +0000)
PWG3/muon/AddTaskTreeBuilder.C [new file with mode: 0644]
PWG3/muon/AliAnalysisTaskMuonTreeBuilder.cxx
PWG3/muon/AliAnalysisTaskMuonTreeBuilder.h

diff --git a/PWG3/muon/AddTaskTreeBuilder.C b/PWG3/muon/AddTaskTreeBuilder.C
new file mode 100644 (file)
index 0000000..56ad23c
--- /dev/null
@@ -0,0 +1,37 @@
+AliAnalysisTaskMuonTreeBuilder *AddTaskTreeBuilder(Bool_t ismc=kFALSE, Int_t run_num=0){
+  printf("Inside add task\n");
+  // Get the pointer to the existing analysis manager via the static access method
+  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+  if (!mgr) {
+    ::Error("AddTaskMuonTreeBuilder", "No analysis manager to connect to");
+    return NULL;
+  }   
+
+  // MC handler if needed
+  if(ismc){
+  AliMCEventHandler *mcH = (AliMCEventHandler*)mgr->GetMCtruthEventHandler();
+  if (!mcH) {
+    ::Error("AddTaskTreeBuilder", "No MC handler connected");
+    return NULL;
+  }    
+  }
+
+  // The task
+  AliAnalysisTaskMuonTreeBuilder *task = new AliAnalysisTaskMuonTreeBuilder("AliAnalysisTaskMuonTreeBuilder");
+  if(ismc) task->SetIsMC(kTRUE);
+
+  //outputs
+//   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0",TList::Class(),AliAnalysisManager::kOutputContainer,"final02.root");
+//   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("ctree0",TTree::Class(),AliAnalysisManager::kOutputContainer,"final02.root");
+  char outname[30];
+  sprintf(outname,"TreeRUN%d.root",run_num);
+  AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("ctree0",TTree::Class(),AliAnalysisManager::kOutputContainer,outname);
+
+  // Adding the task to the analysis manager
+  mgr->AddTask(task);
+  mgr->ConnectInput(task,0,mgr->GetCommonInputContainer());
+  mgr->ConnectOutput(task,1,coutput1);
+//   mgr->ConnectOutput(task,2,coutput2);
+
+  return task;
+}
index 80928a04db152c96223ec96299276c4591660d1f..91b1c387794121da93ffa03710d2c48eb846648f 100644 (file)
@@ -27,6 +27,7 @@
 //     Analysis task for muon-dimuon analysis
 //     Works for real and MC events
 //     Works with the corresponding AddTask macro
+//     Includes a flag for physics selection
 //     
 //     author: L. Bianchi - Universita' & INFN Torino
 
@@ -415,12 +416,12 @@ Double_t AliAnalysisTaskMuonTreeBuilder::Rap(Double_t e, Double_t pz) const
 {
 // calculate rapidity
     Double_t rap;
-    if(e!=pz){
+    if(e>TMath::Abs(pz)){
        rap = 0.5*TMath::Log((e+pz)/(e-pz));
        return rap;
     }
     else{
-       rap = -200;
+       rap = 666.;
        return rap;
     }
 }
@@ -440,47 +441,64 @@ Double_t AliAnalysisTaskMuonTreeBuilder::Phideg(Double_t phi) const
 Double_t AliAnalysisTaskMuonTreeBuilder::CostCS(Double_t px1, Double_t py1, Double_t pz1, Double_t e1,
 Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 {
-// calculate costCS
-  TLorentzVector PMu1CM, PMu2CM, PProjCM, PTargCM, PDimuCM; // In the CM. frame
-  TLorentzVector PMu1Dimu, PMu2Dimu, PProjDimu, PTargDimu; // In the dimuon rest frame
+// Cosine of the theta decay angle (mu+) in the Collins-Soper frame
+
+  TLorentzVector pMu1CM, pMu2CM, pProjCM, pTargCM, pDimuCM; // In the CM. frame
+  TLorentzVector pMu1Dimu, pMu2Dimu, pProjDimu, pTargDimu; // In the dimuon rest frame
   TVector3 beta,zaxisCS;
   Double_t mp=0.93827231;
   //
   // --- Fill the Lorentz vector for projectile and target in the CM frame
   //
-  PProjCM.SetPxPyPzE(0.,0.,-fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
-  PTargCM.SetPxPyPzE(0.,0.,fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
+  pProjCM.SetPxPyPzE(0.,0.,-fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
+  pTargCM.SetPxPyPzE(0.,0.,fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
   //
   // --- Get the muons parameters in the CM frame 
   //
-  PMu1CM.SetPxPyPzE(px1,py1,pz1,e1);
-  PMu2CM.SetPxPyPzE(px2,py2,pz2,e2);
+  pMu1CM.SetPxPyPzE(px1,py1,pz1,e1);
+  pMu2CM.SetPxPyPzE(px2,py2,pz2,e2);
   //
   // --- Obtain the dimuon parameters in the CM frame
   //
-  PDimuCM=PMu1CM+PMu2CM;
+  pDimuCM=pMu1CM+pMu2CM;
   //
   // --- Translate the dimuon parameters in the dimuon rest frame
   //
-  beta=(-1./PDimuCM.E())*PDimuCM.Vect();
-  PMu1Dimu=PMu1CM;
-  PMu2Dimu=PMu2CM;
-  PProjDimu=PProjCM;
-  PTargDimu=PTargCM;
-  PMu1Dimu.Boost(beta);
-  PMu2Dimu.Boost(beta);
-  PProjDimu.Boost(beta);
-  PTargDimu.Boost(beta);
-  //
+  beta=(-1./pDimuCM.E())*pDimuCM.Vect();
+  if(beta.Mag()>=1) return 666.;
+  pMu1Dimu=pMu1CM;
+  pMu2Dimu=pMu2CM;
+  pProjDimu=pProjCM;
+  pTargDimu=pTargCM;
+  pMu1Dimu.Boost(beta);
+  pMu2Dimu.Boost(beta);
+  pProjDimu.Boost(beta);
+  pTargDimu.Boost(beta);
+  
+  //Debugging part -------------------------------------
+  Double_t debugProj[4]={0.,0.,0.,0.};
+  Double_t debugTarg[4]={0.,0.,0.,0.};
+  Double_t debugMu1[4]={0.,0.,0.,0.};
+  Double_t debugMu2[4]={0.,0.,0.,0.};
+  pMu1Dimu.GetXYZT(debugMu1);
+  pMu2Dimu.GetXYZT(debugMu2);
+  pProjDimu.GetXYZT(debugProj);
+  pTargDimu.GetXYZT(debugTarg);
+  if (debugProj[0]!=debugProj[0] ||debugProj[1]!=debugProj[1] || debugProj[2]!=debugProj[2] ||debugProj[3]!=debugProj[3]) return 666; 
+  if (debugTarg[0]!=debugTarg[0] ||debugTarg[1]!=debugTarg[1] || debugTarg[2]!=debugTarg[2] ||debugTarg[3]!=debugTarg[3]) return 666; 
+  if (debugMu1[0]!=debugMu1[0] ||debugMu1[1]!=debugMu1[1] || debugMu1[2]!=debugMu1[2] ||debugMu1[3]!=debugMu1[3]) return 666; 
+  if (debugMu2[0]!=debugMu2[0] ||debugMu2[1]!=debugMu2[1] || debugMu2[2]!=debugMu2[2] ||debugMu2[3]!=debugMu2[3]) return 666; 
+  //----------------------------------------------------
+
   // --- Determine the z axis for the CS angle 
-  //
-  zaxisCS=(((PProjDimu.Vect()).Unit())-((PTargDimu.Vect()).Unit())).Unit();
-  //                                
+  zaxisCS=(((pProjDimu.Vect()).Unit())-((pTargDimu.Vect()).Unit())).Unit();
+                                    
   // --- Determine the CS angle (angle between mu+ and the z axis defined above)
   Double_t cost;
   
-  if(charge1>0) {cost = zaxisCS.Dot((PMu1Dimu.Vect()).Unit());}
-  else {cost = zaxisCS.Dot((PMu2Dimu.Vect()).Unit());}
+  if(charge1>0) {cost = zaxisCS.Dot((pMu1Dimu.Vect()).Unit());}
+  else {cost = zaxisCS.Dot((pMu2Dimu.Vect()).Unit());}
+  
   return cost;
 }
 
@@ -488,37 +506,47 @@ Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 Double_t AliAnalysisTaskMuonTreeBuilder::CostHE(Double_t px1, Double_t py1, Double_t pz1, Double_t e1,
 Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 {
-// calculate costHE
-  TLorentzVector PMu1CM, PMu2CM, PDimuCM; // In the CM frame 
-  TLorentzVector PMu1Dimu, PMu2Dimu; // In the dimuon rest frame
+// Cosine of the theta decay angle (mu+) in the Helicity frame
+  
+  TLorentzVector pMu1CM, pMu2CM, pDimuCM; // In the CM frame 
+  TLorentzVector pMu1Dimu, pMu2Dimu; // In the dimuon rest frame
   TVector3 beta,zaxisCS;
   //
   // --- Get the muons parameters in the CM frame
   //
-  PMu1CM.SetPxPyPzE(px1,py1,pz1,e1);
-  PMu2CM.SetPxPyPzE(px2,py2,pz2,e2);
+  pMu1CM.SetPxPyPzE(px1,py1,pz1,e1);
+  pMu2CM.SetPxPyPzE(px2,py2,pz2,e2);
   //
   // --- Obtain the dimuon parameters in the CM frame
   //
-  PDimuCM=PMu1CM+PMu2CM;
+  pDimuCM=pMu1CM+pMu2CM;
   //
   // --- Translate the muon parameters in the dimuon rest frame
   //
-  beta=(-1./PDimuCM.E())*PDimuCM.Vect();
-  PMu1Dimu=PMu1CM;
-  PMu2Dimu=PMu2CM;
-  PMu1Dimu.Boost(beta);
-  PMu2Dimu.Boost(beta);
-  //
+  beta=(-1./pDimuCM.E())*pDimuCM.Vect();
+  if(beta.Mag()>=1) return 666.;
+  pMu1Dimu=pMu1CM;
+  pMu2Dimu=pMu2CM;
+  pMu1Dimu.Boost(beta);
+  pMu2Dimu.Boost(beta);
+  
+  //Debugging part -------------------------------------
+  Double_t debugMu1[4]={0.,0.,0.,0.};
+  Double_t debugMu2[4]={0.,0.,0.,0.};
+  pMu1Dimu.GetXYZT(debugMu1);
+  pMu2Dimu.GetXYZT(debugMu2);
+  if (debugMu1[0]!=debugMu1[0] ||debugMu1[1]!=debugMu1[1] || debugMu1[2]!=debugMu1[2] ||debugMu1[3]!=debugMu1[3]) return 666; 
+  if (debugMu2[0]!=debugMu2[0] ||debugMu2[1]!=debugMu2[1] || debugMu2[2]!=debugMu2[2] ||debugMu2[3]!=debugMu2[3]) return 666; 
+  //----------------------------------------------------
   // --- Determine the z axis for the calculation of the polarization angle (i.e. the direction of the dimuon in the CM system)
-  //
   TVector3 zaxis;
-  zaxis=(PDimuCM.Vect()).Unit();
+  zaxis=(pDimuCM.Vect()).Unit();
   
   // --- Calculation of the polarization angle (angle between mu+ and the z axis defined above)
   Double_t cost;
-  if(charge1>0) {cost = zaxis.Dot((PMu1Dimu.Vect()).Unit());} 
-  else {cost = zaxis.Dot((PMu2Dimu.Vect()).Unit());} 
+  if(charge1>0) {cost = zaxis.Dot((pMu1Dimu.Vect()).Unit());} 
+  else {cost = zaxis.Dot((pMu2Dimu.Vect()).Unit());} 
   return cost;
 }
 
@@ -526,51 +554,64 @@ Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 Double_t AliAnalysisTaskMuonTreeBuilder::PhiCS(Double_t px1, Double_t py1, Double_t pz1, Double_t e1,
 Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 {
-// calculate phiCS
-   TLorentzVector PMu1CM, PMu2CM, PProjCM, PTargCM, PDimuCM; // In the CM frame
-   TLorentzVector PMu1Dimu, PMu2Dimu, PProjDimu, PTargDimu; // In the dimuon rest frame
+// Phi decay angle (mu+) in the Collins-Soper frame
+
+   TLorentzVector pMu1CM, pMu2CM, pProjCM, pTargCM, pDimuCM; // In the CM frame
+   TLorentzVector pMu1Dimu, pMu2Dimu, pProjDimu, pTargDimu; // In the dimuon rest frame
    TVector3 beta,yaxisCS, xaxisCS, zaxisCS;
    Double_t mp=0.93827231;
-   //
+   
    // --- Fill the Lorentz vector for projectile and target in the CM frame
-   //
-   PProjCM.SetPxPyPzE(0.,0.,-fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
-   PTargCM.SetPxPyPzE(0.,0.,fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
-   //
+   pProjCM.SetPxPyPzE(0.,0.,-fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
+   pTargCM.SetPxPyPzE(0.,0.,fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); 
+   
    // --- Get the muons parameters in the CM frame 
-   //
-   PMu1CM.SetPxPyPzE(px1,py1,pz1,e1);
-   PMu2CM.SetPxPyPzE(px2,py2,pz2,e2);
-   //
+   pMu1CM.SetPxPyPzE(px1,py1,pz1,e1);
+   pMu2CM.SetPxPyPzE(px2,py2,pz2,e2);
+   
    // --- Obtain the dimuon parameters in the CM frame
-   //
-   PDimuCM=PMu1CM+PMu2CM;
-   //
+   pDimuCM=pMu1CM+pMu2CM;
+   
    // --- Translate the dimuon parameters in the dimuon rest frame
-   //
-   beta=(-1./PDimuCM.E())*PDimuCM.Vect();
-   PMu1Dimu=PMu1CM;
-   PMu2Dimu=PMu2CM;
-   PProjDimu=PProjCM;
-   PTargDimu=PTargCM;
-   PMu1Dimu.Boost(beta);
-   PMu2Dimu.Boost(beta);
-   PProjDimu.Boost(beta);
-   PTargDimu.Boost(beta);
-   //
+   beta=(-1./pDimuCM.E())*pDimuCM.Vect();
+   if(beta.Mag()>=1) return 666.;
+   pMu1Dimu=pMu1CM;
+   pMu2Dimu=pMu2CM;
+   pProjDimu=pProjCM;
+   pTargDimu=pTargCM;
+   pMu1Dimu.Boost(beta);
+   pMu2Dimu.Boost(beta);
+   pProjDimu.Boost(beta);
+   pTargDimu.Boost(beta);
+
+   //Debugging part -------------------------------------
+   Double_t debugProj[4]={0.,0.,0.,0.};
+   Double_t debugTarg[4]={0.,0.,0.,0.};
+   Double_t debugMu1[4]={0.,0.,0.,0.};
+   Double_t debugMu2[4]={0.,0.,0.,0.};
+   pMu1Dimu.GetXYZT(debugMu1);
+   pMu2Dimu.GetXYZT(debugMu2);
+   pProjDimu.GetXYZT(debugProj);
+   pTargDimu.GetXYZT(debugTarg);
+   if (debugProj[0]!=debugProj[0] ||debugProj[1]!=debugProj[1] || debugProj[2]!=debugProj[2] ||debugProj[3]!=debugProj[3]) return 666; 
+   if (debugTarg[0]!=debugTarg[0] ||debugTarg[1]!=debugTarg[1] || debugTarg[2]!=debugTarg[2] ||debugTarg[3]!=debugTarg[3]) return 666; 
+   if (debugMu1[0]!=debugMu1[0] ||debugMu1[1]!=debugMu1[1] || debugMu1[2]!=debugMu1[2] ||debugMu1[3]!=debugMu1[3]) return 666; 
+   if (debugMu2[0]!=debugMu2[0] ||debugMu2[1]!=debugMu2[1] || debugMu2[2]!=debugMu2[2] ||debugMu2[3]!=debugMu2[3]) return 666; 
+   //----------------------------------------------------
+
    // --- Determine the z axis for the CS angle 
-   //
-   zaxisCS=(((PProjDimu.Vect()).Unit())-((PTargDimu.Vect()).Unit())).Unit();
-   yaxisCS=(((PProjDimu.Vect()).Unit()).Cross((PTargDimu.Vect()).Unit())).Unit();
+   zaxisCS=(((pProjDimu.Vect()).Unit())-((pTargDimu.Vect()).Unit())).Unit();
+   yaxisCS=(((pProjDimu.Vect()).Unit()).Cross((pTargDimu.Vect()).Unit())).Unit();
    xaxisCS=(yaxisCS.Cross(zaxisCS)).Unit();
  
    Double_t phi=0.;
    if(charge1>0) {
-       phi = TMath::ATan2((PMu1Dimu.Vect()).Dot(yaxisCS),((PMu1Dimu.Vect()).Dot(xaxisCS)));
+       phi = TMath::ATan2((pMu1Dimu.Vect()).Dot(yaxisCS),((pMu1Dimu.Vect()).Dot(xaxisCS)));
    } else {
-       phi = TMath::ATan2((PMu2Dimu.Vect()).Dot(yaxisCS),((PMu2Dimu.Vect()).Dot(xaxisCS)));
+       phi = TMath::ATan2((pMu2Dimu.Vect()).Dot(yaxisCS),((pMu2Dimu.Vect()).Dot(xaxisCS)));
    }
    if (phi>TMath::Pi()) phi=phi-TMath::Pi();
+   
    return phi;
 }
 
@@ -578,50 +619,61 @@ Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 Double_t AliAnalysisTaskMuonTreeBuilder::PhiHE(Double_t px1, Double_t py1, Double_t pz1, Double_t e1,
 Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2)
 {
-// calculate phiHE
-  TLorentzVector PMu1Lab, PMu2Lab, PProjLab, PTargLab, PDimuLab; // In the lab. frame 
-  TLorentzVector PMu1Dimu, PMu2Dimu, PProjDimu, PTargDimu; // In the dimuon rest frame
+// Phi decay angle (mu+) in the Helicity frame
+  TLorentzVector pMu1Lab, pMu2Lab, pProjLab, pTargLab, pDimuLab; // In the lab. frame 
+  TLorentzVector pMu1Dimu, pMu2Dimu, pProjDimu, pTargDimu; // In the dimuon rest frame
   TVector3 beta,xaxis, yaxis,zaxis;
   Double_t mp=0.93827231;
 
-  //
   // --- Get the muons parameters in the LAB frame
-  //
-  PMu1Lab.SetPxPyPzE(px1,py1,pz1,e1);
-  PMu2Lab.SetPxPyPzE(px2,py2,pz2,e2);
-  //
+  pMu1Lab.SetPxPyPzE(px1,py1,pz1,e1);
+  pMu2Lab.SetPxPyPzE(px2,py2,pz2,e2);
+  
   // --- Obtain the dimuon parameters in the LAB frame
-  //
-  PDimuLab=PMu1Lab+PMu2Lab;
+  pDimuLab=pMu1Lab+pMu2Lab;
+  zaxis=(pDimuLab.Vect()).Unit();
   
-  zaxis=(PDimuLab.Vect()).Unit();
-  //
   // --- Translate the muon parameters in the dimuon rest frame
-  //
-  beta=(-1./PDimuLab.E())*PDimuLab.Vect();
+  beta=(-1./pDimuLab.E())*pDimuLab.Vect();
+  if(beta.Mag()>=1.) return 666.;
 
-  PProjLab.SetPxPyPzE(0.,0.,-fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); // proiettile
-  PTargLab.SetPxPyPzE(0.,0.,fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); // bersaglio
+  pProjLab.SetPxPyPzE(0.,0.,-fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); // proiettile
+  pTargLab.SetPxPyPzE(0.,0.,fBeamEnergy,TMath::Sqrt(fBeamEnergy*fBeamEnergy+mp*mp)); // bersaglio
 
-  PProjDimu=PProjLab;
-  PTargDimu=PTargLab;
+  pProjDimu=pProjLab;
+  pTargDimu=pTargLab;
 
-  PProjDimu.Boost(beta);
-  PTargDimu.Boost(beta);
+  pProjDimu.Boost(beta);
+  pTargDimu.Boost(beta);
   
-  yaxis=((PProjDimu.Vect()).Cross(PTargDimu.Vect())).Unit();
+  yaxis=((pProjDimu.Vect()).Cross(pTargDimu.Vect())).Unit();
   xaxis=(yaxis.Cross(zaxis)).Unit();
   
-  PMu1Dimu=PMu1Lab;
-  PMu2Dimu=PMu2Lab;
-  PMu1Dimu.Boost(beta);
-  PMu2Dimu.Boost(beta);
+  pMu1Dimu=pMu1Lab;
+  pMu2Dimu=pMu2Lab;
+  pMu1Dimu.Boost(beta);
+  pMu2Dimu.Boost(beta);
+  
+  //Debugging part -------------------------------------
+  Double_t debugProj[4]={0.,0.,0.,0.};
+  Double_t debugTarg[4]={0.,0.,0.,0.};
+  Double_t debugMu1[4]={0.,0.,0.,0.};
+  Double_t debugMu2[4]={0.,0.,0.,0.};
+  pMu1Dimu.GetXYZT(debugMu1);
+  pMu2Dimu.GetXYZT(debugMu2);
+  pProjDimu.GetXYZT(debugProj);
+  pTargDimu.GetXYZT(debugTarg);
+  if (debugProj[0]!=debugProj[0] ||debugProj[1]!=debugProj[1] || debugProj[2]!=debugProj[2] ||debugProj[3]!=debugProj[3]) return 666; 
+  if (debugTarg[0]!=debugTarg[0] ||debugTarg[1]!=debugTarg[1] || debugTarg[2]!=debugTarg[2] ||debugTarg[3]!=debugTarg[3]) return 666; 
+  if (debugMu1[0]!=debugMu1[0] ||debugMu1[1]!=debugMu1[1] || debugMu1[2]!=debugMu1[2] ||debugMu1[3]!=debugMu1[3]) return 666; 
+  if (debugMu2[0]!=debugMu2[0] ||debugMu2[1]!=debugMu2[1] || debugMu2[2]!=debugMu2[2] ||debugMu2[3]!=debugMu2[3]) return 666; 
+  //----------------------------------------------------
   
   Double_t phi=0.;
    if(charge1 > 0) {
-      phi = TMath::ATan2((PMu1Dimu.Vect()).Dot(yaxis),(PMu1Dimu.Vect()).Dot(xaxis));
+      phi = TMath::ATan2((pMu1Dimu.Vect()).Dot(yaxis),(pMu1Dimu.Vect()).Dot(xaxis));
      } else { 
-      phi = TMath::ATan2((PMu2Dimu.Vect()).Dot(yaxis),(PMu2Dimu.Vect()).Dot(xaxis));
+      phi = TMath::ATan2((pMu2Dimu.Vect()).Dot(yaxis),(pMu2Dimu.Vect()).Dot(xaxis));
    }  
    return phi;
 }
index a1548e0f9a7524aabdbbc9ae4e354edb8878105e..e27669195416721e1404c8d4771659c6bc98e826 100644 (file)
@@ -11,7 +11,6 @@
 class TH1I;
 class TParticle ;
 class TLorentzVector ;
-//class TVector3 ;
 class TFile ;
 class AliStack ;
 class AliESDtrack;
@@ -41,56 +40,56 @@ class AliAnalysisTaskMuonTreeBuilder : public AliAnalysisTaskSE {
   
  protected:
   
-  Double_t             fNevt            ;    // event counter
-  Double_t             fBeamEnergy      ;    // Energy of the beam (required for the CS angle)
-  TList                *fOutput         ;    // output
-  TTree                *fOutputTree     ;    //! tree output
+  Double_t             fNevt            ;      // event counter
+  Double_t             fBeamEnergy      ;      // Energy of the beam (required for the CS angle)
+  TList                *fOutput         ;      // output
+  TTree                *fOutputTree     ;      //! tree output
   
-  Bool_t       fIsMC;
+  Bool_t       fIsMC;                          // if MC truth has to be read
   
-  Bool_t       fIsSelected;
-  char         fTrigClass[100];
+  Bool_t       fIsSelected;                    // physics selection flag
+  char         fTrigClass[100];                // fired trigger classes
 
-  Int_t                fNumMuonTracks  ;               // variables for single mu
-  Int_t                fNumSPDTracklets;
-  Int_t                fNumContributors;
-  Double_t     fVertex[3];
-  Double_t     fpT[10];
-  Double_t     fE[10];
-  Double_t     fpx[10];
-  Double_t     fpy[10];
-  Double_t     fpz[10];
-  Double_t     fpxUncorr[10];
-  Double_t     fpyUncorr[10];
-  Double_t     fpzUncorr[10];
-  Double_t     fy[10];
-  Double_t     feta[10];
-  Double_t     fphi[10];
-  Int_t                fMatchTrig[10];
-  Double_t     fTrackChi2[10];
-  Double_t     fMatchTrigChi2[10];
-  Double_t     fDCA[10];
-  Short_t      fCharge[10];
-  Int_t                fMuFamily[10];
-  Double_t     fRAtAbsEnd[10];
+  Int_t                fNumMuonTracks  ;               // muon tracks in the event
+  Int_t                fNumSPDTracklets;               // spd tracklets
+  Int_t                fNumContributors;               // n contributors
+  Double_t     fVertex[3];                     // x,y,z vertex
+  Double_t     fpT[10];                        // single mu pT
+  Double_t     fE[10];                         // single mu E
+  Double_t     fpx[10];                        // single mu px
+  Double_t     fpy[10];                        // single mu py
+  Double_t     fpz[10];                        // single mu pz
+  Double_t     fpxUncorr[10];                  // single mu px uncorrected
+  Double_t     fpyUncorr[10];                  // single mu py uncorrected
+  Double_t     fpzUncorr[10];                  // single mu pz uncorrected
+  Double_t     fy[10];                         // single mu y
+  Double_t     feta[10];                       // single mu eta
+  Double_t     fphi[10];                       // single mu phi
+  Int_t                fMatchTrig[10];                 // single mu match trigger
+  Double_t     fTrackChi2[10];                 // single mu chi2 track
+  Double_t     fMatchTrigChi2[10];             // single mu chi2 of match trigger
+  Double_t     fDCA[10];                       // single mu DCA
+  Short_t      fCharge[10];                    // single mu charge
+  Int_t                fMuFamily[10];                  // single mu provenience
+  Double_t     fRAtAbsEnd[10];                 // single mu distance from beam center at end abs
 
-  Int_t                fNumDimuons;                    // variables for dimuons
-  Int_t                fDimuonConstituent[45][2];
-  Double_t     fpTdimuon[45];
-  Double_t     fpxdimuon[45];
-  Double_t     fpydimuon[45];
-  Double_t     fpzdimuon[45];
-  Double_t     fydimuon[45];
-  Double_t     fiMassdimuon[45];
-  Double_t     fcostCS[45];
-  Double_t     fcostHE[45];
-  Double_t     fphiCS[45];
-  Double_t     fphiHE[45];
+  Int_t                fNumDimuons;                    // dimuons in the event
+  Int_t                fDimuonConstituent[45][2];      // reference to single mus
+  Double_t     fpTdimuon[45];                  // dimuon pT
+  Double_t     fpxdimuon[45];                  // dimuon px
+  Double_t     fpydimuon[45];                  // dimuon py
+  Double_t     fpzdimuon[45];                  // dimuon pz
+  Double_t     fydimuon[45];                   // dimuon y
+  Double_t     fiMassdimuon[45];               // dimuon invariant mass
+  Double_t     fcostCS[45];                    // dimuon cos theta Collins-Soper
+  Double_t     fcostHE[45];                    // dimuon cos theta Helicity
+  Double_t     fphiCS[45];                     // dimuon phi Collins-Soper
+  Double_t     fphiHE[45];                     // dimuon phi Helicity
   
   //Int_t              fIsPrimary[10];
-  Int_t                fPDG[10];
-  Int_t                fPDGmother[10];
-  Int_t                fPDGdimu[45];
+  Int_t                fPDG[10];                       // PDG single mu
+  Int_t                fPDGmother[10];                 // PDG mother single mu
+  Int_t                fPDGdimu[45];                   // PDG dimu
 
 
 
@@ -100,14 +99,8 @@ class AliAnalysisTaskMuonTreeBuilder : public AliAnalysisTaskSE {
   Int_t          FindMuFamily(AliMCParticle* mcTrack, AliMCEvent* mcEvent) const;
   Double_t Imass  (Double_t e1, Double_t px1, Double_t py1, Double_t pz1, Double_t e2, Double_t px2, Double_t py2, Double_t pz2) const;
   Double_t Rap   (Double_t e, Double_t pz) const;
-//   Double_t Imass(Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t) const;
-//   Double_t Rap(Double_t, Double_t) const;
   Double_t Phideg(Double_t phi) const;
   
-//   Double_t CostCS (Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t);
-//   Double_t CostHE (Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t);
-//   Double_t PhiCS  (Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t);
-//   Double_t PhiHE  (Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t);
   Double_t CostCS (Double_t px1, Double_t py1, Double_t pz1, Double_t e1, Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2);
   Double_t CostHE (Double_t px1, Double_t py1, Double_t pz1, Double_t e1, Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2);
   Double_t PhiCS  (Double_t px1, Double_t py1, Double_t pz1, Double_t e1, Double_t charge1, Double_t px2, Double_t py2, Double_t pz2, Double_t e2);