]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add copy track task
authormverweij <marta.verweij@cern.ch>
Thu, 4 Sep 2014 13:29:33 +0000 (15:29 +0200)
committermverweij <marta.verweij@cern.ch>
Thu, 4 Sep 2014 13:52:39 +0000 (15:52 +0200)
PWGJE/CMakelibPWGJEEMCALJetTasks.pkg
PWGJE/EMCALJetTasks/AliJetModelCopyTracks.cxx [new file with mode: 0644]
PWGJE/EMCALJetTasks/AliJetModelCopyTracks.h [new file with mode: 0644]
PWGJE/EMCALJetTasks/macros/AddTaskModelCopyTracks.C [new file with mode: 0644]
PWGJE/PWGJEEMCALJetTasksLinkDef.h

index 20b1c1265d77baed1d3a92a7922d76d729033cfc..0ca76b662fc18491f131d6bb97a35516374fd77c 100644 (file)
@@ -44,6 +44,7 @@ set ( SRCS
  EMCALJetTasks/AliJetEmbeddingTask.cxx
  EMCALJetTasks/AliJetFastSimulation.cxx
  EMCALJetTasks/AliJetModelBaseTask.cxx
+ EMCALJetTasks/AliJetModelCopyTracks.cxx
  EMCALJetTasks/AliJetModelMergeBranches.cxx
  EMCALJetTasks/AliJetRandomizerTask.cxx
  EMCALJetTasks/AliJetResponseMaker.cxx
diff --git a/PWGJE/EMCALJetTasks/AliJetModelCopyTracks.cxx b/PWGJE/EMCALJetTasks/AliJetModelCopyTracks.cxx
new file mode 100644 (file)
index 0000000..ab869dd
--- /dev/null
@@ -0,0 +1,151 @@
+// $Id$
+//
+// Jet model task to copy tracks while making small change
+// - make particles massless
+//
+// Author: M. Verweij
+
+#include "AliJetModelCopyTracks.h"
+
+#include <TClonesArray.h>
+#include <TFolder.h>
+#include <TLorentzVector.h>
+#include <TParticle.h>
+#include <TParticlePDG.h>
+#include <TRandom3.h>
+#include <TProfile.h>
+#include <TGrid.h>
+#include <TFile.h>
+#include <TF1.h>
+#include "AliAnalysisManager.h"
+#include "AliEMCALDigit.h"
+#include "AliEMCALGeometry.h"
+#include "AliEMCALRecPoint.h"
+#include "AliGenerator.h"
+#include "AliHeader.h"
+#include "AliLog.h"
+#include "AliPicoTrack.h"
+#include "AliRun.h"
+#include "AliRunLoader.h"
+#include "AliStack.h"
+#include "AliStack.h"
+#include "AliVCluster.h"
+#include "AliVEvent.h"
+
+ClassImp(AliJetModelCopyTracks)
+
+//________________________________________________________________________
+AliJetModelCopyTracks::AliJetModelCopyTracks() : 
+AliAnalysisTaskEmcal("AliJetModelCopyTracks",kTRUE),
+  fTracksOutName(""),
+  fTracksOut(0x0),
+  fParticleMass(kMassive),
+  fHistPtOut(0)
+{
+  // Default constructor.
+  SetMakeGeneralHistograms(kTRUE);
+}
+
+//________________________________________________________________________
+AliJetModelCopyTracks::AliJetModelCopyTracks(const char *name) : 
+  AliAnalysisTaskEmcal(name,kTRUE),
+  fTracksOutName(""),
+  fTracksOut(0x0),
+  fParticleMass(kMassive),
+  fHistPtOut(0)
+{
+  // Standard constructor.
+  SetMakeGeneralHistograms(kTRUE);
+}
+
+//________________________________________________________________________
+AliJetModelCopyTracks::~AliJetModelCopyTracks()
+{
+  // Destructor
+
+}
+
+//________________________________________________________________________
+void AliJetModelCopyTracks::ExecOnce() 
+{
+  // Exec only once.
+
+  AliAnalysisTaskEmcal::ExecOnce();
+
+  if (!fTracksOutName.IsNull()) {
+    fTracksOut = new TClonesArray("AliPicoTrack");
+    fTracksOut->SetName(fTracksOutName);
+    if (InputEvent()->FindListObject(fTracksOutName)) {
+      AliFatal(Form("%s: Collection %s is already present in the event!", GetName(), fTracksOutName.Data()));
+      return;
+    }
+    else {
+      InputEvent()->AddObject(fTracksOut);
+    }
+  }
+}
+
+//________________________________________________________________________
+void AliJetModelCopyTracks::UserCreateOutputObjects() 
+{
+  AliAnalysisTaskEmcal::UserCreateOutputObjects();
+
+  const Int_t nBinPt = 100;
+  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] + 1.0;
+    }
+  }
+
+  fHistPtOut = new TH1F("fHistPtOut","fHistPtOut;#it{p}_{T};N",nBinPt,binLimitsPt);
+  fOutput->Add(fHistPtOut);
+
+  PostData(1, fOutput); // Post data for ALL output slots > 0 here.
+}
+
+
+//________________________________________________________________________
+Bool_t AliJetModelCopyTracks::Run() 
+{
+  CopyTracks();
+  return kTRUE;
+}
+
+//________________________________________________________________________
+void AliJetModelCopyTracks::CopyTracks()
+{
+  //Apply toy detector simulation to tracks
+  Int_t nt = 0;
+  const Int_t nTracks = fTracks->GetEntriesFast();
+   for (Int_t i = 0; i < nTracks; ++i) {
+    AliPicoTrack *picotrack = static_cast<AliPicoTrack*>(fTracks->At(i));
+    if (!picotrack)
+      continue;
+
+    Double_t mass = picotrack->M();
+    if(fParticleMass==kMassless) mass = 0.;
+    if(fParticleMass==kPionMass) mass = 0.13957;
+    AliPicoTrack *track = new ((*fTracksOut)[nt]) AliPicoTrack(picotrack->Pt(),
+                                                              picotrack->Eta(),
+                                                              picotrack->Phi(),
+                                                              picotrack->Charge(),
+                                                              picotrack->GetLabel(),
+                                                              AliPicoTrack::GetTrackType(picotrack),
+                                                              picotrack->GetTrackEtaOnEMCal(),
+                                                              picotrack->GetTrackPhiOnEMCal(),
+                                                              picotrack->GetTrackPtOnEMCal(),
+                                                              picotrack->IsEMCAL(),
+                                                              mass); 
+    track->SetBit(TObject::kBitMask,1);
+    fHistPtOut->Fill(track->Pt());
+    nt++;
+   }
+}
+
+
+
+
+
diff --git a/PWGJE/EMCALJetTasks/AliJetModelCopyTracks.h b/PWGJE/EMCALJetTasks/AliJetModelCopyTracks.h
new file mode 100644 (file)
index 0000000..f040e81
--- /dev/null
@@ -0,0 +1,49 @@
+#ifndef ALIJETMODELCOPYTRACKS_H
+#define ALIJETMODELCOPYTRACKS_H
+
+// $Id$
+
+class TClonesArray;
+class TRandom3;
+class AliVParticle;
+class AliPicoTrack;
+
+#include "AliAnalysisTaskEmcal.h"
+
+class AliJetModelCopyTracks : public AliAnalysisTaskEmcal {
+ public:
+  enum ParticleMass {
+    kMassive  = 0,
+    kMassless = 1,
+    kPionMass = 2
+  };
+
+
+  AliJetModelCopyTracks();
+  AliJetModelCopyTracks(const char *name); 
+  virtual ~AliJetModelCopyTracks();
+
+  virtual void           UserCreateOutputObjects();
+
+  void                   SetTracksOutName(const char *n)          { fTracksOutName   = n;    }
+  void                   SetParticleMassType(ParticleMass m)      { fParticleMass    = m;    }
+
+  void                   ExecOnce();
+  Bool_t                 Run();
+
+  void                   CopyTracks();
+
+  TString                fTracksOutName;       // name of output track collection
+  TClonesArray          *fTracksOut;           //!output track collection
+  ParticleMass           fParticleMass;        // particle mass to use
+
+  //Output objects
+  TH1F     *fHistPtOut;                        //!pT spectrum of output particles
+  
+ private:
+  AliJetModelCopyTracks(const AliJetModelCopyTracks&);            // not implemented
+  AliJetModelCopyTracks &operator=(const AliJetModelCopyTracks&); // not implemented
+
+  ClassDef(AliJetModelCopyTracks, 1) // copy tracks class
+};
+#endif
diff --git a/PWGJE/EMCALJetTasks/macros/AddTaskModelCopyTracks.C b/PWGJE/EMCALJetTasks/macros/AddTaskModelCopyTracks.C
new file mode 100644 (file)
index 0000000..b53f4a8
--- /dev/null
@@ -0,0 +1,53 @@
+// $Id$
+
+AliJetModelCopyTracks* AddTaskModelCopyTracks(
+                                            const char     *tracksName1   = "Tracks",
+                                            const char     *tracksName2   = "Tracks2",
+                                            Int_t           massType      = AliJetModelCopyTracks::kMassless,
+                                            const char     *taskName      = "JetModelCopyTracks"
+                                           )
+{  
+  // Get the pointer to the existing analysis manager via the static access method.
+  //==============================================================================
+  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+  if (!mgr)
+  {
+    ::Error("AddTaskModelCopyTracks", "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("AddTaskModelCopyTracks", "This task requires an input event handler");
+    return NULL;
+  }
+  
+  //-------------------------------------------------------
+  // Init the task and do settings
+  //-------------------------------------------------------
+
+  AliJetModelCopyTracks *copyTask = new AliJetModelCopyTracks(taskName);
+  copyTask->SetTracksName(tracksName1);
+  copyTask->SetTracksOutName(tracksName2);
+
+  //-------------------------------------------------------
+  // Final settings, pass to manager and set the containers
+  //-------------------------------------------------------
+  mgr->AddTask(copyTask);
+    
+  // Create containers for input/output
+  mgr->ConnectInput (copyTask, 0, mgr->GetCommonInputContainer() );
+
+  TString contName = taskName;
+  contName += "_histos";
+  TString outputfile = Form("%s",AliAnalysisManager::GetCommonFileName());
+  AliAnalysisDataContainer *outc = mgr->CreateContainer(contName.Data(),
+                                                       TList::Class(),
+                                                       AliAnalysisManager::kOutputContainer,
+                                                       outputfile);
+  mgr->ConnectOutput(copyTask, 1, outc);
+
+  return copyTask;
+}
index 10c7862549592fde59d7c686872b4e8b7943ac02..cb987bf00cdbf9cc18ed069d06a607a9a3754a5b 100644 (file)
@@ -22,6 +22,7 @@
 #pragma link C++ class AliJetEmbeddingFromPYTHIATask+;
 #pragma link C++ class AliJetFastSimulation+;
 #pragma link C++ class AliJetModelBaseTask+;
+#pragma link C++ class AliJetModelCopyTracks+;
 #pragma link C++ class AliJetModelMergeBranches+;
 #pragma link C++ class AliJetRandomizerTask+;
 #pragma link C++ class AliJetConstituentTagCopier+;