#include <EpetraExt_HDF5.h>
| HDF5 (const Epetra_Comm &Comm) | |
| ctor | |
| ~HDF5 () | |
| dtor | |
| void | Create (const string FileName) |
| Creates a new file. | |
| void | Open (const string FileName, int AccessType=H5F_ACC_RDWR) |
| Opens specified file with given access type. | |
| void | Close () |
| Closes the file. | |
| bool | IsOpen () const |
Returns true is a file has already been open using Open()/Create(). | |
| void | CreateGroup (const string &GroupName) |
Creates group GroupName. | |
| bool | IsContained (const string Name) |
Returns true is Name is contained in the database. | |
| void | Write (const string &GroupName, const string &DataSetName, int data) |
Writes an integer in group GroupName using intentified DataSetName. | |
| void | Read (const string &GroupName, const string &DataSetName, int &data) |
Reads an integer from group /GroupName/DataSetName. | |
| void | Write (const string &GroupName, const string &DataSetName, double data) |
Writes a double in group GroupName using intentified DataSetName. | |
| void | Read (const string &GroupName, const string &DataSetName, double &data) |
Reads a double from group /GroupName/DataSetName. | |
| void | Write (const string &GroupName, const string &DataSetName, const string &data) |
Writes a string in group GroupName using intentified DataSetName. | |
| void | Read (const string &GroupName, const string &DataSetName, string &data) |
Reads a string from group /GroupName/DataSetName. | |
| void | Read (const string &GroupName, const string &DataSetName, const int type, const int Length, void *data) |
Reads serial array data, of type type, from group GroupNameusing dataset name DataSetName. | |
| void | Write (const string &GroupName, const string &DataSetName, const int type, const int Length, void *data) |
Writes serial array data, of type type, to group GroupNameusing dataset name DataSetName. | |
| void | WriteComment (const string &GroupName, string Comment) |
Associates string Comment with group GroupName. | |
| void | ReadComment (const string &GroupName, string &Comment) |
Reads the string associated with group GroupName. | |
| void | Write (const string &GroupName, const string &DataSetName, int MySize, int GlobalSize, int type, const void *data) |
Writes distributed array data, of type type, to group GroupNameusing dataset name DataSetName. | |
| void | Read (const string &GroupName, const string &DataSetName, int MySize, int GlobalSize, const int type, void *data) |
Reads distributed array data, of type type, from group GroupNameusing dataset name DataSetName. | |
| void | Write (const string &GroupName, const Epetra_Map &Map) |
Writes a Map to group GroupName. | |
| void | Read (const string &GroupName, Epetra_Map *&Map) |
Reads a map from GroupName. | |
| void | ReadMapProperties (const string &GroupName, int &NumGlobalElements, int &IndexBase, int &NumProc) |
| Reads basic properties of specified Epetra_Map. | |
| void | Read (const string &GroupName, Epetra_BlockMap *&Map) |
Reads a block map from GroupName. | |
| void | Write (const string &GroupName, const Epetra_BlockMap &Map) |
Writes a block map to group GroupName. | |
| void | ReadBlockMapProperties (const string &GroupName, int &NumGlobalElements, int &NumGlobalPoints, int &IndexBase, int &NumProc) |
| Reads basic properties of specified Epetra_BlockMap. | |
| void | Read (const string &GroupName, Epetra_CrsGraph *&Graph) |
Reads a vector from group GroupName, assumes linear distribution. | |
| void | Read (const string &GroupName, const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, Epetra_CrsGraph *&Graph) |
Reads a vector from group GroupName using given map. | |
| void | Write (const string &GroupName, const Epetra_CrsGraph &Graph) |
Writes a distributed vector to group GroupName. | |
| void | ReadCrsGraphProperties (const string &GroupName, int &NumGlobalRows, int &NumGlobalCols, int &NumGlobalNonzeros, int &NumGlobalDiagonals, int &MaxNumIndices) |
| Reads basic properties of specified Epetra_CrsGraph. | |
| void | Write (const string &GroupName, const Epetra_IntVector &x) |
Writes a distributed vector to group GroupName. | |
| void | Read (const string &GroupName, Epetra_IntVector *&X) |
Reads a vector from group GroupName, assumes linear distribution. | |
| void | Read (const string &GroupName, const Epetra_Map &Map, Epetra_IntVector *&X) |
Reads a vector from group GroupName using given map. | |
| void | ReadIntVectorProperties (const string &GroupName, int &GlobalLength) |
| Reads basic properties of specified Epetra_IntVector. | |
| void | Write (const string &GroupName, const Epetra_MultiVector &x) |
Writes a distributed vector to group GroupName. | |
| void | Read (const string &GroupName, Epetra_MultiVector *&X) |
Reads a vector from group GroupName, assumes linear distribution. | |
| void | Read (const string &GroupName, const Epetra_Map &Map, Epetra_MultiVector *&X) |
Reads a vector from group GroupName using given map. | |
| void | ReadMultiVectorProperties (const string &GroupName, int &GlobalLength, int &NumVectors) |
| Reads basic properties of specified Epetra_MultiVector. | |
| void | Write (const string &GroupName, const Epetra_RowMatrix &Matrix) |
Writes a distributed RowMatrix to group GroupName. | |
| void | Read (const string &GroupName, Epetra_CrsMatrix *&A) |
Reads a square matrix from group GroupName, assumes linear distribution. | |
| void | Read (const string &GroupName, const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, Epetra_CrsMatrix *&A) |
Reads a matrix from group GroupName with given range and domain maps. | |
| void | ReadCrsMatrixProperties (const string &GroupName, int &NumGlobalRows, int &NumGlobalCols, int &NumNonzeros, int &NumGlobalDiagonals, int &MaxNumEntries, double &NormOne, double &NormInf) |
| Reads basic properties of specified Epetra_CrsMatrix. | |
| void | Write (const string &GroupName, const Teuchos::ParameterList &List) |
Writes a parameter list to group GroupName. | |
| void | Read (const string &GroupName, Teuchos::ParameterList &List) |
Reads a parameter list from group GroupName. | |
| void | Write (const string &GroupName, const DistArray< int > &array) |
Writes an EpetraExt::DistArray<int> to group GroupName. | |
| void | Read (const string &GroupName, DistArray< int > *&array) |
Reads an EpetraExt::DistArray<int> from group GroupName. | |
| void | Read (const string &GroupName, const Epetra_Map &Map, DistArray< int > *&array) |
Reads an EpetraExt::DistArray<int> from group GroupName. | |
| void | ReadIntDistArrayProperties (const string &GroupName, int &GlobalLength, int &RowSize) |
| Reads the global number of elements and type for a generic handle object. | |
| void | Write (const string &GroupName, const DistArray< double > &array) |
Writes an EpetraExt::DistArray<int> to group GroupName. | |
| void | Read (const string &GroupName, DistArray< double > *&array) |
Reads an EpetraExt::DistArray<int> from group GroupName. | |
| void | Read (const string &GroupName, const Epetra_Map &Map, DistArray< double > *&array) |
Reads an EpetraExt::DistArray<int> from group GroupName. | |
| void | ReadDoubleDistArrayProperties (const string &GroupName, int &GlobalLength, int &RowSize) |
| Reads the global number of elements and type for a generic handle object. | |
| void | Write (const string &GroupName, const Handle &List) |
Writes an Epetra_DistObject to group GroupName. | |
| void | Read (const string &GroupName, Handle &List) |
Reads an Epetra_DistObject from group GroupName. | |
| void | ReadHandleProperties (const string &GroupName, string &Type, int &NumGlobalElements) |
| Reads the global number of elements and type for a generic handle object. |
Introduction
Class HDF5 allows to read and write using the HDF5 parallel binary data format. HDF5 has the following advantages:
The class supports I/O for the following distributed Epetra objects:
The class also supports some non-distributed types:
Some utility methods are provided:
true is the specified group is already contained in the file;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; this class, instead, is only structured around the concept of groups. A group is an entity, like 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.
The HDF5 class has the following limitations:
Errors
When an error occurs, a EpetraExt::Exception is thrown. Method Print() of the Exception class gives a description of what went wrong.
Example of usage
First, one has to create an HDF5 class, then either Open() or Create() the file:
EpetraExt::HDF5 HDF5(Comm); 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);
HDF5.Write("matrix", "quadrature order", 3);
"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 are supposed to have the same value on all processors.
Reading data is a easy as writing. Let us consider how to read an Epetra_CrsMatrix, other Epetra objects having a similar bejavior. Method ReadCrsMatrixProperties() 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. The actual reading is:
Epetra_CrsMatrix* NewMatrix = 0; HDF5.Read("matrix", NewMatrix);
In this case, NewMatrix is based on a linear map. If the DomainMap() and RangeMap() are known and non-trivial, one can use
HDF5.Read("matrix", DomainMap, RangeMap, NewMatrix);
Reading meta-data 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
Reading HDF5 files from MATLAB is very, since the built-in functions hdf5read, hdf5write, hdf5info. For example, to read the above Matrix from MATLAB, one can do:
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);
double() is required by sparse, which does not accept int32 data.
To dump on file matlab.h5 a MATLAB matrix (in this case, A), one can proceed 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');
__type__ specification, that must reflect the Epetra class name.
To dump on file matlab.h5 a MATLAB dense array (in this case, x), one can proceed 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');
To write a Map from MATLAB follows a very similar pattern. The following example shows how to define a map that can be used with 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');
Definition at line 288 of file EpetraExt_HDF5.h.
| EpetraExt::HDF5::HDF5 | ( | const Epetra_Comm & | Comm | ) |
| EpetraExt::HDF5::~HDF5 | ( | ) | [inline] |
| void EpetraExt::HDF5::Create | ( | const string | FileName | ) |
Creates a new file.
| void EpetraExt::HDF5::Open | ( | const string | FileName, | |
| int | AccessType = H5F_ACC_RDWR | |||
| ) |
Opens specified file with given access type.
| void EpetraExt::HDF5::Close | ( | ) | [inline] |
| bool EpetraExt::HDF5::IsOpen | ( | ) | const [inline] |
Returns true is a file has already been open using Open()/Create().
Definition at line 319 of file EpetraExt_HDF5.h.
| void EpetraExt::HDF5::CreateGroup | ( | const string & | GroupName | ) | [inline] |
| bool EpetraExt::HDF5::IsContained | ( | const string | Name | ) |
Returns true is Name is contained in the database.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| int | data | |||
| ) |
Writes an integer in group GroupName using intentified DataSetName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| int & | data | |||
| ) |
Reads an integer from group /GroupName/DataSetName.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| double | data | |||
| ) |
Writes a double in group GroupName using intentified DataSetName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| double & | data | |||
| ) |
Reads a double from group /GroupName/DataSetName.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| const string & | data | |||
| ) |
Writes a string in group GroupName using intentified DataSetName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| string & | data | |||
| ) |
Reads a string from group /GroupName/DataSetName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| const int | type, | |||
| const int | Length, | |||
| void * | data | |||
| ) |
Reads serial array data, of type type, from group GroupNameusing dataset name DataSetName.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| const int | type, | |||
| const int | Length, | |||
| void * | data | |||
| ) |
Writes serial array data, of type type, to group GroupNameusing dataset name DataSetName.
| void EpetraExt::HDF5::WriteComment | ( | const string & | GroupName, | |
| string | Comment | |||
| ) | [inline] |
| void EpetraExt::HDF5::ReadComment | ( | const string & | GroupName, | |
| string & | Comment | |||
| ) | [inline] |
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| int | MySize, | |||
| int | GlobalSize, | |||
| int | type, | |||
| const void * | data | |||
| ) |
Writes distributed array data, of type type, to group GroupNameusing dataset name DataSetName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const string & | DataSetName, | |||
| int | MySize, | |||
| int | GlobalSize, | |||
| const int | type, | |||
| void * | data | |||
| ) |
Reads distributed array data, of type type, from group GroupNameusing dataset name DataSetName.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Epetra_Map & | Map | |||
| ) |
Writes a Map to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Epetra_Map *& | Map | |||
| ) |
Reads a map from GroupName.
| void EpetraExt::HDF5::ReadMapProperties | ( | const string & | GroupName, | |
| int & | NumGlobalElements, | |||
| int & | IndexBase, | |||
| int & | NumProc | |||
| ) |
Reads basic properties of specified Epetra_Map.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Epetra_BlockMap *& | Map | |||
| ) |
Reads a block map from GroupName.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Epetra_BlockMap & | Map | |||
| ) |
Writes a block map to group GroupName.
| void EpetraExt::HDF5::ReadBlockMapProperties | ( | const string & | GroupName, | |
| int & | NumGlobalElements, | |||
| int & | NumGlobalPoints, | |||
| int & | IndexBase, | |||
| int & | NumProc | |||
| ) |
Reads basic properties of specified Epetra_BlockMap.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Epetra_CrsGraph *& | Graph | |||
| ) |
Reads a vector from group GroupName, assumes linear distribution.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const Epetra_Map & | DomainMap, | |||
| const Epetra_Map & | RangeMap, | |||
| Epetra_CrsGraph *& | Graph | |||
| ) |
Reads a vector from group GroupName using given map.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Epetra_CrsGraph & | Graph | |||
| ) |
Writes a distributed vector to group GroupName.
| void EpetraExt::HDF5::ReadCrsGraphProperties | ( | const string & | GroupName, | |
| int & | NumGlobalRows, | |||
| int & | NumGlobalCols, | |||
| int & | NumGlobalNonzeros, | |||
| int & | NumGlobalDiagonals, | |||
| int & | MaxNumIndices | |||
| ) |
Reads basic properties of specified Epetra_CrsGraph.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Epetra_IntVector & | x | |||
| ) |
Writes a distributed vector to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Epetra_IntVector *& | X | |||
| ) |
Reads a vector from group GroupName, assumes linear distribution.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const Epetra_Map & | Map, | |||
| Epetra_IntVector *& | X | |||
| ) |
Reads a vector from group GroupName using given map.
| void EpetraExt::HDF5::ReadIntVectorProperties | ( | const string & | GroupName, | |
| int & | GlobalLength | |||
| ) |
Reads basic properties of specified Epetra_IntVector.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Epetra_MultiVector & | x | |||
| ) |
Writes a distributed vector to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Epetra_MultiVector *& | X | |||
| ) |
Reads a vector from group GroupName, assumes linear distribution.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const Epetra_Map & | Map, | |||
| Epetra_MultiVector *& | X | |||
| ) |
Reads a vector from group GroupName using given map.
| void EpetraExt::HDF5::ReadMultiVectorProperties | ( | const string & | GroupName, | |
| int & | GlobalLength, | |||
| int & | NumVectors | |||
| ) |
Reads basic properties of specified Epetra_MultiVector.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Epetra_RowMatrix & | Matrix | |||
| ) |
Writes a distributed RowMatrix to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Epetra_CrsMatrix *& | A | |||
| ) |
Reads a square matrix from group GroupName, assumes linear distribution.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const Epetra_Map & | DomainMap, | |||
| const Epetra_Map & | RangeMap, | |||
| Epetra_CrsMatrix *& | A | |||
| ) |
Reads a matrix from group GroupName with given range and domain maps.
| void EpetraExt::HDF5::ReadCrsMatrixProperties | ( | const string & | GroupName, | |
| int & | NumGlobalRows, | |||
| int & | NumGlobalCols, | |||
| int & | NumNonzeros, | |||
| int & | NumGlobalDiagonals, | |||
| int & | MaxNumEntries, | |||
| double & | NormOne, | |||
| double & | NormInf | |||
| ) |
Reads basic properties of specified Epetra_CrsMatrix.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Teuchos::ParameterList & | List | |||
| ) |
Writes a parameter list to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Teuchos::ParameterList & | List | |||
| ) |
Reads a parameter list from group GroupName.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const DistArray< int > & | array | |||
| ) |
Writes an EpetraExt::DistArray<int> to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| DistArray< int > *& | array | |||
| ) |
Reads an EpetraExt::DistArray<int> from group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const Epetra_Map & | Map, | |||
| DistArray< int > *& | array | |||
| ) |
Reads an EpetraExt::DistArray<int> from group GroupName.
| void EpetraExt::HDF5::ReadIntDistArrayProperties | ( | const string & | GroupName, | |
| int & | GlobalLength, | |||
| int & | RowSize | |||
| ) |
Reads the global number of elements and type for a generic handle object.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const DistArray< double > & | array | |||
| ) |
Writes an EpetraExt::DistArray<int> to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| DistArray< double > *& | array | |||
| ) |
Reads an EpetraExt::DistArray<int> from group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| const Epetra_Map & | Map, | |||
| DistArray< double > *& | array | |||
| ) |
Reads an EpetraExt::DistArray<int> from group GroupName.
| void EpetraExt::HDF5::ReadDoubleDistArrayProperties | ( | const string & | GroupName, | |
| int & | GlobalLength, | |||
| int & | RowSize | |||
| ) |
Reads the global number of elements and type for a generic handle object.
| void EpetraExt::HDF5::Write | ( | const string & | GroupName, | |
| const Handle & | List | |||
| ) |
Writes an Epetra_DistObject to group GroupName.
| void EpetraExt::HDF5::Read | ( | const string & | GroupName, | |
| Handle & | List | |||
| ) |
Reads an Epetra_DistObject from group GroupName.
| void EpetraExt::HDF5::ReadHandleProperties | ( | const string & | GroupName, | |
| string & | Type, | |||
| int & | NumGlobalElements | |||
| ) |
Reads the global number of elements and type for a generic handle object.
1.4.7