Logo ROOT   6.16/01
Reference Guide
xmlreadfile.C File Reference

Detailed Description

View in nbviewer Open in SWAN Example to read and parse any xml file, supported by TXMLEngine class The input file, produced by xmlnewfile.C macro is used If you need full xml syntax support, use TXMLParser instead

#include "TXMLEngine.h"
#include <stdio.h>
void DisplayNode(TXMLEngine &xml, XMLNodePointer_t node, Int_t level)
{
// this function display all accessible information about xml node and its children
printf("%*c node: %s\n", level, ' ', xml.GetNodeName(node));
// display namespace
XMLNsPointer_t ns = xml.GetNS(node);
if (ns != 0)
printf("%*c namespace: %s refer: %s\n", level + 2, ' ', xml.GetNSName(ns), xml.GetNSReference(ns));
// display attributes
XMLAttrPointer_t attr = xml.GetFirstAttr(node);
while (attr != 0) {
printf("%*c attr: %s value: %s\n", level + 2, ' ', xml.GetAttrName(attr), xml.GetAttrValue(attr));
attr = xml.GetNextAttr(attr);
}
// display content (if exists)
const char *content = xml.GetNodeContent(node);
if (content != 0)
printf("%*c cont: %s\n", level + 2, ' ', content);
// display all child nodes
XMLNodePointer_t child = xml.GetChild(node);
while (child != 0) {
DisplayNode(xml, child, level + 2);
child = xml.GetNext(child);
}
}
void xmlreadfile(const char* filename = "example.xml")
{
// First create engine
// Now try to parse xml file
// Only file with restricted xml syntax are supported
XMLDocPointer_t xmldoc = xml.ParseFile(filename);
if (!xmldoc) return;
// take access to main node
XMLNodePointer_t mainnode = xml.DocGetRootElement(xmldoc);
// display recursively all nodes and subnodes
DisplayNode(xml, mainnode, 1);
// Release memory before exit
xml.FreeDoc(xmldoc);
}
int Int_t
Definition: RtypesCore.h:41
void * XMLNodePointer_t
Definition: TXMLEngine.h:17
void * XMLNsPointer_t
Definition: TXMLEngine.h:18
void * XMLDocPointer_t
Definition: TXMLEngine.h:20
void * XMLAttrPointer_t
Definition: TXMLEngine.h:19
const char * GetNSName(XMLNsPointer_t ns)
return name id of namespace
Definition: TXMLEngine.cxx:770
XMLNodePointer_t GetChild(XMLNodePointer_t xmlnode, Bool_t realnode=kTRUE)
returns first child of xmlnode
void FreeDoc(XMLDocPointer_t xmldoc)
frees allocated document data and deletes document itself
XMLNodePointer_t DocGetRootElement(XMLDocPointer_t xmldoc)
returns root node of document
XMLAttrPointer_t GetNextAttr(XMLAttrPointer_t xmlattr)
return next attribute in the list
Definition: TXMLEngine.cxx:673
const char * GetNSReference(XMLNsPointer_t ns)
return reference id of namespace
Definition: TXMLEngine.cxx:783
const char * GetNodeContent(XMLNodePointer_t xmlnode)
get contents (if any) of xmlnode
XMLNsPointer_t GetNS(XMLNodePointer_t xmlnode)
return namespace attribute (if exists)
Definition: TXMLEngine.cxx:758
const char * GetAttrName(XMLAttrPointer_t xmlattr)
return name of the attribute
Definition: TXMLEngine.cxx:684
XMLAttrPointer_t GetFirstAttr(XMLNodePointer_t xmlnode)
return first attribute in the list, namespace (if exists) will be skipped
Definition: TXMLEngine.cxx:657
const char * GetNodeName(XMLNodePointer_t xmlnode)
returns name of xmlnode
XMLDocPointer_t ParseFile(const char *filename, Int_t maxbuf=100000)
Parses content of file and tries to produce xml structures.
const char * GetAttrValue(XMLAttrPointer_t xmlattr)
return value of attribute
Definition: TXMLEngine.cxx:695
XMLNodePointer_t GetNext(XMLNodePointer_t xmlnode, Bool_t realnode=kTRUE)
return next to xmlnode node if realnode==kTRUE, any special nodes in between will be skipped
static constexpr double ns
Author
Sergey Linev

Definition in file xmlreadfile.C.