Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooAbsRealLValue.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofitcore:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 /**
18 \file RooAbsRealLValue.cxx
19 \class RooAbsRealLValue
20 \ingroup Roofitcore
21 
22 RooAbsRealLValue is the common abstract base class for objects that represent a
23 real value that may appear on the left hand side of an equation ('lvalue').
24 Each implementation must provide a setVal() member to allow direct modification
25 of the value. RooAbsRealLValue may be derived, but its functional relation
26 to other RooAbsArg must be invertible
27 
28 This class has methods that export the defined range of the lvalue,
29 but doesn't hold its values because these limits may be derived
30 from limits of client object. The range serve as integration
31 range when interpreted as a observable and a boundaries when
32 interpreted as a parameter.
33 **/
34 
35 #include "RooFit.h"
36 
37 #include <math.h>
38 #include "Riostream.h"
39 #include "TObjString.h"
40 #include "TTree.h"
41 #include "TH1.h"
42 #include "TH2.h"
43 #include "TH3.h"
44 #include "RooAbsRealLValue.h"
45 #include "RooStreamParser.h"
46 #include "RooRandom.h"
47 #include "RooPlot.h"
48 #include "RooArgList.h"
49 #include "RooAbsBinning.h"
50 #include "RooBinning.h"
51 #include "RooUniformBinning.h"
52 #include "RooCmdConfig.h"
53 #include "RooTreeData.h"
54 #include "RooRealVar.h"
55 #include "RooMsgService.h"
56 
57 
58 
59 using namespace std;
60 
61 ClassImp(RooAbsRealLValue);
62 
63 ////////////////////////////////////////////////////////////////////////////////
64 /// Constructor
65 
66 RooAbsRealLValue::RooAbsRealLValue(const char *name, const char *title, const char *unit) :
67  RooAbsReal(name, title, 0, 0, unit)
68 {
69 }
70 
71 
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 /// Copy constructor
75 
76 RooAbsRealLValue::RooAbsRealLValue(const RooAbsRealLValue& other, const char* name) :
77  RooAbsReal(other,name), RooAbsLValue(other)
78 {
79 }
80 
81 
82 ////////////////////////////////////////////////////////////////////////////////
83 /// Destructor
84 
85 RooAbsRealLValue::~RooAbsRealLValue()
86 {
87 }
88 
89 
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Return kTRUE if the input value is within our fit range. Otherwise, return
93 /// kFALSE and write a clipped value into clippedValPtr if it is non-zero.
94 
95 Bool_t RooAbsRealLValue::inRange(Double_t value, const char* rangeName, Double_t* clippedValPtr) const
96 {
97  // Double_t range = getMax() - getMin() ; // ok for +/-INIFINITY
98  Double_t clippedValue(value);
99  Bool_t isInRange(kTRUE) ;
100 
101  const RooAbsBinning& binning = getBinning(rangeName) ;
102  Double_t min = binning.lowBound() ;
103  Double_t max = binning.highBound() ;
104 
105  // test this value against our upper fit limit
106  if(!RooNumber::isInfinite(max) && value > (max+1e-6)) {
107  if (clippedValPtr) {
108 // coutW(InputArguments) << "RooAbsRealLValue::inFitRange(" << GetName() << "): value " << value
109 // << " rounded down to max limit " << getMax(rangeName) << endl ;
110  }
111  clippedValue = max;
112  isInRange = kFALSE ;
113  }
114  // test this value against our lower fit limit
115  if(!RooNumber::isInfinite(min) && value < min-1e-6) {
116  if (clippedValPtr) {
117 // coutW(InputArguments) << "RooAbsRealLValue::inFitRange(" << GetName() << "): value " << value
118 // << " rounded up to min limit " << getMin(rangeName) << endl;
119  }
120  clippedValue = min ;
121  isInRange = kFALSE ;
122  }
123 
124  if (clippedValPtr) *clippedValPtr=clippedValue ;
125 
126  return isInRange ;
127 }
128 
129 
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// Check if given value is valid
133 
134 Bool_t RooAbsRealLValue::isValidReal(Double_t value, Bool_t verbose) const
135 {
136  if (!inRange(value,0)) {
137  if (verbose)
138  coutI(InputArguments) << "RooRealVar::isValid(" << GetName() << "): value " << value
139  << " out of range (" << getMin() << " - " << getMax() << ")" << endl ;
140  return kFALSE ;
141  }
142  return kTRUE ;
143 }
144 
145 
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 /// Read object contents from given stream
149 
150 Bool_t RooAbsRealLValue::readFromStream(istream& /*is*/, Bool_t /*compact*/, Bool_t /*verbose*/)
151 {
152  return kTRUE ;
153 }
154 
155 
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// Write object contents to given stream
159 
160 void RooAbsRealLValue::writeToStream(ostream& /*os*/, Bool_t /*compact*/) const
161 {
162 }
163 
164 
165 
166 ////////////////////////////////////////////////////////////////////////////////
167 /// Assignment operator from a Double_t
168 
169 RooAbsArg& RooAbsRealLValue::operator=(Double_t newValue)
170 {
171  Double_t clipValue ;
172  // Clip
173  inRange(newValue,0,&clipValue) ;
174  setVal(clipValue) ;
175 
176  return *this ;
177 }
178 
179 
180 ////////////////////////////////////////////////////////////////////////////////
181 /// Assignment operator from other RooAbsReal
182 
183 RooAbsArg& RooAbsRealLValue::operator=(const RooAbsReal& arg)
184 {
185  return operator=(arg.getVal()) ;
186 }
187 
188 
189 
190 ////////////////////////////////////////////////////////////////////////////////
191 /// Create a new RooPlot on the heap with a drawing frame initialized for this
192 /// object, but no plot contents. Use x.frame() as the first argument to a
193 /// y.plotOn(...) method, for example. The caller is responsible for deleting
194 /// the returned object.
195 ///
196 /// <table>
197 /// <tr><th> Optional arguments <th>
198 /// <tr><td> Range(double lo, double hi) <td> Make plot frame for the specified range
199 /// <tr><td> Range(const char* name) <td> Make plot frame for range with the specified name
200 /// <tr><td> Bins(Int_t nbins) <td> Set default binning for datasets to specified number of bins
201 /// <tr><td> AutoRange(const RooAbsData& data, double margin) <td> Specifies range so that all points in given data set fit
202 /// inside the range with given margin.
203 /// <tr><td> AutoSymRange(const RooAbsData& data, double margin) <td> Specifies range so that all points in given data set fit
204 /// inside the range and center of range coincides with mean of distribution in given dataset.
205 /// <tr><td> Name(const char* name) <td> Give specified name to RooPlot object
206 /// <tr><td> Title(const char* title) <td> Give specified title to RooPlot object
207 /// </table>
208 ///
209 RooPlot* RooAbsRealLValue::frame(const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
210  const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
211 {
212  RooLinkedList cmdList ;
213  cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
214  cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
215  cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
216  cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
217 
218  return frame(cmdList) ;
219 }
220 
221 
222 
223 ////////////////////////////////////////////////////////////////////////////////
224 /// Back-end function for named argument frame() method
225 
226 RooPlot* RooAbsRealLValue::frame(const RooLinkedList& cmdList) const
227 {
228  // Define configuration for this method
229  RooCmdConfig pc(Form("RooAbsRealLValue::frame(%s)",GetName())) ;
230  pc.defineDouble("min","Range",0,getMin()) ;
231  pc.defineDouble("max","Range",1,getMax()) ;
232  pc.defineInt("nbins","Bins",0,getBins()) ;
233  pc.defineString("rangeName","RangeWithName",0,"") ;
234  pc.defineString("name","Name",0,"") ;
235  pc.defineString("title","Title",0,"") ;
236  pc.defineMutex("Range","RangeWithName","AutoRange") ;
237  pc.defineObject("rangeData","AutoRange",0,0) ;
238  pc.defineDouble("rangeMargin","AutoRange",0,0.1) ;
239  pc.defineInt("rangeSym","AutoRange",0,0) ;
240 
241  // Process & check varargs
242  pc.process(cmdList) ;
243  if (!pc.ok(kTRUE)) {
244  return 0 ;
245  }
246 
247  // Extract values from named arguments
248  Double_t xmin(getMin()),xmax(getMax()) ;
249  if (pc.hasProcessed("Range")) {
250  xmin = pc.getDouble("min") ;
251  xmax = pc.getDouble("max") ;
252  if (xmin==xmax) {
253  xmin = getMin() ;
254  xmax = getMax() ;
255  }
256  } else if (pc.hasProcessed("RangeWithName")) {
257  const char* rangeName=pc.getString("rangeName",0,kTRUE) ;
258  xmin = getMin(rangeName) ;
259  xmax = getMax(rangeName) ;
260  } else if (pc.hasProcessed("AutoRange")) {
261  RooTreeData* rangeData = static_cast<RooTreeData*>(pc.getObject("rangeData")) ;
262  rangeData->getRange((RooRealVar&)*this,xmin,xmax) ;
263  if (pc.getInt("rangeSym")==0) {
264  // Regular mode: range is from xmin to xmax with given extra margin
265  Double_t margin = pc.getDouble("rangeMargin")*(xmax-xmin) ;
266  xmin -= margin ;
267  xmax += margin ;
268  if (xmin<getMin()) xmin = getMin() ;
269  if (xmin>getMax()) xmax = getMax() ;
270  } else {
271  // Symmetric mode: range is centered at mean of distribution with enough width to include
272  // both lowest and highest point with margin
273  Double_t dmean = rangeData->moment((RooRealVar&)*this,1) ;
274  Double_t ddelta = ((xmax-dmean)>(dmean-xmin)?(xmax-dmean):(dmean-xmin))*(1+pc.getDouble("rangeMargin")) ;
275  xmin = dmean-ddelta ;
276  xmax = dmean+ddelta ;
277  if (xmin<getMin()) xmin = getMin() ;
278  if (xmin>getMax()) xmax = getMax() ;
279  }
280  } else {
281  xmin = getMin() ;
282  xmax = getMax() ;
283  }
284 
285  Int_t nbins = pc.getInt("nbins") ;
286  const char* name = pc.getString("name",0,kTRUE) ;
287  const char* title = pc.getString("title",0,kTRUE) ;
288 
289  RooPlot* theFrame = new RooPlot(*this,xmin,xmax,nbins) ;
290 
291  if (name) {
292  theFrame->SetName(name) ;
293  }
294  if (title) {
295  theFrame->SetTitle(title) ;
296  }
297 
298  return theFrame ;
299 }
300 
301 
302 
303 ////////////////////////////////////////////////////////////////////////////////
304 /// Create a new RooPlot on the heap with a drawing frame initialized for this
305 /// object, but no plot contents. Use x.frame() as the first argument to a
306 /// y.plotOn(...) method, for example. The caller is responsible for deleting
307 /// the returned object.
308 
309 RooPlot *RooAbsRealLValue::frame(Double_t xlo, Double_t xhi, Int_t nbins) const
310 {
311  return new RooPlot(*this,xlo,xhi,nbins);
312 }
313 
314 
315 
316 ////////////////////////////////////////////////////////////////////////////////
317 /// Create a new RooPlot on the heap with a drawing frame initialized for this
318 /// object, but no plot contents. Use x.frame() as the first argument to a
319 /// y.plotOn(...) method, for example. The caller is responsible for deleting
320 /// the returned object.
321 
322 RooPlot *RooAbsRealLValue::frame(Double_t xlo, Double_t xhi) const
323 {
324  return new RooPlot(*this,xlo,xhi,getBins());
325 }
326 
327 
328 
329 ////////////////////////////////////////////////////////////////////////////////
330 /// Create a new RooPlot on the heap with a drawing frame initialized for this
331 /// object, but no plot contents. Use x.frame() as the first argument to a
332 /// y.plotOn(...) method, for example. The caller is responsible for deleting
333 /// the returned object.
334 ///
335 /// The current fit range may not be open ended or empty.
336 
337 RooPlot *RooAbsRealLValue::frame(Int_t nbins) const
338 {
339  // Plot range of variable may not be infinite or empty
340  if (getMin()==getMax()) {
341  coutE(InputArguments) << "RooAbsRealLValue::frame(" << GetName() << ") ERROR: empty fit range, must specify plot range" << endl ;
342  return 0 ;
343  }
344  if (RooNumber::isInfinite(getMin())||RooNumber::isInfinite(getMax())) {
345  coutE(InputArguments) << "RooAbsRealLValue::frame(" << GetName() << ") ERROR: open ended fit range, must specify plot range" << endl ;
346  return 0 ;
347  }
348 
349  return new RooPlot(*this,getMin(),getMax(),nbins);
350 }
351 
352 
353 
354 ////////////////////////////////////////////////////////////////////////////////
355 /// Create a new RooPlot on the heap with a drawing frame initialized for this
356 /// object, but no plot contents. Use x.frame() as the first argument to a
357 /// y.plotOn(...) method, for example. The caller is responsible for deleting
358 /// the returned object.
359 ///
360 /// The current fit range may not be open ended or empty.
361 
362 RooPlot *RooAbsRealLValue::frame() const
363 {
364  // Plot range of variable may not be infinite or empty
365  if (getMin()==getMax()) {
366  coutE(InputArguments) << "RooAbsRealLValue::frame(" << GetName() << ") ERROR: empty fit range, must specify plot range" << endl ;
367  return 0 ;
368  }
369  if (RooNumber::isInfinite(getMin())||RooNumber::isInfinite(getMax())) {
370  coutE(InputArguments) << "RooAbsRealLValue::frame(" << GetName() << ") ERROR: open ended fit range, must specify plot range" << endl ;
371  return 0 ;
372  }
373 
374  return new RooPlot(*this,getMin(),getMax(),getBins());
375 }
376 
377 
378 
379 ////////////////////////////////////////////////////////////////////////////////
380 /// Copy cache of another RooAbsArg to our cache
381 
382 void RooAbsRealLValue::copyCache(const RooAbsArg* source, Bool_t valueOnly, Bool_t setValDirty)
383 {
384  RooAbsReal::copyCache(source,valueOnly,setValDirty) ;
385  setVal(_value) ; // force back-propagation
386 }
387 
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// Structure printing
391 
392 void RooAbsRealLValue::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
393 {
394  RooAbsReal::printMultiline(os,contents,verbose,indent);
395  os << indent << "--- RooAbsRealLValue ---" << endl;
396  TString unit(_unit);
397  if(!unit.IsNull()) unit.Prepend(' ');
398  os << indent << " Fit range is [ ";
399  if(hasMin()) {
400  os << getMin() << unit << " , ";
401  }
402  else {
403  os << "-INF , ";
404  }
405  if(hasMax()) {
406  os << getMax() << unit << " ]" << endl;
407  }
408  else {
409  os << "+INF ]" << endl;
410  }
411 }
412 
413 
414 
415 ////////////////////////////////////////////////////////////////////////////////
416 /// Set a new value sampled from a uniform distribution over the fit range.
417 /// Prints a warning and does nothing if the fit range is not finite.
418 
419 void RooAbsRealLValue::randomize(const char* rangeName)
420 {
421  RooAbsBinning& binning = getBinning(rangeName) ;
422  Double_t min = binning.lowBound() ;
423  Double_t max = binning.highBound() ;
424 
425  if(!RooNumber::isInfinite(min) && !RooNumber::isInfinite(max)) {
426  setValFast(min + RooRandom::uniform()*(max-min));
427  }
428  else {
429  coutE(Generation) << fName << "::" << ClassName() << ":randomize: fails with unbounded fit range" << endl;
430  }
431 }
432 
433 
434 
435 ////////////////////////////////////////////////////////////////////////////////
436 /// Set value to center of bin 'ibin' of binning 'rangeName' (or of
437 /// default binning if no range is specified)
438 
439 void RooAbsRealLValue::setBin(Int_t ibin, const char* rangeName)
440 {
441  // Check range of plot bin index
442  if (ibin<0 || ibin>=numBins(rangeName)) {
443  coutE(InputArguments) << "RooAbsRealLValue::setBin(" << GetName() << ") ERROR: bin index " << ibin
444  << " is out of range (0," << getBins(rangeName)-1 << ")" << endl ;
445  return ;
446  }
447 
448  // Set value to center of requested bin
449  setVal(getBinning(rangeName).binCenter(ibin)) ;
450 }
451 
452 
453 
454 
455 
456 ////////////////////////////////////////////////////////////////////////////////
457 /// Set value to center of bin 'ibin' of binning 'binning'
458 
459 void RooAbsRealLValue::setBin(Int_t ibin, const RooAbsBinning& binning)
460 {
461  // Set value to center of requested bin
462  setVal(binning.binCenter(ibin)) ;
463 }
464 
465 
466 
467 
468 
469 ////////////////////////////////////////////////////////////////////////////////
470 /// Set a new value sampled from a uniform distribution over the fit range.
471 /// Prints a warning and does nothing if the fit range is not finite.
472 
473 void RooAbsRealLValue::randomize(const RooAbsBinning& binning)
474 {
475  Double_t range= binning.highBound() - binning.lowBound() ;
476  setVal(binning.lowBound() + RooRandom::uniform()*range);
477 }
478 
479 
480 
481 
482 
483 ////////////////////////////////////////////////////////////////////////////////
484 /// Set value to center of bin 'ibin' of binning 'rangeName' (or of
485 /// default binning if no range is specified)
486 
487 void RooAbsRealLValue::setBinFast(Int_t ibin, const RooAbsBinning& binning)
488 {
489  // Set value to center of requested bin
490  setValFast(binning.binCenter(ibin)) ;
491 }
492 
493 
494 
495 ////////////////////////////////////////////////////////////////////////////////
496 /// Check if fit range is usable as plot range, i.e. it is neither
497 /// open ended, nor empty
498 
499 Bool_t RooAbsRealLValue::fitRangeOKForPlotting() const
500 {
501  return (hasMin() && hasMax() && (getMin()!=getMax())) ;
502 }
503 
504 
505 
506 ////////////////////////////////////////////////////////////////////////////////
507 /// Check if current value is inside range with given name
508 
509 Bool_t RooAbsRealLValue::inRange(const char* name) const
510 {
511  Double_t val = getVal() ;
512  Double_t epsilon = 1e-8 * fabs(val) ;
513  return (val >= getMin(name)-epsilon && val <= getMax(name)+epsilon) ;
514 }
515 
516 
517 
518 ////////////////////////////////////////////////////////////////////////////////
519 
520 TH1* RooAbsRealLValue::createHistogram(const char *name, const RooCmdArg& arg1, const RooCmdArg& arg2,
521  const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
522  const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
523 
524  // Create an empty ROOT histogram TH1,TH2 or TH3 suitabe to store information represent by the RooAbsRealLValue
525  //
526  // This function accepts the following arguments
527  //
528  // name -- Name of the ROOT histogram
529  //
530  // Binning(const char* name) -- Apply binning with given name to x axis of histogram
531  // Binning(RooAbsBinning& binning) -- Apply specified binning to x axis of histogram
532  // Binning(int_t nbins) -- Apply specified binning to x axis of histogram
533  // Binning(int_t nbins, double lo, double hi) -- Apply specified binning to x axis of histogram
534  // ConditionalObservables(const RooArgSet& set) -- Do not normalized PDF over following observables when projecting PDF into histogram
535  //
536  // YVar(const RooAbsRealLValue& var,...) -- Observable to be mapped on y axis of ROOT histogram
537  // ZVar(const RooAbsRealLValue& var,...) -- Observable to be mapped on z axis of ROOT histogram
538  //
539  // The YVar() and ZVar() arguments can be supplied with optional Binning() arguments to control the binning of the Y and Z axes, e.g.
540  // createHistogram("histo",x,Binning(-1,1,20), YVar(y,Binning(-1,1,30)), ZVar(z,Binning("zbinning")))
541  //
542  // The caller takes ownership of the returned histogram
543 {
544  RooLinkedList l ;
545  l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
546  l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
547  l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
548  l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
549 
550  return createHistogram(name,l) ;
551 }
552 
553 
554 
555 ////////////////////////////////////////////////////////////////////////////////
556 /// Create empty 1,2 or 3D histogram
557 /// Arguments recognized
558 ///
559 /// YVar() -- RooRealVar defining Y dimension with optional range/binning
560 /// ZVar() -- RooRealVar defining Z dimension with optional range/binning
561 /// AxisLabel() -- Vertical axis label
562 /// Binning() -- Range/Binning specification of X axis
563 
564 TH1* RooAbsRealLValue::createHistogram(const char *name, const RooLinkedList& cmdList) const
565 {
566  // Define configuration for this method
567  RooCmdConfig pc(Form("RooAbsRealLValue::createHistogram(%s)",GetName())) ;
568 
569  pc.defineObject("xbinning","Binning",0,0) ;
570  pc.defineString("xbinningName","BinningName",0,"") ;
571  pc.defineInt("nxbins","BinningSpec",0) ;
572  pc.defineDouble("xlo","BinningSpec",0,0) ;
573  pc.defineDouble("xhi","BinningSpec",1,0) ;
574 
575  pc.defineObject("yvar","YVar",0,0) ;
576  pc.defineObject("ybinning","YVar::Binning",0,0) ;
577  pc.defineString("ybinningName","YVar::BinningName",0,"") ;
578  pc.defineInt("nybins","YVar::BinningSpec",0) ;
579  pc.defineDouble("ylo","YVar::BinningSpec",0,0) ;
580  pc.defineDouble("yhi","YVar::BinningSpec",1,0) ;
581 
582  pc.defineObject("zvar","ZVar",0,0) ;
583  pc.defineObject("zbinning","ZVar::Binning",0,0) ;
584  pc.defineString("zbinningName","ZVar::BinningName",0,"") ;
585  pc.defineInt("nzbins","ZVar::BinningSpec",0) ;
586  pc.defineDouble("zlo","ZVar::BinningSpec",0,0) ;
587  pc.defineDouble("zhi","ZVar::BinningSpec",1,0) ;
588 
589  pc.defineString("axisLabel","AxisLabel",0,"Events") ;
590 
591  pc.defineDependency("ZVar","YVar") ;
592 
593  // Process & check varargs
594  pc.process(cmdList) ;
595  if (!pc.ok(kTRUE)) {
596  return 0 ;
597  }
598 
599  // Initialize arrays for call to implementation version of createHistogram
600  const char* axisLabel = pc.getString("axisLabel") ;
601  const RooAbsBinning* binning[3] ;
602  Bool_t ownBinning[3] = { kFALSE, kFALSE, kFALSE } ;
603  RooArgList vars ;
604 
605  // Prepare X dimension
606  vars.add(*this) ;
607  if (pc.hasProcessed("Binning")) {
608  binning[0] = static_cast<RooAbsBinning*>(pc.getObject("xbinning")) ;
609  } else if (pc.hasProcessed("BinningName")) {
610  binning[0] = &getBinning(pc.getString("xbinningName",0,kTRUE)) ;
611  } else if (pc.hasProcessed("BinningSpec")) {
612  Double_t xlo = pc.getDouble("xlo") ;
613  Double_t xhi = pc.getDouble("xhi") ;
614  binning[0] = new RooUniformBinning((xlo==xhi)?getMin():xlo,(xlo==xhi)?getMax():xhi,pc.getInt("nxbins")) ;
615  ownBinning[0] = kTRUE ;
616  } else {
617  binning[0] = &getBinning() ;
618  }
619 
620  if (pc.hasProcessed("YVar")) {
621  RooAbsRealLValue& yvar = *static_cast<RooAbsRealLValue*>(pc.getObject("yvar")) ;
622  vars.add(yvar) ;
623  if (pc.hasProcessed("YVar::Binning")) {
624  binning[1] = static_cast<RooAbsBinning*>(pc.getObject("ybinning")) ;
625  } else if (pc.hasProcessed("YVar::BinningName")) {
626  binning[1] = &yvar.getBinning(pc.getString("ybinningName",0,kTRUE)) ;
627  } else if (pc.hasProcessed("YVar::BinningSpec")) {
628  Double_t ylo = pc.getDouble("ylo") ;
629  Double_t yhi = pc.getDouble("yhi") ;
630  binning[1] = new RooUniformBinning((ylo==yhi)?yvar.getMin():ylo,(ylo==yhi)?yvar.getMax():yhi,pc.getInt("nybins")) ;
631  ownBinning[1] = kTRUE ;
632  } else {
633  binning[1] = &yvar.getBinning() ;
634  }
635  }
636 
637  if (pc.hasProcessed("ZVar")) {
638  RooAbsRealLValue& zvar = *static_cast<RooAbsRealLValue*>(pc.getObject("zvar")) ;
639  vars.add(zvar) ;
640  if (pc.hasProcessed("ZVar::Binning")) {
641  binning[2] = static_cast<RooAbsBinning*>(pc.getObject("zbinning")) ;
642  } else if (pc.hasProcessed("ZVar::BinningName")) {
643  binning[2] = &zvar.getBinning(pc.getString("zbinningName",0,kTRUE)) ;
644  } else if (pc.hasProcessed("ZVar::BinningSpec")) {
645  Double_t zlo = pc.getDouble("zlo") ;
646  Double_t zhi = pc.getDouble("zhi") ;
647  binning[2] = new RooUniformBinning((zlo==zhi)?zvar.getMin():zlo,(zlo==zhi)?zvar.getMax():zhi,pc.getInt("nzbins")) ;
648  ownBinning[2] = kTRUE ;
649  } else {
650  binning[2] = &zvar.getBinning() ;
651  }
652  }
653 
654 
655  TH1* ret = createHistogram(name, vars, axisLabel, binning) ;
656 
657  if (ownBinning[0]) delete binning[0] ;
658  if (ownBinning[1]) delete binning[1] ;
659  if (ownBinning[2]) delete binning[2] ;
660 
661  return ret ;
662 }
663 
664 
665 
666 ////////////////////////////////////////////////////////////////////////////////
667 /// Create an empty 1D-histogram with appropriate scale and labels for this variable.
668 /// This method uses the default plot range which can be changed using the
669 /// setPlotMin(),setPlotMax() methods, and the default binning which can be
670 /// changed with setPlotBins(). The caller takes ownership of the returned
671 /// object and is responsible for deleting it.
672 
673 TH1F *RooAbsRealLValue::createHistogram(const char *name, const char *yAxisLabel) const
674 {
675  // Check if the fit range is usable as plot range
676  if (!fitRangeOKForPlotting()) {
677  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
678  << ") ERROR: fit range empty or open ended, must explicitly specify range" << endl ;
679  return 0 ;
680  }
681 
682  RooArgList list(*this) ;
683  Double_t xlo = getMin() ;
684  Double_t xhi = getMax() ;
685  Int_t nbins = getBins() ;
686 
687  // coverity[ARRAY_VS_SINGLETON]
688  return (TH1F*)createHistogram(name, list, yAxisLabel, &xlo, &xhi, &nbins);
689 }
690 
691 
692 
693 ////////////////////////////////////////////////////////////////////////////////
694 /// Create an empty 1D-histogram with appropriate scale and labels for this variable.
695 /// This method uses the default plot range which can be changed using the
696 /// setPlotMin(),setPlotMax() methods, and the default binning which can be
697 /// changed with setPlotBins(). The caller takes ownership of the returned
698 /// object and is responsible for deleting it.
699 
700 TH1F *RooAbsRealLValue::createHistogram(const char *name, const char *yAxisLabel, Double_t xlo, Double_t xhi, Int_t nBins) const
701 {
702  RooArgList list(*this) ;
703 
704  // coverity[ARRAY_VS_SINGLETON]
705  return (TH1F*)createHistogram(name, list, yAxisLabel, &xlo, &xhi, &nBins);
706 }
707 
708 
709 
710 ////////////////////////////////////////////////////////////////////////////////
711 /// Create an empty 1D-histogram with appropriate scale and labels for this variable.
712 
713 TH1F *RooAbsRealLValue::createHistogram(const char *name, const char *yAxisLabel, const RooAbsBinning& bins) const
714 {
715  RooArgList list(*this) ;
716  const RooAbsBinning* pbins = &bins ;
717 
718  // coverity[ARRAY_VS_SINGLETON]
719  return (TH1F*)createHistogram(name, list, yAxisLabel, &pbins);
720 }
721 
722 
723 
724 ////////////////////////////////////////////////////////////////////////////////
725 /// Create an empty 2D-histogram with appropriate scale and labels for this variable (x)
726 /// and the specified y variable. This method uses the default plot ranges for x and y which
727 /// can be changed using the setPlotMin(),setPlotMax() methods, and the default binning which
728 /// can be changed with setPlotBins(). The caller takes ownership of the returned object
729 /// and is responsible for deleting it.
730 
731 TH2F *RooAbsRealLValue::createHistogram(const char *name, const RooAbsRealLValue &yvar, const char *zAxisLabel,
732  Double_t* xlo, Double_t* xhi, Int_t* nBins) const
733 {
734  if ((!xlo && xhi) || (xlo && !xhi)) {
735  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
736  << ") ERROR must specify either no range, or both limits" << endl ;
737  return 0 ;
738  }
739 
740  Double_t xlo_fit[2] ;
741  Double_t xhi_fit[2] ;
742  Int_t nbins_fit[2] ;
743 
744  Double_t *xlo2 = xlo;
745  Double_t *xhi2 = xhi;
746  Int_t *nBins2 = nBins;
747 
748  if (!xlo2) {
749 
750  if (!fitRangeOKForPlotting()) {
751  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
752  << ") ERROR: fit range empty or open ended, must explicitly specify range" << endl ;
753  return 0 ;
754  }
755  if (!yvar.fitRangeOKForPlotting()) {
756  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
757  << ") ERROR: fit range of " << yvar.GetName() << " empty or open ended, must explicitly specify range" << endl ;
758  return 0 ;
759  }
760 
761  xlo_fit[0] = getMin() ;
762  xhi_fit[0] = getMax() ;
763 
764  xlo_fit[1] = yvar.getMin() ;
765  xhi_fit[1] = yvar.getMax() ;
766 
767  xlo2 = xlo_fit ;
768  xhi2 = xhi_fit ;
769  }
770 
771  if (!nBins2) {
772  nbins_fit[0] = getBins() ;
773  nbins_fit[1] = yvar.getBins() ;
774  nBins2 = nbins_fit ;
775  }
776 
777 
778  RooArgList list(*this,yvar) ;
779  // coverity[OVERRUN_STATIC]
780  return (TH2F*)createHistogram(name, list, zAxisLabel, xlo2, xhi2, nBins2);
781 }
782 
783 
784 
785 ////////////////////////////////////////////////////////////////////////////////
786 /// Create an empty 2D-histogram with appropriate scale and labels for this variable (x)
787 /// and the specified y variable.
788 
789 TH2F *RooAbsRealLValue::createHistogram(const char *name, const RooAbsRealLValue &yvar,
790  const char *zAxisLabel, const RooAbsBinning** bins) const
791 {
792  RooArgList list(*this,yvar) ;
793  return (TH2F*)createHistogram(name, list, zAxisLabel, bins);
794 }
795 
796 
797 
798 ////////////////////////////////////////////////////////////////////////////////
799 /// Create an empty 3D-histogram with appropriate scale and labels for this variable (x)
800 /// and the specified y,z variables. This method uses the default plot ranges for x,y,z which
801 /// can be changed using the setPlotMin(),setPlotMax() methods, and the default binning which
802 /// can be changed with setPlotBins(). The caller takes ownership of the returned object
803 /// and is responsible for deleting it.
804 
805 TH3F *RooAbsRealLValue::createHistogram(const char *name, const RooAbsRealLValue &yvar, const RooAbsRealLValue &zvar,
806  const char *tAxisLabel, Double_t* xlo, Double_t* xhi, Int_t* nBins) const
807 {
808  if ((!xlo && xhi) || (xlo && !xhi)) {
809  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
810  << ") ERROR must specify either no range, or both limits" << endl ;
811  return 0 ;
812  }
813 
814  Double_t xlo_fit[3] ;
815  Double_t xhi_fit[3] ;
816  Int_t nbins_fit[3] ;
817 
818  Double_t *xlo2 = xlo;
819  Double_t *xhi2 = xhi;
820  Int_t* nBins2 = nBins;
821  if (!xlo2) {
822 
823  if (!fitRangeOKForPlotting()) {
824  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
825  << ") ERROR: fit range empty or open ended, must explicitly specify range" << endl ;
826  return 0 ;
827  }
828  if (!yvar.fitRangeOKForPlotting()) {
829  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
830  << ") ERROR: fit range of " << yvar.GetName() << " empty or open ended, must explicitly specify range" << endl ;
831  return 0 ;
832  }
833  if (!zvar.fitRangeOKForPlotting()) {
834  coutE(InputArguments) << "RooAbsRealLValue::createHistogram(" << GetName()
835  << ") ERROR: fit range of " << zvar.GetName() << " empty or open ended, must explicitly specify range" << endl ;
836  return 0 ;
837  }
838 
839  xlo_fit[0] = getMin() ;
840  xhi_fit[0] = getMax() ;
841 
842  xlo_fit[1] = yvar.getMin() ;
843  xhi_fit[1] = yvar.getMax() ;
844 
845  xlo_fit[2] = zvar.getMin() ;
846  xhi_fit[2] = zvar.getMax() ;
847 
848  xlo2 = xlo_fit ;
849  xhi2 = xhi_fit ;
850  }
851 
852  if (!nBins2) {
853  nbins_fit[0] = getBins() ;
854  nbins_fit[1] = yvar.getBins() ;
855  nbins_fit[2] = zvar.getBins() ;
856  nBins2 = nbins_fit ;
857  }
858 
859  RooArgList list(*this,yvar,zvar) ;
860  return (TH3F*)createHistogram(name, list, tAxisLabel, xlo2, xhi2, nBins2);
861 }
862 
863 
864 TH3F *RooAbsRealLValue::createHistogram(const char *name, const RooAbsRealLValue &yvar, const RooAbsRealLValue &zvar,
865  const char* tAxisLabel, const RooAbsBinning** bins) const
866 {
867  // Create an empty 3D-histogram with appropriate scale and labels for this variable (x)
868  // and the specified y,z variables.
869 
870  RooArgList list(*this,yvar,zvar) ;
871  return (TH3F*)createHistogram(name, list, tAxisLabel, bins);
872 }
873 
874 
875 
876 
877 ////////////////////////////////////////////////////////////////////////////////
878 /// Create 1-, 2- or 3-d ROOT histogram with labels taken
879 /// from the variables in 'vars' and the with range and binning
880 /// specified in xlo,xhi and nBins. The dimensions of the arrays xlo,xhi,
881 /// nBins should match the number of objects in vars.
882 
883 TH1 *RooAbsRealLValue::createHistogram(const char *name, RooArgList &vars, const char *tAxisLabel,
884  Double_t* xlo, Double_t* xhi, Int_t* nBins)
885 {
886  const RooAbsBinning* bin[3] ;
887  Int_t ndim = vars.getSize() ;
888  bin[0] = new RooUniformBinning(xlo[0],xhi[0],nBins[0]) ;
889  bin[1] = (ndim>1) ? new RooUniformBinning(xlo[1],xhi[1],nBins[1]) : 0 ;
890  bin[2] = (ndim>2) ? new RooUniformBinning(xlo[2],xhi[2],nBins[2]) : 0 ;
891 
892  TH1* ret = createHistogram(name,vars,tAxisLabel,bin) ;
893 
894  if (bin[0]) delete bin[0] ;
895  if (bin[1]) delete bin[1] ;
896  if (bin[2]) delete bin[2] ;
897  return ret ;
898 }
899 
900 
901 
902 ////////////////////////////////////////////////////////////////////////////////
903 /// Create a 1,2, or 3D-histogram with appropriate scale and labels.
904 /// Binning and ranges are taken from the variables themselves and can be changed by
905 /// calling their setPlotMin/Max() and setPlotBins() methods. A histogram can be filled
906 /// using RooAbsReal::fillHistogram() or RooTreeData::fillHistogram().
907 /// The caller takes ownership of the returned object and is responsible for deleting it.
908 
909 TH1 *RooAbsRealLValue::createHistogram(const char *name, RooArgList &vars, const char *tAxisLabel, const RooAbsBinning** bins)
910 {
911  // Check that we have 1-3 vars
912  Int_t dim= vars.getSize();
913  if(dim < 1 || dim > 3) {
914  oocoutE((TObject*)0,InputArguments) << "RooAbsReal::createHistogram: dimension not supported: " << dim << endl;
915  return 0;
916  }
917 
918  // Check that all variables are AbsReals and prepare a name of the form <name>_<var1>_...
919  TString histName(name);
920  histName.Append("_");
921  const RooAbsRealLValue *xyz[3];
922 
923  Int_t index;
924  for(index= 0; index < dim; index++) {
925  const RooAbsArg *arg= vars.at(index);
926  xyz[index]= dynamic_cast<const RooAbsRealLValue*>(arg);
927  if(!xyz[index]) {
928  oocoutE((TObject*)0,InputArguments) << "RooAbsRealLValue::createHistogram: variable is not real lvalue: " << arg->GetName() << endl;
929  return 0;
930  }
931  histName.Append("_");
932  histName.Append(arg->GetName());
933  }
934  TString histTitle(histName);
935  histTitle.Prepend("Histogram of ");
936 
937  // Create the histogram
938  TH1 *histogram = 0;
939  switch(dim) {
940  case 1:
941  if (bins[0]->isUniform()) {
942  histogram= new TH1F(histName.Data(), histTitle.Data(),
943  bins[0]->numBins(),bins[0]->lowBound(),bins[0]->highBound());
944  } else {
945  histogram= new TH1F(histName.Data(), histTitle.Data(),
946  bins[0]->numBins(),bins[0]->array());
947  }
948  break;
949  case 2:
950  if (bins[0]->isUniform() && bins[1]->isUniform()) {
951  histogram= new TH2F(histName.Data(), histTitle.Data(),
952  bins[0]->numBins(),bins[0]->lowBound(),bins[0]->highBound(),
953  bins[1]->numBins(),bins[1]->lowBound(),bins[1]->highBound());
954  } else {
955  histogram= new TH2F(histName.Data(), histTitle.Data(),
956  bins[0]->numBins(),bins[0]->array(),
957  bins[1]->numBins(),bins[1]->array());
958  }
959  break;
960  case 3:
961  if (bins[0]->isUniform() && bins[1]->isUniform() && bins[2]->isUniform()) {
962  histogram= new TH3F(histName.Data(), histTitle.Data(),
963  bins[0]->numBins(),bins[0]->lowBound(),bins[0]->highBound(),
964  bins[1]->numBins(),bins[1]->lowBound(),bins[1]->highBound(),
965  bins[2]->numBins(),bins[2]->lowBound(),bins[2]->highBound()) ;
966  } else {
967  histogram= new TH3F(histName.Data(), histTitle.Data(),
968  bins[0]->numBins(),bins[0]->array(),
969  bins[1]->numBins(),bins[1]->array(),
970  bins[2]->numBins(),bins[2]->array()) ;
971  }
972  break;
973  }
974  if(!histogram) {
975  oocoutE((TObject*)0,InputArguments) << "RooAbsReal::createHistogram: unable to create a new histogram" << endl;
976  return 0;
977  }
978 
979  // Set the histogram coordinate axis labels from the titles of each variable, adding units if necessary.
980  for(index= 0; index < dim; index++) {
981  TString axisTitle(xyz[index]->getTitle(kTRUE));
982  switch(index) {
983  case 0:
984  histogram->SetXTitle(axisTitle.Data());
985  break;
986  case 1:
987  histogram->SetYTitle(axisTitle.Data());
988  break;
989  case 2:
990  histogram->SetZTitle(axisTitle.Data());
991  break;
992  default:
993  assert(0);
994  break;
995  }
996  }
997 
998  // Set the t-axis title if given one
999  if((0 != tAxisLabel) && (0 != strlen(tAxisLabel))) {
1000  TString axisTitle(tAxisLabel);
1001  axisTitle.Append(" / ( ");
1002  for(Int_t index2= 0; index2 < dim; index2++) {
1003  Double_t delta= bins[index2]->averageBinWidth() ; // xyz[index2]->getBins();
1004  if(index2 > 0) axisTitle.Append(" x ");
1005  axisTitle.Append(Form("%g",delta));
1006  if(strlen(xyz[index2]->getUnit())) {
1007  axisTitle.Append(" ");
1008  axisTitle.Append(xyz[index2]->getUnit());
1009  }
1010  }
1011  axisTitle.Append(" )");
1012  switch(dim) {
1013  case 1:
1014  histogram->SetYTitle(axisTitle.Data());
1015  break;
1016  case 2:
1017  histogram->SetZTitle(axisTitle.Data());
1018  break;
1019  case 3:
1020  // not supported in TH1
1021  break;
1022  default:
1023  assert(0);
1024  break;
1025  }
1026  }
1027 
1028  return histogram;
1029 }
1030 
1031 
1032 Bool_t RooAbsRealLValue::isJacobianOK(const RooArgSet&) const
1033 {
1034  // Interface function to indicate that this lvalue
1035  // has a unit or constant jacobian terms with respect to
1036  // the observable passed as argument. This default implementation
1037  // always returns true (i.e. jacobian is constant)
1038  return kTRUE ;
1039 }