ROOT  6.06/09
Reference Guide
TRootSniffer.cxx
Go to the documentation of this file.
1 // $Id$
2 // Author: Sergey Linev 22/12/2013
3 
4 #include "TRootSniffer.h"
5 
6 #include "TH1.h"
7 #include "TGraph.h"
8 #include "TProfile.h"
9 #include "TCanvas.h"
10 #include "TFile.h"
11 #include "TKey.h"
12 #include "TList.h"
13 #include "TMemFile.h"
14 #include "TStreamerInfo.h"
15 #include "TBufferFile.h"
16 #include "TBufferJSON.h"
17 #include "TBufferXML.h"
18 #include "TROOT.h"
19 #include "TTimer.h"
20 #include "TFolder.h"
21 #include "TTree.h"
22 #include "TBranch.h"
23 #include "TLeaf.h"
24 #include "TClass.h"
25 #include "TMethod.h"
26 #include "TFunction.h"
27 #include "TMethodArg.h"
28 #include "TMethodCall.h"
29 #include "TRealData.h"
30 #include "TDataMember.h"
31 #include "TDataType.h"
32 #include "TBaseClass.h"
33 #include "TObjString.h"
34 #include "TUrl.h"
35 #include "TImage.h"
36 #include "RZip.h"
37 #include "RVersion.h"
38 #include "TRootSnifferStore.h"
39 #include "THttpCallArg.h"
40 
41 #include <stdlib.h>
42 #include <vector>
43 #include <string.h>
44 
45 const char *item_prop_kind = "_kind";
46 const char *item_prop_more = "_more";
47 const char *item_prop_title = "_title";
48 const char *item_prop_hidden = "_hidden";
49 const char *item_prop_typename = "_typename";
50 const char *item_prop_arraydim = "_arraydim";
51 const char *item_prop_realname = "_realname"; // real object name
52 const char *item_prop_user = "_username";
53 const char *item_prop_autoload = "_autoload";
54 const char *item_prop_rootversion = "_root_version";
55 
56 
57 // ============================================================================
58 
59 //////////////////////////////////////////////////////////////////////////
60 // //
61 // TRootSnifferScanRec //
62 // //
63 // Structure used to scan hierarchies of ROOT objects //
64 // Represents single level of hierarchy //
65 // //
66 //////////////////////////////////////////////////////////////////////////
67 
68 ////////////////////////////////////////////////////////////////////////////////
69 /// constructor
70 
72  fParent(0),
73  fMask(0),
74  fSearchPath(0),
75  fLevel(0),
76  fItemName(),
77  fItemsNames(),
78  fRestriction(0),
79  fStore(0),
80  fHasMore(kFALSE),
81  fNodeStarted(kFALSE),
82  fNumFields(0),
83  fNumChilds(0)
84 {
86 }
87 
88 ////////////////////////////////////////////////////////////////////////////////
89 /// destructor
90 
92 {
93  CloseNode();
94 }
95 
96 ////////////////////////////////////////////////////////////////////////////////
97 /// record field for current element
98 
99 void TRootSnifferScanRec::SetField(const char *name, const char *value, Bool_t with_quotes)
100 {
101  if (CanSetFields()) fStore->SetField(fLevel, name, value, with_quotes);
102  fNumFields++;
103 }
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 /// indicates that new child for current element will be started
107 
109 {
111  fNumChilds++;
112 }
113 
114 ////////////////////////////////////////////////////////////////////////////////
115 /// constructs item name from object name
116 /// if special symbols like '/', '#', ':', '&', '?' are used in object name
117 /// they will be replaced with '_'.
118 /// To avoid item name duplication, additional id number can be appended
119 
120 void TRootSnifferScanRec::MakeItemName(const char *objname, TString &itemname)
121 {
122  std::string nnn = objname;
123 
124  size_t pos;
125 
126  // replace all special symbols which can make problem to navigate in hierarchy
127  while ((pos = nnn.find_first_of("- []<>#:&?/\'\"\\")) != std::string::npos)
128  nnn.replace(pos, 1, "_");
129 
130  itemname = nnn.c_str();
131  Int_t cnt = 0;
132 
133  while (fItemsNames.FindObject(itemname.Data())) {
134  itemname.Form("%s_%d", nnn.c_str(), cnt++);
135  }
136 
137  fItemsNames.Add(new TObjString(itemname.Data()));
138 }
139 
140 ////////////////////////////////////////////////////////////////////////////////
141 /// Produce full name, including all parents
142 
144 {
145  if (!prnt) prnt = fParent;
146 
147  if (prnt) {
148  prnt->BuildFullName(buf);
149 
150  buf.Append("/");
151  buf.Append(fItemName);
152  }
153 }
154 
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// creates new node with specified name
158 /// if special symbols like "[]&<>" are used, node name
159 /// will be replaced by default name like "extra_item_N" and
160 /// original node name will be recorded as "_original_name" field
161 /// Optionally, object name can be recorded as "_realname" field
162 
163 void TRootSnifferScanRec::CreateNode(const char *_node_name)
164 {
165  if (!CanSetFields()) return;
166 
168 
170 
171  if (fStore) fStore->CreateNode(fLevel, _node_name);
172 }
173 
174 ////////////////////////////////////////////////////////////////////////////////
175 /// close started node
176 
178 {
179  if (fStore && fNodeStarted) {
182  }
183 }
184 
185 ////////////////////////////////////////////////////////////////////////////////
186 /// set root class name as node kind
187 /// in addition, path to master item (streamer info) specified
188 /// Such master item required to correctly unstream data on JavaScript
189 
191 {
192  if ((cl != 0) && CanSetFields())
193  SetField(item_prop_kind, TString::Format("ROOT.%s", cl->GetName()));
194 }
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 /// returns true if scanning is done
198 /// Can happen when searched element is found
199 
201 {
202  if (fStore == 0)
203  return kFALSE;
204 
205  if ((fMask & kSearch) && fStore->GetResPtr())
206  return kTRUE;
207 
208  if ((fMask & kCheckChilds) && fStore->GetResPtr() &&
209  (fStore->GetResNumChilds() >= 0))
210  return kTRUE;
211 
212  return kFALSE;
213 }
214 
215 ////////////////////////////////////////////////////////////////////////////////
216 /// Checks if result will be accepted.
217 /// Used to verify if sniffer should read object from the file
218 
220 {
221  if (Done()) return kFALSE;
222 
223  // only when doing search, result will be propagated
224  if ((fMask & (kSearch | kCheckChilds)) == 0) return kFALSE;
225 
226  // only when full search path is scanned
227  if (fSearchPath != 0) return kFALSE;
228 
229  if (fStore == 0) return kFALSE;
230 
231  return kTRUE;
232 }
233 
234 ////////////////////////////////////////////////////////////////////////////////
235 /// set results of scanning
236 /// when member should be specified, use SetFoundResult instead
237 
239 {
240  if (member==0) return SetFoundResult(obj, cl);
241 
242  fStore->Error("SetResult", "When member specified, pointer on object (not member) should be provided; use SetFoundResult");
243  return kFALSE;
244 }
245 
246 ////////////////////////////////////////////////////////////////////////////////
247 /// set results of scanning
248 /// when member specified, obj is pointer on object to which member belongs
249 
251 {
252  if (Done()) return kTRUE;
253 
254  if (!IsReadyForResult()) return kFALSE;
255 
256  fStore->SetResult(obj, cl, member, fNumChilds, fRestriction);
257 
258  return Done();
259 }
260 
261 
262 ////////////////////////////////////////////////////////////////////////////////
263 /// returns current depth of scanned hierarchy
264 
266 {
267  Int_t cnt = 0;
268  const TRootSnifferScanRec *rec = this;
269  while (rec->fParent) {
270  rec = rec->fParent;
271  cnt++;
272  }
273 
274  return cnt;
275 }
276 
277 ////////////////////////////////////////////////////////////////////////////////
278 /// returns true if current item can be expanded - means one could explore
279 /// objects members
280 
282 {
283  if (fMask & (kExpand | kSearch | kCheckChilds)) return kTRUE;
284 
285  if (!fHasMore) return kFALSE;
286 
287  // if parent has expand mask, allow to expand item
288  if (fParent && (fParent->fMask & kExpand)) return kTRUE;
289 
290  return kFALSE;
291 }
292 
293 ////////////////////////////////////////////////////////////////////////////////
294 /// returns read-only flag for current item
295 /// Depends from default value and current restrictions
296 
298 {
299  if (fRestriction==0) return dflt;
300 
301  return fRestriction!=2;
302 }
303 
304 ////////////////////////////////////////////////////////////////////////////////
305 /// Method verifies if new level of hierarchy
306 /// should be started with provided object.
307 /// If required, all necessary nodes and fields will be created
308 /// Used when different collection kinds should be scanned
309 
311  const char *obj_name, TRootSniffer* sniffer)
312 {
313  if (super.Done()) return kFALSE;
314 
315  if ((obj != 0) && (obj_name == 0)) obj_name = obj->GetName();
316 
317  // exclude zero names
318  if ((obj_name == 0) || (*obj_name == 0)) return kFALSE;
319 
320  const char *full_name = 0;
321 
322  // remove slashes from file names
323  if (obj && obj->InheritsFrom(TDirectoryFile::Class())) {
324  const char *slash = strrchr(obj_name, '/');
325  if (slash != 0) {
326  full_name = obj_name;
327  obj_name = slash + 1;
328  if (*obj_name == 0) obj_name = "file";
329  }
330  }
331 
332  super.MakeItemName(obj_name, fItemName);
333 
334  if (sniffer && sniffer->HasRestriction(fItemName.Data())) {
335  // check restriction more precisely
336  TString fullname;
337  BuildFullName(fullname, &super);
338  fRestriction = sniffer->CheckRestriction(fullname.Data());
339  if (fRestriction<0) return kFALSE;
340  }
341 
342  fParent = &super;
343  fLevel = super.fLevel;
344  fStore = super.fStore;
345  fSearchPath = super.fSearchPath;
346  fMask = super.fMask & kActions;
347  if (fRestriction==0) fRestriction = super.fRestriction; // get restriction from parent
348  Bool_t topelement(kFALSE);
349 
350  if (fMask & kScan) {
351  // if scanning only fields, ignore all childs
352  if (super.ScanOnlyFields()) return kFALSE;
353  // only when doing scan, increment level, used for text formatting
354  fLevel++;
355  } else {
356  if (fSearchPath == 0) return kFALSE;
357 
358  if (strncmp(fSearchPath, fItemName.Data(), fItemName.Length()) != 0)
359  return kFALSE;
360 
361  const char *separ = fSearchPath + fItemName.Length();
362 
363  Bool_t isslash = kFALSE;
364  while (*separ == '/') {
365  separ++;
366  isslash = kTRUE;
367  }
368 
369  if (*separ == 0) {
370  fSearchPath = 0;
371  if (fMask & kExpand) {
372  topelement = kTRUE;
373  fMask = (fMask & kOnlyFields) | kScan;
374  fHasMore = (fMask & kOnlyFields) == 0;
375  }
376  } else {
377  if (!isslash) return kFALSE;
378  fSearchPath = separ;
379  }
380  }
381 
383 
384  if ((obj_name != 0) && (fItemName != obj_name))
385  SetField(item_prop_realname, obj_name);
386 
387  if (full_name != 0)
388  SetField("_fullname", full_name);
389 
390  if (topelement)
392 
393  if (topelement && sniffer->GetAutoLoad())
395 
396  return kTRUE;
397 }
398 
399 
400 // ====================================================================
401 
402 //////////////////////////////////////////////////////////////////////////
403 // //
404 // TRootSniffer //
405 // //
406 // Sniffer of ROOT objects, data provider for THttpServer //
407 // Provides methods to scan different structures like folders, //
408 // directories, files, trees, collections //
409 // Can locate objects (or its data member) per name //
410 // Can be extended to application-specific classes //
411 // //
412 //////////////////////////////////////////////////////////////////////////
413 
415 
416 ////////////////////////////////////////////////////////////////////////////////
417 /// constructor
418 
419 TRootSniffer::TRootSniffer(const char *name, const char *objpath) :
420  TNamed(name, "sniffer of root objects"),
421  fObjectsPath(objpath),
422  fMemFile(0),
423  fSinfo(0),
424  fReadOnly(kTRUE),
425  fScanGlobalDir(kTRUE),
426  fCurrentArg(0),
427  fCurrentRestrict(0),
428  fCurrentAllowedMethods(0),
429  fRestrictions(),
430  fAutoLoad()
431 {
432  fRestrictions.SetOwner(kTRUE);
433 }
434 
435 ////////////////////////////////////////////////////////////////////////////////
436 /// destructor
437 
439 {
440  if (fSinfo) {
441  delete fSinfo;
442  fSinfo = 0;
443  }
444 
445  if (fMemFile) {
446  delete fMemFile;
447  fMemFile = 0;
448  }
449 }
450 
451 ////////////////////////////////////////////////////////////////////////////////
452 /// set current http arguments, which then used in different process methods
453 /// For instance, if user authorized with some user name,
454 /// depending from restrictions some objects will be invisible
455 /// or user get full access to the element
456 
458 {
459  fCurrentArg = arg;
460  fCurrentRestrict = 0;
462 }
463 
464 ////////////////////////////////////////////////////////////////////////////////
465 /// Restrict access to the specified location
466 ///
467 /// Hides or provides read-only access to different parts of the hierarchy
468 /// Restriction done base on user-name specified with http requests
469 /// Options can be specified in URL style (separated with &)
470 /// Following parameters can be specified:
471 /// visible = [all|user(s)] - make item visible for all users or only specified user
472 /// hidden = [all|user(s)] - make item hidden from all users or only specified user
473 /// readonly = [all|user(s)] - make item read-only for all users or only specified user
474 /// allow = [all|user(s)] - make full access for all users or only specified user
475 /// allow_method = method(s) - allow method(s) execution even when readonly flag specified for the object
476 /// Like make command seen by all but can be executed only by admin
477 /// sniff->Restrict("/CmdReset","allow=admin");
478 /// Or fully hide command from guest account
479 /// sniff->Restrict("/CmdRebin","hidden=guest");
480 
481 void TRootSniffer::Restrict(const char* path, const char* options)
482 {
483  const char* rslash = strrchr(path,'/');
484  if (rslash) rslash++;
485  if ((rslash==0) || (*rslash==0)) rslash = path;
486 
487  fRestrictions.Add(new TNamed(rslash, TString::Format("%s%s%s", path,"%%%",options).Data()));
488 }
489 
490 ////////////////////////////////////////////////////////////////////////////////
491 /// When specified, _autoload attribute will be always add
492 /// to top element of h.json/h.hml requests
493 /// Used to instruct browser automatically load special code
494 
495 void TRootSniffer::SetAutoLoad(const char* scripts)
496 {
497  fAutoLoad = scripts!=0 ? scripts : "";
498 }
499 
500 ////////////////////////////////////////////////////////////////////////////////
501 /// return name of configured autoload scripts (or 0)
502 
503 const char* TRootSniffer::GetAutoLoad() const
504 {
505  return fAutoLoad.Length() > 0 ? fAutoLoad.Data() : 0;
506 }
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Made fast check if item with specified name is in restriction list
510 /// If returns true, requires precise check with CheckRestriction() method
511 
512 Bool_t TRootSniffer::HasRestriction(const char* item_name)
513 {
514  if ((item_name==0) || (*item_name==0) || (fCurrentArg==0)) return kFALSE;
515 
516  return fRestrictions.FindObject(item_name)!=0;
517 }
518 
519 ////////////////////////////////////////////////////////////////////////////////
520 /// return 2 when option match to current user name
521 /// return 1 when option==all
522 /// return 0 when option does not match user name
523 
525 {
526  const char* username = fCurrentArg ? fCurrentArg->GetUserName() : 0;
527 
528  if ((username==0) || (option == 0) || (*option==0)) return 0;
529 
530  if (strcmp(option,"all") == 0) return 1;
531 
532  if (strcmp(username, option) == 0) return 2;
533 
534  if (strstr(option, username) == 0) return -1;
535 
536  TObjArray* arr = TString(option).Tokenize(",");
537 
538  Bool_t find = arr->FindObject(username) != 0;
539 
540  delete arr;
541 
542  return find ? 2 : -1;
543 }
544 
545 ////////////////////////////////////////////////////////////////////////////////
546 /// Checked if restriction is applied to the item
547 /// full_item_name should have full path to the item
548 ///
549 /// Returns -1 - object invisible, cannot be accessed or listed
550 /// 0 - no explicit restrictions, use default
551 /// 1 - read-only access
552 /// 2 - full access
553 
554 Int_t TRootSniffer::CheckRestriction(const char* full_item_name)
555 {
556  if ((full_item_name==0) || (*full_item_name==0)) return 0;
557 
558  const char* item_name = strrchr(full_item_name,'/');
559  if (item_name) item_name++;
560  if ((item_name==0) || (*item_name==0)) item_name = full_item_name;
561 
562  TString pattern1 = TString("*/") + item_name + "%%%";
563  TString pattern2 = TString(full_item_name) + "%%%";
564 
565  const char* options = 0;
567  TObject* obj;
568 
569  while ((obj = iter()) != 0) {
570  const char* title = obj->GetTitle();
571 
572  if (strstr(title,pattern1.Data())==title) { options = title + pattern1.Length(); break; }
573  if (strstr(title,pattern2.Data())==title) { options = title + pattern2.Length(); break; }
574  }
575 
576  if (options==0) return 0;
577 
578  TUrl url;
579  url.SetOptions(options);
580  url.ParseOptions();
581 
582  Int_t can_see = WithCurrentUserName(url.GetValueFromOptions("visible")) -
584 
585  Int_t can_access = WithCurrentUserName(url.GetValueFromOptions("allow")) -
586  WithCurrentUserName(url.GetValueFromOptions("readonly"));
587 
588  if (can_access > 0) return 2; // first of all, if access enabled, provide it
589  if (can_see < 0) return -1; // if object to be hidden, do it
590 
591  const char* methods = url.GetValueFromOptions("allow_method");
592  if (methods!=0) fCurrentAllowedMethods = methods;
593 
594  if (can_access < 0) return 1; // read-only access
595 
596  return 0; // default behavior
597 }
598 
599 ////////////////////////////////////////////////////////////////////////////////
600 /// scan object data members
601 /// some members like enum or static members will be excluded
602 
604 {
605  if ((cl == 0) || (ptr == 0) || rec.Done()) return;
606 
607  // ensure that real class data (including parents) exists
608  if (!(cl->Property() & kIsAbstract)) cl->BuildRealData();
609 
610  // scan only real data
611  TObject *obj = 0;
613  while ((obj = iter()) != 0) {
614  TRealData *rdata = dynamic_cast<TRealData *>(obj);
615  if ((rdata == 0) || strchr(rdata->GetName(),'.')) continue;
616 
617  TDataMember *member = rdata->GetDataMember();
618  // exclude enum or static variables
619  if ((member == 0) || (member->Property() & (kIsStatic | kIsEnum | kIsUnion))) continue;
620  char *member_ptr = ptr + rdata->GetThisOffset();
621 
622  if (member->IsaPointer()) member_ptr = *((char **) member_ptr);
623 
624  TRootSnifferScanRec chld;
625 
626  if (chld.GoInside(rec, member, 0, this)) {
627 
628  TClass *mcl = (member->IsBasic() || member->IsSTLContainer()) ? 0 :
629  gROOT->GetClass(member->GetTypeName());
630 
631  Int_t coll_offset = mcl ? mcl->GetBaseClassOffset(TCollection::Class()) : -1;
632  if (coll_offset >= 0) {
633  chld.SetField(item_prop_more, "true", kFALSE);
634  chld.fHasMore = kTRUE;
635  }
636 
637  if (chld.SetFoundResult(ptr, cl, member)) break;
638 
639  const char *title = member->GetTitle();
640  if ((title != 0) && (strlen(title) != 0))
641  chld.SetField(item_prop_title, title);
642 
643  if (member->GetTypeName())
644  chld.SetField(item_prop_typename, member->GetTypeName());
645 
646  if (member->GetArrayDim() > 0) {
647  // store array dimensions in form [N1,N2,N3,...]
648  TString dim("[");
649  for (Int_t n = 0; n < member->GetArrayDim(); n++) {
650  if (n > 0) dim.Append(",");
651  dim.Append(TString::Format("%d", member->GetMaxIndex(n)));
652  }
653  dim.Append("]");
654  chld.SetField(item_prop_arraydim, dim, kFALSE);
655  } else
656  if (member->GetArrayIndex()!=0) {
657  TRealData *idata = cl->GetRealData(member->GetArrayIndex());
658  TDataMember *imember = (idata!=0) ? idata->GetDataMember() : 0;
659  if ((imember!=0) && (strcmp(imember->GetTrueTypeName(),"int")==0)) {
660  Int_t arraylen = *((int *) (ptr + idata->GetThisOffset()));
661  chld.SetField(item_prop_arraydim, TString::Format("[%d]", arraylen), kFALSE);
662  }
663  }
664 
665  chld.SetRootClass(mcl);
666 
667  if (chld.CanExpandItem()) {
668  if (coll_offset >= 0) {
669  // chld.SetField("#members", "true", kFALSE);
670  ScanCollection(chld, (TCollection *)(member_ptr + coll_offset));
671  }
672  }
673 
674  if (chld.SetFoundResult(ptr, cl, member)) break;
675  }
676  }
677 }
678 
679 ////////////////////////////////////////////////////////////////////////////////
680 /// scans object properties
681 /// here such fields as _autoload or _icon properties depending on class or object name could be assigned
682 /// By default properties, coded in the Class title are scanned. Example:
683 /// ClassDef(UserClassName, 1) // class comments *SNIFF* _field1=value _field2="string value"
684 /// Here *SNIFF* mark is important. After it all expressions like field=value are parsed
685 /// One could use double quotes to code string values with spaces.
686 /// Fields separated from each other with spaces
687 
689 {
690  TClass* cl = obj ? obj->IsA() : 0;
691 
692  const char* pos = strstr(cl ? cl->GetTitle() : "", "*SNIFF*");
693  if (pos==0) return;
694 
695  pos += 7;
696  while (*pos != 0) {
697  if (*pos == ' ') { pos++; continue; }
698  // first locate identifier
699  const char* pos0 = pos;
700  while ((*pos != 0) && (*pos != '=')) pos++;
701  if (*pos == 0) return;
702  TString name(pos0, pos-pos0);
703  pos++;
704  Bool_t quotes = (*pos == '\"');
705  if (quotes) pos++;
706  pos0 = pos;
707  // then value with or without quotes
708  while ((*pos != 0) && (*pos != (quotes ? '\"' : ' '))) pos++;
709  TString value(pos0, pos-pos0);
710  rec.SetField(name, value);
711  if (quotes) pos++;
712  pos++;
713  }
714 }
715 
716 ////////////////////////////////////////////////////////////////////////////////
717 /// scans object childs (if any)
718 /// here one scans collection, branches, trees and so on
719 
721 {
722  if (obj->InheritsFrom(TFolder::Class())) {
723  ScanCollection(rec, ((TFolder *) obj)->GetListOfFolders());
724  } else if (obj->InheritsFrom(TDirectory::Class())) {
725  TDirectory *dir = (TDirectory *) obj;
726  ScanCollection(rec, dir->GetList(), 0, dir->GetListOfKeys());
727  } else if (obj->InheritsFrom(TTree::Class())) {
728  if (!rec.IsReadOnly(fReadOnly)) {
729  rec.SetField("_player", "JSROOT.drawTreePlayer");
730  rec.SetField("_prereq", "jq2d");
731  }
732  ScanCollection(rec, ((TTree *) obj)->GetListOfLeaves());
733  } else if (obj->InheritsFrom(TBranch::Class())) {
734  ScanCollection(rec, ((TBranch *) obj)->GetListOfLeaves());
735  } else if (rec.CanExpandItem()) {
736  ScanObjectMembers(rec, obj->IsA(), (char *) obj);
737  }
738 }
739 
740 ////////////////////////////////////////////////////////////////////////////////
741 /// scan collection content
742 
744  const char *foldername, TCollection *keys_lst)
745 {
746  if (((lst == 0) || (lst->GetSize() == 0)) && ((keys_lst == 0) || (keys_lst->GetSize() == 0))) return;
747 
748  TRootSnifferScanRec folderrec;
749  if (foldername) {
750  if (!folderrec.GoInside(rec, 0, foldername, this)) return;
751  }
752 
753  TRootSnifferScanRec &master = foldername ? folderrec : rec;
754 
755  if (lst != 0) {
756  TIter iter(lst);
757  TObject *next = iter();
758  Bool_t isany = kFALSE;
759 
760  while (next!=0) {
761  if (IsItemField(next)) {
762  // special case - in the beginning one could have items for master folder
763  if (!isany && (next->GetName() != 0) && ((*(next->GetName()) == '_') || master.ScanOnlyFields()))
764  master.SetField(next->GetName(), next->GetTitle());
765  next = iter();
766  continue;
767  }
768 
769  isany = kTRUE;
770  TObject* obj = next;
771 
772  TRootSnifferScanRec chld;
773  if (!chld.GoInside(master, obj, 0, this)) { next = iter(); continue; }
774 
775  if (chld.SetResult(obj, obj->IsA())) return;
776 
777  Bool_t has_kind(kFALSE), has_title(kFALSE);
778 
779  ScanObjectProperties(chld, obj);
780  // now properties, coded as TNamed objects, placed after object in the hierarchy
781  while ((next = iter()) != 0) {
782  if (!IsItemField(next)) break;
783  if ((next->GetName() != 0) && ((*(next->GetName()) == '_') || chld.ScanOnlyFields())) {
784  // only fields starting with _ are stored
785  chld.SetField(next->GetName(), next->GetTitle());
786  if (strcmp(next->GetName(), item_prop_kind)==0) has_kind = kTRUE;
787  if (strcmp(next->GetName(), item_prop_title)==0) has_title = kTRUE;
788  }
789  }
790 
791  if (!has_kind) chld.SetRootClass(obj->IsA());
792  if (!has_title && (obj->GetTitle()!=0)) chld.SetField(item_prop_title, obj->GetTitle());
793 
794  ScanObjectChilds(chld, obj);
795 
796  if (chld.SetResult(obj, obj->IsA())) return;
797  }
798  }
799 
800  if (keys_lst != 0) {
801  TIter iter(keys_lst);
802  TObject *kobj(0);
803 
804  while ((kobj = iter()) != 0) {
805  TKey *key = dynamic_cast<TKey *>(kobj);
806  if (key == 0) continue;
807  TObject *obj = (lst == 0) ? 0 : lst->FindObject(key->GetName());
808 
809  // even object with the name exists, it should also match with class name
810  if ((obj!=0) && (strcmp(obj->ClassName(),key->GetClassName())!=0)) obj = 0;
811 
812  // if object of that name and of that class already in the list, ignore appropriate key
813  if ((obj != 0) && (master.fMask & TRootSnifferScanRec::kScan)) continue;
814 
815  Bool_t iskey = kFALSE;
816  // if object not exists, provide key itself for the scan
817  if (obj == 0) { obj = key; iskey = kTRUE; }
818 
819  TRootSnifferScanRec chld;
820  TString fullname = TString::Format("%s;%d", key->GetName(), key->GetCycle());
821 
822  if (chld.GoInside(master, obj, fullname.Data(), this)) {
823 
824  if (!chld.IsReadOnly(fReadOnly) && iskey && chld.IsReadyForResult()) {
825  TObject *keyobj = key->ReadObj();
826  if (keyobj != 0)
827  if (chld.SetResult(keyobj, keyobj->IsA())) return;
828  }
829 
830  if (chld.SetResult(obj, obj->IsA())) return;
831 
832  TClass *obj_class = obj->IsA();
833 
834  ScanObjectProperties(chld, obj);
835 
836  if (obj->GetTitle()!=0) chld.SetField(item_prop_title, obj->GetTitle());
837 
838  // special handling of TKey class - in non-readonly mode
839  // sniffer allowed to fetch objects
840  if (!chld.IsReadOnly(fReadOnly) && iskey) {
841  if (strcmp(key->GetClassName(), "TDirectoryFile") == 0) {
842  if (chld.fLevel == 0) {
843  TDirectory *dir = dynamic_cast<TDirectory *>(key->ReadObj());
844  if (dir != 0) {
845  obj = dir;
846  obj_class = dir->IsA();
847  }
848  } else {
849  chld.SetField(item_prop_more, "true", kFALSE);
850  chld.fHasMore = kTRUE;
851  }
852  } else {
853  obj_class = TClass::GetClass(key->GetClassName());
854  if (obj_class && obj_class->InheritsFrom(TTree::Class())) {
855  if (rec.CanExpandItem()) {
856  // it is requested to expand tree element - read it
857  obj = key->ReadObj();
858  if (obj) obj_class = obj->IsA();
859  } else {
860  rec.SetField("_player", "JSROOT.drawTreePlayerKey");
861  rec.SetField("_prereq", "jq2d");
862  // rec.SetField("_more", "true"); // one could allow to extend
863  }
864  }
865  }
866  }
867 
868  rec.SetRootClass(obj_class);
869 
870  ScanObjectChilds(chld, obj);
871 
872  // here we should know how many childs are accumulated
873  if (chld.SetResult(obj, obj_class)) return;
874  }
875  }
876  }
877 }
878 
879 ////////////////////////////////////////////////////////////////////////////////
880 /// scan complete ROOT objects hierarchy
881 /// For the moment it includes objects in gROOT directory
882 /// and list of canvases and files
883 /// Also all registered objects are included.
884 /// One could reimplement this method to provide alternative
885 /// scan methods or to extend some collection kinds
886 
888 {
889  rec.SetField(item_prop_kind, "ROOT.Session");
892 
893  // should be on the top while //root/http folder could have properties for itself
894  TFolder *topf = dynamic_cast<TFolder *>(gROOT->FindObject("//root/http"));
895  if (topf) {
896  rec.SetField(item_prop_title, topf->GetTitle());
897  ScanCollection(rec, topf->GetListOfFolders());
898  }
899 
900  {
901  TRootSnifferScanRec chld;
902  if (chld.GoInside(rec, 0, "StreamerInfo", this)) {
903  chld.SetField(item_prop_kind, "ROOT.TStreamerInfoList");
904  chld.SetField(item_prop_title, "List of streamer infos for binary I/O");
905  chld.SetField(item_prop_hidden, "true");
906  chld.SetField("_after_request", "JSROOT.MarkAsStreamerInfo");
907  }
908  }
909 
910  if (IsScanGlobalDir()) {
911  ScanCollection(rec, gROOT->GetList());
912 
913  ScanCollection(rec, gROOT->GetListOfCanvases(), "Canvases");
914 
915  ScanCollection(rec, gROOT->GetListOfFiles(), "Files");
916  }
917 }
918 
919 ////////////////////////////////////////////////////////////////////////////////
920 /// return true if object can be drawn
921 
923 {
924  if (cl == 0) return kFALSE;
925  if (cl->InheritsFrom(TH1::Class())) return kTRUE;
926  if (cl->InheritsFrom(TGraph::Class())) return kTRUE;
927  if (cl->InheritsFrom(TCanvas::Class())) return kTRUE;
928  if (cl->InheritsFrom(TProfile::Class())) return kTRUE;
929  return kFALSE;
930 }
931 
932 ////////////////////////////////////////////////////////////////////////////////
933 /// scan ROOT hierarchy with provided store object
934 
935 void TRootSniffer::ScanHierarchy(const char *topname, const char *path,
936  TRootSnifferStore *store, Bool_t only_fields)
937 {
939  rec.fSearchPath = path;
940  if (rec.fSearchPath) {
941  while(*rec.fSearchPath == '/') rec.fSearchPath++;
942  if (*rec.fSearchPath == 0) rec.fSearchPath = 0;
943  }
944 
945  // if path non-empty, we should find item first and than start scanning
947  if (only_fields) rec.fMask |= TRootSnifferScanRec::kOnlyFields;
948 
949  rec.fStore = store;
950 
951  rec.CreateNode(topname);
952 
953  if (rec.fSearchPath == 0)
955 
956  if ((rec.fSearchPath == 0) && (GetAutoLoad() != 0))
958 
959  ScanRoot(rec);
960 
961  rec.CloseNode();
962 }
963 
964 ////////////////////////////////////////////////////////////////////////////////
965 /// Search element with specified path
966 /// Returns pointer on element
967 /// Optionally one could obtain element class, member description
968 /// and number of childs. When chld!=0, not only element is searched,
969 /// but also number of childs are counted. When member!=0, any object
970 /// will be scanned for its data members (disregard of extra options)
971 
972 void *TRootSniffer::FindInHierarchy(const char *path, TClass **cl,
973  TDataMember **member, Int_t *chld)
974 {
975  if (IsStreamerInfoItem(path)) {
976  // special handling for streamer info
977  CreateMemFile();
978  if (cl && fSinfo) *cl = fSinfo->IsA();
979  return fSinfo;
980  }
981 
982  TRootSnifferStore store;
983 
985  rec.fSearchPath = path;
987  if (*rec.fSearchPath == '/') rec.fSearchPath++;
988  rec.fStore = &store;
989 
990  ScanRoot(rec);
991 
992  TDataMember *res_member = store.GetResMember();
993  TClass *res_cl = store.GetResClass();
994  void *res = store.GetResPtr();
995 
996  if ((res_member!=0) && (res_cl!=0) && (member==0)) {
997  res_cl = (res_member->IsBasic() || res_member->IsSTLContainer()) ? 0 :
998  gROOT->GetClass(res_member->GetTypeName());
999  TRealData *rdata = res_cl->GetRealData(res_member->GetName());
1000  if (rdata) {
1001  res = (char *) res + rdata->GetThisOffset();
1002  if (res_member->IsaPointer()) res = *((char **) res);
1003  } else {
1004  res = 0; // should never happen
1005  }
1006  }
1007 
1008  if (cl) *cl = res_cl;
1009  if (member) *member = res_member;
1010  if (chld) *chld = store.GetResNumChilds();
1011 
1012  // remember current restriction
1013  fCurrentRestrict = store.GetResRestrict();
1014 
1015  return res;
1016 }
1017 
1018 ////////////////////////////////////////////////////////////////////////////////
1019 /// Search element in hierarchy, derived from TObject
1020 
1022 {
1023  TClass *cl(0);
1024 
1025  void *obj = FindInHierarchy(path, &cl);
1026 
1027  return (cl != 0) && (cl->GetBaseClassOffset(TObject::Class()) == 0) ? (TObject *) obj : 0;
1028 }
1029 
1030 ////////////////////////////////////////////////////////////////////////////////
1031 /// Returns hash value for streamer infos
1032 /// At the moment - just number of items in streamer infos list.
1033 
1035 {
1036  return fSinfo ? fSinfo->GetSize() : 0;
1037 }
1038 
1039 ////////////////////////////////////////////////////////////////////////////////
1040 /// Get hash function for specified item
1041 /// used to detect any changes in the specified object
1042 
1043 ULong_t TRootSniffer::GetItemHash(const char *itemname)
1044 {
1045  if (IsStreamerInfoItem(itemname)) return GetStreamerInfoHash();
1046 
1047  TObject *obj = FindTObjectInHierarchy(itemname);
1048 
1049  return obj == 0 ? 0 : TString::Hash(obj, obj->IsA()->Size());
1050 }
1051 
1052 ////////////////////////////////////////////////////////////////////////////////
1053 /// Method verifies if object can be drawn
1054 
1056 {
1057  TClass *obj_cl(0);
1058  void *res = FindInHierarchy(path, &obj_cl);
1059  return (res != 0) && IsDrawableClass(obj_cl);
1060 }
1061 
1062 ////////////////////////////////////////////////////////////////////////////////
1063 /// Method returns true when object has childs or
1064 /// one could try to expand item
1065 
1067 {
1068  TClass *obj_cl(0);
1069  Int_t obj_chld(-1);
1070  void *res = FindInHierarchy(path, &obj_cl, 0, &obj_chld);
1071  return (res != 0) && (obj_chld > 0);
1072 }
1073 
1074 ////////////////////////////////////////////////////////////////////////////////
1075 /// Creates TMemFile instance, which used for objects streaming
1076 /// One could not use TBufferFile directly,
1077 /// while one also require streamer infos list
1078 
1080 {
1081  if (fMemFile != 0) return;
1082 
1083  TDirectory *olddir = gDirectory;
1084  gDirectory = 0;
1085  TFile *oldfile = gFile;
1086  gFile = 0;
1087 
1088  fMemFile = new TMemFile("dummy.file", "RECREATE");
1089  gROOT->GetListOfFiles()->Remove(fMemFile);
1090 
1091  TH1F *d = new TH1F("d", "d", 10, 0, 10);
1092  fMemFile->WriteObject(d, "h1");
1093  delete d;
1094 
1095  TGraph *gr = new TGraph(10);
1096  gr->SetName("abc");
1097  // // gr->SetDrawOptions("AC*");
1098  fMemFile->WriteObject(gr, "gr1");
1099  delete gr;
1100 
1102 
1103  // make primary list of streamer infos
1104  TList *l = new TList();
1105 
1106  l->Add(gROOT->GetListOfStreamerInfo()->FindObject("TGraph"));
1107  l->Add(gROOT->GetListOfStreamerInfo()->FindObject("TH1F"));
1108  l->Add(gROOT->GetListOfStreamerInfo()->FindObject("TH1"));
1109  l->Add(gROOT->GetListOfStreamerInfo()->FindObject("TNamed"));
1110  l->Add(gROOT->GetListOfStreamerInfo()->FindObject("TObject"));
1111 
1112  fMemFile->WriteObject(l, "ll");
1113  delete l;
1114 
1116 
1118 
1119  gDirectory = olddir;
1120  gFile = oldfile;
1121 }
1122 
1123 ////////////////////////////////////////////////////////////////////////////////
1124 /// produce JSON data for specified item
1125 /// For object conversion TBufferJSON is used
1126 
1127 Bool_t TRootSniffer::ProduceJson(const char *path, const char *options,
1128  TString &res)
1129 {
1130  if ((path == 0) || (*path == 0)) return kFALSE;
1131 
1132  if (*path == '/') path++;
1133 
1134  TUrl url;
1135  url.SetOptions(options);
1136  url.ParseOptions();
1137  Int_t compact = -1;
1138  if (url.GetValueFromOptions("compact"))
1139  compact = url.GetIntValueFromOptions("compact");
1140 
1141  TClass *obj_cl(0);
1142  TDataMember *member(0);
1143  void *obj_ptr = FindInHierarchy(path, &obj_cl, &member);
1144  if ((obj_ptr == 0) || ((obj_cl == 0) && (member == 0))) return kFALSE;
1145 
1146  res = TBufferJSON::ConvertToJSON(obj_ptr, obj_cl, compact >= 0 ? compact : 0, member ? member->GetName() : 0);
1147 
1148  return res.Length() > 0;
1149 }
1150 
1151 ////////////////////////////////////////////////////////////////////////////////
1152 /// execute command marked as _kind=='Command'
1153 
1154 Bool_t TRootSniffer::ExecuteCmd(const char *path, const char *options,
1155  TString &res)
1156 {
1157  TFolder *parent(0);
1158  TObject *obj = GetItem(path, parent, kFALSE, kFALSE);
1159 
1160  const char *kind = GetItemField(parent, obj, item_prop_kind);
1161  if ((kind == 0) || (strcmp(kind, "Command") != 0)) {
1162  if (gDebug > 0) Info("ExecuteCmd", "Entry %s is not a command", path);
1163  res = "false";
1164  return kTRUE;
1165  }
1166 
1167  const char *cmethod = GetItemField(parent, obj, "method");
1168  if ((cmethod==0) || (strlen(cmethod)==0)) {
1169  if (gDebug > 0) Info("ExecuteCmd", "Entry %s do not defines method for execution", path);
1170  res = "false";
1171  return kTRUE;
1172  }
1173 
1174  // if read-only specified for the command, it is not allowed for execution
1175  if (fRestrictions.GetLast()>=0) {
1176  FindInHierarchy(path); // one need to call method to check access rights
1177  if (fCurrentRestrict==1) {
1178  if (gDebug > 0) Info("ExecuteCmd", "Entry %s not allowed for specified user", path);
1179  res = "false";
1180  return kTRUE;
1181  }
1182  }
1183 
1184  TString method = cmethod;
1185 
1186  const char *cnumargs = GetItemField(parent, obj, "_numargs");
1187  Int_t numargs = cnumargs ? TString(cnumargs).Atoi() : 0;
1188  if (numargs > 0) {
1189  TUrl url;
1190  url.SetOptions(options);
1191  url.ParseOptions();
1192 
1193  for (Int_t n=0; n<numargs;n++) {
1194  TString argname = TString::Format("arg%d", n+1);
1195  const char* argvalue = url.GetValueFromOptions(argname);
1196  if (argvalue==0) {
1197  if (gDebug > 0) Info("ExecuteCmd", "For command %s argument %s not specified in options %s", path, argname.Data(), options);
1198  res = "false";
1199  return kTRUE;
1200  }
1201 
1202  TString svalue = DecodeUrlOptionValue(argvalue, kTRUE);
1203  argname = TString("%") + argname + TString("%");
1204  method.ReplaceAll(argname, svalue);
1205  }
1206  }
1207 
1208  if (gDebug > 0) Info("ExecuteCmd", "Executing command %s method:%s", path, method.Data());
1209 
1210  TObject *item_obj = 0;
1211  Ssiz_t separ = method.Index("/->");
1212 
1213  if (method.Index("this->") == 0) {
1214  // if command name started with this-> means method of sniffer will be executed
1215  item_obj = this;
1216  separ = 3;
1217  } else
1218  if (separ != kNPOS) {
1219  item_obj = FindTObjectInHierarchy(TString(method.Data(), separ).Data());
1220  }
1221 
1222  if (item_obj != 0) {
1223  method = TString::Format("((%s*)%lu)->%s", item_obj->ClassName(), (long unsigned) item_obj, method.Data() + separ + 3);
1224  if (gDebug > 2) Info("ExecuteCmd", "Executing %s", method.Data());
1225  }
1226 
1227  Long_t v = gROOT->ProcessLineSync(method.Data());
1228 
1229  res.Form("%ld", v);
1230 
1231  return kTRUE;
1232 }
1233 
1234 ////////////////////////////////////////////////////////////////////////////////
1235 /// produce JSON/XML for specified item
1236 /// contrary to h.json request, only fields for specified item are stored
1237 
1238 Bool_t TRootSniffer::ProduceItem(const char *path, const char *options, TString &res, Bool_t asjson)
1239 {
1240  if (asjson) {
1241  TRootSnifferStoreJson store(res, strstr(options, "compact")!=0);
1242  ScanHierarchy("top", path, &store, kTRUE);
1243  } else {
1244  TRootSnifferStoreXml store(res, strstr(options, "compact")!=0);
1245  ScanHierarchy("top", path, &store, kTRUE);
1246  }
1247  return res.Length() > 0;
1248 }
1249 
1250 
1251 ////////////////////////////////////////////////////////////////////////////////
1252 /// produce XML data for specified item
1253 /// For object conversion TBufferXML is used
1254 
1255 Bool_t TRootSniffer::ProduceXml(const char *path, const char * /*options*/,
1256  TString &res)
1257 {
1258  if ((path == 0) || (*path == 0)) return kFALSE;
1259 
1260  if (*path == '/') path++;
1261 
1262  TClass *obj_cl(0);
1263  void *obj_ptr = FindInHierarchy(path, &obj_cl);
1264  if ((obj_ptr == 0) || (obj_cl == 0)) return kFALSE;
1265 
1266  res = TBufferXML::ConvertToXML(obj_ptr, obj_cl);
1267 
1268  return res.Length() > 0;
1269 }
1270 
1271 ////////////////////////////////////////////////////////////////////////////////
1272 /// method replaces all kind of special symbols, which could appear in URL options
1273 
1275 {
1276  if ((value == 0) || (strlen(value) == 0)) return TString();
1277 
1278  TString res = value;
1279 
1280  res.ReplaceAll("%27", "\'");
1281  res.ReplaceAll("%22", "\"");
1282  res.ReplaceAll("%3E", ">");
1283  res.ReplaceAll("%3C", "<");
1284  res.ReplaceAll("%20", " ");
1285  res.ReplaceAll("%5B", "[");
1286  res.ReplaceAll("%5D", "]");
1287 
1288  if (remove_quotes && (res.Length() > 1) &&
1289  ((res[0] == '\'') || (res[0] == '\"')) && (res[0] == res[res.Length() - 1])) {
1290  res.Remove(res.Length() - 1);
1291  res.Remove(0, 1);
1292  }
1293 
1294  return res;
1295 }
1296 
1297 ////////////////////////////////////////////////////////////////////////////////
1298 /// execute command for specified object
1299 /// options include method and extra list of parameters
1300 /// sniffer should be not-readonly to allow execution of the commands
1301 /// reskind defines kind of result 0 - debug, 1 - json, 2 - binary
1302 
1303 Bool_t TRootSniffer::ProduceExe(const char *path, const char *options, Int_t reskind, TString *res_str, void **res_ptr, Long_t *res_length)
1304 {
1305  TString *debug = (reskind == 0) ? res_str : 0;
1306 
1307  if ((path == 0) || (*path == 0)) {
1308  if (debug) debug->Append("Item name not specified\n");
1309  return debug != 0;
1310  }
1311 
1312  if (*path == '/') path++;
1313 
1314  TClass *obj_cl(0);
1315  void *obj_ptr = FindInHierarchy(path, &obj_cl);
1316  if (debug) debug->Append(TString::Format("Item:%s found:%s\n", path, obj_ptr ? "true" : "false"));
1317  if ((obj_ptr == 0) || (obj_cl == 0)) return debug != 0;
1318 
1319  TUrl url;
1320  url.SetOptions(options);
1321 
1322  const char *method_name = url.GetValueFromOptions("method");
1323  TString prototype = DecodeUrlOptionValue(url.GetValueFromOptions("prototype"), kTRUE);
1324  TString funcname = DecodeUrlOptionValue(url.GetValueFromOptions("func"), kTRUE);
1325  TMethod *method = 0;
1326  TFunction *func = 0;
1327  if (method_name != 0) {
1328  if (prototype.Length() == 0) {
1329  if (debug) debug->Append(TString::Format("Search for any method with name \'%s\'\n", method_name));
1330  method = obj_cl->GetMethodAllAny(method_name);
1331  } else {
1332  if (debug) debug->Append(TString::Format("Search for method \'%s\' with prototype \'%s\'\n", method_name, prototype.Data()));
1333  method = obj_cl->GetMethodWithPrototype(method_name, prototype);
1334  }
1335  }
1336 
1337  if (method != 0) {
1338  if (debug) debug->Append(TString::Format("Method: %s\n", method->GetPrototype()));
1339  } else {
1340  if (funcname.Length() > 0) {
1341  if (prototype.Length() == 0) {
1342  if (debug) debug->Append(TString::Format("Search for any function with name \'%s\'\n", funcname.Data()));
1343  func = gROOT->GetGlobalFunction(funcname);
1344  } else {
1345  if (debug) debug->Append(TString::Format("Search for function \'%s\' with prototype \'%s\'\n", funcname.Data(), prototype.Data()));
1346  func = gROOT->GetGlobalFunctionWithPrototype(funcname, prototype);
1347  }
1348  }
1349 
1350  if (func != 0) {
1351  if (debug) debug->Append(TString::Format("Function: %s\n", func->GetPrototype()));
1352  }
1353  }
1354 
1355  if ((method == 0) && (func==0)) {
1356  if (debug) debug->Append("Method not found\n");
1357  return debug != 0;
1358  }
1359 
1360  if ((fReadOnly && (fCurrentRestrict == 0)) || (fCurrentRestrict == 1)) {
1361  if ((method!=0) && (fCurrentAllowedMethods.Index(method_name) == kNPOS)) {
1362  if (debug) debug->Append("Server runs in read-only mode, method cannot be executed\n");
1363  return debug != 0;
1364  } else
1365  if ((func!=0) && (fCurrentAllowedMethods.Index(funcname) == kNPOS)) {
1366  if (debug) debug->Append("Server runs in read-only mode, function cannot be executed\n");
1367  return debug != 0;
1368  } else {
1369  if (debug) debug->Append("For that special method server allows access even read-only mode is specified\n");
1370  }
1371  }
1372 
1373  TList *args = method ? method->GetListOfMethodArgs() : func->GetListOfMethodArgs();
1374 
1375  TList garbage;
1376  garbage.SetOwner(kTRUE); // use as garbage collection
1377  TObject *post_obj = 0; // object reconstructed from post request
1378  TString call_args;
1379 
1380  TIter next(args);
1381  TMethodArg *arg = 0;
1382  while ((arg = (TMethodArg *) next()) != 0) {
1383 
1384  if ((strcmp(arg->GetName(), "rest_url_opt") == 0) &&
1385  (strcmp(arg->GetFullTypeName(), "const char*") == 0) && (args->GetSize() == 1)) {
1386  // very special case - function requires list of options after method=argument
1387 
1388  const char *pos = strstr(options, "method=");
1389  if ((pos == 0) || (strlen(pos) < strlen(method_name) + 8)) return debug != 0;
1390  call_args.Form("\"%s\"", pos + strlen(method_name) + 8);
1391  break;
1392  }
1393 
1394  TString sval;
1395  const char *val = url.GetValueFromOptions(arg->GetName());
1396  if (val) {
1397  sval = DecodeUrlOptionValue(val, kFALSE);
1398  val = sval.Data();
1399  }
1400 
1401  if ((val!=0) && (strcmp(val,"_this_")==0)) {
1402  // special case - object itself is used as argument
1403  sval.Form("(%s*)0x%lx", obj_cl->GetName(), (long unsigned) obj_ptr);
1404  val = sval.Data();
1405  } else
1406  if ((val!=0) && (fCurrentArg!=0) && (fCurrentArg->GetPostData()!=0)) {
1407  // process several arguments which are specific for post requests
1408  if (strcmp(val,"_post_object_xml_")==0) {
1409  // post data has extra 0 at the end and can be used as null-terminated string
1410  post_obj = TBufferXML::ConvertFromXML((const char*) fCurrentArg->GetPostData());
1411  if (post_obj == 0) {
1412  sval = "0";
1413  } else {
1414  sval.Form("(%s*)0x%lx", post_obj->ClassName(), (long unsigned) post_obj);
1415  if (url.HasOption("_destroy_post_")) garbage.Add(post_obj);
1416  }
1417  val = sval.Data();
1418  } else
1419  if ((strcmp(val,"_post_object_")==0) && url.HasOption("_post_class_")) {
1420  TString clname = url.GetValueFromOptions("_post_class_");
1421  TClass* arg_cl = gROOT->GetClass(clname, kTRUE, kTRUE);
1422  if ((arg_cl!=0) && (arg_cl->GetBaseClassOffset(TObject::Class()) == 0) && (post_obj==0)) {
1423  post_obj = (TObject*) arg_cl->New();
1424  if (post_obj==0) {
1425  if (debug) debug->Append(TString::Format("Fail to create object of class %s\n", clname.Data()));
1426  } else {
1427  if (debug) debug->Append(TString::Format("Reconstruct object of class %s from POST data\n", clname.Data()));
1429  buf.MapObject(post_obj, arg_cl);
1430  post_obj->Streamer(buf);
1431  if (url.HasOption("_destroy_post_")) garbage.Add(post_obj);
1432  }
1433  }
1434  sval.Form("(%s*)0x%lx", clname.Data(), (long unsigned) post_obj);
1435  val = sval.Data();
1436  } else
1437  if (strcmp(val,"_post_data_")==0) {
1438  sval.Form("(void*)0x%lx", (long unsigned) *res_ptr);
1439  val = sval.Data();
1440  } else
1441  if (strcmp(val,"_post_length_")==0) {
1442  sval.Form("%ld", (long) *res_length);
1443  val = sval.Data();
1444  }
1445  }
1446 
1447  if (val == 0) val = arg->GetDefault();
1448 
1449  if (debug) debug->Append(TString::Format(" Argument:%s Type:%s Value:%s \n", arg->GetName(), arg->GetFullTypeName(), val ? val : "<missed>"));
1450  if (val == 0) return debug != 0;
1451 
1452  if (call_args.Length() > 0) call_args += ", ";
1453 
1454  if ((strcmp(arg->GetFullTypeName(), "const char*") == 0) || (strcmp(arg->GetFullTypeName(), "Option_t*") == 0)) {
1455  int len = strlen(val);
1456  if ((strlen(val) < 2) || (*val != '\"') || (val[len - 1] != '\"'))
1457  call_args.Append(TString::Format("\"%s\"", val));
1458  else
1459  call_args.Append(val);
1460  } else {
1461  call_args.Append(val);
1462  }
1463  }
1464 
1465  TMethodCall *call = 0;
1466 
1467  if (method!=0) {
1468  call = new TMethodCall(obj_cl, method_name, call_args.Data());
1469  if (debug) debug->Append(TString::Format("Calling obj->%s(%s);\n", method_name, call_args.Data()));
1470 
1471  } else {
1472  call = new TMethodCall(funcname.Data(), call_args.Data());
1473  if (debug) debug->Append(TString::Format("Calling %s(%s);\n", funcname.Data(), call_args.Data()));
1474  }
1475 
1476  garbage.Add(call);
1477 
1478  if (!call->IsValid()) {
1479  if (debug) debug->Append("Fail: invalid TMethodCall\n");
1480  return debug != 0;
1481  }
1482 
1483  Int_t compact = 0;
1484  if (url.GetValueFromOptions("compact"))
1485  compact = url.GetIntValueFromOptions("compact");
1486 
1487  TString res = "null";
1488  void *ret_obj = 0;
1489  TClass *ret_cl = 0;
1490  TBufferFile *resbuf = 0;
1491  if (reskind==2) {
1492  resbuf = new TBufferFile(TBuffer::kWrite, 10000);
1493  garbage.Add(resbuf);
1494  }
1495 
1496  switch (call->ReturnType()) {
1497  case TMethodCall::kLong: {
1498  Long_t l(0);
1499  if (method)
1500  call->Execute(obj_ptr, l);
1501  else
1502  call->Execute(l);
1503  if (resbuf) resbuf->WriteLong(l);
1504  else res.Form("%ld", l);
1505  break;
1506  }
1507  case TMethodCall::kDouble : {
1508  Double_t d(0.);
1509  if (method)
1510  call->Execute(obj_ptr, d);
1511  else
1512  call->Execute(d);
1513  if (resbuf) resbuf->WriteDouble(d);
1514  else res.Form(TBufferJSON::GetFloatFormat(), d);
1515  break;
1516  }
1517  case TMethodCall::kString : {
1518  char *txt(0);
1519  if (method)
1520  call->Execute(obj_ptr, &txt);
1521  else
1522  call->Execute(0, &txt); // here 0 is artificial, there is no proper signature
1523  if (txt != 0) {
1524  if (resbuf) resbuf->WriteString(txt);
1525  else res.Form("\"%s\"", txt);
1526  }
1527  break;
1528  }
1529  case TMethodCall::kOther : {
1530  std::string ret_kind = func ? func->GetReturnTypeNormalizedName() : method->GetReturnTypeNormalizedName();
1531  if ((ret_kind.length() > 0) && (ret_kind[ret_kind.length() - 1] == '*')) {
1532  ret_kind.resize(ret_kind.length() - 1);
1533  ret_cl = gROOT->GetClass(ret_kind.c_str(), kTRUE, kTRUE);
1534  }
1535 
1536  if (ret_cl != 0) {
1537  Long_t l(0);
1538  if (method)
1539  call->Execute(obj_ptr, l);
1540  else
1541  call->Execute(l);
1542  if (l != 0) ret_obj = (void *) l;
1543  } else {
1544  if (method)
1545  call->Execute(obj_ptr);
1546  else
1547  call->Execute();
1548  }
1549 
1550  break;
1551  }
1552  case TMethodCall::kNone : {
1553  if (method)
1554  call->Execute(obj_ptr);
1555  else
1556  call->Execute();
1557  break;
1558  }
1559  }
1560 
1561  const char *_ret_object_ = url.GetValueFromOptions("_ret_object_");
1562  if (_ret_object_ != 0) {
1563  TObject *obj = 0;
1564  if (gDirectory != 0) obj = gDirectory->Get(_ret_object_);
1565  if (debug) debug->Append(TString::Format("Return object %s found %s\n", _ret_object_, obj ? "true" : "false"));
1566 
1567  if (obj == 0) {
1568  res = "null";
1569  } else {
1570  ret_obj = obj;
1571  ret_cl = obj->IsA();
1572  }
1573  }
1574 
1575  if ((ret_obj != 0) && (ret_cl != 0)) {
1576  if ((resbuf != 0) && (ret_cl->GetBaseClassOffset(TObject::Class()) == 0)) {
1577  TObject *obj = (TObject *) ret_obj;
1578  resbuf->MapObject(obj);
1579  obj->Streamer(*resbuf);
1580  if (fCurrentArg) fCurrentArg->SetExtraHeader("RootClassName", ret_cl->GetName());
1581  } else {
1582  res = TBufferJSON::ConvertToJSON(ret_obj, ret_cl, compact);
1583  }
1584  }
1585 
1586  if ((resbuf != 0) && (resbuf->Length() > 0) && (res_ptr != 0) && (res_length != 0)) {
1587  *res_ptr = malloc(resbuf->Length());
1588  memcpy(*res_ptr, resbuf->Buffer(), resbuf->Length());
1589  *res_length = resbuf->Length();
1590  }
1591 
1592  if (debug) debug->Append(TString::Format("Result = %s\n", res.Data()));
1593 
1594  if ((reskind == 1) && res_str) *res_str = res;
1595 
1596  if (url.HasOption("_destroy_result_") && (ret_obj != 0) && (ret_cl != 0)) {
1597  ret_cl->Destructor(ret_obj);
1598  if (debug) debug->Append("Destroy result object at the end\n");
1599  }
1600 
1601  // delete all garbage objects, but should be also done with any return
1602  garbage.Delete();
1603 
1604  return kTRUE;
1605 }
1606 
1607 ////////////////////////////////////////////////////////////////////////////////
1608 /// Process several requests, packing all results into binary or JSON buffer
1609 /// Input parameters should be coded in the POST block and has
1610 /// individual request relative to current path, separated with '\n' symbol like
1611 /// item1/root.bin\n
1612 /// item2/exe.bin?method=GetList\n
1613 /// item3/exe.bin?method=GetTitle\n
1614 /// Request requires 'number' URL option which contains number of requested items
1615 ///
1616 /// In case of binary request output buffer looks like:
1617 /// 4bytes length + payload, 4bytes length + payload, ...
1618 /// In case of JSON request output is array with results for each item
1619 /// multi.json request do not support binary requests for the items
1620 
1621 Bool_t TRootSniffer::ProduceMulti(const char *path, const char *options, void *&ptr, Long_t &length, TString &str, Bool_t asjson)
1622 {
1623  if ((fCurrentArg==0) || (fCurrentArg->GetPostDataLength()<=0) || (fCurrentArg->GetPostData()==0)) return kFALSE;
1624 
1625  const char* args = (const char*) fCurrentArg->GetPostData();
1626  const char* ends = args + fCurrentArg->GetPostDataLength();
1627 
1628  TUrl url;
1629  url.SetOptions(options);
1630 
1631  Int_t number = 0;
1632  if (url.GetValueFromOptions("number"))
1633  number = url.GetIntValueFromOptions("number");
1634 
1635  // binary buffers required only for binary requests, json output can be produced as is
1636  std::vector<void*> mem;
1637  std::vector<Long_t> memlen;
1638 
1639  if (asjson) str = "[";
1640 
1641  for (Int_t n=0;n<number;n++) {
1642  const char* next = args;
1643  while ((next < ends) && (*next != '\n')) next++;
1644  if (next==ends) {
1645  Error("ProduceMulti", "Not enough arguments in POST block");
1646  break;
1647  }
1648 
1649  TString file1(args, next - args);
1650  args = next + 1;
1651 
1652  TString path1, opt1;
1653 
1654  // extract options
1655  Int_t pos = file1.First('?');
1656  if (pos != kNPOS) {
1657  opt1 = file1(pos+1, file1.Length() - pos);
1658  file1.Resize(pos);
1659  }
1660 
1661  // extract extra path
1662  pos = file1.Last('/');
1663  if (pos != kNPOS) {
1664  path1 = file1(0, pos);
1665  file1.Remove(0, pos+1);
1666  }
1667 
1668  if ((path!=0) && (*path!=0)) path1 = TString(path) + "/" + path1;
1669 
1670  void* ptr1 = 0;
1671  Long_t len1 = 0;
1672  TString str1;
1673 
1674  // produce next item request
1675  Produce(path1, file1, opt1, ptr1, len1, str1);
1676 
1677  if (asjson) {
1678  if (n>0) str.Append(", ");
1679  if (ptr1!=0) { str.Append("\"<non-supported binary>\""); free(ptr1); } else
1680  if (str1.Length()>0) str.Append(str1); else str.Append("null");
1681  } else {
1682  if ((str1.Length()>0) && (ptr1==0)) {
1683  len1 = str1.Length();
1684  ptr1 = malloc(len1);
1685  memcpy(ptr1, str1.Data(), len1);
1686  }
1687  mem.push_back(ptr1);
1688  memlen.push_back(len1);
1689  }
1690  }
1691 
1692  if (asjson) {
1693  str.Append("]");
1694  } else {
1695  length = 0;
1696  for (unsigned n=0;n<mem.size();n++) {
1697  length += 4 + memlen[n];
1698  }
1699  ptr = malloc(length);
1700  char* curr = (char*) ptr;
1701  for (unsigned n=0;n<mem.size();n++) {
1702  Long_t l = memlen[n];
1703  *curr++ = (char) (l & 0xff); l = l >> 8;
1704  *curr++ = (char) (l & 0xff); l = l >> 8;
1705  *curr++ = (char) (l & 0xff); l = l >> 8;
1706  *curr++ = (char) (l & 0xff);
1707  if ((mem[n]!=0) && (memlen[n]>0)) memcpy(curr, mem[n], memlen[n]);
1708  curr+=memlen[n];
1709  }
1710  }
1711 
1712  for (unsigned n=0;n<mem.size();n++) free(mem[n]);
1713 
1714  return kTRUE;
1715 }
1716 
1717 
1718 ////////////////////////////////////////////////////////////////////////////////
1719 /// Return true if it is streamer info item name
1720 
1722 {
1723  if ((itemname == 0) || (*itemname == 0)) return kFALSE;
1724 
1725  return (strcmp(itemname, "StreamerInfo") == 0) || (strcmp(itemname, "StreamerInfo/") == 0);
1726 }
1727 
1728 ////////////////////////////////////////////////////////////////////////////////
1729 /// produce binary data for specified item
1730 /// if "zipped" option specified in query, buffer will be compressed
1731 
1732 Bool_t TRootSniffer::ProduceBinary(const char *path, const char * /*query*/, void *&ptr,
1733  Long_t &length)
1734 {
1735  if ((path == 0) || (*path == 0)) return kFALSE;
1736 
1737  if (*path == '/') path++;
1738 
1739  TClass *obj_cl(0);
1740  void *obj_ptr = FindInHierarchy(path, &obj_cl);
1741  if ((obj_ptr == 0) || (obj_cl == 0)) return kFALSE;
1742 
1743  if (obj_cl->GetBaseClassOffset(TObject::Class()) != 0) {
1744  Info("ProduceBinary", "Non-TObject class not supported");
1745  return kFALSE;
1746  }
1747 
1748  // ensure that memfile exists
1749  CreateMemFile();
1750 
1751  TDirectory *olddir = gDirectory;
1752  gDirectory = 0;
1753  TFile *oldfile = gFile;
1754  gFile = 0;
1755 
1756  TObject *obj = (TObject *) obj_ptr;
1757 
1758  TBufferFile *sbuf = new TBufferFile(TBuffer::kWrite, 100000);
1759  sbuf->SetParent(fMemFile);
1760  sbuf->MapObject(obj);
1761  obj->Streamer(*sbuf);
1762  if (fCurrentArg) fCurrentArg->SetExtraHeader("RootClassName", obj_cl->GetName());
1763 
1764  // produce actual version of streamer info
1765  delete fSinfo;
1768 
1769  gDirectory = olddir;
1770  gFile = oldfile;
1771 
1772  ptr = malloc(sbuf->Length());
1773  memcpy(ptr, sbuf->Buffer(), sbuf->Length());
1774  length = sbuf->Length();
1775 
1776  delete sbuf;
1777 
1778  return kTRUE;
1779 }
1780 
1781 ////////////////////////////////////////////////////////////////////////////////
1782 /// Method to produce image from specified object
1783 ///
1784 /// Parameters:
1785 /// kind - image kind TImage::kPng, TImage::kJpeg, TImage::kGif
1786 /// path - path to object
1787 /// options - extra options
1788 ///
1789 /// By default, image 300x200 is produced
1790 /// In options string one could provide following parameters:
1791 /// w - image width
1792 /// h - image height
1793 /// opt - draw options
1794 /// For instance:
1795 /// http://localhost:8080/Files/hsimple.root/hpx/get.png?w=500&h=500&opt=lego1
1796 ///
1797 /// Return is memory with produced image
1798 /// Memory must be released by user with free(ptr) call
1799 
1800 Bool_t TRootSniffer::ProduceImage(Int_t kind, const char *path,
1801  const char *options, void *&ptr,
1802  Long_t &length)
1803 {
1804  ptr = 0;
1805  length = 0;
1806 
1807  if ((path == 0) || (*path == 0)) return kFALSE;
1808  if (*path == '/') path++;
1809 
1810  TClass *obj_cl(0);
1811  void *obj_ptr = FindInHierarchy(path, &obj_cl);
1812  if ((obj_ptr == 0) || (obj_cl == 0)) return kFALSE;
1813 
1814  if (obj_cl->GetBaseClassOffset(TObject::Class()) != 0) {
1815  Error("TRootSniffer", "Only derived from TObject classes can be drawn");
1816  return kFALSE;
1817  }
1818 
1819  TObject *obj = (TObject *) obj_ptr;
1820 
1821  TImage *img = TImage::Create();
1822  if (img == 0) return kFALSE;
1823 
1824  if (obj->InheritsFrom(TPad::Class())) {
1825 
1826  if (gDebug > 1)
1827  Info("TRootSniffer", "Crate IMAGE directly from pad");
1828  img->FromPad((TPad *) obj);
1829  } else if (IsDrawableClass(obj->IsA())) {
1830 
1831  if (gDebug > 1)
1832  Info("TRootSniffer", "Crate IMAGE from object %s", obj->GetName());
1833 
1834  Int_t width(300), height(200);
1835  TString drawopt = "";
1836 
1837  if ((options != 0) && (*options != 0)) {
1838  TUrl url;
1839  url.SetOptions(options);
1840  url.ParseOptions();
1841  Int_t w = url.GetIntValueFromOptions("w");
1842  if (w > 10) width = w;
1843  Int_t h = url.GetIntValueFromOptions("h");
1844  if (h > 10) height = h;
1845  const char *opt = url.GetValueFromOptions("opt");
1846  if (opt != 0) drawopt = opt;
1847  }
1848 
1849  Bool_t isbatch = gROOT->IsBatch();
1850  TVirtualPad *save_gPad = gPad;
1851 
1852  if (!isbatch) gROOT->SetBatch(kTRUE);
1853 
1854  TCanvas *c1 = new TCanvas("__online_draw_canvas__", "title", width, height);
1855  obj->Draw(drawopt.Data());
1856  img->FromPad(c1);
1857  delete c1;
1858 
1859  if (!isbatch) gROOT->SetBatch(kFALSE);
1860  gPad = save_gPad;
1861 
1862  } else {
1863  delete img;
1864  return kFALSE;
1865  }
1866 
1867  TImage *im = TImage::Create();
1868  im->Append(img);
1869 
1870  char *png_buffer(0);
1871  int size(0);
1872 
1873  im->GetImageBuffer(&png_buffer, &size, (TImage::EImageFileTypes) kind);
1874 
1875  if ((png_buffer != 0) && (size > 0)) {
1876  ptr = malloc(size);
1877  length = size;
1878  memcpy(ptr, png_buffer, length);
1879  }
1880 
1881  delete [] png_buffer;
1882  delete im;
1883 
1884  return ptr != 0;
1885 }
1886 
1887 ////////////////////////////////////////////////////////////////////////////////
1888 /// Method produce different kind of data out of object
1889 /// Parameter 'path' specifies object or object member
1890 /// Supported 'file' (case sensitive):
1891 /// "root.bin" - binary data
1892 /// "root.png" - png image
1893 /// "root.jpeg" - jpeg image
1894 /// "root.gif" - gif image
1895 /// "root.xml" - xml representation
1896 /// "root.json" - json representation
1897 /// "exe.json" - method execution with json reply
1898 /// "exe.bin" - method execution with binary reply
1899 /// "exe.txt" - method execution with debug output
1900 /// "cmd.json" - execution of registered commands
1901 /// Result returned either as string or binary buffer,
1902 /// which should be released with free() call
1903 
1904 Bool_t TRootSniffer::Produce(const char *path, const char *file,
1905  const char *options, void *&ptr, Long_t &length, TString &str)
1906 {
1907  if ((file == 0) || (*file == 0)) return kFALSE;
1908 
1909  if (strcmp(file, "root.bin") == 0)
1910  return ProduceBinary(path, options, ptr, length);
1911 
1912  if (strcmp(file, "root.png") == 0)
1913  return ProduceImage(TImage::kPng, path, options, ptr, length);
1914 
1915  if (strcmp(file, "root.jpeg") == 0)
1916  return ProduceImage(TImage::kJpeg, path, options, ptr, length);
1917 
1918  if (strcmp(file, "root.gif") == 0)
1919  return ProduceImage(TImage::kGif, path, options, ptr, length);
1920 
1921  if (strcmp(file, "exe.bin") == 0)
1922  return ProduceExe(path, options, 2, 0, &ptr, &length);
1923 
1924  if (strcmp(file, "root.xml") == 0)
1925  return ProduceXml(path, options, str);
1926 
1927  if (strcmp(file, "root.json") == 0)
1928  return ProduceJson(path, options, str);
1929 
1930  // used for debugging
1931  if (strcmp(file, "exe.txt") == 0)
1932  return ProduceExe(path, options, 0, &str);
1933 
1934  if (strcmp(file, "exe.json") == 0)
1935  return ProduceExe(path, options, 1, &str);
1936 
1937  if (strcmp(file, "cmd.json") == 0)
1938  return ExecuteCmd(path, options, str);
1939 
1940  if (strcmp(file, "item.json") == 0)
1941  return ProduceItem(path, options, str, kTRUE);
1942 
1943  if (strcmp(file, "item.xml") == 0)
1944  return ProduceItem(path, options, str, kFALSE);
1945 
1946  if (strcmp(file, "multi.bin") == 0)
1947  return ProduceMulti(path, options, ptr, length, str, kFALSE);
1948 
1949  if (strcmp(file, "multi.json") == 0)
1950  return ProduceMulti(path, options, ptr, length, str, kTRUE);
1951 
1952  return kFALSE;
1953 }
1954 
1955 ////////////////////////////////////////////////////////////////////////////////
1956 /// return item from the subfolders structure
1957 
1958 TObject *TRootSniffer::GetItem(const char *fullname, TFolder *&parent, Bool_t force, Bool_t within_objects)
1959 {
1960  TFolder *topf = gROOT->GetRootFolder();
1961 
1962  if (topf == 0) {
1963  Error("RegisterObject", "Not found top ROOT folder!!!");
1964  return 0;
1965  }
1966 
1967  TFolder *httpfold = dynamic_cast<TFolder *>(topf->FindObject("http"));
1968  if (httpfold == 0) {
1969  if (!force) return 0;
1970  httpfold = topf->AddFolder("http", "ROOT http server");
1971  httpfold->SetBit(kCanDelete);
1972  // register top folder in list of cleanups
1973  gROOT->GetListOfCleanups()->Add(httpfold);
1974  }
1975 
1976  parent = httpfold;
1977  TObject *obj = httpfold;
1978 
1979  if (fullname==0) return httpfold;
1980 
1981  // when full path started not with slash, "Objects" subfolder is appended
1982  TString path = fullname;
1983  if (within_objects && ((path.Length()==0) || (path[0]!='/')))
1984  path = fObjectsPath + "/" + path;
1985 
1986  TString tok;
1987  Ssiz_t from(0);
1988 
1989  while (path.Tokenize(tok,from,"/")) {
1990  if (tok.Length()==0) continue;
1991 
1992  TFolder *fold = dynamic_cast<TFolder *> (obj);
1993  if (fold == 0) return 0;
1994 
1995  TIter iter(fold->GetListOfFolders());
1996  while ((obj = iter()) != 0) {
1997  if (IsItemField(obj)) continue;
1998  if (tok.CompareTo(obj->GetName())==0) break;
1999  }
2000 
2001  if (obj == 0) {
2002  if (!force) return 0;
2003  obj = fold->AddFolder(tok, "sub-folder");
2004  obj->SetBit(kCanDelete);
2005  }
2006 
2007  parent = fold;
2008  }
2009 
2010  return obj;
2011 }
2012 
2013 ////////////////////////////////////////////////////////////////////////////////
2014 /// creates subfolder where objects can be registered
2015 
2016 TFolder *TRootSniffer::GetSubFolder(const char *subfolder, Bool_t force)
2017 {
2018  TFolder *parent = 0;
2019 
2020  return dynamic_cast<TFolder *> (GetItem(subfolder, parent, force));
2021 }
2022 
2023 ////////////////////////////////////////////////////////////////////////////////
2024 /// Register object in subfolder structure
2025 /// subfolder parameter can have many levels like:
2026 ///
2027 /// TRootSniffer* sniff = new TRootSniffer("sniff");
2028 /// sniff->RegisterObject("my/sub/subfolder", h1);
2029 ///
2030 /// Such objects can be later found in "Objects" folder of sniffer like
2031 ///
2032 /// h1 = sniff->FindTObjectInHierarchy("/Objects/my/sub/subfolder/h1");
2033 ///
2034 /// If subfolder name starts with '/', object will be registered starting from top folder.
2035 ///
2036 /// One could provide additional fields for registered objects
2037 /// For instance, setting "_more" field to true let browser
2038 /// explore objects members. For instance:
2039 ///
2040 /// TEvent* ev = new TEvent("ev");
2041 /// sniff->RegisterObject("Events", ev);
2042 /// sniff->SetItemField("Events/ev", "_more", "true");
2043 
2045 {
2046  TFolder *f = GetSubFolder(subfolder, kTRUE);
2047  if (f == 0) return kFALSE;
2048 
2049  // If object will be destroyed, it will be removed from the folders automatically
2050  obj->SetBit(kMustCleanup);
2051 
2052  f->Add(obj);
2053 
2054  return kTRUE;
2055 }
2056 
2057 ////////////////////////////////////////////////////////////////////////////////
2058 /// unregister (remove) object from folders structures
2059 /// folder itself will remain even when it will be empty
2060 
2062 {
2063  if (obj == 0) return kTRUE;
2064 
2065  TFolder *topf = dynamic_cast<TFolder *>(gROOT->FindObject("//root/http"));
2066 
2067  if (topf == 0) {
2068  Error("UnregisterObject", "Not found //root/http folder!!!");
2069  return kFALSE;
2070  }
2071 
2072  // TODO - probably we should remove all set properties as well
2073  if (topf) topf->RecursiveRemove(obj);
2074 
2075  return kTRUE;
2076 }
2077 
2078 ////////////////////////////////////////////////////////////////////////////////
2079 /// create item element
2080 
2081 Bool_t TRootSniffer::CreateItem(const char *fullname, const char *title)
2082 {
2083  TFolder *f = GetSubFolder(fullname, kTRUE);
2084  if (f == 0) return kFALSE;
2085 
2086  if (title) f->SetTitle(title);
2087 
2088  return kTRUE;
2089 }
2090 
2091 ////////////////////////////////////////////////////////////////////////////////
2092 /// return true when object is TNamed with kItemField bit set
2093 /// such objects used to keep field values for item
2094 
2096 {
2097  return (obj!=0) && (obj->IsA() == TNamed::Class()) && obj->TestBit(kItemField);
2098 }
2099 
2100 ////////////////////////////////////////////////////////////////////////////////
2101 /// set or get field for the child
2102 /// each field coded as TNamed object, placed after chld in the parent hierarchy
2103 
2105  const char *name, const char *value, TNamed **only_get)
2106 {
2107  if (parent==0) return kFALSE;
2108 
2109  if (chld==0) {
2110  Info("SetField", "Should be special case for top folder, support later");
2111  return kFALSE;
2112  }
2113 
2114  TIter iter(parent->GetListOfFolders());
2115 
2116  TObject* obj = 0;
2117  Bool_t find(kFALSE), last_find(kFALSE);
2118  // this is special case of top folder - fields are on very top
2119  if (parent == chld) { last_find = find = kTRUE; }
2120  TNamed* curr = 0;
2121  while ((obj = iter()) != 0) {
2122  if (IsItemField(obj)) {
2123  if (last_find && (obj->GetName()!=0) && !strcmp(name, obj->GetName())) curr = (TNamed*) obj;
2124  } else {
2125  last_find = (obj == chld);
2126  if (last_find) find = kTRUE;
2127  if (find && !last_find) break; // no need to continue
2128  }
2129  }
2130 
2131  // object must be in childs list
2132  if (!find) return kFALSE;
2133 
2134  if (only_get!=0) {
2135  *only_get = curr;
2136  return curr!=0;
2137  }
2138 
2139  if (curr!=0) {
2140  if (value!=0) { curr->SetTitle(value); }
2141  else { parent->Remove(curr); delete curr; }
2142  return kTRUE;
2143  }
2144 
2145  curr = new TNamed(name, value);
2146  curr->SetBit(kItemField);
2147 
2148  if (last_find) {
2149  // object is on last place, therefore just add property
2150  parent->Add(curr);
2151  return kTRUE;
2152  }
2153 
2154  // only here we do dynamic cast to the TList to use AddAfter
2155  TList *lst = dynamic_cast<TList *> (parent->GetListOfFolders());
2156  if (lst==0) {
2157  Error("SetField", "Fail cast to TList");
2158  return kFALSE;
2159  }
2160 
2161  if (parent==chld)
2162  lst->AddFirst(curr);
2163  else
2164  lst->AddAfter(chld, curr);
2165 
2166  return kTRUE;
2167 }
2168 
2169 ////////////////////////////////////////////////////////////////////////////////
2170 /// set field for specified item
2171 
2172 Bool_t TRootSniffer::SetItemField(const char *fullname, const char *name, const char *value)
2173 {
2174  if ((fullname==0) || (name==0)) return kFALSE;
2175 
2176  TFolder *parent(0);
2177  TObject *obj = GetItem(fullname, parent);
2178 
2179  if ((parent==0) || (obj==0)) return kFALSE;
2180 
2181  if (strcmp(name, item_prop_title)==0) {
2182  TNamed* n = dynamic_cast<TNamed*> (obj);
2183  if (n!=0) { n->SetTitle(value); return kTRUE; }
2184  }
2185 
2186  return AccessField(parent, obj, name, value);
2187 }
2188 
2189 ////////////////////////////////////////////////////////////////////////////////
2190 /// return field for specified item
2191 
2192 const char *TRootSniffer::GetItemField(TFolder *parent, TObject *obj, const char *name)
2193 {
2194  if ((parent==0) || (obj==0) || (name==0)) return 0;
2195 
2196  TNamed *field(0);
2197 
2198  if (!AccessField(parent, obj, name, 0, &field)) return 0;
2199 
2200  return field ? field->GetTitle() : 0;
2201 }
2202 
2203 
2204 ////////////////////////////////////////////////////////////////////////////////
2205 /// return field for specified item
2206 
2207 const char *TRootSniffer::GetItemField(const char *fullname, const char *name)
2208 {
2209  if (fullname==0) return 0;
2210 
2211  TFolder *parent(0);
2212  TObject *obj = GetItem(fullname, parent);
2213 
2214  return GetItemField(parent, obj, name);
2215 }
2216 
2217 
2218 ////////////////////////////////////////////////////////////////////////////////
2219 /// Register command which can be executed from web interface
2220 ///
2221 /// As method one typically specifies string, which is executed with
2222 /// gROOT->ProcessLine() method. For instance
2223 /// serv->RegisterCommand("Invoke","InvokeFunction()");
2224 ///
2225 /// Or one could specify any method of the object which is already registered
2226 /// to the server. For instance:
2227 /// serv->Register("/", hpx);
2228 /// serv->RegisterCommand("/ResetHPX", "/hpx/->Reset()");
2229 /// Here symbols '/->' separates item name from method to be executed
2230 ///
2231 /// One could specify additional arguments in the command with
2232 /// syntax like %arg1%, %arg2% and so on. For example:
2233 /// serv->RegisterCommand("/ResetHPX", "/hpx/->SetTitle(\"%arg1%\")");
2234 /// serv->RegisterCommand("/RebinHPXPY", "/hpxpy/->Rebin2D(%arg1%,%arg2%)");
2235 /// Such parameter(s) will be requested when command clicked in the browser.
2236 ///
2237 /// Once command is registered, one could specify icon which will appear in the browser:
2238 /// serv->SetIcon("/ResetHPX", "rootsys/icons/ed_execute.png");
2239 ///
2240 /// One also can set extra property '_fastcmd', that command appear as
2241 /// tool button on the top of the browser tree:
2242 /// serv->SetItemField("/ResetHPX", "_fastcmd", "true");
2243 /// Or it is equivalent to specifying extra argument when register command:
2244 /// serv->RegisterCommand("/ResetHPX", "/hpx/->Reset()", "button;rootsys/icons/ed_delete.png");
2245 
2246 Bool_t TRootSniffer::RegisterCommand(const char *cmdname, const char *method, const char *icon)
2247 {
2248  CreateItem(cmdname, Form("command %s", method));
2249  SetItemField(cmdname, "_kind", "Command");
2250  if (icon != 0) {
2251  if (strncmp(icon, "button;", 7) == 0) {
2252  SetItemField(cmdname, "_fastcmd", "true");
2253  icon += 7;
2254  }
2255  if (*icon != 0) SetItemField(cmdname, "_icon", icon);
2256  }
2257  SetItemField(cmdname, "method", method);
2258  Int_t numargs = 0;
2259  do {
2260  TString nextname = TString::Format("%sarg%d%s","%",numargs+1,"%");
2261  if (strstr(method,nextname.Data())==0) break;
2262  numargs++;
2263  } while (numargs<100);
2264  if (numargs>0)
2265  SetItemField(cmdname, "_numargs", TString::Format("%d", numargs));
2266 
2267  return kTRUE;
2268 }
Bool_t ProduceImage(Int_t kind, const char *path, const char *options, void *&ptr, Long_t &length)
Method to produce image from specified object.
virtual void Append(const TImage *, const char *="+", const char *="#00000000")
Definition: TImage.h:191
TString fItemName
current level of hierarchy
Definition: TRootSniffer.h:43
Bool_t RegisterObject(const char *subfolder, TObject *obj)
Register object in subfolder structure subfolder parameter can have many levels like: ...
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
const char * item_prop_user
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:864
Bool_t SetResult(void *obj, TClass *cl, TDataMember *member=0)
Obsolete, use SetFoundResult instead.
A TFolder object is a collection of objects and folders.
Definition: TFolder.h:32
virtual void WriteString(const char *s)
Write string to I/O buffer.
void * GetResPtr() const
An array of TObjects.
Definition: TObjArray.h:39
Bool_t HasOption(const char *key) const
Returns true if the given key appears in the URL options list.
Definition: TUrl.cxx:670
Bool_t CanExploreItem(const char *path)
Method returns true when object has childs or one could try to expand item.
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:404
EImageFileTypes
Definition: TImage.h:52
ULong_t GetItemHash(const char *itemname)
Get hash function for specified item used to detect any changes in the specified object.
The concrete implementation of TBuffer for writing/reading to/from a ROOT file or socket...
Definition: TBufferFile.h:51
Bool_t Produce(const char *path, const char *file, const char *options, void *&ptr, Long_t &length, TString &str)
Method produce different kind of data out of object Parameter 'path' specifies object or object membe...
Bool_t IsStreamerInfoItem(const char *itemname)
Return true if it is streamer info item name.
UInt_t fMask
pointer on parent record
Definition: TRootSniffer.h:40
Bool_t ProduceBinary(const char *path, const char *options, void *&ptr, Long_t &length)
produce binary data for specified item if "zipped" option specified in query, buffer will be compress...
Int_t Depth() const
Returns depth of hierarchy.
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:487
virtual TList * GetListOfKeys() const
Definition: TDirectory.h:155
void BeforeNextChild()
indicates that new child for current element will be started
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
ClassImp(TSeqCollection) Int_t TSeqCollection TIter next(this)
Return index of object in collection.
Ssiz_t Length() const
Definition: TString.h:390
Collectable string class.
Definition: TObjString.h:32
ClassImp(TRootSniffer) TRootSniffer
constructor
TMemFile * fMemFile
default path for registered objects
Definition: TRootSniffer.h:122
static const EReturnType kOther
Definition: TMethodCall.h:50
TCanvas * c1
Definition: legend1.C:2
Int_t CheckRestriction(const char *item_name)
Checked if restriction is applied to the item full_item_name should have full path to the item...
All ROOT classes may have RTTI (run time type identification) support added.
Definition: TDataMember.h:33
Storage of hierarchy scan in TRootSniffer in JSON format.
This class represents a WWW compatible URL.
Definition: TUrl.h:41
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
static byte * ptr1
Definition: gifdecode.c:16
#define gDirectory
Definition: TDirectory.h:218
TDataMember * GetResMember() const
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
virtual TList * GetList() const
Definition: TDirectory.h:154
TH1 * h
Definition: legend2.C:5
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:892
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
virtual void * FindInHierarchy(const char *path, TClass **cl=0, TDataMember **member=0, Int_t *chld=0)
Search element with specified path Returns pointer on element Optionally one could obtain element cla...
Bool_t GoInside(TRootSnifferScanRec &super, TObject *obj, const char *obj_name=0, TRootSniffer *sniffer=0)
Method verifies if new level of hierarchy should be started with provided object. ...
virtual void Remove(TObject *obj)
Remove object from this folder. obj must be a TObject or a TFolder.
Definition: TFolder.cxx:464
Short_t GetCycle() const
Return cycle number associated to this key.
Definition: TKey.cxx:561
virtual void AddFirst(TObject *obj)
Add object at the beginning of the list.
Definition: TList.cxx:92
Long_t Property() const
Set TObject::fBits and fStreamerType to cache information about the class.
Definition: TClass.cxx:5652
TString DecodeUrlOptionValue(const char *value, Bool_t remove_quotes=kTRUE)
method replaces all kind of special symbols, which could appear in URL options
THist< 1, float > TH1F
Definition: THist.h:315
Bool_t CanSetFields() const
return true when fields could be set to the hierarchy item
Definition: TRootSniffer.h:61
#define gROOT
Definition: TROOT.h:340
Bool_t HasRestriction(const char *item_name)
Made fast check if item with specified name is in restriction list If returns true, requires precise check with CheckRestriction() method.
static const EReturnType kLong
Definition: TMethodCall.h:47
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
Bool_t SetFoundResult(void *obj, TClass *cl, TDataMember *member=0)
Set found element with class and datamember (optional)
virtual const char * GetClassName() const
Definition: TKey.h:77
Bool_t Done() const
Method indicates that scanning can be interrupted while result is set.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Bool_t IsaPointer() const
Return true if data member is a pointer.
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:254
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:496
void SetParent(TObject *parent)
Set parent owning this buffer.
Definition: TBuffer.cxx:237
void CreateMemFile()
Creates TMemFile instance, which used for objects streaming One could not use TBufferFile directly...
virtual void BeforeNextChild(Int_t, Int_t, Int_t)
void Restrict(const char *path, const char *options)
Restrict access to the specified location.
Each ROOT method (see TMethod) has a linked list of its arguments.
Definition: TMethodArg.h:33
An abstract interface to image processing library.
Definition: TImage.h:45
std::string GetReturnTypeNormalizedName() const
Get the normalized name of the return type.
Definition: TFunction.cxx:154
void MapObject(const TObject *obj, UInt_t offset=1)
Add object to the fMap container.
Bool_t UnregisterObject(TObject *obj)
unregister (remove) object from folders structures folder itself will remain even when it will be emp...
#define ROOT_VERSION_CODE
Definition: RVersion.h:21
Int_t fRestriction
list of created items names, need to avoid duplication
Definition: TRootSniffer.h:45
mask for actions, only actions copied to child rec
Definition: TRootSniffer.h:35
Bool_t ProduceItem(const char *path, const char *options, TString &res, Bool_t asjson=kTRUE)
produce JSON/XML for specified item contrary to h.json request, only fields for specified item are st...
Int_t GetMaxIndex(Int_t dim) const
Return maximum index for array dimension "dim".
Abstract interface for storage of hierarchy scan in TRootSniffer.
TMethod * GetMethodAllAny(const char *method)
Return pointer to method without looking at parameters.
Definition: TClass.cxx:4118
TString fCurrentAllowedMethods
current restriction for last-found object
Definition: TRootSniffer.h:128
Bool_t CreateItem(const char *fullname, const char *title)
create item element
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
const char * GetArrayIndex() const
If the data member is pointer and has a valid array size in its comments GetArrayIndex returns a stri...
void CloseNode()
close started node
Bool_t ProduceXml(const char *path, const char *options, TString &res)
produce XML data for specified item For object conversion TBufferXML is used
A TMemFile is like a normal TFile except that it reads and writes only from memory.
Definition: TMemFile.h:19
Bool_t CanExpandItem()
Returns true when item can be expanded.
TList * GetListOfRealData() const
Definition: TClass.h:405
static const EReturnType kString
Definition: TMethodCall.h:49
Bool_t IsScanGlobalDir() const
Definition: TRootSniffer.h:201
const char * Data() const
Definition: TString.h:349
const char * item_prop_rootversion
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
Definition: TObjArray.cxx:395
Int_t GetBaseClassOffset(const TClass *toBase, void *address=0, bool isDerivedObject=true)
Definition: TClass.cxx:2705
const char * GetAutoLoad() const
return name of configured autoload scripts (or 0)
Bool_t IsReadOnly(Bool_t dflt=kTRUE)
Returns read-only flag for current item.
const char * GetTrueTypeName() const
Get full type description of data member, e,g.: "class TDirectory*".
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2334
virtual void ScanObjectChilds(TRootSnifferScanRec &rec, TObject *obj)
scans object childs (if any) here one scans collection, branches, trees and so on ...
Vc_ALWAYS_INLINE void free(T *p)
Frees memory that was allocated with Vc::malloc.
Definition: memory.h:94
Bool_t SetItemField(const char *fullname, const char *name, const char *value)
set field for specified item
void Class()
Definition: Class.C:29
THttpCallArg * fCurrentArg
when enabled (default), scan gROOT for histograms, canvases, open files
Definition: TRootSniffer.h:126
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
Bool_t ExecuteCmd(const char *path, const char *options, TString &res)
execute command marked as _kind=='Command'
Bool_t fReadOnly
last produced streamer info
Definition: TRootSniffer.h:124
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition: TClass.cxx:4683
const char * item_prop_kind
Bool_t ScanOnlyFields() const
Definition: TRootSniffer.h:66
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
Bool_t CanDrawItem(const char *path)
Method verifies if object can be drawn.
const char * item_prop_typename
TFolder * GetSubFolder(const char *foldername, Bool_t force=kFALSE)
creates subfolder where objects can be registered
Int_t GetResRestrict() const
TString & Append(const char *cs)
Definition: TString.h:492
std::vector< std::vector< double > > Data
virtual void Add(TObject *obj)
Add object to this folder. obj must be a TObject or a TFolder.
Definition: TFolder.cxx:168
Int_t GetIntValueFromOptions(const char *key) const
Return a value for a given key from the URL options as an Int_t, a missing key returns -1...
Definition: TUrl.cxx:659
const char * fSearchPath
defines operation kind
Definition: TRootSniffer.h:41
const char * GetDefault() const
Get default value of method argument.
Definition: TMethodArg.cxx:58
static const EReturnType kDouble
Definition: TMethodCall.h:48
Int_t Atoi() const
Return integer value of string.
Definition: TString.cxx:1964
virtual void RecursiveRemove(TObject *obj)
Recursively remove object from a folder.
Definition: TFolder.cxx:456
char * Buffer() const
Definition: TBuffer.h:91
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:30
void * GetPostData() const
Definition: THttpCallArg.h:159
Method or function calling interface.
Definition: TMethodCall.h:41
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
void SetRootClass(TClass *cl)
Mark item with ROOT class and correspondent streamer info.
TFolder * AddFolder(const char *name, const char *title, TCollection *collection=0)
Create a new folder and add it to the list of folders of this folder, return a pointer to the created...
Definition: TFolder.cxx:184
Int_t fNumChilds
number of fields
Definition: TRootSniffer.h:51
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:59
virtual ~TRootSniffer()
destructor
Storage of hierarchy scan in TRootSniffer in XML format.
A doubly linked list.
Definition: TList.h:47
Int_t WithCurrentUserName(const char *option)
return 2 when option match to current user name return 1 when option==all return 0 when option does n...
search for specified item (only objects and collections)
Definition: TRootSniffer.h:32
TRealData * GetRealData(const char *name) const
Return pointer to TRealData element with name "name".
Definition: TClass.cxx:3275
void BuildRealData(void *pointer=0, Bool_t isTransient=kFALSE)
Build a full list of persistent data members.
Definition: TClass.cxx:1937
expand of specified item - allowed to scan object members
Definition: TRootSniffer.h:31
const char * item_prop_arraydim
const char * GetTypeName() const
Get type of data member, e,g.: "class TDirectory*" -> "TDirectory".
const char * GetValueFromOptions(const char *key) const
Return a value for a given key from the URL options.
Definition: TUrl.cxx:647
Double_t length(const TVector2 &v)
Definition: CsgOps.cxx:347
const char * item_prop_hidden
void SetResult(void *_res, TClass *_rescl, TDataMember *_resmemb, Int_t _res_chld, Int_t restr=0)
set pointer on found element, class and number of childs
void SetExtraHeader(const char *name, const char *value)
Definition: THttpCallArg.h:291
SVector< double, 2 > v
Definition: Dict.h:5
virtual void WriteStreamerInfo()
Write the list of TStreamerInfo as a single object in this file The class Streamer description for al...
Definition: TFile.cxx:3561
Int_t WriteObject(const T *obj, const char *name, Option_t *option="", Int_t bufsize=0)
Definition: TDirectory.h:201
normal scan of hierarchy
Definition: TRootSniffer.h:30
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:187
Int_t fLevel
current path searched
Definition: TRootSniffer.h:42
const char * item_prop_realname
Collection abstract base class.
Definition: TCollection.h:48
void Destructor(void *obj, Bool_t dtorOnly=kFALSE)
Explicitly call destructor for object.
Definition: TClass.cxx:5040
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2321
const char * item_prop_title
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
The most important graphics class in the ROOT system.
Definition: TPad.h:46
TString fObjectsPath
Definition: TRootSniffer.h:121
Bool_t RegisterCommand(const char *cmdname, const char *method, const char *icon)
Register command which can be executed from web interface.
char * Form(const char *fmt,...)
Int_t fNumFields
indicate if node was started
Definition: TRootSniffer.h:50
TLine * l
Definition: textangle.C:4
The TRealData class manages the effective list of all data members for a given class.
Definition: TRealData.h:34
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:81
virtual void ScanRoot(TRootSnifferScanRec &rec)
scan complete ROOT objects hierarchy For the moment it includes objects in gROOT directory and list o...
const char * GetUserName() const
Definition: THttpCallArg.h:187
static TString ConvertToXML(const TObject *obj, Bool_t GenericLayout=kFALSE, Bool_t UseNamespaces=kFALSE)
Converts object, inherited from TObject class, to XML string fmt contains configuration of XML layout...
Definition: TBufferXML.cxx:166
virtual TList * GetStreamerInfoList()
Read the list of TStreamerInfo objects written to this file.
Definition: TFile.cxx:1319
virtual ~TRootSnifferScanRec()
destructor
Bool_t AccessField(TFolder *parent, TObject *item, const char *name, const char *value, TNamed **only_get=0)
set or get field for the child each field coded as TNamed object, placed after chld in the parent hie...
void ScanHierarchy(const char *topname, const char *path, TRootSnifferStore *store, Bool_t only_fields=kFALSE)
Method scans normal objects, registered in ROOT.
TString fAutoLoad
list of restrictions for different locations
Definition: TRootSniffer.h:130
void ScanObjectMembers(TRootSnifferScanRec &rec, TClass *cl, char *ptr)
scripts names, which are add as _autoload parameter to h.json request
TGraphErrors * gr
Definition: legend1.C:25
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition: TString.cxx:2240
static TImage * Create()
TRootSnifferScanRec()
number of childs
Bool_t IsReadyForResult() const
Checks if result will be accepted.
virtual Int_t GetLast() const
Returns index of last object in collection.
TString & Remove(Ssiz_t pos)
Definition: TString.h:616
long Long_t
Definition: RtypesCore.h:50
int Ssiz_t
Definition: RtypesCore.h:63
The Canvas class.
Definition: TCanvas.h:48
TList fItemsNames
name of current item
Definition: TRootSniffer.h:44
static TObject * ConvertFromXML(const char *str, Bool_t GenericLayout=kFALSE, Bool_t UseNamespaces=kFALSE)
Read object from XML, produced by ConvertToXML() method.
Definition: TBufferXML.cxx:201
virtual Int_t GetSize() const
Definition: TCollection.h:95
static const char * GetFloatFormat()
return current printf format for float/double members, default "%e"
virtual const char * GetPrototype() const
Returns the prototype of a function as defined by CINT, or 0 in case of error.
Definition: TFunction.cxx:245
double f(double x)
Bool_t ProduceExe(const char *path, const char *options, Int_t reskind, TString *ret_str, void **ret_ptr=0, Long_t *ret_length=0)
execute command for specified object options include method and extra list of parameters sniffer shou...
TCollection * GetListOfFolders() const
Definition: TFolder.h:57
virtual void WriteLong(Long_t l)
Definition: TBufferFile.h:378
TList * GetListOfMethodArgs()
Return list containing the TMethodArgs of a TFunction.
Definition: TFunction.cxx:126
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:415
double Double_t
Definition: RtypesCore.h:55
Bool_t ProduceMulti(const char *path, const char *options, void *&ptr, Long_t &length, TString &str, Bool_t asjson=kTRUE)
Process several requests, packing all results into binary or JSON buffer Input parameters should be c...
TRootSnifferStore * fStore
restriction 0 - default, 1 - read-only, 2 - full access
Definition: TRootSniffer.h:47
Describe directory structure in memory.
Definition: TDirectory.h:41
virtual void FromPad(TVirtualPad *, Int_t=0, Int_t=0, UInt_t=0, UInt_t=0)
Definition: TImage.h:138
virtual void GetImageBuffer(char **, int *, EImageFileTypes=TImage::kPng)
Definition: TImage.h:257
TNamed()
Definition: TNamed.h:40
unsigned long ULong_t
Definition: RtypesCore.h:51
double func(double *x, double *p)
Definition: stressTF1.cxx:213
virtual void WriteDouble(Double_t d)
Definition: TBufferFile.h:413
TList * fSinfo
file used to manage streamer infos
Definition: TRootSniffer.h:123
static TString ConvertToJSON(const TObject *obj, Int_t compact=0, const char *member_name=0)
converts object, inherited from TObject class, to JSON string
TRootSnifferScanRec * fParent
Definition: TRootSniffer.h:39
virtual void CreateNode(Int_t, const char *)
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition: TClass.cxx:2881
TCanvas * slash()
Definition: slash.C:1
virtual void AddAfter(const TObject *after, TObject *obj)
Insert object after object after in the list.
Definition: TList.cxx:220
Bool_t IsItemField(TObject *obj) const
return true when object is TNamed with kItemField bit set such objects used to keep field values for ...
#define name(a, b)
Definition: linkTestLib0.cpp:5
Bool_t ProduceJson(const char *path, const char *options, TString &res)
produce JSON data for specified item For object conversion TBufferJSON is used
Mother of all ROOT objects.
Definition: TObject.h:58
Global functions class (global functions are obtained from CINT).
Definition: TFunction.h:30
Int_t GetArrayDim() const
Return number of array dimensions.
virtual void SetField(Int_t, const char *, const char *, Bool_t)
Int_t IsSTLContainer()
The return type is defined in TDictionary (kVector, kList, etc.)
void SetField(const char *name, const char *value, Bool_t with_quotes=kTRUE)
Set item field only when creating is specified.
void SetAutoLoad(const char *scripts="")
When specified, _autoload attribute will be always add to top element of h.json/h.hml requests Used to instruct browser automatically load special code.
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition: TKey.cxx:727
const char * GetFullTypeName() const
Get full type description of method argument, e.g.: "class TDirectory*".
Definition: TMethodArg.cxx:75
virtual void Add(TObject *obj)
Definition: TList.h:81
const Ssiz_t kNPOS
Definition: Rtypes.h:115
void ScanCollection(TRootSnifferScanRec &rec, TCollection *lst, const char *foldername=0, TCollection *keys_lst=0)
scan collection content
Int_t Length() const
Definition: TBuffer.h:94
void Execute(const char *, const char *, int *=0)
Execute method on this object with the given parameter string, e.g.
Definition: TMethodCall.h:68
Each ROOT class (see TClass) has a linked list of methods.
Definition: TMethod.h:40
virtual void CloseNode(Int_t, Int_t)
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
void SetOptions(const char *opt)
Definition: TUrl.h:96
void ParseOptions() const
Parse URL options into a key/value map.
Definition: TUrl.cxx:616
Bool_t fHasMore
object to store results
Definition: TRootSniffer.h:48
const char * item_prop_more
#define gPad
Definition: TVirtualPad.h:288
R__EXTERN Int_t gDebug
Definition: Rtypes.h:128
bool debug
Bool_t IsBasic() const
Return true if data member is a basic type, e.g. char, int, long...
Long_t Property() const
Get property description word. For meaning of bits see EProperty.
Bool_t fNodeStarted
indicates that potentially there are more items can be found
Definition: TRootSniffer.h:49
virtual TObject * FindObject(const char *name) const
Search object identified by name in the tree of folders inside this folder.
Definition: TFolder.cxx:308
A TTree object has a header with a name and a title.
Definition: TTree.h:94
const char * item_prop_autoload
Long_t GetPostDataLength() const
Definition: THttpCallArg.h:166
TDataMember * GetDataMember() const
Definition: TRealData.h:57
Bool_t InheritsFrom(const char *cl) const
Return kTRUE if this class inherits from a class with name "classname".
Definition: TClass.cxx:4579
static const EReturnType kNone
Definition: TMethodCall.h:53
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
virtual void ScanObjectProperties(TRootSnifferScanRec &rec, TObject *obj)
scans object properties here such fields as _autoload or _icon properties depending on class or objec...
A TTree is a list of TBranches.
Definition: TBranch.h:58
void MakeItemName(const char *objname, TString &itemname)
Construct item name, using object name as basis.
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t GetResNumChilds() const
Long_t GetThisOffset() const
Definition: TRealData.h:59
virtual const char * GetName() const
Returns name of object.
Definition: TRealData.h:56
void SetCurrentCallArg(THttpCallArg *arg)
set current http arguments, which then used in different process methods For instance, if user authorized with some user name, depending from restrictions some objects will be invisible or user get full access to the element
TMethod * GetMethodWithPrototype(const char *method, const char *proto, Bool_t objectIsConst=kFALSE, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Find the method with a given prototype.
Definition: TClass.cxx:4194
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:459
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
TObject * obj
float value
Definition: math.cpp:443
Vc_ALWAYS_INLINE_L T *Vc_ALWAYS_INLINE_R malloc(size_t n)
Allocates memory on the Heap with alignment and padding suitable for vectorized access.
Definition: memory.h:67
void BuildFullName(TString &buf, TRootSnifferScanRec *prnt=0)
Produces full name for the current item.
Int_t fCurrentRestrict
current http arguments (if any)
Definition: TRootSniffer.h:127
const Int_t n
Definition: legend1.C:16
if set, only fields for specified item will be set (but all fields)
Definition: TRootSniffer.h:34
const char * GetItemField(TFolder *parent, TObject *item, const char *name)
return field for specified item
Bool_t IsValid() const
Return true if the method call has been properly initialized and is usable.
Ssiz_t First(char c) const
Find first occurrence of a character c.
Definition: TString.cxx:466
const char * cnt
Definition: TXMLSetup.cxx:75
EReturnType ReturnType()
Returns the return type of the method.
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Compare a string to char *cs2.
Definition: TString.cxx:385
void CreateNode(const char *_node_name)
Starts new node, must be closed at the end.
TObject * GetItem(const char *fullname, TFolder *&parent, Bool_t force=kFALSE, Bool_t within_objects=kTRUE)
return item from the subfolders structure
void Resize(Ssiz_t n)
Resize the string. Truncate or add blanks as necessary.
Definition: TString.cxx:1058
check if there childs, very similar to search
Definition: TRootSniffer.h:33
TList fRestrictions
list of allowed methods, extracted when analyzed object restrictions
Definition: TRootSniffer.h:129
ULong_t GetStreamerInfoHash()
Returns hash value for streamer infos At the moment - just number of items in streamer infos list...
TObject * FindTObjectInHierarchy(const char *path)
Search element in hierarchy, derived from TObject.
TClass * GetResClass() const
#define gFile
Definition: TFile.h:307
static Bool_t IsDrawableClass(TClass *cl)
return true if object can be drawn
UInt_t Hash(ECaseCompare cmp=kExact) const
Return hash value.
Definition: TString.cxx:605