StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
filterStrangeMuDst.C
1 // $Id: filterStrangeMuDst.C,v 2.1 2003/02/10 16:50:08 genevb Exp $
2 // $Log: filterStrangeMuDst.C,v $
3 // Revision 2.1 2003/02/10 16:50:08 genevb
4 // simple updates
5 //
6 // Revision 2.0 2000/06/09 22:12:56 genevb
7 // Updated for version 2 of Strangeness mico DST package
8 //
9 // Revision 1.2 2000/04/12 16:16:55 genevb
10 // Remove unnecessary library loads
11 //
12 // Revision 1.1 2000/04/05 20:27:58 genevb
13 // Introduction of macro to filter strangeness micro DSTs
14 //
15 //
16 //======================================================
17 // owner: Gene Van Buren, UCLA
18 // what it does: uses StStrangeMuDstMaker to read a micro DST
19 // and filter it to a sub-micro DST
20 //
21 // - This example shows how one might select only those V0s
22 // whose mass under the lambda hypothesis falls in a range
23 // near that of that lambda.
24 // - Another possibility is to select entire events, and the
25 // lines commented out during the event loop with totalV0Pt
26 // show how one would filter out those events whose total
27 // Pt from V0's is greater than 350 GeV/c.
28 //======================================================
29 
30 
31 TStopwatch clock;
32 
33 void load() {
34  gSystem->Load("St_base");
35  gSystem->Load("StChain");
36  gSystem->Load("StUtilities");
37  gSystem->Load("StarClassLibrary");
38  gSystem->Load("StEvent");
39  gSystem->Load("StStrangeMuDstMaker");
40 }
41 
42 void run() {
43 
44  // Set number of events to analyse
45  const Int_t Nevents = 10000; // go to EOF
46 
47  StChain chain("myChain");
48 
49  // The maker for the new micro DST must be constructed _before_ the
50  // maker to read the old micro DST. This is because the copying is
51  // done during chain.Clear(), and the new maker's Clear() must be
52  // called to do the copying before the old maker's Clear() is called,
53  // erasing the event.
54 
55  StStrangeMuDstMaker strangeNewDst("strangeNewMuDst");
56  strangeNewDst.SetWrite("newEvMuDst.root"); //
57  strangeNewDst.DoV0(); // Selects V0 vertices for new micro-DST
58  strangeNewDst.DoMc(); // Keep MC info if it is available
59 
60  StStrangeMuDstMaker strangeOldDst("strangeOldMuDst");
61  strangeOldDst.SetRead(); // Sets "read" mode (using default file names)
62  // DoV0() and DoMc() are autmatically called for the old maker by the new.
63 
64  // Now we tell the new maker that it will create a sub-DST of the old one.
65  strangeNewDst.SubDst(&strangeOldDst);
66 
67  // Next, any additional cuts that are being made should be added to
68  // the cuts information in the new DST.
69  strangeNewDst.Cuts().Add("Lambda mass","1.11 < mass < 1.12");
70  // strangeNewDst.Cuts().Add("Large Pt total","Pt total > 350 GeV/c");
71 
72  clock.Start(kTRUE);
73 
74  // Do init
75  Int_t ierr = chain.Init();
76  if( ierr ) { chain.Fatal(ierr,"on init"); return; }
77 
78  // Loop over events
79  for( Int_t i=0; i<Nevents; i++ ) {
80  if( chain.Make() ) break;
81 
82 // Float_t totalV0Pt = 0.;
83  for( Int_t j=0; j<strangeOldDst.GetNV0(); j++ ) {
84  StV0MuDst *v0m = strangeOldDst.GetV0(j);
85  Float_t ml = v0m->massLambda();
86  if ((ml>1.11) && (ml<1.12)) strangeNewDst.SelectV0(j);
87 // totalV0Pt += v0m->ptV0();
88  }
89 // if (totalV0Pt > 350.) strangeNewDst.SelectEvent();
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 
105 void filterStrangeMuDst() {
106  load();
107  run();
108 }
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