EpetraExt  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
EpetraExt::HDF5 Class Reference

class HDF5: A class for storing Epetra objects in parallel binary files More...

#include <EpetraExt_HDF5.h>

 HDF5 (const Epetra_Comm &Comm)
 Constructor. More...
 
 ~HDF5 ()
 Destructor. More...
 
void Create (const std::string FileName)
 Create a new file. More...
 
void Open (const std::string FileName, int AccessType=H5F_ACC_RDWR)
 Open specified file with given access type. More...
 
void Close ()
 Close the file. More...
 
void Flush ()
 Flush the content to the file. More...
 
bool IsOpen () const
 Return true if a file has already been opened using Open()/Create() More...
 
void CreateGroup (const std::string &GroupName)
 Create group GroupName. More...
 
bool IsContained (std::string Name, std::string GroupName="")
 Return true if Name is contained in the database. More...
 
void Write (const std::string &GroupName, const std::string &DataSetName, int data)
 Write an integer in group GroupName using the given DataSetName. More...
 
void Read (const std::string &GroupName, const std::string &DataSetName, int &data)
 Read an integer from group /GroupName/DataSetName. More...
 
void Write (const std::string &GroupName, const std::string &DataSetName, double data)
 Write a double in group GroupName using the given DataSetName. More...
 
void Read (const std::string &GroupName, const std::string &DataSetName, double &data)
 Read a double from group /GroupName/DataSetName. More...
 
void Write (const std::string &GroupName, const std::string &DataSetName, const std::string &data)
 Write a string in group GroupName using the given DataSetName. More...
 
void Read (const std::string &GroupName, const std::string &DataSetName, std::string &data)
 Read a string from group /GroupName/DataSetName. More...
 
void Read (const std::string &GroupName, const std::string &DataSetName, const hid_t type, const int Length, void *data)
 Read the serial array data, of type type, from group GroupName, using the dataset name DataSetName. More...
 
void Write (const std::string &GroupName, const std::string &DataSetName, const hid_t type, const int Length, const void *data)
 Write the serial array data, of type type, to group GroupName, using the dataset name DataSetName. More...
 
void WriteComment (const std::string &GroupName, std::string Comment)
 Associate string Comment with group GroupName. More...
 
void ReadComment (const std::string &GroupName, std::string &Comment)
 Read the string associated with group GroupName. More...
 
void Write (const std::string &GroupName, const std::string &DataSetName, int MySize, int GlobalSize, hid_t type, const void *data)
 Write the distributed array data, of type type, to group GroupName, using dataset name DataSetName. More...
 
void Read (const std::string &GroupName, const std::string &DataSetName, int MySize, int GlobalSize, const hid_t type, void *data)
 Read the distributed array data, of type type, from group GroupName, using dataset name DataSetName. More...
 
void Write (const std::string &GroupName, const Epetra_Map &Map)
 Write a Map to group GroupName. More...
 
void Read (const std::string &GroupName, Epetra_Map *&Map)
 Read a map from GroupName. More...
 
void ReadMapProperties (const std::string &GroupName, int &NumGlobalElements, int &IndexBase, int &NumProc)
 Read basic properties of specified Epetra_Map. More...
 
void Read (const std::string &GroupName, Epetra_BlockMap *&Map)
 Read a block map from GroupName. More...
 
void Write (const std::string &GroupName, const Epetra_BlockMap &Map)
 Write a block map to group GroupName. More...
 
void ReadBlockMapProperties (const std::string &GroupName, int &NumGlobalElements, int &NumGlobalPoints, int &IndexBase, int &NumProc)
 Read basic properties of specified Epetra_BlockMap. More...
 
void Read (const std::string &GroupName, Epetra_CrsGraph *&Graph)
 Read a vector from group GroupName, assuming linear distribution. More...
 
void Read (const std::string &GroupName, const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, Epetra_CrsGraph *&Graph)
 Read a vector from group GroupName using the given map. More...
 
void Write (const std::string &GroupName, const Epetra_CrsGraph &Graph)
 Write a distributed vector to group GroupName. More...
 
void ReadCrsGraphProperties (const std::string &GroupName, int &NumGlobalRows, int &NumGlobalCols, int &NumGlobalNonzeros, int &NumGlobalDiagonals, int &MaxNumIndices)
 Read basic properties of specified Epetra_CrsGraph. More...
 
