StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StMuAnalysisMaker.cxx
1 
16 //
17 // Include header files.
18 #include "TFile.h"
19 #include "StMessMgr.h"
20 #include "TH1.h"
21 
22 #include "StMuDSTMaker/COMMON/StMuTrack.h"
23 #include "StMuDSTMaker/COMMON/StMuDst.h"
24 #include "StMuDSTMaker/COMMON/StMuDebug.h"
25 #include "StMuAnalysisMaker.h"
26 //
27 // Prototype
28 void muEventInfo(StMuEvent&, const int&);
29 
30 ClassImp(StMuAnalysisMaker)
31 
32 // The constructor. Initialize data members here.
33 StMuAnalysisMaker::StMuAnalysisMaker(const Char_t *name) : StMaker(name)
34 { mEventCounter = 0; mFile = 0; }
35 
36 StMuAnalysisMaker::~StMuAnalysisMaker() { /* noop */ }
37 //
38 // Called once at the beginning.
39 Int_t StMuAnalysisMaker::Init()
40 {
41  // Output file and histogram booking
42  mFileName = "muAnalysis.root";
43 
44  mGlobalPt = new TH1D("globalPt","globalPt",100,0.,3.);
45  mPrimaryPt = new TH1D("primaryPt","primaryPt",100,0.,3.);
46  mL3Pt = new TH1D("l3Pt","l3Pt",100,0.,3.);
47  mRefMult = new TH1D("refMult","refMult",100,0.,100.);
48 
49  return StMaker::Init();
50 }
51 //
52 // Called every event after Make().
53 void StMuAnalysisMaker::Clear(Option_t *opt)
54 {
56 }
57 //
58 // Called once at the end.
60 {
61 // Summarize the run.
62  cout << "StMuAnalysisMaker::Finish()\n";
63  cout << "\tProcessed " << mEventCounter << " events." << endl;
64 //
65 // Output histograms
66  mFile = new TFile(mFileName.c_str(), "RECREATE");
67  cout << "\tHistograms will be stored in file '"
68  << mFileName.c_str() << "'" << endl;
69 
70  mGlobalPt->Write();
71  mPrimaryPt->Write();
72  mL3Pt->Write();
73  mRefMult->Write();
74 //
75 // Write histos to file and close it.
76  if( mFile){
77  mFile->Write();
78  mFile->Close();
79  }
80 
81  return kStOK;
82 }
83 //
84 // This method is called every event.
86 {
87  mEventCounter++; // increase counter
88 
89  DEBUGVALUE2(mEventCounter);
90 // Get MuDst
91  StMuDst* mu;
92  mu = (StMuDst*) GetInputDS("MuDst");
93  DEBUGVALUE2(mu);
94 
95  if (!mu){
96  gMessMgr->Warning() << "StMuAnalysisMaker::Make : No MuDst" << endm;
97  return kStOK; // if no event, we're done
98  }
99 //
100 // Check StMuEvent branch
101  StMuEvent* muEvent;
102  muEvent = (StMuEvent*) mu->event();
103  if(muEvent) {
104  int refMult = muEvent->refMult();
105  mRefMult->Fill(refMult);
106  }
107 //
108 // Check track branches
109  StMuTrack* muTrack;
110  int nTracks;
111  nTracks= mu->globalTracks()->GetEntries();
112  printf("Global track # = %d\n",nTracks);
113  for (int l=0; l<nTracks; l++) {
114  muTrack = (StMuTrack*) mu->globalTracks(l);
115  if(muTrack) if (accept(muTrack)) mGlobalPt->Fill(muTrack->pt());
116  }
117  nTracks= mu->primaryTracks()->GetEntries();
118  printf("Primary track # = %d\n",nTracks);
119  for (int l=0; l<nTracks; l++) {
120  muTrack = (StMuTrack*) mu->primaryTracks(l);
121  if(muTrack) if (accept(muTrack)) mPrimaryPt->Fill(muTrack->pt());
122  }
123  nTracks= mu->l3Tracks()->GetEntries();
124  printf("L3 track # = %d\n",nTracks);
125  for (int l=0; l<nTracks; l++) {
126  muTrack = (StMuTrack*) mu->l3Tracks(l);
127  if(muTrack) if (accept(muTrack)) mL3Pt->Fill(muTrack->pt());
128  }
129 //
130 // Printout information of StMuEvent
131  if(muEvent) muEventInfo(*muEvent, mEventCounter);
132 
133  return kStOK;
134 }
135 //
136 // A simple track filter
137 bool StMuAnalysisMaker::accept(StMuTrack* track)
138 {
139 // check for positive flags.
140  return track && track->flag() >= 0;
141 }
static TObjArray * globalTracks()
returns pointer to the global tracks list
Definition: StMuDst.h:303
void Clear(Option_t *option="")
User defined functions.
Double_t pt() const
Returns pT at point of dca to primary vertex.
Definition: StMuTrack.h:256
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StMaker.cxx:634
static TClonesArray * l3Tracks()
returns pointer to the l3Tracks list
Definition: StMuDst.h:307
short flag() const
Returns flag, (see StEvent manual for type information)
Definition: StMuTrack.h:230
static TObjArray * primaryTracks()
returns pointer to a list of tracks belonging to the selected primary vertex
Definition: StMuDst.h:301
unsigned short refMult(int vtx_id=-1)
Reference multiplicity of charged particles as defined in StEventUtilities/StuRefMult.hh for vertex vtx_id (-1 is default index from StMuDst)
Definition: StMuEvent.cxx:195
static StMuEvent * event()
returns pointer to current StMuEvent (class holding the event wise information, e.g. event number, run number)
Definition: StMuDst.h:320
Definition: Stypes.h:40
A typical Analysis Class for MuDst.