]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
-add cut qa class
authorjbook <jbook@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 25 Jun 2013 15:04:38 +0000 (15:04 +0000)
committerjbook <jbook@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 25 Jun 2013 15:04:38 +0000 (15:04 +0000)
-protection for missing V0 branch (nanoAODs)

PWGDQ/CMakelibPWGDQdielectron.pkg
PWGDQ/PWGDQdielectronLinkDef.h
PWGDQ/dielectron/AliAnalysisTaskMultiDielectron.cxx
PWGDQ/dielectron/AliDielectron.cxx
PWGDQ/dielectron/AliDielectron.h
PWGDQ/dielectron/AliDielectronCutQA.cxx [new file with mode: 0644]
PWGDQ/dielectron/AliDielectronCutQA.h [new file with mode: 0644]
PWGDQ/dielectron/AliDielectronV0Cuts.cxx

index 9c56c41cb5584aaad208e7838c2d8718130fe570..c2b8485da43bccc9c8b97f289800a23a5f62d5dd 100644 (file)
@@ -56,6 +56,7 @@ set ( SRCS
       dielectron/AliDielectronTrackRotator.cxx
       dielectron/AliDielectronPID.cxx
       dielectron/AliDielectronCutGroup.cxx
+      dielectron/AliDielectronCutQA.cxx
       dielectron/AliDielectronEventCuts.cxx
       dielectron/AliDielectronHelper.cxx
       dielectron/AliDielectronBtoJPSItoEleCDFfitFCN.cxx
index 2092187dd2e19f46ca1f633f6e481ff4b5554bc5..0c798ab8284bd30d8f7d47c520cb7f26298f3ccd 100644 (file)
@@ -36,6 +36,7 @@
 #pragma link C++ class AliDielectronTrackRotator+;
 #pragma link C++ class AliDielectronPID+;
 #pragma link C++ class AliDielectronCutGroup+;
+#pragma link C++ class AliDielectronCutQA+;
 #pragma link C++ class AliDielectronEventCuts+;
 #pragma link C++ class AliDielectronHelper+;
 #pragma link C++ class AliDielectronBtoJPSItoEleCDFfitFCN+;
index 1349032ace6237a73d8d1facaa614a7147720be3..78943c0d94943e09c2448721413be7aaf94983ef 100644 (file)
@@ -111,6 +111,9 @@ AliAnalysisTaskMultiDielectron::~AliAnalysisTaskMultiDielectron()
   fListHistos.SetOwner(kFALSE);
   fListCF.SetOwner(kFALSE);
   
+  // try to reduce memory issues
+  if(fEventStat)       { delete fEventStat;       fEventStat=0; }
+  if(fTriggerAnalysis) { delete fTriggerAnalysis; fTriggerAnalysis=0; }
 }
 //_________________________________________________________________________________
 void AliAnalysisTaskMultiDielectron::UserCreateOutputObjects()
@@ -129,9 +132,10 @@ void AliAnalysisTaskMultiDielectron::UserCreateOutputObjects()
   AliDielectron *die=0;
   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
     die->Init();
-    if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
-    if (die->GetHistogramArray()) fListHistos.Add(const_cast<TObjArray*>(die->GetHistogramArray()));
-    if (die->GetCFManagerPair()) fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
+    if (die->GetHistogramList())    fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
+    if (die->GetHistogramArray())   fListHistos.Add(const_cast<TObjArray*>(die->GetHistogramArray()));
+    if (die->GetQAHistArray())      fListHistos.Add(const_cast<TObjArray*>(die->GetQAHistArray()));
+    if (die->GetCFManagerPair())    fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
   }
 
   Int_t cuts=fListDielectron.GetEntries();