void Write (const std::string &GroupName, const Epetra_IntVector &x)
 Write a distributed vector to group GroupName. More...
 
void Read (const std::string &GroupName, Epetra_IntVector *&X)
 Read a vector from group GroupName, assuming linear distribution. More...
 
void Read (const std::string &GroupName, const Epetra_Map &Map, Epetra_IntVector *&X)
 Read a vector from group GroupName using the given map. More...
 
void ReadIntVectorProperties (const std::string &GroupName, int &GlobalLength)
 Read basic properties of specified Epetra_IntVector. More...
 
void Write (const std::string &GroupName, const Epetra_MultiVector &x, bool writeTranspose=false)
 Write a distributed vector (or its transpose) to group GroupName. More...
 
void Read (const std::string &GroupName, Epetra_MultiVector *&X, bool writeTranspose=false, const int &indexBase=0)
 Read a vector (or its transpose) from group GroupName. More...
 
void Read (const std::string &GroupName, const Epetra_Map &Map, Epetra_MultiVector *&X, bool writeTranspose=false)
 Read a vector from group GroupName using the given map. More...
 
void ReadMultiVectorProperties (const std::string &GroupName, int &GlobalLength, int &NumVectors)
 Read basic properties of specified Epetra_MultiVector. More...
 
void Write (const std::string &GroupName, const Epetra_RowMatrix &Matrix)
 Write a distributed RowMatrix to group GroupName. More...
 
void Read (const std::string &GroupName, Epetra_CrsMatrix *&A)
 Read a square matrix from group GroupName, assuming linear distribution. More...
 
void Read (const std::string &GroupName, const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, Epetra_CrsMatrix *&A)
 Read a matrix from group GroupName with given range and domain maps. More...
 
void ReadCrsMatrixProperties (const std::string &GroupName, int &NumGlobalRows, int &NumGlobalCols, int &NumNonzeros, int &NumGlobalDiagonals, int &MaxNumEntries, double &NormOne, double &NormInf)
 Read basic properties of specified Epetra_CrsMatrix. More...
 
void Write (const std::string &GroupName, const Teuchos::ParameterList &List)
 Write a parameter list to group GroupName. More...
 
void Read (const std::string &GroupName, Teuchos::ParameterList &List)
 Read a parameter list from group GroupName. More...
 
void Write (const std::string &GroupName, const DistArray< int > &array)
 Write an EpetraExt::DistArray<int> to group GroupName. More...
 
void Read (const std::string &GroupName, DistArray< int > *&array)
 Read an EpetraExt::DistArray<int> from group GroupName. More...
 
void Read (const std::string &GroupName, const Epetra_Map &Map, DistArray< int > *&array)
 Read an EpetraExt::DistArray<int> from group GroupName. More...
 
void ReadIntDistArrayProperties (const std::string &GroupName, int &GlobalLength, int &RowSize)
 Read the global number of elements and type for a generic handle object. More...
 
void Write (const std::string &GroupName, const DistArray< double > &array)
 Write an EpetraExt::DistArray<int> to group GroupName. More...
 
void Read (const std::string &GroupName, DistArray< double > *&array)
 Read an EpetraExt::DistArray<int> from group GroupName. More...
 
void Read (const std::string &GroupName, const Epetra_Map &Map, DistArray< double > *&array)
 Read an EpetraExt::DistArray<int> from group GroupName. More...
 
void ReadDoubleDistArrayProperties (const std::string &GroupName, int &GlobalLength, int &RowSize)
 Read the global number of elements and type for a generic handle object. More...
 
void Write (const std::string &GroupName, const Handle &List)
 Write an Epetra_DistObject to group GroupName. More...
 
void Read (const std::string &GroupName, Handle &List)
 Read an Epetra_DistObject from group GroupName. More...
 
void ReadHandleProperties (const std::string &GroupName, std::string &Type, int &NumGlobalElements)
 Read the global number of elements and type for a generic handle object. More...
 

Detailed Description

class HDF5: A class for storing Epetra objects in parallel binary files

Introduction

The HDF5 class reads and writes data using the HDF5 parallel binary data format. HDF5 has the following advantages:

In order to use HDF5, make sure to set the CMake Boolean option EpetraExt_USING_HDF5 when building Trilinos. You may have to tell CMake where to find the HDF5 library and header files, if they are not installed in their default locations.

