Xdmf3 C++ API

From XdmfWeb
Jump to navigationJump to search

XdmfLogo1.gif

XDMF API

While use of the XDMF API is not necessary to produce or consume valid XDMF datasets, there are many convenience features that make it attractive.

All XDMF Objects are derived from XdmfItem.

XdmfDomain

To understand access to XDMF data, understanding of the XdmfDOM is critical. XDMF uses the libxml2 library to parse and generate XML documents. The Domain is an in-memory tree structure that represents the XML file. Nodes of the tree are added, deleted, queried and serialized.

A Domain can be created inside the code via the New function:

shared_ptr<XdmfDomain> domain = XdmfDomain::New();

A reader can be used to produce a domain from an existing file:

shared_ptr<XdmfReader> reader = XdmfReader::New();
shared_ptr<XdmfDomain> domain =
 shared_dynamic_cast<XdmfDomain>(reader->read("MyXMLFile.xmf"));

Once the XML has been parsed into the Domain, the tree can be navigated and modified. The domain is the root of the Xdmf tree and by accessing its children via XdmfDomain::getUnstructuredGrid(unsigned int) and similar functions.

shared_ptr<XdmfUnstructuredGrid> ungrid = domain->getUnstructuredGrid(0);

The same grid can also be fetched by name.

shared_ptr<XdmfUnstructuredGrid> ungrid = domain->getUnstructuredGrid("Unstructured Grid");

Once you have the child object, that child has its own functions that can be used to modify or retrieve attributes or children.

std::string originalName = ungrid->getName();
ungrid->setName("New Name");

When finished building and modifying the XdmfTree it can then be written out to file with an XdmfWriter visitor.

shared_ptr<XdmfWriter> writer = XdmfWriter::New("outfile.xmf");
domain->accept(writer);

XdmfArray

The XdmfArray class is a self describing data structure. It contains descriptions for its internal number type, precision, and shape (rank and dimensions). Many XDMF classes require an XdmfArray as input to methods. The following C++ example demonstrates creating an interlaced XYZ array from three separate variables:

float           *x, *y, *z;
unsigned int total;
std::vector<unsigned int> dims;
dims.push_back(10);
dims.push_back(20);
dims.push_back(30);
shared_ptr<XdmfArray> xyz = XdmfArray::New();
total = 10 * 20 * 30;
xyz->initialize(XdmfArrayType::Float64(), dims);    // KDim, JDim, IDim
xyz->insert(0, x, 0, total, 3, 1);  //       insert<T>(unsigned int index, T *Values,
                                 //                 unsigned int NumberOfValues, unsigned int ArrayStride=1,
                                 //                 unsigned int ValuesStride=1)
xyz->insert(1, y, 0, total, 3, 1);
xyz->insert(2, z, 0, total, 3, 1);


XdmfTime

When modifying the Xdmf Tree, XdmfTime objects can be added to grids in order to specify timestamps.

shared_ptr<XdmfTime> time = XdmfTime::New(5.0);
ungrid->setTime(time);


XdmfHDF

In XDMF, Light data is stored in XML while the Heavy data is typically stored in an HDF5 file. Passing an XdmfHDF5Writer to an array will write its contents to the specified hdf5 file. A controller object is created to handle the data description and allows the array to read the data back in if needed. The default write mode creates a new hdf5 dataspace for each written array. HDF5 write modes include:

  • Default
  • Overwrite
  • Append
  • Hyperslab
double * Geometry = [-1.75, -1.25, 0, -1.25, -1.25, 0, -0.75, ..];
unsigned int * Connectivity = [3, 2, 5, 1, ..];
int Values = [100, 200, 300, ..];
shared_ptr<XdmfHDF5Writer> heavywriter = XdmfHDF5Writer->New("example.h5");
// Geometry
shared_ptr<XdmfGeometry> GeometryArray = XdmfGeometry::New();
GeometryArray->insert(0, Geometry, geometryNumValues);
GeometryArray->accept(heavywriter);
// Coneectivity
shared_ptr<XdmfTopology> ConnectivityArray = XdmfTopology::New();
ConnectivityArray->insert(0, Connectivity, connectivityNumValues);
ConnectivityArray->accept(heavywriter);
// Values
shared_ptr<XdmfArray> ValueArray = XdmfArray::New();
ValueArray->insert(0, Values, numValues);
ValueArray->accept(heavywriter);

When reading, if reading from an XML file where the arrays have the hdf5 data sets specified, controllers will automatically be created for all arrays with a heavy data description. To access the data in heavy data XdmfArray::read() must be called. If it is required to link heavy data manually a heavy data controller can be created manually and added to an array.

std::string newPath = "File path to hdf5 file goes here";
std::string newSetPath = "path to the set goes here";
shared_ptr<const XdmfArrayType> readType = XdmfArrayType::Int32();
std::vector<unsigned int> readStarts;
//Three dimensions, all starting at index 0
readStarts.push_back(0);
readStarts.push_back(0);
readStarts.push_back(0);
std::vector<unsigned int> readStrides;
//Three dimensions, no skipping between index
readStrides.push_back(1);
readStrides.push_back(1);
readStrides.push_back(1);
std::vector<unsigned int> readCounts;
//Three dimensions, reading 10 values from each
readCounts.push_back(10);
readCounts.push_back(10);
readCounts.push_back(10);
std::vector<unsigned int> readDataSize;
//Three dimensins, each with 20 maximum values
readDataSize.push_back(20);
readDataSize.push_back(20);
readDataSize.push_back(20);
shared_ptr<XdmfHDF5Controller> controller = XdmfHDF5Controller::New(
               newPath,
               newSetPath,
               readType,
               readStarts,
               readStrides,
               readCounts,
               readDataSize);
array->insert(controller);
array->read();

Reading XDMF

Putting all of this together, assume Points.xmf is a valid XDMF XML file with a single Grid. Here is a C++ example to read and print values.

shared_ptr<XdmfReader> reader = XdmfReader->New();
shared_ptr<XdmfDomain> domain =
  shared_dynamic_cast<XdmfDomain>(reader->.read("Points.xmf"));
shared_ptr<XdmfUnstructuredGrid> grid = domain->getUnstructuredGrid(0);
shared_ptr<XdmfTopology> topology = grid->getTopology();
printf("Values = %s\n", topology->getValuesString().c_str());
shared_ptr<XdmfGeometry> geometry = grid->getGeometry();
printf("Geo Type =  %s #  Points = %u\n", geometry->getType()->getName().c_str(), geometry->getNumberOfPoints());
printf("Points = %s\n", geometry->getValuesString());

Writing XDMF

Using the insert() and set() methods, an XDMF dataset can be generated programmatically as well. An object tree mirroring the XML is built as Xdmf objects are attached to each other. Since this is handled via shared pointers the objects are not duplicated if there are multiple instances of that object in the tree and the user does not have to worry about memory cleanup.

shared_ptr<XdmfDomain>  root = XdmfDomain::New();
# Information
shared_ptr<XdmfInformation>  i = XdmfInformation::New(); # Arbitrary Name=Value Facility
i->setName("Time");
i->setValue(0.0123);
root->insert(i); # XdmfDomain is the root of the tree
                 # insert adds the information we just created as a leaf
# Dimensions
shared_ptr<XdmfArray> dimensions = XdmfArray::New();
dimensions->pushBack((unsigned int)10);
dimensions->pushBack((unsigned int)20);
dimensions->pushBack((unsigned int)30);
# Origin
shared_ptr<XdmfArray> origin = XdmfArray::New();
origin->pushBack((unsigned int)1);
origin->pushBack((unsigned int)2);
origin->pushBack((unsigned int)3);
# Brick Size
shared_ptr<XdmfArray> brick = XdmfArray::New();
brick->pushBack((unsigned int)1);
brick->pushBack((unsigned int)2);
brick->pushBack((unsigned int)3);
# Grid
shared_ptr<XdmfRegularGrid> g = XdmfRegularGrid::New(brick, dimensions, origin);
g->setName("Structured Grid");
root->insert(g);
# Attribute
shared_ptr<XdmfAttribute> attr = XdmfAttribute::New();
attr->setName("Pressure");
attr->setAttributeCenter(XdmfAttributeCenter::Cell());
attr->setAttributeType(XdmfAttributeType::Scalar());
attr->initialize(XdmfArrayType::Float64(), 10 * 20 * 30);
double * pressureVals = [0, 1, 2, 3, 4, 5, ...];
attr->insert(0, pressureVals, numPressureVals);
# The heavy data set name is determined by the writer if not set
g->insert(attr);
# Generate a writer
shared_ptr<XdmfWriter> writer = XdmfWriter::New("output.xmf");
domain->accept(writer)