Mesh Oriented datABase  (version 5.4.1)
Array-based unstructured mesh datastructure
WriteHDF5.hpp
Go to the documentation of this file.
00001 /**
00002  * MOAB, a Mesh-Oriented datABase, is a software component for creating,
00003  * storing and accessing finite element mesh data.
00004  *
00005  * Copyright 2004 Sandia Corporation.  Under the terms of Contract
00006  * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government
00007  * retains certain rights in this software.
00008  *
00009  * This library is free software; you can redistribute it and/or
00010  * modify it under the terms of the GNU Lesser General Public
00011  * License as published by the Free Software Foundation; either
00012  * version 2.1 of the License, or (at your option) any later version.
00013  *
00014  */
00015 
00016 #ifndef WRITE_HDF5_HPP
00017 #define WRITE_HDF5_HPP
00018 
00019 #include <list>
00020 #include "moab/MOABConfig.h"
00021 #ifdef MOAB_HAVE_MPI  // include this before HDF5 headers to avoid conflicts
00022 #include "moab_mpi.h"
00023 #endif
00024 #include "moab_mpe.h"
00025 #include "mhdf.h"
00026 #include "moab/Forward.hpp"
00027 #include "moab/Range.hpp"
00028 #include "moab/WriterIface.hpp"
00029 #include "moab/RangeMap.hpp"
00030 #include "moab/WriteUtilIface.hpp"
00031 #include "DebugOutput.hpp"
00032 #include "HDF5Common.hpp"
00033 
00034 namespace moab
00035 {
00036 
00037 class IODebugTrack;
00038 
00039 /* If this define is not set, node->entity adjacencies will not be written */
00040 #undef MB_H5M_WRITE_NODE_ADJACENCIES
00041 
00042 /**
00043  * \brief  Write mesh database to MOAB's native HDF5-based file format.
00044  * \author Jason Kraftcheck
00045  * \date   01 April 2004
00046  */
00047 class WriteHDF5 : public WriterIface
00048 {
00049 
00050   public:
00051     static WriterIface* factory( Interface* );
00052 
00053     WriteHDF5( Interface* iface );
00054 
00055     virtual ~WriteHDF5();
00056 
00057     /** Export specified meshsets to file
00058      * \param filename     The filename to export.
00059      * \param export_sets  Array of handles to sets to export, or NULL to export all.
00060      * \param export_set_count Length of <code>export_sets</code> array.
00061      */
00062     ErrorCode write_file( const char* filename,
00063                           const bool overwrite,
00064                           const FileOptions& opts,
00065                           const EntityHandle* export_sets,
00066                           const int export_set_count,
00067                           const std::vector< std::string >& qa_records,
00068                           const Tag* tag_list = NULL,
00069                           int num_tags        = 0,
00070                           int user_dimension  = 3 );
00071 
00072     /** The type to use for entity IDs w/in the file.
00073      *
00074      * NOTE:  If this is changed, the value of id_type
00075      *        MUST be changed accordingly.
00076      */
00077     typedef EntityHandle wid_t;  // change the name,
00078     //  to avoid conflicts to /usr/include/x86_64-linux-gnu/bits/types.h : id_t , which is unsigned
00079     //  int
00080 
00081     /** HDF5 type corresponding to type of wid_t */
00082     static const hid_t id_type;
00083 
00084     struct ExportType
00085     {
00086         //! The type of the entities in the range
00087         EntityType type;
00088         //! The number of nodes per entity - not used for nodes and sets
00089         int num_nodes;
00090 
00091         virtual ~ExportType() {}
00092 
00093         bool operator==( const ExportType& t ) const
00094         {
00095             return t.type == type && t.num_nodes == num_nodes;
00096         }
00097         bool operator!=( const ExportType& t ) const
00098         {
00099             return t.type != type || t.num_nodes != num_nodes;
00100         }
00101         bool operator<( const ExportType& t ) const
00102         {
00103             return type < t.type || ( type == t.type && num_nodes < t.num_nodes );
00104         }
00105     };
00106 
00107     //! Range of entities, grouped by type, to export
00108     struct ExportSet : public ExportType
00109     {
00110         //! The range of entities.
00111         Range range;
00112         //! The first Id allocated by the mhdf library.  Entities in range have sequential IDs.
00113         wid_t first_id;
00114         //! The offset at which to begin writing this processor's data.
00115         //! Always zero except for parallel IO.
00116         long offset;
00117         //! Offset for adjacency data.  Always zero except for parallel IO
00118         long adj_offset;
00119         //! If doing parallel IO, largest number of entities to write
00120         //! for any processor (needed to do collective IO).  Zero if unused.
00121         long max_num_ents, max_num_adjs;
00122         //! The total number of entities that will be written to the file
00123         //! for this group.  For serial IO, this should always be range.size().
00124         //! For parallel IO, it will be the sum of range size over all processors.
00125         //! For parallel IO, this value is undefined except for on the root
00126         //! processor.
00127         long total_num_ents;
00128 
00129         bool operator<( const ExportType& other ) const
00130         {
00131             return type < other.type || ( type == other.type && num_nodes < other.num_nodes );
00132         }
00133 
00134         bool operator<( std::pair< int, int > other ) const
00135         {
00136             return type < other.first || ( type == other.first && num_nodes < other.second );
00137         }
00138 
00139         bool operator==( const ExportType& other ) const
00140         {
00141             return ( type == other.type && num_nodes == other.num_nodes );
00142         }
00143 
00144         bool operator==( std::pair< int, int > other ) const
00145         {
00146             return ( type == other.first && num_nodes == other.second );
00147         }
00148 
00149         const char* name() const;
00150     };
00151 
00152     //! Tag to write to file.
00153     struct TagDesc
00154     {
00155         //! The tag handle
00156         Tag tag_id;
00157         //! The offset at which to begin writting this processor's data.
00158         //! Always zero except for parallel IO.
00159         wid_t sparse_offset;
00160         //! For variable-length tags, a second offset for the tag data table,
00161         //! separate from the offset used for the ID and Index tables.
00162         //! Always zero except for parallel IO.
00163         wid_t var_data_offset;
00164         //! Write sparse tag data (for serial, is always equal to !range.empty())
00165         bool write_sparse;
00166         //! If doing parallel IO, largest number, over all processes, of entities
00167         //! for which to write tag data.  Zero if unused.
00168         unsigned long max_num_ents;
00169         //! For variable-length tags during parallel IO: the largest number
00170         //! of tag values to be written on by any process, used to calculate
00171         //! the total number of collective writes that all processes must do.
00172         //! Zero for fixed-length tags or if not doing parallel IO.
00173         unsigned long max_num_vals;
00174 
00175         //! List of entity groups for which to write tag data in
00176         //! dense format
00177         std::vector< ExportType > dense_list;
00178 
00179         bool have_dense( const ExportType& type ) const
00180         {
00181             return std::find( dense_list.begin(), dense_list.end(), type ) != dense_list.end();
00182         }
00183 
00184         bool operator<( const TagDesc& ) const;
00185     };
00186 
00187     /** Create attributes holding the HDF5 type handle for the
00188      *  type of a bunch of the default tags.
00189      */
00190     // static ErrorCode register_known_tag_types( Interface* );
00191 
00192     //! Store old HDF5 error handling function
00193     struct HDF5ErrorHandler
00194     {
00195         HDF5_Error_Func_Type func;
00196         void* data;
00197     };
00198 
00199     mhdf_FileHandle file_ptr()
00200     {
00201         return filePtr;
00202     }
00203 
00204     WriteUtilIface* write_util()
00205     {
00206         return writeUtil;
00207     }
00208 
00209   protected:
00210     //! Store old HDF5 error handling function
00211     HDF5ErrorHandler errorHandler;
00212 
00213     /** Function to create the file.  Virtual to allow override
00214      *  for parallel version.
00215      */
00216     virtual ErrorCode parallel_create_file( const char* filename,
00217                                             bool overwrite,
00218                                             const std::vector< std::string >& qa_records,
00219                                             const FileOptions& opts,
00220                                             const Tag* tag_list,
00221                                             int num_tags,
00222                                             int dimension = 3,
00223                                             double* times = 0 );
00224     virtual ErrorCode write_finished();
00225     virtual void debug_barrier_line( int lineno );
00226 
00227     //! Gather tags
00228     ErrorCode gather_tags( const Tag* user_tag_list, int user_tag_list_length );
00229 
00230     /** Check if tag values for a given ExportSet should be written in dense format
00231      *
00232      *\param ents        ExportSet to consider
00233      *\param all_tagged  Range containing all the entities in ents.range for
00234      *                   which an explicit tag value is stored.  Range may
00235      *                   also contain entities not in ents.range, but may
00236      *                   not contain entities in ents.range for which no tag
00237      *                   value is stored.
00238      *\param prefer_dense If true, will return true if at least 2/3 of the
00239      *                   entities are tagged.  This should not be passed as
00240      *                   true if the tag does not have a default value, as
00241      *                   tag values must be stored for all entities in the
00242      *                   ExportSet for dense-formatted data.
00243      */
00244     bool check_dense_format_tag( const ExportSet& ents, const Range& all_tagged, bool prefer_dense );
00245 
00246     /** Helper function for create-file
00247      *
00248      * Calculate the sum of the number of non-set adjacencies
00249      * of all entities in the passed range.
00250      */
00251     ErrorCode count_adjacencies( const Range& elements, wid_t& result );
00252 
00253   public:  // make these public so helper classes in WriteHDF5Parallel can use them
00254     /** Helper function for create-file
00255      *
00256      * Create zero-ed tables where element connectivity and
00257      * adjacency data will be stored.
00258      */
00259     ErrorCode create_elem_table( const ExportSet& block, long num_ents, long& first_id_out );
00260 
00261     /** Helper function for create-file
00262      *
00263      * Create zero-ed table where set descriptions will be written
00264      */
00265     ErrorCode create_set_meta( long num_sets, long& first_id_out );
00266 
00267   protected:
00268     /** Helper function for create-file
00269      *
00270      * Calculate total length of set contents and child tables.
00271      */
00272     ErrorCode count_set_size( const Range& sets,
00273                               long& contents_length_out,
00274                               long& children_length_out,
00275                               long& parents_length_out );
00276 
00277     //! Get information about a meshset
00278     ErrorCode get_set_info( EntityHandle set,
00279                             long& num_entities,
00280                             long& num_children,
00281                             long& num_parents,
00282                             unsigned long& flags );
00283 
00284     /** Helper function for create-file
00285      *
00286      * Create zero-ed tables where set data will be written.
00287      */
00288     ErrorCode create_set_tables( long contents_length, long children_length, long parents_length );
00289 
00290     //! Write exodus-type QA info
00291     ErrorCode write_qa( const std::vector< std::string >& list );
00292 
00293     //!\brief Get tagged entities for which to write tag values
00294     ErrorCode get_num_sparse_tagged_entities( const TagDesc& tag, size_t& count );
00295     //!\brief Get tagged entities for which to write tag values
00296     ErrorCode get_sparse_tagged_entities( const TagDesc& tag, Range& range );
00297     //!\brief Get entities that will be written to file
00298     void get_write_entities( Range& range );
00299 
00300     //! The size of the data buffer (<code>dataBuffer</code>).
00301     size_t bufferSize;
00302     //! A memory buffer to use for all I/O operations.
00303     char* dataBuffer;
00304 
00305     //! Interface pointer passed to constructor
00306     Interface* iFace;
00307     //! Cached pointer to writeUtil interface.
00308     WriteUtilIface* writeUtil;
00309 
00310     //! The file handle from the mhdf library
00311     mhdf_FileHandle filePtr;
00312 
00313     //! Map from entity handles to file IDs
00314     RangeMap< EntityHandle, wid_t > idMap;
00315 
00316     //! The list elements to export.
00317     std::list< ExportSet > exportList;
00318     //! The list of nodes to export
00319     ExportSet nodeSet;
00320     //! The list of sets to export
00321     ExportSet setSet;
00322 
00323     const ExportSet* find( const ExportType& type ) const
00324     {
00325         if( type.type == MBVERTEX )
00326             return &nodeSet;
00327         else if( type.type == MBENTITYSET )
00328             return &setSet;
00329         else
00330         {
00331             std::list< ExportSet >::const_iterator it;
00332             it = std::find( exportList.begin(), exportList.end(), type );
00333             return it == exportList.end() ? 0 : &*it;
00334         }
00335     }
00336 
00337     //! Offset into set contents table (zero except for parallel)
00338     unsigned long setContentsOffset;
00339     //! Offset into set children table (zero except for parallel)
00340     unsigned long setChildrenOffset, setParentsOffset;
00341     //! The largest number of values to write
00342     //! for any processor (needed to do collective IO).
00343     long maxNumSetContents, maxNumSetChildren, maxNumSetParents;
00344     //! Flags idicating if set data should be written.
00345     //! For the normal (non-parallel) case, these values
00346     //! will depend only on whether or not there is any
00347     //! data to be written.  For parallel-meshes, opening
00348     //! the data table is collective so the values must
00349     //! depend on whether or not any processor has meshsets
00350     //! to be written.
00351     bool writeSets, writeSetContents, writeSetChildren, writeSetParents;
00352 
00353     //! Struct describing a set for which the contained and linked entity
00354     //! lists are something other than the local values.  Used to store
00355     //! data for shared sets owned by this process when writing in parallel.
00356     struct SpecialSetData
00357     {
00358         EntityHandle setHandle;
00359         unsigned setFlags;
00360         std::vector< wid_t > contentIds;
00361         std::vector< wid_t > childIds;
00362         std::vector< wid_t > parentIds;
00363     };
00364     struct SpecSetLess
00365     {
00366         bool operator()( const SpecialSetData& a, SpecialSetData b ) const
00367         {
00368             return a.setHandle < b.setHandle;
00369         }
00370     };
00371 
00372     //! Array of special/shared sets, in order of handle value.
00373     std::vector< SpecialSetData > specialSets;
00374     const SpecialSetData* find_set_data( EntityHandle h ) const
00375     {
00376         return const_cast< WriteHDF5* >( this )->find_set_data( h );
00377     }
00378     SpecialSetData* find_set_data( EntityHandle h );
00379 
00380     //! The list of tags to export
00381     std::list< TagDesc > tagList;
00382 
00383     //! True if doing parallel write
00384     bool parallelWrite;
00385     //! True if using collective IO calls for parallel write
00386     bool collectiveIO;
00387     //! True if writing dense-formatted tag data
00388     bool writeTagDense;
00389 
00390     //! Property set to pass to H5Dwrite calls.
00391     //! For serial, should be H5P_DEFAULTS.
00392     //! For parallel, may request collective IO.
00393     hid_t writeProp;
00394 
00395     //! Utility to log debug output
00396     DebugOutput dbgOut;
00397 
00398     static MPEState topState;
00399     static MPEState subState;
00400 
00401     //! Look for overlapping and/or missing writes
00402     bool debugTrack;
00403 
00404     void print_id_map() const;
00405     void print_id_map( std::ostream& str, const char* prefix = "" ) const;
00406 
00407     /** Helper function for create-file
00408      *
00409      * Write tag meta-info and create zero-ed table where
00410      * tag values will be written.
00411      *\param num_entities  Number of entities for which to write tag data.
00412      *\param var_len_total For variable-length tags, the total number of values
00413      *                     in the data table.
00414      */
00415     ErrorCode create_tag( const TagDesc& tag_data, unsigned long num_entities, unsigned long var_len_total );
00416 
00417     /**\brief add entities to idMap */
00418     ErrorCode assign_ids( const Range& entities, wid_t first_id );
00419 
00420     /** Get possibly compacted list of IDs for passed entities
00421      *
00422      * For the passed range of entities, determine if IDs
00423      * can be compacted and write IDs to passed list.
00424      *
00425      * If the IDs are not compacted, the output list will contain
00426      * a simple ordered list of IDs.
00427      *
00428      * If IDs are compacted, the output list will contain
00429      * {start,count} pairs.
00430      *
00431      * If the ID list is compacted, ranged_list will be 'true'.
00432      * Otherwise it will be 'false'.
00433      */
00434     ErrorCode range_to_blocked_list( const Range& input_range,
00435                                      std::vector< wid_t >& output_id_list,
00436                                      bool& ranged_list );
00437 
00438     /** Get possibly compacted list of IDs for passed entities
00439      *
00440      * For the passed range of entities, determine if IDs
00441      * can be compacted and write IDs to passed list.
00442      *
00443      * If the IDs are not compacted, the output list will contain
00444      * a simple ordered list of IDs.
00445      *
00446      * If IDs are compacted, the output list will contain
00447      * {start,count} pairs.
00448      *
00449      * If the ID list is compacted, ranged_list will be 'true'.
00450      * Otherwise it will be 'false'.
00451      */
00452     ErrorCode range_to_blocked_list( const EntityHandle* input_ranges,
00453                                      size_t num_input_ranges,
00454                                      std::vector< wid_t >& output_id_list,
00455                                      bool& ranged_list );
00456 
00457     ErrorCode range_to_id_list( const Range& input_range, wid_t* array );
00458     //! Get IDs for entities
00459     ErrorCode vector_to_id_list( const std::vector< EntityHandle >& input,
00460                                  std::vector< wid_t >& output,
00461                                  bool remove_non_written = false );
00462     //! Get IDs for entities
00463     ErrorCode vector_to_id_list( const EntityHandle* input, wid_t* output, size_t num_entities );
00464     //! Get IDs for entities
00465     ErrorCode vector_to_id_list( const EntityHandle* input,
00466                                  size_t input_len,
00467                                  wid_t* output,
00468                                  size_t& output_len,
00469                                  bool remove_non_written );
00470 
00471     /** When writing tags containing EntityHandles to file, need to convert tag
00472      *  data from EntityHandles to file IDs.  This function does that.
00473      *
00474      * If the handle is not valid or does not correspond to an entity that will
00475      * be written to the file, the file ID is set to zero.
00476      *\param data  The data buffer.  As input, an array of EntityHandles.  As
00477      *             output an array of file IDS, where the size of each integral
00478      *             file ID is the same as the size of EntityHandle.
00479      *\param count The number of handles in the buffer.
00480      *\return true if at least one of the handles is valid and will be written to
00481      *             the file or at least one of the handles is NULL (zero). false
00482      *             otherwise
00483      */
00484     bool convert_handle_tag( EntityHandle* data, size_t count ) const;
00485     bool convert_handle_tag( const EntityHandle* source, EntityHandle* dest, size_t count ) const;
00486 
00487     /** Get IDs of adjacent entities.
00488      *
00489      * For all entities adjacent to the passed entity, if the
00490      * adjacent entity is to be exported (ID is not zero), append
00491      * the ID to the passed list.
00492      */
00493     ErrorCode get_adjacencies( EntityHandle entity, std::vector< wid_t >& adj );
00494 
00495     //! get sum of lengths of tag values (as number of type) for
00496     //! variable length tag data.
00497     ErrorCode get_tag_data_length( const TagDesc& tag_info, const Range& range, unsigned long& result );
00498 
00499   private:
00500     //! Do the actual work of write_file.  Separated from write_file
00501     //! for easier resource cleanup.
00502     ErrorCode write_file_impl( const char* filename,
00503                                const bool overwrite,
00504                                const FileOptions& opts,
00505                                const EntityHandle* export_sets,
00506                                const int export_set_count,
00507                                const std::vector< std::string >& qa_records,
00508                                const Tag* tag_list,
00509                                int num_tags,
00510                                int user_dimension = 3 );
00511 
00512     ErrorCode init();
00513 
00514     ErrorCode serial_create_file( const char* filename,
00515                                   bool overwrite,
00516                                   const std::vector< std::string >& qa_records,
00517                                   const Tag* tag_list,
00518                                   int num_tags,
00519                                   int dimension = 3 );
00520 
00521     /** Get all mesh to export from given list of sets.
00522      *
00523      * Populate exportSets, nodeSet and setSet with lists of
00524      * entities to write.
00525      *
00526      * \param export_sets  The list of meshsets to export
00527      */
00528     ErrorCode gather_mesh_info( const std::vector< EntityHandle >& export_sets );
00529 
00530     //! Same as gather_mesh_info, except for entire mesh
00531     ErrorCode gather_all_mesh();
00532 
00533     //! Initialize internal data structures from gathered mesh
00534     ErrorCode initialize_mesh( const Range entities_by_dim[5] );
00535 
00536     /** Write out the nodes.
00537      *
00538      * Note: Assigns IDs to nodes.
00539      */
00540     ErrorCode write_nodes();
00541 
00542     /** Write out element connectivity.
00543      *
00544      * Write connectivity for passed set of elements.
00545      *
00546      * Note: Assigns element IDs.
00547      * Note: Must do write_nodes first so node IDs get assigned.
00548      */
00549     ErrorCode write_elems( ExportSet& elemset );
00550 
00551     /** Write out meshsets
00552      *
00553      * Write passed set of meshsets, including parent/child relations.
00554      *
00555      * Note: Must have written nodes and element connectivity
00556      *       so entities have assigned IDs.
00557      */
00558     ErrorCode write_sets( double* times );
00559 
00560     /** Write set contents/parents/children lists
00561      *
00562      *\param which_data Which set data to write (contents, parents, or children)
00563      *\param handle     HDF5 handle for data set in which to write data
00564      *\param track      Debugging tool
00565      *\param ranged     Will be populated with handles of sets for which
00566      *                  contents were written in a range-compacted format.
00567      *                  (mhdf_SET_RANGE_BIT).  Should be null for parents/children.
00568      *\param null_stripped Will be populated with handles of sets for which
00569      *                  invalid or null handles were stripped from the contents
00570      *                  list.  This is only done for unordered sets.  This argument
00571      *                  should be null if writing parents/children because those
00572      *                  lists are always ordered.
00573      *\param set_sizes  Will be populated with the length of the data written
00574      *                  for those sets for which the handles were added to
00575      *                  either \c ranged or \c null_stripped.  Values are
00576      *                  in handle order.
00577      */
00578     ErrorCode write_set_data( const WriteUtilIface::EntityListType which_data,
00579                               const hid_t handle,
00580                               IODebugTrack& track,
00581                               Range* ranged                  = 0,
00582                               Range* null_stripped           = 0,
00583                               std::vector< long >* set_sizes = 0 );
00584 
00585     /** Write adjacency info for passed set of elements
00586      *
00587      * Note: Must have written element connectivity so elements
00588      *       have IDs assigned.
00589      */
00590     ErrorCode write_adjacencies( const ExportSet& export_set );
00591 
00592     /** Write tag information and data.
00593      *
00594      * Note: Must have already written nodes, elem connectivity and
00595      *       sets so that entities have IDs assigned.
00596      */
00597 
00598     //! Write tag for all entities.
00599     ErrorCode write_tag( const TagDesc& tag_data, double* times );
00600 
00601     //! Get element connectivity
00602     ErrorCode get_connectivity( Range::const_iterator begin,
00603                                 Range::const_iterator end,
00604                                 int nodes_per_element,
00605                                 wid_t* id_data_out );
00606 
00607     //! Get size data for tag
00608     //!\param tag       MOAB tag ID
00609     //!\param moab_type Output: DataType for tag
00610     //!\param num_bytes Output: MOAB tag size (bits for bit tags).
00611     //!                         MB_VARIABLE_LENGTH for variable-length tags.
00612     //!\param elem_size Output: Size of of the base data type of the
00613     //!                         tag data (e.g. sizeof(double) if
00614     //!                         moab_type == MB_TYPE_DOUBLE).
00615     //!                         One for bit and opaque tags.
00616     //!\param array_size Output: The number of valeus of size elem_size
00617     //!                          for each tag.  Always 1 for opaque data.
00618     //!                          Nubmer of bits for bit tags.
00619     //!\param file_type Output: mhdf type enumeration
00620     //!\param hdf_type  Output: Handle to HDF5 type object.  Caller is
00621     //!                         responsible for releasing this object
00622     //!                         (calling H5Tclose).
00623     ErrorCode get_tag_size( Tag tag,
00624                             DataType& moab_type,
00625                             int& num_bytes,
00626                             int& elem_size,
00627                             int& file_size,
00628                             mhdf_TagDataType& file_type,
00629                             hid_t& hdf_type );
00630 
00631     //! Write ID table for sparse tag
00632     ErrorCode write_sparse_ids( const TagDesc& tag_data,
00633                                 const Range& range,
00634                                 hid_t table_handle,
00635                                 size_t table_size,
00636                                 const char* name = 0 );
00637 
00638     //! Write fixed-length tag data in sparse format
00639     ErrorCode write_sparse_tag( const TagDesc& tag_data,
00640                                 const std::string& tag_name,
00641                                 DataType tag_data_type,
00642                                 hid_t hdf5_data_type,
00643                                 int hdf5_type_size );
00644 
00645     //! Write end index data_set for a variable-length tag
00646     ErrorCode write_var_len_indices( const TagDesc& tag_data,
00647                                      const Range& range,
00648                                      hid_t idx_table,
00649                                      size_t table_size,
00650                                      int type_size,
00651                                      const char* name = 0 );
00652 
00653     //! Write tag value data_set for a variable-length tag
00654     ErrorCode write_var_len_data( const TagDesc& tag_data,
00655                                   const Range& range,
00656                                   hid_t table,
00657                                   size_t table_size,
00658                                   bool handle_tag,
00659                                   hid_t hdf_type,
00660                                   int type_size,
00661                                   const char* name = 0 );
00662 
00663     //! Write varialbe-length tag data
00664     ErrorCode write_var_len_tag( const TagDesc& tag_info,
00665                                  const std::string& tag_name,
00666                                  DataType tag_data_type,
00667                                  hid_t hdf5_type,
00668                                  int hdf5_type_size );
00669 
00670     //! Write dense-formatted tag data
00671     ErrorCode write_dense_tag( const TagDesc& tag_data,
00672                                const ExportSet& elem_data,
00673                                const std::string& tag_name,
00674                                DataType tag_data_type,
00675                                hid_t hdf5_data_type,
00676                                int hdf5_type_size );
00677 
00678     //! Write data for fixed-size tag
00679     ErrorCode write_tag_values( Tag tag_id,
00680                                 hid_t data_table,
00681                                 unsigned long data_offset,
00682                                 const Range& range,
00683                                 DataType tag_data_type,
00684                                 hid_t hdf5_data_type,
00685                                 int hdf5_type_size,
00686                                 unsigned long max_num_ents,
00687                                 IODebugTrack& debug_track );
00688 
00689   protected:
00690     enum TimingValues
00691     {
00692         TOTAL_TIME = 0,
00693         GATHER_TIME,
00694         CREATE_TIME,
00695         CREATE_NODE_TIME,
00696         NEGOTIATE_TYPES_TIME,
00697         CREATE_ELEM_TIME,
00698         FILEID_EXCHANGE_TIME,
00699         CREATE_ADJ_TIME,
00700         CREATE_SET_TIME,
00701         SHARED_SET_IDS,
00702         SHARED_SET_CONTENTS,
00703         SET_OFFSET_TIME,
00704         CREATE_TAG_TIME,
00705         COORD_TIME,
00706         CONN_TIME,
00707         SET_TIME,
00708         SET_META,
00709         SET_CONTENT,
00710         SET_PARENT,
00711         SET_CHILD,
00712         ADJ_TIME,
00713         TAG_TIME,
00714         DENSE_TAG_TIME,
00715         SPARSE_TAG_TIME,
00716         VARLEN_TAG_TIME,
00717         NUM_TIMES
00718     };
00719 
00720     virtual void print_times( const double times[NUM_TIMES] ) const;
00721 };
00722 
00723 }  // namespace moab
00724 
00725 #endif
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines