KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVGVList.cpp
Go to the documentation of this file.
1 // Author: Daniel Cussol
2 //
3 // 17/02/2004
4 #include "Riostream.h"
5 #include "KVGVList.h"
6 
7 #include <KVEvent.h>
8 
10 
11 
12 //_________________________________________________________________
13 
14 
16 void KVGVList::init_KVGVList(void)
17 {
18  fNbBranch = 0;
19  fNbIBranch = 0;
20 }
21 
22 
23 
25 
27 {
28  init_KVGVList();
29 }
30 
31 
32 
33 
35 
37 {
38  init_KVGVList();
39  a.Copy(*this);
40 }
41 
42 
43 
48 
49 void KVGVList::Init(void)
50 {
51  // Initialisation of all global variables in list
52  // As this method may be called several times we ensure that
53  // variables are only initialised once
54 
55  this->R__FOR_EACH(KVVarGlob, ListInit)();
56 }
57 
58 
59 
62 
63 void KVGVList::Reset(void)
64 {
65  // Reset all variables before treating an event
66  this->R__FOR_EACH(KVVarGlob, Reset)();
67  fAbortEventAnalysis = false;
68 }
69 
70 
71 
76 
77 void KVGVList::Fill(const KVNucleus* c)
78 {
79  // Calls KVVarGlob::Fill(KVNucleus*) method of all one-body variables in the list
80  // for all KVNucleus satisfying the KVParticleCondition given to
81  // KVVarGlob::SetSelection() (if no selection given, all nuclei are used).
82 
83  fVG1.R__FOR_EACH(KVVarGlob, Fill)(c);
84 }
85 
86 
87 
88 
93 
94 void KVGVList::Fill2(const KVNucleus* c1, const KVNucleus* c2)
95 {
96  // Calls KVVarGlob::Fill(KVNucleus*,KVNucleus*) method of all two-body variables in the list
97  // for all pairs of KVNucleus (c1,c2) satisfying the KVParticleCondition given to
98  // KVVarGlob::SetSelection() (if no selection given, all nuclei are used).
99 
100  fVG2.R__FOR_EACH(KVVarGlob, Fill2)(c1, c2);
101 }
102 
103 
104 
107 
108 void KVGVList::FillN(const KVEvent* r)
109 {
110  // Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list
111  fVGN.R__FOR_EACH(KVVarGlob, FillN)(r);
112 }
113 
114 
115 
118 
120 {
121  // Calculate all 1-body observables after filling
122  TIter it(&fVG1);
123  KVVarGlob* vg;
124  while ((vg = (KVVarGlob*)it())) {
125  vg->Calculate();
126 #ifdef USING_ROOT6
127  if (!vg->TestEventSelection()) {
128  fAbortEventAnalysis = true;
129  break;
130  }
131 #endif
132  }
133 }
134 
135 
136 
139 
141 {
142  // Calculate all 2-body observables after filling
143  TIter it(&fVG2);
144  KVVarGlob* vg;
145  while ((vg = (KVVarGlob*)it())) {
146  vg->Calculate();
147 #ifdef USING_ROOT6
148  if (!vg->TestEventSelection()) {
149  fAbortEventAnalysis = true;
150  break;
151  }
152 #endif
153  }
154 }
155 
156 
157 
160 
162 {
163  // Calculate all N-body observables after filling
164  TIter it(&fVGN);
165  KVVarGlob* vg;
166  while ((vg = (KVVarGlob*)it())) {
167  vg->Calculate();
168 #ifdef USING_ROOT6
169  if (!vg->TestEventSelection()) {
170  fAbortEventAnalysis = true;
171  break;
172  }
173 #endif
174  }
175 }
176 
177 
178 
188 
190 {
191  // This method will calculate all global variables defined in the list for the event 'e'.
192  //
193  // Note that for 2-body variables, the Fill2 method will be called for each distinct pair of particles
194  // in the event, plus each pair made up of a particle and itself.
195  //
196  // For each variable for which an event selection condition was set (see KVVarGlob::SetEventSelection())
197  // the condition is tested as soon as the variable is calculated. If the condition is not satisfied,
198  // calculation of the other variables is abandonded and method AbortEventAnalysis() returns kTRUE.
199 
200  Reset();
201 
202  TIter it(this);
203  KVVarGlob* vg;
204  while ((vg = (KVVarGlob*)it())) {
205 
206  if (vg->IsGlobalVariable()) {
207  if (vg->IsNBody()) vg->FillN(e);
208  else {
209  for (KVEvent::Iterator it1 = OKEventIterator(*e).begin(); it1 != KVEvent::Iterator::End(); ++it1) {
210  if (vg->IsTwoBody()) {
211  for (KVEvent::Iterator it2(it1); it2 != KVEvent::Iterator::End(); ++it2) {
212  // we use every distinct pair of particles (including identical pairs) in the event
213  vg->Fill2(it1.get_pointer<const KVNucleus>(), it2.get_pointer<const KVNucleus>());
214  }
215  }
216  else {
217  vg->Fill(it1.get_pointer<const KVNucleus>());
218  }
219  }
220  }
221  }
222  vg->Calculate();
223 #ifdef USING_ROOT6
224  if ((fAbortEventAnalysis = !vg->TestEventSelection())) {
225  return;
226  }
227  vg->DefineNewFrame(e);
228 #endif
229  }
230 }
231 
232 
233 
236 
237 KVVarGlob* KVGVList::GetGV(const Char_t* nom) const
238 {
239  //Return pointer to global variable in list with name 'nom'
240 
241  return (KVVarGlob*) FindObject(nom);
242 }
243 
244 
245 
253 
255 {
256  // Overrides KVUniqueNameList::Add(TObject*) so that global variable pointers are sorted
257  // between the 3 lists used for 1-body, 2-body & N-body variables.
258  //
259  // We also print a warning in case the user tries to add a global variable with the
260  // same name as one already in the list. In the case we retain only the first global variable,
261  // any others with the same name are ignored
262 
263  KVUniqueNameList::Add(obj); // add object to main list, check duplicates
264  if (!ObjectAdded()) {
265  Warning("Add", "You tried to add a global variable with the same name as one already in the list");
266  Warning("Add", "Only the first variable added to the list will be used: name=%s class=%s",
267  GetGV(obj->GetName())->GetName(), GetGV(obj->GetName())->ClassName());
268  Warning("Add", "The following global variable (the one you tried to add) will be ignored:");
269  printf("\n");
270  obj->Print();
271  return;
272  }
273  if (obj->InheritsFrom("KVVarGlob")) {
274  // put global variable pointer in appropriate list
275  KVVarGlob* vg = (KVVarGlob*)obj;
276  if (vg->IsOneBody()) fVG1.Add(vg);
277  else if (vg->IsTwoBody()) fVG2.Add(vg);
278  else if (vg->IsNBody()) fVGN.Add(vg);
279  }
280 }
281 
282 
290 
292 {
293  // Overrides KVUniqueNameList::AddFirst(TObject*) so that global variable pointers are sorted
294  // between the 3 lists used for 1-body, 2-body & N-body variables.
295  //
296  // We also print a warning in case the user tries to add a global variable with the
297  // same name as one already in the list. In the case we retain only the first global variable,
298  // any others with the same name are ignored
299 
300  KVUniqueNameList::AddFirst(obj); // add object to main list, check duplicates
301  if (!ObjectAdded()) {
302  Warning("AddFirst", "You tried to add a global variable with the same name as one already in the list");
303  Warning("AddFirst", "Only the first variable added to the list will be used: name=%s class=%s",
304  GetGV(obj->GetName())->GetName(), GetGV(obj->GetName())->ClassName());
305  Warning("AddFirst", "The following global variable (the one you tried to add) will be ignored:");
306  printf("\n");
307  obj->Print();
308  return;
309  }
310  if (obj->InheritsFrom("KVVarGlob")) {
311  // put global variable pointer in appropriate list
312  KVVarGlob* vg = (KVVarGlob*)obj;
313  if (vg->IsOneBody()) fVG1.Add(vg);
314  else if (vg->IsTwoBody()) fVG2.Add(vg);
315  else if (vg->IsNBody()) fVGN.Add(vg);
316  }
317 }
318 
319 
320 
321 
335 
337 {
338  // Create a branch in the TTree for each global variable in the list.
339  // A leaf with the name of each global variable will be created to hold the
340  // value of the variable (result of KVVarGlob::GetValue() method).
341  // For multi-valued global variables we add a branch for each value with name
342  //
343  // `GVname.ValueName`
344  //
345  // To avoid problems with TTree::Draw, 'GVName' will be sanitized i.e. any
346  // mathematical symbols will be replaced by '_'.
347  //
348  // Any variable for which KVVarGlob::SetMaxNumBranches() was called with argument `0`
349  // will not be added to the TTree.
350 
351  if (!tree) return;
352 
353  // Make sure all variables are initialised before proceeding
354  Init();
355 
356  fNbBranch = 0;
357  fNbIBranch = 0;
358  fBranchVar.clear();
359  fIBranchVar.clear();
360 
361  // calculate the number of variables which will be stored in each vector
362  // the space has to be pre-reserved not added by successive "push_back"s, otherwise references to vector elements
363  // are invalidated every time the vector reallocates space
364 
365  TIter next(this);
366  KVVarGlob* ob;
367  while ((ob = (KVVarGlob*)next())) {
368  if (ob->GetNumberOfBranches()) { //skip variables for which KVVarGlob::SetMaxNumBranches(0) was called
369  if (ob->GetNumberOfValues() > 1) {
370  // multi-valued variable
371  for (int i = 0; i < ob->GetNumberOfBranches(); i++) {
372  if (ob->GetValueType(i) == 'I') ++fNbIBranch;
373  else ++fNbBranch;
374  }
375  }
376  else {
377  if (ob->GetValueType(0) == 'I') ++fNbIBranch;
378  else ++fNbBranch;
379  }
380  }
381  }
382 
383  fBranchVar.reserve(fNbBranch);
384  fIBranchVar.reserve(fNbIBranch);
385  fNbBranch = 0;
386  fNbIBranch = 0;
387 
388  next.Reset();
389 
390  while ((ob = (KVVarGlob*)next())) {
391  if (ob->GetNumberOfBranches()) { //skip variables for which KVVarGlob::SetMaxNumBranches(0) was called
392  TString sane_varname = NameSanitizer(ob->GetName());
393  if (ob->GetNumberOfValues() > 1) {
394  // multi-valued variable
395  for (int i = 0; i < ob->GetNumberOfBranches(); i++) {
396  // replace any nasty mathematical symbols which could pose problems
397  // in names of TTree leaves/branches
398  TString sane_name(ob->GetValueName(i));
399  sane_name.ReplaceAll("*", "star");
400  if (ob->GetValueType(i) == 'I') {
401  tree->Branch(Form("%s.%s", sane_varname.Data(), sane_name.Data()), &fIBranchVar[ fNbIBranch ], Form("%s.%s/I", sane_varname.Data(), sane_name.Data()));
402  ++fNbIBranch;
403  }
404  else {
405  tree->Branch(Form("%s.%s", sane_varname.Data(), sane_name.Data()), &fBranchVar[ fNbBranch ], Form("%s.%s/D", sane_varname.Data(), sane_name.Data()));
406  ++fNbBranch;
407  }
408  }
409  }
410  else {
411  if (ob->GetValueType(0) == 'I') {
412  tree->Branch(sane_varname, &fIBranchVar[ fNbIBranch ], Form("%s/I", sane_varname.Data()));
413  ++fNbIBranch;
414  }
415  else {
416  tree->Branch(sane_varname, &fBranchVar[ fNbBranch ], Form("%s/D", sane_varname.Data()));
417  ++fNbBranch;
418  }
419  }
420  }
421  }
422 }
423 
424 
425 
426 
433 
435 {
436  // Use this method ONLY if you first use MakeBranches(TTree*) in order to
437  // automatically create branches for your global variables.
438  // Call this method for each event in order to put the values of the variables
439  // in the branches ready for TTree::Fill to be called (note that TTree::Fill is not
440  // called in this method: you should do it after this).
441 
442  if (!fNbBranch && !fNbIBranch) return; // MakeBranches has not been called
443 
444  int INT_index = 0;
445  int FLT_index = 0;
446  TIter next(this);
447  KVVarGlob* ob;
448  while ((ob = (KVVarGlob*)next())) {
449  if (ob->GetNumberOfBranches()) { //skip variables for which KVVarGlob::SetMaxNumBranches(0) was called
450 
451  if (ob->GetNumberOfValues() > 1) {
452  // multi-valued variable
453  for (int j = 0; j < ob->GetNumberOfBranches(); j++) {
454  if (ob->GetValueType(j) == 'I') {
455  fIBranchVar[ INT_index ] = (Int_t)ob->GetValue(j);
456  ++INT_index;
457  }
458  else {
459  fBranchVar[ FLT_index ] = ob->GetValue(j);
460  ++FLT_index;
461  }
462  }
463  }
464  else {
465  if (ob->GetValueType(0) == 'I') {
466  fIBranchVar[ INT_index ] = (Int_t)ob->GetValue();
467  ++INT_index;
468  }
469  else {
470  fBranchVar[ FLT_index ] = ob->GetValue();
471  ++FLT_index;
472  }
473  }
474 
475  }
476  }
477 }
478 
479 
480 
525 
527 {
528  // Add an event classification object to the list, based on the named global variable
529  // (which must already be in the list).
530  //
531  // Returns a pointer to the object, in order to add either cuts or bins like so:
532  //
533  // #### Using cuts
534  //~~~~{.cpp}
535  // mtot_cuts->AddCut(5);
536  // mtot_cuts->AddCut(10);
537  // mtot_cuts->AddCut(15);
538  // mtot_cuts->AddCut(20);
539  // mtot_cuts->AddCut(25);
540  // mtot_cuts->AddCut(30);
541  //~~~~
542  // This will class events according to:
543  // | mtot | mtot_EC |
544  // |------------|---------|
545  // | <5 | 0 |
546  // | [5, 9] | 1 |
547  // | [10, 14] | 2 |
548  // | [15, 19] | 3 |
549  // | [20, 24] | 4 |
550  // | [25, 29] | 5 |
551  // | >=30 | 6 |
552  //
553  // `mtot_EC` is the default name for an event classification based on `mtot` and will be
554  // used for the branch added to the user's analysis TTree by method MakeBranches().
555  //
556  // #### Using bins
557  // ~~~~{.cpp}
558  // mtot_cuts->AddBin(5,10);
559  // mtot_cuts->AddBin(15,20);
560  // mtot_cuts->AddBin(25,30);
561  //~~~~
562  //This will class events according to the 2 bins with the following numbers:
563  // | mtot | mtot_EC |
564  // |----------|-------------|
565  // | [5, 9] | 0 |
566  // | [15, 19] | 1 |
567  // | [25, 29] | 2 |
568  // | other | -1 |
569  //
570  // Note that in this case, any value outside of a defined bin is unclassified.
571 
572  KVVarGlob* gv = GetGV(varname);
573  if (!gv) {
574  Warning("AddEventClassifier", "Variable %s not found in list. No classification possible.", varname.Data());
575  }
576  KVEventClassifier* ec = new KVEventClassifier(gv);
577  Add(ec);
578  return ec;
579 }
580 
581 
582 
614 
615 KVVarGlob* KVGVList::AddGV(const Char_t* class_name, const Char_t* name)
616 {
617  //Add a global variable to the list.
618  //
619  //`"class_name"` must be the name of a valid class inheriting from KVVarGlob, e.g. any of the default global
620  //variable classes defined as part of the standard KaliVeda package (see the [Global Variables module](group__GlobalVariables.html)),
621  //or the name of a user-defined class (see below).
622  //
623  //`"name"` is a unique name for the new global variable object which will be created and added to the
624  //list of global variables. This name can later be used to retrieve the object (see GetGV()).
625  //
626  //Returns pointer to new global variable object in case more than the usual default initialisation is necessary.
627  //
628  //### User-defined global variables
629  //The user may use her own global variables, without having to add them to the main libraries.
630  //If the given class name is not known, it is assumed to be a user-defined class and we attempt to compile and load
631  //the class from the user's source code. For this to work, the user must:
632  //
633  // 1. add to the ROOT macro path the directory where her class's source code is kept, e.g. in `$HOME/.rootrc` add the following line:
634  //
635  //~~~~~~~~~~~~~~~
636  // +Unix.*.Root.MacroPath: $(HOME)/myVarGlobs
637  //~~~~~~~~~~~~~~~
638  //
639  // 2. for each user-defined class, add a line to $HOME/.kvrootrc to define a "plugin". E.g. for a class called MyNewVarGlob,
640  //
641  //~~~~~~~~~~~~~~~
642  // +Plugin.KVVarGlob: MyNewVarGlob MyNewVarGlob MyNewVarGlob.cpp+ "MyNewVarGlob(const char*)"
643  //~~~~~~~~~~~~~~~
644  //
645  // It is assumed that `MyNewVarGlob.h` and `MyNewVarGlob.cpp` will be found in `$HOME/myVarGlobs` (in this example).
646  // The constructor taking a single character string argument (name of the variable) must be defined in the class.
647 
648  KVVarGlob* vg = prepareGVforAdding(class_name, name);
649  if (vg) Add(vg);
650  return vg;
651 }
652 
653 
654 
656 
657 KVVarGlob* KVGVList::prepareGVforAdding(const Char_t* class_name, const Char_t* name)
658 {
659  KVVarGlob* vg = nullptr;
660  TClass* clas = TClass::GetClass(class_name);
661  if (!clas) {
662  //class not in dictionary - user-defined class ? Look for plugin.
663  TPluginHandler* ph = KVBase::LoadPlugin("KVVarGlob", class_name);
664  if (!ph) {
665  //not found
666  Error("AddGV(const Char_t*,const Char_t*)",
667  "Called with class_name=%s.\nClass is unknown: not in standard libraries, and plugin (user-defined class) not found",
668  class_name);
669  return 0;
670  }
671  else {
672  vg = (KVVarGlob*) ph->ExecPlugin(1, name);
673  }
674  }
675  else if (!clas->InheritsFrom("KVVarGlob")) {
676  Error("AddGV(const Char_t*,const Char_t*)",
677  "%s is not a valid class deriving from KVVarGlob.",
678  class_name);
679  return nullptr;
680  }
681  else {
682  // need to use plugins to create even known classes, in order to call ctor with name
683  TPluginHandler* ph = KVBase::LoadPlugin("KVVarGlob", class_name);
684  if (!ph) {
685  // define plugin handler for known class
686  gPluginMgr->AddHandler("KVVarGlob", class_name, class_name, "KVMultiDetglobvars", Form("%s(const char*)", class_name));
687  ph = KVBase::LoadPlugin("KVVarGlob", class_name);
688  }
689  vg = (KVVarGlob*) ph->ExecPlugin(1, name);
690  }
691  return vg;
692 }
693 
694 
695 
700 
701 KVVarGlob* KVGVList::AddGVFirst(const Char_t* class_name, const Char_t* name)
702 {
703  // Add a global variable at the beginning of the list.
704  //
705  // See AddGV() for details.
706 
707  KVVarGlob* vg = prepareGVforAdding(class_name, name);
708  if (vg) AddFirst(vg);
709  return vg;
710 }
711 
712 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
ROOT::R::TRInterface & r
#define c(i)
#define e(i)
char Char_t
#define R__FOR_EACH(type, proc)
R__EXTERN TPluginManager * gPluginMgr
char * Form(const char *fmt,...)
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:756
Simple class for sorting events according to global variables.
static Iterator End()
Definition: KVEvent.h:347
Base class container for multi-particle events.
Definition: KVEvent.h:176
Manage a list of global variables.
Definition: KVGVList.h:121
void FillBranches()
Definition: KVGVList.cpp:434
KVVarGlob * AddGVFirst(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:701
Int_t fNbIBranch
Definition: KVGVList.h:126
void Reset(void)
Reset all variables before treating an event.
Definition: KVGVList.cpp:63
void Init(void)
Definition: KVGVList.cpp:49
std::vector< Double_t > fBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:123
TList fVG2
two-body variables
Definition: KVGVList.h:145
TString NameSanitizer(const Char_t *s) const
Definition: KVGVList.h:130
Int_t fNbBranch
Definition: KVGVList.h:125
virtual void AddFirst(TObject *obj)
Definition: KVGVList.cpp:291
void Calculate2()
Calculate all 2-body observables after filling.
Definition: KVGVList.cpp:140
void CalculateN()
Calculate all N-body observables after filling.
Definition: KVGVList.cpp:161
void CalculateGlobalVariables(KVEvent *e)
Definition: KVGVList.cpp:189
bool fAbortEventAnalysis
set to false if a global variable fails its own event selection criterion
Definition: KVGVList.h:127
void FillN(const KVEvent *e)
Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list.
Definition: KVGVList.cpp:108
std::vector< Int_t > fIBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:124
TList fVGN
N-body variables.
Definition: KVGVList.h:146
void Calculate()
Calculate all 1-body observables after filling.
Definition: KVGVList.cpp:119
KVVarGlob * prepareGVforAdding(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:657
TList fVG1
one-body variables
Definition: KVGVList.h:144
KVVarGlob * AddGV(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:615
virtual void Add(TObject *obj)
Definition: KVGVList.cpp:254
void Fill(const KVNucleus *c)
Definition: KVGVList.cpp:77
void init_KVGVList(void)
Definition: KVGVList.cpp:16
KVEventClassifier * AddEventClassifier(const TString &varname)
Definition: KVGVList.cpp:526
void MakeBranches(TTree *)
Definition: KVGVList.cpp:336
KVVarGlob * GetGV(const Char_t *nom) const
Return pointer to global variable in list with name 'nom'.
Definition: KVGVList.cpp:237
KVGVList(void)
Definition: KVGVList.cpp:26
void Fill2(const KVNucleus *c1, const KVNucleus *c2)
Definition: KVGVList.cpp:94
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
virtual TObject * FindObject(const char *name) const
Optimised list in which named objects can only be placed once.
Bool_t ObjectAdded() const
virtual void Add(TObject *obj)
virtual void AddFirst(TObject *obj)
Base class for all global variable implementations.
Definition: KVVarGlob.h:217
virtual void Calculate()=0
Bool_t IsOneBody() const
Definition: KVVarGlob.h:351
virtual Int_t GetNumberOfValues() const
Definition: KVVarGlob.h:603
Double_t GetValue(void) const
Definition: KVVarGlob.h:412
void Fill(const KVNucleus *c)
Definition: KVVarGlob.h:388
virtual void FillN(const KVEvent *)
Definition: KVVarGlob.h:406
Int_t GetNumberOfBranches() const
Definition: KVVarGlob.h:609
void Fill2(const KVNucleus *n1, const KVNucleus *n2)
Definition: KVVarGlob.h:396
virtual Bool_t IsGlobalVariable() const
Definition: KVVarGlob.h:369
virtual TString GetValueName(Int_t i) const
Definition: KVVarGlob.h:619
virtual Char_t GetValueType(Int_t) const
Definition: KVVarGlob.h:628
void DefineNewFrame(KVEvent *e) const
Definition: KVVarGlob.h:697
bool TestEventSelection() const
Definition: KVVarGlob.h:676
Bool_t IsNBody() const
Definition: KVVarGlob.h:363
Bool_t IsTwoBody() const
Definition: KVVarGlob.h:357
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
Bool_t InheritsFrom(const char *cl) const
void Reset()
virtual void Add(TObject *obj)
virtual const char * GetName() const
virtual const char * GetName() const
virtual const char * ClassName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Print(Option_t *option="") const
Longptr_t ExecPlugin(int nargs, const T &... params)
void AddHandler(const char *base, const char *regexp, const char *className, const char *pluginName, const char *ctor=0, const char *origin=0)
const char * Data() const
TString & ReplaceAll(const char *s1, const char *s2)
return c1
return c2
KVEvent::Iterator begin() const
Definition: KVEvent.h:509
auto * a