Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TSelEvent.cxx
Go to the documentation of this file.
1 // @(#)root/proof:$Id$
2 // Author: Sangsu Ryu 22/06/2010
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2005, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TSelEvent
13 \ingroup proofbench
14 
15 Selector for PROOF I/O benchmark test.
16 For the I/O benchmark, event files are read in and histograms are filled.
17 For memory clean-up, dedicated files large enough to clean up memory
18 cache on the machine are read in. Or memory clean-up can be
19 accompolished by system call on Linux machine inside SlaveBegin(..)
20 which should be much faster the reading in large files.
21 
22 */
23 
24 #define TSelEvent_cxx
25 
26 #include "TSelEvent.h"
27 #include <TH1F.h>
28 #include <TStyle.h>
29 #include "TParameter.h"
30 #include "TProofBenchTypes.h"
31 #include "TTree.h"
32 #include "TCanvas.h"
33 #include "TFileInfo.h"
34 #include "THashList.h"
35 #include "TClonesArray.h"
36 #include "TRefArray.h"
37 #include <sys/stat.h>
38 #include <fcntl.h>
39 #include <unistd.h>
40 #include "TSystem.h"
41 #include "TROOT.h"
42 
43 ClassImp(TSelEvent);
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// Constructor
47 
48 TSelEvent::TSelEvent(TTree *)
49  : fReadType(0), fDebug(kFALSE), fCHist(0), fPtHist(0),
50  fNTracksHist(0), fEventName(0), fTracks(0), fHighPt(0), fMuons(0),
51  fH(0), b_event_fType(0), b_fEventName(0), b_event_fNtrack(0), b_event_fNseg(0),
52  b_event_fNvertex(0), b_event_fFlag(0), b_event_fTemperature(0),
53  b_event_fMeasures(0), b_event_fMatrix(0), b_fClosestDistance(0),
54  b_event_fEvtHdr(0), b_fTracks(0), b_fHighPt(0), b_fMuons(0),
55  b_event_fLastTrack(0), b_event_fWebHistogram(0), b_fH(0),
56  b_event_fTriggerBits(0), b_event_fIsValid(0)
57 {
58 }
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Constructor
62 
63 TSelEvent::TSelEvent()
64  : fReadType(0), fDebug(kFALSE), fCHist(0), fPtHist(0),
65  fNTracksHist(0), fEventName(0), fTracks(0), fHighPt(0), fMuons(0),
66  fH(0), b_event_fType(0), b_fEventName(0), b_event_fNtrack(0), b_event_fNseg(0),
67  b_event_fNvertex(0), b_event_fFlag(0), b_event_fTemperature(0),
68  b_event_fMeasures(0), b_event_fMatrix(0), b_fClosestDistance(0),
69  b_event_fEvtHdr(0), b_fTracks(0), b_fHighPt(0), b_fMuons(0),
70  b_event_fLastTrack(0), b_event_fWebHistogram(0), b_fH(0),
71  b_event_fTriggerBits(0), b_event_fIsValid(0)
72 {
73 }
74 
75 ////////////////////////////////////////////////////////////////////////////////
76 /// The Begin() function is called at the start of the query.
77 /// When running with PROOF Begin() is only called on the client.
78 /// The tree argument is deprecated (on PROOF 0 is passed).
79 
80 void TSelEvent::Begin(TTree *)
81 {
82  TString option = GetOption();
83 
84  //get parameters
85 
86  Bool_t found_readtype=kFALSE;
87  Bool_t found_debug=kFALSE;
88 
89  TIter nxt(fInput);
90  TString sinput;
91  TObject *obj;
92  while ((obj = nxt())){
93  sinput=obj->GetName();
94  //Info("Begin", "name=%s", sinput.Data());
95  if (sinput.Contains("PROOF_Benchmark_ReadType")){
96  if ((fReadType = dynamic_cast<TPBReadType *>(obj))) found_readtype = kTRUE;
97  continue;
98  }
99  if (sinput.Contains("PROOF_BenchmarkDebug")){
100  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
101  if (a){
102  fDebug= a->GetVal();
103  found_debug=kTRUE;
104  //Info("Begin", "PROOF_BenchmarkDebug=%d", fDebug);
105  }
106  else{
107  Error("Begin", "PROOF_BenchmarkDebug not type TParameter<Int_t>*");
108  }
109  continue;
110  }
111  }
112 
113  if (!found_readtype){
114  fReadType = new TPBReadType(TPBReadType::kReadOpt);
115  Warning("Begin", "PROOF_Benchmark_ReadType not found; using default: %d",
116  fReadType->GetType());
117  }
118  if (!found_debug){
119  Warning("Begin", "PROOF_BenchmarkDebug not found; using default: %d",
120  fDebug);
121  }
122 }
123 
124 ////////////////////////////////////////////////////////////////////////////////
125 /// The SlaveBegin() function is called after the Begin() function.
126 /// When running with PROOF SlaveBegin() is called on each slave server.
127 /// The tree argument is deprecated (on PROOF 0 is passed).
128 
129 void TSelEvent::SlaveBegin(TTree *tree)
130 {
131  Init(tree);
132 
133  TString option = GetOption();
134 
135  Bool_t found_readtype=kFALSE;
136  Bool_t found_debug=kFALSE;
137 
138  //fInput->Print("A");
139  TIter nxt(fInput);
140  TString sinput;
141  TObject *obj;
142  while ((obj = nxt())){
143  sinput=obj->GetName();
144  //Info("SlaveBegin", "name=%s", sinput.Data());
145  if (sinput.Contains("PROOF_Benchmark_ReadType")){
146  if ((fReadType = dynamic_cast<TPBReadType *>(obj))) found_readtype = kTRUE;
147  continue;
148  }
149  if (sinput.Contains("PROOF_BenchmarkDebug")){
150  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
151  if (a){
152  fDebug= a->GetVal();
153  found_debug=kTRUE;
154  //Info("SlaveBegin", "PROOF_BenchmarkDebug=%d", fDebug);
155  }
156  else{
157  Error("SlaveBegin", "PROOF_BenchmarkDebug not type TParameter"
158  "<Int_t>*");
159  }
160  continue;
161  }
162  }
163 
164  if (!found_readtype){
165  fReadType = new TPBReadType(TPBReadType::kReadOpt);
166  Warning("SlaveBegin", "PROOF_Benchmark_ReadType not found; using default: %d",
167  fReadType->GetType());
168  }
169  if (!found_debug){
170  Warning("SlaveBegin", "PROOF_BenchmarkDebug not found; using default: %d",
171  fDebug);
172  }
173 
174  fPtHist = new TH1F("pt_dist","p_{T} Distribution", 100, 0, 5);
175  fPtHist->SetDirectory(0);
176  fPtHist->GetXaxis()->SetTitle("p_{T}");
177  fPtHist->GetYaxis()->SetTitle("dN/p_{T}dp_{T}");
178 
179  fNTracksHist = new TH1F("ntracks_dist","N_{Tracks} per Event"
180  " Distribution", 100, 50, 150);
181  //enable rebinning
182  fNTracksHist->SetCanExtend(TH1::kAllAxes);
183  fNTracksHist->SetDirectory(0);
184  fNTracksHist->GetXaxis()->SetTitle("N_{Tracks}");
185  fNTracksHist->GetYaxis()->SetTitle("N_{Events}");
186 }
187 
188 ////////////////////////////////////////////////////////////////////////////////
189 /// The Process() function is called for each entry in the tree (or possibly
190 /// keyed object in the case of PROOF) to be processed. The entry argument
191 /// specifies which entry in the currently loaded tree is to be processed.
192 /// It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
193 /// to read either all or the required parts of the data. When processing
194 /// keyed objects with PROOF, the object is already loaded and is available
195 /// via the fObject pointer.
196 ///
197 /// This function should contain the "body" of the analysis. It can contain
198 /// simple or elaborate selection criteria, run algorithms on the data
199 /// of the event and typically fill histograms.
200 
201 Bool_t TSelEvent::Process(Long64_t entry)
202 {
203  // WARNING when a selector is used with a TChain, you must use
204  // the pointer to the current TTree to call GetEntry(entry).
205  // The entry is always the local entry number in the current tree.
206  // Assuming that fChain is the pointer to the TChain being processed,
207  // use fChain->GetTree()->GetEntry(entry).
208 
209  if (fReadType->GetType() != TPBReadType::kReadNotSpecified){
210  switch (fReadType->GetType()){
211  case TPBReadType::kReadFull:
212  // Full read
213  fChain->GetTree()->GetEntry(entry);
214  fNTracksHist->Fill(fNtrack);
215 
216  for(Int_t j=0;j<fTracks->GetEntries();j++){
217  Track* curtrack = dynamic_cast<Track*>(fTracks->At(j));
218  fPtHist->Fill(curtrack->GetPt(),1./curtrack->GetPt());
219  }
220  fTracks->Clear("C");
221  break;
222  case TPBReadType::kReadOpt:
223  // Partial read
224  b_event_fNtrack->GetEntry(entry);
225 
226  fNTracksHist->Fill(fNtrack);
227 
228  if (fNtrack>0) {
229  b_fTracks->GetEntry(entry);
230  for(Int_t j=0;j<fTracks->GetEntries();j++){
231  Track* curtrack = dynamic_cast<Track*>(fTracks->At(j));
232  fPtHist->Fill(curtrack->GetPt(),1./curtrack->GetPt());
233  }
234  fTracks->Clear("C");
235  }
236  break;
237  case TPBReadType::kReadNo:
238  // No read
239  break;
240  default:
241  Error("Process", "Read type not supported; %d", fReadType->GetType());
242  return kFALSE;
243  break;
244  }
245  }
246  return kTRUE;
247 }
248 
249 ////////////////////////////////////////////////////////////////////////////////
250 /// The SlaveTerminate() function is called after all entries or objects
251 /// have been processed. When running with PROOF SlaveTerminate() is called
252 /// on each slave server.
253 
254 void TSelEvent::SlaveTerminate()
255 {
256 }
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 /// The Terminate() function is the last function to be called during
260 /// a query. It always runs on the client, it can be used to present
261 /// the results graphically or save the results to file.
262 
263 void TSelEvent::Terminate()
264 {
265 }