]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
from Panos:
authormkrzewic <mkrzewic@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 3 Jun 2013 09:31:02 +0000 (09:31 +0000)
committermkrzewic <mkrzewic@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 3 Jun 2013 09:31:02 +0000 (09:31 +0000)
take care of applying event level cuts in case MC (i.e. generator level) is analysed:
i) external reference multiplicity
ii) impact parameter range from the MC header

PWG/FLOW/Base/AliFlowAnalysisWithQCumulants.cxx
PWG/FLOW/Tasks/AliAnalysisTaskFlowEvent.cxx
PWG/FLOW/Tasks/AliAnalysisTaskFlowStrange.cxx
PWG/FLOW/Tasks/AliAnalysisTaskPIDflowQA.cxx
PWG/FLOW/Tasks/AliAnalysisTaskQAflow.cxx
PWG/FLOW/Tasks/AliFlowEventCuts.cxx
PWG/FLOW/Tasks/AliFlowEventCuts.h

index be74f453f953064c7801b3c525c83fe905104809..6e1db8b224328df859b8b1a0d8d6acbc08bf3c0b 100644 (file)
@@ -366,6 +366,7 @@ void AliFlowAnalysisWithQCumulants::Make(AliFlowEventSimple* anEvent)
  if(fExactNoRPs > 0 && fNumberOfRPsEBE<fExactNoRPs){return;}
  fNumberOfPOIsEBE = anEvent->GetNumberOfPOIs(); // number of POIs (i.e. number of particles of interest)
  fReferenceMultiplicityEBE = anEvent->GetReferenceMultiplicity(); // reference multiplicity for current event
+ //Printf("Reference multiplicity (QC): %.1f",fReferenceMultiplicityEBE);
  Double_t ptEta[2] = {0.,0.}; // 0 = dPt, 1 = dEta
   
  // c) Fill the common control histograms and call the method to fill fAvMultiplicity:
index 3643f8ebb20f18ef960f49747f4594d3bca04a24..e678dcc97c2a2c74c40650b108f5b8b02204e399 100644 (file)
@@ -323,7 +323,8 @@ void AliAnalysisTaskFlowEvent::UserExec(Option_t *)
   if (fAnalysisType == "AUTOMATIC")
   {
     //check event cuts
-    if (InputEvent() && !fCutsEvent->IsSelected(InputEvent())) return;
+    if (InputEvent() && !fCutsEvent->IsSelected(InputEvent(),MCEvent())) 
+      return;
 
     //first attach all possible information to the cuts
     fCutsRP->SetEvent( InputEvent(), MCEvent() );  //attach event
index 190285466fbaa63f3dea77512ba079e12db3c99c..c10479ec53239d7378f426b1d9bea6ffee2ca559 100644 (file)
@@ -450,7 +450,7 @@ void AliAnalysisTaskFlowStrange::MyUserExec(Option_t *)
     Bool_t tESelection = TMath::Abs(tVtxZ-tSPDVtxZ) < 0.5;\r
     Bool_t tDSelection = kTRUE;\r
     if(fUseEventSelection) {\r
-      tDSelection = fCutsEvent->IsSelected(tAOD);\r
+      tDSelection = fCutsEvent->IsSelected(tAOD,0x0);\r
     } else {\r
       if(TMath::Abs(tVtxZ)>10.0) tDSelection = kFALSE; // Cut on VtxZ mandatory!\r
     }\r
index 35c6cfe8995eb90795f344f89b89cab8004fa446..b0724a1e4e3e4703c2d20f3c8e776c575a11d523 100644 (file)
@@ -553,7 +553,7 @@ void  AliAnalysisTaskPIDflowQA::UserExec(Option_t *)
     return;
   }
 
-  if (!(fEventCuts->IsSelected(fESD)))
+  if (!(fEventCuts->IsSelected(fESD,0x0)))
   {
     return;
   }
index 7f2e6f8e90e5b53d5b606bb0472f8e2ef779b058..781cedb37edd60d67d960cfb8e2cbb35b984f7fc 100644 (file)
@@ -222,7 +222,7 @@ void AliAnalysisTaskQAflow::UserExec(Option_t *)
   TH1* hstdspdtrmultB = static_cast<TH1*>(before->At(23));
   TH1* hstdspdtrmultA = static_cast<TH1*>(after->At(23));
 
-  Bool_t passevent = fEventCuts->IsSelected(event);
+  Bool_t passevent = fEventCuts->IsSelected(event,0x0);
   Bool_t isSelectedEventSelection = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kMB);
 
   AliMultiplicity* tracklets = const_cast<AliMultiplicity*>(event->GetMultiplicity());
index f03b7a853886abc569a8e5186e9d744a8c9e3cac..11eefa588ff4c8424446886ec964d9cc471af3ea 100644 (file)
@@ -40,6 +40,8 @@
 #include "AliFlowEventCuts.h"
 #include "AliFlowTrackCuts.h"
 #include "AliTriggerAnalysis.h"
+#include "AliCollisionGeometry.h"
+#include "AliGenEventHeader.h"
 
 ClassImp(AliFlowEventCuts)
 
@@ -84,7 +86,10 @@ AliFlowEventCuts::AliFlowEventCuts():
   fCentralityPercentileMax(100.),
   fCentralityPercentileMin(0.),
   fCutZDCtiming(kFALSE),
-  fTrigAna()
+  fTrigAna(),
+  fCutImpactParameter(kFALSE),
+  fImpactParameterMin(0.0),
+  fImpactParameterMax(100.0)
 {
   //constructor 
 }
@@ -130,7 +135,10 @@ AliFlowEventCuts::AliFlowEventCuts(const char* name, const char* title):
   fCentralityPercentileMax(100.),
   fCentralityPercentileMin(0.),
   fCutZDCtiming(kFALSE),
-  fTrigAna()
+  fTrigAna(),
+  fCutImpactParameter(kFALSE),
+  fImpactParameterMin(0.0),
+  fImpactParameterMax(100.0)
 {
   //constructor 
 }
@@ -176,7 +184,10 @@ AliFlowEventCuts::AliFlowEventCuts(const AliFlowEventCuts& that):
   fCentralityPercentileMax(that.fCentralityPercentileMax),
   fCentralityPercentileMin(that.fCentralityPercentileMin),
   fCutZDCtiming(that.fCutZDCtiming),
-  fTrigAna()
+  fTrigAna(),
+  fCutImpactParameter(that.fCutImpactParameter),
+  fImpactParameterMin(that.fImpactParameterMin),
+  fImpactParameterMax(that.fImpactParameterMax)
 {
   if (that.fQA) DefineHistograms();
   //copy constructor 
@@ -262,15 +273,16 @@ AliFlowEventCuts& AliFlowEventCuts::operator=(const AliFlowEventCuts& that)
 }
 
 //----------------------------------------------------------------------- 
-Bool_t AliFlowEventCuts::IsSelected(TObject* obj)
+Bool_t AliFlowEventCuts::IsSelected(TObject* obj, TObject* objmc)
 {
   //check cuts
   AliVEvent* vevent = dynamic_cast<AliVEvent*>(obj);
-  if (vevent) return PassesCuts(vevent);
+  AliMCEvent* mcevent = dynamic_cast<AliMCEvent*>(objmc);
+  if (vevent) return PassesCuts(vevent,mcevent);;
   return kFALSE;  //when passed wrong type of object
 }
 //----------------------------------------------------------------------- 
-Bool_t AliFlowEventCuts::PassesCuts(AliVEvent *event)
+Bool_t AliFlowEventCuts::PassesCuts(AliVEvent *event, AliMCEvent *mcevent)
 {
   ///check if event passes cuts
   const AliVVertex* pvtx=event->GetPrimaryVertex();
@@ -365,10 +377,10 @@ Bool_t AliFlowEventCuts::PassesCuts(AliVEvent *event)
   }
   if(fCutNumberOfTracks) {if ( event->GetNumberOfTracks() < fNumberOfTracksMin ||
                                event->GetNumberOfTracks() >= fNumberOfTracksMax ) pass=kFALSE;}
-  if(fCutRefMult&&esdevent)
+  if((fCutRefMult&&mcevent)||(fCutRefMult&&esdevent))
   {
     //reference multiplicity still to be defined
-    Double_t refMult = RefMult(event,0x0);
+    Double_t refMult = RefMult(event,mcevent);
     if (refMult < fRefMultMin || refMult >= fRefMultMax )
     {
       pass=kFALSE;
@@ -425,6 +437,19 @@ Bool_t AliFlowEventCuts::PassesCuts(AliVEvent *event)
     if (nselected) meanpt=meanpt/nselected;
     if (meanpt<fMeanPtMin || meanpt >= fMeanPtMax) pass=kFALSE;
   }
+
+  //impact parameter cut
+  if(fCutImpactParameter) {
+    Double_t gImpactParameter = 0.;
+    if(mcevent) {
+      AliCollisionGeometry* headerH = dynamic_cast<AliCollisionGeometry*>(dynamic_cast<AliMCEvent*>(mcevent)->GenEventHeader());
+      if(headerH)
+       gImpactParameter = headerH->ImpactParameter();
+    }
+    if ((gImpactParameter < fImpactParameterMin) || (gImpactParameter >= fImpactParameterMax ))
+      pass=kFALSE;
+  }
+
   if (fQA&&pass) 
   {
     QAafter(1)->Fill(multGlobal,multTPC);
index 2d76801f5566324ac664464ac41c2e7662f05f0f..6ed557e5218d7c7ce93550a1d52484b4be06852d 100644 (file)
@@ -32,9 +32,9 @@ class AliFlowEventCuts : public TNamed {
   AliFlowEventCuts& operator=(const AliFlowEventCuts& someCuts);
   virtual  ~AliFlowEventCuts();
   
-  virtual Bool_t IsSelected(TObject* obj);
+  virtual Bool_t IsSelected(TObject* obj, TObject *objmc);
 
-  Bool_t PassesCuts(AliVEvent* event);
+  Bool_t PassesCuts(AliVEvent* event, AliMCEvent *mcevent);
   
   static AliFlowEventCuts* StandardCuts();
   
@@ -44,6 +44,9 @@ class AliFlowEventCuts : public TNamed {
   void SetRefMultMax(Int_t value) {fRefMultMax=value;fCutRefMult=kTRUE;}
   void SetRefMultMin(Int_t value) {fRefMultMin=value;fCutRefMult=kTRUE;}
   void SetRefMultRange(Int_t min, Int_t max) {fRefMultMin=min;fRefMultMax=max;fCutRefMult=kTRUE;}
+  void SetImpactParameterMax(Double_t value) {fImpactParameterMax=value;fCutImpactParameter=kTRUE;}
+  void SetImpactParameterMin(Double_t value) {fImpactParameterMin=value;fCutImpactParameter=kTRUE;}
+  void SetImpactParameterRange(Double_t min, Double_t max) {fImpactParameterMin=min;fImpactParameterMax=max;fCutImpactParameter=kTRUE;}
   void SetPrimaryVertexXrange(Double_t min, Double_t max)
        {fCutPrimaryVertexX=kTRUE; fPrimaryVertexXmin=min; fPrimaryVertexXmax=max;}
   void SetPrimaryVertexYrange(Double_t min, Double_t max)
@@ -128,6 +131,9 @@ class AliFlowEventCuts : public TNamed {
   Float_t fCentralityPercentileMin; // min centr. perc
   Bool_t fCutZDCtiming;   //cut on ZDC timing
   AliTriggerAnalysis fTrigAna; //trigger analysis object
+  Bool_t fCutImpactParameter; //cut on impact parameter (MC header)
+  Double_t fImpactParameterMin; // min impact parameter
+  Double_t fImpactParameterMax; // max impact parameter
 
   ClassDef(AliFlowEventCuts,4)
 };