00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef EPETRAEXT_HDF5_H
00030 #define EPETRAEXT_HDF5_H
00031
00032 #include "EpetraExt_ConfigDefs.h"
00033 #ifdef HAVE_EPETRAEXT_HDF5
00034
00035 #include "hdf5.h"
00036 class Epetra_Map;
00037 class Epetra_BlockMap;
00038 class Epetra_Comm;
00039 class Epetra_IntVector;
00040 class Epetra_MultiVector;
00041 class Epetra_CrsGraph;
00042 class Epetra_RowMatrix;
00043 class Epetra_CrsMatrix;
00044 class Epetra_VbrMatrix;
00045 namespace Teuchos {
00046 class ParameterList;
00047 }
00048 namespace EpetraExt {
00049 class Handle;
00050 template<class T>
00051 class DistArray;
00052 }
00053
00054 namespace EpetraExt
00055 {
00288 class HDF5
00289 {
00290 public:
00291
00293 HDF5(const Epetra_Comm& Comm);
00294
00296 ~HDF5()
00297 {
00298 if (IsOpen())
00299 Close();
00300 }
00301
00302
00303
00304
00306 void Create(const string FileName);
00307
00309 void Open(const string FileName, int AccessType = H5F_ACC_RDWR);
00310
00312 void Close()
00313 {
00314 H5Fclose(file_id_);
00315 IsOpen_ = false;
00316 }
00317
00319 bool IsOpen() const
00320 {
00321 return(IsOpen_);
00322 }
00323
00325 void CreateGroup(const string& GroupName)
00326 {
00327 hid_t group_id = H5Gcreate(file_id_, GroupName.c_str(), 0);
00328 H5Gclose(group_id);
00329 }
00330
00332 bool IsContained(const string Name);
00333
00334
00335
00336
00338 void Write(const string& GroupName, const string& DataSetName, int data);
00339
00341 void Read(const string& GroupName, const string& DataSetName, int& data);
00342
00344 void Write(const string& GroupName, const string& DataSetName, double data);
00345
00347 void Read(const string& GroupName, const string& DataSetName, double& data);
00348
00350 void Write(const string& GroupName, const string& DataSetName, const string& data);
00351
00353 void Read(const string& GroupName, const string& DataSetName, string& data);
00354
00356 void Read(const string& GroupName, const string& DataSetName,
00357 const int type, const int Length, void* data);
00358
00360 void Write(const string& GroupName, const string& DataSetName,
00361 const int type, const int Length,
00362 void* data);
00363
00365 void WriteComment(const string& GroupName, string Comment)
00366 {
00367 H5Gset_comment(file_id_, GroupName.c_str(), Comment.c_str());
00368 }
00369
00371 void ReadComment(const string& GroupName, string& Comment)
00372 {
00373 char comment[128];
00374 H5Gget_comment(file_id_, GroupName.c_str(), 128, comment);
00375 Comment = comment;
00376 }
00377
00378
00379
00380
00382 void Write(const string& GroupName, const string& DataSetName, int MySize, int GlobalSize, int type, const void* data);
00383
00385 void Read(const string& GroupName, const string& DataSetName,
00386 int MySize, int GlobalSize,
00387 const int type, void* data);
00388
00389
00390
00391
00393 void Write(const string& GroupName, const Epetra_Map& Map);
00394
00396 void Read(const string& GroupName, Epetra_Map*& Map);
00397
00399 void ReadMapProperties(const string& GroupName,
00400 int& NumGlobalElements,
00401 int& IndexBase,
00402 int& NumProc);
00403
00405 void Read(const string& GroupName, Epetra_BlockMap*& Map);
00406
00408 void Write(const string& GroupName, const Epetra_BlockMap& Map);
00409
00411 void ReadBlockMapProperties(const string& GroupName,
00412 int& NumGlobalElements,
00413 int& NumGlobalPoints,
00414 int& IndexBase,
00415 int& NumProc);
00416
00417
00418
00419
00421 void Read(const string& GroupName, Epetra_CrsGraph*& Graph);
00422
00424 void Read(const string& GroupName, const Epetra_Map& DomainMap,
00425 const Epetra_Map& RangeMap, Epetra_CrsGraph*& Graph);
00426
00428 void Write(const string& GroupName, const Epetra_CrsGraph& Graph);
00429
00431 void ReadCrsGraphProperties(const string& GroupName,
00432 int& NumGlobalRows,
00433 int& NumGlobalCols,
00434 int& NumGlobalNonzeros,
00435 int& NumGlobalDiagonals,
00436 int& MaxNumIndices);
00437
00438
00439
00440
00442 void Write(const string& GroupName, const Epetra_IntVector& x);
00443
00445 void Read(const string& GroupName, Epetra_IntVector*& X);
00446
00448 void Read(const string& GroupName, const Epetra_Map& Map, Epetra_IntVector*& X);
00449
00451 void ReadIntVectorProperties(const string& GroupName, int& GlobalLength);
00452
00453
00454
00455
00457 void Write(const string& GroupName, const Epetra_MultiVector& x, bool writeTranspose = false);
00458
00460 void Read(const string& GroupName, Epetra_MultiVector*& X);
00461
00463 void Read(const string& GroupName, const Epetra_Map& Map, Epetra_MultiVector*& X);
00464
00466 void ReadMultiVectorProperties(const string& GroupName,
00467 int& GlobalLength,
00468 int& NumVectors);
00469
00470
00471
00472
00474 void Write(const string& GroupName, const Epetra_RowMatrix& Matrix);
00475
00477 void Read(const string& GroupName, Epetra_CrsMatrix*& A);
00478
00480 void Read(const string& GroupName,
00481 const Epetra_Map& DomainMap,
00482 const Epetra_Map& RangeMap,
00483 Epetra_CrsMatrix*& A);
00484
00486 void ReadCrsMatrixProperties(const string& GroupName,
00487 int& NumGlobalRows,
00488 int& NumGlobalCols,
00489 int& NumNonzeros,
00490 int& NumGlobalDiagonals,
00491 int& MaxNumEntries,
00492 double& NormOne,
00493 double& NormInf);
00494
00495
00496
00497
00499 void Write(const string& GroupName, const Teuchos::ParameterList& List);
00500
00502 void Read(const string& GroupName, Teuchos::ParameterList& List);
00503
00504
00505
00506
00508 void Write(const string& GroupName, const DistArray<int>& array);
00509
00511 void Read(const string& GroupName, DistArray<int>*& array);
00512
00514 void Read(const string& GroupName, const Epetra_Map& Map, DistArray<int>*& array);
00515
00517 void ReadIntDistArrayProperties(const string& GroupName,
00518 int& GlobalLength,
00519 int& RowSize);
00520
00521
00522
00523
00525 void Write(const string& GroupName, const DistArray<double>& array);
00526
00528 void Read(const string& GroupName, DistArray<double>*& array);
00529
00531 void Read(const string& GroupName, const Epetra_Map& Map, DistArray<double>*& array);
00532
00534 void ReadDoubleDistArrayProperties(const string& GroupName,
00535 int& GlobalLength,
00536 int& RowSize);
00537
00538
00539
00540
00542 void Write(const string& GroupName, const Handle& List);
00543
00545 void Read(const string& GroupName, Handle& List);
00546
00548 void ReadHandleProperties(const string& GroupName,
00549 string& Type,
00550 int& NumGlobalElements);
00551
00552
00553 private:
00554
00555
00557 const Epetra_Comm& Comm() const
00558 {
00559 return(Comm_);
00560 }
00561
00563 const Epetra_Comm& Comm_;
00565 string FileName_;
00567 bool IsOpen_;
00568
00570 hid_t file_id_;
00571 hid_t plist_id_;
00572 herr_t status;
00573
00574
00575 };
00576 }
00577 #endif
00578 #endif