Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
16Class for serializing/deserializing object to/from xml.
17
18The simple way to create XML representation is:
19~~~{.cpp}
20 TNamed *obj = new TNamed("name", "title");
21 TString xml = TBufferXML::ToXML(obj);
22~~~
23Produced xml can be decoded into new object:
24~~~{.cpp}
25 TNamed *obj2 = nullptr;
26 TBufferXML::FromXML(obj2, xml);
27~~~
28
29TBufferXML class uses streaming mechanism, provided by ROOT system,
30therefore most of ROOT and user classes can be stored to xml.
31There are limitations for complex objects like TTree, which can not be converted to xml.
32*/
33
34#include "TBufferXML.h"
35
36#include "Compression.h"
37#include "TXMLFile.h"
38#include "TROOT.h"
39#include "TError.h"
40#include "TClass.h"
41#include "TClassTable.h"
42#include "TDataType.h"
43#include "TExMap.h"
44#include "TStreamerInfo.h"
45#include "TStreamerElement.h"
46#include "TMemberStreamer.h"
47#include "TStreamer.h"
48#include "RZip.h"
49#include "snprintf.h"
50
51#include <memory>
52
54
55////////////////////////////////////////////////////////////////////////////////
56/// Creates buffer object to serialize/deserialize data to/from xml.
57/// Mode should be either TBuffer::kRead or TBuffer::kWrite.
58
60{
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Creates buffer object to serialize/deserialize data to/from xml.
65/// This constructor should be used, if data from buffer supposed to be stored in file.
66/// Mode should be either TBuffer::kRead or TBuffer::kWrite.
67
69 : TBufferText(mode, file), TXMLSetup(*file)
70{
71 // this is for the case when StreamerInfo reads elements from
72 // buffer as ReadFastArray. When it checks if size of buffer is
73 // too small and skip reading. Actually, more improved method should
74 // be used here.
75
76 if (XmlFile()) {
77 SetXML(XmlFile()->XML());
80 }
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Destroy xml buffer.
85
87{
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// Returns pointer to TXMLFile object.
92/// Access to file is necessary to produce unique identifier for object references.
93
95{
96 return dynamic_cast<TXMLFile *>(GetParent());
97}
98
99////////////////////////////////////////////////////////////////////////////////
100/// Converts object, inherited from TObject class, to XML string
101/// GenericLayout defines layout choice for XML file
102/// UseNamespaces allow XML namespaces.
103/// See TXMLSetup class for details
104
105TString TBufferXML::ConvertToXML(const TObject *obj, Bool_t GenericLayout, Bool_t UseNamespaces)
106{
107 TClass *clActual = nullptr;
108 void *ptr = (void *)obj;
109
110 if (obj) {
111 clActual = TObject::Class()->GetActualClass(obj);
112 if (!clActual)
113 clActual = TObject::Class();
114 else if (clActual != TObject::Class())
115 ptr = (void *)((Longptr_t)obj - clActual->GetBaseClassOffset(TObject::Class()));
116 }
117
118 return ConvertToXML(ptr, clActual, GenericLayout, UseNamespaces);
119}
120
121////////////////////////////////////////////////////////////////////////////////
122/// Converts any type of object to XML string.
123/// GenericLayout defines layout choice for XML file
124/// UseNamespaces allow XML namespaces.
125/// See TXMLSetup class for details
126
127TString TBufferXML::ConvertToXML(const void *obj, const TClass *cl, Bool_t GenericLayout, Bool_t UseNamespaces)
128{
129 TXMLEngine xml;
130
132 buf.SetXML(&xml);
133 buf.InitMap();
134
136 buf.SetUseNamespaces(UseNamespaces);
137
138 XMLNodePointer_t xmlnode = buf.XmlWriteAny(obj, cl);
139
140 TString res;
141
142 xml.SaveSingleNode(xmlnode, &res);
143
144 xml.FreeNode(xmlnode);
145
146 return res;
147}
148
149////////////////////////////////////////////////////////////////////////////////
150/// Read object from XML, produced by ConvertToXML() method.
151/// If object does not inherit from TObject class, return 0.
152/// GenericLayout and UseNamespaces should be the same as in ConvertToXML()
153
154TObject *TBufferXML::ConvertFromXML(const char *str, Bool_t GenericLayout, Bool_t UseNamespaces)
155{
156 TClass *cl = nullptr;
157 void *obj = ConvertFromXMLAny(str, &cl, GenericLayout, UseNamespaces);
158
159 if (!cl || !obj)
160 return nullptr;
161
162 Int_t delta = cl->GetBaseClassOffset(TObject::Class());
163
164 if (delta < 0) {
165 cl->Destructor(obj);
166 return nullptr;
167 }
168
169 return (TObject *)(((char *)obj) + delta);
170}
171
172////////////////////////////////////////////////////////////////////////////////
173/// Read object of any class from XML, produced by ConvertToXML() method.
174/// If cl!=0, return actual class of object.
175/// GenericLayout and UseNamespaces should be the same as in ConvertToXML()
176
177void *TBufferXML::ConvertFromXMLAny(const char *str, TClass **cl, Bool_t GenericLayout, Bool_t UseNamespaces)
178{
179 TXMLEngine xml;
181
182 buf.SetXML(&xml);
183 buf.InitMap();
184
186 buf.SetUseNamespaces(UseNamespaces);
187
188 XMLNodePointer_t xmlnode = xml.ReadSingleNode(str);
189
190 void *obj = buf.XmlReadAny(xmlnode, nullptr, cl);
191
192 xml.FreeNode(xmlnode);
193
194 return obj;
195}
196
197////////////////////////////////////////////////////////////////////////////////
198/// Convert from XML and check if object derived from specified class
199/// When possible, cast to given class
200
201void *TBufferXML::ConvertFromXMLChecked(const char *xml, const TClass *expectedClass, Bool_t GenericLayout,
202 Bool_t UseNamespaces)
203{
204 TClass *objClass = nullptr;
205 void *res = ConvertFromXMLAny(xml, &objClass, GenericLayout, UseNamespaces);
206
207 if (!res || !objClass)
208 return nullptr;
209
210 if (objClass == expectedClass)
211 return res;
212
213 Int_t offset = objClass->GetBaseClassOffset(expectedClass);
214 if (offset < 0) {
215 ::Error("TBufferXML::ConvertFromXMLChecked", "expected class %s is not base for read class %s",
216 expectedClass->GetName(), objClass->GetName());
217 objClass->Destructor(res);
218 return nullptr;
219 }
220
221 return (char *)res - offset;
222}
223
224////////////////////////////////////////////////////////////////////////////////
225/// Convert object of any class to xml structures
226/// Return pointer on top xml element
227
229{
230 fErrorFlag = 0;
231
232 if (!fXML)
233 return nullptr;
234
235 XMLNodePointer_t res = XmlWriteObject(obj, cl, kTRUE);
236
237 return res;
238}
239
240////////////////////////////////////////////////////////////////////////////////
241/// Recreate object from xml structure.
242/// Return pointer to read object.
243/// if (cl!=0) returns pointer to class of object
244
246{
247 if (!node)
248 return nullptr;
249
250 if (cl)
251 *cl = nullptr;
252
253 fErrorFlag = 0;
254
255 if (!fXML)
256 return nullptr;
257
258 PushStack(node, kTRUE);
259
260 void *res = XmlReadObject(obj, cl);
261
262 PopStack();
263
264 return res;
265}
266
267// TXMLStackObj is used to keep stack of object hierarchy,
268// stored in TBuffer. For example, data for parent class(es)
269// stored in subnodes, but initial object node will be kept.
270
272public:
274 {
275 }
276
278 {
279 if (fIsElemOwner)
280 delete fElem;
281 }
282
284
293};
294
295////////////////////////////////////////////////////////////////////////////////
296/// Add new level to xml stack.
297
299{
300 if (IsReading() && !simple) {
301 current = fXML->GetChild(current);
302 fXML->SkipEmpty(current);
303 }
304
305 fStack.emplace_back(std::make_unique<TXMLStackObj>(current));
306 return fStack.back().get();
307}
308
309////////////////////////////////////////////////////////////////////////////////
310/// Remove one level from xml stack.
311
313{
314 if (fStack.size() > 0)
315 fStack.pop_back();
316 return fStack.size() > 0 ? fStack.back().get() : nullptr;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// Return pointer on current xml node.
321
323{
324 TXMLStackObj *stack = Stack();
325 return stack ? stack->fNode : nullptr;
326}
327
328////////////////////////////////////////////////////////////////////////////////
329/// Shift stack node to next.
330
331void TBufferXML::ShiftStack(const char *errinfo)
332{
333 TXMLStackObj *stack = Stack();
334 if (stack) {
335 fXML->ShiftToNext(stack->fNode);
336 if (gDebug > 4)
337 Info("ShiftStack", "%s to node %s", errinfo, fXML->GetNodeName(stack->fNode));
338 }
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// See comments for function SetCompressionSettings.
343
345{
346 if (algorithm < 0 || algorithm >= ROOT::RCompressionSetting::EAlgorithm::kUndefined)
347 algorithm = 0;
348 if (fCompressLevel < 0) {
350 } else {
351 int level = fCompressLevel % 100;
352 fCompressLevel = 100 * algorithm + level;
353 }
354}
355
356////////////////////////////////////////////////////////////////////////////////
357/// See comments for function SetCompressionSettings.
358
360{
361 if (level < 0)
362 level = 0;
363 if (level > 99)
364 level = 99;
365 if (fCompressLevel < 0) {
366 // if the algorithm is not defined yet use 0 as a default
367 fCompressLevel = level;
368 } else {
369 int algorithm = fCompressLevel / 100;
371 algorithm = 0;
372 fCompressLevel = 100 * algorithm + level;
373 }
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// Used to specify the compression level and algorithm.
378///
379/// See TFile constructor for the details.
380
382{
383 fCompressLevel = settings;
384}
385
386////////////////////////////////////////////////////////////////////////////////
387/// Write binary data block from buffer to xml.
388/// This data can be produced only by direct call of TBuffer::WriteBuf() functions.
389
391{
392 if (!node || (Length() == 0))
393 return;
394
395 const char *src = Buffer();
396 int srcSize = Length();
397
398 char *fZipBuffer = 0;
399
400 Int_t compressionLevel = GetCompressionLevel();
403
404 if ((Length() > 512) && (compressionLevel > 0)) {
405 int zipBufferSize = Length();
406 fZipBuffer = new char[zipBufferSize + 9];
407 int dataSize = Length();
408 int compressedSize = 0;
409 R__zipMultipleAlgorithm(compressionLevel, &dataSize, Buffer(), &zipBufferSize, fZipBuffer, &compressedSize,
410 compressionAlgorithm);
411 if (compressedSize > 0) {
412 src = fZipBuffer;
413 srcSize = compressedSize;
414 } else {
415 delete[] fZipBuffer;
416 fZipBuffer = nullptr;
417 }
418 }
419
420 TString res;
421 char sbuf[500];
422 int block = 0;
423 char *tgt = sbuf;
424 int srcCnt = 0;
425
426 while (srcCnt++ < srcSize) {
427 tgt += sprintf(tgt, " %02x", (unsigned char)*src);
428 src++;
429 if (block++ == 100) {
430 res += sbuf;
431 block = 0;
432 tgt = sbuf;
433 }
434 }
435
436 if (block > 0)
437 res += sbuf;
438
439 XMLNodePointer_t blocknode = fXML->NewChild(node, nullptr, xmlio::XmlBlock, res);
440 fXML->NewIntAttr(blocknode, xmlio::Size, Length());
441
442 if (fZipBuffer) {
443 fXML->NewIntAttr(blocknode, xmlio::Zip, srcSize);
444 delete[] fZipBuffer;
445 }
446}
447
448////////////////////////////////////////////////////////////////////////////////
449/// Read binary block of data from xml.
450
452{
453 if (!blocknode)
454 return;
455
456 Int_t blockSize = fXML->GetIntAttr(blocknode, xmlio::Size);
457 Bool_t blockCompressed = fXML->HasAttr(blocknode, xmlio::Zip);
458 char *fUnzipBuffer = nullptr;
459
460 if (gDebug > 2)
461 Info("XmlReadBlock", "Block size = %d, Length = %d, Compressed = %d", blockSize, Length(), blockCompressed);
462
463 if (blockSize > BufferSize())
464 Expand(blockSize);
465
466 char *tgt = Buffer();
467 Int_t readSize = blockSize;
468
469 TString content = fXML->GetNodeContent(blocknode);
470
471 if (blockCompressed) {
472 Int_t zipSize = fXML->GetIntAttr(blocknode, xmlio::Zip);
473 fUnzipBuffer = new char[zipSize];
474
475 tgt = fUnzipBuffer;
476 readSize = zipSize;
477 }
478
479 char *ptr = (char *)content.Data();
480
481 if (gDebug > 3)
482 Info("XmlReadBlock", "Content %s", ptr);
483
484 for (int i = 0; i < readSize; i++) {
485 while ((*ptr < 48) || ((*ptr > 57) && (*ptr < 97)) || (*ptr > 102))
486 ptr++;
487
488 int b_hi = (*ptr > 57) ? *ptr - 87 : *ptr - 48;
489 ptr++;
490 int b_lo = (*ptr > 57) ? *ptr - 87 : *ptr - 48;
491 ptr++;
492
493 *tgt = b_hi * 16 + b_lo;
494 tgt++;
495
496 if (gDebug > 4)
497 Info("XmlReadBlock", " Buf[%d] = %d", i, b_hi * 16 + b_lo);
498 }
499
500 if (fUnzipBuffer) {
501
502 int srcsize(0), tgtsize(0), unzipRes(0);
503 int status = R__unzip_header(&srcsize, (UChar_t *)fUnzipBuffer, &tgtsize);
504
505 if (status == 0)
506 R__unzip(&readSize, (unsigned char *)fUnzipBuffer, &blockSize, (unsigned char *)Buffer(), &unzipRes);
507
508 if (status != 0 || unzipRes != blockSize)
509 Error("XmlReadBlock", "Decompression error %d", unzipRes);
510 else if (gDebug > 2)
511 Info("XmlReadBlock", "Unzip ok");
512
513 delete[] fUnzipBuffer;
514 }
515}
516
517////////////////////////////////////////////////////////////////////////////////
518/// Add "ptr" attribute to node, if ptr is null or
519/// if ptr is pointer on object, which is already saved in buffer
520/// Automatically add "ref" attribute to node, where referenced object is stored
521
523{
524 if (!node)
525 return kFALSE;
526
527 TString refvalue;
528
529 if (!ptr) {
530 refvalue = xmlio::Null; // null
531 } else {
533 if (!refnode)
534 return kFALSE;
535
536 if (fXML->HasAttr(refnode, xmlio::Ref)) {
537 refvalue = fXML->GetAttr(refnode, xmlio::Ref);
538 } else {
539 refvalue = xmlio::IdBase;
540 if (XmlFile())
541 refvalue += XmlFile()->GetNextRefCounter();
542 else
543 refvalue += GetNextRefCounter();
544 fXML->NewAttr(refnode, nullptr, xmlio::Ref, refvalue.Data());
545 }
546 }
547 if (refvalue.Length() > 0) {
548 fXML->NewAttr(node, nullptr, xmlio::Ptr, refvalue.Data());
549 return kTRUE;
550 }
551
552 return kFALSE;
553}
554
555////////////////////////////////////////////////////////////////////////////////
556/// Searches for "ptr" attribute and returns pointer to object and class,
557/// if "ptr" attribute reference to read object
558
560{
561 cl = nullptr;
562
563 if (!fXML->HasAttr(node, xmlio::Ptr))
564 return kFALSE;
565
566 const char *ptrid = fXML->GetAttr(node, xmlio::Ptr);
567
568 if (!ptrid)
569 return kFALSE;
570
571 // null
572 if (strcmp(ptrid, xmlio::Null) == 0) {
573 ptr = nullptr;
574 return kTRUE;
575 }
576
577 if (strncmp(ptrid, xmlio::IdBase, strlen(xmlio::IdBase)) != 0) {
578 Error("ExtractPointer", "Pointer tag %s not started from %s", ptrid, xmlio::IdBase);
579 return kFALSE;
580 }
581
582 Int_t id = TString(ptrid + strlen(xmlio::IdBase)).Atoi();
583
584 GetMappedObject(id + 1, ptr, cl);
585
586 if (!ptr || !cl)
587 Error("ExtractPointer", "not found ptr %s result %p %s", ptrid, ptr, (cl ? cl->GetName() : "null"));
588
589 return ptr && cl;
590}
591
592////////////////////////////////////////////////////////////////////////////////
593/// Analyze if node has "ref" attribute and register it to object map
594
595void TBufferXML::ExtractReference(XMLNodePointer_t node, const void *ptr, const TClass *cl)
596{
597 if (!node || !ptr)
598 return;
599
600 const char *refid = fXML->GetAttr(node, xmlio::Ref);
601
602 if (!refid)
603 return;
604
605 if (strncmp(refid, xmlio::IdBase, strlen(xmlio::IdBase)) != 0) {
606 Error("ExtractReference", "Reference tag %s not started from %s", refid, xmlio::IdBase);
607 return;
608 }
609
610 Int_t id = TString(refid + strlen(xmlio::IdBase)).Atoi();
611
612 MapObject(ptr, cl, id + 1);
613
614 if (gDebug > 2)
615 Info("ExtractReference", "Find reference %s for object %p class %s", refid, ptr, (cl ? cl->GetName() : "null"));
616}
617
618////////////////////////////////////////////////////////////////////////////////
619/// Check if node has specified name
620
621Bool_t TBufferXML::VerifyNode(XMLNodePointer_t node, const char *name, const char *errinfo)
622{
623 if (!name || !node)
624 return kFALSE;
625
626 if (strcmp(fXML->GetNodeName(node), name) != 0) {
627 if (errinfo) {
628 Error("VerifyNode", "Reading XML file (%s). Get: %s, expects: %s", errinfo, fXML->GetNodeName(node), name);
629 fErrorFlag = 1;
630 }
631 return kFALSE;
632 }
633 return kTRUE;
634}
635
636////////////////////////////////////////////////////////////////////////////////
637/// Check, if stack node has specified name
638
639Bool_t TBufferXML::VerifyStackNode(const char *name, const char *errinfo)
640{
641 return VerifyNode(StackNode(), name, errinfo);
642}
643
644////////////////////////////////////////////////////////////////////////////////
645/// Checks, that attribute of specified name exists and has specified value
646
647Bool_t TBufferXML::VerifyAttr(XMLNodePointer_t node, const char *name, const char *value, const char *errinfo)
648{
649 if (!node || !name || !value)
650 return kFALSE;
651
652 const char *cont = fXML->GetAttr(node, name);
653 if ((!cont || (strcmp(cont, value) != 0))) {
654 if (errinfo) {
655 Error("VerifyAttr", "%s : attr %s = %s, expected: %s", errinfo, name, cont, value);
656 fErrorFlag = 1;
657 }
658 return kFALSE;
659 }
660 return kTRUE;
661}
662
663////////////////////////////////////////////////////////////////////////////////
664/// Checks stack attribute
665
666Bool_t TBufferXML::VerifyStackAttr(const char *name, const char *value, const char *errinfo)
667{
668 return VerifyAttr(StackNode(), name, value, errinfo);
669}
670
671////////////////////////////////////////////////////////////////////////////////
672/// Create item node of specified name
673
675{
676 XMLNodePointer_t node = nullptr;
677 if (GetXmlLayout() == kGeneralized) {
678 node = fXML->NewChild(StackNode(), nullptr, xmlio::Item);
679 fXML->NewAttr(node, nullptr, xmlio::Name, name);
680 } else
681 node = fXML->NewChild(StackNode(), nullptr, name);
682 return node;
683}
684
685////////////////////////////////////////////////////////////////////////////////
686/// Checks, if stack node is item and has specified name
687
688Bool_t TBufferXML::VerifyItemNode(const char *name, const char *errinfo)
689{
690 Bool_t res = kTRUE;
691 if (GetXmlLayout() == kGeneralized)
692 res = VerifyStackNode(xmlio::Item, errinfo) && VerifyStackAttr(xmlio::Name, name, errinfo);
693 else
694 res = VerifyStackNode(name, errinfo);
695 return res;
696}
697
698////////////////////////////////////////////////////////////////////////////////
699/// Create xml node correspondent to TStreamerElement object
700
702{
703 XMLNodePointer_t elemnode = nullptr;
704
705 const char *elemxmlname = XmlGetElementName(elem);
706
707 if (GetXmlLayout() == kGeneralized) {
708 elemnode = fXML->NewChild(StackNode(), nullptr, xmlio::Member);
709 fXML->NewAttr(elemnode, nullptr, xmlio::Name, elemxmlname);
710 } else {
711 // take namesapce for element only if it is not a base class or class name
713 if ((elem->GetType() == TStreamerInfo::kBase) ||
714 ((elem->GetType() == TStreamerInfo::kTNamed) && !strcmp(elem->GetName(), TNamed::Class()->GetName())) ||
715 ((elem->GetType() == TStreamerInfo::kTObject) && !strcmp(elem->GetName(), TObject::Class()->GetName())) ||
716 ((elem->GetType() == TStreamerInfo::kTString) && !strcmp(elem->GetName(), TString::Class()->GetName())))
717 ns = nullptr;
718
719 elemnode = fXML->NewChild(StackNode(), ns, elemxmlname);
720 }
721
722 TXMLStackObj *curr = PushStack(elemnode);
723 curr->fElem = (TStreamerElement *)elem;
724}
725
726////////////////////////////////////////////////////////////////////////////////
727/// Checks if stack node correspond to TStreamerElement object
728
730{
731 const char *elemxmlname = XmlGetElementName(elem);
732
733 if (GetXmlLayout() == kGeneralized) {
735 return kFALSE;
736 if (!VerifyStackAttr(xmlio::Name, elemxmlname))
737 return kFALSE;
738 } else {
739 if (!VerifyStackNode(elemxmlname))
740 return kFALSE;
741 }
742
744
745 TXMLStackObj *curr = PushStack(StackNode()); // set pointer to first data inside element
746 curr->fElem = (TStreamerElement *)elem;
747 return kTRUE;
748}
749
750////////////////////////////////////////////////////////////////////////////////
751/// Write object to buffer
752/// If object was written before, only pointer will be stored
753/// Return pointer to top xml node, representing object
754
755XMLNodePointer_t TBufferXML::XmlWriteObject(const void *obj, const TClass *cl, Bool_t cacheReuse)
756{
757 XMLNodePointer_t objnode = fXML->NewChild(StackNode(), nullptr, xmlio::Object);
758
759 if (!cl)
760 obj = nullptr;
761
762 if (ProcessPointer(obj, objnode))
763 return objnode;
764
765 TString clname = XmlConvertClassName(cl->GetName());
766
767 fXML->NewAttr(objnode, nullptr, xmlio::ObjClass, clname);
768
769 if (cacheReuse)
770 fMap->Add(Void_Hash(obj), (Longptr_t)obj, (Longptr_t)objnode);
771
772 PushStack(objnode);
773
774 ((TClass *)cl)->Streamer((void *)obj, *this);
775
776 PopStack();
777
778 if (gDebug > 1)
779 Info("XmlWriteObject", "Done write for class: %s", cl ? cl->GetName() : "null");
780
781 return objnode;
782}
783
784////////////////////////////////////////////////////////////////////////////////
785/// Read object from the buffer
786
787void *TBufferXML::XmlReadObject(void *obj, TClass **cl)
788{
789 if (cl)
790 *cl = nullptr;
791
792 XMLNodePointer_t objnode = StackNode();
793
794 if (fErrorFlag > 0)
795 return obj;
796
797 if (!objnode)
798 return obj;
799
800 if (!VerifyNode(objnode, xmlio::Object, "XmlReadObjectNew"))
801 return obj;
802
803 TClass *objClass = nullptr;
804
805 if (ExtractPointer(objnode, obj, objClass)) {
806 ShiftStack("readobjptr");
807 if (cl)
808 *cl = objClass;
809 return obj;
810 }
811
812 TString clname = fXML->GetAttr(objnode, xmlio::ObjClass);
813 objClass = XmlDefineClass(clname);
814 if (objClass == TDirectory::Class())
815 objClass = TDirectoryFile::Class();
816
817 if (!objClass) {
818 Error("XmlReadObject", "Cannot find class %s", clname.Data());
819 ShiftStack("readobjerr");
820 return obj;
821 }
822
823 if (gDebug > 1)
824 Info("XmlReadObject", "Reading object of class %s", clname.Data());
825
826 if (!obj)
827 obj = objClass->New();
828
829 ExtractReference(objnode, obj, objClass);
830
831 PushStack(objnode);
832
833 objClass->Streamer((void *)obj, *this);
834
835 PopStack();
836
837 ShiftStack("readobj");
838
839 if (gDebug > 1)
840 Info("XmlReadObject", "Reading object of class %s done", clname.Data());
841
842 if (cl)
843 *cl = objClass;
844
845 return obj;
846}
847
848////////////////////////////////////////////////////////////////////////////////
849/// Function is called from TStreamerInfo WriteBuffer and ReadBuffer functions
850/// and indent new level in xml structure.
851/// This call indicates, that TStreamerInfo functions starts streaming
852/// object data of correspondent class
853
855{
857}
858
859////////////////////////////////////////////////////////////////////////////////
860/// Prepares buffer to stream data of specified class.
861
863{
865
866 if (sinfo)
867 cl = sinfo->GetClass();
868
869 if (!cl)
870 return;
871
872 TString clname = XmlConvertClassName(cl->GetName());
873
874 if (gDebug > 2)
875 Info("IncrementLevel", "Class: %s", clname.Data());
876
877 Bool_t compressClassNode = (fExpectedBaseClass == cl);
878 fExpectedBaseClass = nullptr;
879
880 TXMLStackObj *stack = Stack();
881
882 if (IsWriting()) {
883
884 XMLNodePointer_t classnode = nullptr;
885 if (compressClassNode) {
886 classnode = StackNode();
887 } else {
888 if (GetXmlLayout() == kGeneralized) {
889 classnode = fXML->NewChild(StackNode(), nullptr, xmlio::Class);
890 fXML->NewAttr(classnode, nullptr, "name", clname);
891 } else
892 classnode = fXML->NewChild(StackNode(), nullptr, clname);
893 stack = PushStack(classnode);
894 }
895
896 if (fVersionBuf >= -1) {
897 if (fVersionBuf == -1)
898 fVersionBuf = 1;
900 fVersionBuf = -111;
901 }
902
904 stack->fClassNs = fXML->NewNS(classnode, XmlClassNameSpaceRef(cl), clname);
905
906 } else {
907 if (!compressClassNode) {
908 if (GetXmlLayout() == kGeneralized) {
909 if (!VerifyStackNode(xmlio::Class, "StartInfo"))
910 return;
911 if (!VerifyStackAttr("name", clname, "StartInfo"))
912 return;
913 } else if (!VerifyStackNode(clname, "StartInfo"))
914 return;
915 stack = PushStack(StackNode());
916 }
917 }
918
919 stack->fCompressedClassNode = compressClassNode;
920 stack->fInfo = sinfo;
921 stack->fIsStreamerInfo = kTRUE;
922}
923
924////////////////////////////////////////////////////////////////////////////////
925/// Function is called from TStreamerInfo WriteBuffer and ReadBuffer functions
926/// and decrease level in xml structure.
927
929{
931
933
934 if (gDebug > 2)
935 Info("DecrementLevel", "Class: %s", (info ? info->GetClass()->GetName() : "custom"));
936
937 TXMLStackObj *stack = Stack();
938
939 if (!stack->IsStreamerInfo()) {
941 stack = PopStack(); // remove stack of last element
942 }
943
944 if (stack->fCompressedClassNode) {
945 stack->fInfo = nullptr;
946 stack->fIsStreamerInfo = kFALSE;
948 } else {
949 PopStack(); // back from data of stack info
950 if (IsReading())
951 ShiftStack("declevel"); // shift to next element after streamer info
952 }
953}
954
955////////////////////////////////////////////////////////////////////////////////
956/// Function is called from TStreamerInfo WriteBuffer and ReadBuffer functions
957/// and add/verify next element of xml structure
958/// This calls allows separate data, correspondent to one class member, from another
959
961{
962 WorkWithElement(elem, comptype);
963}
964
965////////////////////////////////////////////////////////////////////////////////
966/// This function is a part of SetStreamerElementNumber method.
967/// It is introduced for reading of data for specified data member of class.
968/// Used also in ReadFastArray methods to resolve problem of compressed data,
969/// when several data members of the same basic type streamed with single ...FastArray call
970
972{
974
976 fExpectedBaseClass = nullptr;
977
978 TXMLStackObj *stack = Stack();
979 if (!stack) {
980 Error("SetStreamerElementNumber", "stack is empty");
981 return;
982 }
983
984 if (!stack->IsStreamerInfo()) { // this is not a first element
986 PopStack(); // go level back
987 if (IsReading())
988 ShiftStack("startelem"); // shift to next element, only for reading
989 stack = Stack();
990 }
991
992 if (!stack) {
993 Error("SetStreamerElementNumber", "Lost of stack");
994 return;
995 }
996
997 if (!elem) {
998 Error("SetStreamerElementNumber", "Problem in Inc/Dec level");
999 return;
1000 }
1001
1002 TStreamerInfo *info = stack->fInfo;
1003
1004 if (!stack->IsStreamerInfo()) {
1005 Error("SetStreamerElementNumber", "Problem in Inc/Dec level");
1006 return;
1007 }
1008 Int_t number = info ? info->GetElements()->IndexOf(elem) : -1;
1009
1010 if (gDebug > 4)
1011 Info("SetStreamerElementNumber", " Next element %s", elem->GetName());
1012
1013 Bool_t isBasicType = (elem->GetType() > 0) && (elem->GetType() < 20);
1014
1016 isBasicType && ((elem->GetType() == comp_type) || (elem->GetType() == comp_type - TStreamerInfo::kConv) ||
1017 (elem->GetType() == comp_type - TStreamerInfo::kSkip));
1018
1019 if ((elem->GetType() == TStreamerInfo::kBase) ||
1020 ((elem->GetType() == TStreamerInfo::kTNamed) && !strcmp(elem->GetName(), TNamed::Class()->GetName())))
1022
1023 if (fExpectedBaseClass && (gDebug > 3))
1024 Info("SetStreamerElementNumber", " Expects base class %s with standard streamer",
1026
1027 if (IsWriting()) {
1028 CreateElemNode(elem);
1029 } else {
1030 if (!VerifyElemNode(elem))
1031 return;
1032 }
1033
1034 stack = Stack();
1035 stack->fElemNumber = number;
1036 stack->fIsElemOwner = (number < 0);
1037}
1038
1039////////////////////////////////////////////////////////////////////////////////
1040/// Should be called at the beginning of custom class streamer.
1041///
1042/// Informs buffer data about class which will be streamed now.
1043/// ClassBegin(), ClassEnd() and ClassMemeber() should be used in
1044/// custom class streamers to specify which kind of data are
1045/// now streamed. Such information is used to correctly
1046/// convert class data to XML. Without that functions calls
1047/// classes with custom streamers cannot be used with TBufferXML
1048
1050{
1051 WorkWithClass(nullptr, cl);
1052}
1053
1054////////////////////////////////////////////////////////////////////////////////
1055/// Should be called at the end of custom streamer
1056/// See TBufferXML::ClassBegin for more details
1057
1059{
1060 DecrementLevel(0);
1061}
1062
1063////////////////////////////////////////////////////////////////////////////////
1064/// Method indicates name and typename of class member,
1065/// which should be now streamed in custom streamer
1066///
1067/// Following combinations are supported:
1068/// -# name = "ClassName", typeName = 0 or typename==ClassName.
1069/// This is a case, when data of parent class "ClassName" should be streamed.
1070/// For instance, if class directly inherited from TObject, custom streamer
1071/// should include following code:
1072/// ~~~{.cpp}
1073/// b.ClassMember("TObject");
1074/// TObject::Streamer(b);
1075/// ~~~
1076/// -# Basic data type
1077/// ~~~{.cpp}
1078/// b.ClassMember("fInt","Int_t");
1079/// b >> fInt;
1080/// ~~~
1081/// -# Array of basic data types
1082/// ~~~{.cpp}
1083/// b.ClassMember("fArr","Int_t", 5);
1084/// b.ReadFastArray(fArr, 5);
1085/// ~~~
1086/// -# Object as data member
1087/// ~~~{.cpp}
1088/// b.ClassMemeber("fName","TString");
1089/// fName.Streamer(b);
1090/// ~~~
1091/// -# Pointer on object as data member
1092/// ~~~{.cpp}
1093/// b.ClassMemeber("fObj","TObject*");
1094/// b.StreamObject(fObj);
1095/// ~~~
1096///
1097/// Arrsize1 and arrsize2 arguments (when specified) indicate first and
1098/// second dimension of array. Can be used for array of basic types.
1099/// See ClassBegin() method for more details.
1100
1101void TBufferXML::ClassMember(const char *name, const char *typeName, Int_t arrsize1, Int_t arrsize2)
1102{
1103 if (!typeName)
1104 typeName = name;
1105
1106 if (!name || (strlen(name) == 0)) {
1107 Error("ClassMember", "Invalid member name");
1108 fErrorFlag = 1;
1109 return;
1110 }
1111
1112 TString tname = typeName;
1113
1114 Int_t typ_id(-1), comp_type(-1);
1115
1116 if (strcmp(typeName, "raw:data") == 0)
1117 typ_id = TStreamerInfo::kMissing;
1118
1119 if (typ_id < 0) {
1120 TDataType *dt = gROOT->GetType(typeName);
1121 if (dt)
1122 if ((dt->GetType() > 0) && (dt->GetType() < 20))
1123 typ_id = dt->GetType();
1124 }
1125
1126 if (typ_id < 0)
1127 if (strcmp(name, typeName) == 0) {
1128 TClass *cl = TClass::GetClass(tname.Data());
1129 if (cl)
1130 typ_id = TStreamerInfo::kBase;
1131 }
1132
1133 if (typ_id < 0) {
1134 Bool_t isptr = kFALSE;
1135 if (tname[tname.Length() - 1] == '*') {
1136 tname.Resize(tname.Length() - 1);
1137 isptr = kTRUE;
1138 }
1139 TClass *cl = TClass::GetClass(tname.Data());
1140 if (!cl) {
1141 Error("ClassMember", "Invalid class specifier %s", typeName);
1142 fErrorFlag = 1;
1143 return;
1144 }
1145
1146 if (cl->IsTObject())
1148 else
1149 typ_id = isptr ? TStreamerInfo::kAnyp : TStreamerInfo::kAny;
1150
1151 if ((cl == TString::Class()) && !isptr)
1152 typ_id = TStreamerInfo::kTString;
1153 }
1154
1155 TStreamerElement *elem = nullptr;
1156
1157 if (typ_id == TStreamerInfo::kMissing) {
1158 elem = new TStreamerElement(name, "title", 0, typ_id, "raw:data");
1159 } else if (typ_id == TStreamerInfo::kBase) {
1160 TClass *cl = TClass::GetClass(tname.Data());
1161 if (cl) {
1162 TStreamerBase *b = new TStreamerBase(tname.Data(), "title", 0);
1163 b->SetBaseVersion(cl->GetClassVersion());
1164 elem = b;
1165 }
1166 } else if ((typ_id > 0) && (typ_id < 20)) {
1167 elem = new TStreamerBasicType(name, "title", 0, typ_id, typeName);
1168 comp_type = typ_id;
1169 } else if ((typ_id == TStreamerInfo::kObject) || (typ_id == TStreamerInfo::kTObject) ||
1170 (typ_id == TStreamerInfo::kTNamed)) {
1171 elem = new TStreamerObject(name, "title", 0, tname.Data());
1172 } else if (typ_id == TStreamerInfo::kObjectp) {
1173 elem = new TStreamerObjectPointer(name, "title", 0, tname.Data());
1174 } else if (typ_id == TStreamerInfo::kAny) {
1175 elem = new TStreamerObjectAny(name, "title", 0, tname.Data());
1176 } else if (typ_id == TStreamerInfo::kAnyp) {
1177 elem = new TStreamerObjectAnyPointer(name, "title", 0, tname.Data());
1178 } else if (typ_id == TStreamerInfo::kTString) {
1179 elem = new TStreamerString(name, "title", 0);
1180 }
1181
1182 if (!elem) {
1183 Error("ClassMember", "Invalid combination name = %s type = %s", name, typeName);
1184 fErrorFlag = 1;
1185 return;
1186 }
1187
1188 if (arrsize1 > 0) {
1189 elem->SetArrayDim(arrsize2 > 0 ? 2 : 1);
1190 elem->SetMaxIndex(0, arrsize1);
1191 if (arrsize2 > 0)
1192 elem->SetMaxIndex(1, arrsize2);
1193 }
1194
1195 // we indicate that there is no streamerinfo
1196 WorkWithElement(elem, comp_type);
1197}
1198
1199////////////////////////////////////////////////////////////////////////////////
1200/// Function is converts TObject and TString structures to more compact representation
1201
1203{
1204 if (GetXmlLayout() == kGeneralized)
1205 return;
1206
1207 const TStreamerElement *elem = Stack()->fElem;
1208 XMLNodePointer_t elemnode = IsWriting() ? Stack()->fNode : Stack(1)->fNode;
1209
1210 if (!elem || !elemnode)
1211 return;
1212
1213 if (elem->GetType() == TStreamerInfo::kTString) {
1214
1215 XMLNodePointer_t node = fXML->GetChild(elemnode);
1216 fXML->SkipEmpty(node);
1217
1218 XMLNodePointer_t nodecharstar(nullptr), nodeuchar(nullptr), nodeint(nullptr), nodestring(nullptr);
1219
1220 while (node) {
1221 const char *name = fXML->GetNodeName(node);
1222 if (strcmp(name, xmlio::String) == 0) {
1223 if (nodestring)
1224 return;
1225 nodestring = node;
1226 } else if (strcmp(name, xmlio::UChar) == 0) {
1227 if (nodeuchar)
1228 return;
1229 nodeuchar = node;
1230 } else if (strcmp(name, xmlio::Int) == 0) {
1231 if (nodeint)
1232 return;
1233 nodeint = node;
1234 } else if (strcmp(name, xmlio::CharStar) == 0) {
1235 if (nodecharstar)
1236 return;
1237 nodecharstar = node;
1238 } else
1239 return; // can not be something else
1240 fXML->ShiftToNext(node);
1241 }
1242
1243 TString str;
1244
1245 if (GetIOVersion() < 3) {
1246 if (!nodeuchar)
1247 return;
1248 if (nodecharstar)
1249 str = fXML->GetAttr(nodecharstar, xmlio::v);
1250 fXML->UnlinkFreeNode(nodeuchar);
1251 fXML->UnlinkFreeNode(nodeint);
1252 fXML->UnlinkFreeNode(nodecharstar);
1253 } else {
1254 if (nodestring)
1255 str = fXML->GetAttr(nodestring, xmlio::v);
1256 fXML->UnlinkFreeNode(nodestring);
1257 }
1258
1259 fXML->NewAttr(elemnode, nullptr, "str", str);
1260 } else if (elem->GetType() == TStreamerInfo::kTObject) {
1261 XMLNodePointer_t node = fXML->GetChild(elemnode);
1262 fXML->SkipEmpty(node);
1263
1264 XMLNodePointer_t vnode = nullptr, idnode = nullptr, bitsnode = nullptr, prnode = nullptr;
1265
1266 while (node) {
1267 const char *name = fXML->GetNodeName(node);
1268
1269 if (strcmp(name, xmlio::OnlyVersion) == 0) {
1270 if (vnode)
1271 return;
1272 vnode = node;
1273 } else if (strcmp(name, xmlio::UInt) == 0) {
1274 if (!idnode)
1275 idnode = node;
1276 else if (!bitsnode)
1277 bitsnode = node;
1278 else
1279 return;
1280 } else if (strcmp(name, xmlio::UShort) == 0) {
1281 if (prnode)
1282 return;
1283 prnode = node;
1284 } else
1285 return;
1286 fXML->ShiftToNext(node);
1287 }
1288
1289 if (!vnode || !idnode || !bitsnode)
1290 return;
1291
1292 TString str = fXML->GetAttr(idnode, xmlio::v);
1293 fXML->NewAttr(elemnode, nullptr, "fUniqueID", str);
1294
1295 str = fXML->GetAttr(bitsnode, xmlio::v);
1296 UInt_t bits;
1297 sscanf(str.Data(), "%u", &bits);
1299
1300 char sbuf[20];
1301 snprintf(sbuf, sizeof(sbuf), "%x", bits);
1302 fXML->NewAttr(elemnode, nullptr, "fBits", sbuf);
1303
1304 if (prnode) {
1305 str = fXML->GetAttr(prnode, xmlio::v);
1306 fXML->NewAttr(elemnode, nullptr, "fProcessID", str);
1307 }
1308
1309 fXML->UnlinkFreeNode(vnode);
1310 fXML->UnlinkFreeNode(idnode);
1311 fXML->UnlinkFreeNode(bitsnode);
1312 fXML->UnlinkFreeNode(prnode);
1313 }
1314}
1315
1316////////////////////////////////////////////////////////////////////////////////
1317/// Function is unpack TObject and TString structures to be able read
1318/// them from custom streamers of this objects
1319
1321{
1322 if (GetXmlLayout() == kGeneralized)
1323 return;
1324 if (!elem || !elemnode)
1325 return;
1326
1327 if (elem->GetType() == TStreamerInfo::kTString) {
1328
1329 if (!fXML->HasAttr(elemnode, "str"))
1330 return;
1331 TString str = fXML->GetAttr(elemnode, "str");
1332 fXML->FreeAttr(elemnode, "str");
1333
1334 if (GetIOVersion() < 3) {
1335 Int_t len = str.Length();
1336 XMLNodePointer_t ucharnode = fXML->NewChild(elemnode, nullptr, xmlio::UChar);
1337 char sbuf[20];
1338 snprintf(sbuf, sizeof(sbuf), "%d", len);
1339 if (len < 255) {
1340 fXML->NewAttr(ucharnode, nullptr, xmlio::v, sbuf);
1341 } else {
1342 fXML->NewAttr(ucharnode, nullptr, xmlio::v, "255");
1343 XMLNodePointer_t intnode = fXML->NewChild(elemnode, nullptr, xmlio::Int);
1344 fXML->NewAttr(intnode, nullptr, xmlio::v, sbuf);
1345 }
1346 if (len > 0) {
1347 XMLNodePointer_t node = fXML->NewChild(elemnode, nullptr, xmlio::CharStar);
1348 fXML->NewAttr(node, nullptr, xmlio::v, str);
1349 }
1350 } else {
1351 XMLNodePointer_t node = fXML->NewChild(elemnode, nullptr, xmlio::String);
1352 fXML->NewAttr(node, nullptr, xmlio::v, str);
1353 }
1354 } else if (elem->GetType() == TStreamerInfo::kTObject) {
1355 if (!fXML->HasAttr(elemnode, "fUniqueID"))
1356 return;
1357 if (!fXML->HasAttr(elemnode, "fBits"))
1358 return;
1359
1360 TString idstr = fXML->GetAttr(elemnode, "fUniqueID");
1361 TString bitsstr = fXML->GetAttr(elemnode, "fBits");
1362 TString prstr = fXML->GetAttr(elemnode, "fProcessID");
1363
1364 fXML->FreeAttr(elemnode, "fUniqueID");
1365 fXML->FreeAttr(elemnode, "fBits");
1366 fXML->FreeAttr(elemnode, "fProcessID");
1367
1368 XMLNodePointer_t node = fXML->NewChild(elemnode, nullptr, xmlio::OnlyVersion);
1369 fXML->NewAttr(node, nullptr, xmlio::v, "1");
1370
1371 node = fXML->NewChild(elemnode, nullptr, xmlio::UInt);
1372 fXML->NewAttr(node, nullptr, xmlio::v, idstr);
1373
1374 UInt_t bits = 0;
1375 sscanf(bitsstr.Data(), "%x", &bits);
1377 char sbuf[20];
1378 snprintf(sbuf, sizeof(sbuf), "%u", bits);
1379
1380 node = fXML->NewChild(elemnode, nullptr, xmlio::UInt);
1381 fXML->NewAttr(node, nullptr, xmlio::v, sbuf);
1382
1383 if (prstr.Length() > 0) {
1384 node = fXML->NewChild(elemnode, nullptr, xmlio::UShort);
1385 fXML->NewAttr(node, nullptr, xmlio::v, prstr.Data());
1386 }
1387 }
1388}
1389
1390////////////////////////////////////////////////////////////////////////////////
1391/// Function is called before any IO operation of TBuffer
1392/// Now is used to store version value if no proper calls are discovered
1393
1395{
1397}
1398
1399////////////////////////////////////////////////////////////////////////////////
1400/// Function to read class from buffer, used in old-style streamers
1401
1403{
1404 const char *clname = nullptr;
1405
1407 clname = XmlReadValue(xmlio::Class);
1408
1409 if (gDebug > 2)
1410 Info("ReadClass", "Try to read class %s", clname ? clname : "---");
1411
1412 return clname ? gROOT->GetClass(clname) : nullptr;
1413}
1414
1415////////////////////////////////////////////////////////////////////////////////
1416/// Function to write class into buffer, used in old-style streamers
1417
1419{
1420 if (gDebug > 2)
1421 Info("WriteClass", "Try to write class %s", cl->GetName());
1422
1424}
1425
1426////////////////////////////////////////////////////////////////////////////////
1427/// Read version value from buffer
1428
1430{
1432
1433 Version_t res = 0;
1434
1435 if (start)
1436 *start = 0;
1437 if (bcnt)
1438 *bcnt = 0;
1439
1442 } else if (fExpectedBaseClass && (fXML->HasAttr(Stack(1)->fNode, xmlio::ClassVersion))) {
1443 res = fXML->GetIntAttr(Stack(1)->fNode, xmlio::ClassVersion);
1444 } else if (fXML->HasAttr(StackNode(), xmlio::ClassVersion)) {
1446 } else {
1447 Error("ReadVersion", "No correspondent tags to read version");
1448 fErrorFlag = 1;
1449 }
1450
1451 if (gDebug > 2)
1452 Info("ReadVersion", "Version = %d", res);
1453
1454 return res;
1455}
1456
1457////////////////////////////////////////////////////////////////////////////////
1458/// Checks buffer, filled by WriteVersion
1459/// if next data is arriving, version should be stored in buffer
1460
1462{
1463 if (IsWriting() && (fVersionBuf >= -100)) {
1464 char sbuf[20];
1465 snprintf(sbuf, sizeof(sbuf), "%d", fVersionBuf);
1467 fVersionBuf = -111;
1468 }
1469}
1470
1471////////////////////////////////////////////////////////////////////////////////
1472/// Copies class version to buffer, but not writes it to xml
1473/// Version will be written with next I/O operation or
1474/// will be added as attribute of class tag, created by IncrementLevel call
1475
1477{
1479
1480 if (fExpectedBaseClass != cl)
1481 fExpectedBaseClass = nullptr;
1482
1484
1485 if (gDebug > 2)
1486 Info("WriteVersion", "Class: %s, version = %d", cl->GetName(), fVersionBuf);
1487
1488 return 0;
1489}
1490
1491////////////////////////////////////////////////////////////////////////////////
1492/// Read object from buffer. Only used from TBuffer
1493
1495{
1497 if (gDebug > 2)
1498 Info("ReadObjectAny", "From node %s", fXML->GetNodeName(StackNode()));
1499 void *res = XmlReadObject(nullptr);
1500 return res;
1501}
1502
1503////////////////////////////////////////////////////////////////////////////////
1504/// Skip any kind of object from buffer
1505/// Actually skip only one node on current level of xml structure
1506
1508{
1509 ShiftStack("skipobjectany");
1510}
1511
1512////////////////////////////////////////////////////////////////////////////////
1513/// Write object to buffer. Only used from TBuffer
1514
1515void TBufferXML::WriteObjectClass(const void *actualObjStart, const TClass *actualClass, Bool_t cacheReuse)
1516{
1518 if (gDebug > 2)
1519 Info("WriteObject", "Class %s", (actualClass ? actualClass->GetName() : " null"));
1520 XmlWriteObject(actualObjStart, actualClass, cacheReuse);
1521}
1522
1523////////////////////////////////////////////////////////////////////////////////
1524/// Template method to read array content
1525
1526template <typename T>
1528{
1529 Int_t indx = 0, cnt, curr;
1530 while (indx < arrsize) {
1531 cnt = 1;
1534 XmlReadBasic(arr[indx]);
1535 curr = indx++;
1536 while (cnt-- > 1)
1537 arr[indx++] = arr[curr];
1538 }
1539}
1540
1541////////////////////////////////////////////////////////////////////////////////
1542/// Template method to read array with size attribute
1543/// If necessary, array is created
1544
1545template <typename T>
1547{
1549 if (!VerifyItemNode(xmlio::Array, is_static ? "ReadStaticArray" : "ReadArray"))
1550 return 0;
1552 if (n <= 0)
1553 return 0;
1554 if (!arr) {
1555 if (is_static)
1556 return 0;
1557 arr = new T[n];
1558 }
1560 XmlReadArrayContent(arr, n);
1561 PopStack();
1562 ShiftStack(is_static ? "readstatarr" : "readarr");
1563 return n;
1564}
1565
1566////////////////////////////////////////////////////////////////////////////////
1567/// Read array of Bool_t from buffer
1568
1570{
1571 return XmlReadArray(b);
1572}
1573
1574////////////////////////////////////////////////////////////////////////////////
1575/// Read array of Char_t from buffer
1576
1578{
1579 return XmlReadArray(c);
1580}
1581
1582////////////////////////////////////////////////////////////////////////////////
1583/// Read array of UChar_t from buffer
1584
1586{
1587 return XmlReadArray(c);
1588}
1589
1590////////////////////////////////////////////////////////////////////////////////
1591/// Read array of Short_t from buffer
1592
1594{
1595 return XmlReadArray(h);
1596}
1597
1598////////////////////////////////////////////////////////////////////////////////
1599/// Read array of UShort_t from buffer
1600
1602{
1603 return XmlReadArray(h);
1604}
1605
1606////////////////////////////////////////////////////////////////////////////////
1607/// Read array of Int_t from buffer
1608
1610{
1611 return XmlReadArray(i);
1612}
1613
1614////////////////////////////////////////////////////////////////////////////////
1615/// Read array of UInt_t from buffer
1616
1618{
1619 return XmlReadArray(i);
1620}
1621
1622////////////////////////////////////////////////////////////////////////////////
1623/// Read array of Long_t from buffer
1624
1626{
1627 return XmlReadArray(l);
1628}
1629
1630////////////////////////////////////////////////////////////////////////////////
1631/// Read array of ULong_t from buffer
1632
1634{
1635 return XmlReadArray(l);
1636}
1637
1638////////////////////////////////////////////////////////////////////////////////
1639/// Read array of Long64_t from buffer
1640
1642{
1643 return XmlReadArray(l);
1644}
1645
1646////////////////////////////////////////////////////////////////////////////////
1647/// Read array of ULong64_t from buffer
1648
1650{
1651 return XmlReadArray(l);
1652}
1653
1654////////////////////////////////////////////////////////////////////////////////
1655/// Read array of Float_t from buffer
1656
1658{
1659 return XmlReadArray(f);
1660}
1661
1662////////////////////////////////////////////////////////////////////////////////
1663/// Read array of Double_t from buffer
1664
1666{
1667 return XmlReadArray(d);
1668}
1669
1670////////////////////////////////////////////////////////////////////////////////
1671/// Read array of Bool_t from buffer
1672
1674{
1675 return XmlReadArray(b, true);
1676}
1677
1678////////////////////////////////////////////////////////////////////////////////
1679/// Read array of Char_t from buffer
1680
1682{
1683 return XmlReadArray(c, true);
1684}
1685
1686////////////////////////////////////////////////////////////////////////////////
1687/// Read array of UChar_t from buffer
1688
1690{
1691 return XmlReadArray(c, true);
1692}
1693
1694////////////////////////////////////////////////////////////////////////////////
1695/// Read array of Short_t from buffer
1696
1698{
1699 return XmlReadArray(h, true);
1700}
1701
1702////////////////////////////////////////////////////////////////////////////////
1703/// Read array of UShort_t from buffer
1704
1706{
1707 return XmlReadArray(h, true);
1708}
1709
1710////////////////////////////////////////////////////////////////////////////////
1711/// Read array of Int_t from buffer
1712
1714{
1715 return XmlReadArray(i, true);
1716}
1717
1718////////////////////////////////////////////////////////////////////////////////
1719/// Read array of UInt_t from buffer
1720
1722{
1723 return XmlReadArray(i, true);
1724}
1725
1726////////////////////////////////////////////////////////////////////////////////
1727/// Read array of Long_t from buffer
1728
1730{
1731 return XmlReadArray(l, true);
1732}
1733
1734////////////////////////////////////////////////////////////////////////////////
1735/// Read array of ULong_t from buffer
1736
1738{
1739 return XmlReadArray(l, true);
1740}
1741
1742////////////////////////////////////////////////////////////////////////////////
1743/// Read array of Long64_t from buffer
1744
1746{
1747 return XmlReadArray(l, true);
1748}
1749
1750////////////////////////////////////////////////////////////////////////////////
1751/// Read array of ULong64_t from buffer
1752
1754{
1755 return XmlReadArray(l, true);
1756}
1757
1758////////////////////////////////////////////////////////////////////////////////
1759/// Read array of Float_t from buffer
1760
1762{
1763 return XmlReadArray(f, true);
1764}
1765
1766////////////////////////////////////////////////////////////////////////////////
1767/// Read array of Double_t from buffer
1768
1770{
1771 return XmlReadArray(d, true);
1772}
1773
1774////////////////////////////////////////////////////////////////////////////////
1775/// Template method to read content of array, which not include size of array
1776/// Also treated situation, when instead of one single array chain
1777/// of several elements should be produced
1778
1779template <typename T>
1781{
1783 if (n <= 0)
1784 return;
1785 if (!VerifyItemNode(xmlio::Array, "ReadFastArray"))
1786 return;
1788 XmlReadArrayContent(arr, n);
1789 PopStack();
1790 ShiftStack("readfastarr");
1791}
1792
1793////////////////////////////////////////////////////////////////////////////////
1794/// Read array of Bool_t from buffer
1795
1797{
1799}
1800
1801////////////////////////////////////////////////////////////////////////////////
1802/// Read array of Char_t from buffer
1803/// if nodename==CharStar, read all array as string
1804
1806{
1807 if ((n > 0) && VerifyItemNode(xmlio::CharStar)) {
1808 const char *buf;
1809 if ((buf = XmlReadValue(xmlio::CharStar))) {
1810 Int_t size = strlen(buf);
1811 if (size < n)
1812 size = n;
1813 memcpy(c, buf, size);
1814 }
1815 } else {
1817 }
1818}
1819
1820////////////////////////////////////////////////////////////////////////////////
1821/// Read array of n characters from the I/O buffer.
1822/// Used only from TLeafC, dummy implementation here
1823
1825{
1826 ReadFastArray(c, n);
1827}
1828
1829////////////////////////////////////////////////////////////////////////////////
1830/// Read array of UChar_t from buffer
1831
1833{
1835}
1836
1837////////////////////////////////////////////////////////////////////////////////
1838/// Read array of Short_t from buffer
1839
1841{
1843}
1844
1845////////////////////////////////////////////////////////////////////////////////
1846/// Read array of UShort_t from buffer
1847
1849{
1851}
1852
1853////////////////////////////////////////////////////////////////////////////////
1854/// Read array of Int_t from buffer
1855
1857{
1858 XmlReadFastArray(i, n);
1859}
1860
1861////////////////////////////////////////////////////////////////////////////////
1862/// Read array of UInt_t from buffer
1863
1865{
1866 XmlReadFastArray(i, n);
1867}
1868
1869////////////////////////////////////////////////////////////////////////////////
1870/// Read array of Long_t from buffer
1871
1873{
1875}
1876
1877////////////////////////////////////////////////////////////////////////////////
1878/// Read array of ULong_t from buffer
1879
1881{
1883}
1884
1885////////////////////////////////////////////////////////////////////////////////
1886/// Read array of Long64_t from buffer
1887
1889{
1891}
1892
1893////////////////////////////////////////////////////////////////////////////////
1894/// Read array of ULong64_t from buffer
1895
1897{
1899}
1900
1901////////////////////////////////////////////////////////////////////////////////
1902/// Read array of Float_t from buffer
1903
1905{
1907}
1908
1909////////////////////////////////////////////////////////////////////////////////
1910/// Read array of Double_t from buffer
1911
1913{
1915}
1916
1917////////////////////////////////////////////////////////////////////////////////
1918/// Read an array of 'n' objects from the I/O buffer.
1919/// Stores the objects read starting at the address 'start'.
1920/// The objects in the array are assume to be of class 'cl'.
1921
1922void TBufferXML::ReadFastArray(void *start, const TClass *cl, Int_t n, TMemberStreamer *streamer,
1923 const TClass *onFileClass)
1924{
1925 if (streamer) {
1926 streamer->SetOnFileClass(onFileClass);
1927 (*streamer)(*this, start, 0);
1928 return;
1929 }
1930
1931 int objectSize = cl->Size();
1932 char *obj = (char *)start;
1933 char *end = obj + n * objectSize;
1934
1935 for (; obj < end; obj += objectSize)
1936 ((TClass *)cl)->Streamer(obj, *this, onFileClass);
1937}
1938
1939////////////////////////////////////////////////////////////////////////////////
1940/// Read an array of 'n' objects from the I/O buffer.
1941///
1942/// The objects read are stored starting at the address '*start'
1943/// The objects in the array are assumed to be of class 'cl' or a derived class.
1944/// 'mode' indicates whether the data member is marked with '->'
1945
1946void TBufferXML::ReadFastArray(void **start, const TClass *cl, Int_t n, Bool_t isPreAlloc, TMemberStreamer *streamer,
1947 const TClass *onFileClass)
1948{
1949
1950 Bool_t oldStyle = kFALSE; // flag used to reproduce old-style I/= actions for kSTLp
1951
1952 if ((GetIOVersion() < 4) && !isPreAlloc) {
1953 TStreamerElement *elem = Stack()->fElem;
1954 if (elem && ((elem->GetType() == TStreamerInfo::kSTLp) ||
1956 oldStyle = kTRUE;
1957 }
1958
1959 if (streamer) {
1960 if (isPreAlloc) {
1961 for (Int_t j = 0; j < n; j++) {
1962 if (!start[j])
1963 start[j] = cl->New();
1964 }
1965 }
1966 streamer->SetOnFileClass(onFileClass);
1967 (*streamer)(*this, (void *)start, oldStyle ? n : 0);
1968 return;
1969 }
1970
1971 if (!isPreAlloc) {
1972
1973 for (Int_t j = 0; j < n; j++) {
1974 if (oldStyle) {
1975 if (!start[j])
1976 start[j] = ((TClass *)cl)->New();
1977 ((TClass *)cl)->Streamer(start[j], *this);
1978 continue;
1979 }
1980 // delete the object or collection
1981 void *old = start[j];
1982 start[j] = ReadObjectAny(cl);
1983 if (old && old != start[j] && TStreamerInfo::CanDelete()
1984 // There are some cases where the user may set up a pointer in the (default)
1985 // constructor but not mark this pointer as transient. Sometime the value
1986 // of this pointer is the address of one of the object with just created
1987 // and the following delete would result in the deletion (possibly of the
1988 // top level object we are goint to return!).
1989 // Eventhough this is a user error, we could prevent the crash by simply
1990 // adding:
1991 // && !CheckObject(start[j],cl)
1992 // However this can increase the read time significantly (10% in the case
1993 // of one TLine pointer in the test/Track and run ./Event 200 0 0 20 30000
1994 //
1995 // If ReadObjectAny returned the same value as we previous had, this means
1996 // that when writing this object (start[j] had already been written and
1997 // is indeed pointing to the same object as the object the user set up
1998 // in the default constructor).
1999 ) {
2000 ((TClass *)cl)->Destructor(old, kFALSE); // call delete and desctructor
2001 }
2002 }
2003
2004 } else {
2005 // case //-> in comment
2006
2007 for (Int_t j = 0; j < n; j++) {
2008 if (!start[j])
2009 start[j] = ((TClass *)cl)->New();
2010 ((TClass *)cl)->Streamer(start[j], *this, onFileClass);
2011 }
2012 }
2013}
2014
2015template <typename T>
2017{
2018 if (fCompressLevel > 0) {
2019 Int_t indx = 0;
2020 while (indx < arrsize) {
2021 XMLNodePointer_t elemnode = XmlWriteBasic(arr[indx]);
2022 Int_t curr = indx++;
2023 while ((indx < arrsize) && (arr[indx] == arr[curr]))
2024 indx++;
2025 if (indx - curr > 1)
2026 fXML->NewIntAttr(elemnode, xmlio::cnt, indx - curr);
2027 }
2028 } else {
2029 for (Int_t indx = 0; indx < arrsize; indx++)
2030 XmlWriteBasic(arr[indx]);
2031 }
2032}
2033
2034////////////////////////////////////////////////////////////////////////////////
2035/// Write array, including it size
2036/// Content may be compressed
2037
2038template <typename T>
2040{
2043 fXML->NewIntAttr(arrnode, xmlio::Size, arrsize);
2044 PushStack(arrnode);
2045 XmlWriteArrayContent(arr, arrsize);
2046 PopStack();
2047}
2048
2049////////////////////////////////////////////////////////////////////////////////
2050/// Write array of Bool_t to buffer
2051
2053{
2054 XmlWriteArray(b, n);
2055}
2056
2057////////////////////////////////////////////////////////////////////////////////
2058/// Write array of Char_t to buffer
2059
2061{
2062 XmlWriteArray(c, n);
2063}
2064
2065////////////////////////////////////////////////////////////////////////////////
2066/// Write array of UChar_t to buffer
2067
2069{
2070 XmlWriteArray(c, n);
2071}
2072
2073////////////////////////////////////////////////////////////////////////////////
2074/// Write array of Short_t to buffer
2075
2077{
2078 XmlWriteArray(h, n);
2079}
2080
2081////////////////////////////////////////////////////////////////////////////////
2082/// Write array of UShort_t to buffer
2083
2085{
2086 XmlWriteArray(h, n);
2087}
2088
2089////////////////////////////////////////////////////////////////////////////////
2090/// Write array of Int_ to buffer
2091
2093{
2094 XmlWriteArray(i, n);
2095}
2096
2097////////////////////////////////////////////////////////////////////////////////
2098/// Write array of UInt_t to buffer
2099
2101{
2102 XmlWriteArray(i, n);
2103}
2104
2105////////////////////////////////////////////////////////////////////////////////
2106/// Write array of Long_t to buffer
2107
2109{
2110 XmlWriteArray(l, n);
2111}
2112
2113////////////////////////////////////////////////////////////////////////////////
2114/// Write array of ULong_t to buffer
2115
2117{
2118 XmlWriteArray(l, n);
2119}
2120
2121////////////////////////////////////////////////////////////////////////////////
2122/// Write array of Long64_t to buffer
2123
2125{
2126 XmlWriteArray(l, n);
2127}
2128
2129////////////////////////////////////////////////////////////////////////////////
2130/// Write array of ULong64_t to buffer
2131
2133{
2134 XmlWriteArray(l, n);
2135}
2136
2137////////////////////////////////////////////////////////////////////////////////
2138/// Write array of Float_t to buffer
2139
2141{
2142 XmlWriteArray(f, n);
2143}
2144
2145////////////////////////////////////////////////////////////////////////////////
2146/// Write array of Double_t to buffer
2147
2149{
2150 XmlWriteArray(d, n);
2151}
2152
2153/////////////////////////////////////////////////////////////////////////////////
2154/// Write array without size attribute
2155/// Also treat situation, when instead of one single array
2156/// chain of several elements should be produced
2157
2158template <typename T>
2160{
2162 if (n <= 0)
2163 return;
2165 PushStack(arrnode);
2166 XmlWriteArrayContent(arr, n);
2167 PopStack();
2168}
2169
2170////////////////////////////////////////////////////////////////////////////////
2171/// Write array of Bool_t to buffer
2172
2174{
2176}
2177
2178////////////////////////////////////////////////////////////////////////////////
2179/// Write array of Char_t to buffer
2180/// If array does not include any special characters,
2181/// it will be reproduced as CharStar node with string as attribute
2182
2184{
2185 Bool_t usedefault = (n == 0);
2186 const Char_t *buf = c;
2187 if (!usedefault)
2188 for (int i = 0; i < n; i++) {
2189 if (*buf < 27) {
2190 usedefault = kTRUE;
2191 break;
2192 }
2193 buf++;
2194 }
2195 if (usedefault) {
2197 } else {
2198 Char_t *buf2 = new Char_t[n + 1];
2199 memcpy(buf2, c, n);
2200 buf2[n] = 0;
2202 delete[] buf2;
2203 }
2204}
2205
2206////////////////////////////////////////////////////////////////////////////////
2207/// Write array of UChar_t to buffer
2208
2210{
2212}
2213
2214////////////////////////////////////////////////////////////////////////////////
2215/// Write array of Short_t to buffer
2216
2218{
2220}
2221
2222////////////////////////////////////////////////////////////////////////////////
2223/// Write array of UShort_t to buffer
2224
2226{
2228}
2229
2230////////////////////////////////////////////////////////////////////////////////
2231/// Write array of Int_t to buffer
2232
2234{
2235 XmlWriteFastArray(i, n);
2236}
2237
2238////////////////////////////////////////////////////////////////////////////////
2239/// Write array of UInt_t to buffer
2240
2242{
2243 XmlWriteFastArray(i, n);
2244}
2245
2246////////////////////////////////////////////////////////////////////////////////
2247/// Write array of Long_t to buffer
2248
2250{
2252}
2253
2254////////////////////////////////////////////////////////////////////////////////
2255/// Write array of ULong_t to buffer
2256
2258{
2260}
2261
2262////////////////////////////////////////////////////////////////////////////////
2263/// Write array of Long64_t to buffer
2264
2266{
2268}
2269
2270////////////////////////////////////////////////////////////////////////////////
2271/// Write array of ULong64_t to buffer
2272
2274{
2276}
2277
2278////////////////////////////////////////////////////////////////////////////////
2279/// Write array of Float_t to buffer
2280
2282{
2284}
2285
2286////////////////////////////////////////////////////////////////////////////////
2287/// Write array of Double_t to buffer
2288
2290{
2292}
2293
2294////////////////////////////////////////////////////////////////////////////////
2295/// Write array of n characters into the I/O buffer.
2296/// Used only by TLeafC, just dummy implementation here
2297
2299{
2300 WriteFastArray(c, n);
2301}
2302
2303////////////////////////////////////////////////////////////////////////////////
2304/// Write an array of object starting at the address 'start' and of length 'n'
2305/// the objects in the array are assumed to be of class 'cl'
2306
2307void TBufferXML::WriteFastArray(void *start, const TClass *cl, Int_t n, TMemberStreamer *streamer)
2308{
2309 if (streamer) {
2310 (*streamer)(*this, start, 0);
2311 return;
2312 }
2313
2314 char *obj = (char *)start;
2315 if (!n)
2316 n = 1;
2317 int size = cl->Size();
2318
2319 for (Int_t j = 0; j < n; j++, obj += size) {
2320 ((TClass *)cl)->Streamer(obj, *this);
2321 }
2322}
2323
2324////////////////////////////////////////////////////////////////////////////////
2325/// Write an array of object starting at the address '*start' and of length 'n'
2326/// the objects in the array are of class 'cl'
2327/// 'isPreAlloc' indicates whether the data member is marked with '->'
2328/// Return:
2329/// - 0: success
2330/// - 2: truncated success (i.e actual class is missing. Only ptrClass saved.)
2331
2332Int_t TBufferXML::WriteFastArray(void **start, const TClass *cl, Int_t n, Bool_t isPreAlloc, TMemberStreamer *streamer)
2333{
2334 // if isPreAlloc is true (data member has a ->) we can assume that the pointer
2335 // is never 0.
2336
2337 Bool_t oldStyle = kFALSE; // flag used to reproduce old-style I/O actions for kSTLp
2338
2339 if ((GetIOVersion() < 4) && !isPreAlloc) {
2340 TStreamerElement *elem = Stack()->fElem;
2341 if (elem && ((elem->GetType() == TStreamerInfo::kSTLp) ||
2343 oldStyle = kTRUE;
2344 }
2345
2346 if (streamer) {
2347 (*streamer)(*this, (void *)start, oldStyle ? n : 0);
2348 return 0;
2349 }
2350
2351 int strInfo = 0;
2352
2353 Int_t res = 0;
2354
2355 if (!isPreAlloc) {
2356
2357 for (Int_t j = 0; j < n; j++) {
2358 // must write StreamerInfo if pointer is null
2359 if (!strInfo && !start[j] && !oldStyle) {
2360 if (cl->Property() & kIsAbstract) {
2361 // Do not try to generate the StreamerInfo for an abstract class
2362 } else {
2363 TStreamerInfo *info = (TStreamerInfo *)((TClass *)cl)->GetStreamerInfo();
2364 ForceWriteInfo(info, kFALSE);
2365 }
2366 }
2367 strInfo = 2003;
2368 if (oldStyle)
2369 ((TClass *)cl)->Streamer(start[j], *this);
2370 else
2371 res |= WriteObjectAny(start[j], cl);
2372 }
2373
2374 } else {
2375 // case //-> in comment
2376
2377 for (Int_t j = 0; j < n; j++) {
2378 if (!start[j])
2379 start[j] = ((TClass *)cl)->New();
2380 ((TClass *)cl)->Streamer(start[j], *this);
2381 }
2382 }
2383 return res;
2384}
2385
2386////////////////////////////////////////////////////////////////////////////////
2387/// Stream object to/from buffer
2388
2389void TBufferXML::StreamObject(void *obj, const TClass *cl, const TClass * /* onfileClass */)
2390{
2391 if (GetIOVersion() < 4) {
2392 TStreamerElement *elem = Stack()->fElem;
2393 if (elem && (elem->GetType() == TStreamerInfo::kTObject)) {
2394 ((TObject *)obj)->TObject::Streamer(*this);
2395 return;
2396 } else if (elem && (elem->GetType() == TStreamerInfo::kTNamed)) {
2397 ((TNamed *)obj)->TNamed::Streamer(*this);
2398 return;
2399 }
2400 }
2401
2403 if (gDebug > 1)
2404 Info("StreamObject", "Class: %s", (cl ? cl->GetName() : "none"));
2405 if (IsReading())
2406 XmlReadObject(obj);
2407 else
2408 XmlWriteObject(obj, cl, kTRUE);
2409}
2410
2411////////////////////////////////////////////////////////////////////////////////
2412/// Reads Bool_t value from buffer
2413
2415{
2417 XmlReadBasic(b);
2418}
2419
2420////////////////////////////////////////////////////////////////////////////////
2421/// Reads Char_t value from buffer
2422
2424{
2426 XmlReadBasic(c);
2427}
2428
2429////////////////////////////////////////////////////////////////////////////////
2430/// Reads UChar_t value from buffer
2431
2433{
2435 XmlReadBasic(c);
2436}
2437
2438////////////////////////////////////////////////////////////////////////////////
2439/// Reads Short_t value from buffer
2440
2442{
2444 XmlReadBasic(h);
2445}
2446
2447////////////////////////////////////////////////////////////////////////////////
2448/// Reads UShort_t value from buffer
2449
2451{
2453 XmlReadBasic(h);
2454}
2455
2456////////////////////////////////////////////////////////////////////////////////
2457/// Reads Int_t value from buffer
2458
2460{
2462 XmlReadBasic(i);
2463}
2464
2465////////////////////////////////////////////////////////////////////////////////
2466/// Reads UInt_t value from buffer
2467
2469{
2471 XmlReadBasic(i);
2472}
2473
2474////////////////////////////////////////////////////////////////////////////////
2475/// Reads Long_t value from buffer
2476
2478{
2480 XmlReadBasic(l);
2481}
2482
2483////////////////////////////////////////////////////////////////////////////////
2484/// Reads ULong_t value from buffer
2485
2487{
2489 XmlReadBasic(l);
2490}
2491
2492////////////////////////////////////////////////////////////////////////////////
2493/// Reads Long64_t value from buffer
2494
2496{
2498 XmlReadBasic(l);
2499}
2500
2501////////////////////////////////////////////////////////////////////////////////
2502/// Reads ULong64_t value from buffer
2503
2505{
2507 XmlReadBasic(l);
2508}
2509
2510////////////////////////////////////////////////////////////////////////////////
2511/// Reads Float_t value from buffer
2512
2514{
2516 XmlReadBasic(f);
2517}
2518
2519////////////////////////////////////////////////////////////////////////////////
2520/// Reads Double_t value from buffer
2521
2523{
2525 XmlReadBasic(d);
2526}
2527
2528////////////////////////////////////////////////////////////////////////////////
2529/// Reads array of characters from buffer
2530
2532{
2534 const char *buf;
2535 if ((buf = XmlReadValue(xmlio::CharStar)))
2536 strcpy(c, buf); // NOLINT unfortunately, size of target buffer cannot be controlled here
2537}
2538
2539////////////////////////////////////////////////////////////////////////////////
2540/// Reads a TString
2541
2543{
2544 if (GetIOVersion() < 3) {
2545 // original TBufferFile method can not be used, while used TString methods are private
2546 // try to reimplement close to the original
2547 Int_t nbig;
2548 UChar_t nwh;
2549 *this >> nwh;
2550 if (nwh == 0) {
2551 s.Resize(0);
2552 } else {
2553 if (nwh == 255)
2554 *this >> nbig;
2555 else
2556 nbig = nwh;
2557
2558 char *data = new char[nbig+1];
2559 data[nbig] = 0;
2560 ReadFastArray(data, nbig);
2561 s = data;
2562 delete[] data;
2563 }
2564 } else {
2566 const char *buf = XmlReadValue(xmlio::String);
2567 if (buf)
2568 s = buf;
2569 }
2570}
2571
2572////////////////////////////////////////////////////////////////////////////////
2573/// Reads a std::string
2574
2575void TBufferXML::ReadStdString(std::string *obj)
2576{
2577 if (GetIOVersion() < 3) {
2578 if (!obj) {
2579 Error("ReadStdString", "The std::string address is nullptr but should not");
2580 return;
2581 }
2582 Int_t nbig;
2583 UChar_t nwh;
2584 *this >> nwh;
2585 if (nwh == 0) {
2586 obj->clear();
2587 } else {
2588 if (obj->size()) {
2589 // Insure that the underlying data storage is not shared
2590 (*obj)[0] = '\0';
2591 }
2592 if (nwh == 255) {
2593 *this >> nbig;
2594 obj->resize(nbig, '\0');
2595 ReadFastArray((char *)obj->data(), nbig);
2596 } else {
2597 obj->resize(nwh, '\0');
2598 ReadFastArray((char *)obj->data(), nwh);
2599 }
2600 }
2601 } else {
2603 const char *buf = XmlReadValue(xmlio::String);
2604 if (buf && obj)
2605 *obj = buf;
2606 }
2607}
2608
2609////////////////////////////////////////////////////////////////////////////////
2610/// Read a char* string
2611
2613{
2614 delete[] s;
2615 s = nullptr;
2616
2617 Int_t nch;
2618 *this >> nch;
2619 if (nch > 0) {
2620 s = new char[nch + 1];
2621 ReadFastArray(s, nch);
2622 s[nch] = 0;
2623 }
2624}
2625
2626////////////////////////////////////////////////////////////////////////////////
2627/// Writes Bool_t value to buffer
2628
2630{
2633}
2634
2635////////////////////////////////////////////////////////////////////////////////
2636/// Writes Char_t value to buffer
2637
2639{
2642}
2643
2644////////////////////////////////////////////////////////////////////////////////
2645/// Writes UChar_t value to buffer
2646
2648{
2651}
2652
2653////////////////////////////////////////////////////////////////////////////////
2654/// Writes Short_t value to buffer
2655
2657{
2660}
2661
2662////////////////////////////////////////////////////////////////////////////////
2663/// Writes UShort_t value to buffer
2664
2666{
2669}
2670
2671////////////////////////////////////////////////////////////////////////////////
2672/// Writes Int_t value to buffer
2673
2675{
2677 XmlWriteBasic(i);
2678}
2679
2680////////////////////////////////////////////////////////////////////////////////
2681/// Writes UInt_t value to buffer
2682
2684{
2686 XmlWriteBasic(i);
2687}
2688
2689////////////////////////////////////////////////////////////////////////////////
2690/// Writes Long_t value to buffer
2691
2693{
2696}
2697
2698////////////////////////////////////////////////////////////////////////////////
2699/// Writes ULong_t value to buffer
2700
2702{
2705}
2706
2707////////////////////////////////////////////////////////////////////////////////
2708/// Writes Long64_t value to buffer
2709
2711{
2714}
2715
2716////////////////////////////////////////////////////////////////////////////////
2717/// Writes ULong64_t value to buffer
2718
2720{
2723}
2724
2725////////////////////////////////////////////////////////////////////////////////
2726/// Writes Float_t value to buffer
2727
2729{
2732}
2733
2734////////////////////////////////////////////////////////////////////////////////
2735/// Writes Double_t value to buffer
2736
2738{
2741}
2742
2743////////////////////////////////////////////////////////////////////////////////
2744/// Writes array of characters to buffer
2745
2747{
2750}
2751
2752////////////////////////////////////////////////////////////////////////////////
2753/// Writes a TString
2754
2756{
2757 if (GetIOVersion() < 3) {
2758 // original TBufferFile method, keep for compatibility
2759 Int_t nbig = s.Length();
2760 UChar_t nwh;
2761 if (nbig > 254) {
2762 nwh = 255;
2763 *this << nwh;
2764 *this << nbig;
2765 } else {
2766 nwh = UChar_t(nbig);
2767 *this << nwh;
2768 }
2769 const char *data = s.Data();
2770 WriteFastArray(data, nbig);
2771 } else {
2774 }
2775}
2776
2777////////////////////////////////////////////////////////////////////////////////
2778/// Writes a std::string
2779
2780void TBufferXML::WriteStdString(const std::string *obj)
2781{
2782 if (GetIOVersion() < 3) {
2783 if (!obj) {
2784 *this << (UChar_t)0;
2785 WriteFastArray("", 0);
2786 return;
2787 }
2788
2789 UChar_t nwh;
2790 Int_t nbig = obj->length();
2791 if (nbig > 254) {
2792 nwh = 255;
2793 *this << nwh;
2794 *this << nbig;
2795 } else {
2796 nwh = UChar_t(nbig);
2797 *this << nwh;
2798 }
2799 WriteFastArray(obj->data(), nbig);
2800 } else {
2802 XmlWriteValue(obj ? obj->c_str() : "", xmlio::String);
2803 }
2804}
2805
2806////////////////////////////////////////////////////////////////////////////////
2807/// Write a char* string
2808
2810{
2811 Int_t nch = 0;
2812 if (s) {
2813 nch = strlen(s);
2814 *this << nch;
2815 WriteFastArray(s, nch);
2816 } else {
2817 *this << nch;
2818 }
2819}
2820
2821////////////////////////////////////////////////////////////////////////////////
2822/// Converts Char_t to string and add xml node to buffer
2823
2825{
2826 char buf[50];
2827 snprintf(buf, sizeof(buf), "%d", value);
2828 return XmlWriteValue(buf, xmlio::Char);
2829}
2830
2831////////////////////////////////////////////////////////////////////////////////
2832/// Converts Short_t to string and add xml node to buffer
2833
2835{
2836 char buf[50];
2837 snprintf(buf, sizeof(buf), "%hd", value);
2838 return XmlWriteValue(buf, xmlio::Short);
2839}
2840
2841////////////////////////////////////////////////////////////////////////////////
2842/// Converts Int_t to string and add xml node to buffer
2843
2845{
2846 char buf[50];
2847 snprintf(buf, sizeof(buf), "%d", value);
2848 return XmlWriteValue(buf, xmlio::Int);
2849}
2850
2851////////////////////////////////////////////////////////////////////////////////
2852/// Converts Long_t to string and add xml node to buffer
2853
2855{
2856 char buf[50];
2857 snprintf(buf, sizeof(buf), "%ld", value);
2858 return XmlWriteValue(buf, xmlio::Long);
2859}
2860
2861////////////////////////////////////////////////////////////////////////////////
2862/// Converts Long64_t to string and add xml node to buffer
2863
2865{
2866 std::string buf = std::to_string(value);
2867 return XmlWriteValue(buf.c_str(), xmlio::Long64);
2868}
2869
2870////////////////////////////////////////////////////////////////////////////////
2871/// Converts Float_t to string and add xml node to buffer
2872
2874{
2875 char buf[200];
2876 ConvertFloat(value, buf, sizeof(buf), kTRUE);
2877 return XmlWriteValue(buf, xmlio::Float);
2878}
2879
2880////////////////////////////////////////////////////////////////////////////////
2881/// Converts Double_t to string and add xml node to buffer
2882
2884{
2885 char buf[1000];
2886 ConvertDouble(value, buf, sizeof(buf), kTRUE);
2887 return XmlWriteValue(buf, xmlio::Double);
2888}
2889
2890////////////////////////////////////////////////////////////////////////////////
2891/// Converts Bool_t to string and add xml node to buffer
2892
2894{
2896}
2897
2898////////////////////////////////////////////////////////////////////////////////
2899/// Converts UChar_t to string and add xml node to buffer
2900
2902{
2903 char buf[50];
2904 snprintf(buf, sizeof(buf), "%u", value);
2905 return XmlWriteValue(buf, xmlio::UChar);
2906}
2907
2908////////////////////////////////////////////////////////////////////////////////
2909/// Converts UShort_t to string and add xml node to buffer
2910
2912{
2913 char buf[50];
2914 snprintf(buf, sizeof(buf), "%hu", value);
2915 return XmlWriteValue(buf, xmlio::UShort);
2916}
2917
2918////////////////////////////////////////////////////////////////////////////////
2919/// Converts UInt_t to string and add xml node to buffer
2920
2922{
2923 char buf[50];
2924 snprintf(buf, sizeof(buf), "%u", value);
2925 return XmlWriteValue(buf, xmlio::UInt);
2926}
2927
2928////////////////////////////////////////////////////////////////////////////////
2929/// Converts ULong_t to string and add xml node to buffer
2930
2932{
2933 char buf[50];
2934 snprintf(buf, sizeof(buf), "%lu", value);
2935 return XmlWriteValue(buf, xmlio::ULong);
2936}
2937
2938////////////////////////////////////////////////////////////////////////////////
2939/// Converts ULong64_t to string and add xml node to buffer
2940
2942{
2943 std::string buf = std::to_string(value);
2944 return XmlWriteValue(buf.c_str(), xmlio::ULong64);
2945}
2946
2947////////////////////////////////////////////////////////////////////////////////
2948/// Create xml node with specified name and adds it to stack node
2949
2950XMLNodePointer_t TBufferXML::XmlWriteValue(const char *value, const char *name)
2951{
2952 XMLNodePointer_t node = nullptr;
2953
2954 if (fCanUseCompact)
2955 node = StackNode();
2956 else
2957 node = CreateItemNode(name);
2958
2959 fXML->NewAttr(node, nullptr, xmlio::v, value);
2960
2962
2963 return node;
2964}
2965
2966////////////////////////////////////////////////////////////////////////////////
2967/// Reads string from current xml node and convert it to Char_t value
2968
2970{
2971 const char *res = XmlReadValue(xmlio::Char);
2972 if (res) {
2973 int n;
2974 sscanf(res, "%d", &n);
2975 value = n;
2976 } else
2977 value = 0;
2978}
2979
2980////////////////////////////////////////////////////////////////////////////////
2981/// Reads string from current xml node and convert it to Short_t value
2982
2984{
2985 const char *res = XmlReadValue(xmlio::Short);
2986 if (res)
2987 sscanf(res, "%hd", &value);
2988 else
2989 value = 0;
2990}
2991
2992////////////////////////////////////////////////////////////////////////////////
2993/// Reads string from current xml node and convert it to Int_t value
2994
2996{
2997 const char *res = XmlReadValue(xmlio::Int);
2998 if (res)
2999 sscanf(res, "%d", &value);
3000 else
3001 value = 0;
3002}
3003
3004////////////////////////////////////////////////////////////////////////////////
3005/// Reads string from current xml node and convert it to Long_t value
3006
3008{
3009 const char *res = XmlReadValue(xmlio::Long);
3010 if (res)
3011 sscanf(res, "%ld", &value);
3012 else
3013 value = 0;
3014}
3015
3016////////////////////////////////////////////////////////////////////////////////
3017/// Reads string from current xml node and convert it to Long64_t value
3018
3020{
3021 const char *res = XmlReadValue(xmlio::Long64);
3022 if (res)
3023 value = (Long64_t)std::stoll(res);
3024 else
3025 value = 0;
3026}
3027
3028////////////////////////////////////////////////////////////////////////////////
3029/// Reads string from current xml node and convert it to Float_t value
3030
3032{
3033 const char *res = XmlReadValue(xmlio::Float);
3034 if (res)
3035 sscanf(res, "%f", &value);
3036 else
3037 value = 0.;
3038}
3039
3040////////////////////////////////////////////////////////////////////////////////
3041/// Reads string from current xml node and convert it to Double_t value
3042
3044{
3045 const char *res = XmlReadValue(xmlio::Double);
3046 if (res)
3047 sscanf(res, "%lf", &value);
3048 else
3049 value = 0.;
3050}
3051
3052////////////////////////////////////////////////////////////////////////////////
3053/// Reads string from current xml node and convert it to Bool_t value
3054
3056{
3057 const char *res = XmlReadValue(xmlio::Bool);
3058 if (res)
3059 value = (strcmp(res, xmlio::True) == 0);
3060 else
3061 value = kFALSE;
3062}
3063
3064////////////////////////////////////////////////////////////////////////////////
3065/// Reads string from current xml node and convert it to UChar_t value
3066
3068{
3069 const char *res = XmlReadValue(xmlio::UChar);
3070 if (res) {
3071 unsigned int n;
3072 sscanf(res, "%ud", &n);
3073 value = n;
3074 } else
3075 value = 0;
3076}
3077
3078////////////////////////////////////////////////////////////////////////////////
3079/// Reads string from current xml node and convert it to UShort_t value
3080
3082{
3083 const char *res = XmlReadValue(xmlio::UShort);
3084 if (res)
3085 sscanf(res, "%hud", &value);
3086 else
3087 value = 0;
3088}
3089
3090////////////////////////////////////////////////////////////////////////////////
3091/// Reads string from current xml node and convert it to UInt_t value
3092
3094{
3095 const char *res = XmlReadValue(xmlio::UInt);
3096 if (res)
3097 sscanf(res, "%u", &value);
3098 else
3099 value = 0;
3100}
3101
3102////////////////////////////////////////////////////////////////////////////////
3103/// Reads string from current xml node and convert it to ULong_t value
3104
3106{
3107 const char *res = XmlReadValue(xmlio::ULong);
3108 if (res)
3109 sscanf(res, "%lu", &value);
3110 else
3111 value = 0;
3112}
3113
3114////////////////////////////////////////////////////////////////////////////////
3115/// Reads string from current xml node and convert it to ULong64_t value
3116
3118{
3119 const char *res = XmlReadValue(xmlio::ULong64);
3120 if (res)
3121 value = (ULong64_t)std::stoull(res);
3122 else
3123 value = 0;
3124}
3125
3126////////////////////////////////////////////////////////////////////////////////
3127/// read string value from current stack node
3128
3129const char *TBufferXML::XmlReadValue(const char *name)
3130{
3131 if (fErrorFlag > 0)
3132 return 0;
3133
3134 Bool_t trysimple = fCanUseCompact;
3136
3137 if (trysimple) {
3138 if (fXML->HasAttr(Stack(1)->fNode, xmlio::v))
3139 fValueBuf = fXML->GetAttr(Stack(1)->fNode, xmlio::v);
3140 else
3141 trysimple = kFALSE;
3142 }
3143
3144 if (!trysimple) {
3145 if (!VerifyItemNode(name, "XmlReadValue"))
3146 return 0;
3148 }
3149
3150 if (gDebug > 4)
3151 Info("XmlReadValue", " Name = %s value = %s", name, fValueBuf.Data());
3152
3153 if (!trysimple)
3154 ShiftStack("readvalue");
3155
3156 return fValueBuf.Data();
3157}
3158
3159////////////////////////////////////////////////////////////////////////////////
3160/// Return current streamer info element
3161
3163{
3164 return Stack()->fInfo;
3165}
#define R__ALWAYS_INLINE
Definition RConfig.hxx:564
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
unsigned short UShort_t
Definition RtypesCore.h:40
long Longptr_t
Definition RtypesCore.h:82
short Version_t
Definition RtypesCore.h:65
unsigned char UChar_t
Definition RtypesCore.h:38
char Char_t
Definition RtypesCore.h:37
const Bool_t kFALSE
Definition RtypesCore.h:101
unsigned long ULong_t
Definition RtypesCore.h:55
long Long_t
Definition RtypesCore.h:54
short Short_t
Definition RtypesCore.h:39
double Double_t
Definition RtypesCore.h:59
long long Long64_t
Definition RtypesCore.h:80
unsigned long long ULong64_t
Definition RtypesCore.h:81
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
@ kIsAbstract
Definition TDictionary.h:71
char name[80]
Definition TGX11.cxx:110
Int_t gDebug
Definition TROOT.cxx:592
#define gROOT
Definition TROOT.h:404
void R__unzip(Int_t *nin, UChar_t *bufin, Int_t *lout, char *bufout, Int_t *nout)
int R__unzip_header(Int_t *nin, UChar_t *bufin, Int_t *lout)
void * XMLNodePointer_t
Definition TXMLEngine.h:17
void * XMLNsPointer_t
Definition TXMLEngine.h:18
#define snprintf
Definition civetweb.c:1540
void InitMap() override
Create the fMap container and initialize them with the null object.
void ForceWriteInfo(TVirtualStreamerInfo *info, Bool_t force) override
force writing the TStreamerInfo to the file
TExMap * fMap
Map containing object,offset pairs for reading/writing.
Definition TBufferIO.h:39
void MapObject(const TObject *obj, UInt_t offset=1) override
Add object to the fMap container.
Long64_t GetObjectTag(const void *obj)
Returns tag for specified object from objects map (if exists) Returns 0 if object not included into o...
static R__ALWAYS_INLINE ULong_t Void_Hash(const void *ptr)
Return hash value for provided object.
Definition TBufferIO.h:53
void GetMappedObject(UInt_t tag, void *&ptr, TClass *&ClassPtr) const override
Retrieve the object stored in the buffer's object map at 'tag' Set ptr and ClassPtr respectively to t...
Int_t WriteObjectAny(const void *obj, const TClass *ptrClass, Bool_t cacheReuse=kTRUE) override
Write object to I/O buffer.
Base class for text-based streamers like TBufferJSON or TBufferXML Special actions list will use meth...
Definition TBufferText.h:20
static const char * ConvertFloat(Float_t v, char *buf, unsigned len, Bool_t not_optimize=kFALSE)
convert float to string with configured format
static const char * ConvertDouble(Double_t v, char *buf, unsigned len, Bool_t not_optimize=kFALSE)
convert float to string with configured format
Class for serializing/deserializing object to/from xml.
Definition TBufferXML.h:33
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...
void WriteLong(Long_t l) final
Writes Long_t value to buffer.
void WriteFastArray(const Bool_t *b, Int_t n) final
Write array of Bool_t to buffer.
void SetXML(TXMLEngine *xml)
Definition TBufferXML.h:229
Int_t GetCompressionSettings() const
Definition TBufferXML.h:348
TXMLStackObj * PushStack(XMLNodePointer_t current, Bool_t simple=kFALSE)
Add new level to xml stack.
Bool_t VerifyAttr(XMLNodePointer_t node, const char *name, const char *value, const char *errinfo=nullptr)
Checks, that attribute of specified name exists and has specified value.
void WorkWithClass(TStreamerInfo *info, const TClass *cl=nullptr)
Prepares buffer to stream data of specified class.
Bool_t VerifyStackNode(const char *name, const char *errinfo=nullptr)
Check, if stack node has specified name.
Int_t GetCompressionAlgorithm() const
Definition TBufferXML.h:336
void ReadUShort(UShort_t &s) final
Reads UShort_t value from buffer.
void ReadUInt(UInt_t &i) final
Reads UInt_t value from buffer.
Int_t fCompressLevel
! Compression level and algorithm
Definition TBufferXML.h:329
void SetCompressionSettings(Int_t settings=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault)
Used to specify the compression level and algorithm.
TString fValueBuf
! Current value buffer
Definition TBufferXML.h:325
void StreamObject(void *obj, const TClass *cl, const TClass *onFileClass=nullptr) final
Stream object to/from buffer.
Int_t GetCompressionLevel() const
Definition TBufferXML.h:342
void WriteStdString(const std::string *s) final
Writes a std::string.
void ReadDouble(Double_t &d) final
Reads Double_t value from buffer.
TBufferXML(TBuffer::EMode mode)
Creates buffer object to serialize/deserialize data to/from xml.
Bool_t VerifyNode(XMLNodePointer_t node, const char *name, const char *errinfo=nullptr)
Check if node has specified name.
void ReadFastArrayString(Char_t *c, Int_t n) final
Read array of n characters from the I/O buffer.
void WriteLong64(Long64_t l) final
Writes Long64_t value to buffer.
void WriteChar(Char_t c) final
Writes Char_t value to buffer.
TXMLEngine * fXML
! instance of TXMLEngine for working with XML structures
Definition TBufferXML.h:322
std::deque< std::unique_ptr< TXMLStackObj > > fStack
! Stack of processed objects
Definition TBufferXML.h:323
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 ...
void WriteFloat(Float_t f) final
Writes Float_t value to buffer.
void ReadTString(TString &s) final
Reads a TString.
void WriteTString(const TString &s) final
Writes a TString.
void IncrementLevel(TVirtualStreamerInfo *) final
Function is called from TStreamerInfo WriteBuffer and ReadBuffer functions and indent new level in xm...
void SetIOVersion(Int_t v)
Definition TBufferXML.h:66
void WriteArray(const Bool_t *b, Int_t n) final
Write array of Bool_t to buffer.
XMLNodePointer_t XmlWriteAny(const void *obj, const TClass *cl)
Convert object of any class to xml structures Return pointer on top xml element.
void XmlReadBasic(Char_t &value)
Reads string from current xml node and convert it to Char_t value.
void ReadULong(ULong_t &l) final
Reads ULong_t value from buffer.
Int_t GetIOVersion() const
Definition TBufferXML.h:65
R__ALWAYS_INLINE void XmlWriteArrayContent(const T *arr, Int_t arrsize)
XMLNodePointer_t XmlWriteValue(const char *value, const char *name)
Create xml node with specified name and adds it to stack node.
void SkipObjectAny() final
Skip any kind of object from buffer Actually skip only one node on current level of xml structure.
Bool_t VerifyStackAttr(const char *name, const char *value, const char *errinfo=nullptr)
Checks stack attribute.
void ReadFloat(Float_t &f) final
Reads Float_t value from buffer.
void ReadULong64(ULong64_t &l) final
Reads ULong64_t value from buffer.
XMLNodePointer_t XmlWriteBasic(Char_t value)
Converts Char_t to string and add xml node to buffer.
void ShiftStack(const char *info=nullptr)
Shift stack node to next.
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 WriteFastArrayString(const Char_t *c, Int_t n) final
Write array of n characters into the I/O buffer.
void ReadShort(Short_t &s) final
Reads Short_t value from buffer.
void BeforeIOoperation()
Function is called before any IO operation of TBuffer Now is used to store version value if no proper...
TClass * ReadClass(const TClass *cl=nullptr, UInt_t *objTag=nullptr) final
Function to read class from buffer, used in old-style streamers.
Int_t fErrorFlag
! Error flag
Definition TBufferXML.h:326
void ReadChar(Char_t &c) final
Reads Char_t value from buffer.
R__ALWAYS_INLINE void XmlReadArrayContent(T *arr, Int_t arrsize)
Template method to read array content.
void ClassEnd(const TClass *) final
Should be called at the end of custom streamer See TBufferXML::ClassBegin for more details.
void ReadLong64(Long64_t &l) final
Reads Long64_t value from buffer.
Version_t fVersionBuf
! Current version buffer
Definition TBufferXML.h:324
void ClassBegin(const TClass *, Version_t=-1) final
Should be called at the beginning of custom class streamer.
void XmlReadBlock(XMLNodePointer_t node)
Read binary block of data from xml.
void * XmlReadObject(void *obj, TClass **cl=nullptr)
Read object from the buffer.
void ReadCharP(Char_t *c) final
Reads array of characters from buffer.
Bool_t fCanUseCompact
! Flag indicate that basic type (like Int_t) can be placed in the same tag
Definition TBufferXML.h:327
R__ALWAYS_INLINE void XmlWriteArray(const T *arr, Int_t arrsize)
Write array, including it size Content may be compressed.
void ReadStdString(std::string *s) final
Reads a std::string.
void ReadBool(Bool_t &b) final
Reads Bool_t value from buffer.
void WriteUShort(UShort_t s) final
Writes UShort_t value to buffer.
void WriteClass(const TClass *cl) final
Function to write class into buffer, used in old-style streamers.
void WriteCharStar(char *s) final
Write a char* string.
const char * XmlReadValue(const char *name)
read string value from current stack node
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...
void ReadFastArray(Bool_t *b, Int_t n) final
Read array of Bool_t from buffer.
void DecrementLevel(TVirtualStreamerInfo *) final
Function is called from TStreamerInfo WriteBuffer and ReadBuffer functions and decrease level in xml ...
R__ALWAYS_INLINE void XmlWriteFastArray(const T *arr, Int_t n)
Write array without size attribute Also treat situation, when instead of one single array chain of se...
void PerformPostProcessing()
Function is converts TObject and TString structures to more compact representation.
void SetStreamerElementNumber(TStreamerElement *elem, Int_t comp_type) final
Function is called from TStreamerInfo WriteBuffer and ReadBuffer functions and add/verify next elemen...
void WriteCharP(const Char_t *c) final
Writes array of characters to buffer.
static TObject * ConvertFromXML(const char *str, Bool_t GenericLayout=kFALSE, Bool_t UseNamespaces=kFALSE)
Read object from XML, produced by ConvertToXML() method.
void ClassMember(const char *name, const char *typeName=nullptr, Int_t arrsize1=-1, Int_t arrsize2=-1) final
Method indicates name and typename of class member, which should be now streamed in custom streamer.
Int_t ReadStaticArray(Bool_t *b) final
Read array of Bool_t from buffer.
void * XmlReadAny(XMLNodePointer_t node, void *obj, TClass **cl)
Recreate object from xml structure.
XMLNodePointer_t StackNode()
Return pointer on current xml node.
void WriteObjectClass(const void *actualObjStart, const TClass *actualClass, Bool_t cacheReuse) final
Write object to buffer. Only used from TBuffer.
R__ALWAYS_INLINE void XmlReadFastArray(T *arr, Int_t n)
Template method to read content of array, which not include size of array Also treated situation,...
Bool_t VerifyItemNode(const char *name, const char *errinfo=nullptr)
Checks, if stack node is item and has specified name.
static void * ConvertFromXMLChecked(const char *xml, const TClass *expectedClass, Bool_t GenericLayout=kFALSE, Bool_t UseNamespaces=kFALSE)
Convert from XML and check if object derived from specified class When possible, cast to given class.
void * ReadObjectAny(const TClass *clCast) final
Read object from buffer. Only used from TBuffer.
void ReadLong(Long_t &l) final
Reads Long_t value from buffer.
void ExtractReference(XMLNodePointer_t node, const void *ptr, const TClass *cl)
Analyze if node has "ref" attribute and register it to object map.
void SetCompressionLevel(Int_t level=ROOT::RCompressionSetting::ELevel::kUseMin)
See comments for function SetCompressionSettings.
Version_t ReadVersion(UInt_t *start=nullptr, UInt_t *bcnt=nullptr, const TClass *cl=nullptr) final
Read version value from buffer.
void ReadInt(Int_t &i) final
Reads Int_t value from buffer.
Int_t ReadArray(Bool_t *&b) final
Read array of Bool_t from buffer.
void SetCompressionAlgorithm(Int_t algorithm=ROOT::RCompressionSetting::EAlgorithm::kUseGlobal)
See comments for function SetCompressionSettings.
void WriteDouble(Double_t d) final
Writes Double_t value to buffer.
void CheckVersionBuf()
Checks buffer, filled by WriteVersion if next data is arriving, version should be stored in buffer.
XMLNodePointer_t XmlWriteObject(const void *obj, const TClass *objClass, Bool_t cacheReuse)
Write object to buffer If object was written before, only pointer will be stored Return pointer to to...
TXMLFile * XmlFile()
Returns pointer to TXMLFile object.
TVirtualStreamerInfo * GetInfo() final
Return current streamer info element.
void WriteUInt(UInt_t i) final
Writes UInt_t value to buffer.
void WorkWithElement(TStreamerElement *elem, Int_t comp_type)
This function is a part of SetStreamerElementNumber method.
virtual ~TBufferXML()
Destroy xml buffer.
void WriteBool(Bool_t b) final
Writes Bool_t value to buffer.
void XmlWriteBlock(XMLNodePointer_t node)
Write binary data block from buffer to xml.
TClass * fExpectedBaseClass
! Pointer to class, which should be stored as parent of current
Definition TBufferXML.h:328
void WriteShort(Short_t s) final
Writes Short_t value to buffer.
void WriteULong64(ULong64_t l) final
Writes ULong64_t value to buffer.
TXMLStackObj * Stack(UInt_t depth=0)
Definition TBufferXML.h:242
XMLNodePointer_t CreateItemNode(const char *name)
Create item node of specified name.
void WriteULong(ULong_t l) final
Writes ULong_t value to buffer.
TXMLStackObj * PopStack()
Remove one level from xml stack.
void CreateElemNode(const TStreamerElement *elem)
Create xml node correspondent to TStreamerElement object.
static void * ConvertFromXMLAny(const char *str, TClass **cl=nullptr, Bool_t GenericLayout=kFALSE, Bool_t UseNamespaces=kFALSE)
Read object of any class from XML, produced by ConvertToXML() method.
Bool_t VerifyElemNode(const TStreamerElement *elem)
Checks if stack node correspond to TStreamerElement object.
R__ALWAYS_INLINE Int_t XmlReadArray(T *&arr, bool is_static=false)
Template method to read array with size attribute If necessary, array is created.
void ReadCharStar(char *&s) final
Read a char* string.
void ReadUChar(UChar_t &c) final
Reads UChar_t value from buffer.
UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE) final
Copies class version to buffer, but not writes it to xml Version will be written with next I/O operat...
void WriteUChar(UChar_t c) final
Writes UChar_t value to buffer.
void WriteInt(Int_t i) final
Writes Int_t value to buffer.
TObject * GetParent() const
Return pointer to parent of this buffer.
Definition TBuffer.cxx:262
void Expand(Int_t newsize, Bool_t copy=kTRUE)
Expand (or shrink) the I/O buffer to newsize bytes.
Definition TBuffer.cxx:223
Int_t BufferSize() const
Definition TBuffer.h:98
@ kWrite
Definition TBuffer.h:73
@ kRead
Definition TBuffer.h:73
Bool_t IsWriting() const
Definition TBuffer.h:87
Bool_t IsReading() const
Definition TBuffer.h:86
Int_t Length() const
Definition TBuffer.h:100
char * Buffer() const
Definition TBuffer.h:96
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:80
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition TClass.cxx:4964
void Destructor(void *obj, Bool_t dtorOnly=kFALSE)
Explicitly call destructor for object.
Definition TClass.cxx:5386
Int_t Size() const
Return size of object of this class.
Definition TClass.cxx:5690
Bool_t IsTObject() const
Return kTRUE is the class inherits from TObject.
Definition TClass.cxx:5924
Int_t GetBaseClassOffset(const TClass *toBase, void *address=0, bool isDerivedObject=true)
Definition TClass.cxx:2789
Long_t Property() const
Returns the properties of the TClass as a bit field stored as a Long_t value.
Definition TClass.cxx:6072
void Streamer(void *obj, TBuffer &b, const TClass *onfile_class=0) const
Definition TClass.h:602
Version_t GetClassVersion() const
Definition TClass.h:417
TClass * GetActualClass(const void *object) const
Return a pointer the the real class of the object.
Definition TClass.cxx:2606
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:2966
Basic data type descriptor (datatype information is obtained from CINT).
Definition TDataType.h:44
Int_t GetType() const
Definition TDataType.h:68
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 SetOnFileClass(const TClass *cl)
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Int_t IndexOf(const TObject *obj) const
Mother of all ROOT objects.
Definition TObject.h:41
@ kIsOnHeap
object is on heap
Definition TObject.h:81
@ kNotDeleted
object has not been deleted
Definition TObject.h:82
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
virtual ~TObject()
TObject destructor.
Definition TObject.cxx:151
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:937
Int_t GetType() const
virtual TClass * GetClassPointer() const
Returns a pointer to the TClass of this element.
virtual void SetArrayDim(Int_t dim)
Set number of array dimensions.
virtual void SetMaxIndex(Int_t dim, Int_t max)
set maximum index for array with dimension dim
Describes a persistent version of a class.
TObjArray * GetElements() const
TClass * GetClass() const
Basic string class.
Definition TString.h:136
Ssiz_t Length() const
Definition TString.h:410
Int_t Atoi() const
Return integer value of string.
Definition TString.cxx:1946
const char * Data() const
Definition TString.h:369
void Resize(Ssiz_t n)
Resize the string. Truncate or add blanks as necessary.
Definition TString.cxx:1120
Abstract Interface class describing Streamer information for one class.
static Bool_t CanDelete()
static function returning true if ReadBuffer can delete object
virtual TClass * GetClass() const =0
XMLNodePointer_t NewChild(XMLNodePointer_t parent, XMLNsPointer_t ns, const char *name, const char *content=nullptr)
create new child element for parent node
XMLNodePointer_t GetChild(XMLNodePointer_t xmlnode, Bool_t realnode=kTRUE)
returns first child of xmlnode
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
void SaveSingleNode(XMLNodePointer_t xmlnode, TString *res, Int_t layout=1)
convert single xmlnode (and its child node) to string if layout<=0, no any spaces or newlines will be...
XMLAttrPointer_t NewIntAttr(XMLNodePointer_t xmlnode, const char *name, Int_t value)
create node attribute with integer value
Bool_t HasAttr(XMLNodePointer_t xmlnode, const char *name)
checks if node has attribute of specified name
XMLNodePointer_t ReadSingleNode(const char *src)
read single xmlnode from provided string
const char * GetNodeContent(XMLNodePointer_t xmlnode)
get contents (if any) of xmlnode
const char * GetNodeName(XMLNodePointer_t xmlnode)
returns name of xmlnode
void FreeAttr(XMLNodePointer_t xmlnode, const char *name)
remove attribute from xmlnode
const char * GetAttr(XMLNodePointer_t xmlnode, const char *name)
returns value of attribute for xmlnode
XMLNsPointer_t NewNS(XMLNodePointer_t xmlnode, const char *reference, const char *name=nullptr)
create namespace attribute for xmlnode.
Int_t GetIntAttr(XMLNodePointer_t node, const char *name)
returns value of attribute as integer
void UnlinkFreeNode(XMLNodePointer_t xmlnode)
combined operation. Unlink node and free used memory
void FreeNode(XMLNodePointer_t xmlnode)
release all memory, allocated from this node and destroys node itself
void SkipEmpty(XMLNodePointer_t &xmlnode)
Skip all current empty nodes and locate on first "true" node.
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 * XmlDefineClass(const char *xmlClassName)
define class for the converted class name, where special symbols were replaced by '_'
const char * XmlClassNameSpaceRef(const TClass *cl)
produce string which used as reference in class namespace definition
EXMLLayout GetXmlLayout() const
Definition TXMLSetup.h:97
virtual void SetUseNamespaces(Bool_t iUseNamespaces=kTRUE)
Definition TXMLSetup.h:105
const char * XmlConvertClassName(const char *name)
convert class name to exclude any special symbols like ':', '<' '>' ',' and spaces
Int_t AtoI(const char *sbuf, Int_t def=0, const char *errinfo=nullptr)
converts string to integer.
const char * XmlGetElementName(const TStreamerElement *el)
return converted name for TStreamerElement
Int_t GetNextRefCounter()
Definition TXMLSetup.h:111
Bool_t IsUseNamespaces() const
Definition TXMLSetup.h:100
@ kSpecialized
Definition TXMLSetup.h:84
@ kGeneralized
Definition TXMLSetup.h:84
virtual void SetXmlLayout(EXMLLayout layout)
Definition TXMLSetup.h:102
Bool_t fIsStreamerInfo
TStreamerInfo * fInfo
TStreamerElement * fElem
Bool_t fCompressedClassNode
XMLNodePointer_t fNode
XMLNsPointer_t fClassNs
Bool_t fIsElemOwner
Bool_t IsStreamerInfo() const
TXMLStackObj(XMLNodePointer_t node)
const Int_t n
Definition legend1.C:16
Definition file.py:1
const char * UChar
Definition TXMLSetup.cxx:89
const char * Ptr
Definition TXMLSetup.cxx:52
const char * Name
Definition TXMLSetup.cxx:67
const char * v
Definition TXMLSetup.cxx:74
const char * Bool
Definition TXMLSetup.cxx:81
const char * Long64
Definition TXMLSetup.cxx:86
const char * False
Definition TXMLSetup.cxx:77
const char * True
Definition TXMLSetup.cxx:76
const char * Int
Definition TXMLSetup.cxx:84
const char * ULong64
Definition TXMLSetup.cxx:93
const char * Member
Definition TXMLSetup.cxx:65
const char * OnlyVersion
Definition TXMLSetup.cxx:51
const char * Long
Definition TXMLSetup.cxx:85
const char * Float
Definition TXMLSetup.cxx:87
const char * Array
Definition TXMLSetup.cxx:80
const char * ClassVersion
Definition TXMLSetup.cxx:49
const char * String
Definition TXMLSetup.cxx:94
const char * Double
Definition TXMLSetup.cxx:88
const char * Object
Definition TXMLSetup.cxx:62
const char * Ref
Definition TXMLSetup.cxx:53
const char * cnt
Definition TXMLSetup.cxx:75
const char * IdBase
Definition TXMLSetup.cxx:55
const char * Size
Definition TXMLSetup.cxx:56
const char * XmlBlock
Definition TXMLSetup.cxx:60
const char * Null
Definition TXMLSetup.cxx:54
const char * Char
Definition TXMLSetup.cxx:82
const char * UShort
Definition TXMLSetup.cxx:90
const char * CharStar
Definition TXMLSetup.cxx:95
const char * UInt
Definition TXMLSetup.cxx:91
const char * ULong
Definition TXMLSetup.cxx:92
const char * Class
Definition TXMLSetup.cxx:64
const char * ObjClass
Definition TXMLSetup.cxx:63
const char * Short
Definition TXMLSetup.cxx:83
const char * Zip
Definition TXMLSetup.cxx:61
const char * Item
Definition TXMLSetup.cxx:66
EValues
Note: this is only temporarily a struct and will become a enum class hence the name.
Definition Compression.h:83
@ kUndefined
Undefined compression algorithm (must be kept the last of the list in case a new algorithm is added).
@ kUseMin
Compression level reserved when we are not sure what to use (1 is for the fastest compression)
Definition Compression.h:68
auto * l
Definition textangle.C:4