00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032 #include <set>
00033 #include <cstring>
00034
00035 #include "Hdf5DataWriter.hpp"
00036
00037 #include "Exception.hpp"
00038 #include "OutputFileHandler.hpp"
00039 #include "PetscTools.hpp"
00040 #include "Version.hpp"
00041
00042 Hdf5DataWriter::Hdf5DataWriter(DistributedVectorFactory& rVectorFactory,
00043 const std::string& rDirectory,
00044 const std::string& rBaseName,
00045 bool cleanDirectory,
00046 bool extendData)
00047 : mrVectorFactory(rVectorFactory),
00048 mDirectory(rDirectory),
00049 mBaseName(rBaseName),
00050 mCleanDirectory(cleanDirectory),
00051 mIsInDefineMode(true),
00052 mIsFixedDimensionSet(false),
00053 mIsUnlimitedDimensionSet(false),
00054 mFileFixedDimensionSize(0u),
00055 mDataFixedDimensionSize(0u),
00056 mLo(mrVectorFactory.GetLow()),
00057 mHi(mrVectorFactory.GetHigh()),
00058 mNumberOwned(0u),
00059 mOffset(0u),
00060 mIsDataComplete(true),
00061 mNeedExtend(false),
00062 mUseMatrixForIncompleteData(false),
00063 mCurrentTimeStep(0),
00064 mSinglePermutation(NULL),
00065 mDoublePermutation(NULL),
00066 mSingleIncompleteOutputMatrix(NULL),
00067 mDoubleIncompleteOutputMatrix(NULL)
00068 {
00069 if (extendData && cleanDirectory)
00070 {
00071 EXCEPTION("You are asking to delete a file and then extend it, change arguments to constructor.");
00072 }
00073
00074 if (extendData)
00075 {
00076
00077 OutputFileHandler output_file_handler(mDirectory, false);
00078 std::string results_dir = output_file_handler.GetOutputDirectoryFullPath();
00079 std::string file_name = results_dir + mBaseName + ".h5";
00080
00081
00082 hid_t property_list_id = H5Pcreate(H5P_FILE_ACCESS);
00083 H5Pset_fapl_mpio(property_list_id, PETSC_COMM_WORLD, MPI_INFO_NULL);
00084
00085
00086 mFileId = H5Fopen(file_name.c_str(), H5F_ACC_RDWR, property_list_id);
00087 H5Pclose(property_list_id);
00088
00089 if (mFileId < 0)
00090 {
00091 mDatasetId = 0;
00092 EXCEPTION("Hdf5DataWriter could not open " + file_name);
00093 }
00094
00095
00096 mDatasetId = H5Dopen(mFileId, "Data");
00097 hid_t variables_dataspace = H5Dget_space(mDatasetId);
00098
00099 hsize_t dataset_max_sizes[DATASET_DIMS];
00100 H5Sget_simple_extent_dims(variables_dataspace, mDatasetDims, dataset_max_sizes);
00101 H5Sclose(variables_dataspace);
00102
00103 if (dataset_max_sizes[0] != H5S_UNLIMITED)
00104 {
00105 H5Dclose(mDatasetId);
00106 H5Fclose(mFileId);
00107 EXCEPTION("Tried to open a datafile for extending which doesn't have an unlimited dimension.");
00108 }
00109 mIsUnlimitedDimensionSet = true;
00110
00111 for (unsigned i=1; i<DATASET_DIMS; i++)
00112 {
00113 assert(mDatasetDims[i] == dataset_max_sizes[i]);
00114 }
00115 mFileFixedDimensionSize = mDatasetDims[1];
00116 mIsFixedDimensionSet = true;
00117 mVariables.reserve(mDatasetDims[2]);
00118
00119
00120 hid_t attribute_id = H5Aopen_name(mDatasetId, "Variable Details");
00121 hid_t attribute_type = H5Aget_type(attribute_id);
00122 hid_t attribute_space = H5Aget_space(attribute_id);
00123 hsize_t attr_dataspace_dim;
00124 H5Sget_simple_extent_dims(attribute_space, &attr_dataspace_dim, NULL);
00125 unsigned num_columns = H5Sget_simple_extent_npoints(attribute_space);
00126 assert(num_columns == mDatasetDims[2]);
00127
00128 char* string_array = (char *)malloc(sizeof(char)*MAX_STRING_SIZE*(int)num_columns);
00129 H5Aread(attribute_id, attribute_type, string_array);
00130
00131 for (unsigned index=0; index<num_columns; index++)
00132 {
00133
00134 std::string column_name_unit(&string_array[MAX_STRING_SIZE*index]);
00135
00136 size_t name_length = column_name_unit.find('(');
00137 size_t unit_length = column_name_unit.find(')') - name_length - 1;
00138
00139 DataWriterVariable var;
00140 var.mVariableName = column_name_unit.substr(0, name_length);
00141 var.mVariableUnits = column_name_unit.substr(name_length+1, unit_length);
00142 mVariables.push_back(var);
00143 }
00144
00145 free(string_array);
00146 H5Tclose(attribute_type);
00147 H5Sclose(attribute_space);
00148 H5Aclose(attribute_id);
00149
00150
00151 mUnlimitedDimensionName = "Time";
00152 mTimeDatasetId = H5Dopen(mFileId, mUnlimitedDimensionName.c_str());
00153 mUnlimitedDimensionUnit = "ms";
00154
00155 hid_t timestep_dataspace = H5Dget_space(mTimeDatasetId);
00156 hsize_t num_timesteps;
00157 H5Sget_simple_extent_dims(timestep_dataspace, &num_timesteps, NULL);
00158 H5Sclose(timestep_dataspace);
00159 mCurrentTimeStep = (long)num_timesteps - 1;
00160
00161
00162 attribute_id = H5Aopen_name(mDatasetId, "IsDataComplete");
00163 if (attribute_id < 0)
00164 {
00165 #define COVERAGE_IGNORE
00166
00167 EXCEPTION("Extending old-format files isn't supported.");
00168 #undef COVERAGE_IGNORE
00169 }
00170 else
00171 {
00172 attribute_type = H5Aget_type(attribute_id);
00173 attribute_space = H5Aget_space(attribute_id);
00174 unsigned is_data_complete;
00175 H5Aread(attribute_id, H5T_NATIVE_UINT, &is_data_complete);
00176 mIsDataComplete = (is_data_complete == 1) ? true : false;
00177 H5Tclose(attribute_type);
00178 H5Sclose(attribute_space);
00179 H5Aclose(attribute_id);
00180 }
00181 if (mIsDataComplete)
00182 {
00183 mNumberOwned = mrVectorFactory.GetLocalOwnership();
00184 mOffset = mLo;
00185 mDataFixedDimensionSize = mFileFixedDimensionSize;
00186 }
00187 else
00188 {
00189
00190 attribute_id = H5Aopen_name(mDatasetId, "NodeMap");
00191 attribute_type = H5Aget_type(attribute_id);
00192 attribute_space = H5Aget_space(attribute_id);
00193
00194 unsigned num_node_indices = H5Sget_simple_extent_npoints(attribute_space);
00195
00196 mIncompleteNodeIndices.clear();
00197 mIncompleteNodeIndices.resize(num_node_indices);
00198 H5Aread(attribute_id, H5T_NATIVE_UINT, &mIncompleteNodeIndices[0]);
00199
00200 H5Tclose(attribute_type);
00201 H5Sclose(attribute_space);
00202 H5Aclose(attribute_id);
00203
00204 mNumberOwned = mrVectorFactory.GetLocalOwnership();
00205 ComputeIncompleteOffset();
00209 mDataFixedDimensionSize = UINT_MAX;
00210 H5Dclose(mDatasetId);
00211 H5Dclose(mTimeDatasetId);
00212 H5Fclose(mFileId);
00213 EXCEPTION("Unable to extend an incomplete data file at present.");
00214 }
00215
00216
00217 mIsInDefineMode = false;
00218 AdvanceAlongUnlimitedDimension();
00219 }
00220 }
00221
00222 Hdf5DataWriter::~Hdf5DataWriter()
00223 {
00224 Close();
00225
00226 if (mSinglePermutation)
00227 {
00228 MatDestroy(mSinglePermutation);
00229 }
00230 if (mDoublePermutation)
00231 {
00232 MatDestroy(mDoublePermutation);
00233 }
00234 if (mSingleIncompleteOutputMatrix)
00235 {
00236 MatDestroy(mSingleIncompleteOutputMatrix);
00237 }
00238 if (mDoubleIncompleteOutputMatrix)
00239 {
00240 MatDestroy(mDoubleIncompleteOutputMatrix);
00241 }
00242 }
00243
00244 void Hdf5DataWriter::DefineFixedDimension(long dimensionSize)
00245 {
00246 if (!mIsInDefineMode)
00247 {
00248 EXCEPTION("Cannot define variables when not in Define mode");
00249 }
00250 if (dimensionSize < 1)
00251 {
00252 EXCEPTION("Fixed dimension must be at least 1 long");
00253 }
00254 if (mIsFixedDimensionSet)
00255 {
00256 EXCEPTION("Fixed dimension already set");
00257 }
00258
00259
00260 mLo = mrVectorFactory.GetLow();
00261 mHi = mrVectorFactory.GetHigh();
00262 mNumberOwned = mrVectorFactory.GetLocalOwnership();
00263 mOffset = mLo;
00264 mFileFixedDimensionSize = dimensionSize;
00265 mDataFixedDimensionSize = dimensionSize;
00266 mIsFixedDimensionSet = true;
00267 }
00268
00269 void Hdf5DataWriter::DefineFixedDimension(const std::vector<unsigned>& rNodesToOuput, long vecSize)
00270 {
00271 unsigned vector_size = rNodesToOuput.size();
00272
00273 for (unsigned index=0; index < vector_size-1; index++)
00274 {
00275 if (rNodesToOuput[index] >= rNodesToOuput[index+1])
00276 {
00277 EXCEPTION("Input should be monotonic increasing");
00278 }
00279 }
00280
00281 if ((int) rNodesToOuput.back() >= vecSize)
00282 {
00283 EXCEPTION("Vector size doesn't match nodes to output");
00284 }
00285
00286 DefineFixedDimension(vecSize);
00287
00288 mFileFixedDimensionSize = vector_size;
00289 mIsDataComplete = false;
00290 mIncompleteNodeIndices = rNodesToOuput;
00291 ComputeIncompleteOffset();
00292 }
00293
00294 void Hdf5DataWriter::DefineFixedDimensionUsingMatrix(const std::vector<unsigned>& rNodesToOuput, long vecSize)
00295 {
00296 unsigned vector_size = rNodesToOuput.size();
00297
00298 for (unsigned index=0; index < vector_size-1; index++)
00299 {
00300 if (rNodesToOuput[index] >= rNodesToOuput[index+1])
00301 {
00302 EXCEPTION("Input should be monotonic increasing");
00303 }
00304 }
00305
00306 if ((int) rNodesToOuput.back() >= vecSize)
00307 {
00308 EXCEPTION("Vector size doesn't match nodes to output");
00309 }
00310
00311 DefineFixedDimension(vecSize);
00312
00313 mFileFixedDimensionSize = vector_size;
00314 mIsDataComplete = false;
00315 mIncompleteNodeIndices = rNodesToOuput;
00316 mUseMatrixForIncompleteData = true;
00317 ComputeIncompleteOffset();
00318
00319
00320 assert(mSingleIncompleteOutputMatrix == NULL);
00321 assert(mDoubleIncompleteOutputMatrix == NULL);
00322 PetscTools::SetupMat(mSingleIncompleteOutputMatrix, mFileFixedDimensionSize, mDataFixedDimensionSize, 2, mNumberOwned, mHi - mLo);
00323 PetscTools::SetupMat(mDoubleIncompleteOutputMatrix, 2*mFileFixedDimensionSize, 2*mDataFixedDimensionSize, 4, 2*mNumberOwned, 2*(mHi - mLo));
00324
00325
00326 for (unsigned row_index = mOffset; row_index < mOffset + mNumberOwned; row_index++)
00327 {
00328
00329 MatSetValue(mSingleIncompleteOutputMatrix, row_index, row_index, 0.0, INSERT_VALUES);
00330
00331 MatSetValue(mSingleIncompleteOutputMatrix, row_index, rNodesToOuput[row_index], 1.0, INSERT_VALUES);
00332
00333 unsigned bi_index=2*row_index;
00334 unsigned perm_index=2*rNodesToOuput[row_index];
00335
00336 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index, bi_index, 0.0, INSERT_VALUES);
00337 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index+1, bi_index+1, 0.0, INSERT_VALUES);
00338
00339 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index, perm_index, 1.0, INSERT_VALUES);
00340 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index+1, perm_index+1, 1.0, INSERT_VALUES);
00341 }
00342 MatAssemblyBegin(mSingleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY);
00343 MatAssemblyBegin(mDoubleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY);
00344 MatAssemblyEnd(mSingleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY);
00345 MatAssemblyEnd(mDoubleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY);
00346
00347
00348
00349 }
00350
00351 void Hdf5DataWriter::ComputeIncompleteOffset()
00352 {
00353 mOffset = 0;
00354 mNumberOwned = 0;
00355 for (unsigned i=0; i<mIncompleteNodeIndices.size(); i++)
00356 {
00357 if (mIncompleteNodeIndices[i] < mLo)
00358 {
00359 mOffset++;
00360 }
00361 else if (mIncompleteNodeIndices[i] < mHi)
00362 {
00363 mNumberOwned++;
00364 }
00365 }
00366 }
00367
00368 int Hdf5DataWriter::DefineVariable(const std::string& rVariableName,
00369 const std::string& rVariableUnits)
00370 {
00371 if (!mIsInDefineMode)
00372 {
00373 EXCEPTION("Cannot define variables when not in Define mode");
00374 }
00375
00376 CheckVariableName(rVariableName);
00377 CheckUnitsName(rVariableUnits);
00378
00379
00380 for (unsigned index=0; index<mVariables.size(); index++)
00381 {
00382 if (mVariables[index].mVariableName == rVariableName)
00383 {
00384 EXCEPTION("Variable name already exists");
00385 }
00386 }
00387
00388 DataWriterVariable new_variable;
00389 new_variable.mVariableName = rVariableName;
00390 new_variable.mVariableUnits = rVariableUnits;
00391 int variable_id;
00392
00393
00394 mVariables.push_back(new_variable);
00395
00396
00397
00398 variable_id = mVariables.size() - 1;
00399
00400 return variable_id;
00401 }
00402
00403 void Hdf5DataWriter::CheckVariableName(const std::string& rName)
00404 {
00405 if (rName.length() == 0)
00406 {
00407 EXCEPTION("Variable name not allowed: may not be blank.");
00408 }
00409 CheckUnitsName(rName);
00410 }
00411
00412 void Hdf5DataWriter::CheckUnitsName(const std::string& rName)
00413 {
00414 for (unsigned i=0; i<rName.length(); i++)
00415 {
00416 if (!isalnum(rName[i]) && !(rName[i]=='_'))
00417 {
00418 std::string error = "Variable name/units '" + rName + "' not allowed: may only contain alphanumeric characters or '_'.";
00419 EXCEPTION(error);
00420 }
00421 }
00422 }
00423
00424 void Hdf5DataWriter::EndDefineMode()
00425 {
00426
00427 if (mVariables.size() < 1)
00428 {
00429 EXCEPTION("Cannot end define mode. No variables have been defined.");
00430 }
00431
00432
00433 if (mIsFixedDimensionSet == false)
00434 {
00435 EXCEPTION("Cannot end define mode. One fixed dimension should be defined.");
00436 }
00437
00438 OutputFileHandler output_file_handler(mDirectory, mCleanDirectory);
00439 std::string results_dir = output_file_handler.GetOutputDirectoryFullPath();
00440 std::string file_name = results_dir + mBaseName + ".h5";
00441
00442
00443 hid_t property_list_id = H5Pcreate(H5P_FILE_ACCESS);
00444 H5Pset_fapl_mpio(property_list_id, PETSC_COMM_WORLD, MPI_INFO_NULL);
00445
00446
00447 mFileId = H5Fcreate(file_name.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, property_list_id);
00448 H5Pclose(property_list_id);
00449 if (mFileId < 0)
00450 {
00451 EXCEPTION("Hdf5DataWriter could not create " + file_name);
00452 }
00453 mIsInDefineMode = false;
00454
00455
00456
00457
00458
00459 mDatasetDims[0] = 1;
00460 mDatasetDims[1] = mFileFixedDimensionSize;
00461 mDatasetDims[2] = mVariables.size();
00462
00463 hsize_t* max_dims = NULL;
00464 hsize_t dataset_max_dims[DATASET_DIMS];
00465
00466 hid_t cparms = H5P_DEFAULT;
00467
00468 if (mIsUnlimitedDimensionSet)
00469 {
00470 dataset_max_dims[0] = H5S_UNLIMITED;
00471 dataset_max_dims[1] = mDatasetDims[1];
00472 dataset_max_dims[2] = mDatasetDims[2];
00473 max_dims = dataset_max_dims;
00474
00475
00476 hsize_t chunk_dims[DATASET_DIMS] ={1, mDatasetDims[1], mDatasetDims[2]};
00477 cparms = H5Pcreate (H5P_DATASET_CREATE);
00478 H5Pset_chunk( cparms, DATASET_DIMS, chunk_dims);
00479 }
00480
00481 hid_t filespace = H5Screate_simple(DATASET_DIMS, mDatasetDims, max_dims);
00482
00483
00484 mDatasetId = H5Dcreate(mFileId, "Data", H5T_NATIVE_DOUBLE, filespace, cparms);
00485 H5Sclose(filespace);
00486
00487
00488 const unsigned MAX_STRING_SIZE = 100;
00489 hsize_t columns[1] = {mVariables.size()};
00490 hid_t colspace = H5Screate_simple(1, columns, NULL);
00491
00492
00493 char* col_data = (char*) malloc(mVariables.size() * sizeof(char) * MAX_STRING_SIZE);
00494
00495 char* col_data_offset = col_data;
00496 for (unsigned var=0; var<mVariables.size(); var++)
00497 {
00498 std::string full_name = mVariables[var].mVariableName + "(" + mVariables[var].mVariableUnits + ")";
00499 strcpy (col_data_offset, full_name.c_str());
00500 col_data_offset += sizeof(char) * MAX_STRING_SIZE;
00501 }
00502
00503
00504 hid_t string_type = H5Tcopy(H5T_C_S1);
00505 H5Tset_size(string_type, MAX_STRING_SIZE );
00506 hid_t attr = H5Acreate(mDatasetId, "Variable Details", string_type, colspace, H5P_DEFAULT);
00507
00508
00509 H5Awrite(attr, string_type, col_data);
00510
00511
00512 free(col_data);
00513 H5Sclose(colspace);
00514 H5Aclose(attr);
00515
00516
00517 columns[0] = 1;
00518 colspace = H5Screate_simple(1, columns, NULL);
00519 attr = H5Acreate(mDatasetId, "IsDataComplete", H5T_NATIVE_UINT, colspace, H5P_DEFAULT);
00520
00521
00522 unsigned is_data_complete = mIsDataComplete ? 1 : 0;
00523 H5Awrite(attr, H5T_NATIVE_UINT, &is_data_complete);
00524
00525 H5Sclose(colspace);
00526 H5Aclose(attr);
00527
00528 if (!mIsDataComplete)
00529 {
00530
00531
00532 columns[0] = mFileFixedDimensionSize;
00533 colspace = H5Screate_simple(1, columns, NULL);
00534 attr = H5Acreate(mDatasetId, "NodeMap", H5T_NATIVE_UINT, colspace, H5P_DEFAULT);
00535
00536
00537 H5Awrite(attr, H5T_NATIVE_UINT, &mIncompleteNodeIndices[0]);
00538
00539 H5Sclose(colspace);
00540 H5Aclose(attr);
00541 }
00542
00543
00544
00545
00546 if (mIsUnlimitedDimensionSet)
00547 {
00548 hsize_t time_dataset_dims[1] = {1};
00549 hsize_t time_dataset_max_dims[1] = {H5S_UNLIMITED};
00550
00551
00552 hsize_t time_chunk_dims[1] ={1};
00553 hid_t time_cparms = H5Pcreate (H5P_DATASET_CREATE);
00554 H5Pset_chunk( time_cparms, 1, time_chunk_dims);
00555
00556 hid_t time_filespace = H5Screate_simple(1, time_dataset_dims, time_dataset_max_dims);
00557
00558
00559 mTimeDatasetId = H5Dcreate(mFileId, mUnlimitedDimensionName.c_str(), H5T_NATIVE_DOUBLE, time_filespace, time_cparms);
00560
00561
00562 hsize_t one = 1;
00563 hid_t one_column_space = H5Screate_simple(1, &one, NULL);
00564
00565
00566 hid_t unit_attr = H5Acreate(mTimeDatasetId, "Unit", string_type, one_column_space, H5P_DEFAULT);
00567
00568
00569 char unit_string[MAX_STRING_SIZE];
00570 strcpy(unit_string, mUnlimitedDimensionUnit.c_str());
00571 H5Awrite(unit_attr, string_type, unit_string);
00572
00573
00574 H5Sclose(one_column_space);
00575 H5Aclose(unit_attr);
00576 H5Sclose(time_filespace);
00577 }
00578
00579
00580
00581
00582
00583 const unsigned MAX_PROVENANCE_STRING_SIZE = 365;
00584 hid_t long_string_type = H5Tcopy(H5T_C_S1);
00585 H5Tset_size(long_string_type, MAX_PROVENANCE_STRING_SIZE );
00586 hsize_t prov_columns[1] = {1};
00587 hid_t provenance_space = H5Screate_simple(1, prov_columns, NULL);
00588 char* provenance_data = (char*) malloc(sizeof(char) * MAX_PROVENANCE_STRING_SIZE);
00589 assert( ChasteBuildInfo::GetProvenanceString().length() < MAX_PROVENANCE_STRING_SIZE);
00590
00591 strcpy(provenance_data, ChasteBuildInfo::GetProvenanceString().c_str());
00592 hid_t prov_attr = H5Acreate(mDatasetId, "Chaste Provenance", long_string_type, provenance_space, H5P_DEFAULT);
00593
00594
00595 H5Awrite(prov_attr, long_string_type, provenance_data);
00596
00597
00598 free(provenance_data);
00599 H5Sclose(provenance_space);
00600 H5Aclose(prov_attr);
00601
00602 }
00603
00604 void Hdf5DataWriter::PutVector(int variableID, Vec petscVector)
00605 {
00606 if (mIsInDefineMode)
00607 {
00608 EXCEPTION("Cannot write data while in define mode.");
00609 }
00610
00611 int vector_size;
00612 VecGetSize(petscVector, &vector_size);
00613
00614 if ((unsigned) vector_size != mDataFixedDimensionSize)
00615 {
00616 EXCEPTION("Vector size doesn't match fixed dimension");
00617 }
00618
00619
00620 PossiblyExtend();
00621
00622 Vec output_petsc_vector;
00623
00624 if (mSinglePermutation == NULL)
00625 {
00626
00627 output_petsc_vector = petscVector;
00628 }
00629 else
00630 {
00631 assert(mIsDataComplete);
00632
00633 VecDuplicate(petscVector, &output_petsc_vector);
00634
00635 MatMult(mSinglePermutation, petscVector, output_petsc_vector);
00636 }
00637
00638 hid_t memspace=0;
00639 if (mNumberOwned !=0)
00640 {
00641 hsize_t v_size[1] = {mNumberOwned};
00642 memspace = H5Screate_simple(1, v_size, NULL);
00643 }
00644
00645 hsize_t count[DATASET_DIMS] = {1, mNumberOwned, 1};
00646 hsize_t offset_dims[DATASET_DIMS] = {mCurrentTimeStep, mOffset, variableID};
00647 hid_t file_dataspace = H5Dget_space(mDatasetId);
00648
00649
00650 hid_t property_list_id = H5Pcreate(H5P_DATASET_XFER);
00651 H5Pset_dxpl_mpio(property_list_id, H5FD_MPIO_COLLECTIVE);
00652
00653 H5Sselect_hyperslab(file_dataspace, H5S_SELECT_SET, offset_dims, NULL, count, NULL);
00654
00655 double* p_petsc_vector;
00656 VecGetArray(output_petsc_vector, &p_petsc_vector);
00657
00658 if (mIsDataComplete)
00659 {
00660 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, file_dataspace, property_list_id, p_petsc_vector);
00661 }
00662 else
00663 {
00664 if (mUseMatrixForIncompleteData)
00665 {
00666
00667 output_petsc_vector = PetscTools::CreateVec(mFileFixedDimensionSize, mNumberOwned);
00668
00669
00670 MatMult(mSingleIncompleteOutputMatrix, petscVector, output_petsc_vector);
00671
00672 double* p_petsc_vector_incomplete;
00673 VecGetArray(output_petsc_vector, &p_petsc_vector_incomplete);
00674 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, file_dataspace, property_list_id, p_petsc_vector_incomplete);
00675 }
00676 else
00677 {
00678
00679 double local_data[mNumberOwned];
00680 for (unsigned i=0; i<mNumberOwned; i++)
00681 {
00682 local_data[i] = p_petsc_vector[ mIncompleteNodeIndices[mOffset+i]-mLo ];
00683
00684 }
00685 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, file_dataspace, property_list_id, local_data);
00686 }
00687 }
00688
00689 VecRestoreArray(output_petsc_vector, &p_petsc_vector);
00690
00691 H5Pclose(property_list_id);
00692 H5Sclose(file_dataspace);
00693 if (mNumberOwned !=0)
00694 {
00695 H5Sclose(memspace);
00696 }
00697
00698 if (petscVector != output_petsc_vector)
00699 {
00700
00701 VecDestroy(output_petsc_vector);
00702 }
00703 }
00704
00705 void Hdf5DataWriter::PutStripedVector(std::vector<int> variableIDs, Vec petscVector)
00706 {
00707 if (mIsInDefineMode)
00708 {
00709 EXCEPTION("Cannot write data while in define mode.");
00710 }
00711
00712 if (variableIDs.size() <= 1)
00713 {
00714 EXCEPTION("The PutStripedVector method requires at least two variables ID. If only one is needed, use PutVector method instead");
00715 }
00716
00717 const int NUM_STRIPES=variableIDs.size();
00718
00719 int firstVariableID=variableIDs[0];
00720
00721
00722 for (unsigned i = 1; i < variableIDs.size(); i++)
00723 {
00724 if (variableIDs[i]-variableIDs[i-1] != 1)
00725 {
00726 EXCEPTION("Columns should be consecutive. Try reordering them.");
00727 }
00728 }
00729
00730 int vector_size;
00731 VecGetSize(petscVector, &vector_size);
00732
00733 if ((unsigned) vector_size != NUM_STRIPES*mDataFixedDimensionSize)
00734 {
00735 EXCEPTION("Vector size doesn't match fixed dimension");
00736 }
00737
00738
00739 PossiblyExtend();
00740
00741 Vec output_petsc_vector;
00742
00743 if (mDoublePermutation == NULL)
00744 {
00745
00746 output_petsc_vector = petscVector;
00747 }
00748 else
00749 {
00750 assert(mIsDataComplete);
00751
00752 VecDuplicate(petscVector, &output_petsc_vector);
00753
00754 MatMult(mDoublePermutation, petscVector, output_petsc_vector);
00755 }
00756
00757 hid_t memspace=0;
00758 if (mNumberOwned !=0)
00759 {
00760 hsize_t v_size[1] = {mNumberOwned*NUM_STRIPES};
00761 memspace = H5Screate_simple(1, v_size, NULL);
00762 }
00763
00764
00765 hsize_t start[DATASET_DIMS] = {mCurrentTimeStep, mOffset, firstVariableID};
00766 hsize_t stride[DATASET_DIMS] = {1, 1, 1};
00767 hsize_t block_size[DATASET_DIMS] = {1, mNumberOwned, 1};
00768 hsize_t number_blocks[DATASET_DIMS] = {1, 1, NUM_STRIPES};
00769
00770 hid_t hyperslab_space = H5Dget_space(mDatasetId);
00771 H5Sselect_hyperslab(hyperslab_space, H5S_SELECT_SET, start, stride, number_blocks, block_size);
00772
00773
00774 hid_t property_list_id = H5Pcreate(H5P_DATASET_XFER);
00775 H5Pset_dxpl_mpio(property_list_id, H5FD_MPIO_COLLECTIVE);
00776
00777 double* p_petsc_vector;
00778 VecGetArray(output_petsc_vector, &p_petsc_vector);
00779
00780 if (mIsDataComplete)
00781 {
00782 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, property_list_id, p_petsc_vector);
00783 }
00784 else
00785 {
00786 if(variableIDs.size() < 3)
00787 {
00788 if (mUseMatrixForIncompleteData)
00789 {
00790
00791 output_petsc_vector = PetscTools::CreateVec(2*mFileFixedDimensionSize, 2*mNumberOwned);
00792
00793
00794 MatMult(mDoubleIncompleteOutputMatrix, petscVector, output_petsc_vector);
00795
00796 double* p_petsc_vector_incomplete;
00797 VecGetArray(output_petsc_vector, &p_petsc_vector_incomplete);
00798
00799 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, property_list_id, p_petsc_vector_incomplete);
00800 }
00801 else
00802 {
00803
00804 double local_data[mNumberOwned*NUM_STRIPES];
00805 for (unsigned i=0; i<mNumberOwned; i++)
00806 {
00808 unsigned local_node_number=mIncompleteNodeIndices[mOffset+i] - mLo;
00809 local_data[NUM_STRIPES*i] = p_petsc_vector[ local_node_number*NUM_STRIPES ];
00810 local_data[NUM_STRIPES*i+1] = p_petsc_vector[ local_node_number*NUM_STRIPES + 1];
00811 }
00812
00813 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, property_list_id, local_data);
00814 }
00815 }
00816 else
00817 {
00818 EXCEPTION("The PutStripedVector functionality for incomplete data is supported for only 2 stripes");
00819 }
00820 }
00821
00822 VecRestoreArray(output_petsc_vector, &p_petsc_vector);
00823
00824 H5Sclose(hyperslab_space);
00825 if (mNumberOwned != 0)
00826 {
00827 H5Sclose(memspace);
00828 }
00829 H5Pclose(property_list_id);
00830
00831 if (petscVector != output_petsc_vector)
00832 {
00833
00834 VecDestroy(output_petsc_vector);
00835 }
00836 }
00837
00838 void Hdf5DataWriter::PutUnlimitedVariable(double value)
00839 {
00840 if (mIsInDefineMode)
00841 {
00842 EXCEPTION("Cannot write data while in define mode.");
00843 }
00844
00845 if (!mIsUnlimitedDimensionSet)
00846 {
00847 EXCEPTION("PutUnlimitedVariable() called but no unlimited dimension has been set");
00848 }
00849
00850
00851 PossiblyExtend();
00852
00853
00854 if (!PetscTools::AmMaster())
00855 {
00856 return;
00857 }
00858
00859 hsize_t size[1] = {1};
00860 hid_t memspace = H5Screate_simple(1, size, NULL);
00861
00862
00863 hsize_t count[1] = {1};
00864 hsize_t offset[1] = {mCurrentTimeStep};
00865 hid_t hyperslab_space = H5Dget_space(mTimeDatasetId);
00866 H5Sselect_hyperslab(hyperslab_space, H5S_SELECT_SET, offset, NULL, count, NULL);
00867
00868 H5Dwrite(mTimeDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, H5P_DEFAULT, &value);
00869
00870 H5Sclose(hyperslab_space);
00871 H5Sclose(memspace);
00872 }
00873
00874 void Hdf5DataWriter::Close()
00875 {
00876 if (mIsInDefineMode)
00877 {
00878 return;
00879 }
00880
00881 H5Dclose(mDatasetId);
00882 if (mIsUnlimitedDimensionSet)
00883 {
00884 H5Dclose(mTimeDatasetId);
00885 }
00886 H5Fclose(mFileId);
00887
00888
00889 mIsInDefineMode = true;
00890 }
00891
00892 void Hdf5DataWriter::DefineUnlimitedDimension(const std::string& rVariableName,
00893 const std::string& rVariableUnits)
00894 {
00895 if (mIsUnlimitedDimensionSet)
00896 {
00897 EXCEPTION("Unlimited dimension already set. Cannot be defined twice");
00898 }
00899
00900 if (!mIsInDefineMode)
00901 {
00902 EXCEPTION("Cannot define variables when not in Define mode");
00903 }
00904
00905 mIsUnlimitedDimensionSet = true;
00906 mUnlimitedDimensionName = rVariableName;
00907 mUnlimitedDimensionUnit = rVariableUnits;
00908 }
00909
00910 void Hdf5DataWriter::AdvanceAlongUnlimitedDimension()
00911 {
00912 if (!mIsUnlimitedDimensionSet)
00913 {
00914 EXCEPTION("Trying to advance along an unlimited dimension without having defined any");
00915 }
00916
00917
00918 mDatasetDims[0]++;
00919 mNeedExtend = true;
00920
00921 mCurrentTimeStep++;
00922 }
00923
00924 void Hdf5DataWriter::PossiblyExtend()
00925 {
00926 if (mNeedExtend)
00927 {
00928 H5Dextend (mDatasetId, mDatasetDims);
00929 H5Dextend (mTimeDatasetId, mDatasetDims);
00930 }
00931 mNeedExtend = false;
00932 }
00933
00934 int Hdf5DataWriter::GetVariableByName(const std::string& rVariableName)
00935 {
00936 int id = -1;
00937
00938 for (unsigned index=0; index<mVariables.size(); index++)
00939 {
00940 if (mVariables[index].mVariableName == rVariableName)
00941 {
00942 id = index;
00943 break;
00944 }
00945 }
00946 if (id == -1)
00947 {
00948 EXCEPTION("Variable does not exist in hdf5 definitions.");
00949 }
00950 return id;
00951 }
00952
00953 bool Hdf5DataWriter::ApplyPermutation(const std::vector<unsigned>& rPermutation)
00954 {
00955 if (!mIsInDefineMode)
00956 {
00957 EXCEPTION("Cannot define permutation when not in Define mode");
00958 }
00959
00960 if (rPermutation.empty())
00961 {
00962 return false;
00963 }
00964
00965 if (rPermutation.size() != mFileFixedDimensionSize ||
00966 rPermutation.size() != mDataFixedDimensionSize)
00967 {
00968 EXCEPTION("Permutation doesn't match the expected problem size");
00969 }
00970
00971 std::set<unsigned> permutation_pigeon_hole;
00972
00973 bool identity_map=true;
00974 for (unsigned i=0; i<mDataFixedDimensionSize; i++)
00975 {
00976 permutation_pigeon_hole.insert(rPermutation[i]);
00977 if (identity_map && i != rPermutation[i])
00978 {
00979 identity_map=false;
00980 }
00981 }
00982 if (identity_map)
00983 {
00984
00985 return false;
00986 }
00987
00988
00989
00990
00991 for (unsigned i=0; i<mDataFixedDimensionSize; i++)
00992 {
00993 if (permutation_pigeon_hole.count(i) != 1u)
00994 {
00995 EXCEPTION("Permutation vector doesn't contain a valid permutation");
00996 }
00997 }
00998
00999 assert(mSinglePermutation == NULL);
01000 assert(mDoublePermutation == NULL);
01001 PetscTools::SetupMat(mSinglePermutation, mDataFixedDimensionSize, mDataFixedDimensionSize, 2, mHi - mLo, mHi - mLo);
01002 PetscTools::SetupMat(mDoublePermutation, 2*mDataFixedDimensionSize, 2*mDataFixedDimensionSize, 4, 2*(mHi - mLo), 2*(mHi - mLo));
01003 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
01004 MatSetOption(mSinglePermutation, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
01005 MatSetOption(mDoublePermutation, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
01006 #else
01007 MatSetOption(mSinglePermutation, MAT_IGNORE_OFF_PROC_ENTRIES);
01008 MatSetOption(mDoublePermutation, MAT_IGNORE_OFF_PROC_ENTRIES);
01009 #endif
01010
01011 for (unsigned row_index=mLo; row_index<mHi; row_index++)
01012 {
01013
01014 MatSetValue(mSinglePermutation, row_index, row_index, 0.0, INSERT_VALUES);
01015
01016 MatSetValue(mSinglePermutation, row_index, rPermutation[row_index], 1.0, INSERT_VALUES);
01017
01018 unsigned bi_index=2*row_index;
01019 unsigned perm_index=2*rPermutation[row_index];
01020
01021 MatSetValue(mDoublePermutation, bi_index, bi_index, 0.0, INSERT_VALUES);
01022 MatSetValue(mDoublePermutation, bi_index+1, bi_index+1, 0.0, INSERT_VALUES);
01023
01024 MatSetValue(mDoublePermutation, bi_index, perm_index, 1.0, INSERT_VALUES);
01025 MatSetValue(mDoublePermutation, bi_index+1, perm_index+1, 1.0, INSERT_VALUES);
01026 }
01027 MatAssemblyBegin(mSinglePermutation, MAT_FINAL_ASSEMBLY);
01028 MatAssemblyBegin(mDoublePermutation, MAT_FINAL_ASSEMBLY);
01029 MatAssemblyEnd(mSinglePermutation, MAT_FINAL_ASSEMBLY);
01030 MatAssemblyEnd(mDoublePermutation, MAT_FINAL_ASSEMBLY);
01031 return true;
01032 }