index 86946203d73873bfbab0f742c14fe8e068dd88c9..e3a778cdc30a9fc932359da8eab5ebe85ea9c517 100644 (file)
@@ -95,6 +95,8 @@ const char* AliDielectron::fgkPairClassNames[11] = {
 //________________________________________________________________
 AliDielectron::AliDielectron() :
   TNamed("AliDielectron","AliDielectron"),
+  fCutQA(kFALSE),
+  fQAmonitor(0x0),
   fEventFilter("EventFilter"),
   fTrackFilter("TrackFilter"),
   fPairPreFilter("PairPreFilter"),
@@ -136,6 +138,8 @@ AliDielectron::AliDielectron() :
 //________________________________________________________________
 AliDielectron::AliDielectron(const char* name, const char* title) :
   TNamed(name,title),
+  fCutQA(kFALSE),
+  fQAmonitor(0x0),
   fEventFilter("EventFilter"),
   fTrackFilter("TrackFilter"),
   fPairPreFilter("PairPreFilter"),
@@ -180,6 +184,7 @@ AliDielectron::~AliDielectron()
   //
   // Default destructor
   //
+  if (fQAmonitor) delete fQAmonitor;
   if (fHistoArray) delete fHistoArray;
   if (fHistos) delete fHistos;
   if (fPairCandidates) delete fPairCandidates;
@@ -219,7 +224,15 @@ void AliDielectron::Init()
     fHistoArray->SetSignalsMC(fSignalsMC);
     fHistoArray->Init();
   }
-} 
+
+  if (fCutQA) {
+    fQAmonitor = new AliDielectronCutQA(Form("QAcuts_%s",GetName()),"QAcuts");
+    fQAmonitor->AddTrackFilter(&fTrackFilter);
+    fQAmonitor->AddPairFilter(&fPairFilter);
+    fQAmonitor->AddEventFilter(&fEventFilter);
+    fQAmonitor->Init();
+  }
+}
 
 //________________________________________________________________
 void AliDielectron::Process(AliVEvent *ev1, AliVEvent *ev2)
@@ -251,6 +264,7 @@ void AliDielectron::Process(AliVEvent *ev1, AliVEvent *ev2)
   }
 
   //in case we have MC load the MC event and process the MC particles
+  // why do not apply the event cuts first ????
   if (AliDielectronMC::Instance()->ConnectMCEvent()){
     ProcessMC(ev1);
   }
@@ -266,8 +280,11 @@ void AliDielectron::Process(AliVEvent *ev1, AliVEvent *ev2)
   UInt_t selectedMask=(1<<fEventFilter.GetCuts()->GetEntries())-1;
 
   //apply event cuts
-    if ((ev1&&fEventFilter.IsSelected(ev1)!=selectedMask) ||
-        (ev2&&fEventFilter.IsSelected(ev2)!=selectedMask)) return;
+  UInt_t cutmask = fEventFilter.IsSelected(ev1);
+  if(fCutQA) fQAmonitor->FillAll(ev1);
+  if(fCutQA) fQAmonitor->Fill(cutmask,ev1);
+  if ((ev1&&cutmask!=selectedMask) ||
+      (ev2&&fEventFilter.IsSelected(ev2)!=selectedMask)) return;
 
   //fill track arrays for the first event
   if (ev1){
@@ -613,12 +630,19 @@ void AliDielectron::FillTrackArrays(AliVEvent * const ev, Int_t eventNr)
     AliVParticle *particle=ev->GetTrack(itrack);
 
     //apply track cuts
-    if (fTrackFilter.IsSelected(particle)!=selectedMask) continue;
+    UInt_t cutmask=fTrackFilter.IsSelected(particle);
+    //fill cut QA
+    if(fCutQA) fQAmonitor->FillAll(particle);
+    if(fCutQA) fQAmonitor->Fill(cutmask,particle);
+
+    if (cutmask!=selectedMask) continue;
 
     //fill selected particle into the corresponding track arrays
     Short_t charge=particle->Charge();
     if (charge>0)      fTracks[eventNr*2].Add(particle);
     else if (charge<0) fTracks[eventNr*2+1].Add(particle);
+
+  
   }
 }
 
@@ -1095,6 +1119,11 @@ void AliDielectron::FillPairArrays(Int_t arr1, Int_t arr2)
       if (fCfManagerPair) fCfManagerPair->Fill(cutMask,candidate);
       //histogram array for the pair
       if (fHistoArray) fHistoArray->Fill(pairIndex,candidate);
+      // cut qa
+      if(pairIndex==AliDielectron::kEv1PM && fCutQA) {
+       fQAmonitor->FillAll(candidate);
+       fQAmonitor->Fill(cutMask,candidate);
+      }
 
       //apply cut
       if (cutMask!=selectedMask) continue;
index d53c1888b836b30a3f3a02fcc8b3e882ff537a4d..116f1918075841b44fc9ef17a7a59661e9d68452 100644 (file)
@@ -24,6 +24,7 @@
 
 #include "AliDielectronHistos.h"
 #include "AliDielectronHF.h"
+#include "AliDielectronCutQA.h"
 
 class AliEventplane;
 class AliVEvent;
@@ -74,6 +75,7 @@ public:
   Int_t GetLeg1Pdg()   const { return fPdgLeg1;   }
   Int_t GetLeg2Pdg()   const { return fPdgLeg2;   }
 
+  void SetCutQA(Bool_t qa=kTRUE) { fCutQA=qa; }
   void SetNoPairing(Bool_t noPairing=kTRUE) { fNoPairing=noPairing; }
   void SetUseKF(Bool_t useKF=kTRUE) { fUseKF=useKF; }
   const TObjArray* GetTrackArray(Int_t i) const {return (i>=0&&i<4)?&fTracks[i]:0;}
@@ -83,6 +85,7 @@ public:
   TObjArray** GetPairArraysPointer() { return &fPairCandidates; }
   void SetHistogramArray(AliDielectronHF * const histoarray) { fHistoArray=histoarray; }
   const TObjArray * GetHistogramArray() const { return fHistoArray?fHistoArray->GetHistArray():0x0; }
+  const TObjArray * GetQAHistArray() const { return fQAmonitor?fQAmonitor->GetQAHistArray():0x0; }
 
   void SetHistogramManager(AliDielectronHistos * const histos) { fHistos=histos; }
   AliDielectronHistos* GetHistoManager() const { return fHistos; }
@@ -130,7 +133,9 @@ public:
   void SaveDebugTree();
 
 private:
-  
+
+  Bool_t fCutQA;                    // monitor cuts
+  AliDielectronCutQA *fQAmonitor;   // monitoring of cuts
   AliAnalysisFilter fEventFilter;    // Event cuts
   AliAnalysisFilter fTrackFilter;    // leg cuts
   AliAnalysisFilter fPairPreFilter;  // pair prefilter cuts
diff --git a/PWGDQ/dielectron/AliDielectronCutQA.cxx b/PWGDQ/dielectron/AliDielectronCutQA.cxx
new file mode 100644 (file)
index 0000000..5f780fa
--- /dev/null
@@ -0,0 +1,247 @@
+/*************************************************************************
+ * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+//////////////////////////////////////////////////////////////////////////
+//                           CutQA                                      //
+//                                                                      //
+/*
+   Allow to monitor how many tracks,pair,events pass the selection criterion 
+   in any of the cuts added to the corresponding filters. All you need to 
+   add to your config is the following:
+
+   dielectron->SetCutQA();
+
+
+*/
+//                                                                      //
+//////////////////////////////////////////////////////////////////////////
+
+#include "AliDielectronCutQA.h"
+
+#include <TList.h>
+#include <TCollection.h>
+
+#include "AliDielectronCutGroup.h"
+#include "AliAnalysisCuts.h"
+
+#include "AliVEvent.h"
+#include "AliVParticle.h"
+#include "AliVTrack.h"
+#include "AliDielectronPair.h"
+
+
+ClassImp(AliDielectronCutQA)
+
+
+AliDielectronCutQA::AliDielectronCutQA() :
+  TNamed(),
+  fQAHistArray(0x0)
+{
+  //
+  // Default constructor
+  //
+  for(Int_t itype=0; itype<kNtypes; itype++) {
+    fCutQA[itype]=0x0;
+    fNCuts[itype]=1;
+    for(Int_t i=0; i<20; i++) {
+      fCutNames[i][itype]="";
+    }
+  }
+  fTypeKeys[kTrack] = "Track";
+  fTypeKeys[kPair] = "Pair";
+  fTypeKeys[kEvent] = "Event";
+
+}
+
+//_____________________________________________________________________
+AliDielectronCutQA::AliDielectronCutQA(const char* name, const char* title) :
+  TNamed(name, title),
+  fQAHistArray(0x0)
+{
+  //
+  // Named Constructor
+  //
+  for(Int_t itype=0; itype<kNtypes; itype++) {
+    fCutQA[itype]=0x0;
+    fNCuts[itype]=1;
+    for(Int_t i=0; i<20; i++) {
+      fCutNames[i][itype]="";
+    }
+  }
+  fTypeKeys[kTrack] = "Track";
+  fTypeKeys[kPair]  = "Pair";
+  fTypeKeys[kEvent] = "Event";
+
+}
+
+//_____________________________________________________________________
+AliDielectronCutQA::~AliDielectronCutQA() 
+{
+  //
+  //Default Destructor
+  //
+  fQAHistArray.Delete();
+}
+
+//_____________________________________________________________________
+void AliDielectronCutQA::Init()
+{
+
+  fQAHistArray.SetName(Form("%s",GetName()));
+
+  // loop over all types
+  for(Int_t itype=0; itype<kNtypes; itype++) {
+    //    printf("\n type: %d\n",itype);
+    fCutNames[0][itype]="no cuts";
+
+    // create histogram based on added cuts
+    fCutQA[itype] = new TH1F(fTypeKeys[itype],
+                            Form("%sQA;cuts;# passed %ss",fTypeKeys[itype],fTypeKeys[itype]),
+                            fNCuts[itype],0,fNCuts[itype]);
+    // loop over all cuts
+    for(Int_t i=0; i<fNCuts[itype]; i++) {
+      fCutQA[itype]->GetXaxis()->SetBinLabel(i+1,fCutNames[i][itype]);
+      //      printf(" %s \n",fCutNames[i][itype]);
+    }
+    fQAHistArray.AddLast(fCutQA[itype]);
+  }
+
+}
+
+//_____________________________________________________________________
+void AliDielectronCutQA::AddTrackFilter(AliAnalysisFilter *trackFilter)
+{
+  //
+  // add track filter cuts to the qa histogram
+  //
+
+
+  TIter listIterator(trackFilter->GetCuts());
+  while (AliAnalysisCuts *thisCut = (AliAnalysisCuts*) listIterator()) {
+    Bool_t addCut=kTRUE;
+
+    // add new cut class to the array
+    if(addCut) {
+      fCutNames[fNCuts[kTrack]][kTrack]=thisCut->GetTitle();
+      //      printf("add cut %s to %d \n",thisCut->GetTitle(),fNCuts[kTrack]);
+      fNCuts[kTrack]++;
+    }
+
+  } // pair filter loop
+
+}
+
+
+//_____________________________________________________________________
+void AliDielectronCutQA::AddPairFilter(AliAnalysisFilter *pairFilter)
+{
+  //
+  // add track filter cuts to the qa histogram
+  //
+
+
+  TIter listIterator(pairFilter->GetCuts());
+  while (AliAnalysisCuts *thisCut = (AliAnalysisCuts*) listIterator()) {
+    Bool_t addCut=kTRUE;
+
+    // add new cut class to the array
+    if(addCut) {
+      fCutNames[fNCuts[kPair]][kPair]=thisCut->GetTitle();
+      //  printf("add cut %s to %d \n",thisCut->GetTitle(),fNCuts[kPair]);
+      fNCuts[kPair]++;
+    }
+
+  } // trk filter loop
+
+}
+
+//_____________________________________________________________________
+void AliDielectronCutQA::AddEventFilter(AliAnalysisFilter *eventFilter)
+{
+  //
+  // add track filter cuts to the qa histogram
+  //
+
+
+  TIter listIterator(eventFilter->GetCuts());
+  while (AliAnalysisCuts *thisCut = (AliAnalysisCuts*) listIterator()) {
+    Bool_t addCut=kTRUE;
+
+    // add new cut class to the array
+    if(addCut) {
+      fCutNames[fNCuts[kEvent]][kEvent]=thisCut->GetTitle();
+      //      printf("add cut %s to %d \n",thisCut->GetTitle(),fNCuts[kEvent]);
+      fNCuts[kEvent]++;
+    }
+
+  } // trk filter loop
+
+}
+
+//_____________________________________________________________________
+void AliDielectronCutQA::Fill(UInt_t mask, TObject *obj)
+{
+  //
+  // fill the corresponding step in the qa histogram
+  //
+
+  UInt_t idx = GetObjIndex(obj);
+
+  Int_t cutstep=1;
+  for (Int_t iCut=0; iCut<fNCuts[idx]-1;++iCut) {
+    //    UInt_t cutMask=1<<iCut;         // for each cut
+    UInt_t cutMask=(1<<(iCut+1))-1; // increasing cut match
+
+    if ((mask&cutMask)==cutMask) {
+      fCutQA[idx]->Fill(cutstep);
+      ++cutstep;
+    }
+
+  }
+
+}
+
+//_____________________________________________________________________
+void AliDielectronCutQA::FillAll(TObject *obj)
+{
+  //
+  // fill the corresponding step in the qa histogram
+  //
+
+  UInt_t idx = GetObjIndex(obj);
+  fCutQA[idx]->Fill(0);
+
+}
+
+//______________________________________________________________________
+UInt_t AliDielectronCutQA::GetObjIndex(TObject *obj)
+{
+  //
+  // return the corresponding idex
+  //
+
+  if(obj->InheritsFrom(AliVTrack::Class())        ) return kTrack;
+  if(obj->InheritsFrom(AliVParticle::Class())     ) return kTrack;
+  if(obj->InheritsFrom(AliDielectronPair::Class())) return kPair;
+  if(obj->InheritsFrom(AliVEvent::Class())        ) return kEvent;
+  printf("FATAL: object type %s not yet supported, please let the author know\n", obj->IsA()->GetName());
+  return -1;
+  //TODO complete
+
+}
+
+
+
+
diff --git a/PWGDQ/dielectron/AliDielectronCutQA.h b/PWGDQ/dielectron/AliDielectronCutQA.h
new file mode 100644 (file)
index 0000000..ff055c1
--- /dev/null
@@ -0,0 +1,65 @@
+#ifndef ALIDIELECTRONCUTQA_H
+#define ALIDIELECTRONCUTQA_H
+
+/* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//#################################################################
+//#                                                               #
+//#             Class AliDielectronCutQA                          #
+//#              Dielectron Group of cuts                         #
+//#                                                               #
+//#  Authors:                                                     #
+//#   Julian    Book,     Uni Ffm / Julian.Book@cern.ch           #
+//#                                                               #
+//#################################################################
+
+#include <TNamed.h>
+#include <AliAnalysisFilter.h>
+#include <TH1F.h>
+#include <TList.h>
+#include <TObjArray.h>
+
+class TCollection;
+
+class AliDielectronCutQA : public TNamed {
+  
+public:
+  enum { kTrack=0, kPair, kEvent, kNtypes };
+
+  AliDielectronCutQA();
+  AliDielectronCutQA(const char*name, const char* title);
+  
+  virtual ~AliDielectronCutQA();
+
+  //Analysis cuts interface
+  //
+  void Init();
+  void AddTrackFilter(     AliAnalysisFilter *trkFilter);
+  /* void AddPrePairFilter(   AliAnalysisFilter *prePairFilter); */
+  /* void AddPrePairLegFilter(AliAnalysisFilter *prePairLegFilter); */
+  void AddPairFilter(      AliAnalysisFilter *pairFilter);
+  void AddEventFilter(     AliAnalysisFilter *eventFilter);
+
+  //  void Fill(AliAnalysisCuts *cut);
+  void Fill(UInt_t mask, TObject *obj);
+  void FillAll(TObject *obj);// { fCutQA->Fill(0); }
+
+  const TObjArray * GetQAHistArray() const { return &fQAHistArray; }
+
+
+private:
+
+  TObjArray fQAHistArray;              //-> array of QA histograms
+  TH1F *fCutQA[kNtypes];               // qa histogram
+  Int_t fNCuts[kNtypes];               // number of cuts
+  const char* fCutNames[20][kNtypes];  // cut names
+  const char* fTypeKeys[kNtypes];      // type names
+
+
+  UInt_t GetObjIndex(TObject *obj);    // return object index
+  
+  ClassDef(AliDielectronCutQA,1) //Group of cuts
+};
+
+#endif
index ad68d85c2235d23f5e0b550e4f2be0432c5a348d..777348df5af1d77428981cc96a342384c1d92e7f 100644 (file)
@@ -225,6 +225,7 @@ void AliDielectronV0Cuts::InitEvent(AliVTrack *trk)
   }
   else if(ev->IsA() == AliAODEvent::Class()) {
     const AliAODEvent *aodEv = static_cast<const AliAODEvent*>(ev);
+    if(!aodEv->GetV0s()) return; // protection for nano AODs
 
     // loop over vertices
     for (Int_t ivertex=0; ivertex<aodEv->GetNumberOfV0s(); ++ivertex){