42 #ifndef EPETRAEXT_HDF5_H 43 #define EPETRAEXT_HDF5_H 46 #ifdef HAVE_EPETRAEXT_HDF5 51 class Epetra_BlockMap;
53 class Epetra_IntVector;
54 class Epetra_MultiVector;
55 class Epetra_CrsGraph;
56 class Epetra_RowMatrix;
57 class Epetra_CrsMatrix;
58 class Epetra_VbrMatrix;
332 HDF5(
const Epetra_Comm& Comm);
345 void Create(
const std::string FileName);
348 void Open(
const std::string FileName,
int AccessType = H5F_ACC_RDWR);
360 H5Fflush(file_id_, H5F_SCOPE_GLOBAL);
372 hid_t group_id = H5Gcreate(file_id_, GroupName.c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
377 bool IsContained(
const std::string Name);
383 void Write(
const std::string& GroupName,
const std::string& DataSetName,
int data);
386 void Read(
const std::string& GroupName,
const std::string& DataSetName,
int& data);
389 void Write(
const std::string& GroupName,
const std::string& DataSetName,
double data);
392 void Read(
const std::string& GroupName,
const std::string& DataSetName,
double& data);
395 void Write(
const std::string& GroupName,
const std::string& DataSetName,
const std::string& data);
398 void Read(
const std::string& GroupName,
const std::string& DataSetName, std::string& data);
401 void Read(
const std::string& GroupName,
const std::string& DataSetName,
402 const int type,
const int Length,
void* data);
405 void Write(
const std::string& GroupName,
const std::string& DataSetName,
406 const int type,
const int Length,
412 H5Gset_comment(file_id_, GroupName.c_str(), Comment.c_str());
416 void ReadComment(
const std::string& GroupName, std::string& Comment)
419 H5Gget_comment(file_id_, GroupName.c_str(), 128, comment);
427 void Write(
const std::string& GroupName,
const std::string& DataSetName,
int MySize,
int GlobalSize,
int type,
const void* data);
430 void Read(
const std::string& GroupName,
const std::string& DataSetName,
431 int MySize,
int GlobalSize,
432 const int type,
void* data);
438 void Write(
const std::string& GroupName,
const Epetra_Map& Map);
441 void Read(
const std::string& GroupName, Epetra_Map*& Map);
444 void ReadMapProperties(
const std::string& GroupName,
445 int& NumGlobalElements,
450 void Read(
const std::string& GroupName, Epetra_BlockMap*& Map);
453 void Write(
const std::string& GroupName,
const Epetra_BlockMap& Map);
456 void ReadBlockMapProperties(
const std::string& GroupName,
457 int& NumGlobalElements,
458 int& NumGlobalPoints,
466 void Read(
const std::string& GroupName, Epetra_CrsGraph*& Graph);
469 void Read(
const std::string& GroupName,
const Epetra_Map& DomainMap,
470 const Epetra_Map& RangeMap, Epetra_CrsGraph*& Graph);
473 void Write(
const std::string& GroupName,
const Epetra_CrsGraph& Graph);
476 void ReadCrsGraphProperties(
const std::string& GroupName,
479 int& NumGlobalNonzeros,
480 int& NumGlobalDiagonals,
487 void Write(
const std::string& GroupName,
const Epetra_IntVector& x);
490 void Read(
const std::string& GroupName, Epetra_IntVector*& X);
493 void Read(
const std::string& GroupName,
const Epetra_Map& Map, Epetra_IntVector*& X);
496 void ReadIntVectorProperties(
const std::string& GroupName,
int& GlobalLength);
504 void Write(
const std::string& GroupName,
const Epetra_MultiVector& x,
bool writeTranspose =
false);
510 void Read(
const std::string& GroupName, Epetra_MultiVector*& X,
511 bool writeTranspose =
false,
const int& indexBase = 0);
516 void Read(
const std::string& GroupName,
const Epetra_Map& Map, Epetra_MultiVector*& X,
517 bool writeTranspose =
false);
520 void ReadMultiVectorProperties(
const std::string& GroupName,
528 void Write(
const std::string& GroupName,
const Epetra_RowMatrix& Matrix);
531 void Read(
const std::string& GroupName, Epetra_CrsMatrix*& A);
534 void Read(
const std::string& GroupName,
535 const Epetra_Map& DomainMap,
536 const Epetra_Map& RangeMap,
537 Epetra_CrsMatrix*& A);
540 void ReadCrsMatrixProperties(
const std::string& GroupName,
544 int& NumGlobalDiagonals,
553 void Write(
const std::string& GroupName,
const Teuchos::ParameterList& List);
556 void Read(
const std::string& GroupName, Teuchos::ParameterList& List);
562 void Write(
const std::string& GroupName,
const DistArray<int>& array);
568 void Read(
const std::string& GroupName,
const Epetra_Map& Map,
DistArray<int>*& array);
571 void ReadIntDistArrayProperties(
const std::string& GroupName,
585 void Read(
const std::string& GroupName,
const Epetra_Map& Map,
DistArray<double>*& array);
588 void ReadDoubleDistArrayProperties(
const std::string& GroupName,
596 void Write(
const std::string& GroupName,
const Handle& List);
599 void Read(
const std::string& GroupName,
Handle& List);
602 void ReadHandleProperties(
const std::string& GroupName,
604 int& NumGlobalElements);
611 const Epetra_Comm& Comm()
const 617 const Epetra_Comm& Comm_;
619 std::string FileName_;
void ReadComment(const std::string &GroupName, std::string &Comment)
Read the string associated with group GroupName.
void CreateGroup(const std::string &GroupName)
Create group GroupName.
bool IsOpen() const
Return true if a file has already been opened using Open()/Create()
void WriteComment(const std::string &GroupName, std::string Comment)
Associate string Comment with group GroupName.
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.
void Flush()
Flush the content to the file.
class HDF5: A class for storing Epetra objects in parallel binary files
void Close()
Close the file.
DistArray<T>: A class to store row-oriented multivectors of type T.