]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWGJE/AliAnalysisTaskJetServices.cxx
Analysis code updated
[u/mrichter/AliRoot.git] / PWGJE / AliAnalysisTaskJetServices.cxx
index 737e2474cb2f4ad272bb1088ae82c62c70c8f5ad..1d00af66124c0a107d2cd7cb6ed136e7dbcf8091 100644 (file)
@@ -86,6 +86,7 @@ AliAnalysisTaskJetServices::AliAnalysisTaskJetServices():
   fFilterMask(0),
   fRPMethod(0),
   fCollisionType(kPbPb),
+  fNTrigger(0),
   fAvgTrials(1),
   fVtxXMean(0),
   fVtxYMean(0),
@@ -101,10 +102,13 @@ AliAnalysisTaskJetServices::AliAnalysisTaskJetServices():
   fRPAngle(0),
   fPsiVZEROA(0),
   fPsiVZEROC(0),
+  fTriggerBit(0x0),
   fRandomizer(0),
   fNonStdFile(""),
+  fTriggerName(0x0),
   fh1Xsec(0x0),
   fh1Trials(0x0),
+  fh1AvgTrials(0x0),
   fh1PtHard(0x0),
   fh1PtHardTrials(0x0),
   fh1SelectionInfoESD(0x0),
@@ -112,6 +116,9 @@ AliAnalysisTaskJetServices::AliAnalysisTaskJetServices():
   fh1CentralityESD(0),
   fh1Centrality(0),
   fh1RP(0),
+  fh2ReducedTrigger(0),
+  fh2CentralityTriggerESD(0),
+  fh2CentralityTrigger(0),
   fh2TriggerCount(0x0),
   fh2ESDTriggerCount(0x0),
   fh2TriggerVtx(0x0),
@@ -153,6 +160,7 @@ AliAnalysisTaskJetServices::AliAnalysisTaskJetServices(const char* name):
   fFilterMask(0),
   fRPMethod(0),
   fCollisionType(kPbPb),
+  fNTrigger(0),
   fAvgTrials(1),
   fVtxXMean(0),
   fVtxYMean(0),
@@ -168,10 +176,13 @@ AliAnalysisTaskJetServices::AliAnalysisTaskJetServices(const char* name):
   fRPAngle(0),
   fPsiVZEROA(0),
   fPsiVZEROC(0),
+  fTriggerBit(0x0),
   fRandomizer(0),
   fNonStdFile(""),
+  fTriggerName(0x0),
   fh1Xsec(0x0),
   fh1Trials(0x0),
+  fh1AvgTrials(0x0),
   fh1PtHard(0x0),
   fh1PtHardTrials(0x0),
   fh1SelectionInfoESD(0x0),
@@ -179,6 +190,9 @@ AliAnalysisTaskJetServices::AliAnalysisTaskJetServices(const char* name):
   fh1CentralityESD(0),
   fh1Centrality(0),
   fh1RP(0),
+  fh2ReducedTrigger(0),
+  fh2CentralityTriggerESD(0),
+  fh2CentralityTrigger(0),
   fh2TriggerCount(0x0),
   fh2ESDTriggerCount(0x0),
   fh2TriggerVtx(0x0),
@@ -238,6 +252,7 @@ Bool_t AliAnalysisTaskJetServices::Notify()
     // construct a poor man average trials 
     Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
     if(ftrials>=nEntries && nEntries>0.)fAvgTrials = ftrials/nEntries;
+    fh1Trials->Fill("#sum{avg ntrials}",ftrials); 
   }  
   return kTRUE;
 }
@@ -267,14 +282,18 @@ void AliAnalysisTaskJetServices::UserCreateOutputObjects()
   fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
   fHistList->Add(fh1Trials);
 
-  const Int_t nBinPt = 100;
+  fh1AvgTrials = new TH1F("fh1AvgTrials","avg trials root file",1,0,1);
+  fh1AvgTrials->GetXaxis()->SetBinLabel(1,"#sum{avg ntrials}");
+  fHistList->Add(fh1AvgTrials);
+
+  const Int_t nBinPt = 125;
   Double_t binLimitsPt[nBinPt+1];
   for(Int_t iPt = 0;iPt <= nBinPt;iPt++){
     if(iPt == 0){
       binLimitsPt[iPt] = 0.0;
     }
     else {// 1.0
-      binLimitsPt[iPt] =  binLimitsPt[iPt-1] + 2.5;
+      binLimitsPt[iPt] =  binLimitsPt[iPt-1] + 2.;
     }
   }
   
@@ -288,6 +307,21 @@ void AliAnalysisTaskJetServices::UserCreateOutputObjects()
   fh1RP = new TH1F("fh1RP","RP;#Psi",440, -1.*TMath::Pi(), 2.*TMath::Pi());
   fHistList->Add(fh1RP);
 
+  fh2ReducedTrigger = new TH2F("fh2ReducedTrigger","red trigger;red trigger;centrality",1<<fNTrigger,-0.5,(1<<fNTrigger)-0.5,100,-0.5,99.5);
+  fHistList->Add(fh2ReducedTrigger);
+
+  fh2CentralityTriggerESD = new TH2F("fh2CentralityTriggerESD",";cent;trigger no",103,-1,102,fNTrigger,-0.5,fNTrigger-0.5);
+  fHistList->Add(fh2CentralityTriggerESD);
+  
+  fh2CentralityTrigger = new TH2F("fh2CentralityTrigger",";cent;trigger no",103,-1,102,fNTrigger,-0.5,fNTrigger-0.5);
+  fHistList->Add(fh2CentralityTrigger);
+
+  for(int i = 0;i<fNTrigger;++i){
+    fh2CentralityTriggerESD->GetYaxis()->SetBinLabel(i+1,fTriggerName[i].Data());
+    fh2CentralityTrigger->GetYaxis()->SetBinLabel(i+1,fTriggerName[i].Data());
+  }
+
+
   fh2TriggerCount = new TH2F("fh2TriggerCount",";Trigger No.;constrained;Count",6,-0.5,5.5,kConstraints,-0.5,kConstraints-0.5); 
   fHistList->Add(fh2TriggerCount);
 
@@ -456,25 +490,6 @@ void AliAnalysisTaskJetServices::UserExec(Option_t */*option*/)
   Bool_t aodVtxValid = false;
   Bool_t aodVtxIn = false;
 
-  if(esd){
-    // trigger analyisis
-    if(!fTriggerAnalysis){
-      fTriggerAnalysis = new AliTriggerAnalysis;
-      fTriggerAnalysis->SetAnalyzeMC(fMC);
-      fTriggerAnalysis->SetSPDGFOThreshhold(1);
-    }
-    //    fTriggerAnalysis->FillTriggerClasses(esd);
-    Bool_t v0A       = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0A);
-    Bool_t v0C       = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0C);
-    Bool_t v0ABG = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0ABG);
-    Bool_t v0CBG = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0CBG);
-    Bool_t spdFO      = fTriggerAnalysis->SPDFiredChips(esd, 0);;
-    if(v0A)fSelectionInfoESD |=  AliAnalysisHelperJetTasks::kV0A;
-    if(v0C)fSelectionInfoESD |=  AliAnalysisHelperJetTasks::kV0C;
-    if(!(v0ABG||v0CBG))fSelectionInfoESD |=  AliAnalysisHelperJetTasks::kNoV0BG;
-    if(spdFO)fSelectionInfoESD |=  AliAnalysisHelperJetTasks::kSPDFO;
-  }
   // Apply additional constraints
   Bool_t esdEventSelected = IsEventSelected(esd);
   Bool_t esdEventPileUp = IsEventPileUp(esd);
@@ -494,6 +509,8 @@ void AliAnalysisTaskJetServices::UserExec(Option_t */*option*/)
   if(physicsSelection) fSelectionInfoESD |=  AliAnalysisHelperJetTasks::kPhysicsSelection;
 
 
