]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
update on output object in AliMuonForwardTrackFinder.cxx
authorauras <auras@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 17 Feb 2012 09:18:01 +0000 (09:18 +0000)
committerauras <auras@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 17 Feb 2012 09:18:01 +0000 (09:18 +0000)
MFT/AliMuonForwardTrackFinder.cxx
MFT/CMakelibMFTsim.pkg
MFT/Config.C
MFT/MFTsimLinkDef.h

index ed02d7572a56420f4af2846e85f31d961c1e45bd..f2172f9d0770eaa1f30fc54760d451367f18b54c 100644 (file)
@@ -491,10 +491,16 @@ Int_t AliMuonForwardTrackFinder::LoadNextTrack() {
   
   // ------------------------------------- ...done!
 
-  if (fMuonTrackReco->GetMCLabel()>=0) fCountRealTracksWithRefMC++;
-  
   fLabelMC = fMuonTrackReco->GetMCLabel();
 
+  Int_t motherPdg=0;
+  if (fLabelMC>=0) {
+    fCountRealTracksWithRefMC++;
+    if (fStack->Particle(fLabelMC)->GetFirstMother() != -1) {
+      motherPdg = fStack->Particle(fStack->Particle(fLabelMC)->GetFirstMother())->GetPdgCode();
+    }
+  }
+
   CheckCurrentMuonTrackable();
 
   if (fMuonTrackReco->GetMatchTrigger()) fCountRealTracksWithRefMC_andTrigger++;
@@ -682,6 +688,7 @@ Int_t AliMuonForwardTrackFinder::LoadNextTrack() {
                                     Double_t(nFinalTracks),
                                     Double_t(fLabelMC>=0),
                                     xVtx, yVtx, zVtx,
+                                    motherPdg,
                                     Double_t(fMuonTrackReco->GetMatchTrigger()),
                                     Double_t(nClustersMC),
                                     Double_t(nGoodClusters),
@@ -742,6 +749,7 @@ Int_t AliMuonForwardTrackFinder::LoadNextTrack() {
                                       Double_t(nFinalTracks),
                                       Double_t(fLabelMC>=0),
                                       xVtx, yVtx, zVtx,
+                                      motherPdg,
                                       Double_t(fMuonTrackReco->GetMatchTrigger()),
                                       Double_t(nClustersMC),
                                       Double_t(nGoodClustersBestCandidate),
@@ -1241,9 +1249,9 @@ void AliMuonForwardTrackFinder::BookHistos() {
     
   }
 
-  fNtuFinalCandidates     = new TNtuple("ntuFinalCandidates",     "Final Candidates (ALL)", "run:event:muonTrack:nFinalCandidates:MCTrackRefExists:xVtx:yVtx:zVtx:triggerMatch:nClustersMC:nGoodClusters:pt:theta:eta:chi2AtPlane0:chi2AtPlane1:chi2AtPlane2:chi2AtPlane3:chi2AtPlane4:chi2AtPlane5:chi2AtPlane6:chi2AtPlane7:chi2AtPlane8");
+  fNtuFinalCandidates     = new TNtuple("ntuFinalCandidates",     "Final Candidates (ALL)", "run:event:muonTrack:nFinalCandidates:MCTrackRefExists:xVtx:yVtx:zVtx:motherPdg:triggerMatch:nClustersMC:nGoodClusters:pt:theta:eta:chi2AtPlane0:chi2AtPlane1:chi2AtPlane2:chi2AtPlane3:chi2AtPlane4:chi2AtPlane5:chi2AtPlane6:chi2AtPlane7:chi2AtPlane8");
 
-  fNtuFinalBestCandidates = new TNtuple("ntuFinalBestCandidates", "Final Best Candidates",  "run:event:muonTrack:nFinalCandidates:MCTrackRefExists:xVtx:yVtx:zVtx:triggerMatch:nClustersMC:nGoodClusters:pt:theta:eta:chi2AtPlane0:chi2AtPlane1:chi2AtPlane2:chi2AtPlane3:chi2AtPlane4:chi2AtPlane5:chi2AtPlane6:chi2AtPlane7:chi2AtPlane8:nClustersAtPlane0:nClustersAtPlane1:nClustersAtPlane2:nClustersAtPlane3:nClustersAtPlane4:nClustersAtPlane5:nClustersAtPlane6:nClustersAtPlane7:nClustersAtPlane8");
+  fNtuFinalBestCandidates = new TNtuple("ntuFinalBestCandidates", "Final Best Candidates",  "run:event:muonTrack:nFinalCandidates:MCTrackRefExists:xVtx:yVtx:zVtx:motherPdg:triggerMatch:nClustersMC:nGoodClusters:pt:theta:eta:chi2AtPlane0:chi2AtPlane1:chi2AtPlane2:chi2AtPlane3:chi2AtPlane4:chi2AtPlane5:chi2AtPlane6:chi2AtPlane7:chi2AtPlane8:nClustersAtPlane0:nClustersAtPlane1:nClustersAtPlane2:nClustersAtPlane3:nClustersAtPlane4:nClustersAtPlane5:nClustersAtPlane6:nClustersAtPlane7:nClustersAtPlane8");
 
 }
 
index ea6256fed39d56c65585632404f831ff439a9436..680067ff8a422d0b622060b79581c6cfaa25cb68 100644 (file)
@@ -25,7 +25,7 @@
 # SHLIBS - Shared Libraries and objects for linking (Executables only)           #
 #--------------------------------------------------------------------------------#
 
-set ( SRCS  AliMFT.cxx AliMFTHit.cxx AliMFTDigitizer.cxx )
+set ( SRCS  AliMFT.cxx AliMFTHit.cxx AliMFTDigitizer.cxx AliGenParamPionsKaons.cxx )
 
 string ( REPLACE ".cxx" ".h" HDRS "${SRCS}" )
 
index 4549300c4f0e30b952ec76a7b3750e5c08878ee5..b7137734857867a53d3cfa9024a75b845ef52b6d 100644 (file)
@@ -42,6 +42,7 @@
 
 enum PDCProc_t {kGenBox,
                kGenMuonLMR,
+               kGenPionKaon,
                kGenCorrHF,
                 kPythia6,
                kPythiaPerugia0, 
@@ -53,6 +54,7 @@ enum PDCProc_t {kGenBox,
 
 const Char_t* pprRunName[] = {"kGenBox",
                              "kGenMuonLMR",
+                             "kGenPionKaon",
                              "kGenCorrHF",
                              "kPythia6",
                              "kPythiaPerugia0", 
@@ -69,7 +71,7 @@ const Char_t* pprField[] = { "kNoField", "k5kG", "kFieldMax" };
 void LoadLibs();
 
 // ----------------------- Generator, field, beam energy,... ------------------------------------------------------------
-static PDCProc_t     proc     = kGenBox;
+static PDCProc_t     proc     = kGenPionKaon;
 static PDCProc_t     signal   = kGenBox;    // only in case kHijing2500Cocktail is the proc
 static Mag_t         mag      = k5kG;
 static Float_t       energy   = 14000.; // energy in CMS
@@ -154,6 +156,7 @@ void Config() {
   else if (proc == kGenBox)               gener = GenBox();
   else if (proc == kGenMuonLMR)           gener = GenMuonLMR();
   else if (proc == kGenCorrHF)            gener = GenCorrHF();
+  else if (proc == kGenPionKaon)          gener = GenParamPionKaon();
 
   // Size of the interaction diamond
   Float_t sigmaz  = 5.4 / TMath::Sqrt(2.);     // [cm]
@@ -225,7 +228,7 @@ void Config() {
   }
   if (iPIPE) {
     //    AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
-    AliPIPE *PIPE = new AliPIPEv4("PIPE", "Beam Pipe", 1.98, 0.08);  // cylindre with new adaptator pipe
+    AliPIPE *PIPE = new AliPIPEv4("PIPE", "Beam Pipe", 1.98, 0.08);
   }
   if (iZDC) {
     AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
@@ -294,6 +297,22 @@ AliGenerator* GenMuonLMR() {
 
 //====================================================================================================================================================
 
+AliGenerator* GenParamPionKaon() {
+  
+  AliGenParamPionsKaons *gener = new AliGenParamPionsKaons(100,"$ALICE_ROOT/MFT/PionKaonKinematics.root");
+  gener->SetPtRange(0, 5.);
+  gener->SetPhiRange(0., 360.);
+  gener->SetYRange(-10., 0.);
+  gener->SetOrigin(0.0, 0.0, 0.0);  // vertex position
+  gener->SetSigma(0.0, 0.0, 0.0);   // vertex position smearing
+  //  gener->SetCutOnChild(1);
+
+  return gener;
+
+}
+
+//====================================================================================================================================================
+
 AliGenerator* GenCorrHF() {
   
   AliGenCorrHF *gener = new AliGenCorrHF(1, 4, 14);  // for charm, 1 pair per event
index 56c06ee36cadc43509a56700dab96ca7b4cf985c..b97fe53b003f5dc8d703285f9bed243e53db03cd 100644 (file)
@@ -9,5 +9,6 @@
 #pragma link C++ class  AliMFT+;
 #pragma link C++ class  AliMFTHit+;
 #pragma link C++ class  AliMFTDigitizer+;
+#pragma link C++ class  AliGenParamPionsKaons+;
 
 #endif