MOAB: Mesh Oriented datABase  (version 5.4.1)
Interface.hpp
Go to the documentation of this file.
00001 /** \mainpage The Mesh-Oriented datABase (MOAB)
00002  *
00003  * MOAB is a component for representing and evaluating mesh data.  MOAB can store
00004  * structured and unstructured mesh, consisting of elements in the finite element “zoo”,
00005  * along with polygons and polyhedra.  The functional interface to MOAB is simple, consisting
00006  * of only four fundamental data types.  This data is quite powerful, allowing the representation
00007  * of most types of metadata commonly found on the mesh.  MOAB is optimized for efficiency in
00008  * space and time, based on access to mesh in chunks rather than through individual entities,
00009  * while also versatile enough to support individual entity access.
00010  *
00011  * The MOAB data model consists of the following four fundamental types: mesh interface instance,
00012  * mesh entities (vertex, edge, tri, etc.), sets, and tags.  Entities are addressed through handles
00013  * rather than pointers, to allow the underlying representation of an entity to change without
00014  * changing the handle to that entity.  Sets are arbitrary groupings of mesh entities and other
00015  * sets.  Sets also support parent/child relationships as a relation distinct from sets containing
00016  * other sets.  The directed-graph provided by set parent/child relationships is useful for modeling
00017  * topological relations from a geometric model and other metadata.  Tags are named data which can
00018  * be assigned to the mesh as a whole, individual entities, or sets.  Tags are a mechanism for
00019  * attaching data to individual entities and sets are a mechanism for describing relations between
00020  * entities; the combination of these two mechanisms is a powerful yet simple interface for
00021  * representing metadata or application-specific data.  For example, sets and tags can be used
00022  * together to describe geometric topology, boundary condition, and inter-processor interface
00023  * groupings in a mesh.
00024  *
00025  * MOAB's API is documented in the moab::Interface class.  Questions and comments should be sent to
00026  * moab-dev _at_ mcs.anl.gov.
00027  *
00028  * \ref userguide "User's Guide"
00029  *
00030  * \ref developerguide "Developer's Guide"
00031  *
00032  * \ref metadata "I/O and Meta-Data Storage Conventions in MOAB"
00033  *
00034  * <a href="pages.html">Full List of Documents</a>
00035  */
00036 
00037 #ifndef MOAB_INTERFACE_HPP
00038 #define MOAB_INTERFACE_HPP
00039 
00040 #include "win32_config.h"
00041 
00042 #define MOAB_API_VERSION        1.01
00043 #define MOAB_API_VERSION_STRING "1.01"
00044 
00045 #include "moab/MOABConfig.h"
00046 #include "moab/Forward.hpp"
00047 #include "moab/Range.hpp"
00048 #include "moab/Compiler.hpp"
00049 #include "moab/ErrorHandler.hpp"
00050 
00051 // include files
00052 #include <string>
00053 #include <functional>
00054 #include <typeinfo>
00055 
00056 //! component architecture definitions
00057 #ifdef XPCOM_MB
00058 
00059 #ifndef __gen_nsISupports_h__
00060 #include "nsISupports.h"
00061 #endif
00062 
00063 #ifndef NS_NO_VTABLE
00064 #define NS_NO_VTABLE
00065 #endif
00066 
00067 #define MBINTERFACE_IID_STR "f728830e-1dd1-11b2-9598-fb9f414f2465"
00068 
00069 #define MBINTERFACE_IID                                    \
00070     {                                                      \
00071         0xf728830e, 0x1dd1, 0x11b2,                        \
00072         {                                                  \
00073             0x95, 0x98, 0xfb, 0x9f, 0x41, 0x4f, 0x24, 0x65 \
00074         }                                                  \
00075     }
00076 
00077 #endif
00078 
00079 #include "moab/UnknownInterface.hpp"
00080 #define MB_INTERFACE_VERSION "2.0.0"
00081 namespace moab
00082 {
00083 
00084 static const MBuuid IDD_MBCore = MBuuid( 0x8956e0a, 0xc300, 0x4005, 0xbd, 0xf6, 0xc3, 0x4e, 0xf7, 0x1f, 0x5a, 0x52 );
00085 
00086 /**
00087  * \class Interface Interface.hpp "moab/Interface.hpp"
00088  * \brief Main interface class to MOAB
00089  * \nosubgrouping
00090  */
00091 #if defined( XPCOM_MB )
00092 class NS_NO_VTABLE Interface : public nsISupports
00093 {
00094 #else
00095 class MOAB_EXPORT Interface : public UnknownInterface
00096 {
00097 #endif
00098 
00099   public:
00100 #ifdef XPCOM_MB
00101     NS_DEFINE_STATIC_IID_ACCESSOR( MBINTERFACE_IID )
00102 #endif
00103 
00104     /** \name Interface */
00105 
00106     /**@{*/
00107 
00108     //! constructor
00109     Interface() {}
00110 
00111     //! destructor
00112     virtual ~Interface() {}
00113 
00114     //! return the entity set representing the whole mesh
00115     virtual EntityHandle get_root_set() = 0;
00116 
00117     //! Get a pointer to an internal MOAB interface
00118     //!\return NULL if not found, iterface pointer otherwise
00119     virtual ErrorCode query_interface_type( const std::type_info& iface_type, void*& iface ) = 0;
00120 
00121     //! Get a pointer to an internal MOAB interface
00122     //!\return NULL if not found, iterface pointer otherwise
00123     template < class IFace >
00124     ErrorCode query_interface( IFace*& ptr )
00125     {
00126         void* tmp_ptr;
00127         ErrorCode result = query_interface_type( typeid( IFace ), tmp_ptr );
00128         ptr              = reinterpret_cast< IFace* >( tmp_ptr );
00129         return result;
00130     }
00131 
00132     //! Release reference to MB interface
00133     virtual ErrorCode release_interface_type( const std::type_info& iface_type, void* iface ) = 0;
00134 
00135     template < class IFace >
00136     ErrorCode release_interface( IFace* interface )
00137     {
00138         return release_interface_type( typeid( IFace ), interface );
00139     }
00140 
00141     //! Release reference to MB interface
00142 
00143     //! Returns the major.minor version number of the interface
00144     /**
00145        \param version_string If non-NULL, will be filled in with a string, possibly
00146        containing implementation-specific information
00147     */
00148     virtual float api_version( std::string* version_string = NULL );
00149 
00150     //! Returns the major.minor version number of the implementation
00151     /**
00152        \param version_string If non-NULL, will be filled in with a string, possibly
00153        containing implementation-specific information
00154     */
00155     virtual float impl_version( std::string* version_string = NULL ) = 0;
00156 
00157     /**@}*/
00158 
00159     /** \name Type and id */
00160 
00161     /**@{*/
00162 
00163     //! Returns the entity type of an EntityHandle.
00164     /** Returns the EntityType (ie, MBVERTEX, MBQUAD, MBHEX ) of <em>handle</em>.
00165         \param handle The EntityHandle you want to find the entity type of.
00166         \return type The entity type of <em>handle</em>.
00167 
00168         Example: \code
00169         EntityType type = type_from_handle( handle);
00170         if( type == MBHEX ) ...  \endcode
00171     */
00172     virtual EntityType type_from_handle( const EntityHandle handle ) const = 0;
00173 
00174     //! Returns the id from an EntityHandle.
00175     /** \param handle The EntityHandle you want to find the id of.
00176         \return id Id of <em>handle</em>
00177 
00178         Example: \code
00179         int id = id_from_handle(handle); \endcode
00180     */
00181     virtual EntityID id_from_handle( const EntityHandle handle ) const = 0;
00182 
00183     //! Returns the topological dimension of an entity
00184     /** Returns the topological dimension of an entity.
00185         \param handle The EntityHandle you want to find the dimension of.
00186         \return type The topological dimension of <em>handle</em>.
00187 
00188         Example: \code
00189         int dim = dimension_from_handle( handle);
00190         if( dim == 0 ) ...  \endcode
00191     */
00192     virtual int dimension_from_handle( const EntityHandle handle ) const = 0;
00193 
00194     //! Gets an entity handle from the data base, if it exists, according to type and id.
00195     /** Given an EntiyType and an id, this function gets the existent EntityHandle.
00196         If no such EntityHandle exits, it returns MB_ENTITY_NOT_FOUND
00197         and sets handle to zero.
00198         \param type The type of the EntityHandle to retrieve from the database.
00199         \param id The id of the EntityHandle to retrieve from the database.
00200         \param handle An EntityHandle of type <em>type</em> and <em>id</em>.
00201 
00202         Example: \code
00203         EntityType handle;
00204         ErrorCode error_code = handle_from_id(MBTRI, 204, handle );
00205         if( error_code == MB_ENTITY_NOT_FOUND ) ... \endcode
00206     */
00207     virtual ErrorCode handle_from_id( const EntityType type, const EntityID, EntityHandle& handle ) const = 0;
00208 
00209     /**@}*/
00210 
00211     /** \name Mesh input/output */
00212 
00213     /**@{*/
00214 
00215     //! Loads a mesh file into the database.
00216     /** Loads the file 'file_name'; types of mesh which can be loaded
00217         depend on modules available at MB compile time.  If
00218         active_block_id_list is NULL, all material sets (blocks in the
00219         ExodusII jargon) are loaded.  Individual material sets  can be
00220         loaded by specifying their ids in 'active_block_id_list'.  All
00221         nodes are loaded on first call for a given file.  Subsequent
00222         calls for a file load any material sets not loaded in previous
00223         calls.
00224         \param file_name Name of file to load into database.
00225         \param active_block_id_list Material set/block ids to load.
00226                 If NULL, ALL blocks of <em>file_name</em> are loaded.
00227         \param num_blocks Number of blocks in active_block_id_list
00228 
00229         Example: \code
00230         std::vector<int> active_block_id_list;
00231         int active_block_id_list[] = {1, 4, 10};
00232         load_mesh( "temp.gen", active_block_id_list, 3 );  //load blocks 1, 4, 10
00233         \endcode
00234     */
00235     virtual ErrorCode load_mesh( const char* file_name,
00236                                  const int* active_block_id_list = NULL,
00237                                  const int num_blocks            = 0 ) = 0;
00238 
00239     /**\brief Load or import a file.
00240      *
00241      * Load a MOAB-native file or import data from some other supported
00242      * file format.
00243      *
00244      *\param file_name The location of the file to read.
00245      *\param file_set  If non-null, this argument must be a pointer to
00246      *                 a valid entity set handle.  All entities read from
00247      *                 the file will be added to this set.  File metadata
00248      *                 will be added to tags on the set.
00249      *\param options A list of string options, separated by semicolons (;).
00250      *               See README.IO for more information.  Options are typically
00251      *               format-specific options or parallel options.  If an
00252      *               option value is unrecognized but the file read otherwise
00253      *               succeeded, MB_UNHANDLED_OPTION will be returned.
00254      *\param set_tag_name The name of a tag used to designate the subset
00255      *               of the file to read.  The name must correspond to
00256      *               data in the file that will be instantiated in MOAB
00257      *               as a tag.
00258      *\param set_tag_values If the name specified in 'set_tag_name'
00259      *               corresponds to a tag with a single integer value,
00260      *               the values in this tag can be used to further
00261      *               limit the subset of data written from the file to
00262      *               only those entities or sets that have a value for
00263      *               the tag that is one of the values in this array.
00264      *\param num_set_tag_values The length of set_tag_values.
00265      *
00266      *\Note file_set is passed by pointer rather than by value (where a
00267      *      zero handle value would indicate no set) so as to intentionally
00268      *      break compatibility with the previous version of this function
00269      *      because the behavior with respect to the file set was changed.
00270      *      The file_set is now an input-only argument.  The previous
00271      *      version of this function unconditionally created a set and
00272      *      passed it back to the caller via a non-const reference argument.
00273      */
00274     virtual ErrorCode load_file( const char* file_name,
00275                                  const EntityHandle* file_set = 0,
00276                                  const char* options          = 0,
00277                                  const char* set_tag_name     = 0,
00278                                  const int* set_tag_values    = 0,
00279                                  int num_set_tag_values       = 0 ) = 0;
00280 
00281     //! Writes mesh to a file.
00282     /** Write mesh to file 'file_name'; if output_list is non-NULL, only
00283         material sets contained in that list will be written.
00284         \param file_name Name of file to write.
00285         \param output_list 1d array of material set handles to write; if
00286                            NULL, all sets are written
00287         \param num_sets Number of sets in output_list array
00288 
00289         Example: \code
00290         EntityHandle output_list[] = {meshset1, meshset2, meshset3};
00291         write_mesh( "output_file.gen", output_list, 3 ); \endcode
00292     */
00293     virtual ErrorCode write_mesh( const char* file_name,
00294                                   const EntityHandle* output_list = NULL,
00295                                   const int num_sets              = 0 ) = 0;
00296 
00297     /**\brief Write or export a file.
00298      *
00299      * Write a MOAB-native file or export data to some other supported
00300      * file format.
00301      *
00302      *\param file_name The location of the file to write.
00303      *\param file_type The type of the file.  If this value is NULL,
00304      *                 then file type will be determined using the
00305      *                 file name suffix.
00306      *\param options   A semicolon-separated list of options.
00307      *                 See README.IO for more information.  Typical options
00308      *                 include the file type, parallel options, and options
00309      *                 specific to certain file formats.
00310      *\param output_sets A list of entity sets to write to the file.  If
00311      *                 no sets are sepcified, the default behavior is to
00312      *                 write all data that is supported by the target file
00313      *                 type.
00314      *\param num_output_sets The length of the output_sets array.
00315      *\param tag_list A list of tags for which to write the tag data.  The
00316      *                write may fail if a tag list is specified but the
00317      *                target file type is not capable of representing the
00318      *                data.  If no tags are specified, the default is to
00319      *                write whatever data the target file format supports.
00320      *\param num_tags The length of tag_list.
00321      */
00322     virtual ErrorCode write_file( const char* file_name,
00323                                   const char* file_type           = 0,
00324                                   const char* options             = 0,
00325                                   const EntityHandle* output_sets = 0,
00326                                   int num_output_sets             = 0,
00327                                   const Tag* tag_list             = 0,
00328                                   int num_tags                    = 0 ) = 0;
00329 
00330     /**\brief Write or export a file.
00331      *
00332      * Write a MOAB-native file or export data to some other supported
00333      * file format.
00334      *
00335      *\param file_name The location of the file to write.
00336      *\param file_type The type of the file.  If this value is NULL,
00337      *                 then file type will be determined using the
00338      *                 file name suffix.
00339      *\param options   A semicolon-separated list of options.
00340      *                 See README.IO for more information.  Typical options
00341      *                 include the file type, parallel options, and options
00342      *                 specific to certain file formats.
00343      *\param output_sets A list of entity sets to write to the file.  If
00344      *                 no sets are sepcified, the default behavior is to
00345      *                 write all data that is supported by the target file
00346      *                 type.
00347      *\param tag_list A list of tags for which to write the tag data.  The
00348      *                write may fail if a tag list is specified but the
00349      *                target file type is not capable of representing the
00350      *                data.  If no tags are specified, the default is to
00351      *                write whatever data the target file format supports.
00352      *\param num_tags The length of tag_list.
00353      */
00354     virtual ErrorCode write_file( const char* file_name,
00355                                   const char* file_type,
00356                                   const char* options,
00357                                   const Range& output_sets,
00358                                   const Tag* tag_list = 0,
00359                                   int num_tags        = 0 ) = 0;
00360 
00361     //! Deletes all mesh entities from this MB instance
00362     virtual ErrorCode delete_mesh() = 0;
00363 
00364     /**@}*/
00365 
00366     /** \name Coordinates and dimensions */
00367 
00368     /**@{*/
00369 
00370     //! Get blocked vertex coordinates for all vertices
00371     /** Blocked = all x, then all y, etc.
00372 
00373     Example: \code
00374     std::vector<double> coords;
00375     get_vertex_coordinates(coords);
00376     double xavg = 0;
00377     for (int i = 0; i < coords.size()/3; i++) xavg += coords[i]; \endcode
00378     */
00379     virtual ErrorCode get_vertex_coordinates( std::vector< double >& coords ) const = 0;
00380 
00381     //! get pointers to coordinate data
00382     /** BEWARE, THIS GIVES ACCESS TO MOAB'S INTERNAL STORAGE, USE WITH CAUTION!
00383      * This function returns pointers to MOAB's internal storage for vertex coordinates.
00384      * Access is similar to tag_iterate, see documentation for that function for details
00385      * about arguments and a coding example.
00386      */
00387     virtual ErrorCode coords_iterate( Range::const_iterator iter,
00388                                       /**< Iterator to first entity you want coordinates for */
00389                                       Range::const_iterator end,
00390                                       /**< Iterator to last entity you want coordinates for */
00391                                       double*& xcoords_ptr,
00392                                       /**< Pointer to x coordinate storage for these entities */
00393                                       double*& ycoords_ptr,
00394                                       /**< Pointer to y coordinate storage for these entities */
00395                                       double*& zcoords_ptr,
00396                                       /**< Pointer to z coordinate storage for these entities */
00397                                       int& count
00398                                       /**< Number of entities for which returned pointers are valid/contiguous */
00399                                       ) = 0;
00400 
00401     //! Gets xyz coordinate information for range of vertices
00402     /** Length of 'coords' should be at least 3*<em>entity_handles.size()</em> before making call.
00403         \param entity_handles Range of vertex handles (error if not of type MBVERTEX)
00404         \param coords Array used to return x, y, and z coordinates.
00405 
00406         Example: \code
00407         double coords[3];
00408         get_coords( vertex_handle, coords );
00409         std::cout<<"x = "<<coords[0]<<std::endl;
00410         std::cout<<"y = "<<coords[1]<<std::endl;
00411         std::cout<<"z = "<<coords[2]<<std::endl; \endcode
00412     */
00413     virtual ErrorCode get_coords( const Range& entity_handles, double* coords ) const = 0;
00414 
00415     //! Gets xyz coordinate information for vector of vertices
00416     /** Identical to range-based function, except entity handles are specified using a 1d vector
00417         and vector length.
00418     */
00419     virtual ErrorCode get_coords( const EntityHandle* entity_handles,
00420                                   const int num_entities,
00421                                   double* coords ) const = 0;
00422 
00423     /**\brief Get vertex coordinates in blocks by dimension.
00424      *
00425      * Get the X, Y, and Z coordinates of a group of vertices.
00426      * Coordinates are returned in separate arrays, one for each
00427      * dimension.  Each coordinate array must be of sufficient
00428      * length to hold the coordinate value for each vertex.  Array
00429      * pointers may be NULL if coordinates in the the respective
00430      * dimension are not desired.
00431      *\param entity_handles  The group of vertex handles for which to get the coordiantes.
00432      *\param x_coords        Output: the X coordinate of each vertex.  May be NULL.
00433      *\param y_coords        Output: the Y coordinate of each vertex.  May be NULL.
00434      *\param z_coords        Output: the Z coordinate of each vertex.  May be NULL.
00435      */
00436     virtual ErrorCode get_coords( const Range& entity_handles,
00437                                   double* x_coords,
00438                                   double* y_coords,
00439                                   double* z_coords ) const = 0;
00440 
00441     //! Sets the xyz coordinates for a vector of vertices
00442     /** An error is returned if any entities in the vector are not vertices.
00443         \param entity_handles EntityHandle's to set coordinates of. (Must be of type MBVERTEX)
00444         \param num_entities Number of entities in entity_handles
00445         \param coords Array containing new xyz coordinates.
00446 
00447         Example: \code
00448         double coords[3] = {0.234, -2.52, 12.023};
00449         set_coords( entity_handle, 1, coords ); \endcode
00450     */
00451     virtual ErrorCode set_coords( const EntityHandle* entity_handles,
00452                                   const int num_entities,
00453                                   const double* coords ) = 0;
00454 
00455     //! Sets the xyz coordinates for a vector of vertices
00456     /** An error is returned if any entities in the vector are not vertices.
00457         \param entity_handles EntityHandle's to set coordinates of. (Must be of type MBVERTEX)
00458         \param num_entities Number of entities in entity_handles
00459         \param coords Array containing new xyz coordinates.
00460 
00461         Example: \code
00462         double coords[3] = {0.234, -2.52, 12.023};
00463         set_coords( entity_handle, 1, coords ); \endcode
00464     */
00465     virtual ErrorCode set_coords( Range entity_handles, const double* coords ) = 0;
00466 
00467     //! Get overall geometric dimension
00468     virtual ErrorCode get_dimension( int& dim ) const = 0;
00469 
00470     //! Set overall geometric dimension
00471     /** Returns error if setting to 3 dimensions, mesh has been created, and
00472      *  there are only 2 dimensions on that mesh
00473      */
00474     virtual ErrorCode set_dimension( const int dim ) = 0;
00475 
00476     /**@}*/
00477 
00478     /** \name Connectivity */
00479 
00480     /**@{*/
00481 
00482     //! get pointers to connectivity data
00483     /** BEWARE, THIS GIVES ACCESS TO MOAB'S INTERNAL STORAGE, USE WITH CAUTION!
00484      * This function returns a pointer to MOAB's internal storage for entity connectivity.
00485      * For each contiguous sub-range of entities, those entities are guaranteed to have
00486      * the same number of vertices (since they're in the same ElementSequence).  Count
00487      * is given in terms of entities, not elements of the connectivity array.
00488      * Access is similar to tag_iterate, see documentation for that function for details
00489      * about arguments and a coding example.
00490      */
00491     virtual ErrorCode connect_iterate( Range::const_iterator iter,
00492                                        /**< Iterator to first entity you want coordinates for */
00493                                        Range::const_iterator end,
00494                                        /**< Iterator to last entity you want coordinates for */
00495                                        EntityHandle*& connect,
00496                                        /**< Pointer to connectivity storage for these entities */
00497                                        int& verts_per_entity,
00498                                        /**< Number of vertices per entity in this block of entities */
00499                                        int& count
00500                                        /**< Number of entities for which returned pointers are valid/contiguous */
00501                                        ) = 0;
00502 
00503     //! Get the connectivity array for all entities of the specified entity type
00504     /**  This function returns the connectivity of just the corner vertices, no higher order nodes
00505          \param type The entity type of elements whose connectivity is to be returned
00506          \param connect an STL vector used to return connectivity array (in the form of entity
00507        handles)
00508     */
00509     virtual ErrorCode get_connectivity_by_type( const EntityType type, std::vector< EntityHandle >& connect ) const = 0;
00510 
00511     //! Gets the connectivity for a vector of elements
00512     /** Same as vector-based version except range is returned (unordered!)
00513      */
00514     virtual ErrorCode get_connectivity( const EntityHandle* entity_handles,
00515                                         const int num_handles,
00516                                         Range& connectivity,
00517                                         bool corners_only = false ) const = 0;
00518 
00519     //! Gets the connectivity for elements
00520     /** Same as vector-based version except range is returned (unordered!)
00521      */
00522     virtual ErrorCode get_connectivity( const Range& entity_handles,
00523                                         Range& connectivity,
00524                                         bool corners_only = false ) const = 0;
00525 
00526     //! Gets the connectivity for a vector of elements
00527     /** Corner vertices or all vertices (including higher-order nodes, if any) are returned.
00528         For non-element handles (ie, MB_MeshSets), returns an error. Connectivity data is copied
00529         from the database into the vector.  Connectivity of a vertex is the same vertex.
00530         The nodes in <em>connectivity</em> are properly ordered according to that element's
00531         canonical ordering.
00532         \param entity_handles Vector of element handles to get connectivity of.
00533         \param num_handles Number of entity handles in <em>entity_handles</em>
00534         \param connectivity Vector in which connectivity of <em>entity_handles</em> is returned.
00535         \param corners_only If true, returns only corner vertices, otherwise returns all of them
00536        (including any higher-order vertices) \param offsets If non-NULL, offsets->[i] stores the
00537        index of the start of entity i's connectivity, with the last value in offsets one beyond the
00538        last entry
00539     */
00540     virtual ErrorCode get_connectivity( const EntityHandle* entity_handles,
00541                                         const int num_handles,
00542                                         std::vector< EntityHandle >& connectivity,
00543                                         bool corners_only           = false,
00544                                         std::vector< int >* offsets = NULL ) const = 0;
00545 
00546     //! Gets a pointer to constant connectivity data of <em>entity_handle</em>
00547     /** Sets <em>number_nodes</em> equal to the number of nodes of the <em>
00548         entity_handle </em>.  Faster then the other <em>get_connectivity</em> function because no
00549         data is copied.  The nodes in 'connectivity' are properly ordered according to the
00550         element's canonical ordering.
00551 
00552 
00553           Example: \code
00554           const EntityHandle* conn;
00555           int number_nodes = 0;
00556           get_connectivity( entity_handle, conn, number_nodes ); \endcode
00557 
00558           Example2: \code
00559           std::vector<EntityHandle> sm_storage;
00560           const EntityHandle* conn;
00561           int number_nodes;
00562           get_connectivity( handle, conn, number_nodes, false, &sm_storage );
00563           if (conn == &sm_storage[0])
00564             std::cout << "Structured mesh element" << std::endl;
00565           \endcode
00566 
00567         \param entity_handle EntityHandle to get connectivity of.
00568         \param connectivity Array in which connectivity of <em>entity_handle</em> is returned.
00569         \param num_nodes Number of MeshVertices in array <em>connectivity</em>.
00570         \param corners_only If true, returns only corner vertices, otherwise returns all of them
00571        (including any higher-order vertices) \param storage Some elements (e.g. structured mesh) may
00572        not have an explicit connectivity list.  This function will normally return
00573        MB_NOT_IMPLEMENTED for such elements.  However, if the caller passes in a non-null value for
00574        this argument, space will be allocated in this vector for the connectivity data and the
00575        connectivity pointer will be set to the data in this vector.
00576     */
00577     virtual ErrorCode get_connectivity( const EntityHandle entity_handle,
00578                                         const EntityHandle*& connectivity,
00579                                         int& num_nodes,
00580                                         bool corners_only                    = false,
00581                                         std::vector< EntityHandle >* storage = 0 ) const = 0;
00582 
00583     //! Sets the connectivity for an EntityHandle.  For non-element handles, return an error.
00584     /** Connectivity is stored exactly as it is ordered in vector <em>connectivity</em>.
00585         \param entity_handle EntityHandle to set connectivity of.
00586         \param connect Vector containing new connectivity of <em>entity_handle</em>.
00587         \param num_connect Number of vertices in <em>connect</em>
00588 
00589         Example: \code
00590         EntityHandle conn[] = {node1, node2, node3};
00591         set_connectivity( tri_element, conn, 3 ); \endcode
00592     */
00593     virtual ErrorCode set_connectivity( const EntityHandle entity_handle,
00594                                         EntityHandle* connect,
00595                                         const int num_connect ) = 0;
00596 
00597     /**@}*/
00598 
00599     /** \name Adjacencies */
00600 
00601     /**@{*/
00602 
00603     //! Get the adjacencies associated with a vector of entities to entities of a specfied
00604     //! dimension.
00605     /** \param from_entities Vector of EntityHandle to get adjacencies of.
00606         \param num_entities Number of entities in <em>from_entities</em>
00607         \param to_dimension Dimension of desired adjacencies
00608         \param create_if_missing If true, MB will create any entities of the specfied dimension
00609         which have not yet been created (only useful when <em>to_dimension <
00610        dim(*from_entities)</em>) \param adj_entities STL vector to which adjacent entities are
00611        appended. \param operation_type Enum of INTERSECT or UNION.  Defines whether to take the
00612        intersection or union of the set of adjacencies recovered for the from_entities.
00613 
00614         The adjacent entities in vector <em>adjacencies</em> are not in any particular
00615         order.
00616 
00617         Example: \code
00618         std::vector<EntityHandle> adjacencies, from_entities = {hex1, hex2};
00619           // generate all edges for these two hexes
00620           get_adjacencies( from_entities, 2, 1, true, adjacencies, Interface::UNION);
00621           adjacencies.clear();
00622             // now find the edges common to both hexes
00623             get_adjacencies( from_entities, 2, 1, false, adjacencies, Interface::INTERSECT);
00624             \endcode
00625     */
00626 
00627     virtual ErrorCode get_adjacencies( const EntityHandle* from_entities,
00628                                        const int num_entities,
00629                                        const int to_dimension,
00630                                        const bool create_if_missing,
00631                                        std::vector< EntityHandle >& adj_entities,
00632                                        const int operation_type = Interface::INTERSECT ) = 0;
00633 
00634     //! Get the adjacencies associated with a vector of entities to entities of a specfied
00635     //! dimension.
00636     /** Identical to vector-based get_adjacencies function, except results are returned in a
00637         range instead of a vector.
00638     */
00639     virtual ErrorCode get_adjacencies( const EntityHandle* from_entities,
00640                                        const int num_entities,
00641                                        const int to_dimension,
00642                                        const bool create_if_missing,
00643                                        Range& adj_entities,
00644                                        const int operation_type = Interface::INTERSECT ) = 0;
00645 
00646     //! Get the adjacencies associated with a range of entities to entities of a specfied dimension.
00647     /** Identical to vector-based get_adjacencies function, except "from" entities specified in a
00648         range instead of a vector.
00649     */
00650     virtual ErrorCode get_adjacencies( const Range& from_entities,
00651                                        const int to_dimension,
00652                                        const bool create_if_missing,
00653                                        Range& adj_entities,
00654                                        const int operation_type = Interface::INTERSECT ) = 0;
00655 
00656     //! Adds adjacencies between "from" and "to" entities.
00657     /** \param from_handle Entities on which the adjacencies are placed
00658         \param to_handles Vector of entities referenced by new adjacencies added to
00659        <em>from_handle</em> \param num_handles Number of entities in <em>to_handles</em> \param
00660        both_ways If true, add the adjacency information in both directions; if false, adjacencies
00661        are added only to <em>from_handle</em>
00662     */
00663     virtual ErrorCode add_adjacencies( const EntityHandle from_handle,
00664                                        const EntityHandle* to_handles,
00665                                        const int num_handles,
00666                                        bool both_ways ) = 0;
00667 
00668     //! Adds adjacencies; same as vector-based, but with range instead
00669     virtual ErrorCode add_adjacencies( const EntityHandle from_handle, Range& adjacencies, bool both_ways ) = 0;
00670 
00671     //! Removes adjacencies between handles.
00672     /** Adjacencies in both directions are removed.
00673         \param from_handle Entity from which adjacencies are being removed.
00674         \param to_handles Entities to which adjacencies are being removed.
00675         \param num_handles Number of handles in <em>to_handles</em>
00676     */
00677     virtual ErrorCode remove_adjacencies( const EntityHandle from_handle,
00678                                           const EntityHandle* to_handles,
00679                                           const int num_handles ) = 0;
00680 
00681     /**\brief Get a ptr to adjacency lists
00682      * Get a pointer to adjacency lists.  These lists are std::vector<EntityHandle>, which are
00683      * pointed to by adjs[i].  Adjacencies are not guaranteed to be in order of increasing
00684      * dimension.  Only a const version of this function is given, because adjacency data is managed
00685      * more carefully in MOAB and should be treated as read-only by applications.  If adjacencies
00686      * have not yet been initialized, adjs_ptr will be NULL (i.e. adjs_ptr == NULL).  There may also
00687      * be NULL entries for individual entities, i.e. adjs_ptr[i] == NULL. \param iter Iterator to
00688      * beginning of entity range desired \param end End iterator for which adjacencies are requested
00689      * \param adjs_ptr Pointer to pointer to const std::vector<EntityHandle>; each member of that
00690      * array is the vector of adjacencies for this entity \param count Number of entities in the
00691      * contiguous chunk starting from *iter
00692      */
00693     virtual ErrorCode adjacencies_iterate( Range::const_iterator iter,
00694                                            Range::const_iterator end,
00695                                            const std::vector< EntityHandle >**& adjs_ptr,
00696                                            int& count ) = 0;
00697     /**@}*/
00698 
00699     //! Enumerated type used in get_adjacencies() and other functions
00700     enum
00701     {
00702         INTERSECT,
00703         UNION
00704     };
00705 
00706     /** \name Getting entities */
00707 
00708     /**@{*/
00709 
00710     //! Retrieves all entities of a given topological dimension in the database or meshset.
00711     /** Appends entities to list passed in.
00712         \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00713         \param dimension Topological dimension of entities desired.
00714         \param entities Range in which entities of dimension <em>dimension</em> are returned.
00715         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00716                          the contents of meshsets, but not the meshsets themselves if true.
00717 
00718         Example: \code
00719           // get 1d (edge) elements in the entire mesh
00720           Range edges;
00721           get_entities_by_dimension( 0, 1, edges );
00722           \endcode
00723     */
00724     virtual ErrorCode get_entities_by_dimension( const EntityHandle meshset,
00725                                                  const int dimension,
00726                                                  Range& entities,
00727                                                  const bool recursive = false ) const = 0;
00728 
00729     //! Retrieves all entities of a given topological dimension in the database or meshset.
00730     /** Appends entities to list passed in.
00731         \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00732         \param dimension Topological dimension of entities desired.
00733         \param entities Range in which entities of dimension <em>dimension</em> are returned.
00734         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00735                          the contents of meshsets, but not the meshsets themselves if true.
00736 
00737         Example: \code
00738           // get 1d (edge) elements in the entire mesh
00739           Range edges;
00740           get_entities_by_dimension( 0, 1, edges );
00741           \endcode
00742     */
00743     virtual ErrorCode get_entities_by_dimension( const EntityHandle meshset,
00744                                                  const int dimension,
00745                                                  std::vector< EntityHandle >& entities,
00746                                                  const bool recursive = false ) const = 0;
00747 
00748     //! Retrieve all entities of a given type in the database or meshset.
00749     /** Appends entities to list passed in.
00750         \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00751         \param type Type of entities to be returned
00752         \param entities Range in which entities of type <em>type</em> are returned.
00753         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00754                          the contents of meshsets, but not the meshsets themselves.  Specifying
00755                          both recursive=true and type=MBENTITYSET is an error, as it would always
00756                          result in an empty list.
00757 
00758         Example: \code
00759           // get the quadrilateral elements in meshset
00760           Range quads;
00761           get_entities_by_type( meshset, MBQUAD, quads );
00762           \endcode
00763     */
00764     virtual ErrorCode get_entities_by_type( const EntityHandle meshset,
00765                                             const EntityType type,
00766                                             Range& entities,
00767                                             const bool recursive = false ) const = 0;
00768 
00769     //! Retrieve all entities of a given type in the database or meshset.
00770     /** Appends entities to list passed in.
00771         \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00772         \param type Type of entities to be returned
00773         \param entities Range in which entities of type <em>type</em> are returned.
00774         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00775                          the contents of meshsets, but not the meshsets themselves.  Specifying
00776                          both recursive=true and type=MBENTITYSET is an error, as it would always
00777                          result in an empty list.
00778 
00779         Example: \code
00780           // get the quadrilateral elements in meshset
00781           Range quads;
00782           get_entities_by_type( meshset, MBQUAD, quads );
00783           \endcode
00784     */
00785     virtual ErrorCode get_entities_by_type( const EntityHandle meshset,
00786                                             const EntityType type,
00787                                             std::vector< EntityHandle >& entities,
00788                                             const bool recursive = false ) const = 0;
00789 
00790     //! Retrieve entities in the database or meshset which have any or all of the tag(s) and
00791     //! (optionally) value(s) specified.
00792     /** \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00793         \param type Type of entities to be returned
00794         \param tag_handles Vector of tag handles entities must have
00795         \param values Vector of pointers to values of tags in <em>tag_handles</em>
00796         \param num_tags Number of tags and values in <em>tag_handles</em> and <em>values</em>
00797         \param entities Range in which entities are returned.
00798         \param condition Boolean condition, either Interface::UNION or Interface::INTERSECT
00799         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00800                          the contents of meshsets, but not the meshsets themselves.  Specifying
00801                          both recursive=true and type=MBENTITYSET is an error, as it would always
00802                          result in an empty list.
00803 
00804         If Interface::UNION is specified as the condition, entities with <em>any</em> of the tags
00805         and values specified are returned.  If Interface::INTERSECT is specified, only entities with
00806         <em>all</em> of the tags/values are returned.
00807 
00808         If <em>values</em> is NULL, entities with the specified tags and any corresponding values
00809        are returned.  Note that if <em>values</em> is non-NULL, it is a vector of <em>pointers</em>
00810        to tag values.
00811 
00812         Example: \code
00813           // get the dirichlet sets in a mesh
00814           Range dir_sets;
00815           Tag dir_tag;
00816           tag_get_handle(DIRICHLET_SET_TAG_NAME, dir_tag, 1, MB_TYPE_INTEGER);
00817           get_entities_by_type_and_tag(0, MBENTITYSET, &dir_tag, NULL, 1, dir_sets,
00818           Interface::UNION);
00819           \endcode
00820     */
00821     virtual ErrorCode get_entities_by_type_and_tag( const EntityHandle meshset,
00822                                                     const EntityType type,
00823                                                     const Tag* tag_handles,
00824                                                     const void* const* values,
00825                                                     const int num_tags,
00826                                                     Range& entities,
00827                                                     const int condition  = Interface::INTERSECT,
00828                                                     const bool recursive = false ) const = 0;
00829 
00830     //! Returns all entities in the data base or meshset, in a range (order not preserved)
00831     /** Appends entities to list passed in.
00832         \param meshset Meshset whose entities are being queried (zero if query is for the entire
00833        mesh). \param entities Range in which entities are returned. \param recursive If true,
00834        meshsets containing meshsets are queried recursively.  Returns the contents of meshsets, but
00835        not the meshsets themselves if true.
00836 
00837         Example: \code
00838         Range entities;
00839           // get all non-meshset entities in meshset, including in contained meshsets
00840           get_entities_by_handle(meshset, entities, true);
00841           \endcode
00842     */
00843     virtual ErrorCode get_entities_by_handle( const EntityHandle meshset,
00844                                               Range& entities,
00845                                               const bool recursive = false ) const = 0;
00846 
00847     //! Returns all entities in the data base or meshset, in a vector (order preserved)
00848     /** Appends entities to list passed in.
00849         \param meshset Meshset whose entities are being queried (zero if query is for the entire
00850        mesh). \param entities STL vector in which entities are returned. \param recursive If true,
00851        meshsets containing meshsets are queried recursively.  Returns the contents of meshsets, but
00852        not the meshsets themselves if true.
00853 
00854         Example: \code
00855         std::vector<EntityHandle> entities;
00856           // get all non-meshset entities in meshset, including in contained meshsets
00857           get_entities_by_handle(meshset, entities, true);
00858           \endcode
00859     */
00860     virtual ErrorCode get_entities_by_handle( const EntityHandle meshset,
00861                                               std::vector< EntityHandle >& entities,
00862                                               const bool recursive = false ) const = 0;
00863 
00864     //! Return the number of entities of given dimension in the database or meshset
00865     /** \param meshset Meshset whose entities are being queried (zero if query is for the entire
00866        mesh). \param dimension Dimension of entities desired. \param num_entities Number of entities
00867        of the given dimension \param recursive If true, meshsets containing meshsets are queried
00868        recursively.  Returns the contents of meshsets, but not the meshsets themselves if true.
00869     */
00870     virtual ErrorCode get_number_entities_by_dimension( const EntityHandle meshset,
00871                                                         const int dimension,
00872                                                         int& num_entities,
00873                                                         const bool recursive = false ) const = 0;
00874 
00875     //! Retrieve the number of entities of a given type in the database or meshset.
00876     /** Identical to get_entities_by_dimension, except returns number instead of entities
00877         \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00878         \param type Type of entities to be returned
00879         \param num_entities Number of entities of type <em>type</em>
00880         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00881                          the contents of meshsets, but not the meshsets themselves.  Specifying
00882                          both recursive=true and type=MBENTITYSET is an error, as it would always
00883                          result in an empty list.
00884     */
00885     virtual ErrorCode get_number_entities_by_type( const EntityHandle meshset,
00886                                                    const EntityType type,
00887                                                    int& num_entities,
00888                                                    const bool recursive = false ) const = 0;
00889 
00890     //! Retrieve number of entities in the database or meshset which have any or all of the
00891     //! tag(s) and (optionally) value(s) specified.
00892     /** Identical to get_entities_by_type_and_tag, except number instead of entities are returned
00893         \param meshset Meshset whose entities are being queried (zero if query is for entire mesh).
00894         \param type Type of entities to be returned
00895         \param tag_handles Vector of tag handles entities must have
00896         \param values Vector of pointers to values of tags in <em>tag_handles</em>
00897         \param num_tags Number of tags and values in <em>tag_handles</em> and <em>values</em>
00898         \param num_entities Range in which number of entities are returned.
00899         \param recursive If true, meshsets containing meshsets are queried recursively.  Returns
00900                          the contents of meshsets, but not the meshsets themselves.  Specifying
00901                          both recursive=true and type=MBENTITYSET is an error, as it would always
00902                          result in an empty list.
00903     */
00904     virtual ErrorCode get_number_entities_by_type_and_tag( const EntityHandle meshset,
00905                                                            const EntityType type,
00906                                                            const Tag* tag_handles,
00907                                                            const void* const* values,
00908                                                            const int num_tags,
00909                                                            int& num_entities,
00910                                                            const int condition  = Interface::INTERSECT,
00911                                                            const bool recursive = false ) const = 0;
00912 
00913     //! Returns number of entities in the data base or meshset
00914     /** Identical to get-entities_by_handle, except number instead of entities are returned
00915         \param meshset Meshset whose entities are being queried (zero if query is for the entire
00916        mesh). \param num_entities Range in which num_entities are returned. \param recursive If
00917        true, meshsets containing meshsets are queried recursively.  Returns the contents of
00918        meshsets, but not the meshsets themselves if true.
00919     */
00920     virtual ErrorCode get_number_entities_by_handle( const EntityHandle meshset,
00921                                                      int& num_entities,
00922                                                      const bool recursive = false ) const = 0;
00923 
00924     /**@}*/
00925 
00926     /** \name Mesh modification */
00927 
00928     /**@{*/
00929 
00930     //! Create an element based on the type and connectivity.
00931     /** Create a new element in the database.  Vertices composing this element must already exist,
00932         and connectivity must be specified in canonical order for the given element type.  If
00933         connectivity vector is not correct for EntityType <em>type</em> (ie, a vector with
00934         3 vertices is passed in to make an MBQUAD), the function returns MB_FAILURE.
00935         \param type Type of element to create. (MBTET, MBTRI, MBKNIFE, etc.)
00936         \param connectivity 1d vector containing connectivity of element to create.
00937         \param num_vertices Number of vertices in element
00938         \param element_handle Handle representing the newly created element in the database.
00939 
00940         Example: \code
00941         EntityHandle quad_conn[] = {vertex0, vertex1, vertex2, vertex3};
00942         EntityHandle quad_handle = 0;
00943         create_element( MBQUAD, quad_conn, 4, quad_handle ); \endcode
00944     */
00945     virtual ErrorCode create_element( const EntityType type,
00946                                       const EntityHandle* connectivity,
00947                                       const int num_vertices,
00948                                       EntityHandle& element_handle ) = 0;
00949 
00950     //! Creates a vertex with the specified coordinates.
00951     /**
00952        \param coordinates Array that has 3 doubles in it.
00953        \param entity_handle EntityHandle representing the newly created vertex in the database.
00954 
00955        Example: \code
00956        double coordinates[] = {1.034, 23.23, -0.432};
00957        EntityHandle new_handle = 0;
00958        create_vertex( coordinates, entity_handle ); \endcode
00959     */
00960     virtual ErrorCode create_vertex( const double coordinates[3], EntityHandle& entity_handle ) = 0;
00961 
00962     //! Create a set of vertices with the specified coordinates
00963     /**
00964        \param coordinates Array that has 3*n doubles in it.
00965        \param nverts Number of vertices to create
00966        \param entity_handles Range passed back with new vertex handles
00967     */
00968     virtual ErrorCode create_vertices( const double* coordinates, const int nverts, Range& entity_handles ) = 0;
00969 
00970     //! Merge two entities into a single entity
00971     /** Merge two entities into a single entities, with <em>entity_to_keep</em> receiving
00972         adjacencies that were on <em>entity_to_remove</em>.
00973         \param entity_to_keep Entity to be kept after merge
00974         \param entity_to_remove Entity to be merged into <em>entity_to_keep</em>
00975         \param auto_merge If false, <em>entity_to_keep</em> and <em>entity_to_remove</em> must share
00976         the same lower-dimensional entities; if true, MB tries to merge those entities automatically
00977         \param delete_removed_entity If true, <em>entity_to_remove</em> is deleted after merge is
00978        complete
00979     */
00980     virtual ErrorCode merge_entities( EntityHandle entity_to_keep,
00981                                       EntityHandle entity_to_remove,
00982                                       bool auto_merge,
00983                                       bool delete_removed_entity ) = 0;
00984 
00985     //! Removes entities in a vector from the data base.
00986     /** If any of the entities are contained in any meshsets, it is removed from those meshsets
00987         which were created with MESHSET_TRACK_OWNER option bit set.  Tags for <em>entity</em> are
00988         removed as part of this function.
00989         \param entities 1d vector of entities to delete
00990         \param num_entities Number of entities in 1d vector
00991     */
00992     virtual ErrorCode delete_entities( const EntityHandle* entities, const int num_entities ) = 0;
00993 
00994     //! Removes entities in a range from the data base.
00995     /** If any of the entities are contained in any meshsets, it is removed from those meshsets
00996         which were created with MESHSET_TRACK_OWNER option bit set.  Tags for <em>entity</em> are
00997         removed as part of this function.
00998         \param entities Range of entities to delete
00999     */
01000     virtual ErrorCode delete_entities( const Range& entities ) = 0;
01001 
01002     /**@}*/
01003 
01004     /** \name Information */
01005 
01006     /**@{*/
01007 
01008     //! List entities to standard output
01009     /** Lists all data pertaining to entities (i.e. vertex coordinates if vertices, connectivity if
01010         elements, set membership if set).  Useful for debugging, but output can become quite long
01011         for large databases.
01012     */
01013     virtual ErrorCode list_entities( const Range& entities ) const = 0;
01014 
01015     //! List entities, or number of entities in database, to standard output
01016     /** Lists data pertaining to entities to standard output.  If <em>entities</em> is NULL and
01017         <em>num_entities</em> is zero, lists only the number of entities of each type in the
01018         database.  If <em>entities</em> is NULL and <em>num_entities</em> is non-zero, lists all
01019         information for all entities in the database.
01020         \param entities 1d vector of entities to list
01021         \param num_entities Number of entities in 1d vector
01022     */
01023     virtual ErrorCode list_entities( const EntityHandle* entities, const int num_entities ) const = 0;
01024 
01025     //! List a single entity; no header printed
01026     /** Lists a single entity, including its connectivity and its adjacencies.
01027      *  No header is printed, because calling function might print information between header
01028      *  and information printed by this function.
01029      *  \param entity The entity to be listed.
01030      */
01031     virtual ErrorCode list_entity( const EntityHandle entity ) const = 0;
01032 
01033     //! Return information about the last error
01034     /** \param info std::string into which information on the last error is written.
01035      */
01036     virtual ErrorCode get_last_error( std::string& info ) const = 0;
01037 
01038     //! Return string representation of given error code
01039     /** \param code Error code for which string is wanted
01040      */
01041     virtual std::string get_error_string( const ErrorCode code ) const = 0;
01042 
01043     /**\brief Calculate amount of memory used to store MOAB data
01044      *
01045      * This function calculates the amount of memory used to store
01046      * MOAB data.
01047      *
01048      * There are two possible values for each catagory of memory use.
01049      * The exact value and the amortized value.  The exact value is the
01050      * amount of memory used to store the data for the specified entities.
01051      * The amortized value includes the exact value and an amortized
01052      * estimate of the memory consumed in overhead for storing the values
01053      * (indexing structures, access structures, etc.)
01054      *
01055      * Note: If ent_array is NULL and total_amortized_storage is *not* NULL,
01056      *       the total memory used by MOAB for storing data all will be
01057      *       returned in the address pointed to by total_amortized_storage.
01058      *
01059      *\param ent_array Array of entities for which to estimate the memory use.
01060      *                 If NULL, estimate is done for all entities.
01061      *\param num_ents The length of ent_array.  Not used if ent_rray is NULL.
01062      *\param total_(amortized_)storage The sum of the memory entity, adjacency,
01063      *                   and all tag storage.
01064      *\param (amortized_)entity_storage The storage for the entity definitions
01065      *                   (connectivity arrays for elements, coordinates for
01066      *                   vertices, list storage within sets, etc.)
01067      *\param (amortized_)adjacency_storage The storage for adjacency data.
01068      *\param tag_array   An array of tags for which to calculate the memory use.
01069      *\param num_tags    The lenght of tag_array
01070      *\param (amortized_)tag_storage If tag_array is not NULL, then one value
01071      *                   for each tag specifying the memory used for storing
01072      *                   that tag.  If tag_array is NULL and this value is not,
01073      *                   the location at which to store the total memory used
01074      *                   for all tags.
01075      */
01076     virtual void estimated_memory_use( const EntityHandle* ent_array                   = 0,
01077                                        unsigned long num_ents                          = 0,
01078                                        unsigned long long* total_storage               = 0,
01079                                        unsigned long long* total_amortized_storage     = 0,
01080                                        unsigned long long* entity_storage              = 0,
01081                                        unsigned long long* amortized_entity_storage    = 0,
01082                                        unsigned long long* adjacency_storage           = 0,
01083                                        unsigned long long* amortized_adjacency_storage = 0,
01084                                        const Tag* tag_array                            = 0,
01085                                        unsigned num_tags                               = 0,
01086                                        unsigned long long* tag_storage                 = 0,
01087                                        unsigned long long* amortized_tag_storage       = 0 ) = 0;
01088 
01089     /**\brief Calculate amount of memory used to store MOAB data
01090      *
01091      * This function calculates the amount of memory used to store
01092      * MOAB data.
01093      *
01094      * There are two possible values for each catagory of memory use.
01095      * The exact value and the amortized value.  The exact value is the
01096      * amount of memory used to store the data for the specified entities.
01097      * The amortized value includes the exact value and an amortized
01098      * estimate of the memory consumed in overhead for storing the values
01099      * (indexing structures, access structures, etc.)
01100      *
01101      *\param ents        Entities for which to estimate the memory use.
01102      *\param total_(amortized_)storage The sum of the memory entity, adjacency,
01103      *                   and all tag storage.
01104      *\param (amortized_)entity_storage The storage for the entity definitions
01105      *                   (connectivity arrays for elements, coordinates for
01106      *                   vertices, list storage within sets, etc.)
01107      *\param (amortized_)adjacency_storage The storage for adjacency data.
01108      *\param tag_array   An array of tags for which to calculate the memory use.
01109      *\param num_tags    The lenght of tag_array
01110      *\param (amortized_)tag_storage If tag_array is not NULL, then one value
01111      *                   for each tag specifying the memory used for storing
01112      *                   that tag.  If tag_array is NULL and this value is not,
01113      *                   the location at which to store the total memory used
01114      *                   for all tags.
01115      */
01116     virtual void estimated_memory_use( const Range& ents,
01117                                        unsigned long long* total_storage               = 0,
01118                                        unsigned long long* total_amortized_storage     = 0,
01119                                        unsigned long long* entity_storage              = 0,
01120                                        unsigned long long* amortized_entity_storage    = 0,
01121                                        unsigned long long* adjacency_storage           = 0,
01122                                        unsigned long long* amortized_adjacency_storage = 0,
01123                                        const Tag* tag_array                            = 0,
01124                                        unsigned num_tags                               = 0,
01125                                        unsigned long long* tag_storage                 = 0,
01126                                        unsigned long long* amortized_tag_storage       = 0 ) = 0;
01127     /**@}*/
01128 
01129     /** \name Higher-order elements */
01130 
01131     /**@{*/
01132 
01133     //! function object for recieving events from MB of higher order nodes added to entities
01134     class MOAB_EXPORT HONodeAddedRemoved
01135     {
01136       public:
01137         //! Constructor
01138         HONodeAddedRemoved() {}
01139 
01140         //! Destructor
01141         virtual ~HONodeAddedRemoved() {}
01142 
01143         //! node_added called when a node was added to an element's connectivity array
01144         //! note: connectivity array of element may be incomplete (corner nodes will exist always)
01145         /**
01146          * \param node Node being added
01147          * \param element Element node is being added to
01148          */
01149         virtual void node_added( EntityHandle node, EntityHandle element ) = 0;
01150 
01151         //! node_added called when a node was added to an element's connectivity array
01152         //! note: connectivity array of element may be incomplete (corner nodes will exist always)
01153         /**
01154          * \param node Node being removed.
01155          */
01156         virtual void node_removed( EntityHandle node ) = 0;
01157     };
01158 
01159     //! Convert entities to higher-order elements by adding mid nodes
01160     /** This function causes MB to create mid-nodes on all edges, faces, and element interiors
01161         for all entities in <em>meshset</em>.  Higher order nodes appear in an element's
01162        connectivity array according to the algorithm described in the documentation for Mesh.  If
01163         <em>HONodeAddedRemoved</em> function is input, this function is called to notify the
01164        application of nodes being added/removed from the mesh. \param meshset The set of entities
01165        being converted \param mid_edge If true, mid-edge nodes are created \param mid_face If true,
01166        mid-face nodes are created \param mid_region If true, mid-element nodes are created \param
01167        function_object If non-NULL, the node_added or node_removed functions on this object are
01168        called when nodes are added or removed from an entity, respectively
01169     */
01170     virtual ErrorCode convert_entities( const EntityHandle meshset,
01171                                         const bool mid_edge,
01172                                         const bool mid_face,
01173                                         const bool mid_region,
01174                                         HONodeAddedRemoved* function_object = 0 ) = 0;
01175 
01176     //! Returns the side number, in canonical ordering, of <em>child</em> with respect to
01177     //! <em>parent</em>
01178     /** Given a parent and child entity, returns the canonical ordering information side number,
01179        sense, and offset of <em>child</em> with respect to <em>parent</em>.  This function returns
01180         MB_FAILURE if <em>child</em> is not related to <em>parent</em>.  This function does *not*
01181         create adjacencies between <em>parent</em> and <em>child</em>.
01182         \param parent Parent entity to be compared
01183         \param child Child entity to be compared
01184         \param side_number Side number in canonical ordering of <em>child</em> with respect to
01185         <em>parent</em>
01186         \param sense Sense of <em>child</em> with respect to <em>parent</em>, assuming ordering of
01187         <em>child</em> as given by get_connectivity called on <em>child</em>; sense is 1, -1
01188         for forward/reverse sense, resp.
01189         \param offset Offset between first vertex of <em>child</em> and first vertex of side
01190         <em>side_number</em> on <em>parent</em>
01191     */
01192     virtual ErrorCode side_number( const EntityHandle parent,
01193                                    const EntityHandle child,
01194                                    int& side_number,
01195                                    int& sense,
01196                                    int& offset ) const = 0;
01197 
01198     //! Find the higher-order node on a subfacet of an entity
01199     /** Given an entity and the connectivity and type of one of its subfacets, find the
01200         high order node on that subfacet, if any.  The number of vertices in <em>subfacet_conn</em>
01201         is derived from <em>subfacet_type</em> and the canonical numbering for that type.
01202         \param parent_handle The element whose subfacet is being queried
01203         \param subfacet_conn The connectivity of the subfacet being queried
01204         \param subfacet_type The type of subfacet being queried
01205         \param high_order_node If the subfacet has a high-order node defined on
01206        <em>parent_handle</em>, the handle for that node.
01207     */
01208     virtual ErrorCode high_order_node( const EntityHandle parent_handle,
01209                                        const EntityHandle* subfacet_conn,
01210                                        const EntityType subfacet_type,
01211                                        EntityHandle& high_order_node ) const = 0;
01212 
01213     //! Return the handle of the side element of a given dimension and index
01214     /** Given a parent entity and a target dimension and side number, return the handle of
01215         the entity corresponding to that side.  If an entity has not been created to represent
01216         that side, one is not created by this function, and zero is returned in
01217        <em>target_entity</em>. \param source_entity The entity whose side is being queried. \param
01218        dim The topological dimension of the side being queried. \param side_number The canonical
01219        index of the side being queried. \param target_entity The handle of the entity representing
01220        this side, if any.
01221     */
01222     virtual ErrorCode side_element( const EntityHandle source_entity,
01223                                     const int dim,
01224                                     const int side_number,
01225                                     EntityHandle& target_entity ) const = 0;
01226 
01227     /**@}*/
01228 
01229     /** \name Tags */
01230 
01231     /**@{*/
01232 
01233     /**\brief Get a tag handle, possibly creating the tag
01234      *
01235      * Get a handle used to associate application-defined values
01236      * with MOAB entities.  If the tag does not already exist then
01237      * \c flags should contain exactly one of \c MB_TAG_SPARSE,
01238      * \c MB_TAG_DENSE, \c MB_TAG_MESH unless \c type is MB_TYPE_BIT,
01239      * which implies \c MB_TAG_BIT storage.
01240      * .
01241      *\param name          The tag name
01242      *\param size          Tag size as number of values of of data type per entity
01243      *                     (or number of bytes if \c MB_TAG_BYTES is passed in flags).  If \c
01244      *MB_TAG_VARLEN is specified, this value is taken to be the size of the default value if one is
01245      *specified and is otherwise ignored. \param type          The type of the data (used for IO)
01246      *\param tag_handle    Output: the resulting tag handle.
01247      *\param flags         Bitwise OR of values from TagType
01248      *\param default_value Optional default value for tag.
01249      *\param created       Optional returned boolean indicating that the tag
01250      *                     was created.
01251      *\return - \c MB_ALREADY_ALLOCATED     if tag exists and \c MB_TAG_EXCL is specified, or
01252      *default values do not match (and \c MB_TAG_ANY or \c MB_TAG_DFTOK not specified).
01253      *        - \c MB_TAG_NOT_FOUND         if tag does not exist and \c MB_TAG_CREAT is not
01254      *specified
01255      *        - \c MB_INVALID_SIZE          if tag value size is not a multiple of the size of the
01256      *data type (and \c MB_TAG_ANY not specified).
01257      *        - \c MB_TYPE_OUT_OF_RANGE     invalid or inconsistent parameter
01258      *        - \c MB_VARIABLE_DATA_LENGTH  if \c MB_TAG_VARLEN and \c default_value is non-null and
01259      *                                      \c default_value_size is not specified.
01260      *
01261      *\NOTE A call to tag_get_handle that includes a default value will fail
01262      * if the tag already exists with a different default value.  A call without
01263      * a default value will succeed if the tag already exists, regardless of
01264      * whether or not the existing tag has a default value.
01265      *
01266      * Examples:
01267      *
01268      * Retrieve a handle for an existing tag, returning a non-success error
01269      * code if the tag does not exist or does not store 1 integer value per
01270      * entity:
01271      *\code
01272      * Tag git_tag;
01273      * mb.tag_get_handle( GLOBAL_ID_TAG_NAME, 1, MB_TYPE_INTEGER, gid_tag );
01274      * \endcode
01275      * Get the tag handle, or create it as a dense tag if it does not already
01276      * exist:
01277      *\code
01278      * Tag gid_tag;
01279      * mb.tag_get_handle( GLOBAL_ID_TAG_NAME, 1, MB_TYPE_INTEGER, gid_tag, MB_TAG_CREAT|MB_TAG_BIT
01280      *); \endcode Create the tag or *fail* if it already exists: \code Tag gid_tag;
01281      * mb.tag_get_handle( GLOBAL_ID_TAG_NAME, 1, MB_TYPE_INTEGER, gid_tag, MB_TAG_EXCL|MB_TAG_DENSE
01282      *); \endcode Get an existing variable length tag, failing if it does not exist or is not
01283      *variable-length or does not contain double values. \code Tag vtag; mb.tag_get_handle(
01284      *tag_name, 0, MB_TYPE_DOUBLE, vtag ); \endcode Get the same variable-length tag, but create it
01285      *with a default value if it doesn't exist.  Note that if the tag already exists this call will
01286      *return a non-success error code if the existing tag has a different default value. \code Tag
01287      *vtag; const double default_val = M_PI; const int def_val_len = 1; mb.tag_get_handle( tag_name,
01288      *def_val_len, MB_TYPE_DOUBLE, vtag, MB_TAG_SPARSE|MB_TAG_VARLEN|MB_TAG_CREAT, &default_val );
01289      * \endcode
01290      */
01291     virtual ErrorCode tag_get_handle( const char* name,
01292                                       int size,
01293                                       DataType type,
01294                                       Tag& tag_handle,
01295                                       unsigned flags            = 0,
01296                                       const void* default_value = 0,
01297                                       bool* created             = 0 ) = 0;
01298 
01299     /**\brief same as non-const version, except that MB_TAG_CREAT flag is ignored. */
01300     virtual ErrorCode tag_get_handle( const char* name,
01301                                       int size,
01302                                       DataType type,
01303                                       Tag& tag_handle,
01304                                       unsigned flags            = 0,
01305                                       const void* default_value = 0 ) const = 0;
01306 
01307     //! Get the name of a tag corresponding to a handle
01308     /** \param tag_handle Tag you want the name of.
01309         \param tag_name Name string for <em>tag_handle</em>.
01310     */
01311     virtual ErrorCode tag_get_name( const Tag tag_handle, std::string& tag_name ) const = 0;
01312 
01313     /**\brief Gets the tag handle corresponding to a name
01314 
01315         If a tag of that name does not exist, returns MB_TAG_NOT_FOUND
01316         \param tag_name Name of the desired tag.
01317         \param tag_handle Tag handle corresponding to <em>tag_name</em>
01318     */
01319     virtual ErrorCode tag_get_handle( const char* tag_name, Tag& tag_handle ) const = 0;
01320 
01321     //! Get the size of the specified tag in bytes
01322     /** Get the size of the specified tag, in bytes (MB_TAG_SPARSE, MB_TAG_DENSE, MB_TAG_MESH)
01323         \note always returns 1 for bit tags.
01324         \param tag Handle of the desired tag.
01325         \param bytes_per_tag Size of the specified tag
01326         \return - MB_TAG_NOT_FOUND for invalid tag handles
01327                 - MB_VARIABLE_DATA_LENGTH for variable-length tags
01328                 - MB_SUCCESS otherwise
01329     */
01330     virtual ErrorCode tag_get_bytes( const Tag tag, int& bytes_per_tag ) const = 0;
01331 
01332     //! Get the array length of a tag
01333     /** Get the size of the specified tag, in as the number of values of
01334         the basic type (e.g. number of integer values for each tag value if
01335         the data type is MB_TYPE_INTEGER).  Gives number of bits for bit tags
01336         and is the same as \c tag_get_bytes for MB_TYPE_OPAQUE tags.
01337         \param tag Handle of the desired tag.
01338         \param length Size of the specified tag
01339         \return - MB_TAG_NOT_FOUND for invalid tag handles
01340                 - MB_VARIABLE_DATA_LENGTH for variable-length tags
01341                 - MB_SUCCESS otherwise
01342     */
01343     virtual ErrorCode tag_get_length( const Tag tag, int& length ) const = 0;
01344 
01345     //! Get the type of the specified tag
01346     /** Get the type of the specified tag
01347         \param tag Handle of the desired tag.
01348         \param tag_type Type of the specified tag
01349     */
01350     virtual ErrorCode tag_get_type( const Tag tag, TagType& tag_type ) const = 0;
01351 
01352     /** \brief Get data type of tag.
01353      *
01354      * Get the type of the tag data.  The tag is data is assumed to
01355      * be a vector of this type.  If the tag data vetcor contains
01356      * more than one value, then the tag size must be a multiple of
01357      * the size of this type.
01358      * \param tag  The tag
01359      * \param type The type of the specified tag (output).
01360      */
01361     virtual ErrorCode tag_get_data_type( const Tag tag, DataType& type ) const = 0;
01362 
01363     //! Get the default value of the specified tag
01364     /** Get the default value of the specified tag
01365         \param tag Handle of the desired tag.
01366         \param def_value Pointer to memory where default value of the specified tag is written
01367         \return - MB_ENTITY_NOT_FOUND If no default value is set for tag.
01368                 - MB_SUCCESS          If success.
01369                 - MB_FAILURE          If <code>def_val</code> is NULL.
01370                 - MB_TAG_NOT_FOUND    If <code>tag_handle</code> is invalid.
01371     */
01372     virtual ErrorCode tag_get_default_value( const Tag tag, void* def_val ) const             = 0;
01373     virtual ErrorCode tag_get_default_value( Tag tag, const void*& def_val, int& size ) const = 0;
01374 
01375     //! Get handles for all tags defined in the mesh instance
01376     /** Get handles for all tags defined on the mesh instance.
01377         \param tag_handles STL vector of all tags
01378     */
01379     virtual ErrorCode tag_get_tags( std::vector< Tag >& tag_handles ) const = 0;
01380 
01381     //! Get handles for all tags defined on this entity
01382     /** Get handles for all tags defined on this entity; if zero, get all tags defined
01383         on mesh instance
01384         \param entity Entity for which you want tags
01385         \param tag_handles STL vector of all tags defined on <em>entity</em>
01386     */
01387     virtual ErrorCode tag_get_tags_on_entity( const EntityHandle entity, std::vector< Tag >& tag_handles ) const = 0;
01388 
01389     //! Get the value of the indicated tag on the specified entities in the specified vector
01390     /** Get the value of the indicated tag on the specified entities; <em>tag_data</em> must contain
01391         enough space (i.e. tag_size*num_entities bytes) to hold all tag data.  MOAB does
01392        <em>not</em> check whether this space is available before writing to it. \note For bit tags,
01393        tag_data must contain one byte per entity.  For each entity, the corresponding byte will
01394        contain the tag bits in the lower bit positions and zero bits in the higher. \param
01395        tag_handle Tag whose values are being queried. \param entity_handles 1d vector of entity
01396        handles whose tag values are being queried \param num_entities Number of entities in 1d
01397        vector of entity handles \param tag_data Pointer to memory into which tag data will be
01398        written
01399     */
01400     virtual ErrorCode tag_get_data( const Tag tag_handle,
01401                                     const EntityHandle* entity_handles,
01402                                     int num_entities,
01403                                     void* tag_data ) const = 0;
01404 
01405     //! Get the value of the indicated tag on the specified entities in the specified range
01406     /** Identical to previous function, except entities are specified using a range instead of a 1d
01407        vector. \param tag_handle Tag whose values are being queried. \param entity_handles Range of
01408        entity handles whose tag values are being queried \param tag_data Pointer to memory into
01409        which tag data will be written
01410     */
01411     virtual ErrorCode tag_get_data( const Tag tag_handle, const Range& entity_handles, void* tag_data ) const = 0;
01412 
01413     //! Set the value of the indicated tag on the specified entities in the specified vector
01414     /** Set the value of the indicated tag on the specified entities; <em>tag_data</em> contains the
01415         values, <em>one value per entity in <em>entity_handles</em></em>.
01416         \note For bit tags, tag_data must contain one byte per entity.  For each
01417               entity, the tag bits will be read from the lower bits of the
01418               corresponding byte.
01419         \param tag_handle Tag whose values are being set
01420         \param entity_handles 1d vector of entity handles whose tag values are being set
01421         \param num_entities Number of entities in 1d vector of entity handles
01422         \param tag_data Pointer to memory holding tag values to be set, <em>one entry per entity
01423        handle</em>
01424     */
01425     virtual ErrorCode tag_set_data( Tag tag_handle,
01426                                     const EntityHandle* entity_handles,
01427                                     int num_entities,
01428                                     const void* tag_data ) = 0;
01429 
01430     //! Set the value of the indicated tag on the specified entities in the specified range
01431     /** Identical to previous function, except entities are specified using a range instead of a 1d
01432        vector. \param tag_handle Tag whose values are being set \param entity_handles Range of
01433        entity handles whose tag values are being set \param tag_data Pointer to memory holding tag
01434        values to be set, <em>one entry per entity handle</em>
01435     */
01436     virtual ErrorCode tag_set_data( Tag tag_handle, const Range& entity_handles, const void* tag_data ) = 0;
01437 
01438     /**\brief Get pointers to tag data
01439      *
01440      * For a tag, get the values for a list of passed entity handles.
01441      *\note  This function may not be used for bit tags.
01442      *\param tag_handle     The tag
01443      *\param entity_handles An array of entity handles for which to retreive tag values.
01444      *\param num_entities   The length of the 'entity_handles' array.
01445      *\param tag_data       An array of 'const void*'.  Array must be at least
01446      *                      'num_entitities' long.  Array is populated (output)
01447      *                      with pointers to the internal storage for the
01448      *                      tag value corresponding to each entity handle.
01449      *\param tag_sizes      The length of each tag value.  Optional for
01450      *                      fixed-length tags.  Required for variable-length tags.
01451      */
01452     virtual ErrorCode tag_get_by_ptr( const Tag tag_handle,
01453                                       const EntityHandle* entity_handles,
01454                                       int num_entities,
01455                                       const void** tag_data,
01456                                       int* tag_sizes = 0 ) const = 0;
01457 
01458     /**\brief Get pointers to tag data
01459      *
01460      * For a tag, get the values for a list of passed entity handles.
01461      *\note  This function may not be used for bit tags.
01462      *\param tag_handle     The tag
01463      *\param entity_handles The entity handles for which to retreive tag values.
01464      *\param tag_data       An array of 'const void*'.  Array is populated (output)
01465      *                      with pointers to the internal storage for the
01466      *                      tag value corresponding to each entity handle.
01467      *\param tag_sizes      The length of each tag value.  Optional for
01468      *                      fixed-length tags.  Required for variable-length tags.
01469      */
01470     virtual ErrorCode tag_get_by_ptr( const Tag tag_handle,
01471                                       const Range& entity_handles,
01472                                       const void** tag_data,
01473                                       int* tag_sizes = 0 ) const = 0;
01474 
01475     /**\brief Set tag data given an array of pointers to tag values.
01476      *
01477      * For a tag, set the values for a list of passed entity handles.
01478      *\note  This function may not be used for bit tags.
01479      *\param tag_handle     The tag
01480      *\param entity_handles An array of entity handles for which to set tag values.
01481      *\param num_entities   The length of the 'entity_handles' array.
01482      *\param tag_data       An array of 'const void*'.  Array must be at least
01483      *                      'num_entitities' long.  Array is expected to
01484      *                      contain pointers to tag values for the corresponding
01485      *                      EntityHandle in 'entity_handles'.
01486      *\param tag_sizes      The length of each tag value.  Optional for
01487      *                      fixed-length tags.  Required for variable-length tags.
01488      */
01489     virtual ErrorCode tag_set_by_ptr( Tag tag_handle,
01490                                       const EntityHandle* entity_handles,
01491                                       int num_entities,
01492                                       void const* const* tag_data,
01493                                       const int* tag_sizes = 0 ) = 0;
01494 
01495     /**\brief Set tag data given an array of pointers to tag values.
01496      *
01497      * For a tag, set the values for a list of passed entity handles.
01498      *\note  This function may not be used for bit tags.
01499      *\param tag_handle     The tag
01500      *\param entity_handles The entity handles for which to set tag values.
01501      *\param tag_data       An array of 'const void*'.  Array is expected to
01502      *                      contain pointers to tag values for the corresponding
01503      *                      EntityHandle in 'entity_handles'.
01504      *\param tag_sizes      The length of each tag value.  Optional for
01505      *                      fixed-length tags.  Required for variable-length tags.
01506      */
01507     virtual ErrorCode tag_set_by_ptr( Tag tag_handle,
01508                                       const Range& entity_handles,
01509                                       void const* const* tag_data,
01510                                       const int* tag_sizes = 0 ) = 0;
01511 
01512     /**\brief Set tag data given value.
01513      *
01514      * For a tag, set the values for a list of passed entity handles to
01515      * the same, specified value.
01516      *
01517      *\param tag_handle     The tag
01518      *\param entity_handles The entity handles for which to set tag values.
01519      *\param tag_data       A pointer to the tag value.
01520      *\param tag_sizes      For variable-length tags, the length of the
01521      *                      tag value.  This argument will be ignored for
01522      *                      fixed-length tags.
01523      */
01524     virtual ErrorCode tag_clear_data( Tag tag_handle,
01525                                       const Range& entity_handles,
01526                                       const void* value,
01527                                       int value_size = 0 ) = 0;
01528 
01529     /**\brief Set tag data given value.
01530      *
01531      * For a tag, set the values for a list of passed entity handles to
01532      * the same, specified value.
01533      *
01534      *\param tag_handle     The tag
01535      *\param entity_handles The entity handles for which to set tag values.
01536      *\param tag_data       A pointer to the tag value.
01537      *\param tag_sizes      For variable-length tags, the length of the
01538      *                      tag value.  This argument will be ignored for
01539      *                      fixed-length tags.
01540      */
01541     virtual ErrorCode tag_clear_data( Tag tag_handle,
01542                                       const EntityHandle* entity_handles,
01543                                       int num_entity_handles,
01544                                       const void* value,
01545                                       int value_size = 0 ) = 0;
01546 
01547     //! Delete the data of a vector of entity handles and sparse tag
01548     /** Delete the data of a tag on a vector of entity handles.  Only sparse tag data are deleted
01549        with this function; dense tags are deleted by deleting the tag itself using tag_delete.
01550         \param tag_handle Handle of the (sparse) tag being deleted from entity
01551         \param entity_handles 1d vector of entity handles from which the tag is being deleted
01552         \param num_handles Number of entity handles in 1d vector
01553     */
01554     virtual ErrorCode tag_delete_data( Tag tag_handle, const EntityHandle* entity_handles, int num_handles ) = 0;
01555 
01556     //! Delete the data of a range of entity handles and sparse tag
01557     /** Delete the data of a tag on a range of entity handles.  Only sparse tag data are deleted
01558        with this function; dense tags are deleted by deleting the tag itself using tag_delete.
01559         \param tag_handle Handle of the (sparse) tag being deleted from entity
01560         \param entity_range Range of entities from which the tag is being deleted
01561     */
01562     virtual ErrorCode tag_delete_data( Tag tag_handle, const Range& entity_range ) = 0;
01563 
01564     /**\brief Access tag data via direct pointer into contiguous blocks
01565      *
01566      * Iteratively obtain direct access to contiguous blocks of tag
01567      * storage.  This function cannot be used with bit tags because
01568      * of the compressed bit storage.  This function cannot be used
01569      * with variable length tags because it does not provide a mechanism
01570      * to determine the length of the value for each entity.  This
01571      * function may be used with sparse tags, but if it is used, it
01572      * will return data for a single entity at a time.
01573      *
01574      *\param tag_handle  The handle of the tag for which to access data
01575      *\param iter        The first entity for which to return data.
01576      *\param end         One past the last entity for which data is desired.
01577      *\param count       The number of entities for which data was returned
01578      *\param data_ptr    Output: pointer to tag storage.
01579      *\param allocate    If true, space for this tag will be allocated, if not it wont
01580      *
01581      *\Note If this function is called for entities for which no tag value
01582      *      has been set, but for which a default value exists, it will
01583      *      force the allocation of explicit storage for each such entity
01584      *      even though MOAB would normally not explicitly store tag values
01585      *      for such entities.
01586      *
01587      *\Example:
01588      *\code
01589      * Range ents; // range to iterate over
01590      * Tag tag; // tag for which to access data
01591      * int bytes;
01592      * ErrorCode err = mb.tag_get_size( tag, bytes );
01593      * if (err) { ... }
01594      *
01595      * ...
01596      * Range::iterator iter = ents.begin();
01597      * while (iter != ents.end()) {
01598      *   int count;
01599      *    // get contiguous block of tag dat
01600      *   void* ptr;
01601      *   err = mb.tag_iterate( tag, iter, ents.end(), count, ptr );
01602      *   if (err) { ... }
01603      *    // do something with tag data
01604      *   process_Data( ptr, count );
01605      *    // advance to next block of data
01606      *   iter += count;
01607      * }
01608      *\endcode
01609      */
01610     virtual ErrorCode tag_iterate( Tag tag_handle,
01611                                    Range::const_iterator begin,
01612                                    Range::const_iterator end,
01613                                    int& count,
01614                                    void*& data_ptr,
01615                                    bool allocate = true ) = 0;
01616 
01617     //! Remove a tag from the database and delete all of its associated data
01618     /** Deletes a tag and all associated data.
01619      */
01620     virtual ErrorCode tag_delete( Tag tag_handle ) = 0;
01621 
01622     /**@}*/
01623 
01624     /** \name Sets */
01625 
01626     /**@{*/
01627 
01628     //! Create a new mesh set
01629     /** Create a new mesh set.  Meshsets can store entities ordered or unordered. A set can include
01630        entities at most once (MESHSET_SET) or more than once.  Meshsets can optionally track its
01631        members using adjacencies (MESHSET_TRACK_OWNER); if set, entities are deleted from tracking
01632        meshsets before being deleted.  This adds data to mesh entities, which can be expensive.
01633         \param options Options bitmask for the new meshset, possible values defined above
01634         \param ms_handle Handle for the meshset created
01635     */
01636     virtual ErrorCode create_meshset( const unsigned int options, EntityHandle& ms_handle, int start_id = 0 ) = 0;
01637 
01638     //! Empty a vector of mesh set
01639     /** Empty a mesh set.
01640         \param ms_handles 1d vector of handles of sets being emptied
01641         \param num_meshsets Number of entities in 1d vector
01642     */
01643     virtual ErrorCode clear_meshset( const EntityHandle* ms_handles, const int num_meshsets ) = 0;
01644 
01645     //! Empty a range of mesh set
01646     /** Empty a mesh set.
01647         \param ms_handles Range of handles of sets being emptied
01648     */
01649     virtual ErrorCode clear_meshset( const Range& ms_handles ) = 0;
01650 
01651     //! Get the options of a mesh set
01652     /** Get the options of a mesh set.
01653         \param ms_handle Handle for mesh set being queried
01654         \param options Bit mask in which mesh set options are returned
01655     */
01656     virtual ErrorCode get_meshset_options( const EntityHandle ms_handle, unsigned int& options ) const = 0;
01657 
01658     //! Set the options of a mesh set
01659     /** Set the options of a mesh set.
01660         \param ms_handle Handle for meshset whose options are being changed
01661         \param options Bit mask of options to be used
01662     */
01663     virtual ErrorCode set_meshset_options( const EntityHandle ms_handle, const unsigned int options ) = 0;
01664 
01665     //! Subtract meshsets
01666     /** Subtract <em>meshset2</em> from <em>meshset1</em>, placing the results in meshset1.
01667         \param meshset1 Mesh set being subtracted from, also used to pass back result
01668         \param meshset2 Mesh set being subtracted from <em>meshset1</em>
01669     */
01670     virtual ErrorCode subtract_meshset( EntityHandle meshset1, const EntityHandle meshset2 ) = 0;
01671 
01672     //! Intersect meshsets
01673     /** Intersect <em>meshset1</em> with <em>meshset2</em>, placing the results in meshset1.
01674         \param meshset1 Mesh set being intersected, also used to pass back result
01675         \param meshset2 Mesh set being intersected with <em>meshset1</em>
01676     */
01677     virtual ErrorCode intersect_meshset( EntityHandle meshset1, const EntityHandle meshset2 ) = 0;
01678 
01679     //! Unite meshsets
01680     /** Unite <em>meshset1</em> with <em>meshset2</em>, placing the results in meshset1.
01681         \param meshset1 Mesh set being united, also used to pass back result
01682         \param meshset2 Mesh set being united with <em>meshset1</em>
01683     */
01684     virtual ErrorCode unite_meshset( EntityHandle meshset1, const EntityHandle meshset2 ) = 0;
01685 
01686     //! Add to a meshset entities in specified range
01687     /** Add to a meshset entities in specified range.  If <em>meshset</em> has MESHSET_TRACK_OWNER
01688         option set, adjacencies are also added to entities in <em>entities</em>.
01689         \param meshset Mesh set being added to
01690         \param entities Range of entities being added to meshset
01691     */
01692     virtual ErrorCode add_entities( EntityHandle meshset, const Range& entities ) = 0;
01693 
01694     //! Add to a meshset entities in specified vector
01695     /** Add to a meshset entities in specified vector.  If <em>meshset</em> has MESHSET_TRACK_OWNER
01696         option set, adjacencies are also added to entities in <em>entities</em>.
01697         \param meshset Mesh set being added to
01698         \param entities 1d vector of entities being added to meshset
01699         \param num_entities Number of entities in 1d vector
01700     */
01701     virtual ErrorCode add_entities( EntityHandle meshset, const EntityHandle* entities, const int num_entities ) = 0;
01702 
01703     //! Remove from a meshset entities in specified range
01704     /** Remove from a meshset entities in specified range.  If <em>meshset</em> has
01705        MESHSET_TRACK_OWNER option set, adjacencies in entities in <em>entities</em> are updated.
01706         \param meshset Mesh set being removed from
01707         \param entities Range of entities being removed from meshset
01708     */
01709     virtual ErrorCode remove_entities( EntityHandle meshset, const Range& entities ) = 0;
01710 
01711     //! Remove from a meshset entities in specified vector
01712     /** Remove from a meshset entities in specified vector.  If <em>meshset</em> has
01713        MESHSET_TRACK_OWNER option set, adjacencies in entities in <em>entities</em> are updated.
01714         \param meshset Mesh set being removed from
01715         \param entities 1d vector of entities being removed from meshset
01716         \param num_entities Number of entities in 1d vector
01717     */
01718     virtual ErrorCode remove_entities( EntityHandle meshset, const EntityHandle* entities, const int num_entities ) = 0;
01719 
01720     //! Return whether a set contains entities
01721     /** Return whether a set contains entities.  Returns true only
01722      * if ALL entities are contained
01723      * \param meshset Mesh set being queried
01724      * \param entities Entities being queried
01725      * \param num_entities Number of entities
01726      * \return bool If true, all entities are contained in set
01727      */
01728     virtual bool contains_entities( EntityHandle meshset,
01729                                     const EntityHandle* entities,
01730                                     int num_entities,
01731                                     const int operation_type = Interface::INTERSECT ) = 0;
01732 
01733     //! Replace entities in a set with other entities
01734     /** Replace entities in a set with other entities
01735      *
01736      * \note  Behavior is undefined if an entity handle exists in both the
01737      *        old_entities and the new_entities arrays or old_entities
01738      *        contains multiple copies of an entity.
01739      * \note  If an entity occurs multiple times in an ordered set, all
01740      *        occurances will be replaced.
01741      * \note  For list-based sets, if not all handles in old_entities
01742      *        occcur in the set, the corresponding new_entities will not
01743      *        be added and  MB_ENTITY_NOT_FOUND will be returned.
01744      *        For set-based sets, all entities in new_entities wll be
01745      *        added and any contained entities in old_entities will be
01746      *        removed, and the return value will be MB_SUCCESS.
01747      * \param meshset Mesh set being modified
01748      * \param old_entities Entities to replace
01749      * \param new_entities New entities to add
01750      * \param num_entities Number of entities in input arrays
01751      * \return - MB_SUCCESS : all entities in old_entities replaced
01752      *         - MB_ENTITY_NOT_FOUND : one or more entities in new_entities
01753      *                        not added to set because corresponding entity
01754      *                        in old_entities did not occur in the ordered
01755      *                        set.
01756      *         - MB_FAILURE : internal error
01757      */
01758     virtual ErrorCode replace_entities( EntityHandle meshset,
01759                                         const EntityHandle* old_entities,
01760                                         const EntityHandle* new_entities,
01761                                         int num_entities ) = 0;
01762     /**@}*/
01763 
01764     /** \name Set parents/children */
01765 
01766     /**@{*/
01767 
01768     //! Get parent mesh sets of a mesh set
01769     /** If <em>num_hops</em> is 1, only immediate parents are returned.  If <em>num_hops</em> is
01770        zero, all ancenstors are returned.  Otherwise, <em>num_hops</em> specifies the maximum number
01771        of generations to traverse. \param meshset The mesh set whose parents are being queried
01772         \param parents STL vector holding the parents returned by this function
01773         \param num_hops Number of generations to traverse (0 = all)
01774     */
01775     virtual ErrorCode get_parent_meshsets( const EntityHandle meshset,
01776                                            std::vector< EntityHandle >& parents,
01777                                            const int num_hops = 1 ) const = 0;
01778 
01779     //! Get parent mesh sets of a mesh set
01780     /** If <em>num_hops</em> is 1, only immediate parents are returned.  If <em>num_hops</em> is
01781        zero, all ancenstors are returned.  Otherwise, <em>num_hops</em> specifies the maximum number
01782        of generations to traverse. \param meshset The mesh set whose parents are being queried
01783         \param parents Range holding the parents returned by this function
01784         \param num_hops Number of generations to traverse (0 = all)
01785     */
01786     virtual ErrorCode get_parent_meshsets( const EntityHandle meshset,
01787                                            Range& parents,
01788                                            const int num_hops = 1 ) const = 0;
01789 
01790     //! Get child mesh sets of a mesh set
01791     /** If <em>num_hops</em> is 1, only immediate children are returned.  If <em>num_hops</em> is
01792        zero, all ancenstors are returned.  Otherwise, <em>num_hops</em> specifies the maximum number
01793        of generations to traverse. \param meshset The mesh set whose children are being queried
01794         \param children STL vector holding the children returned by this function
01795         \param num_hops Number of generations to traverse (0 = all)
01796     */
01797     virtual ErrorCode get_child_meshsets( const EntityHandle meshset,
01798                                           std::vector< EntityHandle >& children,
01799                                           const int num_hops = 1 ) const = 0;
01800 
01801     //! Get child mesh sets of a mesh set
01802     /** If <em>num_hops</em> is 1, only immediate children are returned.  If <em>num_hops</em> is
01803        zero, all ancenstors are returned.  Otherwise, <em>num_hops</em> specifies the maximum number
01804        of generations to traverse. \param meshset The mesh set whose children are being queried
01805         \param children Range holding the children returned by this function
01806         \param num_hops Number of generations to traverse (0 = all)
01807     */
01808     virtual ErrorCode get_child_meshsets( const EntityHandle meshset,
01809                                           Range& children,
01810                                           const int num_hops = 1 ) const = 0;
01811 
01812     /**\brief Get mesh sets contained in a mesh set
01813      *
01814      * If <em>num_hops</em> is 1, only immediate contents are returned.
01815      * Otherwise a recursive query of all contained sets is performed,
01816      * returning every visted set.  The value of num_hops limites the
01817      * depth of the search, with zero indicating no depth limit.
01818      *
01819      *\param meshset The mesh set whose contents are being queried
01820      *\param contained The result list.
01821      *\param num_hops Number of generations to traverse (0 = all)
01822      */
01823     virtual ErrorCode get_contained_meshsets( const EntityHandle meshset,
01824                                               std::vector< EntityHandle >& contained,
01825                                               const int num_hops = 1 ) const = 0;
01826 
01827     /**\brief Get mesh sets contained in a mesh set
01828      *
01829      * If <em>num_hops</em> is 1, only immediate contents are returned.
01830      * Otherwise a recursive query of all contained sets is performed,
01831      * returning every visted set.  The value of num_hops limites the
01832      * depth of the search, with zero indicating no depth limit.
01833      *
01834      *\param meshset The mesh set whose contents are being queried
01835      *\param contained The result list.
01836      *\param num_hops Number of generations to traverse (0 = all)
01837      */
01838     virtual ErrorCode get_contained_meshsets( const EntityHandle meshset,
01839                                               Range& contained,
01840                                               const int num_hops = 1 ) const = 0;
01841 
01842     //! Get the number of parent mesh sets of a mesh set
01843     /** Identical to get_parent_meshsets, only number is returned instead of actual parents.
01844         \param meshset The mesh set whose parents are being queried
01845         \param number Number of parents
01846     */
01847     virtual ErrorCode num_parent_meshsets( const EntityHandle meshset, int* number, const int num_hops = 1 ) const = 0;
01848 
01849     //! Get the number of child mesh sets of a mesh set
01850     /** Identical to get_child_meshsets, only number is returned instead of actual children.
01851         \param meshset The mesh set whose children are being queried
01852         \param number Number of children
01853     */
01854     virtual ErrorCode num_child_meshsets( const EntityHandle meshset, int* number, const int num_hops = 1 ) const = 0;
01855 
01856     /**\brief Get the number of mesh sets contained in a mesh set
01857      *
01858      * Return the number of sets that would be returned by get_contained_meshsets
01859      *
01860      *\param meshset  The initial set to begin the query from.
01861      *\param number   (Output) The result count.
01862      *\param num_hops Search depth (0 => unbounded).
01863      */
01864     virtual ErrorCode num_contained_meshsets( const EntityHandle meshset,
01865                                               int* number,
01866                                               const int num_hops = 1 ) const = 0;
01867 
01868     //! Add a parent mesh set to a mesh set
01869     /** Make <em>parent_meshset</em> a new parent of <em>child_meshset</em>.  This function does
01870         <em>not</em> add a corresponding child link to <em>parent_meshset</em>.
01871         \param child_meshset The child mesh set being given a new parent.
01872         \param parent_meshset The parent being added to <em>child_meshset</em>
01873     */
01874     virtual ErrorCode add_parent_meshset( EntityHandle child_meshset, const EntityHandle parent_meshset ) = 0;
01875 
01876     //! Add a parent mesh sets to a mesh set
01877     /** Make <em>parent_meshset</em> a new parent of <em>child_meshset</em>.  This function does
01878         <em>not</em> add a corresponding child link to <em>parent_meshset</em>.
01879         \param child_meshset The child mesh set being given a new parent.
01880         \param parent_meshset The parent being added to <em>child_meshset</em>
01881     */
01882     virtual ErrorCode add_parent_meshsets( EntityHandle child_meshset,
01883                                            const EntityHandle* parent_meshsets,
01884                                            int num_parent_meshsets ) = 0;
01885 
01886     //! Add a child mesh set to a mesh set
01887     /** Make <em>child_meshset</em> a new child of <em>parent_meshset</em>.  This function does
01888         <em>not</em> add a corresponding parent link to <em>child_meshset</em>.
01889         \param parent_meshset The parent mesh set being given a new child.
01890         \param child_meshset The child being added to <em>parent_meshset</em>
01891     */
01892     virtual ErrorCode add_child_meshset( EntityHandle parent_meshset, const EntityHandle child_meshset ) = 0;
01893 
01894     //! Add a child mesh sets to a mesh set
01895     /** Make <em>child_meshset</em> a new child of <em>parent_meshset</em>.  This function does
01896         <em>not</em> add a corresponding parent link to <em>child_meshset</em>.
01897         \param parent_meshset The parent mesh set being given a new child.
01898         \param child_meshset The child being added to <em>parent_meshset</em>
01899     */
01900     virtual ErrorCode add_child_meshsets( EntityHandle parent_meshset,
01901                                           const EntityHandle* child_meshsets,
01902                                           int num_child_meshsets ) = 0;
01903 
01904     //! Add parent and child links between mesh sets
01905     /** Makes <em>child_meshset</em> a new child of <em>parent_meshset</em>, and vica versa.
01906         \param parent The parent mesh set being given a new child, and the new parent
01907         \param child The child being given a new parent, and the new child
01908     */
01909     virtual ErrorCode add_parent_child( EntityHandle parent, EntityHandle child ) = 0;
01910 
01911     //! Remove parent and child links between mesh sets
01912     /** Removes parent/child links between <em>child_meshset</em> and <em>parent_meshset</em>.
01913         \param parent The parent mesh set being removed from <em>child</em>
01914         \param child The child mesh set being removed from <em>parent</em>
01915     */
01916     virtual ErrorCode remove_parent_child( EntityHandle parent, EntityHandle child ) = 0;
01917 
01918     //! Remove a parent mesh set from a mesh set
01919     /** Removes <em>parent_meshset</em> from the parents of <em>child_meshset</em>.  This function
01920        does <em>not</em> remove a corresponding child link from <em>parent_meshset</em>. \param
01921        child_meshset The child mesh whose parent is being removed \param parent_meshset The parent
01922        being removed from <em>meshset</em>
01923     */
01924     virtual ErrorCode remove_parent_meshset( EntityHandle child_meshset, const EntityHandle parent_meshset ) = 0;
01925 
01926     //! Remove a child mesh set from a mesh set
01927     /** Removes <em>child_meshset</em> from the children of <em>parent_meshset</em>.  This function
01928        does <em>not</em> remove a corresponding parent link from <em>child_meshset</em>. \param
01929        parent_meshset The parent mesh set whose child is being removed \param child_meshset The
01930        child being removed from <em>parent_meshset</em>
01931     */
01932     virtual ErrorCode remove_child_meshset( EntityHandle parent_meshset, const EntityHandle child_meshset ) = 0;
01933 
01934     // ! Returns the global id tag; default value is -1
01935     virtual Tag globalId_tag() = 0;
01936 
01937     /**@}*/
01938 
01939     /** \name Set iterators */
01940 
01941     /**@{*/
01942 
01943     /** \brief Create an iterator over the set
01944      * Create a new iterator that iterates over entities with the specified type or dimension.
01945      * Only one of ent_type or dim can be set; use dim=-1 or ent_type=MBMAXTYPE for the other.
01946      * Iterators for list-type (ordered) sets are stable over set modification, unless entity
01947      * removed or deleted is the one at the current position of the iterator.  If the check_valid
01948      * parameter is passed as true, entities are checked for validity before being passed back by
01949      * get_next_entities function (checking entity validity can have a non-negligible cost).
01950      *
01951      * Iterators returned by this function can be deleted using the normal C++ delete function.
01952      * After creating the iterator through this function, further interactions are through methods
01953      * on the SetIterator class.
01954      * \param meshset The entity set associated with this iterator (use 0 for whole instance)
01955      * \param ent_type Entity type associated with this iterator
01956      * \param ent_dim Dimension associated with this iterator
01957      * \param chunk_size Chunk size of the iterator
01958      * \param check_valid If true, entities are checked for validity before being returned
01959      */
01960 
01961     virtual ErrorCode create_set_iterator( EntityHandle meshset,
01962                                            EntityType ent_type,
01963                                            int ent_dim,
01964                                            int chunk_size,
01965                                            bool check_valid,
01966                                            SetIterator*& set_iter ) = 0;
01967     /**@}*/
01968 
01969     // ************************  Interface options controllable by user  ***************
01970 
01971     /** \name Sequence Option controllers */
01972 
01973     /**@{*/
01974 
01975     /** \brief Interface to control memory allocation for sequences
01976      * Provide a factor that controls the size of the sequence that gets allocated.
01977      * This is typically useful in the parallel setting when a-priori, the number of ghost entities
01978      * and the memory required for them within the same sequence as the owned entities are unknown.
01979      * The default factor is 1.0 but this can be appropriately updated at runtime so that we do not
01980      * have broken sequences.
01981      */
01982     virtual double get_sequence_multiplier() const = 0;
01983 
01984     /** \brief Interface to control memory allocation for sequences
01985      * Provide a factor that controls the size of the sequence that gets allocated.
01986      * This is typically useful in the parallel setting when a-priori, the number of ghost entities
01987      * and the memory required for them within the same sequence as the owned entities are unknown.
01988      * The default factor is 1.0 but this can be appropriately updated at runtime so that we do not
01989      * have broken sequences.
01990      *
01991      * \param meshset User specified multiplier (should be greater than 1.0)
01992      */
01993     virtual void set_sequence_multiplier( double factor ) = 0;
01994 
01995     /**@}*/
01996 };
01997 
01998 //! predicate for STL algorithms.  Returns true if the entity handle is
01999 //! of the specified type.  For example, to remove all the tris out of a list
02000 //! of 2D entities retrieved using get_adjacencies you could do
02001 //! std::remove_if(list.begin(), list.end(), type_equals(gMB, MBTRI));
02002 class type_equals
02003 {
02004     // deprecation of unary_function
02005     typedef EntityHandle argument_type;
02006     typedef bool result_type;
02007 
02008   public:
02009     //! interface object
02010     Interface* meshDB;
02011 
02012     //! type corresponding to this predicate
02013     const EntityType test_type;
02014 
02015     //! Constructor
02016     type_equals( Interface* mdb, const EntityType type ) : meshDB( mdb ), test_type( type ) {}
02017 
02018     //! operator predicate
02019     bool operator()( EntityHandle handle ) const
02020     {
02021         return ( meshDB->type_from_handle( handle ) == test_type );
02022     }
02023 };
02024 
02025 //! predicate for STL algorithms.  Returns true if the entity handle is not
02026 //! of the specified type.  For example, to remove all but the tris out of a list
02027 //! of 2D entities retrieved using get_adjacencies you could do
02028 //! std::remove_if(list.begin(), list.end(), type_not_equals(gMB, MBTRI));
02029 class type_not_equals
02030 {
02031     // deprecation of unary_function
02032     typedef EntityHandle argument_type;
02033     typedef bool result_type;
02034 
02035   public:
02036     //! interface object
02037     Interface* meshDB;
02038 
02039     //! type corresponding to this predicate
02040     const EntityType test_type;
02041 
02042     //! Constructor
02043     type_not_equals( Interface* mdb, const EntityType type ) : meshDB( mdb ), test_type( type ) {}
02044 
02045     //! operator predicate
02046     bool operator()( EntityHandle handle ) const
02047     {
02048         return ( meshDB->type_from_handle( handle ) != test_type );
02049     }
02050 };
02051 
02052 inline float Interface::api_version( std::string* version_string )
02053 {
02054     if( NULL != version_string )
02055         *version_string = std::string( "MOAB API version " ) + std::string( MOAB_API_VERSION_STRING );
02056     return MOAB_API_VERSION;
02057 }
02058 
02059 }  // namespace moab
02060 
02061 #endif  // MB_INTERFACE_HPP
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines