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