]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Update from Salvatore
authorloizides <loizides@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 28 Feb 2013 13:57:59 +0000 (13:57 +0000)
committerloizides <loizides@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 28 Feb 2013 13:57:59 +0000 (13:57 +0000)
PWG/EMCAL/AliAnalysisTaskEmcal.cxx
PWG/EMCAL/AliAnalysisTaskEmcal.h
PWGJE/EMCALJetTasks/AliJetConstituentTagCopier.cxx [new file with mode: 0644]
PWGJE/EMCALJetTasks/AliJetConstituentTagCopier.h [new file with mode: 0644]
PWGJE/EMCALJetTasks/macros/AddTaskTagCopier.C [new file with mode: 0644]

index eadce723f17502a67d01859cdf2bd062a4d40a06..b7b7b4af0fe7d7dfb3000b824be9c2cc601746ab 100644 (file)
@@ -259,7 +259,7 @@ Bool_t AliAnalysisTaskEmcal::AcceptCluster(AliVCluster *clus, Bool_t acceptMC) c
 }
 
 //________________________________________________________________________
-Bool_t AliAnalysisTaskEmcal::AcceptTrack(AliVTrack *track, Bool_t acceptMC) const
+Bool_t AliAnalysisTaskEmcal::AcceptTrack(AliVParticle *track, Bool_t acceptMC) const
 {
   // Return true if track is accepted.
 
index 079831578e92624d2a2801929d6137ce19aa7a13..0cb3245eb1b63fa1bc4b81120d5ccf6bf15750d7 100644 (file)
@@ -10,6 +10,7 @@ class AliEmcalParticle;
 class AliMCParticle;
 class AliVCluster;
 class AliVTrack;
+class AliVParticle;
 class AliVCaloCells;
 class TH1F;
 class AliEMCALGeometry;
@@ -63,7 +64,7 @@ class AliAnalysisTaskEmcal : public AliAnalysisTaskSE {
  protected:
   Bool_t                      AcceptCluster(AliVCluster        *clus,  Bool_t acceptMC = kTRUE) const;
   Bool_t                      AcceptEmcalPart(AliEmcalParticle *part,  Bool_t acceptMC = kTRUE) const;
-  Bool_t                      AcceptTrack(AliVTrack            *track, Bool_t acceptMC = kTRUE) const;
+  Bool_t                      AcceptTrack(AliVParticle         *track, Bool_t acceptMC = kTRUE) const;
   virtual void                ExecOnce();
   virtual Bool_t              FillGeneralHistograms();
   virtual Bool_t              FillHistograms()                                     { return kTRUE                 ; }
diff --git a/PWGJE/EMCALJetTasks/AliJetConstituentTagCopier.cxx b/PWGJE/EMCALJetTasks/AliJetConstituentTagCopier.cxx
new file mode 100644 (file)
index 0000000..59b4e2a
--- /dev/null
@@ -0,0 +1,182 @@
+// $Id: AliJetConstituentTagCopier.cxx  $
+//
+// Copy tags from particle level constituent to detector level
+//
+// Author: S. Aiola
+
+#include "AliJetConstituentTagCopier.h"
+
+#include <TClonesArray.h>
+#include <TH1I.h>
+
+#include "AliVCluster.h"
+#include "AliVParticle.h"
+#include "AliEmcalParticle.h"
+#include "AliLog.h"
+
+ClassImp(AliJetConstituentTagCopier)
+
+//________________________________________________________________________
+AliJetConstituentTagCopier::AliJetConstituentTagCopier() : 
+  AliAnalysisTaskEmcal("AliJetConstituentTagCopier", kFALSE),
+  fMCParticlesName(),
+  fMCParticles(0),
+  fMCParticlesMap(0)
+{
+  // Default constructor.
+}
+
+//________________________________________________________________________
+AliJetConstituentTagCopier::AliJetConstituentTagCopier(const char *name) : 
+  AliAnalysisTaskEmcal(name, kFALSE),
+  fMCParticlesName("MCParticles"),
+  fMCParticles(0),
+  fMCParticlesMap(0)
+{
+  // Standard constructor.
+}
+
+//________________________________________________________________________
+AliJetConstituentTagCopier::~AliJetConstituentTagCopier()
+{
+  // Destructor
+}
+
+//________________________________________________________________________
+void AliJetConstituentTagCopier::ExecOnce()
+{
+  // Execute once.
+
+  if (!fMCParticles) {
+    fMCParticles = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fMCParticlesName));
+    if (!fMCParticles) {
+      AliError(Form("%s: Could not retrieve MC particles %s!", GetName(), fMCParticlesName.Data()));
+      return;
+    }
+    else if (!fMCParticles->GetClass()->GetBaseClass("AliVParticle")) {
+      AliError(Form("%s: Collection %s does not contain AliVParticle objects!", GetName(), fMCParticlesName.Data())); 
+      fMCParticles = 0;
+      return;
+    }
+  }
+
+  if (!fMCParticlesMap) {
+    fMCParticlesMap = dynamic_cast<TH1I*>(InputEvent()->FindListObject(fMCParticlesName + "_Map"));
+    // this is needed to map the MC labels with the indexes of the MC particle collection
+      // if teh map is not given, the MC labels are assumed to be consistent with the indexes (which is not the case if AliEmcalMCTrackSelector is used)
+    if (!fMCParticlesMap) {
+      AliWarning(Form("%s: Could not retrieve map for MC particles %s! Will assume MC labels consistent with indexes...", GetName(), fMCParticlesName.Data())); 
+      fMCParticlesMap = new TH1I("tracksMap","tracksMap",9999,0,1);
+      for (Int_t i = 0; i < 9999; i++) {
+       fMCParticlesMap->SetBinContent(i,i);
+      }
+    }
+  }
+
+  AliAnalysisTaskEmcal::ExecOnce();
+}
+
+//________________________________________________________________________
+Bool_t AliJetConstituentTagCopier::Run()
+{
+  if (fTracks) {
+    if (fTracks->GetClass()->GetBaseClass("AliVParticle"))
+      DoTrackLoop(fTracks);
+    else if (fTracks->GetClass()->GetBaseClass("AliEmcalParticle"))
+      DoEmcalParticleLoop(fTracks);
+    else 
+      AliError(Form("%s: Object type not recognized in collection %s. Nothing will be done.", GetName(), fTracks->GetName()));
+  }
+
+  if (fCaloClusters) {
+    if (fCaloClusters->GetClass()->GetBaseClass("AliVCluster"))
+      DoClusterLoop(fCaloClusters);
+    else if (fCaloClusters->GetClass()->GetBaseClass("AliEmcalParticle"))
+      DoEmcalParticleLoop(fCaloClusters);
+    else 
+      AliError(Form("%s: Object type not recognized in collection %s. Nothing will be done.", GetName(), fCaloClusters->GetName()));
+  }
+
+  return kTRUE;
+}
+
+//________________________________________________________________________
+void AliJetConstituentTagCopier::DoClusterLoop(TClonesArray *array)
+{
+  for (Int_t i = 0; i < array->GetEntries(); i++) {
+    AliVCluster *cluster = static_cast<AliVCluster*>(array->At(i));
+    if (!cluster) {
+      AliError(Form("%s: Could not get cluster %d", GetName(), i));
+      continue;
+    }
+    if (!AcceptCluster(cluster))
+      continue;
+    Int_t mcLabel = cluster->GetLabel();
+    if (mcLabel > 0) {
+      Int_t index = fMCParticlesMap->GetBinContent(mcLabel);
+      AliVParticle *part = static_cast<AliVParticle*>(fMCParticles->At(index));
+      if (!part) {
+       AliError(Form("%s: Could not get MC particle %d", GetName(), index));
+       continue;
+      }
+      UInt_t bits = (UInt_t)part->TestBits(TObject::kBitMask);
+      cluster->SetBit(bits);
+    }
+  }
+}
+
+//________________________________________________________________________
+void AliJetConstituentTagCopier::DoTrackLoop(TClonesArray *array)
+{
+  for (Int_t i = 0; i < array->GetEntries(); i++) {
+    AliVParticle *track = static_cast<AliVParticle*>(array->At(i));
+    if (!track) {
+      AliError(Form("%s: Could not get track %d", GetName(), i));
+      continue;
+    }
+    if (!AcceptTrack(track))
+      continue;
+    Int_t mcLabel = track->GetLabel();
+    if (mcLabel != 0) {
+      Int_t index = fMCParticlesMap->GetBinContent(mcLabel);
+      AliVParticle *part = static_cast<AliVParticle*>(fMCParticles->At(index));
+      if (!part) {
+       AliError(Form("%s: Could not get MC particle %d", GetName(), index));
+       continue;
+      }
+      UInt_t bits = (UInt_t)part->TestBits(TObject::kBitMask);
+      track->SetBit(bits);
+    }
+  }
+}
+
+//________________________________________________________________________
+void AliJetConstituentTagCopier::DoEmcalParticleLoop(TClonesArray *array)
+{
+  for (Int_t i = 0; i < array->GetEntries(); i++) {
+    AliEmcalParticle *emcpart = static_cast<AliEmcalParticle*>(array->At(i));
+    if (!emcpart) {
+      AliError(Form("%s: Could not get EmcalParticle %d", GetName(), i));
+      continue;
+    }
+    if (!AcceptEmcalPart(emcpart))
+      continue;
+    AliVCluster *cluster = emcpart->GetCluster();
+    AliVParticle *track = emcpart->GetTrack();
+    Int_t mcLabel = 0;
+    if (cluster)
+      mcLabel = cluster->GetLabel();
+    else if (track)
+      mcLabel = track->GetLabel();
+    if (mcLabel != 0) {
+      Int_t index = fMCParticlesMap->GetBinContent(mcLabel);
+      AliVParticle *part = static_cast<AliVParticle*>(fMCParticles->At(index));
+      if (!part) {
+       AliError(Form("%s: Could not get MC particle %d", GetName(), index));
+       continue;
+      }
+      UInt_t bits = (UInt_t)part->TestBits(TObject::kBitMask);
+      emcpart->SetBit(bits);
+    }
+  }
+}
diff --git a/PWGJE/EMCALJetTasks/AliJetConstituentTagCopier.h b/PWGJE/EMCALJetTasks/AliJetConstituentTagCopier.h
new file mode 100644 (file)
index 0000000..03a11e4
--- /dev/null
@@ -0,0 +1,37 @@
+#ifndef ALIJETCONSTITUENTTAGCOPIER_H
+#define ALIJETCONSTITUENTTAGCOPIER_H
+
+// $Id: AliJetConstituentTagCopier.h  $
+
+#include "AliAnalysisTaskEmcal.h"
+
+class TClonesArray;
+class TString;
+
+class AliJetConstituentTagCopier : public AliAnalysisTaskEmcal {
+ public:
+  AliJetConstituentTagCopier();
+  AliJetConstituentTagCopier(const char *name);
+  virtual ~AliJetConstituentTagCopier();
+
+  void                        SetMCParticlesName(const char *n)      { fMCParticlesName       = n         ; }
+
+ protected:
+  void                        ExecOnce();
+  Bool_t                      Run();
+  void                        DoClusterLoop(TClonesArray *array);
+  void                        DoTrackLoop(TClonesArray *array);
+  void                        DoEmcalParticleLoop(TClonesArray *array);
+
+  TString                     fMCParticlesName;                       // name of MC particle collection
+  TClonesArray               *fMCParticles;                           //!MC particle collection
+  TH1I                       *fMCParticlesMap;                        //!MC particle map
+
+ private:
+  AliJetConstituentTagCopier(const AliJetConstituentTagCopier&);            // not implemented
+  AliJetConstituentTagCopier &operator=(const AliJetConstituentTagCopier&); // not implemented
+
+  ClassDef(AliJetConstituentTagCopier, 1) // Copy tags from particle level constituent to detector level
+};
+
+#endif
diff --git a/PWGJE/EMCALJetTasks/macros/AddTaskTagCopier.C b/PWGJE/EMCALJetTasks/macros/AddTaskTagCopier.C
new file mode 100644 (file)
index 0000000..c019fe3
--- /dev/null
@@ -0,0 +1,73 @@
+// $Id: AddTaskSAQA.C 60163 2013-01-03 09:37:04Z loizides $
+
+AliJetConstituentTagCopier* AddTaskTagCopier(
+  const char *ntracks            = "Tracks",
+  const char *nclusters          = "CaloClusters",
+  const char *nmcparticles       = "MCParticles",
+  Double_t    trackptcut         = 0.15,
+  Double_t    clusptcut          = 0.30,
+  UInt_t      type               = AliAnalysisTaskEmcal::kTPC,
+  const char *taskname           = "AliJetConstituentTagCopier"
+)
+{  
+  // Get the pointer to the existing analysis manager via the static access method.
+  //==============================================================================
+  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+  if (!mgr)
+  {
+    ::Error("AddTaskTagCopier", "No analysis manager to connect to.");
+    return NULL;
+  }  
+  
+  // Check the analysis type using the event handlers connected to the analysis manager.
+  //==============================================================================
+  if (!mgr->GetInputEventHandler())
+  {
+    ::Error("AddTaskTagCopier", "This task requires an input event handler");
+    return NULL;
+  }
+  
+  //-------------------------------------------------------
+  // Init the task and do settings
+  //-------------------------------------------------------
+  
+  TString name(taskname);
+  if (strcmp(ntracks,"")) {
+    name += "_";
+    name += ntracks;
+  }
+  if (strcmp(nclusters,"")) {
+    name += "_";
+    name += nclusters;
+  }
+  if (strcmp(nmcparticles,"")) {
+    name += "_";
+    name += nmcparticles;
+  }
+  if (type == AliAnalysisTaskEmcal::kTPC) 
+    name += "_TPC";
+  else if (type == AliAnalysisTaskEmcal::kEMCAL) 
+    name += "_EMCAL";
+  else if (type == AliAnalysisTaskEmcal::kUser) 
+    name += "_USER";
+
+  AliJetConstituentTagCopier* task = new AliJetConstituentTagCopier(name);
+  task->SetTracksName(ntracks);
+  task->SetClusName(nclusters);
+  task->SetMCParticlesName(nmcparticles);
+  task->SetTrackPtCut(trackptcut);
+  task->SetClusPtCut(clusptcut);
+  task->SetAnaType(type);
+
+  //-------------------------------------------------------
+  // Final settings, pass to manager and set the containers
+  //-------------------------------------------------------
+
+  mgr->AddTask(task);
+  
+  // Create containers for input/output
+  AliAnalysisDataContainer *cinput1  = mgr->GetCommonInputContainer()  ;
+  mgr->ConnectInput(task, 0, cinput1);
+
+  return task;
+}