Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooGenContext.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 RooGenContext.cxx
19\class RooGenContext
20\ingroup Roofitcore
21
22Class RooGenContext implement a universal generator context for all
23RooAbsPdf classes that do not have or need a specialized generator
24context. This generator context queries the input p.d.f which observables
25it can generate internally and delegates generation of those observables
26to the p.d.f if it deems that safe. The other observables are generated
27use a RooAcceptReject sampling technique.
28**/
29
30#include "RooMsgService.h"
31#include "Riostream.h"
32
33#include "RooGenContext.h"
34#include "RooAbsPdf.h"
35#include "RooDataSet.h"
36#include "RooRealIntegral.h"
37#include "RooAcceptReject.h"
38#include "RooRealVar.h"
39#include "RooDataHist.h"
40#include "RooErrorHandler.h"
41#include "RooNumGenConfig.h"
42#include "RooNumGenFactory.h"
43
44#include "TString.h"
45
46
47using namespace std;
48
50 ;
51
52
53
54////////////////////////////////////////////////////////////////////////////////
55/// Initialize a new context for generating events with the specified
56/// variables, using the specified PDF model. A prototype dataset (if provided)
57/// is not cloned and still belongs to the caller. The contents and shape
58/// of this dataset can be changed between calls to generate() as long as the
59/// expected columns to be copied to the generated dataset are present.
60/// Any argument supplied in the forceDirect RooArgSet are always offered
61/// for internal generation to the p.d.f., even if this is deemed unsafe by
62/// the logic of RooGenContext.
63
65 const RooDataSet *prototype, const RooArgSet* auxProto,
66 bool verbose, const RooArgSet* forceDirect) :
67 RooAbsGenContext(model,vars,prototype,auxProto,verbose),
68 _pdfClone(0), _generator(0),
69 _maxVar(0), _updateFMaxPerEvent(0)
70{
71 cxcoutI(Generation) << "RooGenContext::ctor() setting up event generator context for p.d.f. " << model.GetName()
72 << " for generation of observable(s) " << vars ;
73 if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
74 if (auxProto && !auxProto->empty()) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
75 if (forceDirect && !forceDirect->empty()) ccxcoutI(Generation) << " with internal generation forced for observables " << *forceDirect ;
76 ccxcoutI(Generation) << endl ;
77
78
79 // Clone the model and all nodes that it depends on so that this context
80 // is independent of any existing objects.
81 RooArgSet nodes(model,model.GetName());
82 if (nodes.snapshot(_cloneSet, true)) {
83 coutE(Generation) << "RooGenContext::RooGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
85 }
86
87 // Find the clone in the snapshot list
88 _pdfClone = static_cast<RooAbsPdf*>(_cloneSet.find(model.GetName()));
90
91 // Optionally fix RooAddPdf normalizations
92 if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
93 RooArgSet fullNormSet(vars) ;
94 fullNormSet.add(*prototype->get()) ;
96 }
97
98 // Analyze the list of variables to generate...
99 _isValid= true;
100 const RooAbsArg *arg = 0;
101 for(RooAbsArg * tmp : vars) {
102 if(!_isValid) break;
103
104 // is this argument derived?
105 if(!tmp->isFundamental()) {
106 coutE(Generation) << "RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() << "\"" << endl;
107 _isValid= false;
108 continue;
109 }
110 // lookup this argument in the cloned set of PDF dependents
111 arg= static_cast<RooAbsArg const*>(_cloneSet.find(tmp->GetName()));
112 if(0 == arg) {
113 coutI(Generation) << "RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
114 << "\" which will have uniform distribution" << endl;
115 _uniformVars.add(*tmp);
116 }
117 else {
118
119 // does the model depend on this variable directly, ie, like "x" in
120 // f(x) or f(x,g(x,y)) or even f(x,x) ?
121 const RooAbsArg *direct= arg ;
122 if (forceDirect==0 || !forceDirect->find(direct->GetName())) {
123 if (!_pdfClone->isDirectGenSafe(*arg)) {
124 cxcoutD(Generation) << "RooGenContext::ctor() observable " << arg->GetName() << " has been determined to be unsafe for internal generation" << endl;
125 direct=0 ;
126 }
127 }
128
129 // does the model depend indirectly on this variable through an lvalue chain?
130 // otherwise, this variable will have to be generated with accept/reject
131 if(direct) {
132 _directVars.add(*arg);
133 } else {
134 _otherVars.add(*arg);
135 }
136 }
137 }
138
139 if(!isValid()) {
140 coutE(Generation) << "RooGenContext::ctor() constructor failed with errors" << endl;
141 return;
142 }
143
144 // At this point directVars are all variables that are safe to be generated directly
145 // otherVars are all variables that are _not_ safe to be generated directly
146 if (!_directVars.empty()) {
147 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _directVars << " are safe for internal generator (if supported by p.d.f)" << endl ;
148 }
149 if (!_otherVars.empty()) {
150 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _otherVars << " are NOT safe for internal generator (if supported by p.d.f)" << endl ;
151 }
152
153 // If PDF depends on prototype data, direct generator cannot use static initialization
154 // in initGenerator()
155 bool staticInitOK = !_pdfClone->dependsOn(_protoVars) ;
156 if (!staticInitOK) {
157 cxcoutD(Generation) << "RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << endl ;
158 }
159
160 // Can the model generate any of the direct variables itself?
161 RooArgSet generatedVars;
162 if (!_directVars.empty()) {
163 _code= _pdfClone->getGenerator(_directVars,generatedVars,staticInitOK);
164
165 cxcoutD(Generation) << "RooGenContext::ctor() Model indicates that it can internally generate observables "
166 << generatedVars << " with configuration identifier " << _code << endl ;
167 } else {
168 _code = 0 ;
169 }
170
171 // Move variables which cannot be generated into the list to be generated with accept/reject
172 _directVars.remove(generatedVars);
174
175 // Update _directVars to only include variables that will actually be directly generated
177 _directVars.add(generatedVars);
178
179 cxcoutI(Generation) << "RooGenContext::ctor() Context will" ;
180 if (!_directVars.empty()) ccxcoutI(Generation) << " let model internally generate observables " << _directVars ;
181 if (!_directVars.empty() && !_otherVars.empty()) ccxcoutI(Generation) << " and" ;
182 if (!_otherVars.empty()) ccxcoutI(Generation) << " generate variables " << _otherVars << " with accept/reject sampling" ;
183 if (!_uniformVars.empty()) ccxcoutI(Generation) << ", non-observable variables " << _uniformVars << " will be generated with flat distribution" ;
184 ccxcoutI(Generation) << endl ;
185
186 // initialize the accept-reject generator
187 RooArgSet depList;
189 depList.remove(_otherVars);
190
191 TString nname(_pdfClone->GetName()) ;
192 nname.Append("_AccRej") ;
193 TString ntitle(_pdfClone->GetTitle()) ;
194 ntitle.Append(" (Accept/Reject)") ;
195
196
197 RooArgSet protoDeps;
198 model.getObservables(&_protoVars, protoDeps);
199
200
201 if (_protoVars.empty()) {
202
203 // No prototype variables
204
205 if(depList.empty()) {
206 // All variable are generated with accept-reject
207
208 // Check if PDF supports maximum finding
209 Int_t maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
210 if (maxFindCode != 0) {
211 coutI(Generation) << "RooGenContext::ctor() no prototype data provided, all observables are generated with numerically and "
212 << "model supports analytical maximum findin:, can provide analytical pdf maximum to numeric generator" << endl ;
213 double maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(&_theEvent) ;
214 _maxVar = new RooRealVar("funcMax","function maximum",maxVal) ;
215 cxcoutD(Generation) << "RooGenContext::ctor() maximum value returned by RooAbsPdf::maxVal() is " << maxVal << endl ;
216 }
217 }
218
219 if (_otherVars.getSize()>0) {
220 _pdfClone->getVal(&vars) ; // WVE debug
221 _acceptRejectFunc= std::make_unique<RooRealIntegral>(nname,ntitle,*_pdfClone,depList,&vars);
222 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
223 } else {
224 _acceptRejectFunc = nullptr;
225 }
226
227 } else {
228
229 // Generation _with_ prototype variable
230 depList.remove(_protoVars,true,true);
231 std::unique_ptr<RooAbsReal> integ{_pdfClone->createIntegral(depList,vars)};
232 _acceptRejectFunc= std::unique_ptr<RooRealIntegral>{static_cast<RooRealIntegral*>(integ.release())};
233 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
234
235 // Check if PDF supports maximum finding for the entire phase space
236 RooArgSet allVars(_otherVars);
237 allVars.add(_directVars);
238 Int_t maxFindCode = _pdfClone->getMaxVal(allVars) ;
239 if (maxFindCode != 0) {
240 // Special case: PDF supports max-finding in otherVars, no need to scan other+proto space for maximum
241 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
242 << "model supports analytical maximum finding in the full phase space: "
243 << "can provide analytical pdf maximum to numeric generator" << endl ;
244 double maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(&allVars);
245 _maxVar = new RooRealVar("funcMax", "function maximum", maxVal);
246 } else {
247 maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
248 if (maxFindCode != 0) {
249 _updateFMaxPerEvent = maxFindCode ;
250 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
251 << "model supports analytical maximum finding in the variables which are not"
252 << " internally generated. Can provide analytical pdf maximum to numeric "
253 << "generator" << endl;
254 cxcoutD(Generation) << "RooGenContext::ctor() maximum value must be reevaluated for each "
255 << "event with configuration code " << maxFindCode << endl ;
256 _maxVar = new RooRealVar("funcMax","function maximum",1) ;
257 }
258 }
259
260 if (!_maxVar) {
261
262 // Regular case: First find maximum in other+proto space
263 RooArgSet otherAndProto(_otherVars) ;
264
265 otherAndProto.add(protoDeps) ;
266
267 if (_otherVars.getSize()>0) {
268
269 cxcoutD(Generation) << "RooGenContext::ctor() prototype data provided, observables are generated numerically no "
270 << "analytical estimate of maximum function value provided by model, must determine maximum value through initial sampling space "
271 << "of accept/reject observables plus prototype observables: " << otherAndProto << endl ;
272
273 // Calculate maximum in other+proto space if there are any accept/reject generated observables
275 *model.getGeneratorConfig(),_verbose) ;
276// RooAcceptReject maxFinder(*_acceptRejectFunc,otherAndProto,RooNumGenConfig::defaultConfig(),_verbose) ;
277 double max = maxFinder->getFuncMax() ;
278 _maxVar = new RooRealVar("funcMax","function maximum",max) ;
279
280 if (max==0) {
281 oocoutE(nullptr, Generation) << "RooGenContext::ctor(" << model.GetName()
282 << ") ERROR: generating conditional p.d.f. which requires prior knowledge of function maximum, "
283 << "but chosen numeric generator (" << maxFinder->generatorName() << ") does not support maximum finding" << endl ;
284 delete maxFinder ;
285 throw string("RooGenContext::ctor()") ;
286 }
287 delete maxFinder ;
288
289 cxcoutD(Generation) << "RooGenContext::ctor() maximum function value found through initial sampling is " << max << endl ;
290 }
291 }
292
293 }
294
297 cxcoutD(Generation) << "RooGenContext::ctor() creating MC sampling generator " << _generator->generatorName() << " from function for observables " << _otherVars << endl ;
298 //_generator= new RooAcceptReject(*_acceptRejectFunc,_otherVars,RooNumGenConfig::defaultConfig(),_verbose,_maxVar);
299 } else {
300 _generator = 0 ;
301 }
302
304}
305
306
307////////////////////////////////////////////////////////////////////////////////
308/// Destructor.
309
311{
312 // Clean up our accept/reject generator
313 if (_generator) delete _generator;
314 if (_maxVar) delete _maxVar ;
315}
316
317
318
319////////////////////////////////////////////////////////////////////////////////
320/// Attach the cloned model to the event buffer we will be filling.
321
323{
325 if (_acceptRejectFunc) {
326 _acceptRejectFunc->recursiveRedirectServers(args,false) ; // WVE DEBUG
327 }
328
329 // Attach the RooAcceptReject generator the event buffer
330 if (_generator) {
332 }
333
334}
335
336
337
338////////////////////////////////////////////////////////////////////////////////
339/// Perform one-time initialization of the generator context
340
342{
343 for (auto* arg : theEvent) {
344 arg->setOperMode(RooAbsArg::ADirty) ;
345 }
346
347 attach(theEvent) ;
348
349 // Reset the cloned model's error counters.
351
352 // Initialize the PDFs internal generator
353 if (_directVars.getSize()>0) {
354 cxcoutD(Generation) << "RooGenContext::initGenerator() initializing internal generator of model with code " << _code << endl ;
356 }
357}
358
359
360////////////////////////////////////////////////////////////////////////////////
361/// Generate one event. The 'remaining' integer is not used other than
362/// for printing messages
363
365{
366 if(_otherVars.getSize() > 0) {
367 // call the accept-reject generator to generate its variables
368
369 if (_updateFMaxPerEvent!=0) {
371 cxcoutD(Generation) << "RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << endl ;
372 _maxVar->setVal(max) ;
373 }
374
375 if (_generator) {
376 double resampleRatio(1) ;
377 const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
378 if (resampleRatio<1) {
379 coutI(Generation) << "RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
380 << resampleRatio << " due to increased maximum weight" << endl ;
381 resampleData(resampleRatio) ;
382 }
383 if(0 == subEvent) {
384 coutE(Generation) << "RooGenContext::generateEvent ERROR accept/reject generator failed" << endl ;
385 return;
386 }
387 theEvent.assignValueOnly(*subEvent) ;
388 //theEvent= *subEvent;
389
390 }
391 }
392
393 // Use the model's optimized generator, if one is available.
394 // The generator writes directly into our local 'event' since we attached it above.
395 if(_directVars.getSize() > 0) {
397 }
398
399 // Generate uniform variables (non-dependents)
400 for(auto * uniVar : _uniformVars) {
401 auto * arglv = dynamic_cast<RooAbsLValue*>(uniVar);
402 if (!arglv) {
403 coutE(Generation) << "RooGenContext::generateEvent(" << GetName() << ") ERROR: uniform variable " << uniVar->GetName() << " is not an lvalue" << endl ;
405 }
406 arglv->randomize() ;
407 }
408 theEvent.assign(_uniformVars) ;
409
410}
411
412
413////////////////////////////////////////////////////////////////////////////////
414/// Printing interface
415
416void RooGenContext::printMultiline(ostream &os, Int_t content, bool verbose, TString indent) const
417{
418 RooAbsGenContext::printMultiline(os,content,verbose,indent);
419 os << indent << " --- RooGenContext --- " << endl ;
420 os << indent << "Using PDF ";
422
423 if(verbose) {
424 os << indent << "Use PDF generator for " << _directVars << endl ;
425 os << indent << "Use MC sampling generator " << (_generator ? _generator->generatorName() : "<none>") << " for " << _otherVars << endl ;
426 if (_protoVars.getSize()>0) {
427 os << indent << "Prototype observables are " << _protoVars << endl ;
428 }
429 }
430}
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define oocoutE(o, a)
#define coutE(a)
#define ccxcoutI(a)
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
winID h direct
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
bool recursiveRedirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool recurseInNewSet=true)
Recursively replace all servers with the new servers in newSet.
void setOperMode(OperMode mode, bool recurseADirty=true)
Set the operation mode of this node.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
RooAbsCollection & assignValueOnly(const RooAbsCollection &other, bool forceIfSizeOne=false)
Sets the value of any argument in our set that also appears in the other set.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
RooArgSet _theEvent
Pointer to observable event being generated.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Interface for multi-line printing.
RooArgSet _protoVars
Prototype observables.
bool _verbose
Verbose messaging?
bool _isValid
Is context in valid state?
bool isValid() const
void resampleData(double &ratio)
Rescale existing output buffer with given ratio.
Abstract base class for objects that are lvalues, i.e.
Class RooAbsNumGenerator is the abstract base class for MC event generator implementations like RooAc...
virtual std::string const & generatorName() const =0
Return unique name of generator implementation.
virtual double getFuncMax()
void attachParameters(const RooArgSet &vars)
Reattach original parameters to function clone.
virtual const RooArgSet * generateEvent(UInt_t remaining, double &resampleRatio)=0
double getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:255
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
virtual void resetErrorCounters(Int_t resetValue=10)
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), bool force=true)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables std::liste...
virtual double maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
static void softAbort()
Soft abort function that interrupts macro execution but doesn't kill ROOT.
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Printing interface.
void generateEvent(RooArgSet &theEvent, Int_t remaining) override
Generate one event.
std::unique_ptr< RooRealIntegral > _acceptRejectFunc
Projection function to be passed to accept/reject sampler.
RooArgSet _otherVars
List of observables generated internally, randomly, and by accept/reject sampling.
Int_t _code
Internal generation code.
RooArgSet _uniformVars
void initGenerator(const RooArgSet &theEvent) override
Perform one-time initialization of the generator context.
~RooGenContext() override
Destructor.
RooArgSet _directVars
RooAbsPdf * _pdfClone
Clone of input p.d.f.
RooArgSet _cloneSet
Clone of all nodes of input p.d.f.
Int_t _updateFMaxPerEvent
If true, maximum p.d.f value needs to be recalculated for each event.
void attach(const RooArgSet &params) override
Attach the cloned model to the event buffer we will be filling.
RooAbsNumGenerator * _generator
MC sampling generation engine.
RooRealVar * _maxVar
Variable holding maximum value of p.d.f.
RooGenContext(const RooAbsPdf &model, const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false, const RooArgSet *forceDirect=nullptr)
Initialize a new context for generating events with the specified variables, using the specified PDF ...
RooAbsNumGenerator * createSampler(RooAbsReal &func, const RooArgSet &genVars, const RooArgSet &condVars, const RooNumGenConfig &config, bool verbose=false, RooAbsReal *maxFuncVal=nullptr)
Construct a numeric integrator instance that operates on function 'func' and is configured with 'conf...
static RooNumGenFactory & instance()
Static method returning reference to singleton instance of factory.
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setVal(double value) override
Set value of variable to 'value'.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Basic string class.
Definition TString.h:139
TString & Append(const char *cs)
Definition TString.h:576