MOAB: Mesh Oriented datABase
(version 5.4.1)
|
00001 /** 00002 * MOAB, a Mesh-Oriented datABase, is a software component for creating, 00003 * storing and accessing finite element mesh data. 00004 * 00005 * Copyright 2004 Sandia Corporation. Under the terms of Contract 00006 * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government 00007 * retains certain rights in this software. 00008 * 00009 * This library is free software; you can redistribute it and/or 00010 * modify it under the terms of the GNU Lesser General Public 00011 * License as published by the Free Software Foundation; either 00012 * version 2.1 of the License, or (at your option) any later version. 00013 * 00014 */ 00015 00016 #ifndef MOAB_IMPL_GENERAL_HPP 00017 #define MOAB_IMPL_GENERAL_HPP 00018 00019 #include "moab/Interface.hpp" 00020 #include "moab/ReaderIface.hpp" 00021 #include <map> 00022 #include <list> 00023 00024 namespace moab 00025 { 00026 00027 class WriteUtil; 00028 class ReadUtil; 00029 class ScdInterface; 00030 class AEntityFactory; 00031 class SequenceManager; 00032 class Error; 00033 class HomCoord; 00034 class ReaderWriterSet; 00035 class EntitySequence; 00036 class FileOptions; 00037 class SetIterator; 00038 00039 #ifdef MOAB_HAVE_AHF 00040 class HalfFacetRep; 00041 #endif 00042 00043 /**\class Core 00044 * \brief Implementation of MOAB Interface 00045 * Implementation of the MOAB Interface class. You shouldn't call functions directly 00046 * on an object of type Core (use Interface instead), unless you really have to access 00047 * non-API functionality. 00048 */ 00049 class Core : public Interface 00050 { 00051 00052 public: 00053 friend class SetIterator; 00054 00055 //! constructor 00056 Core(); 00057 00058 //! destructor 00059 ~Core(); 00060 00061 //! Get a pointer to an internal MOAB interface 00062 //!\return NULL if not found, iterface pointer otherwise 00063 virtual ErrorCode query_interface_type( const std::type_info& interface_type, void*& ptr ); 00064 00065 //! Release reference to MB interface 00066 virtual ErrorCode release_interface_type( const std::type_info& interface_type, void* iface ); 00067 00068 virtual int QueryInterface( const MBuuid& uuid, UnknownInterface** iface ); 00069 00070 //! Returns the major.minor version number of the implementation 00071 /** 00072 \param iface_name If non-NULL, will be filled in with a string, possibly 00073 containing implementation-specific information 00074 */ 00075 virtual float impl_version( std::string* version_string = NULL ); 00076 00077 //! get the type from a handle, returns type 00078 virtual EntityType type_from_handle( const EntityHandle handle ) const; 00079 00080 //! get the id from a handle, returns id 00081 virtual EntityID id_from_handle( const EntityHandle handle ) const; 00082 00083 //! get a handle from an id and type 00084 virtual ErrorCode handle_from_id( const EntityType type, const EntityID id, EntityHandle& handle ) const; 00085 00086 virtual int dimension_from_handle( const EntityHandle ) const; 00087 00088 //! load mesh from data in file 00089 //! NOTE: if there is mesh already present, the new mesh will be added 00090 virtual ErrorCode load_mesh( const char* file_name, 00091 const int* active_block_id_list = NULL, 00092 const int num_blocks = 0 ); 00093 00094 /**Load or import a file. */ 00095 virtual ErrorCode load_file( const char* file_name, 00096 const EntityHandle* file_set = 0, 00097 const char* options = 0, 00098 const char* set_tag_name = 0, 00099 const int* set_tag_vals = 0, 00100 int num_set_tag_vals = 0 ); 00101 00102 /**Load or import a file. */ 00103 ErrorCode serial_load_file( const char* file_name, 00104 const EntityHandle* file_set, 00105 const FileOptions& opts, 00106 const ReaderIface::SubsetList* subsets = 0, 00107 const Tag* file_id_tag = 0 ); 00108 00109 ErrorCode serial_read_tag( const char* file_name, 00110 const char* tag_name, 00111 const FileOptions& opts, 00112 std::vector< int >& tag_vals, 00113 const ReaderIface::SubsetList* subsets = 0 ); 00114 00115 virtual ErrorCode write_mesh( const char* file_name, 00116 const EntityHandle* output_list = NULL, 00117 const int num_sets = 0 ); 00118 /** Write or export a file. */ 00119 virtual ErrorCode write_file( const char* file_name, 00120 const char* file_type = 0, 00121 const char* options = 0, 00122 const EntityHandle* output_sets = 0, 00123 int num_output_sets = 0, 00124 const Tag* tag_list = 0, 00125 int num_tags = 0 ); 00126 00127 /** Write or export a file */ 00128 virtual ErrorCode write_file( const char* file_name, 00129 const char* file_type, 00130 const char* options, 00131 const Range& output_sets, 00132 const Tag* tag_list = 0, 00133 int num_tags = 0 ); 00134 00135 //! deletes all mesh entities from this datastore 00136 virtual ErrorCode delete_mesh(); 00137 00138 //! get overall geometric dimension 00139 virtual ErrorCode get_dimension( int& dim ) const; 00140 00141 //! set overall geometric dimension 00142 /** Returns error if setting to 3 dimensions, mesh has been created, and 00143 * there are only 2 dimensions on that mesh 00144 */ 00145 virtual ErrorCode set_dimension( const int dim ); 00146 00147 //! get blocked vertex coordinates for all vertices 00148 /** Blocked = all x, then all y, etc. 00149 */ 00150 virtual ErrorCode get_vertex_coordinates( std::vector< double >& coords ) const; 00151 00152 //! get pointers to coordinate data 00153 virtual ErrorCode coords_iterate( Range::const_iterator iter, 00154 Range::const_iterator end, 00155 double*& xcoords_ptr, 00156 double*& ycoords_ptr, 00157 double*& zcoords_ptr, 00158 int& count ); 00159 00160 //! get the coordinate information for this handle if it is of type Vertex 00161 //! otherwise, return an error 00162 virtual ErrorCode get_coords( const Range& entities, double* coords ) const; 00163 00164 virtual ErrorCode get_coords( const EntityHandle* entities, const int num_entities, double* coords ) const; 00165 00166 virtual ErrorCode get_coords( const EntityHandle entity_handle, 00167 const double*& x, 00168 const double*& y, 00169 const double*& z ) const; 00170 00171 virtual ErrorCode get_coords( const Range& entities, double* x_coords, double* y_coords, double* z_coords ) const; 00172 00173 //! set the coordinate information for this handle if it is of type Vertex 00174 //! otherwise, return an error 00175 virtual ErrorCode set_coords( const EntityHandle* entity_handles, const int num_entities, const double* coords ); 00176 00177 //! set the coordinate information for this handle if it is of type Vertex 00178 //! otherwise, return an error 00179 virtual ErrorCode set_coords( Range entity_handles, const double* coords ); 00180 00181 //! get global connectivity array for specified entity type 00182 /** Assumes just vertices, no higher order nodes 00183 */ 00184 virtual ErrorCode get_connectivity_by_type( const EntityType type, std::vector< EntityHandle >& connect ) const; 00185 00186 //! get pointer to connectivity data 00187 virtual ErrorCode connect_iterate( Range::const_iterator iter, 00188 Range::const_iterator end, 00189 EntityHandle*& connect, 00190 int& verts_per_entity, 00191 int& count ); 00192 00193 //! Gets the connectivity for an element EntityHandle. 00194 /** For non-element handles (ie, MeshSets), 00195 * returns an error. Connectivity data is copied from the database into the vector 00196 * <em>connectivity</em>. The nodes in <em>connectivity</em> are properly ordered. 00197 * \param entity_handle EntityHandle to get connectivity of. 00198 * \param connectivity Vector in which connectivity of <em>entity_handle</em> is returned. 00199 * Should contain MeshVertices. 00200 * \param corners_only If true, returns only corner vertices, otherwise returns all of them 00201 * (including any higher-order vertices) 00202 * 00203 * Example: \code 00204 * std::vector<EntityHandle> conn; 00205 * get_connectivity( entity_handle, conn ); \endcode 00206 */ 00207 virtual ErrorCode get_connectivity( const EntityHandle* entity_handles, 00208 const int num_handles, 00209 std::vector< EntityHandle >& connectivity, 00210 bool corners_only = false, 00211 std::vector< int >* offsets = NULL ) const; 00212 00213 //! Gets the connectivity for a vector of elements 00214 /** Same as vector-based version except range is returned (unordered!) 00215 */ 00216 virtual ErrorCode get_connectivity( const EntityHandle* entity_handles, 00217 const int num_handles, 00218 Range& connectivity, 00219 bool corners_only = false ) const; 00220 00221 //! Gets the connectivity for elements 00222 /** Same as vector-based version except range is returned (unordered!) 00223 */ 00224 virtual ErrorCode get_connectivity( const Range& from_entities, 00225 Range& adj_entities, 00226 bool corners_only = false ) const; 00227 00228 //! Gets a pointer to constant connectivity data of <em>entity_handle</em> 00229 /** Sets <em>number_nodes</em> equal to the number of nodes of the <em> 00230 entity_handle </em>. Faster then the other <em>get_connectivity</em> function. 00231 The nodes in 'connectivity' are properly ordered. 00232 \param entity_handle EntityHandle to get connectivity of. 00233 \param connectivity Array in which connectivity of <em>entity_handle</em> is returned. 00234 Should contain MBVERTEX's. 00235 \param num_nodes Number of MeshVertices in array <em>connectivity</em>. 00236 00237 Example: \code 00238 const EntityHandle* conn; 00239 int number_nodes = 0; 00240 get_connectivity( entity_handle, conn, number_nodes ); \endcode 00241 00242 Example2: \code 00243 std::vector<EntityHandle> sm_storage; 00244 const EntityHandle* conn; 00245 int number_nodes; 00246 get_connectivity( handle, conn, number_nodes, false, &sm_storage ); 00247 if (conn == &sm_storage[0]) 00248 std::cout << "Structured mesh element" << std::endl; 00249 \endcode 00250 */ 00251 virtual ErrorCode get_connectivity( const EntityHandle entity_handle, 00252 const EntityHandle*& connectivity, 00253 int& number_nodes, 00254 bool corners_only = false, 00255 std::vector< EntityHandle >* storage = 0 ) const; 00256 00257 //! Sets the connectivity for an EntityHandle. For non-element handles, return an error. 00258 /** Connectivity is stored exactly as it is ordered in vector <em>connectivity</em>. 00259 \param entity_handle EntityHandle to set connectivity of. 00260 \param connect Vector containing new connectivity of <em>entity_handle</em>. 00261 \param num_connect Number of vertices in <em>connect</em> 00262 00263 Example: \code 00264 std::vector<EntityHandle> conn(3); 00265 conn[0] = node1; 00266 conn[1] = node2; 00267 conn[2] = node3; 00268 set_connectivity( entity_handle, conn, 3 ); \endcode */ 00269 virtual ErrorCode set_connectivity( const EntityHandle entity_handle, 00270 EntityHandle* connect, 00271 const int num_connect ); 00272 00273 //! get the adjacencies associated with a set of entities 00274 /** \param from_entities vector of EntityHandle to get adjacencies of. 00275 \param to_dimension Dimension of desired adjacency information. 00276 \param adj_entities Vector in which adjacent EntityHandles are returned. 00277 \param operation_type enum of INTERSECT or UNION. Defines whether to take 00278 the intersection or union of the set of adjacencies recovered for the from_entities. 00279 00280 The adjacent entities in vector <em>adjacencies</em> are not in any particular 00281 order. 00282 00283 Example: \code 00284 // get the set of edges that are adjacent to all entities in the from_entities list 00285 std::vector<EntityHandle> from_entities = {hex1, hex2}; 00286 std::vector<EntityHandle> adjacencies; 00287 get_adjacencies( from_entities, MB_1D_ENTITY, adjacencies ); 00288 \endcode */ 00289 00290 virtual ErrorCode get_adjacencies( const EntityHandle* from_entities, 00291 const int num_entities, 00292 const int to_dimension, 00293 const bool create_if_missing, 00294 std::vector< EntityHandle >& adj_entities, 00295 const int operation_type = Interface::INTERSECT ); 00296 00297 virtual ErrorCode get_adjacencies( const EntityHandle* from_entities, 00298 const int num_entities, 00299 const int to_dimension, 00300 const bool create_if_missing, 00301 Range& adj_entities, 00302 const int operation_type = Interface::INTERSECT ); 00303 00304 virtual ErrorCode get_adjacencies( const Range& from_entities, 00305 const int to_dimension, 00306 const bool create_if_missing, 00307 Range& adj_entities, 00308 const int operation_type = Interface::INTERSECT ); 00309 00310 /**\brief Get a ptr to adjacency lists 00311 * Get a pointer to adjacency lists. These lists are std::vector<EntityHandle>, which are 00312 * pointed to by adjs[i]. Adjacencies are not guaranteed to be in order of increasing 00313 * dimension. Only a const version of this function is given, because adjacency data is managed 00314 * more carefully in MOAB and should be treated as read-only by applications. If adjacencies 00315 * have not yet been initialized, adjs_ptr will be NULL (i.e. adjs_ptr == NULL). There may also 00316 * be NULL entries for individual entities, i.e. adjs_ptr[i] == NULL. \param iter Iterator to 00317 * beginning of entity range desired \param end End iterator for which adjacencies are requested 00318 * \param adjs_ptr Pointer to pointer to const std::vector<EntityHandle>; each member of that 00319 * array is the vector of adjacencies for this entity \param count Number of entities in the 00320 * contiguous chunk starting from *iter 00321 */ 00322 ErrorCode adjacencies_iterate( Range::const_iterator iter, 00323 Range::const_iterator end, 00324 const std::vector< EntityHandle >**& adjs_ptr, 00325 int& count ); 00326 00327 /**\brief Get all vertices for input entities 00328 * 00329 * Special case of get_adjacencies where to_dimension == 0 00330 * and operation_type == Interface::UNION. 00331 *\Note This is not a variation of get_connectivity because 00332 * the behavior is different for polyhedra. 00333 */ 00334 virtual ErrorCode get_vertices( const Range& from_entities, Range& vertices ); 00335 00336 //! Adds adjacencies 00337 /** \param from_handle entities 00338 \param both_ways add the adjacency information to both the 00339 to_handle and and the from_from :handle 00340 00341 Example: \code 00342 */ 00343 virtual ErrorCode add_adjacencies( const EntityHandle entity_handle, 00344 const EntityHandle* adjacencies, 00345 const int num_handles, 00346 bool both_ways ); 00347 00348 //! Adds adjacencies; same as vector-based, but with range instead 00349 virtual ErrorCode add_adjacencies( const EntityHandle entity_handle, Range& adjacencies, bool both_ways ); 00350 00351 //! Removes adjacencies 00352 /** \param handle EntityHandle to get adjacencies of. 00353 00354 Example: \code 00355 */ 00356 virtual ErrorCode remove_adjacencies( const EntityHandle entity_handle, 00357 const EntityHandle* adjacencies, 00358 const int num_handles ); 00359 00360 //! Retrieves all entities in the database of given dimension. 00361 /** \param dimension Dimension of entities desired. 00362 \param entities Range in which entities of dimension <em>dimension</em> are returned. 00363 00364 Example: \code 00365 int dimension = 2; 00366 Range entities; 00367 get_entities_by_dimension( dimension, entities ); //get 2D EntityHandles in the database 00368 \endcode */ 00369 virtual ErrorCode get_entities_by_dimension( const EntityHandle meshset, 00370 const int dimension, 00371 Range& entities, 00372 const bool recursive = false ) const; 00373 00374 virtual ErrorCode get_entities_by_dimension( const EntityHandle meshset, 00375 const int dimension, 00376 std::vector< EntityHandle >& entities, 00377 const bool recursive = false ) const; 00378 00379 //! Retrieves all entities in the data base of given type. 00380 /** \param type EntityType of entities desired (ie, MBHEX, MBEDGE, MBTRI, etc ) 00381 \param entities Range in which entities of EntityType <em>type</em> are returned. 00382 00383 Example: \code 00384 EntityType type = MBTET; 00385 Range entities; 00386 get_entities_by_dimension( type, entities ); //get MBTET type EntityHandles in the database 00387 \endcode */ 00388 virtual ErrorCode get_entities_by_type( const EntityHandle meshset, 00389 const EntityType type, 00390 Range& entities, 00391 const bool recursive = false ) const; 00392 00393 //! Retrieves all entities in the data base of given type. 00394 /** \param type EntityType of entities desired (ie, MBHEX, MBEDGE, MBTRI, etc ) 00395 \param entities Range in which entities of EntityType <em>type</em> are returned. 00396 00397 Example: \code 00398 EntityType type = MBTET; 00399 Range entities; 00400 get_entities_by_dimension( type, entities ); //get MBTET type EntityHandles in the database 00401 \endcode */ 00402 virtual ErrorCode get_entities_by_type( const EntityHandle meshset, 00403 const EntityType type, 00404 std::vector< EntityHandle >& entities, 00405 const bool recursive = false ) const; 00406 00407 virtual ErrorCode get_entities_by_type_and_tag( const EntityHandle meshset, 00408 const EntityType type, 00409 const Tag* tags, 00410 const void* const* values, 00411 const int num_tags, 00412 Range& entities, 00413 const int condition = Interface::INTERSECT, 00414 const bool recursive = false ) const; 00415 00416 //! Retrieves all entities in the data base 00417 /** \param entities Range in which entities of EntityType <em>type</em> are returned. 00418 00419 Example: \code 00420 Range entities; 00421 get_entities( entities ); //get MBTET type EntityHandles in the database 00422 \endcode */ 00423 virtual ErrorCode get_entities_by_handle( const EntityHandle meshset, 00424 Range& entities, 00425 const bool recursive = false ) const; 00426 00427 //! Retrieves all entities in the data base 00428 /** \param entities Range in which entities of EntityType <em>type</em> are returned. 00429 00430 Example: \code 00431 Range entities; 00432 get_entities( entities ); //get MBTET type EntityHandles in the database 00433 \endcode */ 00434 virtual ErrorCode get_entities_by_handle( const EntityHandle meshset, 00435 std::vector< EntityHandle >& entities, 00436 const bool recursive = false ) const; 00437 00438 //! Retrieves all entities in the database of given dimension. 00439 /** \param dimension Dimension of entities desired. 00440 \param entities Range in which entities of dimension <em>dimension</em> are returned. 00441 00442 Example: \code 00443 int dimension = 2; 00444 Range entities; 00445 get_entities_by_dimension( dimension, entities ); //get 2D EntityHandles in the database 00446 \endcode */ 00447 virtual ErrorCode get_number_entities_by_dimension( const EntityHandle meshset, 00448 const int dimension, 00449 int& number, 00450 const bool recursive = false ) const; 00451 00452 //! Retrieves all entities in the data base of given type. 00453 /** \param type EntityType of entities desired (ie, MBHEX, MBEDGE, MBTRI, etc ) 00454 \param entities Range in which entities of EntityType <em>type</em> are returned. 00455 00456 Example: \code 00457 EntityType type = MBTET; 00458 Range entities; 00459 get_entities_by_dimension( type, entities ); //get MBTET type EntityHandles in the database 00460 \endcode */ 00461 virtual ErrorCode get_number_entities_by_type( const EntityHandle meshset, 00462 const EntityType type, 00463 int& num_entities, 00464 const bool recursive = false ) const; 00465 00466 virtual ErrorCode get_number_entities_by_type_and_tag( const EntityHandle meshset, 00467 const EntityType type, 00468 const Tag* tag_handles, 00469 const void* const* values, 00470 const int num_tags, 00471 int& num_entities, 00472 const int condition = Interface::INTERSECT, 00473 const bool recursive = false ) const; 00474 00475 //! Retrieves all entities in the data base 00476 /** \param entities Range in which entities of EntityType <em>type</em> are returned. 00477 00478 Example: \code 00479 Range entities; 00480 get_entities( entities ); //get MBTET type EntityHandles in the database 00481 \endcode */ 00482 virtual ErrorCode get_number_entities_by_handle( const EntityHandle meshset, 00483 int& num_entities, 00484 const bool recursive = false ) const; 00485 00486 //! Creates an element based on the type and connectivity. 00487 /** If connectivity vector is not correct for EntityType <em>type</em> (ie, a vector with 00488 3 vertices is passed in to make an MBQUAD), the function returns MB_FAILURE. 00489 \param type Type of element to create. (MBTET, MBTRI, MBKNIFE, etc.) 00490 \param connectivity Vector containing connectivity of element to create. 00491 \param handle EntityHandle representing the newly created element in the database. 00492 00493 Example: \code 00494 EntityType type = MBQUAD; 00495 std::vector<EntityHandle> connectivity(4); 00496 quad_conn[0] = vertex0; 00497 quad_conn[1] = vertex1; 00498 quad_conn[2] = vertex2; 00499 quad_conn[3] = vertex3; 00500 EntityHandle element_handle; 00501 create_element( type, connectivity, element_handle ); \endcode */ 00502 virtual ErrorCode create_element( const EntityType type, 00503 const EntityHandle* connectivity, 00504 const int num_nodes, 00505 EntityHandle& element_handle ); 00506 00507 //! Creates a vertex based on coordinates. 00508 /** 00509 \param coordinates Array that has 3 doubles in it. 00510 \param entity_handle EntityHandle representing the newly created vertex in the database. 00511 00512 Example: \code 00513 double *coordinates = double[3]; 00514 coordinates[0] = 1.034; 00515 coordinates[1] = 23.23; 00516 coordinates[2] = -0.432; 00517 EntityHandle entity_handle = 0; 00518 create_vertex( coordinates, entity_handle ); \endcode */ 00519 virtual ErrorCode create_vertex( const double coords[3], EntityHandle& entity_handle ); 00520 00521 //! Create a set of vertices with the specified coordinates 00522 /** 00523 \param coordinates Array that has 3*n doubles in it. 00524 \param nverts Number of vertices to create 00525 \param entity_handles Range passed back with new vertex handles 00526 */ 00527 virtual ErrorCode create_vertices( const double* coordinates, const int nverts, Range& entity_handles ); 00528 00529 //! merges two entities 00530 virtual ErrorCode merge_entities( EntityHandle entity_to_keep, 00531 EntityHandle entity_to_remove, 00532 bool auto_merge, 00533 bool delete_removed_entity ); 00534 00535 //! Removes entities in a vector from the data base. 00536 /** If any of the entities are contained in any meshsets, it is removed from those meshsets 00537 which were created with MESHSET_TRACK_OWNER option bit set. Tags for <em>entity<\em> are 00538 removed as part of this function. 00539 \param entities 1d vector of entities to delete 00540 \param num_entities Number of entities in 1d vector 00541 */ 00542 virtual ErrorCode delete_entities( const EntityHandle* entities, const int num_entities ); 00543 00544 //! Removes entities in a range from the data base. 00545 /** If any of the entities are contained in any meshsets, it is removed from those meshsets 00546 which were created with MESHSET_TRACK_OWNER option bit set. Tags for <em>entity<\em> are 00547 removed as part of this function. 00548 \param entities Range of entities to delete 00549 */ 00550 virtual ErrorCode delete_entities( const Range& range ); 00551 00552 virtual ErrorCode list_entities( const Range& temp_range ) const; 00553 00554 virtual ErrorCode list_entities( const EntityHandle* entities, const int num_entities ) const; 00555 00556 virtual ErrorCode list_entity( const EntityHandle entity ) const; 00557 00558 typedef unsigned long long type_memstorage; 00559 00560 //! function object for recieving events from MB of higher order nodes 00561 //! added to entities 00562 class HONodeAddedRemoved 00563 { 00564 public: 00565 HONodeAddedRemoved() {} 00566 virtual ~HONodeAddedRemoved() {} 00567 //! node_added called when a node was added to an element's connectivity array 00568 //! note: connectivity array of element may be incomplete (corner nodes will exist always) 00569 virtual void node_added( EntityHandle node, EntityHandle element ); 00570 virtual void node_removed( EntityHandle node ); 00571 }; 00572 00573 virtual ErrorCode convert_entities( const EntityHandle meshset, 00574 const bool mid_side, 00575 const bool mid_face, 00576 const bool mid_volume, 00577 Interface::HONodeAddedRemoved* function_object = 0 ); 00578 00579 //! function to get the side number given two elements; returns 00580 //! MB_FAILURE if child not related to parent; does *not* create adjacencies 00581 //! between parent and child 00582 virtual ErrorCode side_number( const EntityHandle parent, 00583 const EntityHandle child, 00584 int& sd_number, 00585 int& sense, 00586 int& offset ) const; 00587 00588 //! given an entity and the connectivity and type of one of its subfacets, find the 00589 //! high order node on that subfacet, if any 00590 virtual ErrorCode high_order_node( const EntityHandle parent_handle, 00591 const EntityHandle* subfacet_conn, 00592 const EntityType subfacet_type, 00593 EntityHandle& hon ) const; 00594 00595 //! given an entity and a target dimension & side number, get that entity 00596 virtual ErrorCode side_element( const EntityHandle source_entity, 00597 const int dim, 00598 const int sd_number, 00599 EntityHandle& target_entity ) const; 00600 00601 //-------------------------Tag Stuff-------------------------------------// 00602 00603 /**\brief Get a tag handle, possibly creating the tag 00604 * 00605 * Get a handle used to associate application-defined values 00606 * with MOAB entities. If the tag does not already exist then 00607 * \c flags should contain exactly one of \c MB_TAG_SPARSE, 00608 * \c MB_TAG_DENSE, \c MB_TAG_MESH unless \c type is MB_TYPE_BIT, 00609 * which implies \c MB_TAG_BIT storage. 00610 * . 00611 *\param name The tag name 00612 *\param size Tag size as number of values of of data type per entity 00613 * (or number of bytes if \c MB_TAG_BYTES is passed in flags). If \c 00614 *MB_TAG_VARLEN is specified, this value is taken to be the size of the default value if one is 00615 *specified and is otherwise ignored. \param type The type of the data (used for IO) 00616 *\param tag_handle Output: the resulting tag handle. 00617 *\param flags Bitwise OR of values from \c TagType 00618 *\param default_value Optional default value for tag. 00619 *\param created Optional returned boolean indicating that the tag 00620 * was created. 00621 *\return - \c MB_ALREADY_ALLOCATED if tag exists and \c MB_TAG_EXCL is specified, or 00622 *default values do not match (and \c MB_TAG_ANY or \c MB_TAG_DFTOK not specified). 00623 * - \c MB_TAG_NOT_FOUND if tag does not exist and \c MB_TAG_CREAT is not 00624 *specified 00625 * - \c MB_INVALID_SIZE if tag value size is not a multiple of the size of the 00626 *data type (and \c MB_TAG_ANY not specified). 00627 * - \c MB_TYPE_OUT_OF_RANGE invalid or inconsistent parameter 00628 * - \c MB_VARIABLE_DATA_LENGTH if \c MB_TAG_VARLEN and \c default_value is non-null and 00629 * \c default_value_size is not specified. 00630 */ 00631 virtual ErrorCode tag_get_handle( const char* name, 00632 int size, 00633 DataType type, 00634 Tag& tag_handle, 00635 unsigned flags = 0, 00636 const void* default_value = 0, 00637 bool* created = 0 ); 00638 00639 /**\brief same as non-const version, except that TAG_CREAT flag is ignored. */ 00640 virtual ErrorCode tag_get_handle( const char* name, 00641 int size, 00642 DataType type, 00643 Tag& tag_handle, 00644 unsigned flags = 0, 00645 const void* default_value = 0 ) const; 00646 00647 //! Gets the tag name string of the tag_handle. 00648 /** \param tag_handle Tag you want the name of. 00649 \param tag_name Name string of <em>tag_handle</em>. 00650 00651 Example: \code 00652 Tag tag_handle = 0; 00653 std::string tag_name = "my_special_tag"; 00654 tag_get_name( tag_handle, tag_name ); //gets the Tag from the tag's name string 00655 \endcode */ 00656 virtual ErrorCode tag_get_name( const Tag tag_handle, std::string& tag_name ) const; 00657 00658 /**\brief Gets the tag handle corresponding to a name 00659 * 00660 * If a tag of that name does not exist, returns MB_TAG_NOT_FOUND 00661 * \param tag_name Name of the desired tag. 00662 * \param tag_handle Tag handle corresponding to <em>tag_name</em> 00663 */ 00664 virtual ErrorCode tag_get_handle( const char* tag_name, Tag& tag_handle ) const; 00665 00666 //! Get handles for all tags defined on this entity 00667 virtual ErrorCode tag_get_tags_on_entity( const EntityHandle entity, std::vector< Tag >& tag_handles ) const; 00668 00669 //! Get the size of the specified tag in bytes 00670 virtual ErrorCode tag_get_bytes( const Tag tag, int& tag_size ) const; 00671 00672 //! Get the array length of a tag 00673 virtual ErrorCode tag_get_length( const Tag tag, int& tag_size ) const; 00674 00675 //! Get the default value of the specified tag 00676 virtual ErrorCode tag_get_default_value( const Tag tag, void* def_val ) const; 00677 virtual ErrorCode tag_get_default_value( Tag tag, const void*& ptr, int& size ) const; 00678 00679 //! get type of tag (sparse, dense, etc.; 0 = dense, 1 = sparse, 2 = bit, 3 = mesh) 00680 virtual ErrorCode tag_get_type( const Tag, TagType& tag_type ) const; 00681 00682 /** \brief Get data type of tag. 00683 * 00684 * Get the type of the tag data. The tag is data is assumed to 00685 * be a vector of this type. If the tag data vetcor contains 00686 * more than one value, then the tag size must be a multiple of 00687 * the size of this type. 00688 * \param tag The tag 00689 * \param type The type of the specified tag (output). 00690 */ 00691 virtual ErrorCode tag_get_data_type( const Tag handle, DataType& type ) const; 00692 00693 //! get handles for all tags defined 00694 virtual ErrorCode tag_get_tags( std::vector< Tag >& tag_handles ) const; 00695 00696 virtual ErrorCode tag_get_data( const Tag tag_handle, 00697 const EntityHandle* entity_handles, 00698 const int num_entities, 00699 void* tag_data ) const; 00700 00701 virtual ErrorCode tag_get_data( const Tag tag_handle, const Range& entity_handles, void* tag_data ) const; 00702 00703 //! Sets the data of a given EntityHandle and Tag. 00704 /** If the <em>tag_handle</em> and the entity type of <em>entity_handle</em> are not 00705 compatible, data of <em>entity_handle</em> never existed and MB_FAILURE 00706 is returned. 00707 \param tag_handle Tag indicating what data is to be set. 00708 \param entity_handle EntityHandle on which to set tag's data. 00709 \param tag_data Data to set the <em>entity_handle</em>'s tag data to. 00710 00711 Example: \code 00712 int tag_data = 1004; 00713 tag_set_data( tag_handle, entity_handle, &tag_data ); \endcode */ 00714 virtual ErrorCode tag_set_data( Tag tag_handle, 00715 const EntityHandle* entity_handles, 00716 int num_entities, 00717 const void* tag_data ); 00718 00719 virtual ErrorCode tag_set_data( Tag tag_handle, const Range& entity_handles, const void* tag_data ); 00720 00721 /**\brief Get pointers to tag data 00722 * 00723 * For a tag, get the values for a list of passed entity handles. 00724 *\note This function may not be used for bit tags. 00725 *\param tag_handle The tag 00726 *\param entity_handles An array of entity handles for which to retreive tag values. 00727 *\param num_entities The length of the 'entity_handles' array. 00728 *\param tag_data An array of 'const void*'. Array must be at least 00729 * 'num_entitities' long. Array is populated (output) 00730 * with pointers to the internal storage for the 00731 * tag value corresponding to each entity handle. 00732 *\param tag_sizes The length of each tag value. Optional for 00733 * fixed-length tags. Required for variable-length tags. 00734 */ 00735 virtual ErrorCode tag_get_by_ptr( const Tag tag_handle, 00736 const EntityHandle* entity_handles, 00737 int num_entities, 00738 const void** tag_data, 00739 int* tag_sizes = 0 ) const; 00740 00741 /**\brief Get pointers to tag data 00742 * 00743 * For a tag, get the values for a list of passed entity handles. 00744 *\note This function may not be used for bit tags. 00745 *\param tag_handle The tag 00746 *\param entity_handles The entity handles for which to retreive tag values. 00747 *\param tag_data An array of 'const void*'. Array is populated (output) 00748 * with pointers to the internal storage for the 00749 * tag value corresponding to each entity handle. 00750 *\param tag_sizes The length of each tag value. Optional for 00751 * fixed-length tags. Required for variable-length tags. 00752 */ 00753 virtual ErrorCode tag_get_by_ptr( const Tag tag_handle, 00754 const Range& entity_handles, 00755 const void** tag_data, 00756 int* tag_sizes = 0 ) const; 00757 00758 /**\brief Set tag data given an array of pointers to tag values. 00759 * 00760 * For a tag, set the values for a list of passed entity handles. 00761 *\note This function may not be used for bit tags. 00762 *\param tag_handle The tag 00763 *\param entity_handles An array of entity handles for which to set tag values. 00764 *\param num_entities The length of the 'entity_handles' array. 00765 *\param tag_data An array of 'const void*'. Array must be at least 00766 * 'num_entitities' long. Array is expected to 00767 * contain pointers to tag values for the corresponding 00768 * EntityHandle in 'entity_handles'. 00769 *\param tag_sizes The length of each tag value. Optional for 00770 * fixed-length tags. Required for variable-length tags. 00771 */ 00772 virtual ErrorCode tag_set_by_ptr( Tag tag_handle, 00773 const EntityHandle* entity_handles, 00774 int num_entities, 00775 void const* const* tag_data, 00776 const int* tag_sizes = 0 ); 00777 00778 /**\brief Set tag data given an array of pointers to tag values. 00779 * 00780 * For a tag, set the values for a list of passed entity handles. 00781 *\note This function may not be used for bit tags. 00782 *\param tag_handle The tag 00783 *\param entity_handles The entity handles for which to set tag values. 00784 *\param tag_data An array of 'const void*'. Array is expected to 00785 * contain pointers to tag values for the corresponding 00786 * EntityHandle in 'entity_handles'. 00787 *\param tag_sizes The length of each tag value. Optional for 00788 * fixed-length tags. Required for variable-length tags. 00789 */ 00790 virtual ErrorCode tag_set_by_ptr( Tag tag_handle, 00791 const Range& entity_handles, 00792 void const* const* tag_data, 00793 const int* tag_sizes = 0 ); 00794 00795 /**\brief Set tag data given value. 00796 * 00797 * For a tag, set the values for a list of passed entity handles to 00798 * the same, specified value. 00799 * 00800 *\param tag_handle The tag 00801 *\param entity_handles The entity handles for which to set tag values. 00802 *\param tag_data A pointer to the tag value. 00803 *\param tag_sizes For variable-length tags, the lenght of the 00804 * tag value. This argument will be ignored for 00805 * fixed-length tags. 00806 */ 00807 virtual ErrorCode tag_clear_data( Tag tag_handle, 00808 const Range& entity_handles, 00809 const void* tag_data, 00810 int tag_size = 0 ); 00811 00812 /**\brief Set tag data given value. 00813 * 00814 * For a tag, set the values for a list of passed entity handles to 00815 * the same, specified value. 00816 * 00817 *\param tag_handle The tag 00818 *\param entity_handles The entity handles for which to set tag values. 00819 *\param tag_data A pointer to the tag value. 00820 *\param tag_sizes For variable-length tags, the lenght of the 00821 * tag value. This argument will be ignored for 00822 * fixed-length tags. 00823 */ 00824 virtual ErrorCode tag_clear_data( Tag tag_handle, 00825 const EntityHandle* entity_handles, 00826 int num_entities, 00827 const void* tag_data, 00828 int tag_size = 0 ); 00829 00830 //! Delete the data of a vector of entity handles and sparse tag 00831 /** Delete the data of a tag on a vector of entity handles. Only sparse tag data are deleted 00832 with this function; dense tags are deleted by deleting the tag itself using tag_delete. 00833 \param tag_handle Handle of the (sparse) tag being deleted from entity 00834 \param entity_handles 1d vector of entity handles from which the tag is being deleted 00835 \param num_handles Number of entity handles in 1d vector 00836 */ 00837 virtual ErrorCode tag_delete_data( Tag tag_handle, const EntityHandle* entity_handles, int num_entities ); 00838 00839 //! Delete the data of a range of entity handles and sparse tag 00840 /** Delete the data of a tag on a range of entity handles. Only sparse tag data are deleted 00841 with this function; dense tags are deleted by deleting the tag itself using tag_delete. 00842 \param tag_handle Handle of the (sparse) tag being deleted from entity 00843 \param entity_range Range of entities from which the tag is being deleted 00844 */ 00845 virtual ErrorCode tag_delete_data( Tag tag_handle, const Range& entity_handles ); 00846 00847 //! Removes the tag from the database and deletes all of its associated data. 00848 virtual ErrorCode tag_delete( Tag tag_handle ); 00849 00850 /**\brief Access tag data via direct pointer into contiguous blocks 00851 * 00852 * Iteratively obtain direct access to contiguous blocks of tag 00853 * storage. This function cannot be used with bit tags because 00854 * of the compressed bit storage. This function cannot be used 00855 * with variable length tags because it does not provide a mechanism 00856 * to determine the length of the value for each entity. This 00857 * function may be used with sparse tags, but if it is used, it 00858 * will return data for a single entity at a time. 00859 * 00860 *\param tag_handle The handle of the tag for which to access data 00861 *\param iter The first entity for which to return data. 00862 *\param end One past the last entity for which data is desired. 00863 *\param count The number of entities for which data was returned 00864 *\param data_ptr Output: pointer to tag storage. 00865 *\param allocate If true, space for this tag will be allocated, if not it wont 00866 * 00867 *\Note If this function is called for entities for which no tag value 00868 * has been set, but for which a default value exists, it will 00869 * force the allocation of explicit storage for each such entity 00870 * even though MOAB would normally not explicitly store tag values 00871 * for such entities. 00872 * 00873 *\Example: 00874 *\code 00875 * Range ents; // range to iterate over 00876 * Tag tag; // tag for which to access data 00877 * int bytes; 00878 * ErrorCode err = mb.tag_get_size( tag, bytes ); 00879 * if (err) { ... } 00880 * 00881 * ... 00882 * Range::iterator iter = ents.begin(); 00883 * while (iter != ents.end()) { 00884 * int count; 00885 * // get contiguous block of tag dat 00886 * void* ptr; 00887 * err = mb.tag_iterate( tag, iter, ents.end(), count, ptr ); 00888 * if (err) { ... } 00889 * // do something with tag data 00890 * process_Data( ptr, count ); 00891 * // advance to next block of data 00892 * iter += count; 00893 * } 00894 *\endcode 00895 */ 00896 virtual ErrorCode tag_iterate( Tag tag_handle, 00897 Range::const_iterator iter, 00898 Range::const_iterator end, 00899 int& count, 00900 void*& data_ptr, 00901 bool allocate = true ); 00902 00903 //! creates a mesh set 00904 virtual ErrorCode create_meshset( const unsigned int options, EntityHandle& ms_handle, int start_id = 0 ); 00905 00906 //! Empty a vector of mesh set 00907 /** Empty a mesh set. 00908 \param ms_handles 1d vector of handles of sets being emptied 00909 \param num_meshsets Number of entities in 1d vector 00910 */ 00911 virtual ErrorCode clear_meshset( const EntityHandle* ms_handles, const int num_meshsets ); 00912 00913 //! Empty a range of mesh set 00914 /** Empty a mesh set. 00915 \param ms_handles Range of handles of sets being emptied 00916 */ 00917 virtual ErrorCode clear_meshset( const Range& ms_handles ); 00918 00919 //! get the options of a mesh set 00920 virtual ErrorCode get_meshset_options( const EntityHandle ms_handle, unsigned int& options ) const; 00921 00922 //! set the options of a mesh set 00923 virtual ErrorCode set_meshset_options( const EntityHandle ms_handle, const unsigned int options ); 00924 00925 //! subtracts meshset2 from meshset1 - modifies meshset1 00926 virtual ErrorCode subtract_meshset( EntityHandle meshset1, const EntityHandle meshset2 ); 00927 00928 //! intersects meshset2 with meshset1 - modifies meshset1 00929 virtual ErrorCode intersect_meshset( EntityHandle meshset1, const EntityHandle meshset2 ); 00930 00931 //! unites meshset2 with meshset1 - modifies meshset1 00932 virtual ErrorCode unite_meshset( EntityHandle meshset1, const EntityHandle meshset2 ); 00933 00934 //! add entities to meshset 00935 virtual ErrorCode add_entities( EntityHandle meshset, const Range& entities ); 00936 00937 //! add entities to meshset 00938 virtual ErrorCode add_entities( EntityHandle meshset, const EntityHandle* entities, const int num_entities ); 00939 00940 //! remove entities from meshset 00941 virtual ErrorCode remove_entities( EntityHandle meshset, const Range& entities ); 00942 00943 //! remove entities from meshset 00944 virtual ErrorCode remove_entities( EntityHandle meshset, const EntityHandle* entities, const int num_entities ); 00945 00946 //! return true if all entities are contained in set 00947 virtual bool contains_entities( EntityHandle meshset, 00948 const EntityHandle* entities, 00949 int num_entities, 00950 const int operation_type = Interface::INTERSECT ); 00951 00952 //! replace entities 00953 virtual ErrorCode replace_entities( EntityHandle meshset, 00954 const EntityHandle* old_entities, 00955 const EntityHandle* new_entities, 00956 int num_entities ); 00957 00958 //------MeshSet Parent/Child functions------ 00959 00960 //! get parent meshsets 00961 virtual ErrorCode get_parent_meshsets( const EntityHandle meshset, 00962 std::vector< EntityHandle >& parents, 00963 const int num_hops = 1 ) const; 00964 00965 //! get parent meshsets 00966 virtual ErrorCode get_parent_meshsets( const EntityHandle meshset, Range& parents, const int num_hops = 1 ) const; 00967 00968 //! get child meshsets 00969 virtual ErrorCode get_child_meshsets( const EntityHandle meshset, 00970 std::vector< EntityHandle >& children, 00971 const int num_hops = 1 ) const; 00972 00973 //! get child meshsets 00974 virtual ErrorCode get_child_meshsets( const EntityHandle meshset, Range& children, const int num_hops = 1 ) const; 00975 00976 //! get contained meshsets 00977 virtual ErrorCode get_contained_meshsets( const EntityHandle meshset, 00978 std::vector< EntityHandle >& children, 00979 const int num_hops = 1 ) const; 00980 00981 //! get contained meshsets 00982 virtual ErrorCode get_contained_meshsets( const EntityHandle meshset, 00983 Range& children, 00984 const int num_hops = 1 ) const; 00985 00986 //! gets number of parent meshsets 00987 virtual ErrorCode num_parent_meshsets( const EntityHandle meshset, int* number, const int num_hops = 1 ) const; 00988 00989 //! gets number of child meshsets 00990 virtual ErrorCode num_child_meshsets( const EntityHandle meshset, int* number, const int num_hops = 1 ) const; 00991 00992 //! gets number of contained meshsets 00993 virtual ErrorCode num_contained_meshsets( const EntityHandle meshset, int* number, const int num_hops = 1 ) const; 00994 00995 //! add a parent meshset 00996 virtual ErrorCode add_parent_meshset( EntityHandle meshset, const EntityHandle parent_meshset ); 00997 00998 //! add parent meshsets 00999 virtual ErrorCode add_parent_meshsets( EntityHandle meshset, const EntityHandle* parents, int count ); 01000 01001 //! add a child meshset 01002 virtual ErrorCode add_child_meshset( EntityHandle meshset, const EntityHandle child_meshset ); 01003 01004 //! add parent meshsets 01005 virtual ErrorCode add_child_meshsets( EntityHandle meshset, const EntityHandle* children, int count ); 01006 01007 //! adds 'parent' to child's parent list and adds 'child' to parent's child list 01008 virtual ErrorCode add_parent_child( EntityHandle parent, EntityHandle child ); 01009 01010 //! removes 'parent' to child's parent list and removes 'child' to parent's child list 01011 virtual ErrorCode remove_parent_child( EntityHandle parent, EntityHandle child ); 01012 01013 //! remove parent meshset 01014 virtual ErrorCode remove_parent_meshset( EntityHandle meshset, const EntityHandle parent_meshset ); 01015 01016 //! remove child meshset 01017 virtual ErrorCode remove_child_meshset( EntityHandle meshset, const EntityHandle child_meshset ); 01018 01019 // ************************ tag information *************** 01020 01021 //! return various specific tag handles 01022 Tag material_tag(); 01023 Tag neumannBC_tag(); 01024 Tag dirichletBC_tag(); 01025 Tag globalId_tag(); 01026 Tag geom_dimension_tag(); 01027 01028 //! get/set the number of nodes 01029 // int total_num_nodes() const; 01030 // void total_num_nodes(const int val); 01031 01032 //! get/set the number of elements 01033 // int total_num_elements() const; 01034 // void total_num_elements(const int val); 01035 01036 // ************************ structured sequence information *************** 01037 01038 //! return a reference to the sequence manager 01039 SequenceManager* sequence_manager() 01040 { 01041 return sequenceManager; 01042 } 01043 const SequenceManager* sequence_manager() const 01044 { 01045 return sequenceManager; 01046 } 01047 01048 /// create structured sequence 01049 ErrorCode create_scd_sequence( const HomCoord& coord_min, 01050 const HomCoord& coord_max, 01051 EntityType type, 01052 EntityID start_id_hint, 01053 EntityHandle& first_handle_out, 01054 EntitySequence*& sequence_out ); 01055 01056 ErrorCode add_vsequence( EntitySequence* vert_seq, 01057 EntitySequence* elem_seq, 01058 const HomCoord& p1, 01059 const HomCoord& q1, 01060 const HomCoord& p2, 01061 const HomCoord& q2, 01062 const HomCoord& p3, 01063 const HomCoord& q3, 01064 bool bb_input = false, 01065 const HomCoord* bb_min = NULL, 01066 const HomCoord* bb_max = NULL ); 01067 01068 //! return the a_entity_factory pointer 01069 AEntityFactory* a_entity_factory() 01070 { 01071 return aEntityFactory; 01072 } 01073 const AEntityFactory* a_entity_factory() const 01074 { 01075 return aEntityFactory; 01076 } 01077 01078 #ifdef MOAB_HAVE_AHF 01079 HalfFacetRep* a_half_facet_rep() 01080 { 01081 return ahfRep; 01082 } 01083 const HalfFacetRep* a_half_facet_rep() const 01084 { 01085 return ahfRep; 01086 } 01087 #endif 01088 01089 //! return set of registered IO tools 01090 ReaderWriterSet* reader_writer_set() 01091 { 01092 return readerWriterSet; 01093 } 01094 01095 //-----------------MeshSet Interface Functions------------------// 01096 01097 void print( const EntityHandle handle, const char* prefix, bool first_call = true ) const; 01098 01099 ErrorCode print_entity_tags( std::string indent_prefix, const EntityHandle handle, TagType tp ) const; 01100 01101 virtual ErrorCode get_last_error( std::string& info ) const; 01102 01103 virtual std::string get_error_string( const ErrorCode code ) const; 01104 01105 //! check all adjacencies for consistency 01106 ErrorCode check_adjacencies(); 01107 01108 //! check some adjacencies for consistency 01109 ErrorCode check_adjacencies( const EntityHandle* ents, int num_ents ); 01110 01111 //! return whether the input handle is valid or not 01112 bool is_valid( const EntityHandle this_ent ) const; 01113 01114 /** \brief Create an iterator over the set 01115 * Create a new iterator that iterates over entities with the specified type or dimension. 01116 * Only one of ent_type or dim can be set; use dim=-1 or ent_type=MBMAXTYPE for the other. 01117 * Iterators for list-type (ordered) sets are stable over set modification, unless entity 01118 * removed or deleted is the one at the current position of the iterator. If the check_valid 01119 * parameter is passed as true, entities are checked for validity before being passed back by 01120 * get_next_entities function (checking entity validity can have a non-negligible cost). 01121 * 01122 * Iterators returned by this function can be deleted using the normal C++ delete function. 01123 * After creating the iterator through this function, further interactions are through methods 01124 * on the SetIterator class. 01125 * \param meshset The entity set associated with this iterator (use 0 for whole instance) 01126 * \param ent_type Entity type associated with this iterator 01127 * \param ent_dim Dimension associated with this iterator 01128 * \param chunk_size Chunk size of the iterator 01129 * \param check_valid If true, entities are checked for validity before being returned 01130 */ 01131 virtual ErrorCode create_set_iterator( EntityHandle meshset, 01132 EntityType ent_type, 01133 int ent_dim, 01134 int chunk_size, 01135 bool check_valid, 01136 SetIterator*& set_iter ); 01137 01138 /** \brief Remove the set iterator from the instance's list 01139 * \param set_iter Set iterator being removed 01140 */ 01141 ErrorCode remove_set_iterator( SetIterator* set_iter ); 01142 01143 /** \brief Get all set iterators associated with the set passed in 01144 * \param meshset Meshset for which iterators are requested 01145 * \param set_iters Set iterators for the set 01146 */ 01147 ErrorCode get_set_iterators( EntityHandle meshset, std::vector< SetIterator* >& set_iters ); 01148 01149 //-----------------Memory Functions------------------// 01150 01151 /**\brief Calculate amount of memory used to store MOAB data 01152 * 01153 * This function calculates the amount of memory used to store 01154 * MOAB data. 01155 * 01156 * There are two possible values for each catagory of memory use. 01157 * The exact value and the amortized value. The exact value is the 01158 * amount of memory used to store the data for the specified entities. 01159 * The amortized value includes the exact value and an amortized 01160 * estimate of the memory consumed in overhead for storing the values 01161 * (indexing structures, access structures, etc.) 01162 * 01163 * Note: If ent_array is NULL, the total memory used by MOAB for storing 01164 * data will be returned in the address pointed to by 01165 * total_amortized_storage, if total_amortized_storage is not NULL. 01166 * 01167 *\param ent_array Array of entities for which to estimate the memory use. 01168 * If NULL, estimate is done for all entities. 01169 *\param num_ents The length of ent_array. Not used if ent_rray is NULL. 01170 *\param total_(amortized_)storage The sum of the memory entity, adjacency, and all tag storage. 01171 *\param (amortized_)entity_storage The storage for the entity definitions 01172 * (connectivity arrays for elements, coordinates for vertices, 01173 * list storage within sets, etc.) 01174 *\param (amortized_)adjacency_storage The storage for adjacency data. 01175 *\param tag_array An array of tags for which to calculate the memory use. 01176 *\param num_tags The lenght of tag_array 01177 *\param (amortized_)tag_storage If tag_array is not NULL, then one value 01178 * for each tag specifying the memory used for storing that 01179 * tag. If tag_array is NULL and this value is not, the 01180 * location at which to store the total memory used for 01181 * all tags. 01182 */ 01183 void estimated_memory_use( const EntityHandle* ent_array = 0, 01184 unsigned long num_ents = 0, 01185 unsigned long long* total_storage = 0, 01186 unsigned long long* total_amortized_storage = 0, 01187 unsigned long long* entity_storage = 0, 01188 unsigned long long* amortized_entity_storage = 0, 01189 unsigned long long* adjacency_storage = 0, 01190 unsigned long long* amortized_adjacency_storage = 0, 01191 const Tag* tag_array = 0, 01192 unsigned num_tags = 0, 01193 unsigned long long* tag_storage = 0, 01194 unsigned long long* amortized_tag_storage = 0 ); 01195 01196 /**\brief Calculate amount of memory used to store MOAB data 01197 * 01198 * This function calculates the amount of memory used to store 01199 * MOAB data. 01200 * 01201 * There are two possible values for each catagory of memory use. 01202 * The exact value and the amortized value. The exact value is the 01203 * amount of memory used to store the data for the specified entities. 01204 * The amortized value includes the exact value and an amortized 01205 * estimate of the memory consumed in overhead for storing the values 01206 * (indexing structures, access structures, etc.) 01207 * 01208 *\param ents Entities for which to estimate the memory use. 01209 *\param total_(amortized_)storage The sum of the memory entity, adjacency, and all tag storage. 01210 *\param (amortized_)entity_storage The storage for the entity definitions 01211 * (connectivity arrays for elements, coordinates for vertices, 01212 * list storage within sets, etc.) 01213 *\param (amortized_)adjacency_storage The storage for adjacency data. 01214 *\param tag_array An array of tags for which to calculate the memory use. 01215 *\param num_tags The lenght of tag_array 01216 *\param (amortized_)tag_storage If tag_array is not NULL, then one value 01217 * for each tag specifying the memory used for storing that 01218 * tag. If tag_array is NULL and this value is not, the 01219 * location at which to store the total memory used for 01220 * all tags. 01221 */ 01222 void estimated_memory_use( const Range& ents, 01223 unsigned long long* total_storage = 0, 01224 unsigned long long* total_amortized_storage = 0, 01225 unsigned long long* entity_storage = 0, 01226 unsigned long long* amortized_entity_storage = 0, 01227 unsigned long long* adjacency_storage = 0, 01228 unsigned long long* amortized_adjacency_storage = 0, 01229 const Tag* tag_array = 0, 01230 unsigned num_tags = 0, 01231 unsigned long long* tag_storage = 0, 01232 unsigned long long* amortized_tag_storage = 0 ); 01233 01234 void print_database() const; 01235 01236 /** \name Sequence Option controllers */ 01237 01238 /**@{*/ 01239 01240 /** \brief Interface to control memory allocation for sequences 01241 * Provide a factor that controls the size of the sequence that gets allocated. 01242 * This is typically useful in the parallel setting when a-priori, the number of ghost entities 01243 * and the memory required for them within the same sequence as the owned entities are unknown. 01244 * The default factor is 1.0 but this can be appropriately updated at runtime so that we do not 01245 * have broken sequences. 01246 */ 01247 virtual double get_sequence_multiplier() const; 01248 01249 /** \brief Interface to control memory allocation for sequences 01250 * Provide a factor that controls the size of the sequence that gets allocated. 01251 * This is typically useful in the parallel setting when a-priori, the number of ghost entities 01252 * and the memory required for them within the same sequence as the owned entities are unknown. 01253 * The default factor is 1.0 but this can be appropriately updated at runtime so that we do not 01254 * have broken sequences. 01255 * 01256 * \param meshset User specified multiplier (should be greater than 1.0) 01257 */ 01258 virtual void set_sequence_multiplier( double factor ); 01259 01260 /**@}*/ 01261 01262 private: 01263 /**\brief Do not allow copying */ 01264 Core( const Core& copy ); 01265 /**\brief Do not allow copying */ 01266 Core& operator=( const Core& copy ); 01267 01268 void estimated_memory_use_internal( const Range* ents, 01269 unsigned long long* total_storage, 01270 unsigned long long* total_amortized_storage, 01271 unsigned long long* entity_storage, 01272 unsigned long long* amortized_entity_storage, 01273 unsigned long long* adjacency_storage, 01274 unsigned long long* amortized_adjacency_storage, 01275 const Tag* tag_array, 01276 unsigned num_tags, 01277 unsigned long long* tag_storage, 01278 unsigned long long* amortized_tag_storage ); 01279 01280 //! database init and de-init routines 01281 ErrorCode initialize(); 01282 void deinitialize(); 01283 01284 //! return the entity set representing the whole mesh 01285 EntityHandle get_root_set(); 01286 01287 //!\brief Clean up after a file reader returns failure. 01288 //! 01289 //! Delete all entities not contained in initial_entities 01290 //! and all tags not contained in initial_tags. 01291 void clean_up_failed_read( const Range& initial_ents, std::vector< Tag > initial_tags ); 01292 01293 // other interfaces for MB 01294 WriteUtil* mMBWriteUtil; 01295 ReadUtil* mMBReadUtil; 01296 ScdInterface* scdInterface; 01297 01298 //! store the total number of elements defined in this interface 01299 // int totalNumElements; 01300 01301 //! store the total number of nodes defined in this interface 01302 // int totalNumNodes; 01303 01304 //! the overall geometric dimension of this mesh 01305 int geometricDimension; 01306 01307 Tag materialTag; 01308 Tag neumannBCTag; 01309 Tag dirichletBCTag; 01310 Tag geomDimensionTag; 01311 Tag globalIdTag; 01312 01313 //! tag server for this interface 01314 std::list< TagInfo* > tagList; 01315 inline bool valid_tag_handle( const TagInfo* t ) const 01316 { 01317 return std::find( tagList.begin(), tagList.end(), t ) != tagList.end(); 01318 } 01319 01320 SequenceManager* sequenceManager; 01321 01322 AEntityFactory* aEntityFactory; 01323 01324 ReaderWriterSet* readerWriterSet; 01325 01326 Error* mError; 01327 bool mpiFinalize; 01328 int writeMPELog; 01329 bool initErrorHandlerInCore; 01330 01331 //! list of iterators 01332 std::vector< SetIterator* > setIterators; 01333 01334 #ifdef MOAB_HAVE_AHF 01335 HalfFacetRep* ahfRep; 01336 bool mesh_modified; 01337 #endif 01338 }; 01339 01340 } // namespace moab 01341 01342 #endif // MOAB_IMPL_GENERAL_HPP