StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
readStrangeMuDst.C
1 // $Id: readStrangeMuDst.C,v 2.0 2000/06/09 22:13:18 genevb Exp $
2 // $Log: readStrangeMuDst.C,v $
3 // Revision 2.0 2000/06/09 22:13:18 genevb
4 // Updated for version 2 of Strangeness mico DST package
5 //
6 // Revision 1.4 2000/04/12 16:16:55 genevb
7 // Remove unnecessary library loads
8 //
9 // Revision 1.3 2000/04/05 14:09:23 genevb
10 // Changed member functions again
11 //
12 // Revision 1.2 2000/03/29 20:58:58 genevb
13 // Modified StV0MuDst member functions
14 //
15 // Revision 1.1 2000/03/29 00:24:54 genevb
16 // Introduction of macro to use StStrangeMuDstMaker
17 //
18 //
19 //======================================================
20 // owner: Gene Van Buren, UCLA
21 // what it does: uses StStrangeMuDstMaker to read a micro DST
22 // and draw some histograms from it
23 //======================================================
24 
25 
26 TStopwatch clock;
27 
28 void load() {
29  gSystem->Load("St_base");
30  gSystem->Load("StChain");
31  gSystem->Load("StUtilities");
32  gSystem->Load("StarClassLibrary");
33  gSystem->Load("StEvent");
34  gSystem->Load("StStrangeMuDstMaker");
35 
36  // Create new canvas
37  TCanvas* c1 =
38  new TCanvas("c1","Getting Started with StV0MiniDst",0,0,600,600);
39 
40  c1->cd();
41  p1 = new TPad("p1","1st Pad",0.01,0.51,0.49,0.99,10,0,0);
42  p2 = new TPad("p2","2nd Pad",0.51,0.51,0.99,0.99,10,0,0);
43  p3 = new TPad("p3","3rd Pad",0.01,0.01,0.49,0.49,10,0,0);
44  p4 = new TPad("p4","4th Pad",0.51,0.01,0.99,0.49,10,0,0);
45 
46  c1->cd();
47  p1->Draw();
48  p2->Draw();
49  p3->Draw();
50  p4->Draw();
51 }
52 
53 void run() {
54  // Define some histograms
55  hX = new TH1F("mX","X coordinate",100,-50,50);
56  hY = new TH1F("mY","Y coordinate",100,-50,50);
57  hZ = new TH1F("mZ","Z coordinate",100,-100,100);
58  hMassLambda = new TH1F("mMassLambda","Lambda Mass",100,1.08,1.2);
59 
60  // Set number of events to analyse
61  const Int_t Nevents = 10000; // go to EOF
62 
63  StChain chain("myChain");
64  StStrangeMuDstMaker strangeDst("strangeMuDst");
65  strangeDst.DoV0(); // Selects V0 vertices for micro-DST
66 // strangeDst.DoXi(); // Selects Xi vertices for micro-DST
67 // strangeDst.DoMc(); // Reads in Monte Carlo information if available
68  strangeDst.SetRead(); // Sets "read" mode (using default file names)
69 
70  clock.Start(kTRUE);
71 
72  // Do init
73  Int_t ierr = chain.Init();
74  if( ierr ) { chain.Fatal(ierr,"on init"); return; }
75 
76  gMessMgr->Info("Here are the cuts used to create this strangeness DST:");
77  strangeDst.Cuts().List();
78 
79  // Loop over events
80  for( Int_t i=0; i<Nevents; i++ ) {
81  if( chain.Make(i) ) break;
82 
83  for( Int_t j=0; j<strangeDst.GetNV0(); j++ ) {
84  StV0MuDst *v0m = strangeDst.GetV0(j);
85  hX->Fill(v0m->decayVertexV0X());
86  hY->Fill(v0m->decayVertexV0Y());
87  hZ->Fill(v0m->decayVertexV0Z());
88  hMassLambda->Fill(v0m->massLambda());
89  }
90 
91  if( i != Nevents) chain.Clear();
92  printf("*** Finished processing event %d\n",i);
93  }
94 
95  // Finish
96  if( Nevents >= 1 ) {
97  chain.Finish();
98  }
99 
100  // Stop the stopwatch
101  clock.Stop();
102  clock.Print();
103 
104  p1->cd();
105  hX->Draw();
106  p2->cd();
107  hY->Draw();
108  p3->cd();
109  hZ->Draw();
110  p4->cd();
111  hMassLambda->Draw();
112 }
113 
114 void readStrangeMuDst() {
115  load();
116  run();
117 }
Float_t massLambda()
Mass assuming lambda hypothesis.
Definition: StV0I.hh:405
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
virtual Int_t Finish()
Definition: StChain.cxx:85
virtual Int_t Make()
Definition: StChain.cxx:110
Float_t decayVertexV0Z() const
Coordinates of decay vertex.
Definition: StV0MuDst.hh:113