The class supports input (I), output (O), or both (I/O) for the following distributed Epetra objects:

The class also supports some non-distributed types:

The HDF5 class assumes that non-distributed data types have the same value on all processors.

This class also provides utility methods:

By using these methods, as well as the other methods to write non-distributed types, one can read and write any serial or distributed object.

Data Model

The HDF5 library itself can be used to define very general data formats. Our HDF5 class, instead, is structured around the concept of groups. A group is an entity, for example a scalar value, an Epetra_Map, or a Teuchos::ParameterList. Within each group, different datasets describe the content of the group. For example, an Epetra_MultiVector is specified by datasets NumVectors and Values, which contain the number of vectors and the numerical values, respectively. The comment of each group is a character string that must match the class name.

Our HDF5 class has the following limitations:

Errors

When an error occurs, a EpetraExt::Exception is thrown. The Print() method of the Exception class returns a description of what went wrong.

Example of usage

First, one must create an HDF5 class, then either Open() or Create() the file:

// Comm is an Epetra_Comm communicator wrapper object.
HDF5.Create("myfile.h5");

Writing commands might be as follows:

Epetra_Map* Map = <create map here>
Epetra_Map* BlockMap = <create block map here>
Epetra_RowMatrix* Matrix = <create matrix here>
Epetra_MultiVector* LHS = <...>
Epetra_MultiVector* RHS = <...>
// write a map, whose group name contains the number of processors
HDF5.Write("map-" + toString(Comm.NumProc()), *Map);
HDF5.Write("matrix", *Matrix);
HDF5.Write("LHS", LHS);
HDF5.Write("RHS", RHS);

To write a Teuchos::ParameterList, simply do

HDF5.Write("EpetraExtList", EpetraExtList);

The file can contain metadata as well. Each metadatum is defined by a group and a dataset name. A group may contain more than one dataset, and may be a new group or an already existing group. For example, to specify the numerical quadrature formula used to assemble the matrix, do as follows:

HDF5.Write("matrix", "quadrature order", 3);

Alternatively, datasets may be assigned to a new group, let's say "my parameters":

HDF5.Write("my parameters", "latitude", 12);
HDF5.Write("my parameters", "longitude", 67);
HDF5.Write("my parameters", "angle", 12.3);
vector<int> iarray(3);
iarray[0] = 0, iarray[1] = 1; iarray[2] = 2;
HDF5.Write("my parameters", "int array", H5T_NATIVE_INT, 3, &iarray[0]);
vector<double> darray(3);
darray[0] = 0.1, darray[1] = 1.1; darray[2] = 2.1;
HDF5.Write("my parameters", "double array", H5T_NATIVE_DOUBLE, 3, &darray[0]);

Note that all non-distributed datasets must have the same value on all processors.

Reading data is as easy as writing. Let us consider how to read an Epetra_CrsMatrix. Other Epetra objects having a similar behavior. The ReadCrsMatrixProperties() method can be used to query for some matrix properties without reading the whole matrix:

int NumGlobalRows, NumGlobalCols, NumGlobalNonzeros;
int NumGlobalDiagonals, MaxNumEntries;
double NormOne, NormInf;
ReadCrsMatrixProperties(GroupName, NumGlobalRows, NumGlobalCols,
NumGlobalNonzeros, NumGlobalDiagonals, MaxNumEntries,
NormOne, NormInf);

The above call is not required, and can be skipped. In order to read the Epetra_CrsMatrix, do as follows:

Epetra_CrsMatrix* NewMatrix = NULL;
HDF5.Read("matrix", NewMatrix);

In this case, NewMatrix is based on a linear map. If the matrix's DomainMap() and RangeMap() are known and non-trivial, one can use

HDF5.Read("matrix", DomainMap, RangeMap, NewMatrix);

Reading metadata looks like:

HDF5.Read("my parameters", "latitude", new_latitude);
HDF5.Read("my parameters", "longitude", new_longitude);
HDF5.Read("my parameters", "int array", H5T_NATIVE_INT, 3, &new_iarray[0]);
HDF5.Read("my parameters", "double array", H5T_NATIVE_DOUBLE, 3, &new_darray[0]);

To analyze the content of the file, one can use "h5dump filename.h5" or "h5dump filename.h5 -H".

MATLAB Interface

MATLAB provides built-in functions to read, write, and query HDF5 files: hdf5read, hdf5write, and hdf5info, respectively. For example, to read the above Epetra_CrsMatrix into MATLAB as a sparse matrix, do as follows:

NumGlobalRows = double(hdf5read('myfile.h5', '/matrix/NumGlobalRows/'));
NumGlobalCols = double(hdf5read('myfile.h5', '/matrix/NumGlobalCols/'));
ROW = double(hdf5read('myfile.h5', '/matrix/ROW/'));
COL = double(hdf5read('myfile.h5', '/matrix/COL/'));
VAL = hdf5read('myfile.h5', '/matrix/VAL/');
A = sparse(ROW + 1, COL + 1, VAL, NumGlobalRows, NumGlobalCols);

The use of double() is required by Matlab's sparse function, since it does not accept int32 data.

To dump a MATLAB matrix (in this case, A) to the file matlab.h5, do as follows:

n = 10;
A = speye(n, n);
[ROW,COL,VAL] = find(A);
hdf5write('matlab.h5', '/speye/__type__', 'Epetra_RowMatrix');
hdf5write('matlab.h5', '/speye/NumGlobalRows', int32(n), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/NumGlobalCols', int32(n), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/NumGlobalNonzeros', int32(n), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/NumGlobalDiagonals', int32(n), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/MaxNumEntries', int32(1), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/NormOne', 1.0, 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/NormInf', 1.0, 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/ROW', int32(ROW - 1), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/COL', int32(COL - 1), 'WriteMode', 'append');
hdf5write('matlab.h5', '/speye/VAL', VAL, 'WriteMode', 'append');

Note that the type specification must reflect the Epetra class name.

To dump a MATLAB dense array (in this case, x) to the file matlab.h5, do as follows:

n = 10;
x = [zeros(n,1), rand(n, 1)]';
hdf5write('matlab.h5', '/x/__type__', 'Epetra_MultiVector');
hdf5write('matlab.h5', '/x/GlobalLength',int32(n), 'WriteMode', 'append');
hdf5write('matlab.h5', '/x/NumVectors', int32(2), 'WriteMode', 'append');
hdf5write('matlab.h5', '/x/Values', x, 'WriteMode', 'append');

Note that MATLAB vectors must be stored as row vectors.

You can also write a Map from MATLAB for use by Epetra. The following example shows how to define an Epetra_Map that distributes data over two processors:

IndexBase = 0;
NumMyElements = [5 5];
n = 10;
MyGlobalElements = [5 6 7 8 9 0 1 2 3 4];
hdf5write('matlab.h5', '/map-2/__type__', 'Epetra_Map');
hdf5write('matlab.h5', '/map-2/NumGlobalElements', int32(n), 'WriteMode', 'append');
hdf5write('matlab.h5', '/map-2/IndexBase', int32(IndexBase), 'WriteMode', 'append');
hdf5write('matlab.h5', '/map-2/NumProc', int32(2), 'WriteMode', 'append');
hdf5write('matlab.h5', '/map-2/NumMyElements', int32(NumMyElements), 'WriteMode', 'append');
hdf5write('matlab.h5', '/map-2/MyGlobalElements', int32(MyGlobalElements), 'WriteMode', 'append');
Author
Marzio Sala, D-INFK/ETHZ
Date
Last updated on 16-Mar-06. Documentation revised by Mark Hoemmen (05 Oct 2011) for spelling, grammar, and clarity.
Todo:
  • all distributed objects are assumed in local state, this is not necessary (just easier)
  • Epetra_VbrMatrix input has to be done, now it is considered as an Epetra_RowMatrix

Definition at line 326 of file EpetraExt_HDF5.h.

Constructor & Destructor Documentation

EpetraExt::HDF5::HDF5 ( const Epetra_Comm Comm)

Constructor.

Definition at line 272 of file EpetraExt_HDF5.cpp.

EpetraExt::HDF5::~HDF5 ( )
inline

Destructor.

Definition at line 335 of file EpetraExt_HDF5.h.

Member Function Documentation

void EpetraExt::HDF5::Create ( const std::string  FileName)

Create a new file.

Definition at line 278 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Open ( const std::string  FileName,
int  AccessType = H5F_ACC_RDWR 
)

Open specified file with given access type.

Definition at line 365 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Close ( )
inline

Close the file.

Definition at line 351 of file EpetraExt_HDF5.h.

void EpetraExt::HDF5::Flush ( )
inline

