Logo ROOT   6.10/09
Reference Guide
TBufferXML.cxx
Go to the documentation of this file.
1 // @(#)root/:$Id: 5400e36954e1dc109fcfc306242c30234beb7312 $
2 // Author: Sergey Linev, Rene Brun 10.05.2004
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2004, 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 /**
13 \class TBufferXML
14 \ingroup IO
15 
16 Class for serializing/deserializing object to/from xml.
17 
18 It redefines most of TBuffer class function to convert simple types,
19 array of simple types and objects to/from xml.
20 Instead of writing a binary data it creates a set of xml structures as
21 nodes and attributes
22 TBufferXML class uses streaming mechanism, provided by ROOT system,
23 therefore most of ROOT and user classes can be stored to xml. There are
24 limitations for complex objects like TTree, which can not be yet converted to xml.
25 */
26 
27 
28 #include "TBufferXML.h"
29 #include "Compression.h"
30 #include "TXMLFile.h"
31 
32 #include "TObjArray.h"
33 #include "TROOT.h"
34 #include "TClass.h"
35 #include "TClassTable.h"
36 #include "TDataType.h"
37 #include "TExMap.h"
38 #include "TMethodCall.h"
39 #include "TStreamerInfo.h"
40 #include "TStreamerElement.h"
41 #include "TProcessID.h"
42 #include "TFile.h"
43 #include "TMemberStreamer.h"
44 #include "TStreamer.h"
45 #include "TStreamerInfoActions.h"
46 #include "RZip.h"
47 
48 #ifdef R__VISUAL_CPLUSPLUS
49 #define FLong64 "%I64d"
50 #define FULong64 "%I64u"
51 #else
52 #define FLong64 "%lld"
53 #define FULong64 "%llu"
54 #endif
55 
57 
58 
59 std::string TBufferXML::fgFloatFmt = "%e";
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 /// Default constructor
63 
65  TBufferFile(),
66  TXMLSetup(),
67  fXML(0),
68  fStack(),
69  fVersionBuf(-111),
70  fObjMap(0),
71  fIdArray(0),
72  fErrorFlag(0),
73  fCanUseCompact(kFALSE),
74  fExpectedChain(kFALSE),
75  fExpectedBaseClass(0),
76  fCompressLevel(0),
77  fIOVersion(3)
78 {
79 }
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 /// Creates buffer object to serailize/deserialize data to/from xml.
83 /// Mode should be either TBuffer::kRead or TBuffer::kWrite.
84 
86  TBufferFile(mode),
87  TXMLSetup(),
88  fXML(0),
89  fStack(),
90  fVersionBuf(-111),
91  fObjMap(0),
92  fIdArray(0),
93  fErrorFlag(0),
97  fCompressLevel(0),
98  fIOVersion(3)
99 {
100  fBufSize = 1000000000;
101 
102  SetParent(0);
105 }
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 /// Creates buffer object to serailize/deserialize data to/from xml.
109 /// This constructor should be used, if data from buffer supposed to be stored in file.
110 /// Mode should be either TBuffer::kRead or TBuffer::kWrite.
111 
113  TBufferFile(mode),
114  TXMLSetup(*file),
115  fXML(0),
116  fStack(),
117  fVersionBuf(-111),
118  fObjMap(0),
119  fIdArray(0),
120  fErrorFlag(0),
124  fCompressLevel(0),
125  fIOVersion(3)
126 {
127  // this is for the case when StreamerInfo reads elements from
128  // buffer as ReadFastArray. When it checks if size of buffer is
129  // too small and skip reading. Actually, more improved method should
130  // be used here.
131  fBufSize = 1000000000;
132 
133  SetParent(file);
136  if (XmlFile()) {
137  SetXML(XmlFile()->XML());
140  }
141 }
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 /// Destroy xml buffer.
145 
147 {
148  if (fObjMap) delete fObjMap;
149  if (fIdArray) delete fIdArray;
150  fStack.Delete();
151 }
152 
153 ////////////////////////////////////////////////////////////////////////////////
154 /// Returns pointer to TXMLFile object.
155 /// Access to file is necessary to produce unique identifier for object references.
156 
158 {
159  return dynamic_cast<TXMLFile*>(GetParent());
160 }
161 
162 ////////////////////////////////////////////////////////////////////////////////
163 /// Converts object, inherited from TObject class, to XML string
164 /// GenericLayout defines layout choice for XML file
165 /// UseNamespaces allow XML namespaces.
166 /// See TXMLSetup class for details
167 
168 TString TBufferXML::ConvertToXML(const TObject* obj, Bool_t GenericLayout, Bool_t UseNamespaces)
169 {
170  TClass *clActual = 0;
171  void *ptr = (void *) obj;
172 
173  if (obj!=0) {
174  clActual = TObject::Class()->GetActualClass(obj);
175  if (!clActual) clActual = TObject::Class(); else
176  if (clActual != TObject::Class())
177  ptr = (void *) ((Long_t) obj - clActual->GetBaseClassOffset(TObject::Class()));
178  }
179 
180  return ConvertToXML(ptr, clActual, GenericLayout, UseNamespaces);
181 }
182 
183 ////////////////////////////////////////////////////////////////////////////////
184 /// Converts any type of object to XML string.
185 /// GenericLayout defines layout choice for XML file
186 /// UseNamespaces allow XML namespaces.
187 /// See TXMLSetup class for details
188 
189 TString TBufferXML::ConvertToXML(const void* obj, const TClass* cl, Bool_t GenericLayout, Bool_t UseNamespaces)
190 {
191  TXMLEngine xml;
192 
194  buf.SetXML(&xml);
195 
197  buf.SetUseNamespaces(UseNamespaces);
198 
199  XMLNodePointer_t xmlnode = buf.XmlWriteAny(obj, cl);
200 
201  TString res;
202 
203  xml.SaveSingleNode(xmlnode, &res);
204 
205  xml.FreeNode(xmlnode);
206 
207  return res;
208 }
209 
210 ////////////////////////////////////////////////////////////////////////////////
211 /// Read object from XML, produced by ConvertToXML() method.
212 /// If object does not inherit from TObject class, return 0.
213 /// GenericLayout and UseNamespaces should be the same as in ConvertToXML()
214 
215 TObject* TBufferXML::ConvertFromXML(const char* str, Bool_t GenericLayout, Bool_t UseNamespaces)
216 {
217  TClass* cl = 0;
218  void* obj = ConvertFromXMLAny(str, &cl, GenericLayout, UseNamespaces);
219 
220  if ((cl==0) || (obj==0)) return 0;
221 
223 
224  if (delta<0) {
225  cl->Destructor(obj);
226  return 0;
227  }
228 
229  return (TObject*) ( ( (char*)obj ) + delta );
230 }
231 
232 ////////////////////////////////////////////////////////////////////////////////
233 /// Read object of any class from XML, produced by ConvertToXML() method.
234 /// If cl!=0, return actual class of object.
235 /// GenericLayout and UseNamespaces should be the same as in ConvertToXML()
236 
237 void* TBufferXML::ConvertFromXMLAny(const char* str, TClass** cl, Bool_t GenericLayout, Bool_t UseNamespaces)
238 {
239  TXMLEngine xml;
241 
242  buf.SetXML(&xml);
243 
245  buf.SetUseNamespaces(UseNamespaces);
246 
247  XMLNodePointer_t xmlnode = xml.ReadSingleNode(str);
248 
249  void* obj = buf.XmlReadAny(xmlnode, 0, cl);
250 
251  xml.FreeNode(xmlnode);
252 
253  return obj;
254 }
255 
256 ////////////////////////////////////////////////////////////////////////////////
257 /// Convert object of any class to xml structures
258 /// Return pointer on top xml element
259 
261 {
262  fErrorFlag = 0;
263 
264  if (fXML==0) return 0;
265 
266  XMLNodePointer_t res = XmlWriteObject(obj, cl);
267 
268  return res;
269 }
270 
271 ////////////////////////////////////////////////////////////////////////////////
272 /// Recreate object from xml structure.
273 /// Return pointer to read object.
274 /// if (cl!=0) returns pointer to class of object
275 
276 void* TBufferXML::XmlReadAny(XMLNodePointer_t node, void* obj, TClass** cl)
277 {
278  if (node==0) return 0;
279  if (cl) *cl = 0;
280 
281  fErrorFlag = 0;
282 
283  if (fXML==0) return 0;
284 
285  PushStack(node, kTRUE);
286 
287  void* res = XmlReadObject(obj, cl);
288 
289  PopStack();
290 
291  return res;
292 }
293 
294 ////////////////////////////////////////////////////////////////////////////////
295 /// Convert object into xml structures.
296 /// !!! Should be used only by TBufferXML itself.
297 /// Use ConvertToXML() methods to convert your object to xml
298 /// Redefined here to avoid gcc 3.x warning
299 
301 {
303 }
304 
305 // TXMLStackObj is used to keep stack of object hierarchy,
306 // stored in TBuffer. For instnace, data for parent class(es)
307 // stored in subnodes, but initial object node will be kept.
308 
309 class TXMLStackObj : public TObject {
310  public:
311  TXMLStackObj(XMLNodePointer_t node) :
312  TObject(),
313  fNode(node),
314  fInfo(0),
315  fElem(0),
316  fElemNumber(0),
317  fCompressedClassNode(kFALSE),
318  fClassNs(0),
319  fIsStreamerInfo(kFALSE),
320  fIsElemOwner(kFALSE)
321  {}
322 
323  virtual ~TXMLStackObj()
324  {
325  if (fIsElemOwner) delete fElem;
326  }
327 
328  Bool_t IsStreamerInfo() const { return fIsStreamerInfo; }
329 
330  XMLNodePointer_t fNode;
332  TStreamerElement* fElem;
333  Int_t fElemNumber;
334  Bool_t fCompressedClassNode;
335  XMLNsPointer_t fClassNs;
336  Bool_t fIsStreamerInfo;
337  Bool_t fIsElemOwner;
338 };
339 
340 ////////////////////////////////////////////////////////////////////////////////
341 /// Add new level to xml stack.
342 
343 TXMLStackObj* TBufferXML::PushStack(XMLNodePointer_t current, Bool_t simple)
344 {
345  if (IsReading() && !simple) {
346  current = fXML->GetChild(current);
347  fXML->SkipEmpty(current);
348  }
349 
350  TXMLStackObj* stack = new TXMLStackObj(current);
351  fStack.Add(stack);
352  return stack;
353 }
354 
355 ////////////////////////////////////////////////////////////////////////////////
356 /// Remove one level from xml stack.
357 
358 TXMLStackObj* TBufferXML::PopStack()
359 {
360  TObject* last = fStack.Last();
361  if (last!=0) {
362  fStack.Remove(last);
363  delete last;
364  fStack.Compress();
365  }
366  return dynamic_cast<TXMLStackObj*> (fStack.Last());
367 }
368 
369 ////////////////////////////////////////////////////////////////////////////////
370 /// Return xml stack object of specified depth.
371 
372 TXMLStackObj* TBufferXML::Stack(Int_t depth)
373 {
374  TXMLStackObj* stack = 0;
375  if (depth<=fStack.GetLast())
376  stack = dynamic_cast<TXMLStackObj*> (fStack.At(fStack.GetLast()-depth));
377  return stack;
378 }
379 
380 ////////////////////////////////////////////////////////////////////////////////
381 /// Return pointer on current xml node.
382 
384 {
385  TXMLStackObj* stack = dynamic_cast<TXMLStackObj*> (fStack.Last());
386  return (stack==0) ? 0 : stack->fNode;
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// Shift stack node to next.
391 
392 void TBufferXML::ShiftStack(const char* errinfo)
393 {
394  TXMLStackObj* stack = dynamic_cast<TXMLStackObj*> (fStack.Last());
395  if (stack) {
396  fXML->ShiftToNext(stack->fNode);
397  if (gDebug>4) Info("ShiftStack","%s to node %s", errinfo, fXML->GetNodeName(stack->fNode));
398  }
399 }
400 
401 ////////////////////////////////////////////////////////////////////////////////
402 /// See comments for function SetCompressionSettings.
403 
405 {
406  if (algorithm < 0 || algorithm >= ROOT::kUndefinedCompressionAlgorithm) algorithm = 0;
407  if (fCompressLevel < 0) {
408  // if the level is not defined yet use 1 as a default
409  fCompressLevel = 100 * algorithm + 1;
410  } else {
411  int level = fCompressLevel % 100;
412  fCompressLevel = 100 * algorithm + level;
413  }
414 }
415 
416 ////////////////////////////////////////////////////////////////////////////////
417 /// See comments for function SetCompressionSettings.
418 
420 {
421  if (level < 0) level = 0;
422  if (level > 99) level = 99;
423  if (fCompressLevel < 0) {
424  // if the algorithm is not defined yet use 0 as a default
425  fCompressLevel = level;
426  } else {
427  int algorithm = fCompressLevel / 100;
428  if (algorithm >= ROOT::kUndefinedCompressionAlgorithm) algorithm = 0;
429  fCompressLevel = 100 * algorithm + level;
430  }
431 }
432 
433 ////////////////////////////////////////////////////////////////////////////////
434 /// Used to specify the compression level and algorithm.
435 ///
436 /// See TFile constructor for the details.
437 
439 {
440  fCompressLevel = settings;
441 }
442 
443 ////////////////////////////////////////////////////////////////////////////////
444 /// Write binary data block from buffer to xml.
445 /// This data can be produced only by direct call of TBuffer::WriteBuf() functions.
446 
448 {
449  if ((node==0) || (Length()==0)) return;
450 
451  const char* src = Buffer();
452  int srcSize = Length();
453 
454  char* fZipBuffer = 0;
455 
456  Int_t compressionLevel = GetCompressionLevel();
457  Int_t compressionAlgorithm = GetCompressionAlgorithm();
458 
459  if ((Length() > 512) && (compressionLevel > 0)) {
460  int zipBufferSize = Length();
461  fZipBuffer = new char[zipBufferSize + 9];
462  int dataSize = Length();
463  int compressedSize = 0;
464  R__zipMultipleAlgorithm(compressionLevel, &dataSize, Buffer(), &zipBufferSize,
465  fZipBuffer, &compressedSize, compressionAlgorithm);
466  if (compressedSize > 0) {
467  src = fZipBuffer;
468  srcSize = compressedSize;
469  } else {
470  delete[] fZipBuffer;
471  fZipBuffer = 0;
472  }
473  }
474 
475  TString res;
476  char sbuf[500];
477  int block = 0;
478  char* tgt = sbuf;
479  int srcCnt = 0;
480 
481  while (srcCnt++<srcSize) {
482  tgt+=sprintf(tgt, " %02x", (unsigned char) *src);
483  src++;
484  if (block++==100) {
485  res += sbuf;
486  block = 0;
487  tgt = sbuf;
488  }
489  }
490 
491  if (block>0) res += sbuf;
492 
493  XMLNodePointer_t blocknode = fXML->NewChild(node, 0, xmlio::XmlBlock, res);
494  fXML->NewIntAttr(blocknode, xmlio::Size, Length());
495 
496  if (fZipBuffer) {
497  fXML->NewIntAttr(blocknode, xmlio::Zip, srcSize);
498  delete[] fZipBuffer;
499  }
500 }
501 
502 ////////////////////////////////////////////////////////////////////////////////
503 /// Read binary block of data from xml.
504 
506 {
507  if (blocknode==0) return;
508 
509  Int_t blockSize = fXML->GetIntAttr(blocknode, xmlio::Size);
510  Bool_t blockCompressed = fXML->HasAttr(blocknode, xmlio::Zip);
511  char* fUnzipBuffer = 0;
512 
513  if (gDebug>2)
514  Info("XmlReadBlock","Block size = %d, Length = %d, Compressed = %d",
515  blockSize, Length(), blockCompressed);
516 
517  if (blockSize>BufferSize()) Expand(blockSize);
518 
519  char* tgt = Buffer();
520  Int_t readSize = blockSize;
521 
522  TString content = fXML->GetNodeContent(blocknode);
523 
524  if (blockCompressed) {
525  Int_t zipSize = fXML->GetIntAttr(blocknode, xmlio::Zip);
526  fUnzipBuffer = new char[zipSize];
527 
528  tgt = fUnzipBuffer;
529  readSize = zipSize;
530  }
531 
532  char* ptr = (char*) content.Data();
533 
534  if (gDebug>3)
535  Info("XmlReadBlock","Content %s", ptr);
536 
537  for (int i=0;i<readSize;i++) {
538  while ((*ptr<48) || ((*ptr>57) && (*ptr<97)) || (*ptr>102)) ptr++;
539 
540  int b_hi = (*ptr>57) ? *ptr-87 : *ptr-48;
541  ptr++;
542  int b_lo = (*ptr>57) ? *ptr-87 : *ptr-48;
543  ptr++;
544 
545  *tgt=b_hi*16+b_lo;
546  tgt++;
547 
548  if (gDebug>4) Info("XmlReadBlock"," Buf[%d] = %d", i, b_hi*16+b_lo);
549  }
550 
551  if (fUnzipBuffer) {
552 
553  int srcsize;
554  int tgtsize;
555  int status = R__unzip_header(&srcsize, (UChar_t*) fUnzipBuffer, &tgtsize);
556 
557  int unzipRes = 0;
558  if (status == 0) {
559  R__unzip(&readSize, (unsigned char*) fUnzipBuffer, &blockSize,
560  (unsigned char*) Buffer(), &unzipRes);
561  }
562  if (status != 0 || unzipRes!=blockSize)
563  Error("XmlReadBlock", "Decompression error %d", unzipRes);
564  else
565  if (gDebug>2) Info("XmlReadBlock","Unzip ok");
566  delete[] fUnzipBuffer;
567  }
568 }
569 
570 ////////////////////////////////////////////////////////////////////////////////
571 /// Add "ptr" attribute to node, if ptr is null or
572 /// if ptr is pointer on object, which is already saved in buffer
573 /// Automatically add "ref" attribute to node, where referenced object is stored
574 
576 {
577  if (node==0) return kFALSE;
578 
579  TString refvalue;
580 
581  if (ptr==0)
582  refvalue = xmlio::Null; //null
583  else {
584  if (fObjMap==0) return kFALSE;
585 
586  ULong_t hash = TString::Hash(&ptr, sizeof(void*));
587 
588  XMLNodePointer_t refnode = (XMLNodePointer_t) (Long_t)fObjMap->GetValue(hash, (Long_t) ptr);
589  if (refnode==0) return kFALSE;
590 
591  if (fXML->HasAttr(refnode, xmlio::Ref))
592  refvalue = fXML->GetAttr(refnode, xmlio::Ref);
593  else {
594  refvalue = xmlio::IdBase;
595  if (XmlFile())
596  refvalue += XmlFile()->GetNextRefCounter();
597  else
598  refvalue += GetNextRefCounter();
599  fXML->NewAttr(refnode, 0, xmlio::Ref, refvalue.Data());
600  }
601  }
602  if (refvalue.Length()>0) {
603  fXML->NewAttr(node, 0, xmlio::Ptr, refvalue.Data());
604  return kTRUE;
605  }
606 
607  return kFALSE;
608 }
609 
610 ////////////////////////////////////////////////////////////////////////////////
611 /// Register pair of object pointer and node, where this object is saved,
612 /// in object map
613 
615 {
616  if ((node==0) || (ptr==0)) return;
617 
618  ULong_t hash = TString::Hash(&ptr, sizeof(void*));
619 
620  if (fObjMap==0) fObjMap = new TExMap();
621 
622  if (fObjMap->GetValue(hash, (Long_t) ptr)==0)
623  fObjMap->Add(hash, (Long_t) ptr, (Long_t) node);
624 }
625 
626 ////////////////////////////////////////////////////////////////////////////////
627 /// Searches for "ptr" attribute and returns pointer to object and class,
628 /// if "ptr" attribute reference to read object
629 
631 {
632  cl = 0;
633 
634  if (!fXML->HasAttr(node,xmlio::Ptr)) return kFALSE;
635 
636  const char* ptrid = fXML->GetAttr(node, xmlio::Ptr);
637 
638  if (ptrid==0) return kFALSE;
639 
640  // null
641  if (strcmp(ptrid, xmlio::Null)==0) {
642  ptr = 0;
643  return kTRUE;
644  }
645 
646  if ((fIdArray==0) || (fObjMap==0)) return kFALSE;
647 
648  TNamed* obj = (TNamed*) fIdArray->FindObject(ptrid);
649  if (obj) {
650  ptr = (void*) (Long_t)fObjMap->GetValue((Long_t) fIdArray->IndexOf(obj));
651  cl = TClass::GetClass(obj->GetTitle());
652  return kTRUE;
653  }
654  return kFALSE;
655 }
656 
657 ////////////////////////////////////////////////////////////////////////////////
658 /// Analyse, if node has "ref" attribute and register it to object map
659 
660 void TBufferXML::ExtractReference(XMLNodePointer_t node, const void* ptr, const TClass* cl)
661 {
662  if ((node==0) || (ptr==0)) return;
663 
664  const char* refid = fXML->GetAttr(node, xmlio::Ref);
665 
666  if (refid==0) return;
667 
668  if (fIdArray==0) {
669  fIdArray = new TObjArray;
671  }
672  TNamed* nid = new TNamed(refid, cl->GetName());
673  fIdArray->Add(nid);
674 
675  if (fObjMap==0) fObjMap = new TExMap();
676 
677  fObjMap->Add((Long_t) fIdArray->IndexOf(nid), (Long_t) ptr);
678 
679  if (gDebug>2)
680  Info("ExtractReference","Find reference %s for object %p", refid, ptr);
681 }
682 
683 ////////////////////////////////////////////////////////////////////////////////
684 /// Check, if node has specified name
685 
686 Bool_t TBufferXML::VerifyNode(XMLNodePointer_t node, const char* name, const char* errinfo)
687 {
688  if ((name==0) || (node==0)) return kFALSE;
689 
690  if (strcmp(fXML->GetNodeName(node), name)!=0) {
691  if (errinfo) {
692  Error("VerifyNode","Reading XML file (%s). Get: %s, expects: %s",
693  errinfo, fXML->GetNodeName(node), name);
694  fErrorFlag = 1;
695  }
696  return kFALSE;
697  }
698  return kTRUE;
699 }
700 
701 ////////////////////////////////////////////////////////////////////////////////
702 /// Check, if stack node has specified name
703 
704 Bool_t TBufferXML::VerifyStackNode(const char* name, const char* errinfo)
705 {
706  return VerifyNode(StackNode(), name, errinfo);
707 }
708 
709 
710 ////////////////////////////////////////////////////////////////////////////////
711 /// Checks, that attribute of specified name exists and has specified value
712 
713 Bool_t TBufferXML::VerifyAttr(XMLNodePointer_t node, const char* name, const char* value, const char* errinfo)
714 {
715  if ((node==0) || (name==0) || (value==0)) return kFALSE;
716  const char* cont = fXML->GetAttr(node, name);
717  if (((cont==0) || (strcmp(cont, value)!=0))) {
718  if (errinfo) {
719  Error("VerifyAttr","%s : attr %s = %s, expected: %s", errinfo, name, cont, value);
720  fErrorFlag = 1;
721  }
722  return kFALSE;
723  }
724  return kTRUE;
725 }
726 
727 ////////////////////////////////////////////////////////////////////////////////
728 /// Checks stack attribute
729 
730 Bool_t TBufferXML::VerifyStackAttr(const char* name, const char* value, const char* errinfo)
731 {
732  return VerifyAttr(StackNode(), name, value, errinfo);
733 }
734 
735 ////////////////////////////////////////////////////////////////////////////////
736 /// Create item node of specified name
737 
739 {
740  XMLNodePointer_t node = 0;
741  if (GetXmlLayout()==kGeneralized) {
742  node = fXML->NewChild(StackNode(), 0, xmlio::Item, 0);
743  fXML->NewAttr(node, 0, xmlio::Name, name);
744  } else
745  node = fXML->NewChild(StackNode(), 0, name, 0);
746  return node;
747 }
748 
749 ////////////////////////////////////////////////////////////////////////////////
750 /// Checks, if stack node is item and has specified name
751 
752 Bool_t TBufferXML::VerifyItemNode(const char* name, const char* errinfo)
753 {
754  Bool_t res = kTRUE;
755  if (GetXmlLayout()==kGeneralized)
756  res = VerifyStackNode(xmlio::Item, errinfo) &&
757  VerifyStackAttr(xmlio::Name, name, errinfo);
758  else
759  res = VerifyStackNode(name, errinfo);
760  return res;
761 }
762 
763 ////////////////////////////////////////////////////////////////////////////////
764 /// Create xml node correspondent to TStreamerElement object
765 
767 {
768  XMLNodePointer_t elemnode = 0;
769 
770  const char* elemxmlname = XmlGetElementName(elem);
771 
772  if (GetXmlLayout()==kGeneralized) {
773  elemnode = fXML->NewChild(StackNode(), 0, xmlio::Member, 0);
774  fXML->NewAttr(elemnode, 0, xmlio::Name, elemxmlname);
775  } else {
776  // take namesapce for element only if it is not a base class or class name
777  XMLNsPointer_t ns = Stack()->fClassNs;
778  if ((elem->GetType()==TStreamerInfo::kBase)
779  || ((elem->GetType()==TStreamerInfo::kTNamed) && !strcmp(elem->GetName(), TNamed::Class()->GetName()))
780  || ((elem->GetType()==TStreamerInfo::kTObject) && !strcmp(elem->GetName(), TObject::Class()->GetName()))
781  || ((elem->GetType()==TStreamerInfo::kTString) && !strcmp(elem->GetName(), TString::Class()->GetName())))
782  ns = 0;
783 
784  elemnode = fXML->NewChild(StackNode(), ns, elemxmlname, 0);
785  }
786 
787  TXMLStackObj* curr = PushStack(elemnode);
788  curr->fElem = (TStreamerElement*)elem;
789 }
790 
791 ////////////////////////////////////////////////////////////////////////////////
792 /// Checks, if stack node correspond to TStreamerElement object
793 
795 {
796  const char* elemxmlname = XmlGetElementName(elem);
797 
798  if (GetXmlLayout()==kGeneralized) {
799  if (!VerifyStackNode(xmlio::Member)) return kFALSE;
800  if (!VerifyStackAttr(xmlio::Name, elemxmlname)) return kFALSE;
801  } else {
802  if (!VerifyStackNode(elemxmlname)) return kFALSE;
803  }
804 
806 
807  TXMLStackObj* curr = PushStack(StackNode()); // set pointer to first data inside element
808  curr->fElem = (TStreamerElement*)elem;
809  return kTRUE;
810 }
811 
812 ////////////////////////////////////////////////////////////////////////////////
813 /// Write object to buffer
814 /// If object was written before, only pointer will be stored
815 /// Return pointer to top xml node, representing object
816 
818 {
820 
821  if (!cl) obj = 0;
822  if (ProcessPointer(obj, objnode)) return objnode;
823 
824  TString clname = XmlConvertClassName(cl->GetName());
825 
826  fXML->NewAttr(objnode, 0, xmlio::ObjClass, clname);
827 
828  RegisterPointer(obj, objnode);
829 
830  PushStack(objnode);
831 
832  ((TClass*)cl)->Streamer((void*)obj, *this);
833 
834  PopStack();
835 
836  if (gDebug>1)
837  Info("XmlWriteObject","Done write for class: %s", cl ? cl->GetName() : "null");
838 
839  return objnode;
840 }
841 
842 ////////////////////////////////////////////////////////////////////////////////
843 /// Read object from the buffer
844 
845 void* TBufferXML::XmlReadObject(void* obj, TClass** cl)
846 {
847  if (cl) *cl = 0;
848 
849  XMLNodePointer_t objnode = StackNode();
850 
851  if (fErrorFlag>0) return obj;
852 
853  if (objnode==0) return obj;
854 
855  if (!VerifyNode(objnode, xmlio::Object, "XmlReadObjectNew")) return obj;
856 
857  TClass* objClass = 0;
858 
859  if (ExtractPointer(objnode, obj, objClass)) {
860  ShiftStack("readobjptr");
861  if (cl) *cl = objClass;
862  return obj;
863  }
864 
865  TString clname = fXML->GetAttr(objnode, xmlio::ObjClass);
866  objClass = XmlDefineClass(clname);
867  if (objClass == TDirectory::Class()) objClass = TDirectoryFile::Class();
868 
869  if (objClass==0) {
870  Error("XmlReadObject", "Cannot find class %s", clname.Data());
871  ShiftStack("readobjerr");
872  return obj;
873  }
874 
875  if (gDebug>1)
876  Info("XmlReadObject", "Reading object of class %s", clname.Data());
877 
878  if (obj==0) obj = objClass->New();
879 
880  ExtractReference(objnode, obj, objClass);
881 
882  PushStack(objnode);
883 
884  objClass->Streamer((void*)obj, *this);
885 
886  PopStack();
887 
888  ShiftStack("readobj");
889 
890  if (gDebug>1)
891  Info("XmlReadObject", "Reading object of class %s done", clname.Data());
892 
893  if (cl) *cl = objClass;
894 
895  return obj;
896 }
897 
898 ////////////////////////////////////////////////////////////////////////////////
899 /// Function is called from TStreamerInfo WriteBuffer and Readbuffer functions
900 /// and indent new level in xml structure.
901 /// This call indicates, that TStreamerInfo functions starts streaming
902 /// object data of correspondent class
903 
905 {
907 }
908 
909 ////////////////////////////////////////////////////////////////////////////////
910 /// Prepares buffer to stream data of specified class.
911 
913 {
916 
917  if (sinfo!=0) cl = sinfo->GetClass();
918 
919  if (cl==0) return;
920 
921  TString clname = XmlConvertClassName(cl->GetName());
922 
923  if (gDebug>2) Info("IncrementLevel","Class: %s", clname.Data());
924 
925  Bool_t compressClassNode = fExpectedBaseClass==cl;
926  fExpectedBaseClass = 0;
927 
928  TXMLStackObj* stack = Stack();
929 
930  if (IsWriting()) {
931 
932  XMLNodePointer_t classnode = 0;
933  if (compressClassNode) {
934  classnode = StackNode();
935  } else {
936  if (GetXmlLayout()==kGeneralized) {
937  classnode = fXML->NewChild(StackNode(), 0, xmlio::Class, 0);
938  fXML->NewAttr(classnode, 0, "name", clname);
939  } else
940  classnode = fXML->NewChild(StackNode(), 0, clname, 0);
941  stack = PushStack(classnode);
942  }
943 
944  if (fVersionBuf>=-1) {
945  if (fVersionBuf == -1) fVersionBuf = 1;
947  fVersionBuf = -111;
948  }
949 
951  stack->fClassNs = fXML->NewNS(classnode, XmlClassNameSpaceRef(cl), clname);
952 
953  } else {
954  if (!compressClassNode) {
955  if (GetXmlLayout()==kGeneralized) {
956  if (!VerifyStackNode(xmlio::Class, "StartInfo")) return;
957  if (!VerifyStackAttr("name", clname, "StartInfo")) return;
958  } else
959  if (!VerifyStackNode(clname, "StartInfo")) return;
960  stack = PushStack(StackNode());
961  }
962  }
963 
964  stack->fCompressedClassNode = compressClassNode;
965  stack->fInfo = sinfo;
966  stack->fIsStreamerInfo = kTRUE;
967 }
968 
969 ////////////////////////////////////////////////////////////////////////////////
970 /// Function is called from TStreamerInfo WriteBuffer and Readbuffer functions
971 /// and decrease level in xml structure.
972 
974 {
975  CheckVersionBuf();
976 
979 
980  if (gDebug>2)
981  Info("DecrementLevel","Class: %s", (info ? info->GetClass()->GetName() : "custom"));
982 
983  TXMLStackObj* stack = Stack();
984 
985  if (!stack->IsStreamerInfo()) {
987  stack = PopStack(); // remove stack of last element
988  }
989 
990  if (stack->fCompressedClassNode) {
991  stack->fInfo = 0;
992  stack->fIsStreamerInfo = kFALSE;
993  stack->fCompressedClassNode = kFALSE;
994  } else {
995  PopStack(); // back from data of stack info
996  if (IsReading()) ShiftStack("declevel"); // shift to next element after streamer info
997  }
998 }
999 
1000 ////////////////////////////////////////////////////////////////////////////////
1001 /// Function is called from TStreamerInfo WriteBuffer and Readbuffer functions
1002 /// and add/verify next element of xml structure
1003 /// This calls allows separate data, correspondent to one class member, from another
1004 
1006 {
1007  WorkWithElement(elem, comptype);
1008 }
1009 
1010 ////////////////////////////////////////////////////////////////////////////////
1011 /// This function is a part of SetStreamerElementNumber method.
1012 /// It is introduced for reading of data for specified data memeber of class.
1013 /// Used also in ReadFastArray methods to resolve problem of compressed data,
1014 /// when several data memebers of the same basic type streamed with single ...FastArray call
1015 
1017 {
1018  CheckVersionBuf();
1019 
1022  fExpectedBaseClass = 0;
1023 
1024  TXMLStackObj* stack = Stack();
1025  if (stack==0) {
1026  Error("SetStreamerElementNumber", "stack is empty");
1027  return;
1028  }
1029 
1030  if (!stack->IsStreamerInfo()) { // this is not a first element
1032  PopStack(); // go level back
1033  if (IsReading()) ShiftStack("startelem"); // shift to next element, only for reading
1034  stack = dynamic_cast<TXMLStackObj*> (fStack.Last());
1035  }
1036 
1037  if (stack==0) {
1038  Error("SetStreamerElementNumber", "Lost of stack");
1039  return;
1040  }
1041 
1042  if (!elem) {
1043  Error("SetStreamerElementNumber", "Problem in Inc/Dec level");
1044  return;
1045  }
1046 
1047  TStreamerInfo* info = stack->fInfo;
1048 
1049  if (!stack->IsStreamerInfo()) {
1050  Error("SetStreamerElementNumber", "Problem in Inc/Dec level");
1051  return;
1052  }
1053  Int_t number = info ? info->GetElements()->IndexOf(elem) : -1;
1054 
1055  if (gDebug>4) Info("SetStreamerElementNumber", " Next element %s", elem->GetName());
1056 
1057  Bool_t isBasicType = (elem->GetType()>0) && (elem->GetType()<20);
1058 
1059  fExpectedChain = isBasicType && (comp_type - elem->GetType() == TStreamerInfo::kOffsetL);
1060 
1061  if (fExpectedChain && (gDebug>3)) {
1062  Info("SetStreamerElementNumber",
1063  " Expects chain for elem %s number %d",
1064  elem->GetName(), number);
1065  }
1066 
1067  fCanUseCompact = isBasicType && ((elem->GetType()==comp_type) ||
1068  (elem->GetType()==comp_type-TStreamerInfo::kConv) ||
1069  (elem->GetType()==comp_type-TStreamerInfo::kSkip));
1070 
1071  if ((elem->GetType()==TStreamerInfo::kBase) ||
1072  ((elem->GetType()==TStreamerInfo::kTNamed) && !strcmp(elem->GetName(), TNamed::Class()->GetName())))
1074 
1075  if (fExpectedBaseClass && (gDebug>3))
1076  Info("SetStreamerElementNumber",
1077  " Expects base class %s with standard streamer",
1079 
1080  if (IsWriting()) {
1081  CreateElemNode(elem);
1082  } else {
1083  if (!VerifyElemNode(elem)) return;
1084  }
1085 
1086  stack = Stack();
1087  stack->fElemNumber = number;
1088  stack->fIsElemOwner = (number<0);
1089 }
1090 
1091 ////////////////////////////////////////////////////////////////////////////////
1092 /// Should be called at the beginning of custom class streamer.
1093 ///
1094 /// Informs buffer data about class which will be streamed now.
1095 /// ClassBegin(), ClassEnd() and ClassMemeber() should be used in
1096 /// custom class streamers to specify which kind of data are
1097 /// now streamed. Such information is used to correctly
1098 /// convert class data to XML. Without that functions calls
1099 /// classes with custom streamers cannot be used with TBufferXML
1100 
1102 {
1103  WorkWithClass(0, cl);
1104 }
1105 
1106 ////////////////////////////////////////////////////////////////////////////////
1107 /// Should be called at the end of custom streamer
1108 /// See TBufferXML::ClassBegin for more details
1109 
1111 {
1112  DecrementLevel(0);
1113 }
1114 
1115 ////////////////////////////////////////////////////////////////////////////////
1116 /// Method indicates name and typename of class member,
1117 /// which should be now streamed in custom streamer
1118 ///
1119 /// Following combinations are supported:
1120 /// -# name = "ClassName", typeName = 0 or typename==ClassName.
1121 /// This is a case, when data of parent class "ClassName" should be streamed.
1122 /// For instance, if class directly inherited from TObject, custom streamer
1123 /// should include following code:
1124 /// ~~~{.cpp}
1125 /// b.ClassMember("TObject");
1126 /// TObject::Streamer(b);
1127 /// ~~~
1128 /// -# Basic data type
1129 /// ~~~{.cpp}
1130 /// b.ClassMember("fInt","Int_t");
1131 /// b >> fInt;
1132 /// ~~~
1133 /// -# Array of basic data types
1134 /// ~~~{.cpp}
1135 /// b.ClassMember("fArr","Int_t", 5);
1136 /// b.ReadFastArray(fArr, 5);
1137 /// ~~~
1138 /// -# Object as data member
1139 /// ~~~{.cpp}
1140 /// b.ClassMemeber("fName","TString");
1141 /// fName.Streamer(b);
1142 /// ~~~
1143 /// -# Pointer on object as data member
1144 /// ~~~{.cpp}
1145 /// b.ClassMemeber("fObj","TObject*");
1146 /// b.StreamObject(fObj);
1147 /// ~~~
1148 ///
1149 /// Arrsize1 and arrsize2 arguments (when specified) indicate first and
1150 /// second dimension of array. Can be used for array of basic types.
1151 /// See ClassBegin() method for more details.
1152 
1153 void TBufferXML::ClassMember(const char* name, const char* typeName, Int_t arrsize1, Int_t arrsize2)
1154 {
1155  if (typeName==0) typeName = name;
1156 
1157  if ((name==0) || (strlen(name)==0)) {
1158  Error("ClassMember","Invalid member name");
1159  fErrorFlag = 1;
1160  return;
1161  }
1162 
1163  TString tname = typeName;
1164 
1165  Int_t typ_id(-1), comp_type(-1);
1166 
1167  if (strcmp(typeName,"raw:data")==0)
1168  typ_id = TStreamerInfo::kMissing;
1169 
1170  if (typ_id<0) {
1171  TDataType *dt = gROOT->GetType(typeName);
1172  if (dt!=0)
1173  if ((dt->GetType()>0) && (dt->GetType()<20))
1174  typ_id = dt->GetType();
1175  }
1176 
1177  if (typ_id<0)
1178  if (strcmp(name, typeName)==0) {
1179  TClass* cl = TClass::GetClass(tname.Data());
1180  if (cl!=0) typ_id = TStreamerInfo::kBase;
1181  }
1182 
1183  if (typ_id<0) {
1184  Bool_t isptr = kFALSE;
1185  if (tname[tname.Length()-1]=='*') {
1186  tname.Resize(tname.Length()-1);
1187  isptr = kTRUE;
1188  }
1189  TClass* cl = TClass::GetClass(tname.Data());
1190  if (cl==0) {
1191  Error("ClassMember","Invalid class specifier %s", typeName);
1192  fErrorFlag = 1;
1193  return;
1194  }
1195 
1196  if (cl->IsTObject())
1198  else
1199  typ_id = isptr ? TStreamerInfo::kAnyp : TStreamerInfo::kAny;
1200 
1201  if ((cl==TString::Class()) && !isptr)
1202  typ_id = TStreamerInfo::kTString;
1203  }
1204 
1205  TStreamerElement* elem = 0;
1206 
1207  if (typ_id == TStreamerInfo::kMissing) {
1208  elem = new TStreamerElement(name,"title",0, typ_id, "raw:data");
1209  } else
1210 
1211  if (typ_id==TStreamerInfo::kBase) {
1212  TClass* cl = TClass::GetClass(tname.Data());
1213  if (cl!=0) {
1214  TStreamerBase* b = new TStreamerBase(tname.Data(), "title", 0);
1215  b->SetBaseVersion(cl->GetClassVersion());
1216  elem = b;
1217  }
1218  } else
1219 
1220  if ((typ_id>0) && (typ_id<20)) {
1221  elem = new TStreamerBasicType(name, "title", 0, typ_id, typeName);
1222  comp_type = typ_id;
1223  } else
1224 
1225  if ((typ_id==TStreamerInfo::kObject) ||
1226  (typ_id==TStreamerInfo::kTObject) ||
1227  (typ_id==TStreamerInfo::kTNamed)) {
1228  elem = new TStreamerObject(name, "title", 0, tname.Data());
1229  } else
1230 
1231  if (typ_id==TStreamerInfo::kObjectp) {
1232  elem = new TStreamerObjectPointer(name, "title", 0, tname.Data());
1233  } else
1234 
1235  if (typ_id==TStreamerInfo::kAny) {
1236  elem = new TStreamerObjectAny(name, "title", 0, tname.Data());
1237  } else
1238 
1239  if (typ_id==TStreamerInfo::kAnyp) {
1240  elem = new TStreamerObjectAnyPointer(name, "title", 0, tname.Data());
1241  } else
1242 
1243  if (typ_id==TStreamerInfo::kTString) {
1244  elem = new TStreamerString(name, "title", 0);
1245  }
1246 
1247  if (elem==0) {
1248  Error("ClassMember","Invalid combination name = %s type = %s", name, typeName);
1249  fErrorFlag = 1;
1250  return;
1251  }
1252 
1253  if (arrsize1>0) {
1254  elem->SetArrayDim(arrsize2>0 ? 2 : 1);
1255  elem->SetMaxIndex(0, arrsize1);
1256  if (arrsize2>0)
1257  elem->SetMaxIndex(1, arrsize2);
1258  }
1259 
1260  // we indicate that there is no streamerinfo
1261  WorkWithElement(elem, comp_type);
1262 }
1263 
1264 ////////////////////////////////////////////////////////////////////////////////
1265 /// Function is converts TObject and TString structures to more compact representation
1266 
1268 {
1269  if (GetXmlLayout()==kGeneralized) return;
1270 
1271  const TStreamerElement* elem = Stack()->fElem;
1272  XMLNodePointer_t elemnode = IsWriting() ? Stack()->fNode : Stack(1)->fNode;
1273 
1274  if ((elem==0) || (elemnode==0)) return;
1275 
1276  if (elem->GetType()==TStreamerInfo::kTString) {
1277 
1278  XMLNodePointer_t node = fXML->GetChild(elemnode);
1279  fXML->SkipEmpty(node);
1280 
1281  XMLNodePointer_t nodecharstar(0), nodeuchar(0), nodeint(0), nodestring(0);
1282 
1283  while (node!=0) {
1284  const char* name = fXML->GetNodeName(node);
1285  if (strcmp(name, xmlio::String)==0) {
1286  if (nodestring) return;
1287  nodestring = node;
1288  } else
1289  if (strcmp(name, xmlio::UChar)==0) {
1290  if (nodeuchar) return;
1291  nodeuchar = node;
1292  } else
1293  if (strcmp(name, xmlio::Int)==0) {
1294  if (nodeint) return;
1295  nodeint = node;
1296  } else
1297  if (strcmp(name, xmlio::CharStar)==0) {
1298  if (nodecharstar!=0) return;
1299  nodecharstar = node;
1300  } else return; // can not be something else
1301  fXML->ShiftToNext(node);
1302  }
1303 
1304  TString str;
1305 
1306  if (GetIOVersion()<3) {
1307  if (nodeuchar==0) return;
1308  if (nodecharstar!=0)
1309  str = fXML->GetAttr(nodecharstar, xmlio::v);
1310  fXML->UnlinkFreeNode(nodeuchar);
1311  fXML->UnlinkFreeNode(nodeint);
1312  fXML->UnlinkFreeNode(nodecharstar);
1313  } else {
1314  if (nodestring!=0)
1315  str = fXML->GetAttr(nodestring, xmlio::v);
1316  fXML->UnlinkFreeNode(nodestring);
1317  }
1318 
1319  fXML->NewAttr(elemnode, 0, "str", str);
1320  } else
1321  if (elem->GetType()==TStreamerInfo::kTObject) {
1322  XMLNodePointer_t node = fXML->GetChild(elemnode);
1323  fXML->SkipEmpty(node);
1324 
1325  XMLNodePointer_t vnode = 0;
1326  XMLNodePointer_t idnode = 0;
1327  XMLNodePointer_t bitsnode = 0;
1328  XMLNodePointer_t prnode = 0;
1329  while (node!=0) {
1330  const char* name = fXML->GetNodeName(node);
1331 
1332  if (strcmp(name, xmlio::OnlyVersion)==0) {
1333  if (vnode) return;
1334  vnode = node;
1335  } else
1336  if (strcmp(name, xmlio::UInt)==0) {
1337  if (idnode==0) idnode = node; else
1338  if (bitsnode==0) bitsnode = node; else return;
1339  } else
1340  if (strcmp(name, xmlio::UShort)==0) {
1341  if (prnode) return;
1342  prnode = node;
1343  } else return;
1344  fXML->ShiftToNext(node);
1345  }
1346 
1347  if ((vnode==0) || (idnode==0) || (bitsnode==0)) return;
1348 
1349  TString str = fXML->GetAttr(idnode,xmlio::v);
1350  fXML->NewAttr(elemnode, 0, "fUniqueID", str);
1351 
1352  str = fXML->GetAttr(bitsnode, xmlio::v);
1353  UInt_t bits;
1354  sscanf(str.Data(),"%u", &bits);
1355 
1356  char sbuf[20];
1357  snprintf(sbuf, sizeof(sbuf), "%x",bits);
1358  fXML->NewAttr(elemnode, 0, "fBits", sbuf);
1359 
1360  if (prnode!=0) {
1361  str = fXML->GetAttr(prnode,xmlio::v);
1362  fXML->NewAttr(elemnode, 0, "fProcessID", str);
1363  }
1364 
1365  fXML->UnlinkFreeNode(vnode);
1366  fXML->UnlinkFreeNode(idnode);
1367  fXML->UnlinkFreeNode(bitsnode);
1368  fXML->UnlinkFreeNode(prnode);
1369  }
1370 }
1371 
1372 ////////////////////////////////////////////////////////////////////////////////
1373 /// Function is unpack TObject and TString structures to be able read
1374 /// them from custom streamers of this objects
1375 
1377 {
1378  if (GetXmlLayout()==kGeneralized) return;
1379  if ((elem==0) || (elemnode==0)) return;
1380 
1381  if (elem->GetType()==TStreamerInfo::kTString) {
1382 
1383  if (!fXML->HasAttr(elemnode,"str")) return;
1384  TString str = fXML->GetAttr(elemnode, "str");
1385  fXML->FreeAttr(elemnode, "str");
1386 
1387  if (GetIOVersion()<3) {
1388  Int_t len = str.Length();
1389  XMLNodePointer_t ucharnode = fXML->NewChild(elemnode, 0, xmlio::UChar,0);
1390  char sbuf[20];
1391  snprintf(sbuf, sizeof(sbuf), "%d", len);
1392  if (len<255) {
1393  fXML->NewAttr(ucharnode,0,xmlio::v,sbuf);
1394  } else {
1395  fXML->NewAttr(ucharnode,0,xmlio::v,"255");
1396  XMLNodePointer_t intnode = fXML->NewChild(elemnode, 0, xmlio::Int, 0);
1397  fXML->NewAttr(intnode, 0, xmlio::v, sbuf);
1398  }
1399  if (len>0) {
1400  XMLNodePointer_t node = fXML->NewChild(elemnode, 0, xmlio::CharStar, 0);
1401  fXML->NewAttr(node, 0, xmlio::v, str);
1402  }
1403  } else {
1404  XMLNodePointer_t node = fXML->NewChild(elemnode, 0, xmlio::String, 0);
1405  fXML->NewAttr(node, 0, xmlio::v, str);
1406  }
1407  } else
1408  if (elem->GetType()==TStreamerInfo::kTObject) {
1409  if (!fXML->HasAttr(elemnode, "fUniqueID")) return;
1410  if (!fXML->HasAttr(elemnode, "fBits")) return;
1411 
1412  TString idstr = fXML->GetAttr(elemnode, "fUniqueID");
1413  TString bitsstr = fXML->GetAttr(elemnode, "fBits");
1414  TString prstr = fXML->GetAttr(elemnode, "fProcessID");
1415 
1416  fXML->FreeAttr(elemnode, "fUniqueID");
1417  fXML->FreeAttr(elemnode, "fBits");
1418  fXML->FreeAttr(elemnode, "fProcessID");
1419 
1420  XMLNodePointer_t node = fXML->NewChild(elemnode, 0, xmlio::OnlyVersion, 0);
1421  fXML->NewAttr(node, 0, xmlio::v, "1");
1422 
1423  node = fXML->NewChild(elemnode, 0, xmlio::UInt, 0);
1424  fXML->NewAttr(node, 0, xmlio::v, idstr);
1425 
1426  UInt_t bits;
1427  sscanf(bitsstr.Data(),"%x", &bits);
1428  char sbuf[20];
1429  snprintf(sbuf, sizeof(sbuf), "%u", bits);
1430 
1431  node = fXML->NewChild(elemnode, 0, xmlio::UInt, 0);
1432  fXML->NewAttr(node, 0, xmlio::v, sbuf);
1433 
1434  if (prstr.Length()>0) {
1435  node = fXML->NewChild(elemnode, 0, xmlio::UShort, 0);
1436  fXML->NewAttr(node, 0, xmlio::v, prstr.Data());
1437  }
1438  }
1439 }
1440 
1441 ////////////////////////////////////////////////////////////////////////////////
1442 /// Function is called before any IO operation of TBuffer
1443 /// Now is used to store version value if no proper calls are discovered
1444 
1446 {
1447  CheckVersionBuf();
1448 }
1449 
1450 ////////////////////////////////////////////////////////////////////////////////
1451 /// Function to read class from buffer, used in old-style streamers
1452 
1454 {
1455  const char* clname = 0;
1456 
1458  clname = XmlReadValue(xmlio::Class);
1459  }
1460 
1461  if (gDebug>2) Info("ReadClass", "Try to read class %s", clname ? clname : "---");
1462 
1463  return clname ? gROOT->GetClass(clname) : 0;
1464 }
1465 
1466 ////////////////////////////////////////////////////////////////////////////////
1467 /// Function to write class into buffer, used in old-style streamers
1468 
1470 {
1471  if (gDebug>2) Info("WriteClass", "Try to write class %s", cl->GetName());
1472 
1474 }
1475 
1476 ////////////////////////////////////////////////////////////////////////////////
1477 /// Suppressed function of TBuffer
1478 
1479 Int_t TBufferXML::CheckByteCount(UInt_t /*r_s */, UInt_t /*r_c*/, const TClass* /*cl*/)
1480 {
1481  return 0;
1482 }
1483 
1484 ////////////////////////////////////////////////////////////////////////////////
1485 /// Suppressed function of TBuffer
1486 
1488 {
1489  return 0;
1490 }
1491 
1492 ////////////////////////////////////////////////////////////////////////////////
1493 /// Suppressed function of TBuffer
1494 
1496 {
1497 }
1498 
1499 ////////////////////////////////////////////////////////////////////////////////
1500 /// Skip class version from I/O buffer.
1501 
1503 {
1504  ReadVersion(0,0,cl);
1505 }
1506 
1507 ////////////////////////////////////////////////////////////////////////////////
1508 /// Read version value from buffer
1509 
1510 Version_t TBufferXML::ReadVersion(UInt_t *start, UInt_t *bcnt, const TClass * /*cl*/)
1511 {
1513 
1514  Version_t res = 0;
1515 
1516  if (start) *start = 0;
1517  if (bcnt) *bcnt = 0;
1518 
1521  } else
1522  if ((fExpectedBaseClass!=0) && (fXML->HasAttr(Stack(1)->fNode, xmlio::ClassVersion))) {
1523  res = fXML->GetIntAttr(Stack(1)->fNode, xmlio::ClassVersion);
1524  } else
1527  } else {
1528  Error("ReadVersion", "No correspondent tags to read version");;
1529  fErrorFlag = 1;
1530  }
1531 
1532  if (gDebug>2) Info("ReadVersion","Version = %d", res);
1533 
1534  return res;
1535 }
1536 
1537 ////////////////////////////////////////////////////////////////////////////////
1538 /// Checks buffer, filled by WriteVersion
1539 /// if next data is arriving, version should be stored in buffer
1540 
1542 {
1543  if (IsWriting() && (fVersionBuf>=-100)) {
1544  char sbuf[20];
1545  snprintf(sbuf, sizeof(sbuf), "%d", fVersionBuf);
1547  fVersionBuf = -111;
1548  }
1549 }
1550 
1551 ////////////////////////////////////////////////////////////////////////////////
1552 /// Copies class version to buffer, but not writes it to xml
1553 /// Version will be written with next I/O operation or
1554 /// will be added as attribute of class tag, created by IncrementLevel call
1555 
1557 {
1559 
1560  if (fExpectedBaseClass!=cl)
1561  fExpectedBaseClass = 0;
1562 
1563  fVersionBuf = cl->GetClassVersion();
1564 
1565  if (gDebug>2)
1566  Info("WriteVersion", "Class: %s, version = %d",
1567  cl->GetName(), fVersionBuf);
1568 
1569  return 0;
1570 }
1571 
1572 ////////////////////////////////////////////////////////////////////////////////
1573 /// Read object from buffer. Only used from TBuffer
1574 
1576 {
1578  if (gDebug>2)
1579  Info("ReadObjectAny","From node %s", fXML->GetNodeName(StackNode()));
1580  void* res = XmlReadObject(0);
1581  return res;
1582 }
1583 
1584 ////////////////////////////////////////////////////////////////////////////////
1585 /// Skip any kind of object from buffer
1586 /// Actually skip only one node on current level of xml structure
1587 
1589 {
1590  ShiftStack("skipobjectany"); \
1591 }
1592 
1593 ////////////////////////////////////////////////////////////////////////////////
1594 /// Write object to buffer. Only used from TBuffer
1595 
1596 void TBufferXML::WriteObjectClass(const void *actualObjStart, const TClass *actualClass)
1597 {
1599  if (gDebug>2)
1600  Info("WriteObject","Class %s", (actualClass ? actualClass->GetName() : " null"));
1601  XmlWriteObject(actualObjStart, actualClass);
1602 }
1603 
1604 // Macro to read content of uncompressed array
1605 #define TXMLReadArrayNoncompress(vname) \
1606 { \
1607  for(Int_t indx=0;indx<n;indx++) \
1608  XmlReadBasic(vname[indx]); \
1609 }
1610 
1611 // macro to read content of array with compression
1612 #define TXMLReadArrayContent(vname, arrsize) \
1613 { \
1614  Int_t indx = 0; \
1615  while(indx<arrsize) { \
1616  Int_t cnt = 1; \
1617  if (fXML->HasAttr(StackNode(), xmlio::cnt)) \
1618  cnt = fXML->GetIntAttr(StackNode(), xmlio::cnt); \
1619  XmlReadBasic(vname[indx]); \
1620  Int_t curr = indx; indx++; \
1621  while(cnt>1) { \
1622  vname[indx] = vname[curr]; \
1623  cnt--; indx++; \
1624  } \
1625  } \
1626 }
1627 
1628 // macro to read array, which include size attribute
1629 #define TBufferXML_ReadArray(tname, vname) \
1630 { \
1631  BeforeIOoperation(); \
1632  if (!VerifyItemNode(xmlio::Array,"ReadArray")) return 0; \
1633  Int_t n = fXML->GetIntAttr(StackNode(), xmlio::Size); \
1634  if (n<=0) return 0; \
1635  if (!vname) vname = new tname[n]; \
1636  PushStack(StackNode()); \
1637  TXMLReadArrayContent(vname, n); \
1638  PopStack(); \
1639  ShiftStack("readarr"); \
1640  return n; \
1641 }
1642 
1643 ////////////////////////////////////////////////////////////////////////////////
1644 /// Read a Float16_t from the buffer
1645 
1647 {
1649  XmlReadBasic(*f);
1650 }
1651 
1652 ////////////////////////////////////////////////////////////////////////////////
1653 /// Read a Double32_t from the buffer
1654 
1656 {
1658  XmlReadBasic(*d);
1659 }
1660 
1661 ////////////////////////////////////////////////////////////////////////////////
1662 /// Read a Double32_t from the buffer when the factor and minimun value have been specified
1663 /// see comments about Double32_t encoding at TBufferFile::WriteDouble32().
1664 /// Currently TBufferXML does not optimize space in this case.
1665 
1666 void TBufferXML::ReadWithFactor(Float_t *ptr, Double_t /* factor */, Double_t /* minvalue */)
1667 {
1669  XmlReadBasic(*ptr);
1670 }
1671 
1672 ////////////////////////////////////////////////////////////////////////////////
1673 /// Read a Float16_t from the buffer when the number of bits is specified (explicitly or not)
1674 /// see comments about Float16_t encoding at TBufferFile::WriteFloat16().
1675 /// Currently TBufferXML does not optimize space in this case.
1676 
1677 void TBufferXML::ReadWithNbits(Float_t *ptr, Int_t /* nbits */)
1678 {
1680  XmlReadBasic(*ptr);
1681 }
1682 
1683 ////////////////////////////////////////////////////////////////////////////////
1684 /// Read a Double32_t from the buffer when the factor and minimun value have been specified
1685 /// see comments about Double32_t encoding at TBufferFile::WriteDouble32().
1686 /// Currently TBufferXML does not optimize space in this case.
1687 
1688 void TBufferXML::ReadWithFactor(Double_t *ptr, Double_t /* factor */, Double_t /* minvalue */)
1689 {
1691  XmlReadBasic(*ptr);
1692 }
1693 
1694 ////////////////////////////////////////////////////////////////////////////////
1695 /// Read a Double32_t from the buffer when the number of bits is specified (explicitly or not)
1696 /// see comments about Double32_t encoding at TBufferFile::WriteDouble32().
1697 /// Currently TBufferXML does not optimize space in this case.
1698 
1700 {
1702  XmlReadBasic(*ptr);
1703 }
1704 
1705 ////////////////////////////////////////////////////////////////////////////////
1706 /// Write a Float16_t to the buffer
1707 
1709 {
1711  XmlWriteBasic(*f);
1712 }
1713 
1714 ////////////////////////////////////////////////////////////////////////////////
1715 /// Write a Double32_t to the buffer
1716 
1718 {
1720  XmlWriteBasic(*d);
1721 }
1722 
1723 ////////////////////////////////////////////////////////////////////////////////
1724 /// Read array of Bool_t from buffer
1725 
1727 {
1729 }
1730 
1731 ////////////////////////////////////////////////////////////////////////////////
1732 /// Read array of Char_t from buffer
1733 
1735 {
1737 }
1738 
1739 ////////////////////////////////////////////////////////////////////////////////
1740 /// Read array of UChar_t from buffer
1741 
1743 {
1745 }
1746 
1747 ////////////////////////////////////////////////////////////////////////////////
1748 /// Read array of Short_t from buffer
1749 
1751 {
1753 }
1754 
1755 ////////////////////////////////////////////////////////////////////////////////
1756 /// Read array of UShort_t from buffer
1757 
1759 {
1761 }
1762 
1763 ////////////////////////////////////////////////////////////////////////////////
1764 /// Read array of Int_t from buffer
1765 
1767 {
1769 }
1770 
1771 ////////////////////////////////////////////////////////////////////////////////
1772 /// Read array of UInt_t from buffer
1773 
1775 {
1777 }
1778 
1779 ////////////////////////////////////////////////////////////////////////////////
1780 /// Read array of Long_t from buffer
1781 
1783 {
1785 }
1786 
1787 ////////////////////////////////////////////////////////////////////////////////
1788 /// Read array of ULong_t from buffer
1789 
1791 {
1793 }
1794 
1795 ////////////////////////////////////////////////////////////////////////////////
1796 /// Read array of Long64_t from buffer
1797 
1799 {
1801 }
1802 
1803 ////////////////////////////////////////////////////////////////////////////////
1804 /// Read array of ULong64_t from buffer
1805 
1807 {
1809 }
1810 
1811 ////////////////////////////////////////////////////////////////////////////////
1812 /// Read array of Float_t from buffer
1813 
1815 {
1817 }
1818 
1819 ////////////////////////////////////////////////////////////////////////////////
1820 /// Read array of Double_t from buffer
1821 
1823 {
1825 }
1826 
1827 ////////////////////////////////////////////////////////////////////////////////
1828 /// Read array of Float16_t from buffer
1829 
1831 {
1833 }
1834 
1835 ////////////////////////////////////////////////////////////////////////////////
1836 /// Read array of Double32_t from buffer
1837 
1839 {
1841 }
1842 
1843 // macro to read array from xml buffer
1844 #define TBufferXML_ReadStaticArray(vname) \
1845 { \
1846  BeforeIOoperation(); \
1847  if (!VerifyItemNode(xmlio::Array,"ReadStaticArray")) return 0; \
1848  Int_t n = fXML->GetIntAttr(StackNode(), xmlio::Size); \
1849  if (n<=0) return 0; \
1850  if (!vname) return 0; \
1851  PushStack(StackNode()); \
1852  TXMLReadArrayContent(vname, n); \
1853  PopStack(); \
1854  ShiftStack("readstatarr"); \
1855  return n; \
1856 }
1857 
1858 ////////////////////////////////////////////////////////////////////////////////
1859 /// Read array of Bool_t from buffer
1860 
1862 {
1864 }
1865 
1866 ////////////////////////////////////////////////////////////////////////////////
1867 /// Read array of Char_t from buffer
1868 
1870 {
1872 }
1873 
1874 ////////////////////////////////////////////////////////////////////////////////
1875 /// Read array of UChar_t from buffer
1876 
1878 {
1880 }
1881 
1882 ////////////////////////////////////////////////////////////////////////////////
1883 /// Read array of Short_t from buffer
1884 
1886 {
1888 }
1889 
1890 ////////////////////////////////////////////////////////////////////////////////
1891 /// Read array of UShort_t from buffer
1892 
1894 {
1896 }
1897 
1898 ////////////////////////////////////////////////////////////////////////////////
1899 /// Read array of Int_t from buffer
1900 
1902 {
1904 }
1905 
1906 ////////////////////////////////////////////////////////////////////////////////
1907 /// Read array of UInt_t from buffer
1908 
1910 {
1912 }
1913 
1914 ////////////////////////////////////////////////////////////////////////////////
1915 /// Read array of Long_t from buffer
1916 
1918 {
1920 }
1921 
1922 ////////////////////////////////////////////////////////////////////////////////
1923 /// Read array of ULong_t from buffer
1924 
1926 {
1928 }
1929 
1930 ////////////////////////////////////////////////////////////////////////////////
1931 /// Read array of Long64_t from buffer
1932 
1934 {
1936 }
1937 
1938 ////////////////////////////////////////////////////////////////////////////////
1939 /// Read array of ULong64_t from buffer
1940 
1942 {
1944 }
1945 
1946 ////////////////////////////////////////////////////////////////////////////////
1947 /// Read array of Float_t from buffer
1948 
1950 {
1952 }
1953 
1954 ////////////////////////////////////////////////////////////////////////////////
1955 /// Read array of Double_t from buffer
1956 
1958 {
1960 }
1961 
1962 ////////////////////////////////////////////////////////////////////////////////
1963 /// Read array of Float16_t from buffer
1964 
1966 {
1968 }
1969 
1970 ////////////////////////////////////////////////////////////////////////////////
1971 /// Read array of Double32_t from buffer
1972 
1974 {
1976 }
1977 
1978 // macro to read content of array, which not include size of array
1979 // macro also treat situation, when instead of one single array chain
1980 // of several elements should be produced
1981 #define TBufferXML_ReadFastArray(vname) \
1982 { \
1983  BeforeIOoperation(); \
1984  if (n<=0) return; \
1985  TStreamerElement* elem = Stack(0)->fElem; \
1986  if ((elem!=0) && (elem->GetType()>TStreamerInfo::kOffsetL) && \
1987  (elem->GetType()<TStreamerInfo::kOffsetP) && \
1988  (elem->GetArrayLength()!=n)) fExpectedChain = kTRUE; \
1989  if (fExpectedChain) { \
1990  fExpectedChain = kFALSE; \
1991  Int_t startnumber = Stack(0)->fElemNumber; \
1992  TStreamerInfo* info = Stack(1)->fInfo; \
1993  Int_t index = 0; \
1994  while (index<n) { \
1995  elem = (TStreamerElement*)info->GetElements()->At(startnumber++); \
1996  if (elem->GetType()<TStreamerInfo::kOffsetL) { \
1997  if (index>0) { PopStack(); ShiftStack("chainreader"); VerifyElemNode(elem); } \
1998  fCanUseCompact = kTRUE; \
1999  XmlReadBasic(vname[index]); \
2000  index++; \
2001  } else { \
2002  if (!VerifyItemNode(xmlio::Array,"ReadFastArray")) return; \
2003  PushStack(StackNode()); \
2004  Int_t elemlen = elem->GetArrayLength(); \
2005  TXMLReadArrayContent((vname+index), elemlen); \
2006  PopStack(); \
2007  ShiftStack("readfastarr"); \
2008  index+=elemlen; \
2009  } \
2010  } \
2011  } else { \
2012  if (!VerifyItemNode(xmlio::Array,"ReadFastArray")) return; \
2013  PushStack(StackNode()); \
2014  TXMLReadArrayContent(vname, n); \
2015  PopStack(); \
2016  ShiftStack("readfastarr"); \
2017  } \
2018 }
2019 
2020 ////////////////////////////////////////////////////////////////////////////////
2021 /// Read array of Bool_t from buffer
2022 
2024 {
2026 }
2027 
2028 ////////////////////////////////////////////////////////////////////////////////
2029 /// Read array of Char_t from buffer
2030 /// if nodename==CharStar, read all array as string
2031 
2033 {
2034  if ((n>0) && VerifyItemNode(xmlio::CharStar)) {
2035  const char* buf;
2036  if ((buf = XmlReadValue(xmlio::CharStar))) {
2037  Int_t size = strlen(buf);
2038  if (size<n) size = n;
2039  memcpy(c, buf, size);
2040  }
2041  } else
2043 }
2044 
2045 ////////////////////////////////////////////////////////////////////////////////
2046 /// Read array of UChar_t from buffer
2047 
2049 {
2051 }
2052 
2053 ////////////////////////////////////////////////////////////////////////////////
2054 /// Read array of Short_t from buffer
2055 
2057 {
2059 }
2060 
2061 ////////////////////////////////////////////////////////////////////////////////
2062 /// Read array of UShort_t from buffer
2063 
2065 {
2067 }
2068 
2069 ////////////////////////////////////////////////////////////////////////////////
2070 /// Read array of Int_t from buffer
2071 
2073 {
2075 }
2076 
2077 ////////////////////////////////////////////////////////////////////////////////
2078 /// Read array of UInt_t from buffer
2079 
2081 {
2083 }
2084 
2085 ////////////////////////////////////////////////////////////////////////////////
2086 /// Read array of Long_t from buffer
2087 
2089 {
2091 }
2092 
2093 ////////////////////////////////////////////////////////////////////////////////
2094 /// Read array of ULong_t from buffer
2095 
2097 {
2099 }
2100 
2101 ////////////////////////////////////////////////////////////////////////////////
2102 /// Read array of Long64_t from buffer
2103 
2105 {
2107 }
2108 
2109 ////////////////////////////////////////////////////////////////////////////////
2110 /// Read array of ULong64_t from buffer
2111 
2113 {
2115 }
2116 
2117 ////////////////////////////////////////////////////////////////////////////////
2118 /// Read array of Float_t from buffer
2119 
2121 {
2123 }
2124 
2125 ////////////////////////////////////////////////////////////////////////////////
2126 /// Read array of Double_t from buffer
2127 
2129 {
2131 }
2132 
2133 ////////////////////////////////////////////////////////////////////////////////
2134 /// Read array of Float16_t from buffer
2135 
2137 {
2139 }
2140 
2141 ////////////////////////////////////////////////////////////////////////////////
2142 /// Read array of Float16_t from buffer
2143 
2145 {
2147 }
2148 
2149 ////////////////////////////////////////////////////////////////////////////////
2150 /// Read array of Float16_t from buffer
2151 
2153 {
2155 }
2156 
2157 ////////////////////////////////////////////////////////////////////////////////
2158 /// Read array of Double32_t from buffer
2159 
2161 {
2163 }
2164 
2165 ////////////////////////////////////////////////////////////////////////////////
2166 /// Read array of Double32_t from buffer
2167 
2168 void TBufferXML::ReadFastArrayWithFactor(Double_t *d, Int_t n, Double_t /* factor */, Double_t /* minvalue */)
2169 {
2171 }
2172 
2173 ////////////////////////////////////////////////////////////////////////////////
2174 /// Read array of Double32_t from buffer
2175 
2177 {
2179 }
2180 
2181 ////////////////////////////////////////////////////////////////////////////////
2182 /// redefined here to avoid warning message from gcc
2183 
2184 void TBufferXML::ReadFastArray(void *start, const TClass *cl, Int_t n, TMemberStreamer *s, const TClass *onFileClass)
2185 {
2186  TBufferFile::ReadFastArray(start, cl, n, s, onFileClass);
2187 }
2188 
2189 ////////////////////////////////////////////////////////////////////////////////
2190 /// redefined here to avoid warning message from gcc
2191 
2192 void TBufferXML::ReadFastArray(void **startp, const TClass *cl, Int_t n, Bool_t isPreAlloc, TMemberStreamer *s, const TClass *onFileClass)
2193 {
2194  TBufferFile::ReadFastArray(startp, cl, n, isPreAlloc, s, onFileClass);
2195 }
2196 
2197 // macro to write content of noncompressed array
2198 #define TXMLWriteArrayNoncompress(vname, arrsize) \
2199 { \
2200  for(Int_t indx=0;indx<arrsize;indx++) \
2201  XmlWriteBasic(vname[indx]); \
2202 }
2203 
2204 // macro to write content of compressed array
2205 #define TXMLWriteArrayCompress(vname, arrsize) \
2206 { \
2207  Int_t indx = 0; \
2208  while(indx<arrsize) { \
2209  XMLNodePointer_t elemnode = XmlWriteBasic(vname[indx]); \
2210  Int_t curr = indx; indx++; \
2211  while ((indx<arrsize) && (vname[indx]==vname[curr])) indx++; \
2212  if (indx-curr > 1) \
2213  fXML->NewIntAttr(elemnode, xmlio::cnt, indx-curr); \
2214  } \
2215 }
2216 
2217 #define TXMLWriteArrayContent(vname, arrsize) \
2218 { \
2219  if (fCompressLevel>0) { \
2220  TXMLWriteArrayCompress(vname, arrsize) \
2221  } else { \
2222  TXMLWriteArrayNoncompress(vname, arrsize) \
2223  } \
2224 }
2225 
2226 // macro to write array, which include size
2227 #define TBufferXML_WriteArray(vname) \
2228 { \
2229  BeforeIOoperation(); \
2230  XMLNodePointer_t arrnode = CreateItemNode(xmlio::Array); \
2231  fXML->NewIntAttr(arrnode, xmlio::Size, n); \
2232  PushStack(arrnode); \
2233  TXMLWriteArrayContent(vname, n); \
2234  PopStack(); \
2235 }
2236 
2237 ////////////////////////////////////////////////////////////////////////////////
2238 /// Write array of Bool_t to buffer
2239 
2241 {
2243 }
2244 
2245 ////////////////////////////////////////////////////////////////////////////////
2246 /// Write array of Char_t to buffer
2247 
2249 {
2251 }
2252 
2253 ////////////////////////////////////////////////////////////////////////////////
2254 /// Write array of UChar_t to buffer
2255 
2257 {
2259 }
2260 
2261 ////////////////////////////////////////////////////////////////////////////////
2262 /// Write array of Short_t to buffer
2263 
2265 {
2267 }
2268 
2269 ////////////////////////////////////////////////////////////////////////////////
2270 /// Write array of UShort_t to buffer
2271 
2273 {
2275 }
2276 
2277 ////////////////////////////////////////////////////////////////////////////////
2278 /// Write array of Int_ to buffer
2279 
2281 {
2283 }
2284 
2285 ////////////////////////////////////////////////////////////////////////////////
2286 /// Write array of UInt_t to buffer
2287 
2289 {
2291 }
2292 
2293 ////////////////////////////////////////////////////////////////////////////////
2294 /// Write array of Long_t to buffer
2295 
2297 {
2299 }
2300 
2301 ////////////////////////////////////////////////////////////////////////////////
2302 /// Write array of ULong_t to buffer
2303 
2305 {
2307 }
2308 
2309 ////////////////////////////////////////////////////////////////////////////////
2310 /// Write array of Long64_t to buffer
2311 
2313 {
2315 }
2316 
2317 ////////////////////////////////////////////////////////////////////////////////
2318 /// Write array of ULong64_t to buffer
2319 
2321 {
2323 }
2324 
2325 ////////////////////////////////////////////////////////////////////////////////
2326 /// Write array of Float_t to buffer
2327 
2329 {
2331 }
2332 
2333 ////////////////////////////////////////////////////////////////////////////////
2334 /// Write array of Double_t to buffer
2335 
2337 {
2339 }
2340 
2341 ////////////////////////////////////////////////////////////////////////////////
2342 /// Write array of Float16_t to buffer
2343 
2345 {
2347 }
2348 
2349 ////////////////////////////////////////////////////////////////////////////////
2350 /// Write array of Double32_t to buffer
2351 
2353 {
2355 }
2356 
2357 // write array without size attribute
2358 // macro also treat situation, when instead of one single array
2359 // chain of several elements should be produced
2360 #define TBufferXML_WriteFastArray(vname) \
2361 { \
2362  BeforeIOoperation(); \
2363  if (n<=0) return; \
2364  TStreamerElement* elem = Stack(0)->fElem; \
2365  if ((elem!=0) && (elem->GetType()>TStreamerInfo::kOffsetL) && \
2366  (elem->GetType()<TStreamerInfo::kOffsetP) && \
2367  (elem->GetArrayLength()!=n)) fExpectedChain = kTRUE; \
2368  if (fExpectedChain) { \
2369  TStreamerInfo* info = Stack(1)->fInfo; \
2370  Int_t startnumber = Stack(0)->fElemNumber; \
2371  fExpectedChain = kFALSE; \
2372  Int_t index = 0; \
2373  while (index<n) { \
2374  elem =(TStreamerElement*)info->GetElements()->At(startnumber++); \
2375  if (elem->GetType()<TStreamerInfo::kOffsetL) { \
2376  if(index>0) { PopStack(); CreateElemNode(elem); } \
2377  fCanUseCompact = kTRUE; \
2378  XmlWriteBasic(vname[index]); \
2379  index++; \
2380  } else { \
2381  XMLNodePointer_t arrnode = CreateItemNode(xmlio::Array); \
2382  Int_t elemlen = elem->GetArrayLength(); \
2383  PushStack(arrnode); \
2384  TXMLWriteArrayContent((vname+index), elemlen); \
2385  index+=elemlen; \
2386  PopStack(); \
2387  } \
2388  } \
2389  } else { \
2390  XMLNodePointer_t arrnode = CreateItemNode(xmlio::Array); \
2391  PushStack(arrnode); \
2392  TXMLWriteArrayContent(vname, n); \
2393  PopStack(); \
2394  } \
2395 }
2396 
2397 ////////////////////////////////////////////////////////////////////////////////
2398 /// Write array of Bool_t to buffer
2399 
2401 {
2403 }
2404 
2405 ////////////////////////////////////////////////////////////////////////////////
2406 /// Write array of Char_t to buffer
2407 /// If array does not include any special characters,
2408 /// it will be reproduced as CharStar node with string as attribute
2409 
2411 {
2412  Bool_t usedefault = (n==0) || fExpectedChain;
2413  const Char_t* buf = c;
2414  if (!usedefault)
2415  for (int i=0;i<n;i++) {
2416  if (*buf < 27) { usedefault = kTRUE; break; }
2417  buf++;
2418  }
2419  if (usedefault) {
2421  } else {
2422  Char_t* buf2 = new Char_t[n+1];
2423  memcpy(buf2, c, n);
2424  buf2[n] = 0;
2426  delete[] buf2;
2427  }
2428 }
2429 
2430 ////////////////////////////////////////////////////////////////////////////////
2431 /// Write array of UChar_t to buffer
2432 
2434 {
2436 }
2437 
2438 ////////////////////////////////////////////////////////////////////////////////
2439 /// Write array of Short_t to buffer
2440 
2442 {
2444 }
2445 
2446 ////////////////////////////////////////////////////////////////////////////////
2447 /// Write array of UShort_t to buffer
2448 
2450 {
2452 }
2453 
2454 ////////////////////////////////////////////////////////////////////////////////
2455 /// Write array of Int_t to buffer
2456 
2458 {
2460 }
2461 
2462 ////////////////////////////////////////////////////////////////////////////////
2463 /// Write array of UInt_t to buffer
2464 
2466 {
2468 }
2469 
2470 ////////////////////////////////////////////////////////////////////////////////
2471 /// Write array of Long_t to buffer
2472 
2474 {
2476 }
2477 
2478 ////////////////////////////////////////////////////////////////////////////////
2479 /// Write array of ULong_t to buffer
2480 
2482 {
2484 }
2485 
2486 ////////////////////////////////////////////////////////////////////////////////
2487 /// Write array of Long64_t to buffer
2488 
2490 {
2492 }
2493 
2494 ////////////////////////////////////////////////////////////////////////////////
2495 /// Write array of ULong64_t to buffer
2496 
2498 {
2500 }
2501 
2502 ////////////////////////////////////////////////////////////////////////////////
2503 /// Write array of Float_t to buffer
2504 
2506 {
2508 }
2509 
2510 ////////////////////////////////////////////////////////////////////////////////
2511 /// Write array of Double_t to buffer
2512 
2514 {
2516 }
2517 
2518 ////////////////////////////////////////////////////////////////////////////////
2519 /// Write array of Float16_t to buffer
2520 
2522 {
2524 }
2525 
2526 ////////////////////////////////////////////////////////////////////////////////
2527 /// Write array of Double32_t to buffer
2528 
2530 {
2532 }
2533 
2534 ////////////////////////////////////////////////////////////////////////////////
2535 /// Recall TBuffer function to avoid gcc warning message
2536 
2537 void TBufferXML::WriteFastArray(void *start, const TClass *cl, Int_t n, TMemberStreamer *s)
2538 {
2539  TBufferFile::WriteFastArray(start, cl, n, s);
2540 }
2541 
2542 ////////////////////////////////////////////////////////////////////////////////
2543 /// Recall TBuffer function to avoid gcc warning message
2544 
2545 Int_t TBufferXML::WriteFastArray(void **startp, const TClass *cl, Int_t n, Bool_t isPreAlloc, TMemberStreamer *s)
2546 {
2547  return TBufferFile::WriteFastArray(startp, cl, n, isPreAlloc, s);
2548 }
2549 
2550 ////////////////////////////////////////////////////////////////////////////////
2551 /// steram object to/from buffer
2552 
2553 void TBufferXML::StreamObject(void *obj, const std::type_info &typeinfo, const TClass* /* onFileClass */ )
2554 {
2555  StreamObject(obj, TClass::GetClass(typeinfo));
2556 }
2557 
2558 ////////////////////////////////////////////////////////////////////////////////
2559 /// steram object to/from buffer
2560 
2561 void TBufferXML::StreamObject(void *obj, const char *className, const TClass* /* onFileClass */ )
2562 {
2563  StreamObject(obj, TClass::GetClass(className));
2564 }
2565 
2567 {
2568  // steram object to/from buffer
2569 
2570  StreamObject(obj, obj ? obj->IsA() : TObject::Class());
2571 }
2572 
2573 ////////////////////////////////////////////////////////////////////////////////
2574 /// Steram object to/from buffer
2575 
2576 void TBufferXML::StreamObject(void *obj, const TClass *cl, const TClass* /* onfileClass */ )
2577 {
2579  if (gDebug>1)
2580  Info("StreamObject","Class: %s", (cl ? cl->GetName() : "none"));
2581  if (IsReading())
2582  XmlReadObject(obj);
2583  else
2584  XmlWriteObject(obj, cl);
2585 }
2586 
2587 // macro for right shift operator for basic type
2588 #define TBufferXML_operatorin(vname) \
2589 { \
2590  BeforeIOoperation(); \
2591  XmlReadBasic(vname); \
2592 }
2593 
2594 ////////////////////////////////////////////////////////////////////////////////
2595 /// Reads Bool_t value from buffer
2596 
2598 {
2600 }
2601 
2602 ////////////////////////////////////////////////////////////////////////////////
2603 /// Reads Char_t value from buffer
2604 
2606 {
2608 }
2609 
2610 ////////////////////////////////////////////////////////////////////////////////
2611 /// Reads UChar_t value from buffer
2612 
2614 {
2616 }
2617 
2618 ////////////////////////////////////////////////////////////////////////////////
2619 /// Reads Short_t value from buffer
2620 
2622 {
2624 }
2625 
2626 ////////////////////////////////////////////////////////////////////////////////
2627 /// Reads UShort_t value from buffer
2628 
2630 {
2632 }
2633 
2634 ////////////////////////////////////////////////////////////////////////////////
2635 /// Reads Int_t value from buffer
2636 
2638 {
2640 }
2641 
2642 ////////////////////////////////////////////////////////////////////////////////
2643 /// Reads UInt_t value from buffer
2644 
2646 {
2648 }
2649 
2650 ////////////////////////////////////////////////////////////////////////////////
2651 /// Reads Long_t value from buffer
2652 
2654 {
2656 }
2657 
2658 ////////////////////////////////////////////////////////////////////////////////
2659 /// Reads ULong_t value from buffer
2660 
2662 {
2664 }
2665 
2666 ////////////////////////////////////////////////////////////////////////////////
2667 /// Reads Long64_t value from buffer
2668 
2670 {
2672 }
2673 
2674 ////////////////////////////////////////////////////////////////////////////////
2675 /// Reads ULong64_t value from buffer
2676 
2678 {
2680 }
2681 
2682 ////////////////////////////////////////////////////////////////////////////////
2683 /// Reads Float_t value from buffer
2684 
2686 {
2688 }
2689 
2690 ////////////////////////////////////////////////////////////////////////////////
2691 /// Reads Double_t value from buffer
2692 
2694 {
2696 }
2697 
2698 ////////////////////////////////////////////////////////////////////////////////
2699 /// Reads array of characters from buffer
2700 
2702 {
2704  const char* buf;
2705  if ((buf = XmlReadValue(xmlio::CharStar)))
2706  strcpy(c, buf);
2707 }
2708 
2709 ////////////////////////////////////////////////////////////////////////////////
2710 /// Reads a TString
2711 
2713 {
2714  if (GetIOVersion()<3) {
2716  } else {
2718  const char* buf;
2719  if ((buf = XmlReadValue(xmlio::String)))
2720  s = buf;
2721  }
2722 }
2723 
2724 ////////////////////////////////////////////////////////////////////////////////
2725 /// Reads a std::string
2726 
2727 void TBufferXML::ReadStdString(std::string *s)
2728 {
2729  if (GetIOVersion()<3) {
2731  } else {
2733  const char* buf;
2734  if ((buf = XmlReadValue(xmlio::String)))
2735  if (s) *s = buf;
2736  }
2737 }
2738 
2739 ////////////////////////////////////////////////////////////////////////////////
2740 /// Read a char* string
2741 
2743 {
2745 }
2746 
2747 
2748 // macro for left shift operator for basic types
2749 #define TBufferXML_operatorout(vname) \
2750 { \
2751  BeforeIOoperation(); \
2752  XmlWriteBasic(vname); \
2753 }
2754 
2755 ////////////////////////////////////////////////////////////////////////////////
2756 /// Writes Bool_t value to buffer
2757 
2759 {
2761 }
2762 
2763 ////////////////////////////////////////////////////////////////////////////////
2764 /// Writes Char_t value to buffer
2765 
2767 {
2769 }
2770 
2771 ////////////////////////////////////////////////////////////////////////////////
2772 /// Writes UChar_t value to buffer
2773 
2775 {
2777 }
2778 
2779 ////////////////////////////////////////////////////////////////////////////////
2780 /// Writes Short_t value to buffer
2781 
2783 {
2785 }
2786 
2787 ////////////////////////////////////////////////////////////////////////////////
2788 /// Writes UShort_t value to buffer
2789 
2791 {
2793 }
2794 
2795 ////////////////////////////////////////////////////////////////////////////////
2796 /// Writes Int_t value to buffer
2797 
2799 {
2801 }
2802 
2803 ////////////////////////////////////////////////////////////////////////////////
2804 /// Writes UInt_t value to buffer
2805 
2807 {
2809 }
2810 
2811 ////////////////////////////////////////////////////////////////////////////////
2812 /// Writes Long_t value to buffer
2813 
2815 {
2817 }
2818 
2819 ////////////////////////////////////////////////////////////////////////////////
2820 /// Writes ULong_t value to buffer
2821 
2823 {
2825 }
2826 
2827 ////////////////////////////////////////////////////////////////////////////////
2828 /// Writes Long64_t value to buffer
2829 
2831 {
2833 }
2834 
2835 ////////////////////////////////////////////////////////////////////////////////
2836 /// Writes ULong64_t value to buffer
2837 
2839 {
2841 }
2842 
2843 ////////////////////////////////////////////////////////////////////////////////
2844 /// Writes Float_t value to buffer
2845 
2847 {
2849 }
2850 
2851 ////////////////////////////////////////////////////////////////////////////////
2852 /// Writes Double_t value to buffer
2853 
2855 {
2857 }
2858 
2859 ////////////////////////////////////////////////////////////////////////////////
2860 /// Writes array of characters to buffer
2861 
2863 {
2866 }
2867 
2868 ////////////////////////////////////////////////////////////////////////////////
2869 /// Writes a TString
2870 
2872 {
2873  if (GetIOVersion()<3) {
2875  } else {
2878  }
2879 }
2880 
2881 ////////////////////////////////////////////////////////////////////////////////
2882 /// Writes a TString
2883 
2884 void TBufferXML::WriteStdString(const std::string *s)
2885 {
2886  if (GetIOVersion()<3) {
2888  } else {
2890  if (s) XmlWriteValue(s->c_str(), xmlio::String);
2891  else XmlWriteValue("", xmlio::String);
2892  }
2893 }
2894 
2895 ////////////////////////////////////////////////////////////////////////////////
2896 /// Write a char* string
2897 
2899 {
2901 }
2902 
2903 
2904 ////////////////////////////////////////////////////////////////////////////////
2905 /// Converts Char_t to string and add xml node to buffer
2906 
2908 {
2909  char buf[50];
2910  snprintf(buf, sizeof(buf), "%d",value);
2911  return XmlWriteValue(buf, xmlio::Char);
2912 }
2913 
2914 ////////////////////////////////////////////////////////////////////////////////
2915 /// Converts Short_t to string and add xml node to buffer
2916 
2918 {
2919  char buf[50];
2920  snprintf(buf, sizeof(buf), "%hd", value);
2921  return XmlWriteValue(buf, xmlio::Short);
2922 }
2923 
2924 ////////////////////////////////////////////////////////////////////////////////
2925 /// Converts Int_t to string and add xml node to buffer
2926 
2928 {
2929  char buf[50];
2930  snprintf(buf, sizeof(buf), "%d", value);
2931  return XmlWriteValue(buf, xmlio::Int);
2932 }
2933 
2934 ////////////////////////////////////////////////////////////////////////////////
2935 /// Converts Long_t to string and add xml node to buffer
2936 
2938 {
2939  char buf[50];
2940  snprintf(buf, sizeof(buf), "%ld", value);
2941  return XmlWriteValue(buf, xmlio::Long);
2942 }
2943 
2944 ////////////////////////////////////////////////////////////////////////////////
2945 /// Converts Long64_t to string and add xml node to buffer
2946 
2948 {
2949  char buf[50];
2950  snprintf(buf, sizeof(buf), FLong64, value);
2951  return XmlWriteValue(buf, xmlio::Long64);
2952 }
2953 
2954 ////////////////////////////////////////////////////////////////////////////////
2955 /// Converts Float_t to string and add xml node to buffer
2956 
2958 {
2959  char buf[200];
2960  snprintf(buf, sizeof(buf), fgFloatFmt.c_str(), value);
2961  return XmlWriteValue(buf, xmlio::Float);
2962 }
2963 
2964 ////////////////////////////////////////////////////////////////////////////////
2965 /// Converts Double_t to string and add xml node to buffer
2966 
2968 {
2969  char buf[1000];
2970  snprintf(buf, sizeof(buf), fgFloatFmt.c_str(), value);
2971  return XmlWriteValue(buf, xmlio::Double);
2972 }
2973 
2974 ////////////////////////////////////////////////////////////////////////////////
2975 /// Converts Bool_t to string and add xml node to buffer
2976 
2978 {
2980 }
2981 
2982 ////////////////////////////////////////////////////////////////////////////////
2983 /// Converts UChar_t to string and add xml node to buffer
2984 
2986 {
2987  char buf[50];
2988  snprintf(buf, sizeof(buf), "%u", value);
2989  return XmlWriteValue(buf, xmlio::UChar);
2990 }
2991 
2992 ////////////////////////////////////////////////////////////////////////////////
2993 /// Converts UShort_t to string and add xml node to buffer
2994 
2996 {
2997  char buf[50];
2998  snprintf(buf, sizeof(buf), "%hu", value);
2999  return XmlWriteValue(buf, xmlio::UShort);
3000 }
3001 
3002 ////////////////////////////////////////////////////////////////////////////////
3003 /// Converts UInt_t to string and add xml node to buffer
3004 
3006 {
3007  char buf[50];
3008  snprintf(buf, sizeof(buf), "%u", value);
3009  return XmlWriteValue(buf, xmlio::UInt);
3010 }
3011 
3012 ////////////////////////////////////////////////////////////////////////////////
3013 /// Converts ULong_t to string and add xml node to buffer
3014 
3016 {
3017  char buf[50];
3018  snprintf(buf, sizeof(buf), "%lu", value);
3019  return XmlWriteValue(buf, xmlio::ULong);
3020 }
3021 
3022 ////////////////////////////////////////////////////////////////////////////////
3023 /// Converts ULong64_t to string and add xml node to buffer
3024 
3026 {
3027  char buf[50];
3028  snprintf(buf, sizeof(buf), FULong64, value);
3029  return XmlWriteValue(buf, xmlio::ULong64);
3030 }
3031 
3032 ////////////////////////////////////////////////////////////////////////////////
3033 /// Create xml node with specified name and adds it to stack node
3034 
3035 XMLNodePointer_t TBufferXML::XmlWriteValue(const char* value, const char* name)
3036 {
3037  XMLNodePointer_t node = 0;
3038 
3039  if (fCanUseCompact)
3040  node = StackNode();
3041  else
3042  node = CreateItemNode(name);
3043 
3044  fXML->NewAttr(node, 0, xmlio::v, value);
3045 
3047 
3048  return node;
3049 }
3050 
3051 ////////////////////////////////////////////////////////////////////////////////
3052 /// Reads string from current xml node and convert it to Char_t value
3053 
3055 {
3056  const char* res = XmlReadValue(xmlio::Char);
3057  if (res) {
3058  int n;
3059  sscanf(res,"%d", &n);
3060  value = n;
3061  } else
3062  value = 0;
3063 }
3064 
3065 ////////////////////////////////////////////////////////////////////////////////
3066 /// Reads string from current xml node and convert it to Short_t value
3067 
3069 {
3070  const char* res = XmlReadValue(xmlio::Short);
3071  if (res)
3072  sscanf(res,"%hd", &value);
3073  else
3074  value = 0;
3075 }
3076 
3077 ////////////////////////////////////////////////////////////////////////////////
3078 /// Reads string from current xml node and convert it to Int_t value
3079 
3081 {
3082  const char* res = XmlReadValue(xmlio::Int);
3083  if (res)
3084  sscanf(res,"%d", &value);
3085  else
3086  value = 0;
3087 }
3088 
3089 ////////////////////////////////////////////////////////////////////////////////
3090 /// Reads string from current xml node and convert it to Long_t value
3091 
3093 {
3094  const char* res = XmlReadValue(xmlio::Long);
3095  if (res)
3096  sscanf(res,"%ld", &value);
3097  else
3098  value = 0;
3099 }
3100 
3101 ////////////////////////////////////////////////////////////////////////////////
3102 /// Reads string from current xml node and convert it to Long64_t value
3103 
3105 {
3106  const char* res = XmlReadValue(xmlio::Long64);
3107  if (res)
3108  sscanf(res, FLong64, &value);
3109  else
3110  value = 0;
3111 }
3112 
3113 ////////////////////////////////////////////////////////////////////////////////
3114 /// Reads string from current xml node and convert it to Float_t value
3115 
3117 {
3118  const char* res = XmlReadValue(xmlio::Float);
3119  if (res)
3120  sscanf(res, "%f", &value);
3121  else
3122  value = 0.;
3123 }
3124 
3125 ////////////////////////////////////////////////////////////////////////////////
3126 /// Reads string from current xml node and convert it to Double_t value
3127 
3129 {
3130  const char* res = XmlReadValue(xmlio::Double);
3131  if (res)
3132  sscanf(res, "%lf", &value);
3133  else
3134  value = 0.;
3135 }
3136 
3137 ////////////////////////////////////////////////////////////////////////////////
3138 /// Reads string from current xml node and convert it to Bool_t value
3139 
3141 {
3142  const char* res = XmlReadValue(xmlio::Bool);
3143  if (res)
3144  value = (strcmp(res, xmlio::True)==0);
3145  else
3146  value = kFALSE;
3147 }
3148 
3149 ////////////////////////////////////////////////////////////////////////////////
3150 /// Reads string from current xml node and convert it to UChar_t value
3151 
3153 {
3154  const char* res = XmlReadValue(xmlio::UChar);
3155  if (res) {
3156  unsigned int n;
3157  sscanf(res,"%ud", &n);
3158  value = n;
3159  } else
3160  value = 0;
3161 }
3162 
3163 ////////////////////////////////////////////////////////////////////////////////
3164 /// Reads string from current xml node and convert it to UShort_t value
3165 
3167 {
3168  const char* res = XmlReadValue(xmlio::UShort);
3169  if (res)
3170  sscanf(res,"%hud", &value);
3171  else
3172  value = 0;
3173 }
3174 
3175 ////////////////////////////////////////////////////////////////////////////////
3176 /// Reads string from current xml node and convert it to UInt_t value
3177 
3179 {
3180  const char* res = XmlReadValue(xmlio::UInt);
3181  if (res)
3182  sscanf(res,"%u", &value);
3183  else
3184  value = 0;
3185 }
3186 
3187 ////////////////////////////////////////////////////////////////////////////////
3188 /// Reads string from current xml node and convert it to ULong_t value
3189 
3191 {
3192  const char* res = XmlReadValue(xmlio::ULong);
3193  if (res)
3194  sscanf(res,"%lu", &value);
3195  else
3196  value = 0;
3197 }
3198 
3199 ////////////////////////////////////////////////////////////////////////////////
3200 /// Reads string from current xml node and convert it to ULong64_t value
3201 
3203 {
3204  const char* res = XmlReadValue(xmlio::ULong64);
3205  if (res)
3206  sscanf(res, FULong64, &value);
3207  else
3208  value = 0;
3209 }
3210 
3211 ////////////////////////////////////////////////////////////////////////////////
3212 /// read string value from current stack node
3213 
3214 const char* TBufferXML::XmlReadValue(const char* name)
3215 {
3216  if (fErrorFlag>0) return 0;
3217 
3218  Bool_t trysimple = fCanUseCompact;
3220 
3221  if (trysimple) {
3222  if (fXML->HasAttr(Stack(1)->fNode,xmlio::v))
3223  fValueBuf = fXML->GetAttr(Stack(1)->fNode, xmlio::v);
3224  else
3225  trysimple = kFALSE;
3226  }
3227 
3228  if (!trysimple) {
3229  if (!VerifyItemNode(name, "XmlReadValue")) return 0;
3231  }
3232 
3233  if (gDebug>4)
3234  Info("XmlReadValue"," Name = %s value = %s", name, fValueBuf.Data());
3235 
3236  if (!trysimple)
3237  ShiftStack("readvalue");
3238 
3239  return fValueBuf.Data();
3240 }
3241 
3242 void TBufferXML::SetFloatFormat(const char* fmt)
3243 {
3244  // Set printf format for float/double members, default "%e"
3245  // This method is not thread-safe as it changes a global state.
3246 
3247  if (!fmt) fgFloatFmt = "%e";
3248  fgFloatFmt = fmt;
3249 }
3250 
3252 {
3253  // return current printf format for float/double members, default "%e"
3254 
3255  return fgFloatFmt.c_str();
3256 }
3257 
3258 ////////////////////////////////////////////////////////////////////////////////
3259 /// Read one collection of objects from the buffer using the StreamerInfoLoopAction.
3260 /// The collection needs to be a split TClonesArray or a split vector of pointers.
3261 
3263 {
3264  TVirtualStreamerInfo *info = sequence.fStreamerInfo;
3265  IncrementLevel(info);
3266 
3267  if (gDebug) {
3268  //loop on all active members
3269  TStreamerInfoActions::ActionContainer_t::const_iterator end = sequence.fActions.end();
3270  for(TStreamerInfoActions::ActionContainer_t::const_iterator iter = sequence.fActions.begin();
3271  iter != end;
3272  ++iter) {
3273  // Idea: Try to remove this function call as it is really needed only for XML streaming.
3274  SetStreamerElementNumber((*iter).fConfiguration->fCompInfo->fElem,(*iter).fConfiguration->fCompInfo->fType);
3275  (*iter).PrintDebug(*this,obj);
3276  (*iter)(*this,obj);
3277  }
3278 
3279  } else {
3280  //loop on all active members
3281  TStreamerInfoActions::ActionContainer_t::const_iterator end = sequence.fActions.end();
3282  for(TStreamerInfoActions::ActionContainer_t::const_iterator iter = sequence.fActions.begin();
3283  iter != end;
3284  ++iter) {
3285  // Idea: Try to remove this function call as it is really needed only for XML streaming.
3286  SetStreamerElementNumber((*iter).fConfiguration->fCompInfo->fElem,(*iter).fConfiguration->fCompInfo->fType);
3287  (*iter)(*this,obj);
3288  }
3289  }
3290 
3291  DecrementLevel(info);
3292  return 0;
3293 }
3294 
3295 ////////////////////////////////////////////////////////////////////////////////
3296 /// Read one collection of objects from the buffer using the StreamerInfoLoopAction.
3297 /// The collection needs to be a split TClonesArray or a split vector of pointers.
3298 
3299 Int_t TBufferXML::ApplySequenceVecPtr(const TStreamerInfoActions::TActionSequence &sequence, void *start_collection, void *end_collection)
3300 {
3301  TVirtualStreamerInfo *info = sequence.fStreamerInfo;
3302  IncrementLevel(info);
3303 
3304  if (gDebug) {
3305  //loop on all active members
3306  TStreamerInfoActions::ActionContainer_t::const_iterator end = sequence.fActions.end();
3307  for(TStreamerInfoActions::ActionContainer_t::const_iterator iter = sequence.fActions.begin();
3308  iter != end;
3309  ++iter) {
3310  // Idea: Try to remove this function call as it is really needed only for XML streaming.
3311  SetStreamerElementNumber((*iter).fConfiguration->fCompInfo->fElem,(*iter).fConfiguration->fCompInfo->fType);
3312  (*iter).PrintDebug(*this,*(char**)start_collection); // Warning: This limits us to TClonesArray and vector of pointers.
3313  (*iter)(*this,start_collection,end_collection);
3314  }
3315 
3316  } else {
3317  //loop on all active members
3318  TStreamerInfoActions::ActionContainer_t::const_iterator end = sequence.fActions.end();
3319  for(TStreamerInfoActions::ActionContainer_t::const_iterator iter = sequence.fActions.begin();
3320  iter != end;
3321  ++iter) {
3322  // Idea: Try to remove this function call as it is really needed only for XML streaming.
3323  SetStreamerElementNumber((*iter).fConfiguration->fCompInfo->fElem,(*iter).fConfiguration->fCompInfo->fType);
3324  (*iter)(*this,start_collection,end_collection);
3325  }
3326  }
3327 
3328  DecrementLevel(info);
3329  return 0;
3330 }
3331 
3332 ////////////////////////////////////////////////////////////////////////////////
3333 /// Read one collection of objects from the buffer using the StreamerInfoLoopAction.
3334 
3335 Int_t TBufferXML::ApplySequence(const TStreamerInfoActions::TActionSequence &sequence, void *start_collection, void *end_collection)
3336 {
3337  TVirtualStreamerInfo *info = sequence.fStreamerInfo;
3338  IncrementLevel(info);
3339 
3341  if (gDebug) {
3342 
3343  // Get the address of the first item for the PrintDebug.
3344  // (Performance is not essential here since we are going to print to
3345  // the screen anyway).
3346  void *arr0 = loopconfig->GetFirstAddress(start_collection,end_collection);
3347  // loop on all active members
3348  TStreamerInfoActions::ActionContainer_t::const_iterator end = sequence.fActions.end();
3349  for(TStreamerInfoActions::ActionContainer_t::const_iterator iter = sequence.fActions.begin();
3350  iter != end;
3351  ++iter) {
3352  // Idea: Try to remove this function call as it is really needed only for XML streaming.
3353  SetStreamerElementNumber((*iter).fConfiguration->fCompInfo->fElem,(*iter).fConfiguration->fCompInfo->fType);
3354  (*iter).PrintDebug(*this,arr0);
3355  (*iter)(*this,start_collection,end_collection,loopconfig);
3356  }
3357 
3358  } else {
3359  //loop on all active members
3360  TStreamerInfoActions::ActionContainer_t::const_iterator end = sequence.fActions.end();
3361  for(TStreamerInfoActions::ActionContainer_t::const_iterator iter = sequence.fActions.begin();
3362  iter != end;
3363  ++iter) {
3364  // Idea: Try to remove this function call as it is really needed only for XML streaming.
3365  SetStreamerElementNumber((*iter).fConfiguration->fCompInfo->fElem,(*iter).fConfiguration->fCompInfo->fType);
3366  (*iter)(*this,start_collection,end_collection,loopconfig);
3367  }
3368  }
3369 
3370  DecrementLevel(info);
3371  return 0;
3372 }
Int_t fCompressLevel
! Compression level and algorithm
Definition: TBufferXML.h:337
const char * ULong
Definition: TXMLSetup.cxx:92
Describe Streamer information for one class version.
Definition: TStreamerInfo.h:43
XMLNodePointer_t CreateItemNode(const char *name)
Create item node of specified name.
Definition: TBufferXML.cxx:738
virtual Int_t ApplySequenceVecPtr(const TStreamerInfoActions::TActionSequence &sequence, void *start_collection, void *end_collection)
Read one collection of objects from the buffer using the StreamerInfoLoopAction.
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
const char * Object
Definition: TXMLSetup.cxx:62
#define TBufferXML_operatorin(vname)
virtual ~TBufferXML()
Destroy xml buffer.
Definition: TBufferXML.cxx:146
#define TBufferXML_operatorout(vname)
void XmlReadBlock(XMLNodePointer_t node)
Read binary block of data from xml.
Definition: TBufferXML.cxx:505
Bool_t IsReading() const
Definition: TBuffer.h:81
const char * Long64
Definition: TXMLSetup.cxx:86
const char * XmlBlock
Definition: TXMLSetup.cxx:60
TXMLStackObj * PushStack(XMLNodePointer_t current, Bool_t simple=kFALSE)
Add new level to xml stack.
Definition: TBufferXML.cxx:343
const char * UInt
Definition: TXMLSetup.cxx:91
virtual Int_t ReadStaticArray(Bool_t *b)
Read array of Bool_t from buffer.
void Add(ULong64_t hash, Long64_t key, Long64_t value)
Add an (key,value) pair to the table. The key should be unique.
Definition: TExMap.cxx:87
virtual void WriteObject(const TObject *obj)
Convert object into xml structures.
Definition: TBufferXML.cxx:300
virtual Int_t ReadStaticArrayFloat16(Float_t *f, TStreamerElement *ele=0)
Read array of Float16_t from buffer.
Bool_t IsWriting() const
Definition: TBuffer.h:82
An array of TObjects.
Definition: TObjArray.h:37
Int_t GetCompressionSettings() const
Definition: TBufferXML.h:358
#define FULong64
Definition: TBufferXML.cxx:53
virtual void ReadFloat16(Float_t *f, TStreamerElement *ele=0)
Read a Float16_t from the buffer.
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:847
const char * Ref
Definition: TXMLSetup.cxx:53
The concrete implementation of TBuffer for writing/reading to/from a ROOT file or socket...
Definition: TBufferFile.h:47
virtual void WriteLong64(Long64_t l)
Writes Long64_t value to buffer.
const char * ObjClass
Definition: TXMLSetup.cxx:63
void SetCompressionAlgorithm(Int_t algorithm=0)
See comments for function SetCompressionSettings.
Definition: TBufferXML.cxx:404
virtual void ReadFastArrayDouble32(Double_t *d, Int_t n, TStreamerElement *ele=0)
Read array of Double32_t from buffer.
long long Long64_t
Definition: RtypesCore.h:69
void WorkWithElement(TStreamerElement *elem, Int_t comp_type)
This function is a part of SetStreamerElementNumber method.
Bool_t ProcessPointer(const void *ptr, XMLNodePointer_t node)
Add "ptr" attribute to node, if ptr is null or if ptr is pointer on object, which is already saved in...
Definition: TBufferXML.cxx:575
void XmlWriteBlock(XMLNodePointer_t node)
Write binary data block from buffer to xml.
Definition: TBufferXML.cxx:447
void FreeAttr(XMLNodePointer_t xmlnode, const char *name)
remove attribute from xmlnode
Definition: TXMLEngine.cxx:526
virtual void ReadWithFactor(Float_t *ptr, Double_t factor, Double_t minvalue)
Read a Double32_t from the buffer when the factor and minimun value have been specified see comments ...
const char * Double
Definition: TXMLSetup.cxx:88
virtual void ReadLong(Long_t &l)
Reads Long_t value from buffer.
short Version_t
Definition: RtypesCore.h:61
virtual void WriteFastArray(const Bool_t *b, Int_t n)
Write array of Bool_t to buffer.
virtual Int_t ReadArray(Bool_t *&b)
Read array of Bool_t from buffer.
TString fValueBuf
Definition: TBufferXML.h:330
virtual TClass * GetClass() const =0
TLoopConfiguration * fLoopConfig
If this is a bundle of memberwise streaming action, this configures the looping.
static void * ConvertFromXMLAny(const char *str, TClass **cl=0, Bool_t GenericLayout=kFALSE, Bool_t UseNamespaces=kFALSE)
Read object of any class from XML, produced by ConvertToXML() method.
Definition: TBufferXML.cxx:237
float Float_t
Definition: RtypesCore.h:53
virtual void WriteLong(Long_t l)
Writes Long_t value to buffer.
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:329
void SetCompressionSettings(Int_t settings=1)
Used to specify the compression level and algorithm.
Definition: TBufferXML.cxx:438
void SetCompressionLevel(Int_t level=1)
See comments for function SetCompressionSettings.
Definition: TBufferXML.cxx:419
const char * Size
Definition: TXMLSetup.cxx:56
TObject * GetParent() const
Return pointer to parent of this buffer.
Definition: TBuffer.cxx:231
EXMLLayout GetXmlLayout() const
Definition: TXMLSetup.h:97
const char * v
Definition: TXMLSetup.cxx:74
#define TBufferXML_WriteArray(vname)
virtual void ReadTString(TString &s)
Read TString from TBuffer.
unsigned short UShort_t
Definition: RtypesCore.h:36
virtual void ClassMember(const char *name, const char *typeName=0, Int_t arrsize1=-1, Int_t arrsize2=-1)
Method indicates name and typename of class member, which should be now streamed in custom streamer...
virtual TClass * GetClassPointer() const
Returns a pointer to the TClass of this element.
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TH1 * h
Definition: legend2.C:5
XMLNodePointer_t XmlWriteObject(const void *obj, const TClass *objClass)
Write object to buffer If object was written before, only pointer will be stored Return pointer to to...
Definition: TBufferXML.cxx:817
virtual void WriteULong64(ULong64_t l)
Writes ULong64_t value to buffer.
const char * ULong64
Definition: TXMLSetup.cxx:93
virtual void ReadFastArray(Bool_t *b, Int_t n)
Read array of Bool_t from buffer.
virtual void SkipVersion(const TClass *cl=0)
Skip class version from I/O buffer.
Version_t fVersionBuf
Definition: TBufferXML.h:325
#define gROOT
Definition: TROOT.h:375
virtual TObject * Remove(TObject *obj)
Remove object from array.
Definition: TObjArray.cxx:653
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)
Suppressed function of TBuffer.
virtual void WriteCharStar(char *s)
Write a char* string.
void RegisterPointer(const void *ptr, XMLNodePointer_t node)
Register pair of object pointer and node, where this object is saved, in object map.
Definition: TBufferXML.cxx:614
Basic string class.
Definition: TString.h:129
virtual Int_t ReadStaticArrayDouble32(Double_t *d, TStreamerElement *ele=0)
Read array of Double32_t from buffer.
virtual void WriteChar(Char_t c)
Writes Char_t value to buffer.
virtual void WriteCharStar(char *s)
Write char* into TBuffer.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual TClass * ReadClass(const TClass *cl=0, UInt_t *objTag=0)
Function to read class from buffer, used in old-style streamers.
void SetParent(TObject *parent)
Set parent owning this buffer.
Definition: TBuffer.cxx:239
XMLNodePointer_t XmlWriteAny(const void *obj, const TClass *cl)
Convert object of any class to xml structures Return pointer on top xml element.
Definition: TBufferXML.cxx:260
virtual void ReadFastArrayWithNbits(Float_t *ptr, Int_t n, Int_t nbits)
Read array of Float16_t from buffer.
const char * Item
Definition: TXMLSetup.cxx:66
virtual void SetArrayDim(Int_t dim)
Set number of array dimensions.
void CreateElemNode(const TStreamerElement *elem)
Create xml node correspondent to TStreamerElement object.
Definition: TBufferXML.cxx:766
TObject * At(Int_t idx) const
Definition: TObjArray.h:165
virtual void WriteTString(const TString &s)
Writes a TString.
const char * Class
Definition: TXMLSetup.cxx:64
virtual void WriteDouble32(Double_t *d, TStreamerElement *ele=0)
Write a Double32_t to the buffer.
virtual void SetMaxIndex(Int_t dim, Int_t max)
set maximum index for array with dimension dim
Float_t delta
const char * False
Definition: TXMLSetup.cxx:77
void UnlinkFreeNode(XMLNodePointer_t xmlnode)
combined operation. Unlink node and free used memory
Definition: TXMLEngine.cxx:921
XMLNodePointer_t XmlWriteValue(const char *value, const char *name)
Create xml node with specified name and adds it to stack node.
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:687
XMLAttrPointer_t NewIntAttr(XMLNodePointer_t xmlnode, const char *name, Int_t value)
create node attribute with integer value
Definition: TXMLEngine.cxx:514
const char * Name
Definition: TXMLSetup.cxx:67
#define TBufferXML_WriteFastArray(vname)
virtual void WriteULong(ULong_t l)
Writes ULong_t value to buffer.
const char * GetNodeContent(XMLNodePointer_t xmlnode)
get contents (if any) of xml node
Definition: TXMLEngine.cxx:938
Int_t Length() const
Definition: TBuffer.h:94
virtual void WriteStdString(const std::string *s)
Writes a TString.
TXMLEngine * fXML
Definition: TBufferXML.h:321
const char * Float
Definition: TXMLSetup.cxx:87
const char * String
Definition: TXMLSetup.cxx:94
virtual void ReadCharP(Char_t *c)
Reads array of characters from buffer.
TXMLFile * XmlFile()
Returns pointer to TXMLFile object.
Definition: TBufferXML.cxx:157
Int_t GetBaseClassOffset(const TClass *toBase, void *address=0, bool isDerivedObject=true)
Definition: TClass.cxx:2709
virtual void WriteInt(Int_t i)
Writes Int_t value to buffer.
Bool_t VerifyAttr(XMLNodePointer_t node, const char *name, const char *value, const char *errinfo=0)
Checks, that attribute of specified name exists and has specified value.
Definition: TBufferXML.cxx:713
TObject * Last() const
Return the object in the last filled slot. Returns 0 if no entries.
Definition: TObjArray.cxx:479
virtual void WriteShort(Short_t s)
Writes Short_t value to buffer.
TExMap * fObjMap
Definition: TBufferXML.h:327
virtual void ReadFastArrayFloat16(Float_t *f, Int_t n, TStreamerElement *ele=0)
Read array of Float16_t from buffer.
virtual void * GetFirstAddress(void *start, const void *end) const =0
virtual void DecrementLevel(TVirtualStreamerInfo *)
Function is called from TStreamerInfo WriteBuffer and Readbuffer functions and decrease level in xml ...
Definition: TBufferXML.cxx:973
Bool_t VerifyNode(XMLNodePointer_t node, const char *name, const char *errinfo=0)
Check, if node has specified name.
Definition: TBufferXML.cxx:686
virtual Int_t ApplySequence(const TStreamerInfoActions::TActionSequence &sequence, void *object)
Read one collection of objects from the buffer using the StreamerInfoLoopAction.
virtual void WriteUChar(UChar_t c)
Writes UChar_t value to buffer.
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)
Read version value from buffer.
void Class()
Definition: Class.C:29
virtual void ReadCharStar(char *&s)
Read char* from TBuffer.
virtual void ReadCharStar(char *&s)
Read a char* string.
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:396
void BeforeIOoperation()
Function is called before any IO operation of TBuffer Now is used to store version value if no proper...
virtual void ReadShort(Short_t &s)
Reads Short_t value from buffer.
virtual void WriteStdString(const std::string *s)
Write std::string to TBuffer.
virtual void StreamObject(void *obj, const std::type_info &typeinfo, const TClass *onFileClass=0)
steram object to/from buffer
char * Buffer() const
Definition: TBuffer.h:91
Int_t AtoI(const char *sbuf, Int_t def=0, const char *errinfo=0)
converts string to integer.
Definition: TXMLSetup.cxx:287
virtual void WriteUShort(UShort_t s)
Writes UShort_t value to buffer.
virtual void WriteUInt(UInt_t i)
Writes UInt_t value to buffer.
XMLNsPointer_t NewNS(XMLNodePointer_t xmlnode, const char *reference, const char *name=0)
create namespace attribute for xmlnode.
Definition: TXMLEngine.cxx:641
int R__unzip_header(Int_t *nin, UChar_t *bufin, Int_t *lout)
void SetIOVersion(Int_t v)
Definition: TBufferXML.h:51
virtual void ReadStdString(std::string *s)
Reads a std::string.
virtual void WriteFastArrayDouble32(const Double_t *d, Int_t n, TStreamerElement *ele=0)
Write array of Double32_t to buffer.
virtual void WriteCharP(const Char_t *c)
Writes array of characters to buffer.
const char * UChar
Definition: TXMLSetup.cxx:89
const char * Char
Definition: TXMLSetup.cxx:82
virtual void ReadStdString(std::string *s)
Read std::string from TBuffer.
Base class of the Configurations for the member wise looping routines.
TObject()
TObject constructor.
Definition: TObject.h:213
virtual void ReadULong64(ULong64_t &l)
Reads ULong64_t value from buffer.
virtual void ReadBool(Bool_t &b)
Reads Bool_t value from buffer.
const char * GetNodeName(XMLNodePointer_t xmlnode)
returns name of xmlnode
Definition: TXMLEngine.cxx:930
XMLNodePointer_t StackNode()
Return pointer on current xml node.
Definition: TBufferXML.cxx:383
void Expand(Int_t newsize, Bool_t copy=kTRUE)
Expand (or shrink) the I/O buffer to newsize bytes.
Definition: TBuffer.cxx:201
virtual void ReadUChar(UChar_t &c)
Reads UChar_t value from buffer.
virtual void WriteTString(const TString &s)
Write TString to TBuffer.
Int_t GetType() const
Definition: TDataType.h:68
virtual void ReadDouble32(Double_t *d, TStreamerElement *ele=0)
Read a Double32_t from the buffer.
TStreamerInfo * fInfo
Pointer to TStreamerInfo object writing/reading the buffer.
Definition: TBufferFile.h:58
Int_t GetLast() const
Return index of last object in array.
Definition: TObjArray.cxx:528
TObjArray fStack
Definition: TBufferXML.h:323
Int_t fIOVersion
! Indicates format of ROOT xml file
Definition: TBufferXML.h:338
virtual void ReadUInt(UInt_t &i)
Reads UInt_t value from buffer.
const char * XmlGetElementName(const TStreamerElement *el)
return converted name for TStreamerElement
Definition: TXMLSetup.cxx:247
TClass * fExpectedBaseClass
! Pointer to class, which should be stored as parent of current
Definition: TBufferXML.h:336
Long64_t GetValue(ULong64_t hash, Long64_t key)
Return the value belonging to specified key and hash value.
Definition: TExMap.cxx:173
void ShiftStack(const char *info=0)
Shift stack node to next.
Definition: TBufferXML.cxx:392
Bool_t VerifyStackAttr(const char *name, const char *value, const char *errinfo=0)
Checks stack attribute.
Definition: TBufferXML.cxx:730
virtual void * ReadObjectAny(const TClass *clCast)
Read object from buffer. Only used from TBuffer.
void ShiftToNext(XMLNodePointer_t &xmlnode, Bool_t realnode=kTRUE)
shifts specified node to next if realnode==kTRUE, any special nodes in between will be skipped ...
TClass * GetClass() const
void SetXML(TXMLEngine *xml)
Definition: TBufferXML.h:248
#define TBufferXML_ReadArray(tname, vname)
virtual void IncrementLevel(TVirtualStreamerInfo *)
Function is called from TStreamerInfo WriteBuffer and Readbuffer functions and indent new level in xm...
Definition: TBufferXML.cxx:904
Basic data type descriptor (datatype information is obtained from CINT).
Definition: TDataType.h:44
virtual Int_t ReadArrayFloat16(Float_t *&f, TStreamerElement *ele=0)
Read array of Float16_t from buffer.
TBufferXML()
Default constructor.
Definition: TBufferXML.cxx:64
virtual void WriteFloat16(Float_t *f, TStreamerElement *ele=0)
Write a Float16_t to the buffer.
void Destructor(void *obj, Bool_t dtorOnly=kFALSE)
Explicitly call destructor for object.
Definition: TClass.cxx:5063
static std::string fgFloatFmt
! Printf argument for floats and doubles, either "%f" or "%e" or "%10f" and so on ...
Definition: TBufferXML.h:340
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
Ssiz_t Length() const
Definition: TString.h:388
#define TBufferXML_ReadStaticArray(vname)
virtual void ClassEnd(const TClass *)
Should be called at the end of custom streamer See TBufferXML::ClassBegin for more details...
short Short_t
Definition: RtypesCore.h:35
TLine * l
Definition: textangle.C:4
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:71
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:168
Int_t GetNextRefCounter()
Definition: TXMLSetup.h:111
Bool_t fCanUseCompact
! Flag indicate that basic type (like Int_t) can be placed in the same tag
Definition: TBufferXML.h:334
virtual void WriteArrayFloat16(const Float_t *f, Int_t n, TStreamerElement *ele=0)
Write array of Float16_t to buffer.
Bool_t VerifyStackNode(const char *name, const char *errinfo=0)
Check, if stack node has specified name.
Definition: TBufferXML.cxx:704
void * XMLNodePointer_t
Definition: TXMLEngine.h:17
Bool_t VerifyElemNode(const TStreamerElement *elem)
Checks, if stack node correspond to TStreamerElement object.
Definition: TBufferXML.cxx:794
void ExtractReference(XMLNodePointer_t node, const void *ptr, const TClass *cl)
Analyse, if node has "ref" attribute and register it to object map.
Definition: TBufferXML.cxx:660
static void SetFloatFormat(const char *fmt="%e")
void PerformPreProcessing(const TStreamerElement *elem, XMLNodePointer_t elemnode)
Function is unpack TObject and TString structures to be able read them from custom streamers of this ...
void Streamer(void *obj, TBuffer &b, const TClass *onfile_class=0) const
Definition: TClass.h:530
virtual void WriteObjectClass(const void *actualObjStart, const TClass *actualClass)
Write object to buffer. Only used from TBuffer.
Bool_t IsUseNamespaces() const
Definition: TXMLSetup.h:100
void SkipEmpty(XMLNodePointer_t &xmlnode)
Skip all current empty nodes and locate on first "true" node.
const char * XmlConvertClassName(const char *name)
convert class name to exclude any special symbols like &#39;:&#39;, &#39;<&#39; &#39;>&#39; &#39;,&#39; and spaces ...
Definition: TXMLSetup.cxx:220
virtual void WriteFloat(Float_t f)
Writes Float_t value to buffer.
void CheckVersionBuf()
Checks buffer, filled by WriteVersion if next data is arriving, version should be stored in buffer...
const char * Ptr
Definition: TXMLSetup.cxx:52
const Bool_t kFALSE
Definition: RtypesCore.h:92
Bool_t HasAttr(XMLNodePointer_t xmlnode, const char *name)
checks if node has attribute of specified name
Definition: TXMLEngine.cxx:446
static const char * GetFloatFormat()
#define FLong64
Definition: TBufferXML.cxx:52
long Long_t
Definition: RtypesCore.h:50
virtual void WriteArrayDouble32(const Double_t *d, Int_t n, TStreamerElement *ele=0)
Write array of Double32_t to buffer.
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:215
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)
Copies class version to buffer, but not writes it to xml Version will be written with next I/O operat...
TObjArray * fIdArray
Definition: TBufferXML.h:328
Version_t GetClassVersion() const
Definition: TClass.h:372
XMLAttrPointer_t NewAttr(XMLNodePointer_t xmlnode, XMLNsPointer_t, const char *name, const char *value)
creates new attribute for xmlnode, namespaces are not supported for attributes
Definition: TXMLEngine.cxx:488
virtual void ReadWithNbits(Float_t *ptr, Int_t nbits)
Read a Float16_t from the buffer when the number of bits is specified (explicitly or not) see comment...
virtual void ReadDouble(Double_t &d)
Reads Double_t value from buffer.
#define ClassImp(name)
Definition: Rtypes.h:336
double f(double x)
double Double_t
Definition: RtypesCore.h:55
virtual void WriteArray(const Bool_t *b, Int_t n)
Write array of Bool_t to buffer.
Bool_t ExtractPointer(XMLNodePointer_t node, void *&ptr, TClass *&cl)
Searches for "ptr" attribute and returns pointer to object and class, if "ptr" attribute reference to...
Definition: TBufferXML.cxx:630
const char * GetAttr(XMLNodePointer_t xmlnode, const char *name)
returns value of attribute for xmlnode
Definition: TXMLEngine.cxx:460
void SaveSingleNode(XMLNodePointer_t xmlnode, TString *res, Int_t layout=1)
convert single xml node (and its child node) to string if layout<=0, no any spaces or newlines will b...
virtual void ReadInt(Int_t &i)
Reads Int_t value from buffer.
virtual void ReadUShort(UShort_t &s)
Reads UShort_t value from buffer.
unsigned long long ULong64_t
Definition: RtypesCore.h:70
Int_t IndexOf(const TObject *obj) const
Definition: TObjArray.cxx:552
unsigned long ULong_t
Definition: RtypesCore.h:51
virtual void ReadULong(ULong_t &l)
Reads ULong_t value from buffer.
Bool_t fExpectedChain
! Flag to resolve situation when several elements of same basic type stored as FastArray ...
Definition: TBufferXML.h:335
void WorkWithClass(TStreamerInfo *info, const TClass *cl=0)
Prepares buffer to stream data of specified class.
Definition: TBufferXML.cxx:912
Int_t GetCompressionLevel() const
Definition: TBufferXML.h:352
const char * XmlClassNameSpaceRef(const TClass *cl)
produce string which used as reference in class namespace definition
Definition: TXMLSetup.cxx:234
virtual void ReadTString(TString &s)
Reads a TString.
Int_t BufferSize() const
Definition: TBuffer.h:92
virtual void WriteObject(const TObject *obj)
Write object to I/O buffer.
virtual void WriteDouble(Double_t d)
Writes Double_t value to buffer.
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:2885
void FreeNode(XMLNodePointer_t xmlnode)
release all memory, allocated fro this node and destroyes node itself
Definition: TXMLEngine.cxx:893
void R__unzip(Int_t *nin, UChar_t *bufin, Int_t *lout, char *bufout, Int_t *nout)
Bool_t VerifyItemNode(const char *name, const char *errinfo=0)
Checks, if stack node is item and has specified name.
Definition: TBufferXML.cxx:752
TXMLStackObj * Stack(Int_t depth=0)
Return xml stack object of specified depth.
Definition: TBufferXML.cxx:372
virtual void WriteBool(Bool_t b)
Writes Bool_t value to buffer.
Mother of all ROOT objects.
Definition: TObject.h:37
TObjArray * GetElements() const
void * XMLNsPointer_t
Definition: TXMLEngine.h:18
XMLNodePointer_t XmlWriteBasic(Char_t value)
Converts Char_t to string and add xml node to buffer.
char Char_t
Definition: RtypesCore.h:29
const char * Null
Definition: TXMLSetup.cxx:54
Bool_t IsTObject() const
Return kTRUE is the class inherits from TObject.
Definition: TClass.cxx:5575
const char * CharStar
Definition: TXMLSetup.cxx:95
Class for serializing/deserializing object to/from xml.
Definition: TBufferXML.h:34
virtual void ReadLong64(Long64_t &l)
Reads Long64_t value from buffer.
const char * Member
Definition: TXMLSetup.cxx:65
XMLNodePointer_t GetChild(XMLNodePointer_t xmlnode, Bool_t realnode=kTRUE)
returns first child of xml node
Definition: TXMLEngine.cxx:993
void * XmlReadObject(void *obj, TClass **cl=0)
Read object from the buffer.
Definition: TBufferXML.cxx:845
virtual void SetUseNamespaces(Bool_t iUseNamespaces=kTRUE)
Definition: TXMLSetup.h:105
Definition: file.py:1
Int_t fErrorFlag
Definition: TBufferXML.h:332
const char * Int
Definition: TXMLSetup.cxx:84
Int_t GetIOVersion() const
Definition: TBufferXML.h:50
Int_t GetCompressionAlgorithm() const
Definition: TBufferXML.h:346
virtual void SetByteCount(UInt_t cntpos, Bool_t packInVersion=kFALSE)
Suppressed function of TBuffer.
XMLNodePointer_t NewChild(XMLNodePointer_t parent, XMLNsPointer_t ns, const char *name, const char *content=0)
create new child element for parent node
Definition: TXMLEngine.cxx:614
virtual void WriteFastArray(const Bool_t *b, Int_t n)
Write array of n bools into the I/O buffer.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
virtual void ReadFastArrayWithFactor(Float_t *ptr, Int_t n, Double_t factor, Double_t minvalue)
Read array of Float16_t from buffer.
#define snprintf
Definition: civetweb.c:822
virtual void SkipObjectAny()
Skip any kind of object from buffer Actually skip only one node on current level of xml structure...
Int_t GetIntAttr(XMLNodePointer_t node, const char *name)
returns value of attribute as integer
Definition: TXMLEngine.cxx:475
Int_t fBufSize
Definition: TBuffer.h:47
virtual void WriteClass(const TClass *cl)
Function to write class into buffer, used in old-style streamers.
R__EXTERN Int_t gDebug
Definition: Rtypes.h:83
virtual void SetStreamerElementNumber(TStreamerElement *elem, Int_t comp_type)
Function is called from TStreamerInfo WriteBuffer and Readbuffer functions and add/verify next elemen...
void Add(TObject *obj)
Definition: TObjArray.h:73
void SetBaseVersion(Int_t v)
TClass * XmlDefineClass(const char *xmlClassName)
define class for the converted class name, where special symbols were replaced by &#39;_&#39; ...
Definition: TXMLSetup.cxx:270
const char * Short
Definition: TXMLSetup.cxx:83
const char * OnlyVersion
Definition: TXMLSetup.cxx:51
unsigned char UChar_t
Definition: RtypesCore.h:34
void PerformPostProcessing()
Function is converts TObject and TString structures to more compact representation.
virtual void WriteFastArrayFloat16(const Float_t *d, Int_t n, TStreamerElement *ele=0)
Write array of Float16_t to buffer.
virtual void ReadChar(Char_t &c)
Reads Char_t value from buffer.
virtual void ClassBegin(const TClass *, Version_t=-1)
Should be called at the beginning of custom class streamer.
XMLNodePointer_t ReadSingleNode(const char *src)
read single xml node from provided string
#define TBufferXML_ReadFastArray(vname)
virtual void Compress()
Remove empty slots from array.
Definition: TObjArray.cxx:309
Abstract Interface class describing Streamer information for one class.
const char * ClassVersion
Definition: TXMLSetup.cxx:49
const char * True
Definition: TXMLSetup.cxx:76
virtual Int_t ReadArrayDouble32(Double_t *&d, TStreamerElement *ele=0)
Read array of Double32_t from buffer.
const Bool_t kTRUE
Definition: RtypesCore.h:91
const char * UShort
Definition: TXMLSetup.cxx:90
virtual void SetXmlLayout(EXMLLayout layout)
Definition: TXMLSetup.h:102
Int_t GetType() const
const char * IdBase
Definition: TXMLSetup.cxx:55
const Int_t n
Definition: legend1.C:16
virtual void ReadFloat(Float_t &f)
Reads Float_t value from buffer.
TXMLStackObj * PopStack()
Remove one level from xml stack.
Definition: TBufferXML.cxx:358
void XmlReadBasic(Char_t &value)
Reads string from current xml node and convert it to Char_t value.
const char * Zip
Definition: TXMLSetup.cxx:61
This class stores a (key,value) pair using an external hash.
Definition: TExMap.h:33
const char * Long
Definition: TXMLSetup.cxx:85
virtual void ReadFastArray(Bool_t *b, Int_t n)
Read array of n bools from the I/O buffer.
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
const char * XmlReadValue(const char *name)
read string value from current stack node
const char * Bool
Definition: TXMLSetup.cxx:81
void * XmlReadAny(XMLNodePointer_t node, void *obj, TClass **cl)
Recreate object from xml structure.
Definition: TBufferXML.cxx:276
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition: TClass.cxx:4706
const char * Data() const
Definition: TString.h:347
TVirtualStreamerInfo * fStreamerInfo
StreamerInfo used to derive these actions.