ROOT  6.06/09
Reference Guide
TSelHist.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 //////////////////////////////////////////////////////////////////////////
13 // //
14 // TSelHist //
15 // PROOF selector for CPU-intensive benchmark test. //
16 // Events are generated and 1-D, 2-D, and/or 3-D histograms are filled. //
17 // //
18 //////////////////////////////////////////////////////////////////////////
19 
20 #define TSelHist_cxx
21 
22 #include "TSelHist.h"
23 #include "TProofBenchTypes.h"
24 #include <TCanvas.h>
25 #include <TPaveText.h>
26 #include <TFormula.h>
27 #include <TF1.h>
28 #include <TH1F.h>
29 #include <TH2F.h>
30 #include <TH3F.h>
31 #include <TMath.h>
32 #include <TRandom3.h>
33 #include <TString.h>
34 #include <TStyle.h>
35 #include <TSystem.h>
36 #include <TParameter.h>
37 #include <TROOT.h>
38 
40 
41 ////////////////////////////////////////////////////////////////////////////////
42 ///Constructor
43 
45  : fHistType(0), fNHists(16), fDraw(0), fHist1D(0), fHist2D(0), fHist3D(0),
46  fRandom(0), fCHist1D(0), fCHist2D(0), fCHist3D(0)
47 {
48 }
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 /// Destructor
52 
54 {
55  //if (fRandom) delete fRandom;
57 
58  // Info("TSelHist","destroying ...");
59 
60  if (!fDraw){
61  for (Int_t i=0; i < fNHists; i++) {
62  if (fHist1D && fHist1D[i] && !fOutput->FindObject(fHist1D[i])) {
63  SafeDelete(fHist1D[i]);
64  }
65  if (fHist2D && fHist2D[i] && !fOutput->FindObject(fHist2D[i])) {
66  SafeDelete(fHist2D[i]);
67  }
68  if (fHist3D && fHist3D[i] && !fOutput->FindObject(fHist3D[i])) {
69  SafeDelete(fHist3D[i]);
70  }
71  }
72  }
76 }
77 
78 ////////////////////////////////////////////////////////////////////////////////
79 /// The Begin() function is called at the start of the query.
80 /// When running with PROOF Begin() is only called on the client.
81 /// The tree argument is deprecated (on PROOF 0 is passed).
82 
83 void TSelHist::Begin(TTree * /*tree*/)
84 {
85  TString option = GetOption();
86 
87  Bool_t found_histtype=kFALSE;
88  Bool_t found_nhists=kFALSE;
89  Bool_t found_draw=kFALSE;
90 
91  TIter nxt(fInput);
92  TString sinput;
93  TObject *obj;
94 
95 
96  while ((obj = nxt())){
97  sinput=obj->GetName();
98  //Info("Begin", "object name=%s", sinput.Data());
99  if (sinput.Contains("PROOF_Benchmark_HistType")){
100  if ((fHistType = dynamic_cast<TPBHistType *>(obj))) found_histtype = kTRUE;
101  continue;
102  }
103  if (sinput.Contains("PROOF_BenchmarkNHists")){
104  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
105  if (a){
106  fNHists= a->GetVal();
107  found_nhists=kTRUE;
108  //Info("Begin", "PROOF_BenchmarkNHists=%d", fNHists);
109  }
110  else{
111  Error("Begin", "PROOF_BenchmarkNHists not type TParameter<Int_t>*");
112  }
113  continue;
114  }
115  if (sinput.Contains("PROOF_BenchmarkDraw")){
116  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
117  if (a){
118  fDraw= a->GetVal();
119  found_draw=kTRUE;
120  //Info("Begin", "PROOF_BenchmarkDraw=%d", fDraw);
121  }
122  else{
123  Error("Begin", "PROOF_BenchmarkDraw not type TParameter<Int_t>*");
124  }
125  continue;
126  }
127  }
128 
129  if (!found_histtype){
131  Warning("Begin", "PROOF_Benchmark_HistType not found; using default: %d",
132  (Int_t) fHistType->GetType());
133  }
134  if (!found_nhists){
135  Warning("Begin", "PROOF_BenchmarkNHists not found; using default: %d",
136  fNHists);
137  }
138  if (!found_draw){
139  Warning("Begin", "PROOF_BenchmarkDraw not found; using default: %d",
140  fDraw);
141  }
142 
143  if (fDraw) {
147  }
148 }
149 
150 ////////////////////////////////////////////////////////////////////////////////
151 /// The SlaveBegin() function is called after the Begin() function.
152 /// When running with PROOF SlaveBegin() is called on each slave server.
153 /// The tree argument is deprecated (on PROOF 0 is passed).
154 
155 void TSelHist::SlaveBegin(TTree * /*tree*/)
156 {
157  TString option = GetOption();
158 
159  Bool_t found_histtype=kFALSE;
160  Bool_t found_nhists=kFALSE;
161  Bool_t found_draw=kFALSE;
162 
163  TIter nxt(fInput);
164  TString sinput;
165  TObject *obj;
166 
167  while ((obj = nxt())){
168  sinput=obj->GetName();
169  //Info("SlaveBegin", "object name=%s", sinput.Data());
170  if (sinput.Contains("PROOF_Benchmark_HistType")){
171  if ((fHistType = dynamic_cast<TPBHistType *>(obj))) found_histtype = kTRUE;
172  continue;
173  }
174  if (sinput.Contains("PROOF_BenchmarkNHists")){
175  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
176  if (a){
177  fNHists=a->GetVal();
178  found_nhists=kTRUE;
179  //Info("SlaveBegin", "PROOF_BenchmarkNHists=%d", fNHists);
180  }
181  else{
182  Error("SlaveBegin", "PROOF_BenchmarkNHists not type TParameter"
183  "<Int_t>*");
184  }
185  continue;
186  }
187  if (sinput.Contains("PROOF_BenchmarkDraw")){
188  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
189  if (a){
190  fDraw=a->GetVal();
191  found_draw=kTRUE;
192  //Info("SlaveBegin", "PROOF_BenchmarkDraw=%d", fDraw);
193  }
194  else{
195  Error("SlaveBegin", "PROOF_BenchmarkDraw not type TParameter"
196  "<Int_t>*");
197  }
198  continue;
199  }
200  }
201 
202  if (!found_histtype){
204  Warning("SlaveBegin", "PROOF_Benchmark_HistType not found; using default: %d",
205  fHistType->GetType());
206  }
207  if (!found_nhists){
208  Warning("SlaveBegin", "PROOF_BenchmarkNHists not found; using default: %d",
209  fNHists);
210  }
211  if (!found_draw){
212  Warning("SlaveBegin", "PROOF_BenchmarkDraw not found; using default: %d",
213  fDraw);
214  }
215 
216  // Create the histogram
218  fHist1D = new TH1F*[fNHists];
219  for (Int_t i=0; i < fNHists; i++) {
220  fHist1D[i] = new TH1F(Form("h1d_%d",i), Form("h1d_%d",i), 100, -3., 3.);
222  if (fDraw) fOutput->Add(fHist1D[i]);
223  }
224  }
226  fHist2D = new TH2F*[fNHists];
227  for (Int_t i=0; i < fNHists; i++) {
228  fHist2D[i] = new TH2F(Form("h2d_%d",i), Form("h2d_%d",i), 100, -3., 3.,
229  100, -3., 3.);
231  if (fDraw) fOutput->Add(fHist2D[i]);
232  }
233  }
235  fHist3D = new TH3F*[fNHists];
236  for (Int_t i=0; i < fNHists; i++) {
237  fHist3D[i] = new TH3F(Form("h3d_%d",i), Form("h3d_%d",i), 100, -3., 3.,
238  100, -3., 3., 100, -3., 3.);
240  if (fDraw) fOutput->Add(fHist3D[i]);
241  }
242  }
243  // Set random seed
244  fRandom = new TRandom3(0);
245 }
246 
247 ////////////////////////////////////////////////////////////////////////////////
248 /// The Process() function is called for each entry in the tree (or possibly
249 /// keyed object in the case of PROOF) to be processed. The entry argument
250 /// specifies which entry in the currently loaded tree is to be processed.
251 /// It can be passed to either TSelHist::GetEntry() or TBranch::GetEntry()
252 /// to read either all or the required parts of the data. When processing
253 /// keyed objects with PROOF, the object is already loaded and is available
254 /// via the fObject pointer.
255 ///
256 /// This function should contain the "body" of the analysis. It can contain
257 /// simple or elaborate selection criteria, run algorithms on the data
258 /// of the event and typically fill histograms.
259 ///
260 /// The processing can be stopped by calling Abort().
261 ///
262 /// Use fStatus to set the return value of TTree::Process().
263 ///
264 /// The return value is currently not used.
265 
267 {
268  Double_t x, y, z;
270  for (Int_t i=0; i < fNHists; i++) {
271  if (fRandom && fHist1D[i]) {
272  x = fRandom->Gaus(0.,1.);
273  fHist1D[i]->Fill(x);
274  }
275  }
276  }
278  for (Int_t i=0; i < fNHists; i++) {
279  if (fRandom && fHist2D[i]) {
280  x = fRandom->Gaus(0.,1.);
281  y = fRandom->Gaus(0.,1.);
282  fHist2D[i]->Fill(x, y);
283  }
284  }
285  }
287  for (Int_t i=0; i < fNHists; i++) {
288  if (fRandom && fHist3D[i]) {
289  x = fRandom->Gaus(0.,1.);
290  y = fRandom->Gaus(0.,1.);
291  z = fRandom->Gaus(0.,1.);
292  fHist3D[i]->Fill(x, y, z);
293  }
294  }
295  }
296 
297  return kTRUE;
298 }
299 
300 ////////////////////////////////////////////////////////////////////////////////
301 /// The SlaveTerminate() function is called after all entries or objects
302 /// have been processed. When running with PROOF SlaveTerminate() is called
303 /// on each slave server.
304 
306 {
307 }
308 
309 ////////////////////////////////////////////////////////////////////////////////
310 /// The Terminate() function is the last function to be called during
311 /// a query. It always runs on the client, it can be used to present
312 /// the results graphically or save the results to file.
313 
315 {
316  //
317  // Create a canvas, with 100 pads
318  //
319 
320  if (!fDraw || gROOT->IsBatch()){
321  return;
322  }
323 
325  fCHist1D=dynamic_cast<TCanvas*>(gROOT->FindObject("CHist1D"));
326  if (!fCHist1D){
327  fCHist1D = new TCanvas("CHist1D","Proof TSelHist Canvas (1D)", 200, 10,
328  700,700);
330  nside = (nside*nside < fNHists) ? nside+1 : nside;
331  fCHist1D->Divide(nside,nside,0,0);
332  }
333 
334  for (Int_t i=0; i < fNHists; i++) {
335  fHist1D[i] = dynamic_cast<TH1F *>
336  (fOutput->FindObject(Form("h1d_%d",i)));
337  fCHist1D->cd(i+1);
338  if (fHist1D[i]) fHist1D[i]->Draw();
339  }
340  // Final update
341  fCHist1D->cd();
342  fCHist1D->Update();
343  }
345  fCHist2D=dynamic_cast<TCanvas*>(gROOT->FindObject("CHist2D"));
346  if (!fCHist2D){
347  fCHist2D = new TCanvas("CHist2D","Proof TSelHist Canvas (2D)", 200, 10,
348  700,700);
350  nside = (nside*nside < fNHists) ? nside+1 : nside;
351  fCHist2D->Divide(nside,nside,0,0);
352  }
353  for (Int_t i=0; i < fNHists; i++) {
354  fHist2D[i] = dynamic_cast<TH2F *>
355  (fOutput->FindObject(Form("h2d_%d",i)));
356  fCHist2D->cd(i+1);
357  if (fHist2D[i]) fHist2D[i]->Draw("SURF");
358  }
359  // Final update
360  fCHist2D->cd();
361  fCHist2D->Update();
362  }
363 
365  fCHist3D=dynamic_cast<TCanvas*>(gROOT->FindObject("CHist3D"));
366  if (!fCHist3D){
367  fCHist3D = new TCanvas("CHist3D","Proof TSelHist Canvas (3D)", 200, 10,
368  700,700);
370  nside = (nside*nside < fNHists) ? nside+1 : nside;
371  fCHist3D->Divide(nside,nside,0,0);
372  }
373 
374  fOutput->Print("a");
375  for (Int_t i=0; i < fNHists; i++) {
376  fHist3D[i] = dynamic_cast<TH3F *>
377  (fOutput->FindObject(Form("h3d_%d",i)));
378  fCHist3D->cd(i+1);
379  if (fHist3D[i]) printf("fHist3D[%d] found\n", i);
380  if (fHist3D[i]) fHist3D[i]->Draw();
381  }
382  // Final update
383  fCHist3D->cd();
384  fCHist3D->Update();
385  }
386 
387 }
virtual const char * GetOption() const
Definition: TSelector.h:65
virtual void Begin(TTree *tree)
The Begin() function is called at the start of the query.
Definition: TSelHist.cxx:83
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
TCanvas * fCHist3D
Definition: TSelHist.h:47
Random number generator class based on M.
Definition: TRandom3.h:29
TSelectorList * fOutput
Definition: TSelector.h:50
long long Long64_t
Definition: RtypesCore.h:69
float Float_t
Definition: RtypesCore.h:53
3-D histogram with a float per channel (see TH1 documentation)}
Definition: TH3.h:272
THist< 2, float > TH2F
Definition: THist.h:321
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
Definition: Rtypes.h:61
TObject * FindObject(const char *name) const
Find object using its name.
Definition: THashList.cxx:212
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
THist< 1, float > TH1F
Definition: THist.h:315
#define gROOT
Definition: TROOT.h:340
TCanvas * fCHist1D
Definition: TSelHist.h:45
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void SlaveTerminate()
The SlaveTerminate() function is called after all entries or objects have been processed.
Definition: TSelHist.cxx:305
THist< 3, float > TH3F
Definition: THist.h:327
#define SafeDelete(p)
Definition: RConfig.h:436
Double_t x[n]
Definition: legend1.C:17
Bool_t fDraw
Definition: TSelHist.h:40
virtual void SlaveBegin(TTree *tree)
The SlaveBegin() function is called after the Begin() function.
Definition: TSelHist.cxx:155
TCanvas * fCHist2D
Definition: TSelHist.h:46
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Named parameter, streamable and storable.
Definition: TParameter.h:49
virtual void Terminate()
The Terminate() function is the last function to be called during a query.
Definition: TSelHist.cxx:314
EHistType GetType() const
TH3F ** fHist3D
Definition: TSelHist.h:43
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
ClassImp(TSelHist) TSelHist
Constructor.
Definition: TSelHist.cxx:39
TH2F ** fHist2D
Definition: TSelHist.h:42
char * Form(const char *fmt,...)
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH3.cxx:275
virtual void Print(Option_t *option="") const
Default print for collections, calls Print(option, 1).
virtual Bool_t Process(Long64_t entry)
The Process() function is called for each entry in the tree (or possibly keyed object in the case of ...
Definition: TSelHist.cxx:266
The Canvas class.
Definition: TCanvas.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:415
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Double_t y[n]
Definition: legend1.C:17
virtual ~TSelHist()
Destructor.
Definition: TSelHist.cxx:53
Mother of all ROOT objects.
Definition: TObject.h:58
Int_t fNHists
Definition: TSelHist.h:39
TList * fInput
Current object if processing object (vs. TTree)
Definition: TSelector.h:49
TPBHistType * fHistType
Definition: TSelHist.h:38
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1077
virtual void Add(TObject *obj)
Definition: TList.h:81
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
TH1F ** fHist1D
Definition: TSelHist.h:41
A TTree object has a header with a name and a title.
Definition: TTree.h:94
const AParamType & GetVal() const
Definition: TParameter.h:77
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2179
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:285
TObject * obj
TRandom3 * fRandom
Definition: TSelHist.h:44
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904