StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SpinAnalysisTreeReaderPDSF.C
1 void SpinAnalysisTreeReaderPDSF(const long nevents = 20) {
2  LoadSpinTreeLibs();
3 
4  //create a new reader
5  StSpinTreeReader *reader = new StSpinTreeReader();
6 
7  //add some files to analyze, one at a time or in a text file
8  reader->selectDataset("$STAR/StRoot/StSpinPool/StSpinTree/datasets/run6_pdsf.dataset");
9  //reader->selectFile("./spinAnalyses_6119039.tree.root");
10 
11  //configure the branches you're interested in (default = true)
12  reader->connectJets = true;
13  reader->connectNeutralJets = false;
14  reader->connectChargedPions = true;
15  reader->connectBemcPions = true;
16  reader->connectBemcElectrons = true;
17  reader->connectEemcPions = false; //not added yet
18 
19  //optionally filter events by run and trigger
20  //reader->selectRunlist("$STAR/StRoot/StSpinPool/StSpinTree/filters/run6_jets.runlist");
21  reader->selectRun(7143025);
22 
23  //select events that passed hardware OR software trigger for any trigger in list
24  reader->selectTrigger(137221);
25  reader->selectTrigger(137222);
26  reader->selectTrigger(137611);
27  reader->selectTrigger(137622);
28  reader->selectTrigger(5);
29 
30  //we can change the OR to AND by doing
31  reader->requireDidFire = true;
32  reader->requireShouldFire = true;
33 
34  StJetSkimEvent *ev = reader->event();
35 
36  long entries = reader->GetEntries();
37  if(entries > nevents) entries = nevents;
38  for(int i=0; i<entries; i++) {
39  reader->GetEntry(i);
40 
41  //the basics: run, event, etc.
42  int runId = ev->runId();
43  int eventId = ev->eventId();
44 
45  printf("----------------Reading Event %d of %d----------------\n",i+1,entries);
46  printf("basics: Run = %d, Event = %d\n",runId,eventId);
47 
48  //triggers -- note prescale/threshold access
49  TClonesArray *trigs = ev->triggers();
50  for(int j=0; j<trigs->GetEntries(); j++) {
51  StJetSkimTrig* aTrig = (StJetSkimTrig*)trigs->At(j);
52  int trigId = aTrig->trigId();
53  bool didFire = aTrig->didFire();
54  int shouldFire = aTrig->shouldFire();
55 
56  StJetSkimTrigHeader *header = ev->trigHeader(trigId);
57  if(!header) {
58  printf("ERROR LOADING TRIGGER HEADER FOR %d\n",trigId);
59  break;
60  }
61  float prescale = header->prescale;
62 
63  printf("trigger = %6d prescale = %8.1f didFire = %d shouldFire =% d\n",trigId,prescale,didFire,shouldFire);
64  }
65 
66  //vertices -- bestVert() will return NULL if none were found
67  int nVertices = ev->vertices()->GetEntries();
68  StJetSkimVert *bestVert = ev->bestVert();
69  if(bestVert) {
70  printf("nVertices = %d position of best = %f\n",nVertices,bestVert->position()[2]);
71  }
72  else if(nVertices == 0){
73  printf("no vertices found in this event\n");
74  }
75  else {
76  printf("ERROR LOADING BEST VERTEX IN THIS EVENT\n");
77  }
78 
79  //for more details on event-level quantities see StJetMaker skim event documentation and examples
80 
81  //jets -- could also access using jets pointer from above
82  printf("nJets = %d\n",reader->nJets());
83  for(int j=0; j<reader->nJets(); j++) {
84  StJet* aJet = reader->jet(j);
85  double R = aJet->btowEtSum/aJet->Et();
86  if(bestVert) printf("jet pt=%7.4f jet eta=% 1.4f det eta=% 1.4f E_neu/E_tot=%1.4f\n",aJet->jetPt,aJet->jetEta,aJet->detEta(bestVert->position()[2]),R);
87  else printf("jet pt=%7.4f jet eta=% 1.4f no vertex E_neu/E_tot=%1.4f\n",aJet->jetPt,aJet->jetEta,R);
88  }
89 
90  //charged pions
91  printf("nChargedPions = %d\n",reader->nChargedPions());
92  for(int j=0; j<reader->nChargedPions(); j++) {
93  StChargedPionTrack* aPion = reader->chargedPion(j);
94  printf("pt=%7.4f eta=% 1.4f nSigmaPion=% 1.4f\n",aPion->pt(),aPion->eta(),aPion->nSigmaPion());
95  }
96 
97  //neutral pions
98  printf("nNeutralPions = %d\n",reader->nBemcPions());
99  for(int j=0; j<reader->nBemcPions(); j++) {
100  TPi0Candidate* aPi0 = reader->bemcPion(j);
101  printf("neutral pion mass=%f pt=%f eta=%f\n",aPi0->Mass(),aPi0->Pt(),aPi0->Eta());
102  }
103 
104  cout << "-----------------------------------------------------" << endl;
105  }
106 
107  delete reader;
108 }
109 
110 void LoadSpinTreeLibs() {
111  gSystem->Load("libPhysics");
112  gSystem->Load("St_base");
113  gSystem->Load("StChain");
114  gSystem->Load("St_Tables");
115  gSystem->Load("StEvent");
116  gSystem->Load("StDetectorDbMaker");
117  gSystem->Load("StEmcUtil");
118  gSystem->Load("StStrangeMuDstMaker");
119  gSystem->Load("StMuDSTMaker");
120  gSystem->Load("StSpinDbMaker");
121  gSystem->Load("StEmcTriggerMaker");
122  gSystem->Load("StJetFinder");
123  gSystem->Load("StJetMaker");
124  gSystem->Load("StMcEvent");
125  gSystem->Load("StChargedPionAnalysisMaker");
126  gSystem->Load("StSpinTree");
127 }
float nSigmaPion() const
Returns Craig&#39;s distance to the calculated dE/dx band for pions in units of sigma.
float jetPt
Pt (stored for convenience when drawing TTree)
Definition: StJet.h:125
float jetEta
Eta (stored for convenience when drawing TTree)
Definition: StJet.h:128
float pt() const
Returns pT at point of dca to primary vertex.
float eta() const
Returns pseudo rapidity at point of dca to primary vertex.
float btowEtSum
The summed Et from Barrel towers.
Definition: StJet.h:116
Definition: StJet.h:91