Flush the content to the file.

Definition at line 358 of file EpetraExt_HDF5.h.

bool EpetraExt::HDF5::IsOpen ( ) const
inline

Return true if a file has already been opened using Open()/Create()

Definition at line 364 of file EpetraExt_HDF5.h.

void EpetraExt::HDF5::CreateGroup ( const std::string &  GroupName)
inline

Create group GroupName.

Definition at line 370 of file EpetraExt_HDF5.h.

bool EpetraExt::HDF5::IsContained ( std::string  Name,
std::string  GroupName = "" 
)

Return true if Name is contained in the database.

Definition at line 395 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const std::string &  DataSetName,
int  data 
)

Write an integer in group GroupName using the given DataSetName.

Definition at line 1460 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const std::string &  DataSetName,
int &  data 
)

Read an integer from group /GroupName/DataSetName.

Definition at line 1504 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const std::string &  DataSetName,
double  data 
)

Write a double in group GroupName using the given DataSetName.

Definition at line 1482 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const std::string &  DataSetName,
double &  data 
)

Read a double from group /GroupName/DataSetName.

Definition at line 1526 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const std::string &  DataSetName,
const std::string &  data 
)

Write a string in group GroupName using the given DataSetName.

Definition at line 1548 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const std::string &  DataSetName,
std::string &  data 
)

Read a string from group /GroupName/DataSetName.

Definition at line 1579 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const std::string &  DataSetName,
const hid_t  type,
const int  Length,
void *  data 
)

Read the serial array data, of type type, from group GroupName, using the dataset name DataSetName.

Definition at line 1639 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const std::string &  DataSetName,
const hid_t  type,
const int  Length,
const void *  data 
)

Write the serial array data, of type type, to group GroupName, using the dataset name DataSetName.

Definition at line 1612 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::WriteComment ( const std::string &  GroupName,
std::string  Comment 
)
inline

Associate string Comment with group GroupName.

Definition at line 410 of file EpetraExt_HDF5.h.

void EpetraExt::HDF5::ReadComment ( const std::string &  GroupName,
std::string &  Comment 
)
inline

Read the string associated with group GroupName.

Definition at line 416 of file EpetraExt_HDF5.h.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const std::string &  DataSetName,
int  MySize,
int  GlobalSize,
hid_t  type,
const void *  data 
)

Write the distributed array data, of type type, to group GroupName, using dataset name DataSetName.

Definition at line 1670 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const std::string &  DataSetName,
int  MySize,
int  GlobalSize,
const hid_t  type,
void *  data 
)

Read the distributed array data, of type type, from group GroupName, using dataset name DataSetName.

Definition at line 1720 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Epetra_Map Map 
)

Write a Map to group GroupName.

Definition at line 520 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Epetra_Map *&  Map 
)

Read a map from GroupName.

Definition at line 539 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadMapProperties ( const std::string &  GroupName,
int &  NumGlobalElements,
int &  IndexBase,
int &  NumProc 
)

Read basic properties of specified Epetra_Map.

Definition at line 565 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Epetra_BlockMap *&  Map 
)

Read a block map from GroupName.

Definition at line 459 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Epetra_BlockMap Map 
)

Write a block map to group GroupName.

Definition at line 428 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadBlockMapProperties ( const std::string &  GroupName,
int &  NumGlobalElements,
int &  NumGlobalPoints,
int &  IndexBase,
int &  NumProc 
)

Read basic properties of specified Epetra_BlockMap.

Definition at line 495 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Epetra_CrsGraph *&  Graph 
)

Read a vector from group GroupName, assuming linear distribution.

Definition at line 728 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const Epetra_Map DomainMap,
const Epetra_Map RangeMap,
Epetra_CrsGraph *&  Graph 
)

Read a vector from group GroupName using the given map.

Definition at line 771 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Epetra_CrsGraph Graph 
)

Write a distributed vector to group GroupName.

Definition at line 690 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadCrsGraphProperties ( const std::string &  GroupName,
int &  NumGlobalRows,
int &  NumGlobalCols,
int &  NumGlobalNonzeros,
int &  NumGlobalDiagonals,
int &  MaxNumIndices 
)

Read basic properties of specified Epetra_CrsGraph.

Definition at line 744 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Epetra_IntVector x 
)

Write a distributed vector to group GroupName.

Definition at line 592 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Epetra_IntVector *&  X 
)

Read a vector from group GroupName, assuming linear distribution.

Definition at line 619 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const Epetra_Map Map,
Epetra_IntVector *&  X 
)

Read a vector from group GroupName using the given map.

Definition at line 634 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadIntVectorProperties ( const std::string &  GroupName,
int &  GlobalLength 
)

Read basic properties of specified Epetra_IntVector.

Definition at line 667 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Epetra_MultiVector x,
bool  writeTranspose = false 
)

Write a distributed vector (or its transpose) to group GroupName.

Write the transpose if writeTranspose is true.

Definition at line 1017 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Epetra_MultiVector *&  X,
bool  writeTranspose = false,
const int &  indexBase = 0 
)

Read a vector (or its transpose) from group GroupName.

This method assumes a linear distribution. Read the transpose if writeTranspose is true.

Definition at line 1126 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const Epetra_Map Map,
Epetra_MultiVector *&  X,
bool  writeTranspose = false 
)

Read a vector from group GroupName using the given map.

Read the transpose if writeTranspose is true.

Definition at line 1109 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadMultiVectorProperties ( const std::string &  GroupName,
int &  GlobalLength,
int &  NumVectors 
)

Read basic properties of specified Epetra_MultiVector.

Definition at line 1231 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Epetra_RowMatrix Matrix 
)

Write a distributed RowMatrix to group GroupName.

Definition at line 824 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Epetra_CrsMatrix *&  A 
)

Read a square matrix from group GroupName, assuming linear distribution.

Definition at line 868 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const Epetra_Map DomainMap,
const Epetra_Map RangeMap,
Epetra_CrsMatrix *&  A 
)

Read a matrix from group GroupName with given range and domain maps.

Definition at line 886 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadCrsMatrixProperties ( const std::string &  GroupName,
int &  NumGlobalRows,
int &  NumGlobalCols,
int &  NumNonzeros,
int &  NumGlobalDiagonals,
int &  MaxNumEntries,
double &  NormOne,
double &  NormInf 
)

Read basic properties of specified Epetra_CrsMatrix.

Definition at line 928 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Teuchos::ParameterList &  List 
)

Write a parameter list to group GroupName.

Definition at line 963 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Teuchos::ParameterList &  List 
)

Read a parameter list from group GroupName.

Definition at line 981 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const DistArray< int > &  array 
)

Write an EpetraExt::DistArray<int> to group GroupName.

Definition at line 1256 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
DistArray< int > *&  array 
)

Read an EpetraExt::DistArray<int> from group GroupName.

Definition at line 1318 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const Epetra_Map Map,
DistArray< int > *&  array 
)

Read an EpetraExt::DistArray<int> from group GroupName.

Definition at line 1286 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadIntDistArrayProperties ( const std::string &  GroupName,
int &  GlobalLength,
int &  RowSize 
)

Read the global number of elements and type for a generic handle object.

Definition at line 1333 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const DistArray< double > &  array 
)

Write an EpetraExt::DistArray<int> to group GroupName.

Definition at line 1358 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
DistArray< double > *&  array 
)

Read an EpetraExt::DistArray<int> from group GroupName.

Definition at line 1420 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
const Epetra_Map Map,
DistArray< double > *&  array 
)

Read an EpetraExt::DistArray<int> from group GroupName.

Definition at line 1388 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::ReadDoubleDistArrayProperties ( const std::string &  GroupName,
int &  GlobalLength,
int &  RowSize 
)

Read the global number of elements and type for a generic handle object.

Definition at line 1435 of file EpetraExt_HDF5.cpp.

void EpetraExt::HDF5::Write ( const std::string &  GroupName,
const Handle List 
)

Write an Epetra_DistObject to group GroupName.

Definition at line 71 of file EpetraExt_HDF5_DistObject.cpp.

void EpetraExt::HDF5::Read ( const std::string &  GroupName,
Handle List 
)

Read an Epetra_DistObject from group GroupName.

Definition at line 158 of file EpetraExt_HDF5_DistObject.cpp.

void EpetraExt::HDF5::ReadHandleProperties ( const std::string &  GroupName,
std::string &  Type,
int &  NumGlobalElements 
)

Read the global number of elements and type for a generic handle object.

Definition at line 231 of file EpetraExt_HDF5_DistObject.cpp.


The documentation for this class was generated from the following files: