]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
added option to create flow event from SPD tracklets
authorsnelling <snelling@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 15 Jun 2010 19:20:01 +0000 (19:20 +0000)
committersnelling <snelling@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 15 Jun 2010 19:20:01 +0000 (19:20 +0000)
PWG2/FLOW/AliFlowTasks/AliAnalysisTaskFlowEvent.cxx
PWG2/FLOW/AliFlowTasks/AliAnalysisTaskFlowEvent.h
PWG2/FLOW/AliFlowTasks/AliFlowEvent.cxx
PWG2/FLOW/AliFlowTasks/AliFlowEvent.h
PWG2/FLOW/macros/AddTaskFlow.C

index d291a63d9db7694c14378c6fcbcddc6e1d4e79bc..6253eb1830042e301a7361885829bd9ce07c52f4 100644 (file)
@@ -72,6 +72,7 @@ AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent() :
   AliAnalysisTaskSE(),
   //  fOutputFile(NULL),
   fAnalysisType("ESD"),
+  fRPType("Global"),
   fCFManager1(NULL),
   fCFManager2(NULL),
   fQAInt(NULL),
@@ -102,6 +103,7 @@ AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent(const char *name, Bool_t on,
   AliAnalysisTaskSE(name),
   //  fOutputFile(NULL),
   fAnalysisType("ESD"),
+  fRPType("Global"),
   fCFManager1(NULL),
   fCFManager2(NULL),
   fQAInt(NULL),
@@ -176,6 +178,9 @@ void AliAnalysisTaskFlowEvent::UserExec(Option_t *)
   AliMCEvent*  mcEvent = MCEvent();                              // from TaskSE
   AliESDEvent* myESD = dynamic_cast<AliESDEvent*>(InputEvent()); // from TaskSE
   AliAODEvent* myAOD = dynamic_cast<AliAODEvent*>(InputEvent()); // from TaskSE
+  AliMultiplicity* myTracklets = NULL;
+
+  /*test*/ cout<<"(AliAnalysisTaskFlowEvent::UserExec)fRPType = "<<fRPType<<endl;
 
   // Make the FlowEvent for MC input
   if (fAnalysisType == "MC")
@@ -206,6 +211,7 @@ void AliAnalysisTaskFlowEvent::UserExec(Option_t *)
     //make event
     flowEvent = new AliFlowEvent(mcEvent,fCFManager1,fCFManager2);
   }
+
   // Make the FlowEvent for ESD input
   else if (fAnalysisType == "ESD")
   {
@@ -230,8 +236,16 @@ void AliAnalysisTaskFlowEvent::UserExec(Option_t *)
     }
 
     //make event
-    flowEvent = new AliFlowEvent(myESD,fCFManager1,fCFManager2);
+    if (fRPType == "Global") {
+      flowEvent = new AliFlowEvent(myESD,fCFManager1,fCFManager2);
+    }
+    else if (fRPType == "Tracklet"){
+      flowEvent = new AliFlowEvent(myESD,myTracklets,fCFManager2);
+    }
+    // if monte carlo event get reaction plane from monte carlo (depends on generator)
+    if (mcEvent && mcEvent->GenEventHeader()) flowEvent->SetMCReactionPlaneAngle(mcEvent);
   }
+
   // Make the FlowEvent for ESD input combined with MC info
   else if (fAnalysisType == "ESDMCkineESD" || fAnalysisType == "ESDMCkineMC" )
   {
@@ -292,40 +306,11 @@ void AliAnalysisTaskFlowEvent::UserExec(Option_t *)
     AliWarning("FlowEvent cut on multiplicity"); return;
   }
 
-  //tag subevents
+  //tag subEvents
   flowEvent->TagSubeventsInEta(fMinA,fMaxA,fMinB,fMaxB);
 
-  ////TODO afterburner
-  //if (fbAfterburnerOn && fMyTRandom3) {
-  //  // set the new value of the values using a after burner
-  //  cout << "settings for afterburner in TaskFlowEvent.cxx:" << endl;
-  //  cout << "fCount = " << fCount << endl;
-  //  cout << "fNoOfLoops = " << fNoOfLoops << endl;
-  //  cout << "fEllipticFlowValue = " << fEllipticFlowValue << endl;
-  //  cout << "fSigmaEllipticFlowValue = " << fSigmaEllipticFlowValue << endl;
-  //  cout << "fMultiplicityOflowEvent = " << fMultiplicityOflowEvent << endl;
-  //  cout << "fSigmaMultiplicityOflowEvent = " << fSigmaMultiplicityOflowEvent << endl;
-  //  Double_t xRPangle;
-  //  if (!mcEvent) { xRPangle = TMath::TwoPi()*(fMyTRandom3->Rndm()); }
-  //  else { xRPangle = fRP; }
-  //  Double_t xNewFlowValue = fMyTRandom3->Gaus(fEllipticFlowValue,fSigmaEllipticFlowValue);
-  //  Int_t nNewMultOflowEvent =  TMath::Nint(fMyTRandom3->Gaus(fMultiplicityOflowEvent,fSigmaMultiplicityOflowEvent));
-  //  cout << "xRPangle = " << xRPangle << endl;
-  //  cout << "xNewFlowValue = " << xNewFlowValue << endl;
-  //  cout << "nNewMultOflowEvent = " << nNewMultOflowEvent << endl;
-  //  cout << "settings for after burner" << endl;
-  //  flowEventMaker->SetMCReactionPlaneAngle(xRPangle);
-  //  flowEventMaker->SetNoOfLoops(fNoOfLoops);
-  //  flowEventMaker->SetEllipticFlowValue(xNewFlowValue);
-  //  flowEventMaker->SetMultiplicityOflowEvent(nNewMultOflowEvent);
-  //  //end settings afterburner
-  //}
-
-  // if monte carlo event get reaction plane from monte carlo (depends on generator)
-  if (mcEvent && mcEvent->GenEventHeader()) flowEvent->SetMCReactionPlaneAngle(mcEvent);
-
   //fListHistos->Print();
-  //  fOutputFile->WriteObject(flowEvent,"myFlowEventSimple");
+  //fOutputFile->WriteObject(flowEvent,"myFlowEventSimple");
   PostData(1,flowEvent);
   if (fQA)
   {
index e02f3c36bb83b6b93bffeb25614d36e0921cb63d..cfccaba92f5dda6f369a23baefc7fc6ccc7d3045 100644 (file)
@@ -31,6 +31,9 @@ class AliAnalysisTaskFlowEvent : public AliAnalysisTaskSE {
   void SetAnalysisType(TString type) { this->fAnalysisType = type; }
   TString GetAnalysisType() const    { return this->fAnalysisType; }
 
+  void SetRPType(TString rptype) { this->fRPType = rptype; }
+  TString GetRPType() const    { return this->fRPType; }
+
   void    SetMinMult(Int_t multmin)    {this->fMinMult = multmin; }
   Int_t   GetMinMult() const           {return this->fMinMult; }
   void    SetMaxMult(Int_t multmax)    {this->fMaxMult = multmax; }
@@ -77,6 +80,7 @@ class AliAnalysisTaskFlowEvent : public AliAnalysisTaskSE {
   //  AliESDEvent*  fESD;                   // ESD object
   //  AliAODEvent*  fAOD;                   // AOD object
   TString       fAnalysisType;          // can be MC, ESD or AOD
+  TString       fRPType;                // can be Global or Tracklet
   AliCFManager* fCFManager1;            // correction framework manager
   AliCFManager* fCFManager2;            // correction framework manager
   TList*        fQAInt;                 // QA histogram list
index 0cc88cb498f0c76e2e9a89ee45a318b8a23dae36..883fd496e7af4d820bb3397b62a5b2b820902c04 100644 (file)
@@ -362,3 +362,73 @@ AliFlowEvent::AliFlowEvent( const AliESDEvent* anInput,
   SetMCReactionPlaneAngle(anInputMc);
 }
 
+//-----------------------------------------------------------------------
+AliFlowEvent::AliFlowEvent( const AliESDEvent* anInput,
+                                             const AliMultiplicity* anInputTracklets,
+                                             const AliCFManager* poiCFManager ):
+  AliFlowEventSimple(20)
+{
+
+  //Select the particles of interest from the ESD
+  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
+
+  //loop over tracks
+  for (Int_t itrkN=0; itrkN<iNumberOfInputTracks; itrkN++)
+    {
+      AliESDtrack* pParticle = anInput->GetTrack(itrkN);   //get input particle
+
+      //check if pParticle passes the cuts
+      Bool_t poiOK = kTRUE;
+      if (poiCFManager)
+       {
+         poiOK = ( poiCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle) &&
+                   poiCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,pParticle));
+       }
+      if (!poiOK) continue;
+      
+      //make new AliFLowTrack
+      AliFlowTrack* pTrack = new AliFlowTrack();
+      pTrack->SetPt(pParticle->Pt() );
+      pTrack->SetEta(pParticle->Eta() );
+      pTrack->SetPhi(pParticle->Phi() );
+          
+      //marking the particles used for the particle of interest (POI) selection:
+      if(poiOK && poiCFManager)
+       {
+         pTrack->SetForPOISelection(kTRUE);
+         pTrack->SetSource(AliFlowTrack::kFromESD);
+       }
+
+      AddTrack(pTrack);
+    }//end of while (itrkN < iNumberOfInputTracks)
+
+  //Select the reference particles from the SPD tracklets
+  anInputTracklets = anInput->GetMultiplicity();
+  Int_t multSPD = anInputTracklets->GetNumberOfTracklets();
+
+  cout << "N tracklets: " << multSPD << endl; //for testing
+
+  //loop over tracklets
+  for (Int_t itracklet=0; itracklet<multSPD; ++itracklet) {
+    Float_t thetaTr= anInputTracklets->GetTheta(itracklet);
+    Float_t phiTr= anInputTracklets->GetPhi(itracklet);
+    // calculate eta
+    Float_t etaTr = -TMath::Log(TMath::Tan(thetaTr/2.));
+    
+    //make new AliFLowTrackSimple
+    AliFlowTrack* pTrack = new AliFlowTrack();
+    pTrack->SetPt(0.0);
+    pTrack->SetEta(etaTr);
+    pTrack->SetPhi(phiTr);
+    //marking the particles used for the reference particle (RP) selection:
+    fNumberOfRPs++;
+    pTrack->SetForRPSelection(kTRUE);
+    pTrack->SetSource(AliFlowTrack::kFromTracklet);
+
+    //Add the track to the flowevent
+    AddTrack(pTrack);
+  }
+
+}
+
+
index 80076b512c10d6069e6e8b8ad767daf5c91ec103..6f72d3559a6919e4edbe628182c14b52896483ec 100644 (file)
@@ -49,7 +49,7 @@ public:
                 const AliCFManager* poiCFManager=NULL );  //use CF(2x)
   AliFlowEvent( const AliESDEvent* anInput,
                 const AliMultiplicity* anInputTracklets,
-                const AliCFManager* poiCFManager ){}
+                const AliCFManager* poiCFManager );
 
   void SetMCReactionPlaneAngle(const AliMCEvent* mcEvent);
 
