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