]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
accept also neutral particles + add eta,phi histo
authormverweij <marta.verweij@cern.ch>
Mon, 7 Jul 2014 18:23:34 +0000 (20:23 +0200)
committermverweij <marta.verweij@cern.ch>
Wed, 30 Jul 2014 15:48:51 +0000 (11:48 -0400)
setter for charged only

PWGJE/EMCALJetTasks/AliJetEmbeddingFromGenTask.cxx
PWGJE/EMCALJetTasks/AliJetEmbeddingFromGenTask.h

index bcf971dc5936d023a1a39e15c3692f5130528f2d..6fe2b3ddf77825818ab158ca542ca2311d5eea1e 100644 (file)
@@ -36,7 +36,9 @@ AliJetEmbeddingFromGenTask::AliJetEmbeddingFromGenTask() :
   AliJetModelBaseTask("AliJetEmbeddingFromGenTask"),
   fGen(0),
   fMassless(kFALSE),
+  fChargedOnly(kFALSE),
   fHistPt(0),
+  fHistEtaPhi(0),
   fHistTrials(0),
   fHistXsection(0),
   fHistPtHard(0)
@@ -50,7 +52,9 @@ AliJetEmbeddingFromGenTask::AliJetEmbeddingFromGenTask(const char *name, Bool_t
   AliJetModelBaseTask(name,drawqa),
   fGen(0),
   fMassless(kFALSE),
+  fChargedOnly(kFALSE),
   fHistPt(0),
+  fHistEtaPhi(0),
   fHistTrials(0),
   fHistXsection(0),
   fHistPtHard(0)
@@ -78,6 +82,9 @@ void AliJetEmbeddingFromGenTask::UserCreateOutputObjects()
   fHistPt = new TH1F("fHistpt","fHistPt;#it{p}_{T};N",100,0.,100.);
   fOutput->Add(fHistPt);
 
+  fHistEtaPhi = new TH2F("fHistEtapHI","fHistEtaPhi;#eta;#varphi",100,-3.,3.,100.,0.,TMath::TwoPi());
+  fOutput->Add(fHistEtaPhi);
+
   fHistTrials = new TH1F("fHistTrials", "fHistTrials", 1, 0, 1);
   fHistTrials->GetYaxis()->SetTitle("trials");
   fOutput->Add(fHistTrials);
@@ -143,8 +150,7 @@ void AliJetEmbeddingFromGenTask::Run()
     if (!pdg) 
       continue;
     Int_t c = (Int_t)(TMath::Abs(pdg->Charge()));
-    if (c==0) 
-      continue;
+    if (fChargedOnly && c==0)  continue;
     Double_t pt = part->Pt();
     Double_t eta = part->Eta();
     Double_t phi = part->Phi();
@@ -163,6 +169,7 @@ void AliJetEmbeddingFromGenTask::Run()
     Double_t mass = part->GetMass();
     if(fMassless) mass = 0.;
     fHistPt->Fill(pt);
+    fHistEtaPhi->Fill(eta,phi);
     AddTrack(pt, eta, phi,0,0,0,0,0,0,c,mass);
   }
 
index 15115854c1948a62c5bca271f117adff7d0ca87d..a3b2995e2ac25924df37ad19388426f1b6678a9a 100644 (file)
@@ -19,8 +19,9 @@ class AliJetEmbeddingFromGenTask : public AliJetModelBaseTask {
   void           UserCreateOutputObjects();
   void           FillPythiaHistograms();
 
-  void           SetGen(AliGenerator *gen)      { fGen = gen    ; }
-  void           SetMasslessParticles(Bool_t b) { fMassless = b ; }
+  void           SetGen(AliGenerator *gen)      { fGen = gen       ; }
+  void           SetMasslessParticles(Bool_t b) { fMassless = b    ; }
+  void           SetChargedOnly(Bool_t b)       { fChargedOnly = b ; }
 
  protected:
   Bool_t         ExecOnce();
@@ -28,8 +29,10 @@ class AliJetEmbeddingFromGenTask : public AliJetModelBaseTask {
 
   AliGenerator  *fGen;                    //generator
   Bool_t         fMassless;               //make particles massless
+  Bool_t         fChargedOnly;            //accept only charged particles
 
-  TH1F          *fHistPt;                 //!pT spectrum of embedded particles
+  TH1F          *fHistPt;                 //!pT spectrum of generated particles
+  TH2F          *fHistEtaPhi;             //!eta-phi of generated particles
   TH1F          *fHistTrials;             //!trials from generator
   TProfile      *fHistXsection;           //!x-section from generator
   TH1           *fHistPtHard;             //!pt hard distribution
@@ -38,6 +41,6 @@ class AliJetEmbeddingFromGenTask : public AliJetModelBaseTask {
   AliJetEmbeddingFromGenTask(const AliJetEmbeddingFromGenTask&);            // not implemented
   AliJetEmbeddingFromGenTask &operator=(const AliJetEmbeddingFromGenTask&); // not implemented
 
-  ClassDef(AliJetEmbeddingFromGenTask, 3) // Jet embedding task
+  ClassDef(AliJetEmbeddingFromGenTask, 4) // Jet embedding task
 };
 #endif