index f7730514d5fba2184d2cd8f3ec6dc05abd42e613..4b9f8296657a8aef5a6fd85ffd3835cfbd716b5e 100644 (file)
@@ -8,10 +8,10 @@
 /////////////////////////////////////////////////////////////////////////////////////////////
 
 // Define the range for eta subevents (for SP method)
-Double_t minA = -0.9;
+Double_t minA = -2.0;
 Double_t maxA = -0.1;
 Double_t minB = 0.1;
-Double_t maxB = 0.9;
+Double_t maxB = 2.0;
 
 // use physics selection class
 Bool_t  UsePhysicsSelection = kTRUE;
@@ -43,6 +43,10 @@ const Int_t multmax = 10000;     //used for AliFlowEventSimple (to set the centr
 
 
 //----------For RP selection----------
+// Use Global tracks ("Global") or SPD tracklets ("Tracklet") for the RP selection
+const TString rptype = "Global";
+//const TString rptype = "Tracklet";
+
 //KINEMATICS (on generated and reconstructed tracks)
 Bool_t UseKineforRP =  kTRUE;
 const Double_t ptminRP = 0.0;
@@ -103,8 +107,8 @@ const Int_t  minTrackrefsMuonRP = 0;
 Bool_t UseKineforPOI = kTRUE;
 const Double_t ptminPOI = 0.0;
 const Double_t ptmaxPOI = 10.0;
-const Double_t etaminPOI  = -0.9;
-const Double_t etamaxPOI  = 0.9;
+const Double_t etaminPOI  = -0.5;
+const Double_t etamaxPOI  = 0.5;
 const Int_t    chargePOI = 1;  //not used
 const Bool_t   isChargedPOI = kTRUE;
 
@@ -309,6 +313,7 @@ AliAnalysisTaskFlowEvent* AddTaskFlow(TString type, Bool_t* METHODS, Bool_t QA,
       taskFE->SetEllipticFlowValue(ellipticFlow); }    //TEST
     else {taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",kTRUE); }
     taskFE->SetAnalysisType(type);
+    taskFE->SetRPType(rptype); //only for ESD
     if (UseMultCut) {
       taskFE->SetMinMult(multmin);
       taskFE->SetMaxMult(multmax);
@@ -693,15 +698,6 @@ AliAnalysisTaskFlowEvent* AddTaskFlow(TString type, Bool_t* METHODS, Bool_t QA,
     taskGFC->SetUsePhiWeights(WEIGHTS[0]); 
     taskGFC->SetUsePtWeights(WEIGHTS[1]);
     taskGFC->SetUseEtaWeights(WEIGHTS[2]); 
-    // calculation vs multiplicity:
-    taskGFC->SetCalculateVsMultiplicity(kFALSE);   
-    taskGFC->SetnBinsMult(10000);
-    taskGFC->SetMinMult(0);
-    taskGFC->SetMaxMult(10000);   
-    // tuning of interpolating parameters:
-    taskGFC->SetTuneParameters(kFALSE);
-    Double_t r0[10] = {1.8,1.9,2.0,2.1,2.2,2.3,2.4,2.5,2.6,2.7}; // up to 10 values allowed
-    for(Int_t r=0;r<10;r++) {taskGFC->SetTuningR0(r0[r],r);}
     mgr->AddTask(taskGFC);
   }
   if (QC){