Chaste Release::3.1
|
00001 /* 00002 00003 Copyright (c) 2005-2012, University of Oxford. 00004 All rights reserved. 00005 00006 University of Oxford means the Chancellor, Masters and Scholars of the 00007 University of Oxford, having an administrative office at Wellington 00008 Square, Oxford OX1 2JD, UK. 00009 00010 This file is part of Chaste. 00011 00012 Redistribution and use in source and binary forms, with or without 00013 modification, are permitted provided that the following conditions are met: 00014 * Redistributions of source code must retain the above copyright notice, 00015 this list of conditions and the following disclaimer. 00016 * Redistributions in binary form must reproduce the above copyright notice, 00017 this list of conditions and the following disclaimer in the documentation 00018 and/or other materials provided with the distribution. 00019 * Neither the name of the University of Oxford nor the names of its 00020 contributors may be used to endorse or promote products derived from this 00021 software without specific prior written permission. 00022 00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" 00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE 00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE 00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR 00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE 00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) 00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT 00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT 00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00033 00034 */ 00035 00036 /* 00037 * Implementation file for Hdf5DataWriter class. 00038 * 00039 */ 00040 #include <set> 00041 #include <cstring> //For strcmp etc. Needed in gcc-4.4 00042 00043 #include "Hdf5DataWriter.hpp" 00044 00045 #include "Exception.hpp" 00046 #include "OutputFileHandler.hpp" 00047 #include "PetscTools.hpp" 00048 #include "Version.hpp" 00049 00050 Hdf5DataWriter::Hdf5DataWriter(DistributedVectorFactory& rVectorFactory, 00051 const std::string& rDirectory, 00052 const std::string& rBaseName, 00053 bool cleanDirectory, 00054 bool extendData) 00055 : mrVectorFactory(rVectorFactory), 00056 mDirectory(rDirectory), 00057 mBaseName(rBaseName), 00058 mCleanDirectory(cleanDirectory), 00059 mIsInDefineMode(true), 00060 mIsFixedDimensionSet(false), 00061 mIsUnlimitedDimensionSet(false), 00062 mEstimatedUnlimitedLength(1u), 00063 mFileFixedDimensionSize(0u), 00064 mDataFixedDimensionSize(0u), 00065 mLo(mrVectorFactory.GetLow()), 00066 mHi(mrVectorFactory.GetHigh()), 00067 mNumberOwned(0u), 00068 mOffset(0u), 00069 mIsDataComplete(true), 00070 mNeedExtend(false), 00071 mUseMatrixForIncompleteData(false), 00072 mCurrentTimeStep(0), 00073 mSinglePermutation(NULL), 00074 mDoublePermutation(NULL), 00075 mSingleIncompleteOutputMatrix(NULL), 00076 mDoubleIncompleteOutputMatrix(NULL), 00077 mUseOptimalChunkSizeAlgorithm(true), 00078 mFixedChunkSize(0) 00079 { 00080 if (extendData && cleanDirectory) 00081 { 00082 EXCEPTION("You are asking to delete a file and then extend it, change arguments to constructor."); 00083 } 00084 00085 if (extendData) 00086 { 00087 // Where to find the file 00088 OutputFileHandler output_file_handler(mDirectory, false); 00089 std::string results_dir = output_file_handler.GetOutputDirectoryFullPath(); 00090 std::string file_name = results_dir + mBaseName + ".h5"; 00091 00092 FileFinder h5_file(file_name,RelativeTo::Absolute); 00093 if (!h5_file.Exists()) 00094 { 00095 EXCEPTION("Hdf5DataWriter could not open " + file_name + " , as it does not exist."); 00096 } 00097 00098 // Set up a property list saying how we'll open the file 00099 hid_t property_list_id = H5Pcreate(H5P_FILE_ACCESS); 00100 H5Pset_fapl_mpiposix(property_list_id, PETSC_COMM_WORLD, 0); 00101 00102 // Open the file and free the property list 00103 mFileId = H5Fopen(file_name.c_str(), H5F_ACC_RDWR, property_list_id); 00104 H5Pclose(property_list_id); 00105 00106 if (mFileId < 0) 00107 { 00108 mDatasetId = 0; 00109 EXCEPTION("Hdf5DataWriter could not open " << file_name << 00110 " , H5Fopen error code = " << mFileId); 00111 } 00112 00113 // Open the main dataset, and figure out its size/shape 00114 mDatasetId = H5Dopen(mFileId, "Data"); 00115 hid_t variables_dataspace = H5Dget_space(mDatasetId); 00116 //unsigned variables_dataset_rank = H5Sget_simple_extent_ndims(variables_dataspace); 00117 hsize_t dataset_max_sizes[DATASET_DIMS]; 00118 H5Sget_simple_extent_dims(variables_dataspace, mDatasetDims, dataset_max_sizes); 00119 H5Sclose(variables_dataspace); 00120 00121 // Check that an unlimited dimension is defined 00122 if (dataset_max_sizes[0] != H5S_UNLIMITED) 00123 { 00124 H5Dclose(mDatasetId); 00125 H5Fclose(mFileId); 00126 EXCEPTION("Tried to open a datafile for extending which doesn't have an unlimited dimension."); 00127 } 00128 mIsUnlimitedDimensionSet = true; 00129 00130 // Sanity check other dimension sizes 00131 for (unsigned i=1; i<DATASET_DIMS; i++) // Zero is excluded since it is unlimited 00132 { 00133 assert(mDatasetDims[i] == dataset_max_sizes[i]); 00134 } 00135 mFileFixedDimensionSize = mDatasetDims[1]; 00136 mIsFixedDimensionSet = true; 00137 mVariables.reserve(mDatasetDims[2]); 00138 00139 // Figure out what the variables are 00140 hid_t attribute_id = H5Aopen_name(mDatasetId, "Variable Details"); 00141 hid_t attribute_type = H5Aget_type(attribute_id); 00142 hid_t attribute_space = H5Aget_space(attribute_id); 00143 hsize_t attr_dataspace_dim; 00144 H5Sget_simple_extent_dims(attribute_space, &attr_dataspace_dim, NULL); 00145 unsigned num_columns = H5Sget_simple_extent_npoints(attribute_space); 00146 assert(num_columns == mDatasetDims[2]); // I think... 00147 00148 char* string_array = (char *)malloc(sizeof(char)*MAX_STRING_SIZE*(int)num_columns); 00149 H5Aread(attribute_id, attribute_type, string_array); 00150 00151 // Loop over columns/variables 00152 for (unsigned index=0; index<num_columns; index++) 00153 { 00154 // Read the string from the raw vector 00155 std::string column_name_unit(&string_array[MAX_STRING_SIZE*index]); 00156 00157 // Find location of unit name 00158 size_t name_length = column_name_unit.find('('); 00159 size_t unit_length = column_name_unit.find(')') - name_length - 1; 00160 00161 // Create the variable 00162 DataWriterVariable var; 00163 var.mVariableName = column_name_unit.substr(0, name_length); 00164 var.mVariableUnits = column_name_unit.substr(name_length+1, unit_length); 00165 mVariables.push_back(var); 00166 } 00167 00168 // Free memory, release ids 00169 free(string_array); 00170 H5Tclose(attribute_type); 00171 H5Sclose(attribute_space); 00172 H5Aclose(attribute_id); 00173 00174 // Now deal with time 00175 mUnlimitedDimensionName = "Time"; // Assumed by the reader... 00176 mTimeDatasetId = H5Dopen(mFileId, mUnlimitedDimensionName.c_str()); 00177 mUnlimitedDimensionUnit = "ms"; // Assumed by Chaste... 00178 00179 // How many time steps have been written so far? 00180 hid_t timestep_dataspace = H5Dget_space(mTimeDatasetId); 00181 hsize_t num_timesteps; 00182 H5Sget_simple_extent_dims(timestep_dataspace, &num_timesteps, NULL); 00183 H5Sclose(timestep_dataspace); 00184 mCurrentTimeStep = (long)num_timesteps - 1; 00185 00186 // Incomplete data? 00187 attribute_id = H5Aopen_name(mDatasetId, "IsDataComplete"); 00188 if (attribute_id < 0) 00189 { 00190 #define COVERAGE_IGNORE 00191 // Old format, before we added the attribute. 00192 EXCEPTION("Extending old-format files isn't supported."); 00193 #undef COVERAGE_IGNORE 00194 } 00195 else 00196 { 00197 attribute_type = H5Aget_type(attribute_id); 00198 attribute_space = H5Aget_space(attribute_id); 00199 unsigned is_data_complete; 00200 H5Aread(attribute_id, H5T_NATIVE_UINT, &is_data_complete); 00201 mIsDataComplete = (is_data_complete == 1) ? true : false; 00202 H5Tclose(attribute_type); 00203 H5Sclose(attribute_space); 00204 H5Aclose(attribute_id); 00205 } 00206 if (mIsDataComplete) 00207 { 00208 mNumberOwned = mrVectorFactory.GetLocalOwnership(); 00209 mOffset = mLo; 00210 mDataFixedDimensionSize = mFileFixedDimensionSize; 00211 } 00212 else 00213 { 00214 // Read which nodes appear in the file (mIncompleteNodeIndices) 00215 attribute_id = H5Aopen_name(mDatasetId, "NodeMap"); 00216 attribute_type = H5Aget_type(attribute_id); 00217 attribute_space = H5Aget_space(attribute_id); 00218 00219 // Get the dataset/dataspace dimensions 00220 unsigned num_node_indices = H5Sget_simple_extent_npoints(attribute_space); 00221 00222 // Read data from hyperslab in the file into the hyperslab in memory 00223 mIncompleteNodeIndices.clear(); 00224 mIncompleteNodeIndices.resize(num_node_indices); 00225 H5Aread(attribute_id, H5T_NATIVE_UINT, &mIncompleteNodeIndices[0]); 00226 00227 // Release ids 00228 H5Tclose(attribute_type); 00229 H5Sclose(attribute_space); 00230 H5Aclose(attribute_id); 00231 00232 // Set up what data we can 00233 mNumberOwned = mrVectorFactory.GetLocalOwnership(); 00234 ComputeIncompleteOffset(); 00238 mDataFixedDimensionSize = UINT_MAX; 00239 H5Dclose(mDatasetId); 00240 H5Dclose(mTimeDatasetId); 00241 H5Fclose(mFileId); 00242 EXCEPTION("Unable to extend an incomplete data file at present."); 00243 } 00244 00245 // Done 00246 mIsInDefineMode = false; 00247 AdvanceAlongUnlimitedDimension(); 00248 } 00249 } 00250 00251 Hdf5DataWriter::~Hdf5DataWriter() 00252 { 00253 Close(); 00254 00255 if (mSinglePermutation) 00256 { 00257 PetscTools::Destroy(mSinglePermutation); 00258 } 00259 if (mDoublePermutation) 00260 { 00261 PetscTools::Destroy(mDoublePermutation); 00262 } 00263 if (mSingleIncompleteOutputMatrix) 00264 { 00265 PetscTools::Destroy(mSingleIncompleteOutputMatrix); 00266 } 00267 if (mDoubleIncompleteOutputMatrix) 00268 { 00269 PetscTools::Destroy(mDoubleIncompleteOutputMatrix); 00270 } 00271 } 00272 00273 void Hdf5DataWriter::DefineFixedDimension(long dimensionSize) 00274 { 00275 if (!mIsInDefineMode) 00276 { 00277 EXCEPTION("Cannot define variables when not in Define mode"); 00278 } 00279 if (dimensionSize < 1) 00280 { 00281 EXCEPTION("Fixed dimension must be at least 1 long"); 00282 } 00283 if (mIsFixedDimensionSet) 00284 { 00285 EXCEPTION("Fixed dimension already set"); 00286 } 00287 00288 // Work out the ownership details 00289 mLo = mrVectorFactory.GetLow(); 00290 mHi = mrVectorFactory.GetHigh(); 00291 mNumberOwned = mrVectorFactory.GetLocalOwnership(); 00292 mOffset = mLo; 00293 mFileFixedDimensionSize = dimensionSize; 00294 mDataFixedDimensionSize = dimensionSize; 00295 mIsFixedDimensionSet = true; 00296 } 00297 00298 void Hdf5DataWriter::DefineFixedDimension(const std::vector<unsigned>& rNodesToOuput, long vecSize) 00299 { 00300 unsigned vector_size = rNodesToOuput.size(); 00301 00302 for (unsigned index=0; index < vector_size-1; index++) 00303 { 00304 if (rNodesToOuput[index] >= rNodesToOuput[index+1]) 00305 { 00306 EXCEPTION("Input should be monotonic increasing"); 00307 } 00308 } 00309 00310 if ((int) rNodesToOuput.back() >= vecSize) 00311 { 00312 EXCEPTION("Vector size doesn't match nodes to output"); 00313 } 00314 00315 DefineFixedDimension(vecSize); 00316 00317 mFileFixedDimensionSize = vector_size; 00318 mIsDataComplete = false; 00319 mIncompleteNodeIndices = rNodesToOuput; 00320 ComputeIncompleteOffset(); 00321 } 00322 00323 void Hdf5DataWriter::DefineFixedDimensionUsingMatrix(const std::vector<unsigned>& rNodesToOuput, long vecSize) 00324 { 00325 unsigned vector_size = rNodesToOuput.size(); 00326 00327 for (unsigned index=0; index < vector_size-1; index++) 00328 { 00329 if (rNodesToOuput[index] >= rNodesToOuput[index+1]) 00330 { 00331 EXCEPTION("Input should be monotonic increasing"); 00332 } 00333 } 00334 00335 if ((int) rNodesToOuput.back() >= vecSize) 00336 { 00337 EXCEPTION("Vector size doesn't match nodes to output"); 00338 } 00339 00340 DefineFixedDimension(vecSize); 00341 00342 mFileFixedDimensionSize = vector_size; 00343 mIsDataComplete = false; 00344 mIncompleteNodeIndices = rNodesToOuput; 00345 mUseMatrixForIncompleteData = true; 00346 ComputeIncompleteOffset(); 00347 00348 // Make sure we've not done it already 00349 assert(mSingleIncompleteOutputMatrix == NULL); 00350 assert(mDoubleIncompleteOutputMatrix == NULL); 00351 PetscTools::SetupMat(mSingleIncompleteOutputMatrix, mFileFixedDimensionSize, mDataFixedDimensionSize, 2, mNumberOwned, mHi - mLo); 00352 PetscTools::SetupMat(mDoubleIncompleteOutputMatrix, 2*mFileFixedDimensionSize, 2*mDataFixedDimensionSize, 4, 2*mNumberOwned, 2*(mHi - mLo)); 00353 00354 // Only do local rows 00355 for (unsigned row_index = mOffset; row_index < mOffset + mNumberOwned; row_index++) 00356 { 00357 // Put zero on the diagonal 00358 MatSetValue(mSingleIncompleteOutputMatrix, row_index, row_index, 0.0, INSERT_VALUES); 00359 00360 // Put one at (i,j) 00361 MatSetValue(mSingleIncompleteOutputMatrix, row_index, rNodesToOuput[row_index], 1.0, INSERT_VALUES); 00362 00363 unsigned bi_index = 2*row_index; 00364 unsigned perm_index = 2*rNodesToOuput[row_index]; 00365 00366 // Put zeroes on the diagonal 00367 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index, bi_index, 0.0, INSERT_VALUES); 00368 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index+1, bi_index+1, 0.0, INSERT_VALUES); 00369 00370 // Put ones at (i,j) 00371 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index, perm_index, 1.0, INSERT_VALUES); 00372 MatSetValue(mDoubleIncompleteOutputMatrix, bi_index+1, perm_index+1, 1.0, INSERT_VALUES); 00373 } 00374 MatAssemblyBegin(mSingleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY); 00375 MatAssemblyBegin(mDoubleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY); 00376 MatAssemblyEnd(mSingleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY); 00377 MatAssemblyEnd(mDoubleIncompleteOutputMatrix, MAT_FINAL_ASSEMBLY); 00378 00379 // MatView(mSingleIncompleteOutputMatrix, PETSC_VIEWER_STDOUT_WORLD); 00380 } 00381 00382 void Hdf5DataWriter::ComputeIncompleteOffset() 00383 { 00384 mOffset = 0; 00385 mNumberOwned = 0; 00386 for (unsigned i=0; i<mIncompleteNodeIndices.size(); i++) 00387 { 00388 if (mIncompleteNodeIndices[i] < mLo) 00389 { 00390 mOffset++; 00391 } 00392 else if (mIncompleteNodeIndices[i] < mHi) 00393 { 00394 mNumberOwned++; 00395 } 00396 } 00397 } 00398 00399 int Hdf5DataWriter::DefineVariable(const std::string& rVariableName, 00400 const std::string& rVariableUnits) 00401 { 00402 if (!mIsInDefineMode) 00403 { 00404 EXCEPTION("Cannot define variables when not in Define mode"); 00405 } 00406 00407 CheckVariableName(rVariableName); 00408 CheckUnitsName(rVariableUnits); 00409 00410 // Check for the variable being already defined 00411 for (unsigned index=0; index<mVariables.size(); index++) 00412 { 00413 if (mVariables[index].mVariableName == rVariableName) 00414 { 00415 EXCEPTION("Variable name already exists"); 00416 } 00417 } 00418 00419 DataWriterVariable new_variable; 00420 new_variable.mVariableName = rVariableName; 00421 new_variable.mVariableUnits = rVariableUnits; 00422 int variable_id; 00423 00424 // Add the variable to the variable vector 00425 mVariables.push_back(new_variable); 00426 00427 // Use the index of the variable vector as the variable ID. 00428 // This is ok since there is no way to remove variables. 00429 variable_id = mVariables.size() - 1; 00430 00431 return variable_id; 00432 } 00433 00434 void Hdf5DataWriter::CheckVariableName(const std::string& rName) 00435 { 00436 if (rName.length() == 0) 00437 { 00438 EXCEPTION("Variable name not allowed: may not be blank."); 00439 } 00440 CheckUnitsName(rName); 00441 } 00442 00443 void Hdf5DataWriter::CheckUnitsName(const std::string& rName) 00444 { 00445 for (unsigned i=0; i<rName.length(); i++) 00446 { 00447 if (!isalnum(rName[i]) && !(rName[i]=='_')) 00448 { 00449 std::string error = "Variable name/units '" + rName + "' not allowed: may only contain alphanumeric characters or '_'."; 00450 EXCEPTION(error); 00451 } 00452 } 00453 } 00454 00455 void Hdf5DataWriter::EndDefineMode() 00456 { 00457 // Check that at least one variable has been defined 00458 if (mVariables.size() < 1) 00459 { 00460 EXCEPTION("Cannot end define mode. No variables have been defined."); 00461 } 00462 00463 // Check that a fixed dimension has been defined 00464 if (mIsFixedDimensionSet == false) 00465 { 00466 EXCEPTION("Cannot end define mode. One fixed dimension should be defined."); 00467 } 00468 00469 OutputFileHandler output_file_handler(mDirectory, mCleanDirectory); 00470 std::string results_dir = output_file_handler.GetOutputDirectoryFullPath(); 00471 std::string file_name = results_dir + mBaseName + ".h5"; 00472 00473 // Set up a property list saying how we'll open the file 00474 hid_t property_list_id = H5Pcreate(H5P_FILE_ACCESS); 00475 H5Pset_fapl_mpiposix(property_list_id, PETSC_COMM_WORLD, 0); 00476 00477 // Create a file (collectively) and free the property list 00478 mFileId = H5Fcreate(file_name.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, property_list_id); 00479 H5Pclose(property_list_id); 00480 if (mFileId < 0) 00481 { 00482 EXCEPTION("Hdf5DataWriter could not create " << file_name << 00483 " , H5Fcreate error code = " << mFileId); 00484 } 00485 mIsInDefineMode = false; 00486 00487 /* 00488 * Create "Data" dataset 00489 */ 00490 00491 // Create the dataspace for the dataset. 00492 mDatasetDims[0] = 1; // While developing we got a non-documented "only the first dimension can be extendible" error. 00493 mDatasetDims[1] = mFileFixedDimensionSize; 00494 mDatasetDims[2] = mVariables.size(); 00495 00496 hsize_t* max_dims = NULL; 00497 hsize_t dataset_max_dims[DATASET_DIMS]; // dataset max dimensions 00498 00499 hid_t cparms = H5P_DEFAULT; 00500 00501 if (mIsUnlimitedDimensionSet) 00502 { 00503 dataset_max_dims[0] = H5S_UNLIMITED; 00504 dataset_max_dims[1] = mDatasetDims[1]; 00505 dataset_max_dims[2] = mDatasetDims[2]; 00506 max_dims = dataset_max_dims; 00507 00508 hsize_t chunk_size; 00509 if (mUseOptimalChunkSizeAlgorithm) 00510 { 00511 /* 00512 * Modify dataset creation properties to enable chunking. We don't want 00513 * more than 100 chunks, as performance degrades significantly if there 00514 * are too many, where "too many" appears to be about 1000. HDF5's caching 00515 * won't apply if the chunks are too large, but this seems to have less of 00516 * an impact. 00517 */ 00518 chunk_size = mEstimatedUnlimitedLength/100; 00519 if (chunk_size < 100) 00520 { 00521 chunk_size = 100; 00522 } 00523 } 00524 else 00525 { 00526 chunk_size = mFixedChunkSize; 00527 } 00528 /* 00529 * If the size of a chunk in bytes is bigger than 4GB then there may be problems. 00530 * HDF5 1.6.x does not check for this - but there may be snags further down the line. 00531 * HDF5 1.8.x does more error checking and produces std::cout errors and a file with 00532 * no data in it. 00533 */ 00534 if (chunk_size * mDatasetDims[1] * mDatasetDims[2] > (uint64_t)0xffffffff) 00535 { 00536 /* 00537 * Note: this exception can be avoided by altering the lines above 00538 * where chunk_size is set (at a loss of efficiency). 00539 */ 00540 mIsInDefineMode = true; // To stop things that would be created below from being deleted on Close() 00541 H5Fclose(mFileId); // This is the one thing which we have made 00542 EXCEPTION("HDF5 may be writing more than 4GB to disk at any time and would fail. It may be possible to tune the Chaste code to get around this"); 00543 } 00544 00545 hsize_t chunk_dims[DATASET_DIMS] = {chunk_size, mDatasetDims[1], mDatasetDims[2]}; 00546 cparms = H5Pcreate (H5P_DATASET_CREATE); 00547 00548 H5Pset_chunk( cparms, DATASET_DIMS, chunk_dims); 00549 } 00550 00551 hid_t filespace = H5Screate_simple(DATASET_DIMS, mDatasetDims, max_dims); 00552 00553 // Create the dataset and close filespace 00554 mDatasetId = H5Dcreate(mFileId, "Data", H5T_NATIVE_DOUBLE, filespace, cparms); 00555 H5Sclose(filespace); 00556 00557 // Create dataspace for the name, unit attribute 00558 const unsigned MAX_STRING_SIZE = 100; 00559 hsize_t columns[1] = {mVariables.size()}; 00560 hid_t colspace = H5Screate_simple(1, columns, NULL); 00561 00562 // Create attribute for variable names 00563 char* col_data = (char*) malloc(mVariables.size() * sizeof(char) * MAX_STRING_SIZE); 00564 00565 char* col_data_offset = col_data; 00566 for (unsigned var=0; var<mVariables.size(); var++) 00567 { 00568 std::string full_name = mVariables[var].mVariableName + "(" + mVariables[var].mVariableUnits + ")"; 00569 strcpy (col_data_offset, full_name.c_str()); 00570 col_data_offset += sizeof(char) * MAX_STRING_SIZE; 00571 } 00572 00573 // Create the type 'string' 00574 hid_t string_type = H5Tcopy(H5T_C_S1); 00575 H5Tset_size(string_type, MAX_STRING_SIZE ); 00576 hid_t attr = H5Acreate(mDatasetId, "Variable Details", string_type, colspace, H5P_DEFAULT); 00577 00578 // Write to the attribute 00579 H5Awrite(attr, string_type, col_data); 00580 00581 // Close dataspace & attribute 00582 free(col_data); 00583 H5Sclose(colspace); 00584 H5Aclose(attr); 00585 00586 // Create "boolean" attribute telling the data to be incomplete or not 00587 columns[0] = 1; 00588 colspace = H5Screate_simple(1, columns, NULL); 00589 attr = H5Acreate(mDatasetId, "IsDataComplete", H5T_NATIVE_UINT, colspace, H5P_DEFAULT); 00590 00591 // Write to the attribute - note that the native boolean is not predictable 00592 unsigned is_data_complete = mIsDataComplete ? 1 : 0; 00593 H5Awrite(attr, H5T_NATIVE_UINT, &is_data_complete); 00594 00595 H5Sclose(colspace); 00596 H5Aclose(attr); 00597 00598 if (!mIsDataComplete) 00599 { 00600 // We need to write a map 00601 // Create "unsigned" attribute with the map 00602 columns[0] = mFileFixedDimensionSize; 00603 colspace = H5Screate_simple(1, columns, NULL); 00604 attr = H5Acreate(mDatasetId, "NodeMap", H5T_NATIVE_UINT, colspace, H5P_DEFAULT); 00605 00606 // Write to the attribute 00607 H5Awrite(attr, H5T_NATIVE_UINT, &mIncompleteNodeIndices[0]); 00608 00609 H5Sclose(colspace); 00610 H5Aclose(attr); 00611 } 00612 00613 /* 00614 * Create "Time" dataset 00615 */ 00616 if (mIsUnlimitedDimensionSet) 00617 { 00618 hsize_t time_dataset_dims[1] = {1}; 00619 hsize_t time_dataset_max_dims[1] = {H5S_UNLIMITED}; 00620 00621 // Modify dataset creation properties to enable chunking. 00622 hsize_t time_chunk_dims[1] ={1}; 00623 hid_t time_cparms = H5Pcreate (H5P_DATASET_CREATE); 00624 H5Pset_chunk( time_cparms, 1, time_chunk_dims); 00625 00626 hid_t time_filespace = H5Screate_simple(1, time_dataset_dims, time_dataset_max_dims); 00627 00628 // Create the dataset 00629 mTimeDatasetId = H5Dcreate(mFileId, mUnlimitedDimensionName.c_str(), H5T_NATIVE_DOUBLE, time_filespace, time_cparms); 00630 00631 // Create the dataspace for the attribute 00632 hsize_t one = 1; 00633 hid_t one_column_space = H5Screate_simple(1, &one, NULL); 00634 00635 // Create an attribute for the time unit 00636 hid_t unit_attr = H5Acreate(mTimeDatasetId, "Unit", string_type, one_column_space, H5P_DEFAULT); 00637 00638 // Copy the unit to a string MAX_STRING_SIZE long and write it 00639 char unit_string[MAX_STRING_SIZE]; 00640 strcpy(unit_string, mUnlimitedDimensionUnit.c_str()); 00641 H5Awrite(unit_attr, string_type, unit_string); 00642 00643 // Close the filespace 00644 H5Sclose(one_column_space); 00645 H5Aclose(unit_attr); 00646 H5Sclose(time_filespace); 00647 } 00648 00649 /* 00650 * Create the provenance attribute 00651 */ 00652 00653 // Create a longer type for 'string' 00654 const unsigned MAX_PROVENANCE_STRING_SIZE = 1023; 00655 hid_t long_string_type = H5Tcopy(H5T_C_S1); 00656 H5Tset_size(long_string_type, MAX_PROVENANCE_STRING_SIZE ); 00657 hsize_t prov_columns[1] = {1}; 00658 hid_t provenance_space = H5Screate_simple(1, prov_columns, NULL); 00659 char* provenance_data = (char*) malloc(sizeof(char) * MAX_PROVENANCE_STRING_SIZE); 00660 assert( ChasteBuildInfo::GetProvenanceString().length() < MAX_PROVENANCE_STRING_SIZE); 00661 00662 strcpy(provenance_data, ChasteBuildInfo::GetProvenanceString().c_str()); 00663 hid_t prov_attr = H5Acreate(mDatasetId, "Chaste Provenance", long_string_type, provenance_space, H5P_DEFAULT); 00664 00665 // Write to the attribute 00666 H5Awrite(prov_attr, long_string_type, provenance_data); 00667 00668 // Close dataspace & attribute 00669 free(provenance_data); 00670 H5Sclose(provenance_space); 00671 H5Aclose(prov_attr); 00672 } 00673 00674 void Hdf5DataWriter::PutVector(int variableID, Vec petscVector) 00675 { 00676 if (mIsInDefineMode) 00677 { 00678 EXCEPTION("Cannot write data while in define mode."); 00679 } 00680 00681 int vector_size; 00682 VecGetSize(petscVector, &vector_size); 00683 00684 if ((unsigned) vector_size != mDataFixedDimensionSize) 00685 { 00686 EXCEPTION("Vector size doesn't match fixed dimension"); 00687 } 00688 00689 // Make sure that everything is actually extended to the correct dimension. 00690 PossiblyExtend(); 00691 00692 Vec output_petsc_vector; 00693 00694 // Decide what to write 00695 if (mSinglePermutation == NULL) 00696 { 00697 // No permutation - just write 00698 output_petsc_vector = petscVector; 00699 } 00700 else 00701 { 00702 assert(mIsDataComplete); 00703 // Make a vector with the same pattern (doesn't copy the data) 00704 VecDuplicate(petscVector, &output_petsc_vector); 00705 00706 // Apply the permutation matrix 00707 MatMult(mSinglePermutation, petscVector, output_petsc_vector); 00708 } 00709 00710 // Define a dataset in memory for this process 00711 hid_t memspace=0; 00712 if (mNumberOwned != 0) 00713 { 00714 hsize_t v_size[1] = {mNumberOwned}; 00715 memspace = H5Screate_simple(1, v_size, NULL); 00716 } 00717 00718 // Select hyperslab in the file 00719 hsize_t count[DATASET_DIMS] = {1, mNumberOwned, 1}; 00720 hsize_t offset_dims[DATASET_DIMS] = {mCurrentTimeStep, mOffset, variableID}; 00721 hid_t file_dataspace = H5Dget_space(mDatasetId); 00722 00723 // Create property list for collective dataset 00724 hid_t property_list_id = H5Pcreate(H5P_DATASET_XFER); 00725 H5Pset_dxpl_mpio(property_list_id, H5FD_MPIO_COLLECTIVE); 00726 00727 H5Sselect_hyperslab(file_dataspace, H5S_SELECT_SET, offset_dims, NULL, count, NULL); 00728 00729 double* p_petsc_vector; 00730 VecGetArray(output_petsc_vector, &p_petsc_vector); 00731 00732 if (mIsDataComplete) 00733 { 00734 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, file_dataspace, property_list_id, p_petsc_vector); 00735 } 00736 else 00737 { 00738 if (mUseMatrixForIncompleteData) 00739 { 00740 // Make a vector of the required size 00741 output_petsc_vector = PetscTools::CreateVec(mFileFixedDimensionSize, mNumberOwned); 00742 00743 // Fill the vector by multiplying complete data by incomplete output matrix 00744 MatMult(mSingleIncompleteOutputMatrix, petscVector, output_petsc_vector); 00745 00746 double* p_petsc_vector_incomplete; 00747 VecGetArray(output_petsc_vector, &p_petsc_vector_incomplete); 00748 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, file_dataspace, property_list_id, p_petsc_vector_incomplete); 00749 } 00750 else 00751 { 00752 // Make a local copy of the data you own 00753 double local_data[mNumberOwned]; 00754 for (unsigned i=0; i<mNumberOwned; i++) 00755 { 00756 local_data[i] = p_petsc_vector[ mIncompleteNodeIndices[mOffset+i]-mLo ]; 00757 00758 } 00759 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, file_dataspace, property_list_id, local_data); 00760 } 00761 } 00762 00763 VecRestoreArray(output_petsc_vector, &p_petsc_vector); 00764 00765 H5Pclose(property_list_id); 00766 H5Sclose(file_dataspace); 00767 if (mNumberOwned !=0) 00768 { 00769 H5Sclose(memspace); 00770 } 00771 00772 if (petscVector != output_petsc_vector) 00773 { 00774 // Free local vector 00775 PetscTools::Destroy(output_petsc_vector); 00776 } 00777 } 00778 00779 void Hdf5DataWriter::PutStripedVector(std::vector<int> variableIDs, Vec petscVector) 00780 { 00781 if (mIsInDefineMode) 00782 { 00783 EXCEPTION("Cannot write data while in define mode."); 00784 } 00785 00786 if (variableIDs.size() <= 1) 00787 { 00788 EXCEPTION("The PutStripedVector method requires at least two variables ID. If only one is needed, use PutVector method instead"); 00789 } 00790 00791 const int NUM_STRIPES=variableIDs.size(); 00792 00793 int firstVariableID=variableIDs[0]; 00794 00795 // Currently the method only works with consecutive columns, can be extended if needed 00796 for (unsigned i = 1; i < variableIDs.size(); i++) 00797 { 00798 if (variableIDs[i]-variableIDs[i-1] != 1) 00799 { 00800 EXCEPTION("Columns should be consecutive. Try reordering them."); 00801 } 00802 } 00803 00804 int vector_size; 00805 VecGetSize(petscVector, &vector_size); 00806 00807 if ((unsigned) vector_size != NUM_STRIPES*mDataFixedDimensionSize) 00808 { 00809 EXCEPTION("Vector size doesn't match fixed dimension"); 00810 } 00811 00812 // Make sure that everything is actually extended to the correct dimension 00813 PossiblyExtend(); 00814 00815 Vec output_petsc_vector; 00816 00817 // Decide what to write 00818 if (mDoublePermutation == NULL) 00819 { 00820 // No permutation - just write 00821 output_petsc_vector = petscVector; 00822 } 00823 else 00824 { 00825 assert(mIsDataComplete); 00826 // Make a vector with the same pattern (doesn't copy the data) 00827 VecDuplicate(petscVector, &output_petsc_vector); 00828 00829 // Apply the permutation matrix 00830 MatMult(mDoublePermutation, petscVector, output_petsc_vector); 00831 } 00832 // Define a dataset in memory for this process 00833 hid_t memspace=0; 00834 if (mNumberOwned !=0) 00835 { 00836 hsize_t v_size[1] = {mNumberOwned*NUM_STRIPES}; 00837 memspace = H5Screate_simple(1, v_size, NULL); 00838 } 00839 00840 // Select hyperslab in the file 00841 hsize_t start[DATASET_DIMS] = {mCurrentTimeStep, mOffset, firstVariableID}; 00842 hsize_t stride[DATASET_DIMS] = {1, 1, 1};//we are imposing contiguous variables, hence the stride is 1 (3rd component) 00843 hsize_t block_size[DATASET_DIMS] = {1, mNumberOwned, 1}; 00844 hsize_t number_blocks[DATASET_DIMS] = {1, 1, NUM_STRIPES}; 00845 00846 hid_t hyperslab_space = H5Dget_space(mDatasetId); 00847 H5Sselect_hyperslab(hyperslab_space, H5S_SELECT_SET, start, stride, number_blocks, block_size); 00848 00849 // Create property list for collective dataset write, and write! Finally. 00850 hid_t property_list_id = H5Pcreate(H5P_DATASET_XFER); 00851 H5Pset_dxpl_mpio(property_list_id, H5FD_MPIO_COLLECTIVE); 00852 00853 double* p_petsc_vector; 00854 VecGetArray(output_petsc_vector, &p_petsc_vector); 00855 00856 if (mIsDataComplete) 00857 { 00858 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, property_list_id, p_petsc_vector); 00859 } 00860 else 00861 { 00862 if (variableIDs.size() < 3) // incomplete data and striped vector is supported only for NUM_STRIPES=2...for the moment 00863 { 00864 if (mUseMatrixForIncompleteData) 00865 { 00866 // Make a vector of the required size 00867 output_petsc_vector = PetscTools::CreateVec(2*mFileFixedDimensionSize, 2*mNumberOwned); 00868 00869 // Fill the vector by multiplying complete data by incomplete output matrix 00870 MatMult(mDoubleIncompleteOutputMatrix, petscVector, output_petsc_vector); 00871 00872 double* p_petsc_vector_incomplete; 00873 VecGetArray(output_petsc_vector, &p_petsc_vector_incomplete); 00874 00875 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, property_list_id, p_petsc_vector_incomplete); 00876 } 00877 else 00878 { 00879 // Make a local copy of the data you own 00880 double local_data[mNumberOwned*NUM_STRIPES]; 00881 for (unsigned i=0; i<mNumberOwned; i++) 00882 { 00883 unsigned local_node_number = mIncompleteNodeIndices[mOffset+i] - mLo; 00884 local_data[NUM_STRIPES*i] = p_petsc_vector[ local_node_number*NUM_STRIPES ]; 00885 local_data[NUM_STRIPES*i+1] = p_petsc_vector[ local_node_number*NUM_STRIPES + 1]; 00886 } 00887 00888 H5Dwrite(mDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, property_list_id, local_data); 00889 } 00890 } 00891 else 00892 { 00893 EXCEPTION("The PutStripedVector functionality for incomplete data is supported for only 2 stripes"); 00894 } 00895 } 00896 00897 VecRestoreArray(output_petsc_vector, &p_petsc_vector); 00898 00899 H5Sclose(hyperslab_space); 00900 if (mNumberOwned != 0) 00901 { 00902 H5Sclose(memspace); 00903 } 00904 H5Pclose(property_list_id); 00905 00906 if (petscVector != output_petsc_vector) 00907 { 00908 // Free local vector 00909 PetscTools::Destroy(output_petsc_vector); 00910 } 00911 } 00912 00913 void Hdf5DataWriter::PutUnlimitedVariable(double value) 00914 { 00915 if (mIsInDefineMode) 00916 { 00917 EXCEPTION("Cannot write data while in define mode."); 00918 } 00919 00920 if (!mIsUnlimitedDimensionSet) 00921 { 00922 EXCEPTION("PutUnlimitedVariable() called but no unlimited dimension has been set"); 00923 } 00924 00925 // Make sure that everything is actually extended to the correct dimension. 00926 PossiblyExtend(); 00927 00928 // This data is only written by the master 00929 if (!PetscTools::AmMaster()) 00930 { 00931 return; 00932 } 00933 00934 hsize_t size[1] = {1}; 00935 hid_t memspace = H5Screate_simple(1, size, NULL); 00936 00937 // Select hyperslab in the file. 00938 hsize_t count[1] = {1}; 00939 hsize_t offset[1] = {mCurrentTimeStep}; 00940 hid_t hyperslab_space = H5Dget_space(mTimeDatasetId); 00941 H5Sselect_hyperslab(hyperslab_space, H5S_SELECT_SET, offset, NULL, count, NULL); 00942 00943 H5Dwrite(mTimeDatasetId, H5T_NATIVE_DOUBLE, memspace, hyperslab_space, H5P_DEFAULT, &value); 00944 00945 H5Sclose(hyperslab_space); 00946 H5Sclose(memspace); 00947 } 00948 00949 void Hdf5DataWriter::Close() 00950 { 00951 if (mIsInDefineMode) 00952 { 00953 return; // Nothing to do... 00954 } 00955 00956 H5Dclose(mDatasetId); 00957 if (mIsUnlimitedDimensionSet) 00958 { 00959 H5Dclose(mTimeDatasetId); 00960 } 00961 H5Fclose(mFileId); 00962 00963 // Cope with being called twice (e.g. if a user calls Close then the destructor) 00964 mIsInDefineMode = true; 00965 } 00966 00967 void Hdf5DataWriter::DefineUnlimitedDimension(const std::string& rVariableName, 00968 const std::string& rVariableUnits, 00969 unsigned estimatedLength) 00970 { 00971 if (mIsUnlimitedDimensionSet) 00972 { 00973 EXCEPTION("Unlimited dimension already set. Cannot be defined twice"); 00974 } 00975 00976 if (!mIsInDefineMode) 00977 { 00978 EXCEPTION("Cannot define variables when not in Define mode"); 00979 } 00980 00981 mIsUnlimitedDimensionSet = true; 00982 mUnlimitedDimensionName = rVariableName; 00983 mUnlimitedDimensionUnit = rVariableUnits; 00984 mEstimatedUnlimitedLength = estimatedLength; 00985 } 00986 00987 void Hdf5DataWriter::AdvanceAlongUnlimitedDimension() 00988 { 00989 if (!mIsUnlimitedDimensionSet) 00990 { 00991 EXCEPTION("Trying to advance along an unlimited dimension without having defined any"); 00992 } 00993 00994 mCurrentTimeStep++; 00995 00996 /* 00997 * Extend the dataset. 00998 * 00999 * If the user provided an estimate for the length of the 01000 * unlimited dimension, allocate all that space. 01001 */ 01002 if ( mEstimatedUnlimitedLength > mDatasetDims[0] ) 01003 { 01004 mDatasetDims[0] = mEstimatedUnlimitedLength; 01005 mNeedExtend = true; 01006 } 01007 01008 // If you are beyond the user estimate increment step by step 01009 if ( mCurrentTimeStep >= (long) mEstimatedUnlimitedLength ) 01010 { 01011 mDatasetDims[0]++; 01012 mNeedExtend = true; 01013 } 01014 } 01015 01016 void Hdf5DataWriter::PossiblyExtend() 01017 { 01018 if (mNeedExtend) 01019 { 01020 H5Dextend (mDatasetId, mDatasetDims); 01021 H5Dextend (mTimeDatasetId, mDatasetDims); 01022 } 01023 mNeedExtend = false; 01024 } 01025 01026 int Hdf5DataWriter::GetVariableByName(const std::string& rVariableName) 01027 { 01028 int id = -1; 01029 01030 // Check for the variable name in the existing variables 01031 for (unsigned index=0; index<mVariables.size(); index++) 01032 { 01033 if (mVariables[index].mVariableName == rVariableName) 01034 { 01035 id = index; 01036 break; 01037 } 01038 } 01039 if (id == -1) 01040 { 01041 EXCEPTION("Variable does not exist in hdf5 definitions."); 01042 } 01043 return id; 01044 } 01045 01046 bool Hdf5DataWriter::ApplyPermutation(const std::vector<unsigned>& rPermutation) 01047 { 01048 if (!mIsInDefineMode) 01049 { 01050 EXCEPTION("Cannot define permutation when not in Define mode"); 01051 } 01052 01053 if (rPermutation.empty()) 01054 { 01055 return false; 01056 } 01057 01058 if (rPermutation.size() != mFileFixedDimensionSize || 01059 rPermutation.size() != mDataFixedDimensionSize) 01060 { 01061 EXCEPTION("Permutation doesn't match the expected problem size"); 01062 } 01063 01064 // Permutation checker 01065 std::set<unsigned> permutation_pigeon_hole; 01066 01067 // Fill up the pigeon holes 01068 bool identity_map = true; 01069 for (unsigned i=0; i<mDataFixedDimensionSize; i++) 01070 { 01071 permutation_pigeon_hole.insert(rPermutation[i]); 01072 if (identity_map && i != rPermutation[i]) 01073 { 01074 identity_map = false; 01075 } 01076 } 01077 if (identity_map) 01078 { 01079 // Do nothing 01080 return false; 01081 } 01082 01083 /* 01084 * The pigeon-hole principle states that each index appears exactly once 01085 * so if any don't appear then either one appears twice or something out of 01086 * scope has appeared. 01087 */ 01088 for (unsigned i=0; i<mDataFixedDimensionSize; i++) 01089 { 01090 if (permutation_pigeon_hole.count(i) != 1u) 01091 { 01092 EXCEPTION("Permutation vector doesn't contain a valid permutation"); 01093 } 01094 } 01095 // Make sure we've not done it already 01096 assert(mSinglePermutation == NULL); 01097 assert(mDoublePermutation == NULL); 01098 PetscTools::SetupMat(mSinglePermutation, mDataFixedDimensionSize, mDataFixedDimensionSize, 2, mHi - mLo, mHi - mLo); 01099 PetscTools::SetupMat(mDoublePermutation, 2*mDataFixedDimensionSize, 2*mDataFixedDimensionSize, 4, 2*(mHi - mLo), 2*(mHi - mLo)); 01100 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x 01101 MatSetOption(mSinglePermutation, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE); 01102 MatSetOption(mDoublePermutation, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE); 01103 #else 01104 MatSetOption(mSinglePermutation, MAT_IGNORE_OFF_PROC_ENTRIES); 01105 MatSetOption(mDoublePermutation, MAT_IGNORE_OFF_PROC_ENTRIES); 01106 #endif 01107 // Only do local rows 01108 for (unsigned row_index=mLo; row_index<mHi; row_index++) 01109 { 01110 // Put zero on the diagonal 01111 MatSetValue(mSinglePermutation, row_index, row_index, 0.0, INSERT_VALUES); 01112 01113 // Put one at (i,j) 01114 MatSetValue(mSinglePermutation, row_index, rPermutation[row_index], 1.0, INSERT_VALUES); 01115 01116 unsigned bi_index = 2*row_index; 01117 unsigned perm_index = 2*rPermutation[row_index]; 01118 01119 // Put zeroes on the diagonal 01120 MatSetValue(mDoublePermutation, bi_index, bi_index, 0.0, INSERT_VALUES); 01121 MatSetValue(mDoublePermutation, bi_index+1, bi_index+1, 0.0, INSERT_VALUES); 01122 01123 // Put ones at (i,j) 01124 MatSetValue(mDoublePermutation, bi_index, perm_index, 1.0, INSERT_VALUES); 01125 MatSetValue(mDoublePermutation, bi_index+1, perm_index+1, 1.0, INSERT_VALUES); 01126 } 01127 MatAssemblyBegin(mSinglePermutation, MAT_FINAL_ASSEMBLY); 01128 MatAssemblyBegin(mDoublePermutation, MAT_FINAL_ASSEMBLY); 01129 MatAssemblyEnd(mSinglePermutation, MAT_FINAL_ASSEMBLY); 01130 MatAssemblyEnd(mDoublePermutation, MAT_FINAL_ASSEMBLY); 01131 return true; 01132 } 01133 01134 void Hdf5DataWriter::SetFixedChunkSize(unsigned chunkSize) 01135 { 01136 assert(mIsInDefineMode); 01137 01138 mUseOptimalChunkSizeAlgorithm = false; 01139 mFixedChunkSize = chunkSize; 01140 }