Mesh Oriented datABase  (version 5.4.1)
Array-based unstructured mesh datastructure
WriteUtilIface.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 MOAB_WRITE_UTIL_IFACE_HPP
00017 #define MOAB_WRITE_UTIL_IFACE_HPP
00018 
00019 #include "moab/Range.hpp"
00020 #include "moab/Compiler.hpp"
00021 #include <vector>
00022 
00023 namespace moab
00024 {
00025 
00026 //! Interface implemented in MOAB which provides memory for mesh reading utilities
00027 class WriteUtilIface
00028 {
00029   public:
00030     //! Constructor
00031     WriteUtilIface() {}
00032 
00033     //! Destructor
00034     virtual ~WriteUtilIface() {}
00035 
00036     //! Check if the specified file already exists.
00037     //! Returns MB_SUCCESS if file does not exist, MB_ALREADY_ALLOCATED
00038     //! if file does exist, or MB_FAILURE for some other error condition.
00039     virtual ErrorCode check_doesnt_exist( const char* file_name ) = 0;
00040 
00041     //! Gather all entities in the mesh, or in the sets specified
00042     virtual ErrorCode gather_entities(
00043         Range& all_ents,                     /**< range in which entities are returned */
00044         const EntityHandle* ent_sets = NULL, /**< entity sets whose contents are to be gathered */
00045         int num_sets                 = 0 /**< number of sets in list */ ) = 0;
00046 
00047     //! Given information about the nodes to be written, and pointers to memory
00048     //! to which coordinates will be written, writes coordinate data there, and
00049     //! also assigns global ids to nodes & writes to a tag, if a tag is specified
00050     //! \param num_arrays Number of coordinate arrays requested
00051     //! \param num_nodes Number of nodes to be written
00052     //! \param entities Range of nodes to be written
00053     //! \param node_id_tag Tag used to write ids to nodes
00054     //! \param start_node_id Starting value for node ids
00055     //! \param arrays Pointers to memory where coordinate data will be written
00056     //! \return status Return status
00057     virtual ErrorCode get_node_coords( const int num_arrays,
00058                                        const int num_nodes,
00059                                        const Range& entities,
00060                                        Tag node_id_tag,
00061                                        const int start_node_id,
00062                                        std::vector< double* >& arrays ) = 0;
00063 
00064     /** Get array of coordinate values for nodes
00065      *
00066      * Given a range of node handles, retrieve a single or multiple coordinate
00067      * value(s) for each.
00068      *
00069      * Failure conditions:
00070      *  - invalid entity handles (not vertices, non-existent entity, etc.)
00071      *  - range is empty (<code>iter == end</code>)
00072      *  - <code>output_array</code> is null
00073      *  - insufficient space in <code>output_array</code>
00074      *
00075      *\param which_array  The coordinate to retrieve (0-&gt;X, 1-&gt;Y, 2-&gt;Z, -1-&gt;all)
00076      *\param begin        The first node handle.
00077      *\param end          One past the last node handle.
00078      *\param output_size  The size of <code>output_array</code>.
00079      *\param output_array The memory in which to write the node coordinates.
00080      *\author Jason Kraftcheck
00081      */
00082     virtual ErrorCode get_node_coords( const int which_array,
00083                                        Range::const_iterator begin,
00084                                        const Range::const_iterator& end,
00085                                        const size_t output_size,
00086                                        double* const output_array ) = 0;
00087 
00088     //! Given information about elements to be written and a pointer to memory
00089     //! where connectivity for those elements should be written, writes connectivity
00090     //! to that memory; uses node ids stored in a tag during call to <em>get_node_coords</em>
00091     //! function
00092     //! \param num_elements Number of elements to be written
00093     //! \param verts_per_element Number of vertices per element
00094     //! \param node_id_tag Tag used to store node ids
00095     //! \param entities Range of elements to be written
00096     //! \param element_id_tag Tag which should be used to store element ids
00097     //! \param start_element_id Starting value for element ids
00098     //! \param array Pointer to memory where connectivity data will be written
00099     //! \return status Return status
00100     virtual ErrorCode get_element_connect( const int num_elements,
00101                                            const int verts_per_element,
00102                                            Tag node_id_tag,
00103                                            const Range& entities,
00104                                            Tag element_id_tag,
00105                                            int start_element_id,
00106                                            int* array,
00107                                            bool add_sizes = false ) = 0;
00108 
00109     /** Get connectivity for elements
00110      *
00111      * Get the connectivity list for a range of elements.
00112      *
00113      * Failure cases:
00114      *  - Passed range is empty (<code>begin == end</code>).
00115      *  - <code>vertices_per_elem</code> is less than one
00116      *  - <code>element_array</code> is null.
00117      *  - The range contains invalid handles (non-existent entities,
00118      *      not an element, etc.)
00119      *  - Retrieving ID tag for an entity failed.
00120      *  - Insufficient space in passed array.
00121      *
00122      *\param begin        The first element handle
00123      *\param end          One past the last element handle
00124      *\param vertices_per_elem Number of vertices to retrieve for each
00125      *                    element. If the element has more vertices, the
00126      *                    element connectivity will be truncated. If
00127      *                    <code>vertices_per_elem</code> is greater than the
00128      *                    number of nodes for an element, the data will be
00129      *                    padded with zeros.
00130      *\param node_id_tag  A tag with integer values.
00131      *\param array_size   The length of <code>element_array</code>
00132      *\param element_array The memory location at which to store the
00133      *                    connectivity list.
00134      *\param add_sizes If true, writes size of connect array before connectivity in array
00135      *\author Jason Kraftcheck
00136      */
00137     virtual ErrorCode get_element_connect( Range::const_iterator begin,
00138                                            const Range::const_iterator& end,
00139                                            const int vertices_per_elem,
00140                                            Tag node_id_tag,
00141                                            const size_t array_size,
00142                                            int* const element_array,
00143                                            bool add_sizes = false ) = 0;
00144 
00145     /** Get connectivity for elements
00146      *
00147      * Get the connectivity list for a range of elements.
00148      *
00149      * Failure cases:
00150      *  - Passed range is empty (<code>begin == end</code>).
00151      *  - <code>vertices_per_elem</code> is less than one
00152      *  - <code>element_array</code> is null.
00153      *  - The range contains invalid handles (non-existent entities,
00154      *      not an element, etc.)
00155      *  - Insufficient space in passed array.
00156      *
00157      *\param begin        The first element handle
00158      *\param end          One past the last element handle
00159      *\param vertices_per_elem Number of vertices to retrieve for each
00160      *                    element. If the element has more vertices, the
00161      *                    element connectivity will be truncated. If
00162      *                    <code>vertices_per_elem</code> is greater than the
00163      *                    number of nodes for an element, the data will be
00164      *                    padded with zeros.
00165      *\param array_size   The length of <code>element_array</code>
00166      *\param element_array The memory location at which to store the
00167      *                    connectivity list.
00168      *\author Jason Kraftcheck
00169      */
00170     virtual ErrorCode get_element_connect( Range::const_iterator begin,
00171                                            const Range::const_iterator& end,
00172                                            const int vertices_per_elem,
00173                                            const size_t array_size,
00174                                            EntityHandle* const element_array ) = 0;
00175 
00176     /** Get poly (polygon or polyhedron) connectivity size
00177      *\param begin  First iterator in range of poly
00178      *\param end    One past last in range of poly.
00179      *\param connectivity_size  The length of the connectivity list
00180      *              For the specified range of polyhedra.
00181      *\author Jason Kraftcheck
00182      */
00183     virtual ErrorCode get_poly_connect_size( Range::const_iterator begin,
00184                                              const Range::const_iterator& end,
00185                                              int& connectivity_size ) = 0;
00186 
00187     /** Get poly (polygon or polyhedron) connectivity.
00188      *
00189      * Connectivity is returned in two arrays. The first is
00190      * an array of global IDs that is the concatenation of the
00191      * connectivity for the entire range of polys. The second
00192      * is the last index of the connectivity data for each poly
00193      * in the global ID array.
00194      *
00195      * This function will add as many polys as possible to the
00196      * passed arrays given the sizes of those arrays. It will
00197      * then pass back position at which it stopped and the sizes
00198      * of the data written to the arrays.
00199      *
00200      * Failure cases:
00201      *  - Passed range is empty (<code>begin == end</code>).
00202      *  - <code>element_array</code> or <code>index_array</code> is null.
00203      *  - The range contains invalid handles (non-existent entities,
00204      *      not an poly, etc.)
00205      *  - Retrieving ID tag for an entity failed.
00206      *
00207      *\param iter               As input, the first element handle.
00208      *                          As output, one past the last element handle
00209      *                          for which data was written to the arrays.
00210      *\param end                The iterator at which to stop.
00211      *\param node_id_tag        A tag with integer values.
00212      *\param element_array_len  As input, length of <code>element_array</code>.
00213      *                          As output, the number of entries written in that
00214      *                          array.
00215      *\param element_array      The memory location at which to store the
00216      *                          connectivity list.
00217      *\param index_array_len    As input, the length of <code>index_array</code>.
00218      *                          As output, the number of entries written in that
00219      *                          array.
00220      *\param index_array        The memory location at which to store offsets.
00221      *\param index_offset       Value to offset (add to) index values. As output
00222      *                          the input value plus the amount of data
00223      *                          written to the element array. (The value you
00224      *                          presumably want to pass to the next call.)
00225      *\author Jason Kraftcheck
00226      */
00227     virtual ErrorCode get_poly_connect( Range::const_iterator& iter,
00228                                         const Range::const_iterator& end,
00229                                         const Tag node_id_tag,
00230                                         size_t& element_array_len,
00231                                         int* const element_array,
00232                                         size_t& index_array_len,
00233                                         int* const index_array,
00234                                         int& index_offset ) = 0;
00235 
00236     //! Given elements to be written, gather all the nodes which define those elements
00237     //! \param elements Range of elements to be written
00238     //! \param node_bit_mark_tag Bit tag to use to identify nodes
00239     //! \param nodes Range of nodes gathered from elements (returned)
00240     //! \return status Return status
00241     virtual ErrorCode gather_nodes_from_elements( const Range& elements,
00242                                                   const Tag node_bit_mark_tag,
00243                                                   Range& nodes ) = 0;
00244 
00245     //! assign ids to input entities starting with start_id, written to id_tag
00246     //! if id_tag is zero, assigns to GLOBAL_ID_TAG_NAME
00247     //! \param elements Entities to be written
00248     //! \param id_tag Tag used to store entity id
00249     //! \param start_id Starting value for entity ids
00250     //! \return status Return status
00251     virtual ErrorCode assign_ids( Range& elements, Tag id_tag, const int start_id ) = 0;
00252 
00253     /** Get explicit adjacencies
00254      *
00255      * Get explicit adjacences stored in database.
00256      * Does not create any explicit adjacencies or search for
00257      * implicit ones.
00258      *
00259      *\param entity  The entity to retrieve adjacencies for.
00260      *\param id_tag  The global ID tag
00261      *\param adj     The output list of global IDs of adjacent entities.
00262      *\author Jason Kraftcheck
00263      */
00264     virtual ErrorCode get_adjacencies( EntityHandle entity, Tag id_tag, std::vector< int >& adj ) = 0;
00265 
00266     virtual ErrorCode get_adjacencies( EntityHandle entity, const EntityHandle*& adj_array, int& num_adj ) = 0;
00267 
00268     /**\brief Re-order outgoing element connectivity
00269      *
00270      * Permute the connectivity of each element such that the node
00271      * order is that of the target file format rather than that of MBCN.
00272      *\param order The permutation to use. Must be an array of 'node_per_elem'
00273      *             integers and be a permutation of the values [0..node_per_elem-1].
00274      *             Such that for a single element:
00275      *             target_conn[i] == mbcn_conn[order[i]]
00276      *\param conn  The connectivity array to re-order
00277      *\param num_elem  The number of elements in the connectivity array
00278      *\param node_per_elem The number of nodes in each element's connectivity list.
00279      */
00280     template < typename T >
00281     static inline void reorder( const int* order, T* conn, int num_elem, int node_per_elem );
00282 
00283     /**\brief Get list of tags to write.
00284      *
00285      * Get the list of tags to write to the file, possibly using
00286      * an optional user-specified tag list. This function consolidates
00287      * some common code for file writers to use to figure out what
00288      * tag data to write to the file. It provides the following features:
00289      *  o filter list based on user-specified array of tag handles
00290      *  o filter internal tags (those for which the name is prefixed with
00291      *    two underscore characters)
00292      *  o filter anonymous tags
00293      *  o optionally filter variable-length tags.
00294      *
00295      *\author Jason Kraftcheck
00296      *\param result_list List of tag handles for which to write data
00297      *\param user_tag_list Optional array of tag handles passed by user
00298      *                     to write to file.
00299      *\param include_variable_length_tags If false, return only fixed-length
00300      *                                    tags.
00301      */
00302     virtual ErrorCode get_tag_list( std::vector< Tag >& result_list,
00303                                     const Tag* user_tag_list          = 0,
00304                                     int user_tag_list_length          = 0,
00305                                     bool include_variable_length_tags = true ) = 0;
00306 
00307     enum EntityListType
00308     {
00309         CONTENTS    = 0,
00310         CHILDREN    = 1,
00311         PARENTS     = 2,
00312         TOPOLOGICAL = 1
00313     };
00314 
00315     /*\brief Get pointers to internal storage of entity data
00316      *
00317      * Get pointers to element connectivity or set content storage.
00318      *\param query_begin Start of range of entities for which to return results
00319      *\param query_end   End of range of entities for which to return results.
00320      *\param output_pointer_array Result list of pointers. Points to either
00321      *          element connectivity or set contents. Note: set contents
00322      *          may be in range-compacted format.
00323      *\param lengths Optional per-entity length of list. If passed, then
00324      *          always set, for each entity, to the number of values in the
00325      *          array passed back in \c output_pointer_array
00326      *\param relation If entity is entity set, which set data to return
00327      *          (contents array, parent array, or child array). If
00328      *          entity is an element, then CONTENTS for complete connectivity
00329      *          or TOPOLOGICAL for only corner vertices.
00330      *\param flags Optional per-entity flag values. If passed, then
00331      *          always set to zero for elements and set to set creation
00332      *          flags for entity sets.
00333      *\return MB_STRUCTURED_MESH if one or more input elements are stored as
00334      *          structured mesh and therefore do not have explicit storage.
00335      *        MB_TYPE_OUT_OF_RANGE if called for vertices.
00336      */
00337     virtual ErrorCode get_entity_list_pointers( Range::const_iterator query_begin,
00338                                                 Range::const_iterator query_end,
00339                                                 EntityHandle const** output_pointer_array,
00340                                                 EntityListType relation = CONTENTS,
00341                                                 int* lengths            = 0,
00342                                                 unsigned char* flags    = 0 ) = 0;
00343 
00344     /*\brief Get pointers to internal storage of entity data
00345      *
00346      * Get pointers to element connectivity or set content storage.
00347      *\param entities Pointer to list of entities for which to return results
00348      *\param num_entities Number of entities in list
00349      *\param output_pointer_array Result list of pointers. Points to either
00350      *          element connectivity or set contents. Note: set contents
00351      *          may be in range-compacted format.
00352      *\param lengths Optional per-entity length of list. If passed, then
00353      *          always set, for each entity, to the number of values in the
00354      *          array passed back in \c output_pointer_array
00355      *\param relation If entity is entity set, which set data to return
00356      *          (contents array, parent array, or child array). If
00357      *          entity is an element, then CONTENTS for complete connectivity
00358      *          or TOPOLOGICAL for only corner vertices.
00359      *\param flags Optional per-entity flag values. If passed, then
00360      *          always set to zero for elements and set to set creation
00361      *          flags for entity sets.
00362      *\return MB_STRUCTURED_MESH if one or more input elements are stored as
00363      *          structured mesh and therefore do not have explicit storage.
00364      *        MB_TYPE_OUT_OF_RANGE if called for vertices.
00365      */
00366     virtual ErrorCode get_entity_list_pointers( EntityHandle const* entities,
00367                                                 int num_entities,
00368                                                 EntityHandle const** output_pointer_array,
00369                                                 EntityListType relation = CONTENTS,
00370                                                 int* lengths            = 0,
00371                                                 unsigned char* flags    = 0 ) = 0;
00372 };
00373 
00374 template < typename T >
00375 inline void WriteUtilIface::reorder( const int* order, T* conn, int num_elem, int node_per_elem )
00376 {
00377     std::vector< T > elem( node_per_elem );
00378     T* const end = conn + num_elem * node_per_elem;
00379     while( conn != end )
00380     {
00381         std::copy( conn, conn + node_per_elem, elem.begin() );
00382         for( int j = 0; j < node_per_elem; ++j, ++conn )
00383             *conn = elem[order[j]];
00384     }
00385 }
00386 
00387 }  // namespace moab
00388 
00389 #endif
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines