00001 #include "EpetraExt_ConfigDefs.h"
00002
00003
00004 #ifdef HAVE_EPETRAEXT_HDF5
00005
00006
00007 #ifdef HAVE_MPI
00008 # include "Epetra_MpiComm.h"
00009 # include "mpi.h"
00010 #else
00011 # include "Epetra_SerialComm.h"
00012 #endif
00013
00014 #include "Teuchos_ParameterList.hpp"
00015 #include "Teuchos_RefCountPtr.hpp"
00016 #include "Epetra_Map.h"
00017 #include "Epetra_BlockMap.h"
00018 #include "Epetra_CrsGraph.h"
00019 #include "Epetra_FECrsGraph.h"
00020 #include "Epetra_RowMatrix.h"
00021 #include "Epetra_CrsMatrix.h"
00022 #include "Epetra_FECrsMatrix.h"
00023 #include "Epetra_IntVector.h"
00024 #include "Epetra_MultiVector.h"
00025 #include "Epetra_Import.h"
00026 #include "EpetraExt_Exception.h"
00027 #include "EpetraExt_Utils.h"
00028 #include "EpetraExt_HDF5.h"
00029 #include "EpetraExt_DistArray.h"
00030
00031 #define CHECK_HID(hid_t) \
00032 { if (hid_t < 0) \
00033 throw(EpetraExt::Exception(__FILE__, __LINE__, \
00034 "hid_t is negative")); }
00035
00036 #define CHECK_STATUS(status) \
00037 { if (status < 0) \
00038 throw(EpetraExt::Exception(__FILE__, __LINE__, \
00039 "function H5Giterater returned a negative value")); }
00040
00041
00042
00043 struct FindDataset_t
00044 {
00045 std::string name;
00046 bool found;
00047 };
00048
00049 static herr_t FindDataset(hid_t loc_id, const char *name, void *opdata)
00050 {
00051 std::string& token = ((FindDataset_t*)opdata)->name;
00052 if (token == name)
00053 ((FindDataset_t*)opdata)->found = true;
00054
00055 return(0);
00056 }
00057
00058
00059
00060 static void WriteParameterListRecursive(const Teuchos::ParameterList& params,
00061 hid_t group_id)
00062 {
00063 Teuchos::ParameterList::ConstIterator it = params.begin();
00064 for (; it != params.end(); ++ it)
00065 {
00066 std::string key(params.name(it));
00067 if (params.isSublist(key))
00068 {
00069
00070
00071
00072 hid_t child_group_id = H5Gcreate(group_id, key.c_str(), 0);
00073 WriteParameterListRecursive(params.sublist(key), child_group_id);
00074 H5Gclose(child_group_id);
00075 }
00076 else
00077 {
00078
00079
00080
00081
00082
00083 herr_t status;
00084 hsize_t one = 1;
00085 hid_t dataspace_id, dataset_id;
00086 bool found = false;
00087
00088
00089 if (params.isType<std::string>(key))
00090 {
00091 std::string value = params.get<std::string>(key);
00092 hsize_t len = value.size() + 1;
00093 dataspace_id = H5Screate_simple(1, &len, NULL);
00094 dataset_id = H5Dcreate(group_id, key.c_str(), H5T_C_S1, dataspace_id, H5P_DEFAULT);
00095 status = H5Dwrite(dataset_id, H5T_C_S1, H5S_ALL, H5S_ALL,
00096 H5P_DEFAULT, value.c_str());
00097 CHECK_STATUS(status);
00098 status = H5Dclose(dataset_id);
00099 CHECK_STATUS(status);
00100 status = H5Sclose(dataspace_id);
00101 CHECK_STATUS(status);
00102 found = true;
00103 }
00104
00105 if (params.isType<bool>(key))
00106 {
00107
00108 unsigned short value = params.get<bool>(key) ? 1 : 0;
00109 dataspace_id = H5Screate_simple(1, &one, NULL);
00110 dataset_id = H5Dcreate(group_id, key.c_str(), H5T_NATIVE_USHORT, dataspace_id, H5P_DEFAULT);
00111 status = H5Dwrite(dataset_id, H5T_NATIVE_USHORT, H5S_ALL, H5S_ALL,
00112 H5P_DEFAULT, &value);
00113 CHECK_STATUS(status);
00114 status = H5Dclose(dataset_id);
00115 CHECK_STATUS(status);
00116 status = H5Sclose(dataspace_id);
00117 CHECK_STATUS(status);
00118 found = true;
00119 }
00120
00121 if (params.isType<int>(key))
00122 {
00123 int value = params.get<int>(key);
00124 dataspace_id = H5Screate_simple(1, &one, NULL);
00125 dataset_id = H5Dcreate(group_id, key.c_str(), H5T_NATIVE_INT, dataspace_id, H5P_DEFAULT);
00126 status = H5Dwrite(dataset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL,
00127 H5P_DEFAULT, &value);
00128 CHECK_STATUS(status);
00129 status = H5Dclose(dataset_id);
00130 CHECK_STATUS(status);
00131 status = H5Sclose(dataspace_id);
00132 CHECK_STATUS(status);
00133 found = true;
00134 }
00135
00136 if (params.isType<double>(key))
00137 {
00138 double value = params.get<double>(key);
00139 dataspace_id = H5Screate_simple(1, &one, NULL);
00140 dataset_id = H5Dcreate(group_id, key.c_str(), H5T_NATIVE_DOUBLE, dataspace_id, H5P_DEFAULT);
00141 status = H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL,
00142 H5P_DEFAULT, &value);
00143 CHECK_STATUS(status);
00144 status = H5Dclose(dataset_id);
00145 CHECK_STATUS(status);
00146 status = H5Sclose(dataspace_id);
00147 CHECK_STATUS(status);
00148 found = true;
00149 }
00150
00151 if (!found)
00152 {
00153 throw(EpetraExt::Exception(__FILE__, __LINE__,
00154 "type for parameter " + key + " not supported"));
00155 }
00156 }
00157 }
00158 }
00159
00160
00161
00162
00163 static herr_t f_operator(hid_t loc_id, const char *name, void *opdata)
00164 {
00165 H5G_stat_t statbuf;
00166 hid_t dataset_id, space_id, type_id;
00167 Teuchos::ParameterList* sublist;
00168 Teuchos::ParameterList* params = (Teuchos::ParameterList*)opdata;
00169
00170
00171
00172
00173
00174 H5Gget_objinfo(loc_id, name, 0, &statbuf);
00175 if (strcmp(name, "__type__") == 0)
00176 return(0);
00177
00178 switch (statbuf.type) {
00179 case H5G_GROUP:
00180 sublist = ¶ms->sublist(name);
00181 H5Giterate(loc_id, name , NULL, f_operator, sublist);
00182 break;
00183 case H5G_DATASET:
00184 hsize_t len;
00185 dataset_id = H5Dopen(loc_id, name);
00186 space_id = H5Dget_space(dataset_id);
00187 if (H5Sget_simple_extent_ndims(space_id) != 1)
00188 throw(EpetraExt::Exception(__FILE__, __LINE__,
00189 "dimensionality of parameters must be 1."));
00190 H5Sget_simple_extent_dims(space_id, &len, NULL);
00191 type_id = H5Dget_type(dataset_id);
00192 if (H5Tequal(type_id, H5T_NATIVE_DOUBLE) > 0) {
00193 double value;
00194 H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, &value);
00195 params->set(name, value);
00196 } else if (H5Tequal(type_id, H5T_NATIVE_INT) > 0) {
00197 int value;
00198 H5Dread(dataset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &value);
00199 params->set(name, value);
00200 } else if (H5Tequal(type_id, H5T_C_S1) > 0) {
00201 char* buf = new char[len];
00202 H5Dread(dataset_id, H5T_C_S1, H5S_ALL, H5S_ALL, H5P_DEFAULT, buf);
00203 params->set(name, std::string(buf));
00204 delete[] buf;
00205 } else if (H5Tequal(type_id, H5T_NATIVE_USHORT) > 0) {
00206 unsigned short value;
00207 H5Dread(dataset_id, H5T_NATIVE_USHORT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &value);
00208 params->set(name, value != 0 ? true : false);
00209 } else {
00210 throw(EpetraExt::Exception(__FILE__, __LINE__,
00211 "unsupported datatype"));
00212 }
00213 H5Tclose(type_id);
00214 H5Sclose(space_id);
00215 H5Dclose(dataset_id);
00216 break;
00217 default:
00218 throw(EpetraExt::Exception(__FILE__, __LINE__,
00219 "unsupported datatype"));
00220 }
00221 return 0;
00222 }
00223
00224
00225 EpetraExt::HDF5::HDF5(const Epetra_Comm& Comm) :
00226 Comm_(Comm),
00227 IsOpen_(false)
00228 {}
00229
00230
00231 void EpetraExt::HDF5::Create(const std::string FileName)
00232 {
00233 if (IsOpen())
00234 throw(Exception(__FILE__, __LINE__,
00235 "an HDF5 is already open, first close the current one",
00236 "using method Close(), then open/create a new one"));
00237
00238 FileName_ = FileName;
00239
00240
00241 plist_id_ = H5Pcreate(H5P_FILE_ACCESS);
00242 #ifdef HAVE_MPI
00243
00244 H5Pset_fapl_mpio(plist_id_, MPI_COMM_WORLD, MPI_INFO_NULL);
00245 #endif
00246
00247 #if 0
00248 unsigned int boh = H5Z_FILTER_MAX;
00249 H5Pset_filter(plist_id_, H5Z_FILTER_DEFLATE, H5Z_FILTER_MAX, 0, &boh);
00250 #endif
00251
00252
00253 file_id_ = H5Fcreate(FileName.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT,
00254 plist_id_);
00255 H5Pclose(plist_id_);
00256
00257 IsOpen_ = true;
00258 }
00259
00260
00261 void EpetraExt::HDF5::Open(const std::string FileName, int AccessType)
00262 {
00263 if (IsOpen())
00264 throw(Exception(__FILE__, __LINE__,
00265 "an HDF5 is already open, first close the current one",
00266 "using method Close(), then open/create a new one"));
00267
00268 FileName_ = FileName;
00269
00270
00271 file_id_ = H5Fopen(FileName.c_str(), AccessType, H5P_DEFAULT);
00272
00273 #ifdef HAVE_MPI
00274
00275
00276
00277 #endif
00278
00279 IsOpen_ = true;
00280 }
00281
00282
00283 bool EpetraExt::HDF5::IsContained(std::string Name)
00284 {
00285 if (!IsOpen())
00286 throw(Exception(__FILE__, __LINE__, "no file open yet"));
00287
00288 FindDataset_t data;
00289 data.name = Name;
00290 data.found = false;
00291
00292 int idx_f = H5Giterate(file_id_, "/", NULL, FindDataset, (void*)&data);
00293
00294 return(data.found);
00295 }
00296
00297
00298
00299
00300
00301
00302 void EpetraExt::HDF5::Write(const std::string& Name, const Epetra_BlockMap& BlockMap)
00303 {
00304 int NumMyPoints = BlockMap.NumMyPoints();
00305 int NumMyElements = BlockMap.NumMyElements();
00306 int NumGlobalElements = BlockMap.NumGlobalElements();
00307 int NumGlobalPoints = BlockMap.NumGlobalPoints();
00308 int* MyGlobalElements = BlockMap.MyGlobalElements();
00309 int* ElementSizeList = BlockMap.ElementSizeList();
00310
00311 std::vector<int> NumMyElements_v(Comm_.NumProc());
00312 Comm_.GatherAll(&NumMyElements, &NumMyElements_v[0], 1);
00313
00314 std::vector<int> NumMyPoints_v(Comm_.NumProc());
00315 Comm_.GatherAll(&NumMyPoints, &NumMyPoints_v[0], 1);
00316
00317 Write(Name, "MyGlobalElements", NumMyElements, NumGlobalElements,
00318 H5T_NATIVE_INT, MyGlobalElements);
00319 Write(Name, "ElementSizeList", NumMyElements, NumGlobalElements,
00320 H5T_NATIVE_INT, ElementSizeList);
00321 Write(Name, "NumMyPoints", H5T_NATIVE_INT, Comm_.NumProc(), &NumMyPoints_v[0]);
00322
00323
00324 Write(Name, "NumProc", Comm_.NumProc());
00325
00326 Write(Name, "NumGlobalPoints", 1, Comm_.NumProc(), H5T_NATIVE_INT, &NumGlobalPoints);
00327 Write(Name, "NumGlobalElements", 1, Comm_.NumProc(), H5T_NATIVE_INT, &NumGlobalElements);
00328 Write(Name, "IndexBase", BlockMap.IndexBase());
00329 Write(Name, "__type__", "Epetra_BlockMap");
00330 }
00331
00332
00333 void EpetraExt::HDF5::Read(const std::string& GroupName, Epetra_BlockMap*& BlockMap)
00334 {
00335 int NumGlobalElements, NumGlobalPoints, IndexBase, NumProc;
00336
00337 ReadBlockMapProperties(GroupName, NumGlobalElements, NumGlobalPoints,
00338 IndexBase, NumProc);
00339
00340 std::vector<int> NumMyPoints_v(Comm_.NumProc());
00341 std::vector<int> NumMyElements_v(Comm_.NumProc());
00342
00343 Read(GroupName, "NumMyElements", H5T_NATIVE_INT, Comm_.NumProc(), &NumMyElements_v[0]);
00344 Read(GroupName, "NumMyPoints", H5T_NATIVE_INT, Comm_.NumProc(), &NumMyPoints_v[0]);
00345 int NumMyElements = NumMyElements_v[Comm_.MyPID()];
00346 int NumMyPoints = NumMyPoints_v[Comm_.MyPID()];
00347
00348 if (NumProc != Comm_.NumProc())
00349 throw(Exception(__FILE__, __LINE__,
00350 "requested map not compatible with current number of",
00351 "processors, " + toString(Comm_.NumProc()) +
00352 " vs. " + toString(NumProc)));
00353
00354 std::vector<int> MyGlobalElements(NumMyElements);
00355 std::vector<int> ElementSizeList(NumMyElements);
00356
00357 Read(GroupName, "MyGlobalElements", NumMyElements, NumGlobalElements,
00358 H5T_NATIVE_INT, &MyGlobalElements[0]);
00359
00360 Read(GroupName, "ElementSizeList", NumMyElements, NumGlobalElements,
00361 H5T_NATIVE_INT, &ElementSizeList[0]);
00362
00363 BlockMap = new Epetra_BlockMap(NumGlobalElements, NumMyElements,
00364 &MyGlobalElements[0], &ElementSizeList[0],
00365 IndexBase, Comm_);
00366 }
00367
00368
00369 void EpetraExt::HDF5::ReadBlockMapProperties(const std::string& GroupName,
00370 int& NumGlobalElements,
00371 int& NumGlobalPoints,
00372 int& IndexBase,
00373 int& NumProc)
00374 {
00375 if (!IsContained(GroupName))
00376 throw(Exception(__FILE__, __LINE__,
00377 "requested group " + GroupName + " not found"));
00378
00379 std::string Label;
00380 Read(GroupName, "__type__", Label);
00381
00382 if (Label != "Epetra_BlockMap")
00383 throw(Exception(__FILE__, __LINE__,
00384 "requested group " + GroupName + " is not an Epetra_BlockMap",
00385 "__type__ = " + Label));
00386
00387 Read(GroupName, "NumGlobalElements", NumGlobalElements);
00388 Read(GroupName, "NumGlobalPoints", NumGlobalPoints);
00389 Read(GroupName, "IndexBase", IndexBase);
00390 Read(GroupName, "NumProc", NumProc);
00391 }
00392
00393
00394 void EpetraExt::HDF5::Write(const std::string& Name, const Epetra_Map& Map)
00395 {
00396 int MySize = Map.NumMyElements();
00397 int GlobalSize = Map.NumGlobalElements();
00398 int* MyGlobalElements = Map.MyGlobalElements();
00399
00400 std::vector<int> NumMyElements(Comm_.NumProc());
00401 Comm_.GatherAll(&MySize, &NumMyElements[0], 1);
00402
00403 Write(Name, "MyGlobalElements", MySize, GlobalSize,
00404 H5T_NATIVE_INT, MyGlobalElements);
00405 Write(Name, "NumMyElements", H5T_NATIVE_INT, Comm_.NumProc(), &NumMyElements[0]);
00406 Write(Name, "NumGlobalElements", 1, Comm_.NumProc(), H5T_NATIVE_INT, &GlobalSize);
00407 Write(Name, "NumProc", Comm_.NumProc());
00408 Write(Name, "IndexBase", Map.IndexBase());
00409 Write(Name, "__type__", "Epetra_Map");
00410 }
00411
00412
00413 void EpetraExt::HDF5::Read(const std::string& GroupName, Epetra_Map*& Map)
00414 {
00415 int NumGlobalElements, IndexBase, NumProc;
00416
00417 ReadMapProperties(GroupName, NumGlobalElements, IndexBase, NumProc);
00418
00419 std::vector<int> NumMyElements_v(Comm_.NumProc());
00420
00421 Read(GroupName, "NumMyElements", H5T_NATIVE_INT, Comm_.NumProc(), &NumMyElements_v[0]);
00422 int NumMyElements = NumMyElements_v[Comm_.MyPID()];
00423
00424 if (NumProc != Comm_.NumProc())
00425 throw(Exception(__FILE__, __LINE__,
00426 "requested map not compatible with current number of",
00427 "processors, " + toString(Comm_.NumProc()) +
00428 " vs. " + toString(NumProc)));
00429
00430 std::vector<int> MyGlobalElements(NumMyElements);
00431
00432 Read(GroupName, "MyGlobalElements", NumMyElements, NumGlobalElements,
00433 H5T_NATIVE_INT, &MyGlobalElements[0]);
00434
00435 Map = new Epetra_Map(-1, NumMyElements, &MyGlobalElements[0], IndexBase, Comm_);
00436 }
00437
00438
00439 void EpetraExt::HDF5::ReadMapProperties(const std::string& GroupName,
00440 int& NumGlobalElements,
00441 int& IndexBase,
00442 int& NumProc)
00443 {
00444 if (!IsContained(GroupName))
00445 throw(Exception(__FILE__, __LINE__,
00446 "requested group " + GroupName + " not found"));
00447
00448 std::string Label;
00449 Read(GroupName, "__type__", Label);
00450
00451 if (Label != "Epetra_Map")
00452 throw(Exception(__FILE__, __LINE__,
00453 "requested group " + GroupName + " is not an Epetra_Map",
00454 "__type__ = " + Label));
00455
00456 Read(GroupName, "NumGlobalElements", NumGlobalElements);
00457 Read(GroupName, "IndexBase", IndexBase);
00458 Read(GroupName, "NumProc", NumProc);
00459 }
00460
00461
00462
00463
00464
00465
00466 void EpetraExt::HDF5::Write(const std::string& Name, const Epetra_IntVector& x)
00467 {
00468 if (x.Map().LinearMap())
00469 {
00470 Write(Name, "GlobalLength", x.GlobalLength());
00471 Write(Name, "Values", x.Map().NumMyElements(), x.Map().NumGlobalElements(),
00472 H5T_NATIVE_INT, x.Values());
00473 }
00474 else
00475 {
00476
00477
00478 const Epetra_BlockMap& OriginalMap = x.Map();
00479 Epetra_Map LinearMap(OriginalMap.NumGlobalElements(), OriginalMap.IndexBase(), Comm_);
00480 Epetra_Import Importer(LinearMap, OriginalMap);
00481
00482 Epetra_IntVector LinearX(LinearMap);
00483 LinearX.Import(x, Importer, Insert);
00484
00485 Write(Name, "GlobalLength", x.GlobalLength());
00486 Write(Name, "Values", LinearMap.NumMyElements(), LinearMap.NumGlobalElements(),
00487 H5T_NATIVE_INT, LinearX.Values());
00488 }
00489 Write(Name, "__type__", "Epetra_IntVector");
00490 }
00491
00492
00493 void EpetraExt::HDF5::Read(const std::string& GroupName, Epetra_IntVector*& X)
00494 {
00495
00496 int GlobalLength;
00497 ReadIntVectorProperties(GroupName, GlobalLength);
00498
00499
00500 Epetra_Map LinearMap(GlobalLength, 0, Comm_);
00501 X = new Epetra_IntVector(LinearMap);
00502
00503 Read(GroupName, "Values", LinearMap.NumMyElements(),
00504 LinearMap.NumGlobalElements(), H5T_NATIVE_INT, X->Values());
00505 }
00506
00507
00508 void EpetraExt::HDF5::Read(const std::string& GroupName, const Epetra_Map& Map,
00509 Epetra_IntVector*& X)
00510 {
00511
00512 int GlobalLength;
00513 ReadIntVectorProperties(GroupName, GlobalLength);
00514
00515 if (Map.LinearMap())
00516 {
00517 X = new Epetra_IntVector(Map);
00518
00519 Read(GroupName, "Values", Map.NumMyElements(), Map.NumGlobalElements(),
00520 H5T_NATIVE_INT, X->Values());
00521
00522 }
00523 else
00524 {
00525
00526
00527 Epetra_Map LinearMap(GlobalLength, Map.IndexBase(), Comm_);
00528 Epetra_IntVector LinearX(LinearMap);
00529
00530 Read(GroupName, "Values", LinearMap.NumMyElements(),
00531 LinearMap.NumGlobalElements(),
00532 H5T_NATIVE_INT, LinearX.Values());
00533
00534 Epetra_Import Importer(Map, LinearMap);
00535 X = new Epetra_IntVector(Map);
00536 X->Import(LinearX, Importer, Insert);
00537 }
00538 }
00539
00540
00541 void EpetraExt::HDF5::ReadIntVectorProperties(const std::string& GroupName,
00542 int& GlobalLength)
00543 {
00544 if (!IsContained(GroupName))
00545 throw(Exception(__FILE__, __LINE__,
00546 "requested group " + GroupName + " not found"));
00547
00548 std::string Label;
00549 Read(GroupName, "__type__", Label);
00550
00551 if (Label != "Epetra_IntVector")
00552 throw(Exception(__FILE__, __LINE__,
00553 "requested group " + GroupName + " is not an Epetra_IntVector",
00554 "__type__ = " + Label));
00555
00556 Read(GroupName, "GlobalLength", GlobalLength);
00557 }
00558
00559
00560
00561
00562
00563
00564 void EpetraExt::HDF5::Write(const std::string& Name, const Epetra_CrsGraph& Graph)
00565 {
00566 if (!Graph.Filled())
00567 throw(Exception(__FILE__, __LINE__,
00568 "input Epetra_CrsGraph is not FillComplete()'d"));
00569
00570
00571 int MySize = Graph.NumMyNonzeros();
00572 int GlobalSize = Graph.NumGlobalNonzeros();
00573 std::vector<int> ROW(MySize);
00574 std::vector<int> COL(MySize);
00575
00576 int count = 0;
00577 int* RowIndices;
00578 int NumEntries;
00579
00580 for (int i = 0; i < Graph.NumMyRows(); ++i)
00581 {
00582 Graph.ExtractMyRowView(i, NumEntries, RowIndices);
00583 for (int j = 0; j < NumEntries; ++j)
00584 {
00585 ROW[count] = Graph.GRID(i);
00586 COL[count] = Graph.GCID(RowIndices[j]);
00587 ++count;
00588 }
00589 }
00590
00591 Write(Name, "ROW", MySize, GlobalSize, H5T_NATIVE_INT, &ROW[0]);
00592 Write(Name, "COL", MySize, GlobalSize, H5T_NATIVE_INT, &COL[0]);
00593 Write(Name, "MaxNumIndices", Graph.MaxNumIndices());
00594 Write(Name, "NumGlobalRows", Graph.NumGlobalRows());
00595 Write(Name, "NumGlobalCols", Graph.NumGlobalCols());
00596 Write(Name, "NumGlobalNonzeros", Graph.NumGlobalNonzeros());
00597 Write(Name, "NumGlobalDiagonals", Graph.NumGlobalDiagonals());
00598 Write(Name, "__type__", "Epetra_CrsGraph");
00599 }
00600
00601
00602 void EpetraExt::HDF5::Read(const std::string& GroupName, Epetra_CrsGraph*& Graph)
00603 {
00604 int NumGlobalRows, NumGlobalCols;
00605 int NumGlobalNonzeros, NumGlobalDiagonals, MaxNumIndices;
00606
00607 ReadCrsGraphProperties(GroupName, NumGlobalRows,
00608 NumGlobalCols, NumGlobalNonzeros,
00609 NumGlobalDiagonals, MaxNumIndices);
00610
00611 Epetra_Map RangeMap(NumGlobalRows, 0, Comm_);
00612 Epetra_Map DomainMap(NumGlobalCols, 0, Comm_);
00613
00614 Read(GroupName, DomainMap, RangeMap, Graph);
00615 }
00616
00617
00618 void EpetraExt::HDF5::ReadCrsGraphProperties(const std::string& GroupName,
00619 int& NumGlobalRows,
00620 int& NumGlobalCols,
00621 int& NumGlobalNonzeros,
00622 int& NumGlobalDiagonals,
00623 int& MaxNumIndices)
00624 {
00625 if (!IsContained(GroupName))
00626 throw(Exception(__FILE__, __LINE__,
00627 "requested group " + GroupName + " not found"));
00628
00629 std::string Label;
00630 Read(GroupName, "__type__", Label);
00631
00632 if (Label != "Epetra_CrsGraph")
00633 throw(Exception(__FILE__, __LINE__,
00634 "requested group " + GroupName + " is not an Epetra_CrsGraph",
00635 "__type__ = " + Label));
00636
00637 Read(GroupName, "NumGlobalRows", NumGlobalRows);
00638 Read(GroupName, "NumGlobalCols", NumGlobalCols);
00639 Read(GroupName, "NumGlobalNonzeros", NumGlobalNonzeros);
00640 Read(GroupName, "NumGlobalDiagonals", NumGlobalDiagonals);
00641 Read(GroupName, "MaxNumIndices", MaxNumIndices);
00642 }
00643
00644
00645 void EpetraExt::HDF5::Read(const std::string& GroupName, const Epetra_Map& DomainMap,
00646 const Epetra_Map& RangeMap, Epetra_CrsGraph*& Graph)
00647 {
00648 if (!IsContained(GroupName))
00649 throw(Exception(__FILE__, __LINE__,
00650 "requested group " + GroupName + " not found in database"));
00651
00652 std::string Label;
00653 Read(GroupName, "__type__", Label);
00654
00655 if (Label != "Epetra_CrsGraph")
00656 throw(Exception(__FILE__, __LINE__,
00657 "requested group " + GroupName + " is not an Epetra_CrsGraph",
00658 "__type__ = " + Label));
00659
00660 int NumGlobalRows, NumGlobalCols, NumGlobalNonzeros;
00661 Read(GroupName, "NumGlobalNonzeros", NumGlobalNonzeros);
00662 Read(GroupName, "NumGlobalRows", NumGlobalRows);
00663 Read(GroupName, "NumGlobalCols", NumGlobalCols);
00664
00665
00666 int NumMyNonzeros = NumGlobalNonzeros / Comm_.NumProc();
00667 if (Comm_.MyPID() == 0)
00668 NumMyNonzeros += NumGlobalNonzeros % Comm_.NumProc();
00669
00670 std::vector<int> ROW(NumMyNonzeros);
00671 Read(GroupName, "ROW", NumMyNonzeros, NumGlobalNonzeros, H5T_NATIVE_INT, &ROW[0]);
00672
00673 std::vector<int> COL(NumMyNonzeros);
00674 Read(GroupName, "COL", NumMyNonzeros, NumGlobalNonzeros, H5T_NATIVE_INT, &COL[0]);
00675
00676 Epetra_FECrsGraph* Graph2 = new Epetra_FECrsGraph(Copy, DomainMap, 0);
00677
00678 for (int i = 0; i < NumMyNonzeros; )
00679 {
00680 int count = 1;
00681 while (ROW[i + count] == ROW[i])
00682 ++count;
00683
00684 Graph2->InsertGlobalIndices(1, &ROW[i], count, &COL[i]);
00685 i += count;
00686 }
00687
00688 Graph2->FillComplete(DomainMap, RangeMap);
00689
00690 Graph = Graph2;
00691 }
00692
00693
00694
00695
00696
00697
00698 void EpetraExt::HDF5::Write(const std::string& Name, const Epetra_RowMatrix& Matrix)
00699 {
00700 if (!Matrix.Filled())
00701 throw(Exception(__FILE__, __LINE__,
00702 "input Epetra_RowMatrix is not FillComplete()'d"));
00703
00704 int MySize = Matrix.NumMyNonzeros();
00705 int GlobalSize = Matrix.NumGlobalNonzeros();
00706 std::vector<int> ROW(MySize);
00707 std::vector<int> COL(MySize);
00708 std::vector<double> VAL(MySize);
00709
00710 int count = 0;
00711 int Length = Matrix.MaxNumEntries();
00712 std::vector<int> Indices(Length);
00713 std::vector<double> Values(Length);
00714 int NumEntries;
00715
00716 for (int i = 0; i < Matrix.NumMyRows(); ++i)
00717 {
00718 Matrix.ExtractMyRowCopy(i, Length, NumEntries, &Values[0], &Indices[0]);
00719 for (int j = 0; j < NumEntries; ++j)
00720 {
00721 ROW[count] = Matrix.RowMatrixRowMap().GID(i);
00722 COL[count] = Matrix.RowMatrixColMap().GID(Indices[j]);
00723 VAL[count] = Values[j];
00724 ++count;
00725 }
00726 }
00727
00728 Write(Name, "ROW", MySize, GlobalSize, H5T_NATIVE_INT, &ROW[0]);
00729 Write(Name, "COL", MySize, GlobalSize, H5T_NATIVE_INT, &COL[0]);
00730 Write(Name, "VAL", MySize, GlobalSize, H5T_NATIVE_DOUBLE, &VAL[0]);
00731 Write(Name, "NumGlobalRows", Matrix.NumGlobalRows());
00732 Write(Name, "NumGlobalCols", Matrix.NumGlobalCols());
00733 Write(Name, "NumGlobalNonzeros", Matrix.NumGlobalNonzeros());
00734 Write(Name, "NumGlobalDiagonals", Matrix.NumGlobalDiagonals());
00735 Write(Name, "MaxNumEntries", Matrix.MaxNumEntries());
00736 Write(Name, "NormOne", Matrix.NormOne());
00737 Write(Name, "NormInf", Matrix.NormInf());
00738 Write(Name, "__type__", "Epetra_RowMatrix");
00739 }
00740
00741
00742 void EpetraExt::HDF5::Read(const std::string& GroupName, Epetra_CrsMatrix*& A)
00743 {
00744 int NumGlobalRows, NumGlobalCols, NumGlobalNonzeros;
00745 int NumGlobalDiagonals, MaxNumEntries;
00746 double NormOne, NormInf;
00747
00748 ReadCrsMatrixProperties(GroupName, NumGlobalRows, NumGlobalCols,
00749 NumGlobalNonzeros, NumGlobalDiagonals, MaxNumEntries,
00750 NormOne, NormInf);
00751
00752
00753 Epetra_Map RangeMap(NumGlobalRows, 0, Comm_);
00754 Epetra_Map DomainMap(NumGlobalCols, 0, Comm_);
00755
00756 Read(GroupName, DomainMap, RangeMap, A);
00757 }
00758
00759
00760 void EpetraExt::HDF5::Read(const std::string& GroupName, const Epetra_Map& DomainMap,
00761 const Epetra_Map& RangeMap, Epetra_CrsMatrix*& A)
00762 {
00763 int NumGlobalRows, NumGlobalCols, NumGlobalNonzeros;
00764 int NumGlobalDiagonals, MaxNumEntries;
00765 double NormOne, NormInf;
00766
00767 ReadCrsMatrixProperties(GroupName, NumGlobalRows, NumGlobalCols,
00768 NumGlobalNonzeros, NumGlobalDiagonals, MaxNumEntries,
00769 NormOne, NormInf);
00770
00771 int NumMyNonzeros = NumGlobalNonzeros / Comm_.NumProc();
00772 if (Comm_.MyPID() == 0)
00773 NumMyNonzeros += NumGlobalNonzeros % Comm_.NumProc();
00774
00775 std::vector<int> ROW(NumMyNonzeros);
00776 Read(GroupName, "ROW", NumMyNonzeros, NumGlobalNonzeros, H5T_NATIVE_INT, &ROW[0]);
00777
00778 std::vector<int> COL(NumMyNonzeros);
00779 Read(GroupName, "COL", NumMyNonzeros, NumGlobalNonzeros, H5T_NATIVE_INT, &COL[0]);
00780
00781 std::vector<double> VAL(NumMyNonzeros);
00782 Read(GroupName, "VAL", NumMyNonzeros, NumGlobalNonzeros, H5T_NATIVE_DOUBLE, &VAL[0]);
00783
00784 Epetra_FECrsMatrix* A2 = new Epetra_FECrsMatrix(Copy, DomainMap, 0);
00785
00786 for (int i = 0; i < NumMyNonzeros; )
00787 {
00788 int count = 1;
00789 while (ROW[i + count] == ROW[i])
00790 ++count;
00791
00792 A2->InsertGlobalValues(1, &ROW[i], count, &COL[i], &VAL[i]);
00793 i += count;
00794 }
00795
00796 A2->FillComplete(DomainMap, RangeMap);
00797
00798 A = A2;
00799 }
00800
00801
00802 void EpetraExt::HDF5::ReadCrsMatrixProperties(const std::string& GroupName,
00803 int& NumGlobalRows,
00804 int& NumGlobalCols,
00805 int& NumGlobalNonzeros,
00806 int& NumGlobalDiagonals,
00807 int& MaxNumEntries,
00808 double& NormOne,
00809 double& NormInf)
00810 {
00811 if (!IsContained(GroupName))
00812 throw(Exception(__FILE__, __LINE__,
00813 "requested group " + GroupName + " not found"));
00814
00815 std::string Label;
00816 Read(GroupName, "__type__", Label);
00817
00818 if (Label != "Epetra_RowMatrix")
00819 throw(Exception(__FILE__, __LINE__,
00820 "requested group " + GroupName + " is not an Epetra_RowMatrix",
00821 "__type__ = " + Label));
00822
00823 Read(GroupName, "NumGlobalRows", NumGlobalRows);
00824 Read(GroupName, "NumGlobalCols", NumGlobalCols);
00825 Read(GroupName, "NumGlobalNonzeros", NumGlobalNonzeros);
00826 Read(GroupName, "NumGlobalDiagonals", NumGlobalDiagonals);
00827 Read(GroupName, "MaxNumEntries", MaxNumEntries);
00828 Read(GroupName, "NormOne", NormOne);
00829 Read(GroupName, "NormInf", NormInf);
00830 }
00831
00832
00833
00834
00835
00836
00837 void EpetraExt::HDF5::Write(const std::string& GroupName, const Teuchos::ParameterList& params)
00838 {
00839 if (!IsOpen())
00840 throw(Exception(__FILE__, __LINE__, "no file open yet"));
00841
00842 hid_t group_id = H5Gcreate(file_id_, GroupName.c_str(), 0);
00843
00844
00845 WriteParameterListRecursive(params, group_id);
00846
00847
00848 status = H5Gclose(group_id);
00849 CHECK_STATUS(status);
00850
00851 Write(GroupName, "__type__", "Teuchos::ParameterList");
00852 }
00853
00854
00855 void EpetraExt::HDF5::Read(const std::string& GroupName, Teuchos::ParameterList& params)
00856 {
00857 if (!IsOpen())
00858 throw(Exception(__FILE__, __LINE__, "no file open yet"));
00859
00860 std::string Label;
00861 Read(GroupName, "__type__", Label);
00862
00863 if (Label != "Teuchos::ParameterList")
00864 throw(Exception(__FILE__, __LINE__,
00865 "requested group " + GroupName + " is not a Teuchos::ParameterList",
00866 "__type__ = " + Label));
00867
00868
00869 hid_t group_id;
00870 herr_t status;
00871
00872
00873 group_id = H5Gopen(file_id_, GroupName.c_str());
00874 CHECK_HID(group_id);
00875
00876
00877 std::string xxx = "/" + GroupName;
00878 status = H5Giterate(group_id, xxx.c_str() , NULL, f_operator, ¶ms);
00879 CHECK_STATUS(status);
00880
00881
00882 status = H5Gclose(group_id);
00883 CHECK_STATUS(status);
00884 }
00885
00886
00887
00888
00889
00890
00891 void EpetraExt::HDF5::Write(const std::string& GroupName, const Epetra_MultiVector& X, bool writeTranspose)
00892 {
00893
00894 if (!IsOpen())
00895 throw(Exception(__FILE__, __LINE__, "no file open yet"));
00896
00897 hid_t group_id, dset_id;
00898 hid_t filespace_id, memspace_id;
00899 herr_t status;
00900
00901
00902 Teuchos::RefCountPtr<Epetra_MultiVector> LinearX;
00903
00904 if (Comm().NumProc() == 1 || X.Map().LinearMap())
00905 LinearX = Teuchos::rcp(const_cast<Epetra_MultiVector*>(&X), false);
00906 else
00907 {
00908 Epetra_Map LinearMap(X.GlobalLength(), X.Map().IndexBase(), Comm_);
00909 LinearX = Teuchos::rcp(new Epetra_MultiVector(LinearMap, X.NumVectors()));
00910 Epetra_Import Importer(LinearMap, X.Map());
00911 LinearX->Import(X, Importer, Insert);
00912 }
00913
00914 int NumVectors = X.NumVectors();
00915 int GlobalLength = X.GlobalLength();
00916
00917
00918
00919
00920 int indexT(0);
00921 if (writeTranspose) indexT = 1;
00922
00923 hsize_t q_dimsf[] = {GlobalLength, GlobalLength};
00924 q_dimsf[indexT] = NumVectors;
00925
00926 filespace_id = H5Screate_simple(2, q_dimsf, NULL);
00927
00928 if (!IsContained(GroupName))
00929 CreateGroup(GroupName);
00930
00931 group_id = H5Gopen(file_id_, GroupName.c_str());
00932
00933
00934 dset_id = H5Dcreate(group_id, "Values", H5T_NATIVE_DOUBLE, filespace_id, H5P_DEFAULT);
00935
00936
00937 plist_id_ = H5Pcreate(H5P_DATASET_XFER);
00938 #ifdef HAVE_MPI
00939 H5Pset_dxpl_mpio(plist_id_, H5FD_MPIO_COLLECTIVE);
00940 #endif
00941
00942
00943
00944 hsize_t offset[] = {LinearX->Map().GID(0)-X.Map().IndexBase(),
00945 LinearX->Map().GID(0)-X.Map().IndexBase()};
00946 hsize_t stride[] = {1, 1};
00947 hsize_t count[] = {LinearX->MyLength(),
00948 LinearX->MyLength()};
00949 hsize_t block[] = {1, 1};
00950
00951
00952 for (int n(0); n < NumVectors; ++n)
00953 {
00954
00955 offset[indexT] = n;
00956 count [indexT] = 1;
00957
00958 filespace_id = H5Dget_space(dset_id);
00959 H5Sselect_hyperslab(filespace_id, H5S_SELECT_SET, offset, stride,
00960 count, block);
00961
00962
00963 hsize_t dimsm[] = {LinearX->MyLength()};
00964 memspace_id = H5Screate_simple(1, dimsm, NULL);
00965
00966
00967 status = H5Dwrite(dset_id, H5T_NATIVE_DOUBLE, memspace_id, filespace_id,
00968 plist_id_, LinearX->operator[](n));
00969 CHECK_STATUS(status);
00970 }
00971 H5Gclose(group_id);
00972 H5Sclose(memspace_id);
00973 H5Sclose(filespace_id);
00974 H5Dclose(dset_id);
00975 H5Pclose(plist_id_);
00976
00977 Write(GroupName, "GlobalLength", GlobalLength);
00978 Write(GroupName, "NumVectors", NumVectors);
00979 Write(GroupName, "__type__", "Epetra_MultiVector");
00980 }
00981
00982
00983 void EpetraExt::HDF5::Read(const std::string& GroupName, const Epetra_Map& Map,
00984 Epetra_MultiVector*& X)
00985 {
00986
00987 Epetra_MultiVector* LinearX;
00988 Read(GroupName, LinearX);
00989
00990
00991 Epetra_Import Importer(Map, X->Map());
00992 X = new Epetra_MultiVector(Map, LinearX->NumVectors());
00993 X->Import(*LinearX, Importer, Insert);
00994
00995
00996 delete LinearX;
00997 }
00998
00999
01000 void EpetraExt::HDF5::Read(const std::string& GroupName, Epetra_MultiVector*& LinearX)
01001 {
01002 int GlobalLength, NumVectors;
01003
01004 ReadMultiVectorProperties(GroupName, GlobalLength, NumVectors);
01005
01006 hid_t group_id;
01007
01008
01009 hsize_t q_dimsf[] = {NumVectors, GlobalLength};
01010 hid_t filespace_id = H5Screate_simple(2, q_dimsf, NULL);
01011
01012 if (!IsContained(GroupName))
01013 throw(Exception(__FILE__, __LINE__,
01014 "requested group " + GroupName + " not found"));
01015
01016 group_id = H5Gopen(file_id_, GroupName.c_str());
01017
01018
01019 hid_t dset_id = H5Dopen(group_id, "Values");
01020
01021
01022 plist_id_ = H5Pcreate(H5P_DATASET_XFER);
01023 #ifdef HAVE_MPI
01024 H5Pset_dxpl_mpio(plist_id_, H5FD_MPIO_COLLECTIVE);
01025 #endif
01026 H5Pclose(plist_id_);
01027
01028 Epetra_Map LinearMap(GlobalLength, 0, Comm());
01029 LinearX = new Epetra_MultiVector(LinearMap, NumVectors);
01030
01031
01032 hsize_t offset[] = {0, LinearMap.GID(0)};
01033 hsize_t stride[] = {1, 1};
01034 hsize_t count[] = {NumVectors, 1};
01035 hsize_t block[] = {1, LinearX->MyLength()};
01036
01037 filespace_id = H5Dget_space(dset_id);
01038 H5Sselect_hyperslab(filespace_id, H5S_SELECT_SET, offset, stride,
01039 count, block);
01040
01041
01042 hsize_t dimsm[] = {NumVectors * LinearX->MyLength()};
01043 hid_t memspace_id = H5Screate_simple(1, dimsm, NULL);
01044
01045
01046 CHECK_STATUS(H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace_id, filespace_id,
01047 H5P_DEFAULT, LinearX->Values()));
01048
01049 CHECK_STATUS(H5Gclose(group_id));
01050 CHECK_STATUS(H5Sclose(memspace_id));
01051 CHECK_STATUS(H5Sclose(filespace_id));
01052 CHECK_STATUS(H5Dclose(dset_id));
01053 }
01054
01055
01056 void EpetraExt::HDF5::ReadMultiVectorProperties(const std::string& GroupName,
01057 int& GlobalLength,
01058 int& NumVectors)
01059 {
01060 if (!IsContained(GroupName))
01061 throw(Exception(__FILE__, __LINE__,
01062 "requested group " + GroupName + " not found"));
01063
01064 std::string Label;
01065 Read(GroupName, "__type__", Label);
01066
01067 if (Label != "Epetra_MultiVector")
01068 throw(Exception(__FILE__, __LINE__,
01069 "requested group " + GroupName + " is not an Epetra_MultiVector",
01070 "__type__ = " + Label));
01071
01072 Read(GroupName, "GlobalLength", GlobalLength);
01073 Read(GroupName, "NumVectors", NumVectors);
01074 }
01075
01076
01077
01078
01079
01080
01081 void EpetraExt::HDF5::Write(const std::string& GroupName, const DistArray<int>& x)
01082 {
01083 if (x.Map().LinearMap())
01084 {
01085 Write(GroupName, "Values", x.Map().NumMyElements() * x.RowSize(),
01086 x.Map().NumGlobalElements() * x.RowSize(),
01087 H5T_NATIVE_INT, x.Values());
01088 }
01089 else
01090 {
01091
01092
01093 const Epetra_BlockMap& OriginalMap = x.Map();
01094 Epetra_Map LinearMap(OriginalMap.NumGlobalElements(), OriginalMap.IndexBase(), Comm_);
01095 Epetra_Import Importer(LinearMap, OriginalMap);
01096
01097 EpetraExt::DistArray<int> LinearX(LinearMap, x.RowSize());
01098 LinearX.Import(x, Importer, Insert);
01099
01100 Write(GroupName, "Values", LinearMap.NumMyElements() * x.RowSize(),
01101 LinearMap.NumGlobalElements() * x.RowSize(),
01102 H5T_NATIVE_INT, LinearX.Values());
01103 }
01104
01105 Write(GroupName, "__type__", "EpetraExt::DistArray<int>");
01106 Write(GroupName, "GlobalLength", x.GlobalLength());
01107 Write(GroupName, "RowSize", x.RowSize());
01108 }
01109
01110
01111 void EpetraExt::HDF5::Read(const std::string& GroupName, const Epetra_Map& Map,
01112 DistArray<int>*& X)
01113 {
01114
01115 int GlobalLength, RowSize;
01116 ReadIntDistArrayProperties(GroupName, GlobalLength, RowSize);
01117
01118 if (Map.LinearMap())
01119 {
01120 X = new EpetraExt::DistArray<int>(Map, RowSize);
01121
01122 Read(GroupName, "Values", Map.NumMyElements() * RowSize,
01123 Map.NumGlobalElements() * RowSize, H5T_NATIVE_INT, X->Values());
01124 }
01125 else
01126 {
01127
01128
01129 Epetra_Map LinearMap(GlobalLength, Map.IndexBase(), Comm_);
01130 EpetraExt::DistArray<int> LinearX(LinearMap, RowSize);
01131
01132 Read(GroupName, "Values", LinearMap.NumMyElements() * RowSize,
01133 LinearMap.NumGlobalElements() * RowSize,
01134 H5T_NATIVE_INT, LinearX.Values());
01135
01136 Epetra_Import Importer(Map, LinearMap);
01137 X = new EpetraExt::DistArray<int>(Map, RowSize);
01138 X->Import(LinearX, Importer, Insert);
01139 }
01140 }
01141
01142
01143 void EpetraExt::HDF5::Read(const std::string& GroupName, DistArray<int>*& X)
01144 {
01145
01146 int GlobalLength, RowSize;
01147 ReadIntDistArrayProperties(GroupName, GlobalLength, RowSize);
01148
01149
01150 Epetra_Map LinearMap(GlobalLength, 0, Comm_);
01151 X = new EpetraExt::DistArray<int>(LinearMap, RowSize);
01152
01153 Read(GroupName, "Values", LinearMap.NumMyElements() * RowSize,
01154 LinearMap.NumGlobalElements() * RowSize, H5T_NATIVE_INT, X->Values());
01155 }
01156
01157
01158 void EpetraExt::HDF5::ReadIntDistArrayProperties(const std::string& GroupName,
01159 int& GlobalLength,
01160 int& RowSize)
01161 {
01162 if (!IsContained(GroupName))
01163 throw(Exception(__FILE__, __LINE__,
01164 "requested group " + GroupName + " not found"));
01165
01166 std::string Label;
01167 Read(GroupName, "__type__", Label);
01168
01169 if (Label != "EpetraExt::DistArray<int>")
01170 throw(Exception(__FILE__, __LINE__,
01171 "requested group " + GroupName + " is not an EpetraExt::DistArray<int>",
01172 "__type__ = " + Label));
01173
01174 Read(GroupName, "GlobalLength", GlobalLength);
01175 Read(GroupName, "RowSize", RowSize);
01176 }
01177
01178
01179
01180
01181
01182
01183 void EpetraExt::HDF5::Write(const std::string& GroupName, const DistArray<double>& x)
01184 {
01185 if (x.Map().LinearMap())
01186 {
01187 Write(GroupName, "Values", x.Map().NumMyElements() * x.RowSize(),
01188 x.Map().NumGlobalElements() * x.RowSize(),
01189 H5T_NATIVE_DOUBLE, x.Values());
01190 }
01191 else
01192 {
01193
01194
01195 const Epetra_BlockMap& OriginalMap = x.Map();
01196 Epetra_Map LinearMap(OriginalMap.NumGlobalElements(), OriginalMap.IndexBase(), Comm_);
01197 Epetra_Import Importer(LinearMap, OriginalMap);
01198
01199 EpetraExt::DistArray<double> LinearX(LinearMap, x.RowSize());
01200 LinearX.Import(x, Importer, Insert);
01201
01202 Write(GroupName, "Values", LinearMap.NumMyElements() * x.RowSize(),
01203 LinearMap.NumGlobalElements() * x.RowSize(),
01204 H5T_NATIVE_DOUBLE, LinearX.Values());
01205 }
01206
01207 Write(GroupName, "__type__", "EpetraExt::DistArray<double>");
01208 Write(GroupName, "GlobalLength", x.GlobalLength());
01209 Write(GroupName, "RowSize", x.RowSize());
01210 }
01211
01212
01213 void EpetraExt::HDF5::Read(const std::string& GroupName, const Epetra_Map& Map,
01214 DistArray<double>*& X)
01215 {
01216
01217 int GlobalLength, RowSize;
01218 ReadDoubleDistArrayProperties(GroupName, GlobalLength, RowSize);
01219
01220 if (Map.LinearMap())
01221 {
01222 X = new EpetraExt::DistArray<double>(Map, RowSize);
01223
01224 Read(GroupName, "Values", Map.NumMyElements() * RowSize,
01225 Map.NumGlobalElements() * RowSize, H5T_NATIVE_DOUBLE, X->Values());
01226 }
01227 else
01228 {
01229
01230
01231 Epetra_Map LinearMap(GlobalLength, Map.IndexBase(), Comm_);
01232 EpetraExt::DistArray<double> LinearX(LinearMap, RowSize);
01233
01234 Read(GroupName, "Values", LinearMap.NumMyElements() * RowSize,
01235 LinearMap.NumGlobalElements() * RowSize,
01236 H5T_NATIVE_DOUBLE, LinearX.Values());
01237
01238 Epetra_Import Importer(Map, LinearMap);
01239 X = new EpetraExt::DistArray<double>(Map, RowSize);
01240 X->Import(LinearX, Importer, Insert);
01241 }
01242 }
01243
01244
01245 void EpetraExt::HDF5::Read(const std::string& GroupName, DistArray<double>*& X)
01246 {
01247
01248 int GlobalLength, RowSize;
01249 ReadDoubleDistArrayProperties(GroupName, GlobalLength, RowSize);
01250
01251
01252 Epetra_Map LinearMap(GlobalLength, 0, Comm_);
01253 X = new EpetraExt::DistArray<double>(LinearMap, RowSize);
01254
01255 Read(GroupName, "Values", LinearMap.NumMyElements() * RowSize,
01256 LinearMap.NumGlobalElements() * RowSize, H5T_NATIVE_DOUBLE, X->Values());
01257 }
01258
01259
01260 void EpetraExt::HDF5::ReadDoubleDistArrayProperties(const std::string& GroupName,
01261 int& GlobalLength,
01262 int& RowSize)
01263 {
01264 if (!IsContained(GroupName))
01265 throw(Exception(__FILE__, __LINE__,
01266 "requested group " + GroupName + " not found"));
01267
01268 std::string Label;
01269 Read(GroupName, "__type__", Label);
01270
01271 if (Label != "EpetraExt::DistArray<double>")
01272 throw(Exception(__FILE__, __LINE__,
01273 "requested group " + GroupName + " is not an EpetraExt::DistArray<double>",
01274 "__type__ = " + Label));
01275
01276 Read(GroupName, "GlobalLength", GlobalLength);
01277 Read(GroupName, "RowSize", RowSize);
01278 }
01279
01280
01281
01282
01283
01284
01285 void EpetraExt::HDF5::Write(const std::string& GroupName, const std::string& DataSetName,
01286 int what)
01287 {
01288 if (!IsContained(GroupName))
01289 CreateGroup(GroupName);
01290
01291 hid_t filespace_id = H5Screate(H5S_SCALAR);
01292 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01293 hid_t dset_id = H5Dcreate(group_id, DataSetName.c_str(), H5T_NATIVE_INT,
01294 filespace_id, H5P_DEFAULT);
01295
01296 herr_t status = H5Dwrite(dset_id, H5T_NATIVE_INT, H5S_ALL, filespace_id,
01297 H5P_DEFAULT, &what);
01298 CHECK_STATUS(status);
01299
01300
01301 H5Dclose(dset_id);
01302 H5Gclose(group_id);
01303 H5Sclose(filespace_id);
01304 }
01305
01306
01307 void EpetraExt::HDF5::Write(const std::string& GroupName, const std::string& DataSetName,
01308 double what)
01309 {
01310 if (!IsContained(GroupName))
01311 CreateGroup(GroupName);
01312
01313 hid_t filespace_id = H5Screate(H5S_SCALAR);
01314 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01315 hid_t dset_id = H5Dcreate(group_id, DataSetName.c_str(), H5T_NATIVE_DOUBLE,
01316 filespace_id, H5P_DEFAULT);
01317
01318 herr_t status = H5Dwrite(dset_id, H5T_NATIVE_DOUBLE, H5S_ALL,
01319 filespace_id, H5P_DEFAULT, &what);
01320 CHECK_STATUS(status);
01321
01322
01323 H5Dclose(dset_id);
01324 H5Gclose(group_id);
01325 H5Sclose(filespace_id);
01326 }
01327
01328
01329 void EpetraExt::HDF5::Read(const std::string& GroupName, const std::string& DataSetName, int& data)
01330 {
01331 if (!IsContained(GroupName))
01332 throw(Exception(__FILE__, __LINE__,
01333 "requested group " + GroupName + " not found"));
01334
01335
01336 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01337
01338 hid_t filespace_id = H5Screate(H5S_SCALAR);
01339 hid_t dset_id = H5Dopen(group_id, DataSetName.c_str());
01340
01341 herr_t status = H5Dread(dset_id, H5T_NATIVE_INT, H5S_ALL, filespace_id,
01342 H5P_DEFAULT, &data);
01343 CHECK_STATUS(status);
01344
01345 H5Sclose(filespace_id);
01346 H5Dclose(dset_id);
01347 H5Gclose(group_id);
01348 }
01349
01350
01351 void EpetraExt::HDF5::Read(const std::string& GroupName, const std::string& DataSetName, double& data)
01352 {
01353 if (!IsContained(GroupName))
01354 throw(Exception(__FILE__, __LINE__,
01355 "requested group " + GroupName + " not found"));
01356
01357
01358 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01359
01360 hid_t filespace_id = H5Screate(H5S_SCALAR);
01361 hid_t dset_id = H5Dopen(group_id, DataSetName.c_str());
01362
01363 herr_t status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, H5S_ALL, filespace_id,
01364 H5P_DEFAULT, &data);
01365 CHECK_STATUS(status);
01366
01367 H5Sclose(filespace_id);
01368 H5Dclose(dset_id);
01369 H5Gclose(group_id);
01370 }
01371
01372
01373 void EpetraExt::HDF5::Write(const std::string& GroupName,
01374 const std::string& DataSetName,
01375 const std::string& data)
01376 {
01377 if (!IsContained(GroupName))
01378 CreateGroup(GroupName);
01379
01380 hsize_t len = 1;
01381
01382 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01383
01384 hid_t dataspace_id = H5Screate_simple(1, &len, NULL);
01385
01386 hid_t atype = H5Tcopy(H5T_C_S1);
01387 H5Tset_size(atype, data.size() + 1);
01388
01389 hid_t dataset_id = H5Dcreate(group_id, DataSetName.c_str(), atype,
01390 dataspace_id, H5P_DEFAULT);
01391 CHECK_STATUS(H5Dwrite(dataset_id, atype, H5S_ALL, H5S_ALL,
01392 H5P_DEFAULT, data.c_str()));
01393
01394 CHECK_STATUS(H5Tclose(atype));
01395
01396 CHECK_STATUS(H5Dclose(dataset_id));
01397
01398 CHECK_STATUS(H5Sclose(dataspace_id));
01399
01400 CHECK_STATUS(H5Gclose(group_id));
01401 }
01402
01403
01404 void EpetraExt::HDF5::Read(const std::string& GroupName,
01405 const std::string& DataSetName,
01406 std::string& data)
01407 {
01408 if (!IsContained(GroupName))
01409 throw(Exception(__FILE__, __LINE__,
01410 "requested group " + GroupName + " not found"));
01411
01412 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01413
01414 hid_t dataset_id = H5Dopen(group_id, DataSetName.c_str());
01415
01416 hid_t datatype_id = H5Dget_type(dataset_id);
01417 size_t typesize_id = H5Tget_size(datatype_id);
01418 H5T_class_t typeclass_id = H5Tget_class(datatype_id);
01419
01420 if(typeclass_id != H5T_STRING)
01421 throw(Exception(__FILE__, __LINE__,
01422 "requested group " + GroupName + " is not a std::string"));
01423 char data2[160];
01424 CHECK_STATUS(H5Dread(dataset_id, datatype_id, H5S_ALL, H5S_ALL,
01425 H5P_DEFAULT, data2));
01426 data = data2;
01427
01428 H5Dclose(dataset_id);
01429 H5Gclose(group_id);
01430 }
01431
01432
01433
01434
01435
01436
01437 void EpetraExt::HDF5::Write(const std::string& GroupName, const std::string& DataSetName,
01438 const int type, const int Length,
01439 void* data)
01440 {
01441 if (!IsContained(GroupName))
01442 CreateGroup(GroupName);
01443
01444 hsize_t dimsf = Length;
01445
01446 hid_t filespace_id = H5Screate_simple(1, &dimsf, NULL);
01447
01448 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01449
01450 hid_t dset_id = H5Dcreate(group_id, DataSetName.c_str(), type,
01451 filespace_id, H5P_DEFAULT);
01452
01453 herr_t status = H5Dwrite(dset_id, type, H5S_ALL, H5S_ALL,
01454 H5P_DEFAULT, data);
01455 CHECK_STATUS(status);
01456
01457
01458 H5Dclose(dset_id);
01459 H5Gclose(group_id);
01460 H5Sclose(filespace_id);
01461 }
01462
01463
01464 void EpetraExt::HDF5::Read(const std::string& GroupName, const std::string& DataSetName,
01465 const int type, const int Length,
01466 void* data)
01467 {
01468 if (!IsContained(GroupName))
01469 throw(Exception(__FILE__, __LINE__,
01470 "requested group " + GroupName + " not found"));
01471
01472 hsize_t dimsf = Length;
01473
01474 hid_t filespace_id = H5Screate_simple(1, &dimsf, NULL);
01475
01476 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01477
01478 hid_t dset_id = H5Dopen(group_id, DataSetName.c_str());
01479
01480 herr_t status = H5Dread(dset_id, type, H5S_ALL, filespace_id,
01481 H5P_DEFAULT, data);
01482 CHECK_STATUS(status);
01483
01484
01485 H5Dclose(dset_id);
01486 H5Gclose(group_id);
01487 H5Sclose(filespace_id);
01488 }
01489
01490
01491
01492
01493
01494
01495 void EpetraExt::HDF5::Write(const std::string& GroupName, const std::string& DataSetName,
01496 int MySize, int GlobalSize, int type, const void* data)
01497 {
01498 int Offset;
01499 Comm_.ScanSum(&MySize, &Offset, 1);
01500 Offset -= MySize;
01501
01502 hsize_t MySize_t = MySize;
01503 hsize_t GlobalSize_t = GlobalSize;
01504 hsize_t Offset_t = Offset;
01505
01506 hid_t filespace_id = H5Screate_simple(1, &GlobalSize_t, NULL);
01507
01508
01509 if (!IsContained(GroupName))
01510 CreateGroup(GroupName);
01511
01512 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01513 hid_t dset_id = H5Dcreate(group_id, DataSetName.c_str(), type, filespace_id,
01514 H5P_DEFAULT);
01515
01516 H5Sclose(filespace_id);
01517
01518
01519
01520
01521 hid_t memspace_id = H5Screate_simple(1, &MySize_t, NULL);
01522
01523
01524 filespace_id = H5Dget_space(dset_id);
01525 H5Sselect_hyperslab(filespace_id, H5S_SELECT_SET, &Offset_t, NULL, &MySize_t, NULL);
01526
01527 plist_id_ = H5Pcreate(H5P_DATASET_XFER);
01528 #ifdef HAVE_MPI
01529 H5Pset_dxpl_mpio(plist_id_, H5FD_MPIO_COLLECTIVE);
01530 #endif
01531
01532 status = H5Dwrite(dset_id, type, memspace_id, filespace_id,
01533 plist_id_, data);
01534 CHECK_STATUS(status);
01535
01536
01537 H5Dclose(dset_id);
01538 H5Gclose(group_id);
01539 H5Sclose(filespace_id);
01540 H5Sclose(memspace_id);
01541 H5Pclose(plist_id_);
01542 }
01543
01544
01545 void EpetraExt::HDF5::Read(const std::string& GroupName, const std::string& DataSetName,
01546 int MySize, int GlobalSize,
01547 const int type, void* data)
01548 {
01549 if (!IsOpen())
01550 throw(Exception(__FILE__, __LINE__, "no file open yet"));
01551
01552 hsize_t MySize_t = MySize;
01553
01554
01555 int itmp;
01556 Comm_.ScanSum(&MySize, &itmp, 1);
01557 hsize_t Offset_t = itmp - MySize;
01558
01559 hid_t group_id = H5Gopen(file_id_, GroupName.c_str());
01560 hid_t dataset_id = H5Dopen(group_id, DataSetName.c_str());
01561
01562
01563
01564 hid_t filespace_id = H5Dget_space(dataset_id);
01565 H5Sselect_hyperslab(filespace_id, H5S_SELECT_SET, &Offset_t, NULL,
01566 &MySize_t, NULL);
01567
01568 hid_t mem_dataspace = H5Screate_simple (1, &MySize_t, NULL);
01569
01570 herr_t status = H5Dread(dataset_id, type, mem_dataspace, filespace_id,
01571 H5P_DEFAULT, data);
01572 CHECK_STATUS(status);
01573
01574 H5Sclose(mem_dataspace);
01575 H5Gclose(group_id);
01576
01577 H5Dclose(dataset_id);
01578
01579 }
01580
01581
01582 #endif // HAVE_EPETRAEXT_HDF5
01583