+
+
   // here we have all selection information, fill histogram
   for(unsigned int i = 1;i<(UInt_t)fh1SelectionInfoESD->GetNbinsX();i++){
     if((i&fSelectionInfoESD)==i)fh1SelectionInfoESD->Fill(i);
@@ -582,6 +599,14 @@ void AliAnalysisTaskJetServices::UserExec(Option_t */*option*/)
        tmpCent = esd->GetCentrality()->GetCentralityPercentile("V0M");
        if(tmpCent<0)tmpCent = 101;
        fh1CentralityESD->Fill(tmpCent);
+       UInt_t ir = 0;
+       for(int it = 0;it<fNTrigger;++it){
+         if(fInputHandler->IsEventSelected()&fTriggerBit[it]){
+           fh2CentralityTriggerESD->Fill(tmpCent,it);
+           ir |= (1<<it);
+         }
+       }
+       fh2ReducedTrigger->Fill(ir,cent);
       }
     }
   }
@@ -628,6 +653,9 @@ void AliAnalysisTaskJetServices::UserExec(Option_t */*option*/)
       fh2TriggerCount->Fill(iCl,kSelected);
       fh2TriggerCount->Fill(0.,kSelected);
       fh1Centrality->Fill(cent);
+      for(int it = 0;it<fNTrigger;++it){
+       if(fInputHandler->IsEventSelected()&fTriggerBit[it])fh2CentralityTrigger->Fill(cent,it);
+      }
       AliAnalysisHelperJetTasks::Selected(kTRUE,kTRUE);// select this event
       if(aodH&&cand&&fFilterAODCollisions&&!esd){
        if(fCentrality<=80&&aodVtxIn){
@@ -666,7 +694,7 @@ void AliAnalysisTaskJetServices::UserExec(Option_t */*option*/)
   Double_t ptHard = 0; 
   Double_t nTrials = 1; // Trials for MC trigger 
 
-  fh1Trials->Fill("#sum{ntrials}",fAvgTrials); 
+  fh1AvgTrials->Fill("#sum{avg ntrials}",fAvgTrials); 
   AliMCEvent* mcEvent = MCEvent();
   //    AliStack *pStack = 0; 
   if(mcEvent){
@@ -766,6 +794,8 @@ AliAnalysisTaskJetServices::~AliAnalysisTaskJetServices(){
   delete fp1CalibRPYA;
   delete fp1CalibRPXC;
   delete fp1CalibRPYC;
+  delete [] fTriggerBit;
+  delete [] fTriggerName;
 
 }
 
@@ -871,15 +901,16 @@ Bool_t  AliAnalysisTaskJetServices::IsVertexValid ( const AliAODVertex* vtx) con
     return kFALSE;
   }
 
-
+  TString vtxName(vtx->GetName());
   if(fDebug){
     Printf(" n contrib %d",vtx->GetNContributors());
+    Printf("vtxName: %s",vtxName.Data());
     vtx->Print();
   }
   
   //  if(vtx->GetNContributors()<3)return kFALSE;
   // do not want tpc only primary vertex
-  TString vtxName(vtx->GetName());
+
   if(vtxName.Contains("TPCVertex"))return kFALSE;
 
   // no dispersion yet...
@@ -1185,3 +1216,30 @@ Int_t  AliAnalysisTaskJetServices::GetListOfTracks(TList *list){
 
 }
 
+void AliAnalysisTaskJetServices::SetNTrigger(Int_t n){
+  if(n>0){
+      fNTrigger = n;
+      delete [] fTriggerName;
+      fTriggerName = new TString [fNTrigger];
+      delete [] fTriggerBit;fTriggerBit = 0;
+      fTriggerBit = new UInt_t [fNTrigger];
+  }
+  else{
+    fNTrigger = 0;
+  }
+}
+
+void AliAnalysisTaskJetServices::SetTrigger(Int_t i,UInt_t it,const char* c){
+  if(i<fNTrigger){
+    Printf("%d",it);
+    Printf("%p",c);
+    Printf("%s",c);    
+    Printf("%p",&fTriggerName[i]);
+
+    fTriggerBit[i] = it;
+    fTriggerName[i] =  c; // placement new
+    //    new(&fTriggerName[i]) TString(c); // placement new
+    Printf("%s",fTriggerName[i].Data());
+    
+  } 
+}