MOAB: Mesh Oriented datABase  (version 5.4.1)
moabUG.h
Go to the documentation of this file.
00001 /*! \page userguide User's Guide 
00002  
00003   \subpage team 
00004  
00005   \subpage contents
00006  
00007   \subpage figures
00008  
00009   \subpage tables
00010  
00011   \subpage differences
00012 
00013   \subpage building
00014 
00015   \page team MOAB team members
00016  <h2>The MOAB Team, including: </h2>
00017  
00018  - Vijay S. Mahadevan (Argonne National Lab)
00019  - Timothy J. Tautges (Siemens, Univ Wisconsin-Madison)
00020  - Iulian Grindeanu (Argonne National Lab) 
00021  - Rajeev Jain (Argonne National Lab)
00022  - Danqing Wu  (Argonne National Lab)
00023  - Navamita Ray (Los Alamos National Lab)
00024  - Jane Hu (Univ Wisconsin-Madison)
00025  - Paul Wilson (Univ Wisconsin-Madison)
00026  - Patrick Shriwise (Argonne National Lab)
00027  - Anthony Scopatz (The University of South Carolina)
00028 
00029 
00030  <h2>Emeritus members:</h2>
00031  
00032  - Jason A. Kraftcheck
00033  - Brandon M. Smith
00034  - Hong-Jun Kim
00035  - Jim Porter
00036  - Xiabing Xu
00037  
00038   \page contents Table of Contents
00039  
00040   \ref introduction  
00041 
00042   \ref interface     
00043 
00044     \ref twoone    
00045 
00046     \ref twotwo     
00047 
00048     \ref twothree       
00049 
00050     \ref twofour   
00051 
00052   \ref api     
00053 
00054   \ref services      
00055 
00056     \ref fourone    
00057 
00058     \ref fourtwo   
00059 
00060     \ref fourthree  
00061 
00062     \ref fourfour      
00063 
00064     \ref fourfive    
00065 
00066     \ref foursix
00067 
00068     \ref fourseven
00069 
00070     \ref foureight
00071 
00072   \ref parallel      
00073 
00074     \ref fiveone    
00075 
00076     \ref fivetwo     
00077 
00078     \ref fivethree    
00079 
00080     \ref fivefour      
00081 
00082   \ref applications   
00083 
00084   \ref implementation         
00085 
00086   \ref pymoab
00087 
00088   \ref representation     
00089 
00090   \ref element    
00091 
00092     \ref tenone  
00093 
00094     \ref tentwo        
00095 
00096     \ref tenthree      
00097 
00098   \ref performance   
00099 
00100   \ref error-handling
00101 
00102   \ref conclusions    
00103 
00104   \ref references 
00105 
00106   \section introduction 1.Introduction
00107 
00108 In scientific computing, systems of partial differential equations (PDEs) are solved on computers.  One of the most widely used methods to solve PDEs numerically is to solve over discrete neighborhoods or “elements” of the domain.  Popular discretization methods include Finite Difference (FD), Finite Element (FE), and Finite Volume (FV).  These methods require the decomposition of the domain into a discretized representation, which is referred to as a “mesh”.  The mesh is one of the fundamental types of data linking the various tools in the analysis process (mesh generation, analysis, visualization, etc.).  Thus, the representation of mesh data and operations on those data play a very important role in PDE-based simulations.
00109  
00110 MOAB is a component for representing and evaluating mesh data. It is part of the Scalable Interfaces for Geometry and Mesh based Applications (SIGMA) toolkit [1]. MOAB can store structured and unstructured mesh, consisting of elements in the finite element “zoo”, along with polygons and polyhedra.  The functional interface to MOAB is simple, consisting of only four fundamental data types.  This data is quite powerful, allowing the representation of most types of metadata commonly found on the mesh.  Internally, MOAB uses array-based storage for fine-grained data, which in many cases provides more efficient access, especially for large portions of mesh and associated data.  MOAB is optimized for efficiency in space and time, based on access to mesh in chunks rather than through individual entities, while also versatile enough to support individual entity access. MOAB is also unique in that it maintains a completely parallel view of the mesh database so that queries to traverse the cells and manipulation of the grid in memory is scalable in parallel.
00111 
00112 The MOAB data model consists of the following four fundamental types: mesh interface instance, mesh entities (vertex, edge, tri, etc.), sets, and tags.  Entities are addressed through handles rather than pointers, to allow the underlying representation of an entity to change without changing the handle to that entity.  Sets are arbitrary groupings of mesh entities and other sets.  Sets also support parent/child relationships as a relation distinct from sets containing other sets.  The directed graph provided by set parent/child relationships is useful for embedding graphs whose nodes include collections of mesh entities; this approach has been used to represent a wide variety of application-specific data, including geometric model topology, processor partitions, and various types of search trees.  Tags are named data which can be assigned to the mesh as a whole, individual entities, or sets.  Tags are a mechanism for attaching data to individual entities, and sets are a mechanism for describing relations between entities; the combination of these two mechanisms is a powerful yet simple interface for representing metadata or application-specific data.
00113 
00114 Various mesh-related tools are provided with MOAB or can be used directly with MOAB.  These tools can be used for mesh format translation (mbconvert), mesh skinning (Skinner class), solution transfer between meshes (MBCoupler tool), ray tracing and other geometric searches (OrientedBoxTreeTool, AdaptiveKDTree), visualization (vtkMOABReader tool), mesh optimization/smoothing (Mesquite interfaces), and relation between mesh and geometric models (the separately-packed Lasso tool).  These tools are described later in this document.
00115 
00116 MOAB is written in the C++ programming language, with applications interacting with MOAB mostly through its moab::Interface class.  All of the MOAB functions and classes are isolated in and accessed through the moab namespace<sup>1</sup>. The remainder of this report gives class and function names without the “moab::” namespace qualification; unless otherwise noted, the namespace qualifier should be added to all class and function names referenced here.  MOAB also implements the iMesh interface, which is specified in C but can be called directly from other languages.  Almost all of the functionality in MOAB can be accessed through the iMesh interface.  MOAB is developed and supported on the Linux and MacOS operating systems, as well as various HPC operating systems.  MOAB can be used on parallel computing systems as well, including both clusters and high-end parallel systems like IBM BG/P and Cray systems.  MOAB is released under a standard LGPL open source software license.
00117 
00118 MOAB is used in several ways in various applications.  MOAB serves as the underlying mesh data representation in several scientific computing applications [2], [3], [4], [5], [6].  MOAB can also be used as a mesh format translator, using readers and writers included in MOAB.  MOAB has also been used as a bridge to couple results in multi-physics analysis and to link these applications with other mesh services [7] in order to consistently solve problems at scale.
00119 
00120 The remainder of this guide is organized as follows.  Section 2, “Getting Started”, provides a few simple examples of using MOAB to perform simple tasks on a mesh.  Section 3 discusses the MOAB data model in more detail, including some aspects of the implementation.  Section 4 summarizes the MOAB function API.  Section 5 describes some of the tools included with MOAB, and the implementation of mesh readers/writers for MOAB.  Section 6 describes how to build MOAB-based applications.  Section 7 contains a brief description of MOAB’s relation to the iMesh mesh interface.  Sections 8 and 9 discuss MOAB's representations of structured and spectral element meshes, respectively.  Section 10 gives helpful hints for accessing MOAB in an efficient manner from applications.  Section 11 gives a conclusion and future plans for MOAB development.  Section 12 gives references cited in this report.
00121 
00122 Several other sources of information about MOAB may also be of interest to readers.  Meta-data conventions define how sets and /or tags are used together to represent various commonly-used simulation constructs in MOAB [8]. MOAB also has an active mailing list [9] to provide guidance and encourage discussions about usage and development of algorithms using MOAB. A separate mailing list [10] for MOAB-related release announcements and future guidance may also be subscribed by users.  Potential users are encouraged to interact with the MOAB team using these mailing lists, and watch the SIGMA website [1] for updates.
00123 
00124 <sup>1</sup> Non-namespaced names are also provided for backward compatibility, with the “MB” prefix added to the class or variable name.
00125 
00126  \ref contents
00127 
00128  \section interface 2.MOAB Data Model
00129 The MOAB data model describes the basic types used in MOAB and the language used to communicate that data to applications.  This chapter describes that data model, along with some of the reasons for some of the design choices in MOAB.
00130 
00131  \ref contents
00132 
00133  \subsection twoone 2.1. MOAB Interface 
00134 
00135 MOAB is written in C++.  The primary interface with applications is through member functions of the abstract base class Interface.  The MOAB library is created by instantiating Core, which implements the Interface API.  Multiple instances of MOAB can exist concurrently in the same application; mesh entities are not shared between these instances<sup>2</sup>.  MOAB is most easily viewed as a database of mesh objects accessed through the instance.  No other assumptions explicitly made about the nature of the mesh stored there; for example, there is no fundamental requirement that elements fill space or do not overlap each other geometrically.
00136  
00137 <sup>2</sup> One exception to this statement is when the parallel interface to MOAB is used; in this case, entity sharing between instances is handled explicitly using message passing.  This is described in more detail in Section 5 of this document.
00138 
00139  \ref contents
00140 
00141  \subsection twotwo 2.2. Mesh Entities
00142 MOAB represents the following topological mesh entities: vertex, edge, triangle, quadrilateral, polygon, tetrahedron, pyramid, prism, knife, hexahedron, polyhedron.  MOAB uses the EntityType enumeration to refer to these entity types (see Table 1).  This enumeration has several special characteristics, chosen intentionally: the types begin with vertex, entity types are grouped by topological dimension, with lower-dimensional entities appearing before higher dimensions; the enumeration includes an entity type for sets (described in the next section); and MBMAXTYPE is included at the end of this enumeration, and can be used to terminate loops over type.  In addition to these defined values, the an increment operator (++) is defined such that variables of type EntityType can be used as iterators in loops.
00143 MOAB refers to entities using “handles”.  Handles are implemented as long integer data types, with the four highest-order bits used to store the entity type (mesh vertex, edge, tri, etc.) and the remaining bits storing the entity id.  This scheme is convenient for applications because:
00144 - Handles sort lexicographically by type and dimension; this can be useful for grouping and iterating over entities by type.
00145 - The type of an entity is indicated by the handle itself, without needing to call a function.
00146 - Entities allocated in sequence will typically have contiguous handles; this characteristic can be used to efficiently store and operate on large lists of handles.
00147 .
00148 
00149 This handle implementation is exposed to applications intentionally, because of optimizations that it enables, and is unlikely to change in future versions.
00150 
00151   \subsection tableone Table 1: Values defined for the EntityType enumerated type.
00152 <table border="1">
00153 <tr>
00154 <td>MBVERTEX = 0</td>
00155 <td>MBPRISM</td>
00156 </tr>
00157 <tr>
00158 <td>MBEDGE</td>
00159 <td>MBKNIFE</td>
00160 </tr>
00161 <tr>
00162 <td>MBTRI</td>
00163 <td>MBHEX</td>
00164 </tr>
00165 <tr>
00166 <td>MBQUAD</td>
00167 <td>MBPOLYHEDRON</td>
00168 </tr>
00169 <tr>
00170 <td>MBPOLYGON</td>
00171 <td>MBENTITYSET</td>
00172 </tr>
00173 <tr>
00174 <td>MBTET</td>
00175 <td>MBMAXTYPE</td>
00176 </tr>
00177 <tr>
00178 <td>MBPYRAMID</td>
00179 <td></td>
00180 </tr>
00181 </table>
00182 
00183 MOAB defines a special class for storing lists of entity handles, named Range.  This class stores handles as a series of (start_handle, end_handle) subrange tuples.  If a list of handles has large contiguous ranges, it can be represented in almost constant size using Range.  Since entities are typically created in groups, e.g. during mesh generation or file import, a high degree of contiguity in handle space is typical.  Range provides an interface similar to C++ STL containers like std::vector, containing iterator data types and functions for initializing and iterating over entity handles stored in the range.  Range also provides functions for efficient Boolean operations like subtraction and intersection.  Most API functions in MOAB come in both range-based and vector-based variants.  By definition, a list of entities stored in an Range is always sorted, and can contain a given entity handle only once.  Range cannot store the handle 0 (zero).
00184 
00185 Typical usage of an Range object would look like:
00186 
00187 \code
00188 using namespace moab;
00189    int my_function(Range &from_range) {
00190           int num_in_range = from_range.size();
00191           Range to_range;
00192           Range::iterator rit;
00193     for (rit = from_range.begin(); rit != from_range.end(); ++rit) {
00194             EntityHandle this_ent = *rit;
00195             to_range.insert(this_ent);
00196           }
00197         }
00198 \endcode
00199 
00200 Here, the range is iterated similar to how std::vector is iterated.
00201 
00202   \ref contents
00203 
00204  \subsection adjacencies 2.2.1. Adjacencies & AEntities 
00205 
00206 The term adjacencies is used to refer to those entities topologically connected to a given entity, e.g. the faces bounded by a given edge or the vertices bounding a given region.  The same term is used for both higher-dimensional (or bounded) and lower-dimensional (or bounding) adjacent entities.  MOAB provides functions for querying adjacent entities by target dimension, using the same functions for higher- and lower-dimension adjacencies.  By default, MOAB stores the minimum data necessary to recover adjacencies between entities.  When a mesh is initially loaded into MOAB, only entity-vertex (i.e. “downward”) adjacencies are stored, in the form of entity connectivity.  When “upward” adjacencies are requested for the first time, e.g. from vertices to regions, MOAB stores all vertex-entity adjacencies explicitly, for all entities in the mesh.  Non-vertex entity to entity adjacencies are never stored, unless explicitly requested by the application.
00207 
00208 In its most fundamental form, a mesh need only be represented by its vertices and the entities of maximal topological dimension.  For example, a hexahedral mesh can be represented as the connectivity of the hex elements and the vertices forming the hexes.  Edges and faces in a 3D mesh need not be explicitly represented.  We refer to such entities as “AEntities”, where ‘A’ refers to “Auxiliary”, “Ancillary”, and a number of other words mostly beginning with ‘A’.  Individual AEntities are created only when requested by applications, either using mesh modification functions or by requesting adjacencies with a special “create if missing” flag passed as “true”.  This reduces the overall memory usage when representing large meshes.  Note entities must be explicitly represented before they can be assigned tag values or added to entity sets (described in following Sections).
00209 
00210 \ref contents
00211 
00212  \subsection twothree 2.3. Entity Sets
00213 Entity sets are also known as "mesh sets", or when the context is clear, not to be confused with std::set, just "sets". Entity sets are used to store arbitrary collections of entities and other sets.  Sets are used for a variety of things in mesh-based applications, from the set of entities discretizing a given geometric model entity to the entities partitioned to a specific processor in a parallel finite element application.  MOAB entity sets can also store parent/child relations with other entity sets, with these relations distinct from contains relations.  Parent/child relations are useful for building directed graphs with graph nodes representing collections of mesh entities; this construct can be used, for example, to represent an interface of mesh faces shared by two distinct collections of mesh regions.  MOAB also defines one special set, the “root set” or the interface itself; all entities are part of this set by definition.  Defining a root set allows the use of a single set of MOAB API functions to query entities in the overall mesh as well as its subsets.
00214 
00215 MOAB entity sets can be one of two distinct types: list-type entity sets preserve the order in which entities are added to the set, and can store a given entity handle multiple times in the same set; set-type sets are always ordered by handle, regardless of the order of addition to the set, and can store a given entity handle only once.  This characteristic is assigned when the set is created, and cannot be changed during the set’s lifetime.
00216 
00217 MOAB provides the option to track or not track entities in a set.  When entities (and sets) are deleted by other operations in MOAB, they will also be removed from containing sets for which tracking has been enabled.  This behavior is assigned when the set is created, and cannot be changed during the set’s lifetime.  The cost of turning tracking on for a given set is sizeof(EntityHandle) for each entity added to the set; MOAB stores containing sets in the same list which stores adjacencies to other entities.
00218 
00219 Using an entity set looks like the following:
00220 \code
00221 using namespace moab;
00222 // load a file using MOAB, putting the loaded mesh into a file set
00223 EntityHandle file_set;
00224 ErrorCode rval = moab->create_meshset(MESHSET_SET, file_set);
00225 rval = moab->load_file(“fname.vtk”, &file_set);
00226 Range set_ents;
00227 // get all the 3D entities in the set
00228 rval = moab->get_entities_by_dimension(file_set, 3, set_ents);
00229 \endcode
00230 
00231 Entity sets are often used in conjunction with tags (described in the next section), and provide a powerful mechanism to store a variety of meta-data with meshes.
00232 
00233 \ref contents
00234 
00235  \subsection twofour 2.4. Tags 
00236 
00237 Applications of a mesh database often need to attach data to mesh entities.  The types of attached data are often not known at compile time, and can vary across individual entities and entity types.  MOAB refers to this attached data as a “tag”.  Tags can be thought of loosely as a variable, which can be given a distinct value for individual entities, entity sets, or for the interface itself.  A tag is referenced using a handle, similarly to how entities are referenced in MOAB.  Each MOAB tag has the following characteristics, which can be queried through the MOAB interface:
00238 - Name
00239 - Size (in bytes)
00240 - Storage type
00241 - Data type (integer, double, opaque, entity handle)
00242 - Handle
00243 .
00244 
00245 The storage type determines how tag values are stored on entities.  
00246 
00247 - Dense: Dense tag values are stored in arrays which match arrays of contiguous entity handles.  Dense tags are more efficient in both storage and memory if large numbers of entities are assigned the same tag.  Storage for a given dense tag is not allocated until a tag value is set on an entity; memory for a given dense tag is allocated for all entities in a given sequence at the same time.
00248 - Sparse: Sparse tags are stored as a list of (entity handle, tag value) tuples, one list per sparse tag, sorted by entity handle.
00249 - Bit: Bit tags are stored similarly to dense tags, but with special handling to allow allocation in bit-size amounts per entity.
00250 .
00251 
00252 MOAB also supports variable-length tags, which can have a different length for each entity they are assigned to.  Variable length tags are stored similarly to sparse tags.
00253 
00254 The data type of a tag can either be one understood at compile time (integer, double, entity handle), in which case the tag value can be saved and restored properly to/from files and between computers of different architecture (MOAB provides a native HDF5-based save/restore format for this purpose; see Section 4.6).  The opaque data type is used for character strings, or for allocating “raw memory” for use by applications (e.g. for storage application-defined structures or other abstract data types).  These tags are saved and restored as raw memory, with no special handling for endian or precision differences.
00255 
00256 An application would use the following code to attach a double-precision tag to vertices in a mesh, e.g. to assign a temperature field to those vertices:
00257 
00258 \code
00259 using namespace moab;
00260 // load a file using MOAB and get the vertices
00261 ErrorCode rval = moab->load_file(“fname.vtk”);
00262 Range verts;
00263 rval = moab->get_entities_by_dimension(0, 0, verts);
00264 // create a tag called “TEMPERATURE”
00265 Tag temperature;
00266 double def_val = -1.0d-300, new_val = 273.0;
00267 rval = moab->tag_create(“TEMPERATURE”, sizeof(double), MB_TAG_DENSE, 
00268                         MB_TYPE_DOUBLE, temperature, &def_val);
00269 // assign a value to vertices
00270 for (Range::iterator vit = verts.begin(); 
00271      vit != verts.end(); ++vit)
00272   rval = moab->tag_set_data(temperature, &(*rit), 1, &new_val);
00273 
00274 \endcode
00275 
00276 The semantic meaning of a tag is determined by applications using it.  However, to promote interoperability between applications, there are a number of tag names reserved by MOAB which are intended to be used by convention.  Mesh readers and writers in MOAB use these tag conventions, and applications can use them as well to access the same data. Ref. [1] maintains an up-to-date list of conventions for meta-data usage in MOAB.
00277 
00278   \ref contents
00279 
00280   \section api 3.MOAB API Design Philosophy and Summary
00281 
00282 This section describes the design philosophy behind MOAB, and summarizes the functions, data types and enumerated variables in the MOAB API.  A complete description of the MOAB API is available in online documentation in the MOAB distribution [11].
00283 
00284 MOAB is designed to operate efficiently on collections of entities.  Entities are often created or referenced in groups (e.g. the mesh faces discretizing a given geometric face, the 3D elements read from a file), with those groups having some form of temporal or spatial locality.  The interface provides special mechanisms for reading data directly into the native storage used in MOAB, and for writing large collections of entities directly from that storage, to avoid data copies.  MOAB applications structured to take advantage of that locality will typically operate more efficiently.
00285 
00286 MOAB has been designed to maximize the flexibility of mesh data which can be represented.  There is no explicit constraint on the geometric structure of meshes represented in MOAB, or on the connectivity between elements.  In particular, MOAB allows the representation of multiple entities with the same exact connectivity; however, in these cases, explicit adjacencies must be used to distinguish adjacencies with AEntities bounding such entities.
00287 
00288 The number of vertices used to represent a given topological entity can vary, depending on analysis needs; this is often the case in FEA.  For example, applications often use “quadratic” or 10-vertex tetrahedral, with vertices at edge midpoints as well as corners.  MOAB does not distinguish these variants by entity type, referring to all variants as “tetrahedra”.  The number of vertices for a given entity is used to distinguish the variants, with canonical numbering conventions used to determine placement of the vertices [12].  This is similar to how such variations are represented in the Exodus [13] and Patran [14] file formats.  In practice, we find that this simplifies coding in applications, since in many cases the handling of entities depends only on the number of corner vertices in the element.  Some MOAB API functions provide a flag which determines whether corner or all vertices are requested.
00289 
00290 The MOAB API is designed to balance complexity and ease of use.  This balance is evident in the following general design characteristics:
00291 
00292 - Entity lists: Lists of entities are passed to and from MOAB in a variety of forms.  Lists output from MOAB are passed as either STL vector or Range data types.  Either of these constructs may be more efficient in both time and memory, depending on the semantics of the data being requested.  Input lists are passed as either Range’s, or as a pointer to EntityHandle and a size.  The latter allows the same function to be used when passing individual entities, without requiring construction of an otherwise unneeded STL vector.
00293 - Entity sets: Most query functions accept an entity set as input.  Applications can pass zero to indicate a request for the whole interface.  Note that this convention applies only to query functions; attempts to add or subtract entities to/from the interface using set-based modification functions, or to add parents or children to the interface set, will fail.  Allowing specification of the interface set in this manner avoids the need for a separate set of API functions to query the database as a whole.
00294 - Implicit Booleans in output lists: A number of query functions in MOAB allow specification of a Boolean operation (Interface::INTERSECT or Interface::UNION).  This operation is applied to the results of the query, often eliminating the need for code the application would need to otherwise implement.  For example, to find the set of vertices shared by a collection of quadrilaterals, the application would pass that list of quadrilaterals to a request for vertex adjacencies, with Interface::INTERSECT passed for the Boolean flag.  The list of vertices returned would be the same as if the application called that function for each individual entity, and computed the intersection of the results over all the quadrilaterals.  Applications may also input non-empty lists to store the results, in which case the intersection is also performed with entities already in the list.  In many cases, this allows optimizations in both time and memory inside the MOAB implementation. 
00295 .
00296 
00297 Since these objectives are at odds with each other, tradeoffs had to be made between them.  Some specific issues that came up are:
00298 
00299 - Using ranges: Where possible, entities can be referenced using either ranges (which allow efficient storage of long lists) or STL vectors (which allow list order to be preserved), in both input and output arguments.
00300 - Entities in sets: Accessing the entities in a set is done using the same functions which access entities in the entire mesh.  The whole mesh is referenced by specifying a set handle of zero<sup>3</sup>.
00301 - Entity vectors on input: Functions which could normally take a single entity as input are specified to take a vector of handles instead.  Single entities are specified by taking the address of that entity handle and specifying a list length of one.  This minimizes the number of functions, while preserving the ability to input single entities.<sup>4</sup>
00302 .
00303 
00304 Table 2 lists basic data types and enumerated variables defined and used by MOAB.  Values of the ErrorCode enumeration are returned from most MOAB functions, and can be compared to those listed in Appendix [ref-appendix].
00305 
00306 MOAB uses several pre-defined tag names to define data commonly found in various mesh-based analyses.  Ref. [8] and the \href{metadata.html}{meta-data conventions guide} provide an up-to-date decription as new conventions emerge for using sets and tags in MOAB applications.
00307 
00308   \subsection tabletwo Table 2: Basic data types and enums defined in MOAB.
00309 
00310 <table border="1">
00311 <tr>
00312 <th>Enum / Type</th>
00313 <th>Description</th>
00314 </tr>
00315 <tr>
00316 <td>ErrorCode</td>
00317 <td>Specific error codes returned from MOAB</td>
00318 </tr>
00319 <tr>
00320 <td>EntityHandle</td>
00321 <td>Type used to represent entity handles</td>
00322 </tr>
00323 <tr>
00324 <td>Tag</td>
00325 <td>Type used to represent tag handles</td>
00326 </tr>
00327 <tr>
00328 <td>TagType</td>
00329 <td>Type used to represent tag storage type</td>
00330 </tr>
00331 <tr>
00332 <td>DataType</td>
00333 <td>Type used to represent tag data type</td>
00334 </tr>
00335 </table>
00336 
00337 Table 3 lists the various groups of functions that comprise the MOAB API.  This is listed here strictly as a reference to the various types of functionality supported by MOAB; for a more detailed description of the scope and syntax of the MOAB API, see the online documentation [11].
00338 
00339   \subsection tablethree Table 3: Groups of functions in MOAB API.  See Ref. [11] for more details.
00340 
00341 <table border="1">
00342 <tr>
00343 <th>Function group</th>
00344 <th>Examples</th>
00345 <th>Description</th>
00346 </tr>
00347 <tr>
00348 <td>Constructor, destructor, interface</td>
00349 <td>Interface, ~Core, query_interface</td>
00350 <td>Construct/destroy interface; get pointer to read/write interface</td>
00351 </tr>
00352 <tr>
00353 <td>Entity query</td>
00354 <td>get_entities_by_dimension, get_entities_by_handle</td>
00355 <td>Get entities by dimension, type, etc.</td>
00356 </tr>
00357 <tr>
00358 <td>Adjacencies</td>
00359 <td>get_adjacencies, set_adjacencies, add_adjacencies</td>
00360 <td>Get topologically adjacent entities; set or add explicit adjacencies</td>
00361 </tr>
00362 <tr>
00363 <td>Vertex coordinates</td>
00364 <td>get_coords, set_coords</td>
00365 <td>Get/set vertex coordinates</td>
00366 </tr>
00367 <tr>
00368 <td>Connectivity</td>
00369 <td>get_connectivity, set_connectivity</td>
00370 <td>Get/set connectivity of non-vertex entities</td>
00371 </tr>
00372 <tr>
00373 <td>Sets</td>
00374 <td>create_meshset, add_entities, add_parent_child</td>
00375 <td>Create and work with entity sets</td>
00376 </tr>
00377 <tr>
00378 <td>Tags</td>
00379 <td>tag_get_data, tag_create</td>
00380 <td>Create, read, write tag data</td>
00381 </tr>
00382 <tr>
00383 <td>Handles</td>
00384 <td>type_from_handle, id_from_handle</td>
00385 <td>Go between handles and types/ids</td>
00386 </tr>
00387 <tr>
00388 <td>File handling</td>
00389 <td>load_mesh, save_mesh</td>
00390 <td>Read/write mesh files</td>
00391 </tr>
00392 <tr>
00393 <td>Geometric dimension</td>
00394 <td>get_dimension, set_dimension</td>
00395 <td>Get/set geometric dimension of mesh</td>
00396 </tr>
00397 <tr>
00398 <td>Mesh modification</td>
00399 <td>create_vertex, delete_entity</td>
00400 <td>Create or delete mesh entities</td>
00401 </tr>
00402 <tr>
00403 <td>Information</td>
00404 <td>list_entities, get_last_error</td>
00405 <td>Get or print certain information</td>
00406 </tr>
00407 <tr>
00408 <td>High-order nodes</td>
00409 <td>high_order_node</td>
00410 <td>Get information on high-order nodes</td>
00411 </tr>
00412 <tr>
00413 <td>Canonical numbering</td>
00414 <td>side_number</td>
00415 <td>Get canonical numbering information</td>
00416 </tr>
00417 </table>
00418 
00419 <sup>3</sup>In iMesh, the whole mesh is specified by a special entity set handle, referred to as the “root set”.
00420 
00421 <sup>4</sup>Note that STL vectors of entity handles can be input in this manner by using &vector[0] and vector.size() for the 1d vector address and size, respectively.
00422 
00423  \ref contents
00424 
00425  \section services 4.Related Mesh Services
00426 
00427 A number of mesh-based services are often used in conjunction with a mesh library.  For example, parallel applications often need to visualize the mesh and associated data.  Other services, like spatial interpolation or finding the faces on the “skin” of a 3D mesh, can be implemented more efficiently using knowledge of specific data structures in MOAB.  Several of these services provided with MOAB are described in this chapter.
00428 
00429  \ref contents
00430 
00431   \subsection fourone 4.1. Visualization
00432 
00433 Visualization is one of the most common needs associated with meshes.  The primary tool used to visualize MOAB meshes is VisIt [15].  Users can download a VisIt version that has the MOAB plugin pre-compiled (versions >= 2.12), and then select the MOAB format to directly visualize native files (default extension h5m) in HDF5 format. The visualization plugin has also been written with scalability in mind and can be launched on several processors to perform visualization and analysis of large datasets resulting from applications.
00434 
00435 There are capabilities in VisIt for viewing and manipulation of tag data and some types of entity sets.  Dense tag data is visualized using the same mechanisms used to view other field data in VisIt, e.g. using a pseudocolor plot; Material sets, neumann and dirichlet sets, along with parallel partition sets can be visualized using the subset capability. 
00436 
00437  \ref contents
00438 
00439   \subsection fourtwo 4.2. Parallel Decomposition
00440 
00441 To support parallel simulation, applications often need to partition a mesh into parts, designed to balance the load and minimize communication between sets.  MOAB includes the mbpart tool for this purpose, constructed on the well-known Zoltan partitioning library [16] and Metis [17].  After computing the partition using Zoltan or Metis, MOAB stores the partition as either tags on individual entities in the partition, or as tagged sets, one set per part.  Since a partition often exhibits locality similar to how the entities were created, storing it as sets (based on Range’s) is often more memory-efficient than an entity tag-based representation.  Figure below shows a couple of partitioned meshes computed with mbpart with -z option for Zoltan and visualized in VisIt.
00442 
00443 
00444  \image html vis_part.png
00445 
00446 
00447 
00448  \ref contents
00449 
00450   \subsection fourthree 4.3. Skinner
00451 
00452 An operation commonly applied to mesh is to compute the outermost “skin” bounding a contiguous block of elements.  This skin consists of elements of one fewer topological dimension, arranged in one or more topological balls on the boundary of the elements.  The Skinner tool computes the skin of a mesh in a memory-efficient manner.  Skinner uses knowledge about whether vertex-entity adjacencies and AEntities exist to minimize memory requirements and searching time required during the skinning process.  This skin can be provided as a single collection of entities, or as sets of entities distinguished by forward and reverse orientation with respect to higher-dimensional entities in the set being skinned.
00453 
00454 The following code fragment shows how Skinner can be used to compute the skin of a range of hex elements:
00455 
00456  \code
00457 using namespace moab;
00458 Range hexes, faces;
00459 ErrorCode rval = moab->get_entities_by_dimension(0, 3, hexes);
00460 Skinner myskinner(moab);
00461 bool verts_too = false;
00462 ErrorCode rval = myskinner.find_skin(hexes, verts_too, faces);
00463 \endcode
00464 
00465 Skinner can also skin a mesh based on geometric topology groupings imported with the mesh.  The geometric topology groupings contain information about the mesh “owned” by each of the entities in the geometric model, e.g. the model vertices, edges, etc.  Links between the mesh sets corresponding to those entities can be inferred directly from the mesh.  Skinning a mesh this way will typically be much faster than doing so on the actual mesh elements, because there is no need to create and destroy interior faces on the mesh.
00466 
00467  \ref contents
00468 
00469   \subsection fourfour 4.4. Tree Decompositions
00470 
00471 MOAB provides several mechanisms for spatial decomposition and searching in a mesh:
00472 
00473 - AdaptiveKDTree: Adaptive KD tree, a space-filling decomposition with axis-aligned splitting planes, enabling fast searching.
00474 - BSPTree: Binary Space Partition tree, with non-axis-aligned partitions, for fast spatial searches with slightly better memory efficiency than KD trees.
00475 - OrientedBoxTreeTool: Oriented Bounding Box tree hierarchy, useful for fast ray-tracing on collections of mesh facets.
00476 .
00477 
00478 These trees have various space and time searching efficiencies.  All are implemented based on entity sets and parent/child relations between those sets, allowing storage of a tree decomposition using MOAB’s native file storage mechanism (see Section 4.6.1).  MOAB’s entity set implementation is specialized for memory efficiency when representing binary trees.  Tree decompositions in MOAB have been used to implement fast ray tracing to support radiation transport [6], [18], solution coupling between meshes [2], [3], [7], and embedded boundary mesh generation [19].  MOAB also includes the DagMC tool to support several Monte Carlo radiation transport solvers by exposing queries related to ray tracing on the underlying geometry.
00479 
00480 The following code fragment shows very basic use of AdaptiveKDTree.  A range of entities is put in the tree; the leaf containing a given point is found, and the entities in that leaf are returned.
00481 \code
00482 using namespace moab;
00483 // create the adaptive kd tree from a range of tets
00484 std::ostringstream options;
00485 options << "MAX_PER_LEAF=1;SPLITS_PER_DIR=1;PLANE_SET=0;";
00486 FileOptions opts(options.str().c_str());
00487 EntityHandle root;
00488 ErrorCode rval;
00489 AdaptiveKDTree tool( &moab );
00490 rval = tool.build_tree( tets, &root, &opts);
00491 
00492 // get the overall bounding box corners
00493 const double * boxmax, * boxmin;
00494 boxmax = myTree.box_max();
00495 boxmin = myTree.box_min;
00496 
00497 // get the tree leaf containing point xyz, and the tets in that leaf
00498 AdaptiveKDTreeIter treeiter;
00499 rval = myTree.point_search(xyz, treeiter);
00500 Range leaf_tets;
00501 rval = moab->get_entities_by_dimension(treeiter.handle(), 3, 
00502                                        leaf_tets, false);
00503 \endcode
00504 
00505 More detailed examples of using the various tree decompositions in MOAB can be found in [ref-treeexamples].
00506 
00507  \ref contents
00508 
00509   \subsection fourfive 4.5. File Reader/Writer Interfaces
00510 
00511 Mesh readers and writers communicate mesh into/out of MOAB from/to disk files.  Reading a mesh often involves importing large sets of data, for example coordinates of all the nodes in the mesh.  Normally, this process would involve reading data from the file into a temporary data buffer, then copying data from there into its destination in MOAB.  To avoid the expense of copying data, MOAB has implemented a reader/writer interface that provides direct access to blocks of memory used to represent mesh.
00512 
00513 The reader interface, declared in ReadUtilIface, is used to request blocks of memory for storing coordinate positions and element connectivity.  The pointers returned from these functions point to the actual memory used to represent those data in MOAB.  Once data is written to that memory, no further copying is done.  This not only saves time, but it also eliminates the need to allocate a large memory buffer for intermediate storage of these data. 
00514 
00515 MOAB allocates memory for nodes and elements (and their corresponding dense tags) in chunks, to avoid frequent allocation/de-allocation of small chunks of memory.  The chunk size used depends on from where the mesh is being created, and can strongly affect the performance (and memory layout) of MOAB.  Since dense tags are allocated at the chunk size, this can also affect overall memory usage in cases where the mesh size is small but the number of dense tags or dense tag size is large.  When creating vertices and elements through the normal MOAB API, default chunk sizes defined in the SequenceManager class are used.  However, most of the file readers in MOAB allocate the exact amount of space necessary to represent the mesh being read.  There are also a few exceptions to this:
00516 
00517 - When compiled in parallel, this space is increased by a factor of 1.5, to allow subsequent creation of ghost vertices/elements in the same chunk as the original mesh.
00518 - The .cub file reader, which creates nodes and elements for individual geometric entities in separate calls, allocates using the default vertex/element sequence sizes, which are defined in the SequenceManager class in MOAB.
00519 .
00520 
00521 Applications calling the reader interface functions directly can specify the allocation chunk size as an optional parameter.
00522 
00523 The reader interface consists of the following functions:
00524 
00525 - get_node_coords: Given the number of vertices requested, the number of geometric dimensions, and a requested start id, allocates a block of vertex handles and returns pointers to coordinate arrays in memory, along with the actual start handle for that block of vertices.
00526 - get_element_connect: Given the number of elements requested, the number of vertices per element, the element type and the requested start id, allocates the block of elements, and returns a pointer to the connectivity array for those elements and the actual start handle for that block.  The number of vertices per element is necessary because those elements may include higher-order nodes, and MOAB stores these as part of the normal connectivity array.
00527 - update_adjacencies: This function takes the start handle for a block of elements and the connectivity of those elements, and updates adjacencies for those elements.  Which adjacencies are updated depends on the options set in AEntityFactory.
00528 .
00529 
00530 The following code fragment illustrates the use of ReadUtilIface to read a mesh directly into MOAB’s native representation.  This code assumes that connectivity is specified in terms of vertex indices, with vertex indices starting from 1.
00531 
00532 \code
00533 // get the read iface from moab
00534 ReadUtilIface *iface;
00535 ErrorCode rval = moab->query_interface(iface);
00536 
00537 // allocate a block of vertex handles and read xyz’s into them
00538 std::vector<double*> arrays;
00539 EntityHandle startv, *starth;
00540 rval = iface->get_node_coords(3, num_nodes, 0, startv, arrays);
00541 for (int i = 0; i < num_nodes; i++)
00542   infile >> arrays[0][i] >> arrays[1][i] >> arrays[2][i];
00543 
00544 // allocate block of hex handles and read connectivity into them
00545 rval = iface->get_element_connect(num_hexes, 8, MBHEX, 0, starth);
00546 for (int i = 0; i < 8*num_hexes; i++)
00547   infile >> starth[i];
00548 
00549 // change connectivity indices to vertex handles
00550 for (int i = 0; i < 8*num_hexes; i++)
00551   starth[i] += startv-1;
00552 \endcode
00553 
00554 The writer interface, declared in WriteUtilIface, provides functions that support writing vertex coordinates and element connectivity to storage locations input by the application.  Assembling these data is a common task for writing mesh, and can be non-trivial when exporting only subsets of a mesh.  The writer interface declares the following functions:
00555 
00556 - get_node_coords: Given already-allocated memory and the number of vertices and dimensions, and a range of vertices, this function writes vertex coordinates to that memory.  If a tag is input, that tag is also written with integer vertex ids, starting with 1, corresponding to the order the vertices appear in that sequence (these ids are used to write the connectivity array in the form of vertex indices).
00557 - get_element_connect: Given a range of elements and the tag holding vertex ids, and a pointer to memory, the connectivity of the specified elements are written to that memory, in terms of the indices referenced by the specified tag.  Again, the number of vertices per element is input, to allow the direct output of higher-order vertices.
00558 - gather_nodes_from_elements: Given a range of elements, this function returns the range of vertices used by those elements.  If a bit-type tag is input, vertices returned are also marked with 0x1 using that tag.  If no tag is input, the implementation of this function uses its own bit tag for marking, to avoid using an n2 algorithm for gathering vertices.
00559 - reorder: Given a permutation vector, this function reorders the connectivity for entities with specified type and number of vertices per entity to match that permutation.  This function is needed for writing connectivity into numbering systems other than that used internally in MOAB.
00560 .
00561 
00562 The following code fragment shows how to use WriteUtilIface to write the vertex coordinates and connectivity indices for a subset of entities.
00563 
00564 \code
00565 using namespace moab;
00566 // get the write iface from moab
00567 WriteUtilIface *iface;
00568 ErrorCode rval = moab->query_interface(iface);
00569 
00570 // get all hexes the model, and choose the first 10 of those
00571 Range tmp_hexes, hexes, verts;
00572 rval = moab->get_entities_by_type(0, MBHEX, tmp_hexes);
00573 for (int i = 0; i < 10; i++) hexes.insert(tmp_hexes[i]);
00574 rval = iface->gather_nodes_from_elements(hexes, 0, verts);
00575 
00576 // assign vertex ids
00577 iface->assign_ids(verts, 0, 1);
00578 
00579 // allocate space for coordinates & write them
00580 std::vector<double*> arrays(3);
00581 for (int i = 0; i < 3; i++) arrays[i] = new double[verts.size()];
00582 iface->get_node_coords(3, verts.size(), verts, 0, 1, arrays);
00583 
00584 // put connect’y in array, in the form of indices into vertex array
00585 std::vector<int> conn(8*hexes.size());
00586 iface->get_element_connect(hexes.size(), 8, 0, hexes, 0, 1, &conn[0]);
00587 \endcode
00588 
00589  \ref contents
00590 
00591   \subsection foursix 4.6. File Readers/Writers Packaged With MOAB
00592 
00593 MOAB has been designed to efficiently represent data and metadata commonly found in finite element mesh files.  Readers and writers are included with MOAB which import/export specific types of metadata in terms of MOAB sets and tags, as described earlier in this document.  The number of readers and writers in MOAB will probably grow over time, and so they are not enumerated here.  See the src/io/README file in the MOAB source distribution for a current list of supported formats.
00594 
00595 Because of its generic support for readers and writers, described in the previous section, MOAB is also a good environment for constructing new mesh readers and writers.  The ReadTemplate and WriteTemplate classes in src/io are useful starting points for constructing new file readers/writers; applications are encouraged to submit their own readers/writers for inclusion in MOAB’s contrib/io directory in the MOAB source. 
00596 
00597 The usefulness of a file reader/writer is determined not only by its ability to read and write nodes and elements, but also in its ability to store the various types of meta-data included with the typical mesh.  MOAB readers and writers are distinguished by their ability to preserve meta-data in meshes that they read and write.  For example, MOAB’s CUB reader imports not only the mesh saved from CUBIT, but also the grouping of mesh entities into sets which reflect the geometric topology of the model used to generate the mesh.  A more detailed description of meta-data conventions used in MOAB’s file readers and writers, and the individual file reader/writer header files in src/io for details about the specific readers and writers are available \href{metadata.html}{here}.
00598 
00599 Three specific file readers in MOAB bear further discussion: MOAB’s native HDF5-based file reader/writer; the CUB reader, used to import mesh and meta-data represented in CUBIT; and the CGM reader, which imports geometric models.  These are described next.
00600 
00601  \ref contents
00602 
00603   \subsection native 4.6.1. Native HD5-Based Reader/Writer
00604 
00605 A mesh database must be able to save and restore the data stored in its data model, at least to the extent to which it understands the semantics of that data.  MOAB defines an HDF5-based file format that can store data embedded in MOAB.  By convention, these files are given an “.h5m” file extension.  When reading or writing large amounts of data, it is recommended to use this file format, as it is the most complete and also the most efficient of the file readers/writers in MOAB. 
00606 
00607   \subsection cub 4.6.2. CUB Reader
00608 
00609 CUBIT is a toolkit for generating tetrahedral and hexahedral finite element meshes from solid model geometry [20].  This tool saves and restores data in a custom “.cub” file, which stores both mesh and geometry (and data relating the two).  The CUB reader in MOAB can import and interpret much of the meta-data information saved in .cub files.  The information read from .cub files, and stored in the MOAB data model, includes:
00610 
00611 - Geometric model entities and topology
00612 - Model entity names and ids
00613 - Groups, element blocks, nodesets, and sidesets, including model entities stored in them
00614 - Mesh scheme and interval size information assigned to model entities
00615 .
00616 
00617 Note that although information about model entities is recovered, MOAB by default does not depend on a solid modeling engine; this information is stored in the form of entity sets and parent/child relations between them.  
00618 
00619  \ref contents
00620 
00621   \subsection cgm 4.6.3. CGM Reader
00622 
00623 The Common Geometry Module (CGM) [21] is a library for representing solid model and other types of solid geometry data.  The CUBIT mesh generation toolkit uses CGM for its geometric modeling support, and CGM can restore geometric models in the exact state in which they were represented in CUBIT.  MOAB contains a CGM reader, which can be enabled with a configure option.  Using this reader, MOAB can read geometric models, and represent their model topology using entity sets linked by parent/child relations.  The mesh in these models comes directly from the modeling engine faceting routines; these are the same facets used to visualize solid models in other graphics engines.  When used in conjunction with the VisIt visualization tool (see Section 4.1), this provides a solution for visualizing geometric models.  The figure below  shows a model imported using MOAB’s CGM reader and visualized with VisIt.
00624 
00625 \image html simple.png
00626 
00627 \ref contents
00628 
00629  \subsection fourseven 4.7. AHF Representation
00630 
00631 Currently, the upward (vertex to entities) adjacencies are created and stored the first time a query requiring the adjacency is performed. Any non-vertex entity to entity adjacencies are performed using boolean operations on vertex-entity adjacencies. Because of this approach, such adjacency queries might become expensive with increasing dimension. We have added an alternative approach for obtaining adjacencies using the Array-based Half-Facet (AHF) representation[22]. The AHF uses sibling half-facets as a core abstraction which are generalizations of the opposite half-edge and half-face data structure for 2D and 3D manifold meshes. The AHF data structure consists of two essential maps: 1) the mapping between all sibling half-facets (sibhfs) and, 2) the mapping from each vertex to some incident half-facet (v2hf). The entire range of adjacencies (higher-, same- and lower-dimension) are computed using these two maps.
00632 
00633 The easiest way to avail this feature is to configure MOAB with " --enable-ahf " option. The adjacency queries can then be performed through calls to the preserved interface function "get_adjacencies" returning the values in a standard vector. Currently, returning adjacent entityhandles in MOAB::Range is not supported for AHF-based queries. There is one key difference between the native MOAB (adjacency-list based) and AHF based adjacency calls using the "get_adjacencies" interface. In native MOAB adjacency calls, the same-dimensional queries return the query entities whereas for AHF it would return the same-dimensional entities connected to the query entities via a lower dimensional facet. Thus the entire range ( higher-dimensional, same-dimensional, lower-dimensional) of adjacencies can be obtained using the same interface. Similar to MOAB's native adjacency lists, the AHF maps are created during the first adjacency call which will make the first adjacency call expensive.
00634 
00635 In the current release, AHF based adjacencies calls do not support the following cases:
00636   - polygon/polyhedral meshes,
00637   - mixed entity type meshes,
00638   - meshsets,
00639   - create_if_missing option set to true, and
00640   - modified meshes.
00641 
00642 The support for these would be gradually added in the next releases. In these cases, any adjacency call would revert back to MOAB's native adjacency list based queries.
00643 
00644 If for some reason, the user does not want to configure MOAB with AHF but would still like to use the AHF-based adjacencies for certain queries, they could use the following three interface functions provided in the HalfFacetRep class which implements the AHF maps and adjacency queries:
00645   - initialize : This function creates all the necessary AHF maps for the input mesh and hence should be called before any adjacency calls are made.
00646   - get_adjacencies: Function for adjacency calls.
00647   - deinitialize: This function deletes all the AHF maps and should be called after all AHF-based adjacency calls have been performed.
00648 
00649 TODO:: Other features to be added
00650   - obtain ring neighborhoods with support for half-rings
00651   - efficient extraction of boundaries
00652 
00653   \ref contents
00654 
00655   \subsection foureight 4.8. Uniform Mesh Refinement
00656   Many applications require a hierarchy of successively refined meshes for a number of purposes such as convergence studies, to use multilevel methods like multigrid, generate large meshes in parallel computing to support increasing mesh sizes with increase in number of processors, etc. Uniform mesh refinement provides a simple and efficient way to generate such hierarchies via successive refinement of the mesh at a previous level. It also provides a natural hierarchy via parent and child type of relationship between entities of meshes at different levels. Generally, the standard nested refinement patterns used are the subdivision schemes from 1 to 4 for 2D (triangles, quads) and 1 to 8 for 3D (tets, hexes) entity types. However, many applications might require degree 3 or more for p-refinements.
00657 
00658  MOAB supports generation of a mesh hierarchy i.e., a sequence of meshes with user specified degrees for each level of refinement, from an initial unstructured mesh with support for higher degrees of refinement (supported degrees are listed later).  Thus MOAB supports multi-degree and multi-level mesh generation via uniform refinement. The following figure shows the initial and most refined mesh for four simple meshes to illustrate the multi-degree capability.
00659 
00660   \image html uref_allEtype.png "Uniform Refinement of 2D and 3D meshes"
00661 
00662   Applications using mesh hierarchies require two types of mesh access: intralevel and interlevel. The intralevel access involves working with the mesh at a particular level whereas interlevel access involves querying across different levels. In order to achieve data locality with reduced cache misses for efficient intralevel mesh access, old vertices in the previous i.e. immediate parent mesh are duplicated in the current level. All the entities thus created for the current level use the new entityhandles of old vertices along with the handles of the new vertices. This design makes mesh at each level of the hierarchy independent of those at previous levels. For each mesh in the hierarchy, a MESHSET is created and all entities of the mesh are added to this MESHSET. Thus the meshes of the hierarchy are accessible via these level-wise MESHSET handles.
00663 
00664  For interlevel queries, separate interface functions are defined to allow queries across different levels. These queries mainly allow obtaining the parent at some specified parent level for a child from later levels and vice versa. The child-parent queries are not restricted to a level difference of one as the internal array-based layout of the memory allows traversing between different levels via index relations.
00665 
00666  The hierarchy generation capability is implemented in NestedRefine class. In Table 4, the user interface functions are briefly described. The main hierarchy generating function takes as input a sequence of refinement degrees to be used for generating mesh at each level from the previous level, the total number of levels in the hierarchy. It returns EntityHandles for the meshsets created for the mesh at each level. The number of levels in the hierarchy is prefixed (by the user) and cannot change during hierarchy generation. An example of how to generate a hierarchy can be found  under examples/UniformRefinement.cpp.
00667 
00668  The next three functions for getting the coordinates, connectivity and adjacencies are standard operations to access and query a mesh. The coordinates and connectivity functions, which are similar to their counterparts in MOAB Interface class, allows one to query the mesh at a specific level via its level index. The reason to provide such similar functions is to increase computational efficiency by utilizing the direct access to memory pointers to EntitySequences created during hierarchy generation under the NestedRefine class instead of calling the standard interfaces where a series of calls have to made before the EntitySequence of the requested entity is located. It is important to note that any calls to the standard interfaces available under Interface class should work.
00669 
00670  The underlying mesh data structure used for uniform refinement is the AHF datastructure implemented in HalfFacetRep class. This direct dependence on HalfFacetRep class removes the necessity to configure MOAB with --enable-ahf flag in order to use the uniform refinement capability. During the creation of an object of the NestedRefine class, it will internally initialize all the relevant AHF maps for the mesh in memory. During mesh hierarchy generation, after the creation of mesh at each level all the relevant AHF maps are updated to allow query over it.
00671 
00672  For interlevel queries, currently three kinds of queries are supported. The child to parent function allows querying for a parent in any of the previous levels (including the initial mesh). The parent to child, on the other hand, returns all its children from a requested child level. These two types of queries are only supported for entities, not vertices. A separate vertex to entities function is provided which returns entities from the previous level that are either incident or contain this vertex.
00673 
00674   \subsection tablefour Table 4: User interface functions NestedRefine class.
00675 
00676 <table border="1">
00677 <tr>
00678 <th>Function group</th>
00679 <th>Function</th>
00680 <th>Description</th>
00681 </tr>
00682 <tr>
00683 <td>Hierarchy generation</td>
00684 <td>generate_mesh_hierarchy</td>
00685 <td>Generate a mesh hierarchy with a given sequence of refinement degree for each level. </td>
00686 </tr>
00687 <tr>
00688 <td>Vertex coordinates</td>
00689 <td>get_coordinates</td>
00690 <td>Get vertex coordinates</td>
00691 </tr>
00692 <tr>
00693 <td>Connectivity</td>
00694 <td>get_connectivity</td>
00695 <td>Get connectivity of non-vertex entities</td>
00696 </tr>
00697 <tr>
00698 <td>Adjacencies</td>
00699 <td>get_adjacencies</td>
00700 <td>Get topologically adjacent entities</td>
00701 </tr>
00702 <tr>
00703 <td>Interlevel Queries</td>
00704 <td>child_to_parent</td>
00705 <td>Get the parent entity of a child from a specified parent level</td>
00706 </tr>
00707 <tr>
00708 <td>Interlevel Queries</td>
00709 <td>parent_to_child</td>
00710 <td>Get all children of the parent from a specified child level</td>
00711 </tr>
00712 <tr>
00713 <td>Interlevel Queries</td>
00714 <td>vertex_to_entities</td>
00715 <td>Get all the entities in the previous level incident on or containing the vertex</td>
00716 </tr>
00717 </table>
00718 
00719   \subsection tablefive Table 5: The refinement degrees currently supported.
00720 
00721 <table border="1">
00722 <tr>
00723 <th>Mesh Dimension</th>
00724 <th>EntityTypes</th>
00725 <th>Degree</th>
00726 <th>Number of children</th>
00727 </tr>
00728 <tr>
00729 <td>1</td>
00730 <td>MBEDGE</td>
00731 <td>2, 3, 5</td>
00732 <td>2, 3, 5 </td>
00733 </tr>
00734 <tr>
00735 <td>2</td>
00736 <td>MBTRI, MBQUAD</td>
00737 <td>2, 3, 5</td>
00738 <td>4, 9, 25 </td>
00739 </tr>
00740 <tr>
00741 <td>3</td>
00742 <td>MBTET, MBHEX</td>
00743 <td>2, 3</td>
00744 <td>8, 27 </td>
00745 </tr>
00746 </table>
00747 
00748 
00749 In Table 5, the currently supported degrees of refinement for each dimension is listed along with the number of children created for each such degree for a single entity. The following figure shows the cpu times(serial run) for generating hierarchies with various degrees of refinement for each dimension and can be used by the user to guide in choosing the degrees of refinement for the hierarchy. For example, if a multilevel hierarchy is required, a degree 2 refinement per level would give a gradually increasing mesh with more number of levels. If a very refined mesh is desired quickly, then a small hierarchy with high-order refinement should be generated.
00750 
00751 \image html uref_timeEtype.png "Mesh sizes Vs. Time"
00752 
00753 
00754   Current support:
00755   - Single dimension(i.e, no curves in surface meshes or curves/surfaces in volume meshes)
00756   - Linear point projection
00757   - Serial
00758 
00759   TODO:
00760    - Mixed-dimensional meshes
00761    - Mixed-entity types
00762    - High-Order point projection
00763    - Parallel
00764 
00765 
00766  \ref contents
00767 
00768   \section parallel 5.Parallel Mesh Representation and Query
00769 
00770 A parallel mesh representation must strike a careful balance between providing an interface to mesh similar to that of a serial mesh, while also allowing the discovery of parallel aspects of the mesh and performance of parallel mesh-based operations efficiently.  MOAB supports a spatial domain-decomposed view of a parallel mesh, where each subdomain is assigned to a processor, lower-dimensional entities on interfaces between subdomains are shared between processors, and ghost entities can be exchanged with neighboring processors.  Locally, each subdomain, along with any locally-represented ghost entities, are accessed through a local MOAB instance.  Parallel aspects of the mesh, e.g. whether entities are shared, on an interface, or ghost entities, are embedded in the same data model (entities, sets, tags, interface) used in the rest of MOAB.  MOAB provides a suite of parallel functions for initializing and communicating with a parallel mesh, along with functions to query the parallel aspects of the mesh.
00771 
00772   \ref contents
00773 
00774   \subsection fiveone 5.1. Nomenclature & Representation
00775 
00776 Before discussing how to access parallel aspects of a mesh, several terms need to be defined:  
00777 
00778 <B>Shared entity:</B> An entity shared by one or several other processors.
00779 
00780 <B>Multi-shared entity:</B> An entity shared by more than two processors.
00781 
00782 <B>Owning Processor:</B> Each shared entity is owned by exactly one processor.  This processor has the right to set tag values on the entity and have those values propagated to any sharing processors.  
00783 
00784 <B>Part:</B> The collection of entities assigned to a given processor.  When reading mesh in parallel, the entities in a Part, along with any lower-dimensional entities adjacent to those, will be read onto the assigned processor.
00785 
00786 <B>Partition:</B> A collection of Parts which take part in parallel collective communication, usually associated with an MPI communicator.
00787 
00788 <B>Interface:</B> A collection of mesh entities bounding entities in multiple parts.  Interface entities are owned by a single processor, but are represented on all parts/processors they bound.
00789 
00790 <B>Ghost entity:</B> A shared, non-interface, non-owned entity.
00791 
00792 <B>Parallel status:</B> A characteristic of entities and sets represented in parallel. The parallel status, or “pstatus”, is represented by a bit field stored in an unsigned character, with bit values as described in Table 6.
00793 
00794   \subsection tablesix Table 6: Bits representing various parallel characteristics of a mesh.  Also listed are enumerated values that can be used in bitmask expressions; these enumerated variables are declared in MBParallelConventions.h.
00795 
00796 <table border="1">
00797 <tr>
00798 <th>Bit</th>
00799 <th>Name</th>
00800 <th>Represents</th>
00801 </tr>
00802 <tr>
00803 <td>0x1</td>
00804 <td>PSTATUS_NOT_OWNED</td>
00805 <td>Not owned by the local processor</td>
00806 </tr>
00807 <tr>
00808 <td>0x2</td>
00809 <td>PSTATUS_SHARED</td>
00810 <td>Shared by exactly two processorstd>
00811 </tr>
00812 <tr>
00813 <td>0x4</td>
00814 <td>PSTATUS_MULTISHARED</td>
00815 <td>Shared by three or more processors</td>
00816 </tr>
00817 <tr>
00818 <td>0x8</td>
00819 <td>PSTATUS_INTERFACE</td>
00820 <td>Part of lower-dimensional interface shared by multiple processors</td>
00821 </tr>
00822 <tr>
00823 <td>0x10</td>
00824 <td>PSTATUS_GHOST</td>
00825 <td>Non-owned, non-interface entities represented locally</td>
00826 </tr>
00827 </table>
00828 
00829 Parallel functionality is described in the following sections and in [26].  First, methods to load a mesh into a parallel representation are described; next, functions for accessing parallel aspects of a mesh are described; functions for communicating mesh and tag data are described.
00830 
00831   \ref contents
00832 
00833   \subsection fivetwo 5.2. Parallel Mesh Initialization
00834 
00835 Parallel mesh is initialized in MOAB in several steps:
00836 
00837 -#  Establish a local mesh on each processor, either by reading the mesh into that representation from disk, or by creating mesh locally through the normal MOAB interface.  
00838 -#  Find vertices, then other entities, shared between processors, based on a globally-consistent vertex numbering stored on the GLOBAL_ID tag.  
00839 -#  Exchange ghost or halo elements within a certain depth of processor interfaces with neighboring processors.  
00840 .
00841 
00842 These steps can be executed by a single call to MOAB’s load_file function, using the procedure described in Section 5.2.1.  Or, they can be executed in smaller increments calling specific functions in MOAB’s ParallelComm class, as described in Section 5.2.2.  Closely related to the latter method is the handling of communicators, described in more detail in Section.
00843 
00844   \ref contents
00845 
00846   \subsection initialization 5.2.1. Parallel Mesh Initialization by Loading a File
00847 
00848 In the file reading approach, a mesh must contain some definition of the partition (the assignment of mesh, usually regions, to processors).  Partitions can be derived from other set structures already on the mesh, or can be computed explicitly for that purpose by tools like mbzoltan (see Section 4.2).  For example, geometric volumes used to generate the mesh, and region-based material type assignments, are both acceptable partitions.  In addition to defining the groupings of regions into parts, the assignment of specific parts to processors can be done implicitly or using additional data stored with the partition.
00849 
00850 MOAB implements several specific methods for loading mesh into a parallel representation:
00851 
00852 - READ_PART: each processor reads only the mesh used by its part(s).
00853 
00854 - READ_DELETE: each processor reads the entire mesh, then deletes the mesh not used by its part(s).
00855 
00856 - BCAST_DELETE: the root processor reads and broadcasts the mesh; each processor then deletes the mesh not used by its part(s).
00857 
00858 The READ_DELETE and BCAST_DELETE methods are supported for all file types MOAB is able to read, while READ_PART is only implemented for MOAB’s native HDF5-based file format and netcdf and pnetcdf files from climate application.
00859 
00860 Various other options control the selection of part sets or other details of the parallel read process.  For example, the application can specify the tags, and optionally tag values, which identify parts, and whether those parts should be distributed according to tag value or using round-robin assignment.
00861 
00862 The options used to specify loading method, the data used to identify parts, and other parameters controlling the parallel read process, are shown in Table 7.
00863   \subsection tableseven Table 7: Options passed to MOAB’s load_file function identifying the partition and other parameters controlling the parallel read of mesh data.  Options and values should appear in option string as “option=val”, with a delimiter (usually “;”) between options.
00864 
00865 <table border="1">
00866 <tr>
00867 <th>Option</th>
00868 <th>Value</th>
00869 <th>Description</th>
00870 </tr>
00871 <tr>
00872 <td>PARTITION</td>
00873 <td><tag_name></td>
00874 <td>Sets with the specified tag name should be used as part sets</td>
00875 </tr>
00876 <tr>
00877 <td>PARTITION_VAL</td>
00878 <td><val1, val2-val3, ...></td>
00879 <td>Integer values to be combined with tag name, with ranges input using val1, val2-val3.  Not meaningful unless PARTITION option is also given.</td>
00880 </tr>
00881 <tr>
00882 <td>PARTITION_DISTRIBUTE</td>
00883 <td>(none)</td>
00884 <td>If present, or values are not input using PARTITION_VAL, sets with tag indicated in PARTITION option are partitioned across processors in round-robin fashion.</td>
00885 </tr>
00886 <tr>
00887 <td>PARALLEL_RESOLVE_SHARED_ENTS</td>
00888 <td><pd.sd></td>
00889 <td>Resolve entities shared between processors, where partition is made up of pd- dimensional entities, and entities of dimension sd and lower should be resolved.</td>
00890 </tr>
00891 <tr>
00892 <td>PARALLEL_GHOSTS</td>
00893 <td><gd.bd.nl[.ad]></td>
00894 <td>Exchange ghost elements at shared inter-processor interfaces.  Ghost elements of dimension gd will be exchanged.  Ghost elements are chosen going through bd-dimensional interface entities.  Number of layers of ghost elements is specified in nl.  If ad is present, lower-dimensional entities bounding exchanged ghost entities will also be exchanged; allowed values for ad are 1 (exchange bounding edges), 2 (faces), or 3 (edges and faces).</td>
00895 </tr>
00896 <td>PARALLEL_COMM</td>
00897 <td><id></td>
00898 <td>Use the ParallelComm with index <id>.  Index for a ParallelComm object can be checked with ParallelComm::get_id(), and a ParallelComm with a given index can be retrieved using ParallelComm::get_pcomm(id).</td>
00899 </tr>
00900 <tr>
00901 <tr>
00902 <td>CPUTIME</td>
00903 <td>(none)</td>
00904 <td>Print cpu time required for each step of parallel read & initialization.</td>
00905 </tr>
00906 <tr>
00907 <td>MPI_IO_RANK</td>
00908 <td><r></td>
00909 <td>If read method requires reading mesh onto a single processor, processor with rank r is used to do that read.</td>
00910 </tr>
00911 </table>
00912 
00913 Several example option strings controlling parallel reading and initialization are:
00914 
00915 <B>“PARALLEL=READ_DELETE; PARTITION=MATERIAL_SET; PARTITION_VAL=100, 200, 600-700”:</B> The whole mesh is read by every processor; this processor keeps mesh in sets assigned the tag whose name is “MATERIAL_SET” and whose value is any one of 100, 200, and 600-700 inclusive.
00916 
00917 <B>“PARALLEL=READ_PART; PARTITION=PARALLEL_PARTITION, PARTITION_VAL=2”:</B> Each processor reads only its mesh; this processor, whose rank is 2, is responsible for elements in a set with the PARALLEL_PARTITION tag whose value is 2.  This would by typical input for a mesh which had already been partitioned with e.g. Zoltan or Parmetis.
00918 
00919 <B>“PARALLEL=BCAST_DELETE; PARTITION=GEOM_DIMENSION, PARTITION_VAL=3, PARTITION_DISTRIBUTE”:</B> The root processor reads the file and broadcasts the whole mesh to all processors.  If a list is constructed with entity sets whose GEOM_DIMENSION tag is 3, i.e. sets corresponding to geometric volumes in the original geometric model, this processor is responsible for all elements with index R+iP, i >= 0 (i.e. a round-robin distribution).
00920 
00921  \ref contents
00922 
00923   \subsection functions 5.2.2. Parallel Mesh Initialization Using Functions
00924 
00925 After creating the local mesh on each processor, an application can call the following functions in ParallelComm to establish information on shared mesh entities.  See the [http://ftp.mcs.anl.gov/pub/fathom/moab-docs/HelloParMOAB_8cpp-example.html example] in the MOAB source tree for a complete example of how this is done from an application.
00926 
00927 - ParallelComm::resolve_shared_entities (collective): Resolves shared entities between processors, based by default on GLOBAL_ID tag values of vertices. For meshes with more than 2 billion vertices, an opaque tag of size 8 should be used, in which a 64 bit integer is casted to that opaque tag.  Various forms are available, based on entities to be evaluated and maximum dimension for which entity sharing should be found.
00928 
00929 - ParallelComm::exchange_ghost_cells (collective): Exchange ghost entities with processors sharing an interface with this processor, based on specified ghost dimension (dimension of ghost entities exchanged), bridge dimension, number of layers, and type of adjacencies to ghost entities.  An entity is sent as a ghost if it is within that number of layers, across entities of the bridge dimension, with entities owned by the receiving processor, or if it is a lower-dimensional entity adjacent to a ghost entity and that option is requested.
00930 .
00931 
00932  \ref contents
00933 
00934   \subsection communicator 5.2.3. Communicator Handling
00935 
00936 The ParallelComm constructor takes arguments for an MPI communicator and a MOAB instance.  The ParallelComm instance stores the MPI communicator, and registers itself with the MOAB instance.  Applications can specify the ParallelComm index to be used for a given file operation, thereby specifying the MPI communicator for that parallel operation.  For example:
00937 
00938 \code
00939 using namespace moab;
00940 // pass a communicator to the constructor, getting back the index
00941 MPI_Comm my_mpicomm;
00942 int pcomm_index;
00943 ParallelComm my_pcomm(moab, my_mpicomm, &pcomm_index);
00944 
00945 // write the pcomm index into a string option
00946 char load_opt[32];
00947 sprintf(load_opt, "PARALLEL=BCAST_DELETE;PARALLEL_COMM=%d", 
00948    pcomm_index);
00949 
00950 // specify that option in a parallel read operation
00951 ErrorCode rval = moab->load_file(load_opt, fname, ...)
00952 \endcode
00953 
00954 In the above code fragment, the ParallelComm instance with index pcomm_index will be used in the parallel file read, so that the operation executes over the specified MPI communicator.  If no ParallelComm instance is specified for a parallel file operation, a default instance will be defined, using MPI_COMM_WORLD.
00955 
00956 Applications needing to retrieve a ParallelComm instance created previously and stored with the MOAB instance, e.g. by a different code component, can do so using a static function on ParallelComm:
00957 
00958 \code
00959 ParallelComm *my_pcomm = ParallelComm::get_pcomm(moab, pcomm_index);
00960 \endcode
00961 
00962 ParallelComm also provides the ParallelComm::get_all_pcomm function, for retrieving all ParallelComm instances stored with a MOAB instance.  For syntax and usage of this function, see the MOAB online documentation for ParallelComm.hpp [11].
00963 
00964  \ref contents
00965 
00966   \subsection fivethree 5.3. Parallel Mesh Query Functions
00967 
00968 Various functions are commonly used in parallel mesh-based applications.  Functions marked as being collective must be called collectively for all processors that are members of the communicator associated with the ParallelComm instance used for the call.
00969 
00970 <B>ParallelComm::get_pstatus:</B>  Get the parallel status for the entity.
00971 
00972 <B>ParallelComm::get_pstatus_entities:</B> Get all entities whose pstatus satisfies (pstatus & val).
00973 
00974 <B>ParallelComm::get_owner:</B> Get the rank of the owning processor for the specified entity.
00975 
00976 <B>ParallelComm::get_owner_handle:</B> Get the rank of the owning processor for the specified entity, and the entity's handle on the owning processor.
00977 
00978 <B>ParallelComm::get_sharing_data:</B> Get the sharing processor(s) and handle(s) for an entity or entities.  Various overloaded versions are available, some with an optional “operation” argument, where Interface::INTERSECT or Interface::UNION can be specified.  This is similar to the operation arguments to Interface::get_adjacencies.
00979 
00980 <B>ParallelComm::get_shared_entities:</B> Get entities shared with the given processor, or with all processors.  This function has optional arguments for specifying dimension, whether interface entities are requested, and whether to return just owned entities.
00981 
00982 <B>ParallelComm::get_interface_procs:</B> Return all processors with whom this processor shares an interface.
00983 
00984 <B>ParallelComm::get_comm_procs:</B> Return all processors with whom this processor communicates.
00985 
00986   \ref contents
00987 
00988   \subsection fivefour 5.4. Parallel Mesh Communication
00989 
00990 Once a parallel mesh has been initialized, applications can call the ParallelComm::exchange_tags function for exchanging tag values between processors.  This function causes the owning processor to send the specified tag values for all shared, owned entities to other processors sharing those entities.  Asynchronous communication is used to hide latency, and only point-to-point communication is used in these calls.
00991 
00992   \ref contents
00993 
00994   \section applications 6.Building MOAB-Based Applications
00995 
00996 There are two primary mechanisms supported by MOAB for building applications, one based on MOAB-defined make variables, and the other based on the use of libtool and autoconf.  Both assume the use of a “make”-based build system.  
00997 
00998 The easiest way to incorporate MOAB into an application’s build process is to include the “moab.make” file into the application’s Makefile, adding the make variables MOAB_INCLUDES and MOAB_LIBS_LINK to application’s compile and link commands, respectively.  MOAB_INCLUDES contains compiler options specifying the location of MOAB include files, and any preprocessor definitions required by MOAB.  MOAB_LIBS_LINK contains both the options telling where libraries can be found, and the link options which incorporate those libraries into the application.  Any libraries depended on by the particular configuration of MOAB are included in that definition, e.g. the HDF5 library.  Using this method to incorporate MOAB is the most straightforward; for example, the following Makefile is used to build one of the example problems packaged with the MOAB source:
00999 \code
01000 include ${MOAB_LIB_DIR}/moab.make
01001 
01002 GetEntities : GetEntities.o
01003     ${CXX} $< ${MOAB_LIBS_LINK} -o $@
01004 
01005 .cpp.o : 
01006     ${CXX} ${MOAB_INCLUDES} -c $<
01007 \endcode
01008 
01009 Here, the MOAB_LIB_DIR environment variable or make argument definition specifies where the MOAB library is installed; this is also the location of the moab.make file.  Once that file has been included, MOAB_INCLUDES and MOAB_LIBS_LINK can be used, as shown.
01010 
01011 Other make variables are defined in the moab.make file which simplify building applications:
01012 
01013 - MOAB_LIBDIR, MOAB_INCLUDEDIR: the directories into which MOAB libraries and include files will be installed, respectively.  Note that some include files are put in a subdirectory named “moab” below that directory, to reflect namespace naming conventions used in MOAB.
01014 
01015 - MOAB_CXXFLAGS, MOAB_CFLAGS, MOAB_LDFLAGS: Options passed to the C++ and C compilers and the linker, respectively.
01016 
01017 - MOAB_CXX, MOAB_CC, MOAB_FC: C++, C, and Fortran compilers specified to MOAB at configure time, respectively.
01018 .
01019 
01020 The second method for incorporating MOAB into an application’s build system is to use autoconf and libtool.  MOAB is configured using these tools, and generates the “.la” files that hold information on library dependencies that can be used in application build systems also based on autoconf and libtool.  Further information on this subject is beyond the scope of this User’s Guide; see the “.la” files as installed by MOAB, and contact the MOAB developer’s mailing list [9] for more details.
01021 
01022   \ref contents
01023 
01024   \section implementation  7.iMesh (ITAPS Mesh Interface) Implementation in MOAB
01025 
01026 iMesh is a common API to mesh data developed as part of the Interoperable Tools for Advanced Petascale Simulations (ITAPS) project [23].  Applications using the iMesh interface can operate on any implementation of that interface, including MOAB.  MOAB-based applications can take advantage of other services implemented on top of iMesh, including the MESQUITE mesh improvement toolkit [24].
01027 
01028 MOAB’s native interface is accessed through the Interface abstract C++ base class.  Wrappers are not provided in other languages; rather, applications wanting to access MOAB from those languages should do so through iMesh.  In most cases, the data models and functionality available through MOAB and iMesh are identical.  However, there are a few differences, subtle and not-so-subtle, between the two:
01029 
01030 <B>SPARSE tags used by default:</B> MOAB’s iMesh implementation creates SPARSE tags by default, because of semantic requirements of other tag-related functions in iMesh.  To create DENSE tags through iMesh, use the iMesh_createTagWithOptions extension function (see below).
01031 
01032 <B>Higher-order elements:</B> ITAPS currently handles higher-order elements (e.g. a 10-node tetrahedron and a 27-node hexahedra) [23].  As described in [sec-entities], MOAB supports higher-order entities by allowing various numbers of vertices to define topological entities like quadrilateral or tetrahedron.  Applications can specify flags to the connectivity and adjacency functions specifying whether corner or all vertices are requested.
01033 
01034 <B>Self-adjacencies:</B> In MOAB’s native interface, entities are always self-adjacent<sup>6</sup>; that is, adjacencies of equal dimension requested from an entity will always include that entity, while from iMesh will not include that entity.
01035 
01036 <B>Option strings:</B> The iMesh specification requires that options in the options string passed to various functions (e.g. iMesh_load) be prepended with the implementation name required to parse them, and delimited with spaces.  Thus, a MOAB-targeted option would appear as “moab:PARALLEL=READ_PART moab:PARTITION=MATERIAL_SET”.
01037 
01038 To provide complete MOAB support from other languages through iMesh, a collection of iMesh extension functions are also available.  A general description of these extensions appears below; for a complete description, see the online documentation for iMesh-extensions.h [11].
01039 
01040 - Recursive get_entities functions: There are many cases where sets include other sets.  MOAB provides iMesh_getEntitiesRec, and other recursive-supporting functions, to get all non-set entities of a given type or topology accessible from input set(s).  Similar functions are available for number of entities of a given type/topology.
01041 
01042 - Get entities by tag, and optionally tag value: It is common to search for entities with a given tag, and possibly tag value(s); functions like iMesh_getEntitiesByTag are provided for this purpose.
01043 
01044 - Options to createTag: To provide more control over the tag type, the iMesh_createTagWithOptions is provided.  The storage type is controlled with the “
01045 
01046 - MBCNType: Canonical numbering evaluations are commonly needed by applications, e.g. to apply boundary conditions locally.  The MBCN package provides these evaluations in terms of entity types defined in MOAB [11]; the getMBCNType is required to translate between iMesh_Topology and MBCN type.
01047 
01048 - Iterator step: Step an iterator a specified number of entities; allows advancement of an iterator without needing to allocate memory to hold the entity handles stepped over.
01049 
01050 - Direct access to tag storage: The Interface::tag_iterate function allows an application get a pointer to the memory used to store a given tag.  For dense tags on contiguous ranges of entities, this provides more efficient access to tags.  The iMesh functionn iMesh_tagIterate provides access to this functionlity.  See examples/TagIterateC.c and examples/TagIterateF.F for examples of how to use this from C and Fortran, respectively. 
01051 .
01052 
01053 As required by the iMesh specification, MOAB generates the “iMesh-Defs.inc” file and installs it with the iMesh and MOAB libraries.  This file defines make variables which can be used to build iMesh-based applications.  The method used here is quite similar to that used for MOAB itself (see Section 6).  In particular, the IMESH_INCLUDES and IMESH_LIBS make variables can be used with application compile and link commands, respectively, with other make variables similar to those provided in moab.make also available.
01054 
01055 Note that using the iMesh interface from Fortran-based applications requires a compiler that supports Cray pointers, along with the pass-by-value (%VAL) extension.  Almost all compilers support those extensions; however, if using the gcc series of compilers, you must use gfortran 4.3 or later.
01056 
01057 <sup>5</sup>There are currently no implementations of this interface.
01058 
01059 <sup>6</sup>iMesh and MOAB both define adjacencies using the topological concept of closure.  Since the closure of an entity includes the entity itself, the d-dimensional entities on the closure of a given entity should include the entity itself.
01060 
01061   \ref contents
01062 
01063   \section pymoab 8.Python Interface (PyMOAB)
01064 
01065 A Python interface to MOAB's essential core functionality and a few other tools has been added as of Version 5.0. The pymoab module can be used to interactively interrogate existing mesh files or prototype MOAB-based algorithms. It can also be connected to other Python applications or modules for generation, manipulation, and visualization of a MOAB mesh and mesh data. Examples of this can be found in the laplaciansmoother.py and yt2moab.py files. Interaction with the PyMOAB interface is intended to be somewhat analagous to interaction with the MOAB C++ API. A simple example of file loading and mesh interrogation can be found in interrogate_mesh.py
01066 
01067 <B>PyMOAB uses NumPy internally </B> to represent data and vertex coordinates though other properly formed data structures can be used to create vertices, set data, etc. Data and vertex coordinates will always be returned in the form of NumPy arrays, however. 
01068 
01069 <B>EntityHandles</B> are represented by Python long integers. Arrays of these values are commonly returned from calls as MOAB Ranges.
01070 
01071 <B>MOAB ErrorCodes</B> - Errors are automatically checked internally by the PyMOAB instance, raising various exceptions depending on the error that occurs. These exceptions can be handled as raised from the functions or acceptable error values returned can also be specified via the exceptions parameter common to all PyMOAB functions in which case no exception will be raised.
01072 
01073 Documentation for PyMOAB functions is provided as part of this User's Guide, but can also be accessed in the Python interpreter by calling `help(<function_or_method>)`.
01074 
01075 Examples of PyMOAB usage can be found in the /examples/python/ directory.
01076 
01077 
01078   \ref contents
01079 
01080   \section representation 9.Structured Mesh Representation
01081 
01082 A structured mesh is defined as a D-dimensional mesh whose interior vertices have 2D connected edges.   Structured mesh can be stored without connectivity, if certain information is kept about the parametric space of each structured block of mesh.  MOAB can represent structured mesh with implicit connectivity, saving approximately 57% of the storage cost compared to an unstructured representation<sup>7</sup>.  Since connectivity must be computed on the fly, these queries execute a bit slower than those for unstructured mesh.  More information on the theory, and design behind MOAB's structured mesh representation can be found in [8].
01083 
01084 Currently, MOAB's structured mesh representation can only be used by creating structured mesh at runtime; that is, structured mesh is saved/restored in an unstructured format in MOAB's HDF5-based native save format.  For more details on how to use MOAB's structured mesh representation, see the scdseq_test.cpp source file in the test/ directory.
01085 
01086 <sup>7</sup> This assumes vertex coordinates are represented explicitly, and that there are approximately the same number of vertices and hexahedra in a structured hex mesh.
01087 
01088  \ref contents
01089 
01090   \section element 10.Spectral Element Meshes
01091 
01092 The Spectral Element Method (SEM) is a high-order method, using a polynomial Legendre interpolation basis with Gauss-Lobatto quadrature points, in contrast to the Lagrange basis used in (linear) finite elements [25].  SEM obtains exponential convergence with decreasing mesh characteristic sizes, and codes implementing this method typically have high floating-point intensity, making the method highly efficient on modern CPUs.  Most Nth-order SEM codes require tensor product cuboid (quad/hex) meshes, with each d-dimensional element containing (N+1)^d degrees of freedom (DOFs).  There are various methods for representing SEM meshes and solution fields on them; this document discusses these methods and the tradeoffs between them.  The mesh parts of this discussion are given in terms of the iMesh mesh interface and its implementation by the MOAB mesh library.
01093 
01094 The figure above shows a two-dimensional 3rd-order SEM mesh consisting of four quadrilaterals.  For this mesh, each quadrilateral has (N+1)^2=16 DOFs, with corner and edge degrees of freedom shared between neighboring quadrilaterals.
01095 
01096   \ref contents
01097 
01098   \subsection tenone 10.1. Representations
01099 
01100 There are various representations of this mesh in a mesh database like MOAB, depending on how DOFs are related to mesh entities and tags on those entities.  We mention several possible representations:
01101 
01102 -# Corner vertices, element-based DOFs: Each quadrilateral is defined by four vertices, ordered in CCW order typical of FE meshes.  DOFs are stored as tags on quadrilaterals, with size (N+1)^2 values, ordered lexicographically (i.e. as a 2D array tag(i,j) with i varying faster than j.)  In the figure above, the connectivity for face 1 would be (1, 4, 16, 13), and DOFs would be ordered (1..16).  Note that in this representation, tag values for DOFs shared by neighboring elements must be set multiple times, since there are as many copies of these DOFs as elements sharing them.
01103 -# High-order FE-like elements: Each DOF is represented by a mesh vertex. Quadrilaterals each have (N+1)^2 vertices, ordered as they would be for high-order finite elements (corner vertices first, then mid-edge and mid-face elements; see [22]).  Mid -face, -edge, and -region vertices for a given edge/face/region would be ordered lexicographically, according to positive direction in a corresponding reference element.  In the figure above, the connectivity array for face 1 would be (1, 4, 16, 13, 2, 3, 8, 12, 14, 15, 5, 9, 6, 7, 10, 11).  DOF values are stored as tags on vertices.  Since DOFs are uniquely associated with vertices and vertices are shared by neighboring elements, tag values only need to be set once.  Full vertex-quadrilateral adjacencies are available, for all vertices.
01104 -# Linear FE-like elements, one vertex per DOF, array with DOF vertices: Each quadrilateral is defined by four (corner) vertices, with additional vertices representing mid-edge and mid-face DOFs.  An additional “DOF array” tag is assigned to each quadrilateral, storing the array of vertices representing the (N+1)^2 DOFs for the quadrilateral, ordered lexicographically.  For the figure above, the connectivity array for face 1 would be (1, 4, 16, 13), and the DOF array would be (1..16), assuming that vertex handles are integers as shown in the figure.  DOF values are stored as tags on vertices, and lexicographically-ordered arrays of DOFs can be retrieved using the DOF array tag as input to the tag_get_data function in MOAB.  Adjacency functions would only be meaningful for corner vertices, but tag values would only need to be set once per DOF.
01105 -# High-order FE-like elements, array with DOF vertices: This is a combination of options 2 and 3.  The advantage would be full vertex-quad adjacency support and direct availability of lexicographically-ordered vertex arrays, at the expense of more memory.
01106 -# Convert to linear mesh: Since a spectral element is a cuboid with higher-order vertices, it can always be converted to N^2 linear cuboids using the high-order vertices as corners of the finer quads/hexes.  This is how readers in ParaView and VisIt typically import spectral meshes (CAM-SE also exports connectivity in this form).
01107 
01108 As a convenience for applications, functions could also be provided for important tasks, like assembling the vertex handles for an entity in lexographic order (useful for option 2 above), and getting an array of tag values in lexicographic order (for option 3 above).
01109 
01110   \ref contents
01111 
01112   \subsection tentwo 10.2. Tradeoffs
01113 
01114 There are various competing tradeoffs in the various representation types.  These include:
01115 
01116 - Adjacencies: being able to retrieve the element(s) using a given (corner or higher-order) vertex.
01117 - Connectivity list: being able to retrieve the connectivity of a given element, consisting of all (corner + higher-order) vertices in the element, usually in lexicographical order.  This is closely linked with being able to access the connectivity list as a const*, i.e. using the list straight from memory without needing to copy it.
01118 - Memory vs. time: There is a memory vs. execution time tradeoff between duplicating interface vertex solution/tag variables in neighboring elements (more memory but more time-efficient and allows direct access to tag storage by applications) versus using vertex-based tags (less memory but requires assembly of variables into lexicographically-ordered arrays, and prevents direct access from applications).
01119 .
01120 
01121 The lower-memory option (storing variables on vertices and assembling into lexicographically-ordered arrays for application use) usually ends up costing more in memory anyway, since applications must allocate their own storage for these arrays.  On the other hand, certain applications will always choose to do that, instead of sharing storage with MOAB for these variables.  In the case where applications do share memory with MOAB, other tools would need to interpret the lexicographically-ordered field arrays specially, instead of simply treating the vertex tags as a point-based field.
01122 
01123   \ref contents
01124 
01125   \subsection tenthree 10.3. MOAB Representation
01126 In choosing the right MOAB representation for spectral meshes, we are trying to balance a) minimal memory usage, b) access to properly-ordered and -aligned tag storage, and c) maximal compatibility with tools likely to use MOAB.  The solution we propose is to use a representation most like option 2) above, with a few optional behaviors based on application requirements.  
01127 
01128 In brief, we propose to represent elements using the linear, FE-ordered connectivity list (containing only corner vertices from the spectral element), with field variables written to either vertices, lexicographically-ordered arrays on elements, or both, and with a lexicographically-ordered array (stored on tag SPECTRAL_VERTICES) of all (corner+higher-order) vertices stored on elements.  In the either/or case, the choice will be evident from the tag size and the entities on which the tag is set.  In the both case, the tag name will have a “-LEX” suffix for the element tags, and the size of the element tag will be (N+1)^2 times that of the vertex-based tag.  Finally, the file set containing the spectral elements (or the root set, if no file set was input to the read) will contain a “SPECTRAL_ORDER” tag whose value is N.  These conventions are described in the “Metadata Information” document distributed with the MOAB source code.
01129 
01130   \ref contents
01131 
01132   \section performance 11.Performance and Using MOAB Efficiently from Applications
01133 
01134 MOAB is designed to operate efficiently on groups of entities and for large meshes.  Applications will be most efficient when they operate on entities in groups, especially groups which are close in their order of creation.  The MOAB API is structured to encourage operations on groups of entities.  Conversely, MOAB will not perform as well as other libraries if there are frequent deletion and creation of entities.  For those types of applications, a mesh library using a C++ object-based representation is more appropriate.  In this section, performance of MOAB when executing a variety of tasks is described, and compared to that of other representations.  Of course, these metrics are based on the particular models and environments where they are run, and may or may not be representative of other application types.
01135 
01136 One useful measure of MOAB performance is in the representation and query of a large mesh.  MOAB includes a performance test, located in the test/perf directory, in which a single rectangular region of hexahedral elements is created then queried; the following steps are performed:
01137 
01138 - Create the vertices and hexes in the mesh
01139 - For each vertex, get the set of connected hexahedra
01140 - For each hex, get the connected vertices, their coordinates, average them, and assign them as a tag on the hexes
01141 .
01142 
01143 This test can be run on your system to determine the runtime and memory performance for these queries in MOAB.
01144 
01145   \ref contents
01146 
01147   \section error-handling 12.Error Handling
01148 
01149 Errors are handled through the routine MBError(). This routine calls MBTraceBackErrorHandler(), the default error handler which tries to print a traceback.
01150 
01151 The arguments to MBTraceBackErrorHandler() are the line number where the error occurred, the function where error was detected, the file in which
01152 the error was detected, the corresponding directory, the error message, and the error type.
01153 
01154 A small set of macros is used to make the error handling lightweight. These macros are used throughout
01155 the MOAB libraries and can be employed by the application programmer as well. When an error is first
01156 detected, one should set it by calling\n
01157 \code
01158 MB_SET_ERR(err_code, err_msg);
01159 \endcode
01160 Note, err_msg can be a string literal, or a C++ style output stream with << operators, such that the error message string is formated, like\n
01161 \code
01162 MB_SET_ERR(MB_FAILURE, "Failed " << n << " times");
01163 \endcode
01164 
01165 The user should check the return codes for all MOAB routines (and possibly user-defined routines as well) with\n
01166 \code
01167 ErrorCode rval = MOABRoutine(...);MB_CHK_ERR(rval);
01168 \endcode
01169 To pass back a new error message (if rval is not MB_SUCCESS), use\n
01170 \code
01171 ErrorCode rval = MOABRoutine(...);MB_CHK_SET_ERR(rval, "User specified error message string (or stream)");
01172 \endcode
01173 If this procedure is followed throughout all of the user’s libraries and codes, any error will by default generate
01174 a clean traceback of the location of the error.
01175 
01176 In addition to the basic macros mentioned above, there are some variations, such as (for more information, refer to src/moab/ErrorHandler.hpp):
01177 - MB_SET_GLB_ERR() to set a globally fatal error (for all processors)
01178 - MB_SET_ERR_RET() for functions that return void type
01179 - MB_SET_ERR_RET_VAL() for functions that return any data type
01180 - MB_SET_ERR_CONT() to continue execution instead of returning from current function
01181 
01182 The error control mechanism are enabled by default if a MOAB Core instance is created. Otherwise, the user needs to call
01183 MBErrorHandler_Init() and MBErrorHandler_Finalize() at the application level in the main function.
01184 For example code on error handling, please refer to examples/TestErrorHandling.cpp, examples/TestErrorHandlingPar.cpp and examples/ErrorHandlingSimulation.cpp.
01185 
01186   \ref contents
01187 
01188   \section conclusions 13.Conclusions and Future Plans
01189 
01190 MOAB, a Mesh-Oriented datABase, provides a simple but powerful data abstraction to structured and unstructured mesh, and makes that abstraction available through a function API.  MOAB provides the mesh representation for the VERDE mesh verification tool, which demonstrates some of the powerful mesh metadata representation capabilities in MOAB.  MOAB includes modules that import mesh in the ExodusII, CUBIT .cub and Vtk file formats, as well as the capability to write mesh to ExodusII, all without licensing restrictions normally found in ExodusII-based applications.  MOAB also has the capability to represent and query structured mesh in a way that optimizes storage space using the parametric space of a structured mesh; see Ref. [17] for details.
01191 
01192 Initial results have demonstrated that the data abstraction provided by MOAB is powerful enough to represent many different kinds of mesh data found in real applications, including geometric topology groupings and relations, boundary condition groupings, and inter-processor interface representation.  Our future plans are to further explore how these abstractions can be used in the design through analysis process.
01193 
01194   \ref contents
01195 
01196   \section references 14.References
01197 
01198 [1] Mahadevan, Vijay S., Iulian Grindeanu, Rajeev Jain, Patrick Shriwise, Navamita Ray, Paul Wilson, Tautges, Timothy J., "SIGMA -- MOAB.", URL: \href{http://sigma.mcs.anl.gov/}{http://sigma.mcs.anl.gov/}
01199 
01200 [2] Mahadevan, Vijay S., Iulian Grindeanu, Robert Jacob, and Jason Sarich. "Improving climate model coupling through a complete mesh representation: a case study with E3SM (v1) and MOAB (v5. x).", Geosci. Model Dev. Discuss., https://doi.org/10.5194/gmd-2018-280, in review, 2018.
01201 
01202 [3] Mahadevan, Vijay S., Elia Merzari, Timothy Tautges, Rajeev Jain, Aleksandr Obabko, Michael Smith, and Paul Fischer. "High-resolution coupled physics solvers for analysing fine-scale nuclear reactor design problems." Philosophical Transactions of the Royal Society A: Mathematical, Physical and Engineering Sciences 372, no. 2021 (2014): 20130381.
01203 
01204 [4] Yan, Mi, Kirk Jordan, Dinesh Kaushik, Michael Perrone, Vipin Sachdeva, Timothy J. Tautges, and John Magerlein. "Coupling a basin modeling and a seismic code using MOAB." Procedia Computer Science 9 (2012): 986-993.
01205 
01206 [5] Jacob, Robert, Jayesh Krishna, Xiabing Xu, Tim Tautges, Iulian Grindeanu, Rob Latham, Kara Peterson et al. "ParNCL and ParGAL: Data-parallel tools for postprocessing of large-scale Earth science data." Procedia Computer Science 18 (2013): 1245-1254.
01207 
01208 [6] Bohm, Tim D., Mohamed E. Sawan, Steve T. Jackson, and Paul PH Wilson. "Detailed nuclear analysis of ITER ELM coils." Fusion Engineering and Design 87, no. 5-6 (2012): 657-661.
01209 
01210 [7] Tautges, Timothy J., and Alvaro Caceres. "Scalable parallel solution coupling for multiphysics reactor simulation." In Journal of Physics: Conference Series, vol. 180, no. 1, p. 012017. IOP Publishing, 2009.
01211 
01212 [8] Tautges, Timothy J. "MOAB-SD: integrated structured and unstructured mesh representation." Engineering With Computers 20, no. 3 (2004): 286-293.
01213 
01214 [9] MOAB Users and Developers Email List., [email protected].
01215 
01216 [10] MOAB Announcement Email List., [email protected].
01217 
01218 [11] MOAB online documentation., http://ftp.mcs.anl.gov/pub/fathom/moab-docs/index.html
01219 
01220 [12] T.J. Tautges, “Canonical numbering systems for finite-element codes,” Communications in Numerical Methods in Engineering,  vol. Online, Mar. 2009.
01221 
01222 [13] L.A. Schoof and V.R. Yarberry, EXODUS II: A Finite Element Data Model,  Albuquerque, NM: Sandia National Laboratories, 1994.
01223 
01224 [14] M. PATRAN, “PATRAN User’s Manual,” 2005.
01225 
01226 [15] VisIt User's Guide.
01227 
01228 [16] Devine, Karen, Erik Boman, Robert Heaphy, Bruce Hendrickson, and Courtenay Vaughan. "Zoltan data management service for parallel dynamic applications." Computing in Science & Engineering 4, no. 2 (2002): 90.
01229 
01230 [17] METIS - Serial Graph Partitioning and Fill-reducing Matrix Ordering. http://glaros.dtc.umn.edu/gkhome/views/metis
01231 
01232 [18] Tautges, Timothy J., Paul PH Wilson, Jason A. Kraftcheck, Brandon M. Smith, and Douglass L. Henderson. "Acceleration techniques for direct use of CAD-based geometries in Monte Carlo radiation transport." (2009).
01233 
01234 [19] Kim, Hong-Jun, and Timothy J. Tautges. "EBMesh: An embedded boundary meshing tool." In Proceedings of the 19th International Meshing Roundtable, pp. 227-242. Springer, Berlin, Heidelberg, 2010.
01235 
01236 [20] G.D. Sjaardema, T.J. Tautges, T.J. Wilson, S.J. Owen, T.D. Blacker, W.J. Bohnhoff, T.L. Edwards, J.R. Hipp, R.R. Lober, and S.A. Mitchell, CUBIT mesh generation environment Volume 1: Users manual, Sandia National Laboratories, May 1994, 1994.
01237 
01238 [21] T.J. Tautges, “CGM: A geometry interface for mesh generation, analysis and other applications,” Engineering with Computers,  vol. 17, 2001, pp. 299-314.
01239 
01240 [22] V. Dyedov, N. Ray, D.Einstein, X. Jiao, T. Tautges, “AHF: Array-based half-facet data structure for mixed-dimensional and non-manifold meshes”, In Proceedings of 22nd International Meshing Roundtable, Orlando, Florida, October 2013.
01241 
01242 [23] Li, Xiaolin. Interoperable Technologies for Advanced Petascale Simulations. No. 41293. Stony Brook University, 2013.
01243 
01244 [24] Knupp, Patrick. "Mesh quality improvement for SciDAC applications." In Journal of Physics: Conference Series, vol. 46, no. 1, p. 458. IOP Publishing, 2006.
01245 
01246 [25] Deville, Michel O., Paul F. Fischer, Paul F. Fischer, and E. H. Mund. High-order methods for incompressible fluid flow. Vol. 9. Cambridge university press, 2002.
01247 
01248 [26] T. J. Tautges, J. A. Kraftcheck, N. Bertram, V. Sachdeva and J. Magerlein, "Mesh Interface Resolution and Ghost Exchange in a Parallel Mesh Representation," 2012 IEEE 26th International Parallel and Distributed Processing Symposium Workshops & PhD Forum, Shanghai, 2012, pp. 1670-1679, doi: 10.1109/IPDPSW.2012.208.
01249   \ref contents
01250 
01251   \page differences Differences Between iMesh and MOAB
01252 
01253   The data models used in MOAB and iMesh are quite similar, but not identical.The most significant differences are the following:
01254 
01255 - Tags: MOAB differentiates between DENSE, SPARSE, and BIT tags, using different storage models for each, while iMesh uses a single tag concept.  iMesh allows application to query whether an entity has been given a tag of a specified type; this query is incompatible with the concept of a DENSE tag with a default value.  Thus, MOAB’s iMesh implementation creates SPARSE tags by default, and tags created and accessed through this interface will use more memory than DENSE tags created through MOAB’s native interface.  To mitigate this problem, MOAB implements an extension of the iMesh_createTag function which allows specification of the tag type (DENSE, SPARSE, etc.) to be created.  See later in this section for more information.
01256 
01257 - Higher-order nodes: ITAPS currently handles higher-order elements (e.g. a 10-node tetrahedron) using a special “Shape” interface.  In this interface, higher-order nodes are only accessible through the AEntities which they resolve.  MOAB’s iMesh implementation provides access to higher-order nodes in the same manner described in Section  , by varying the number of vertices defining each entity.  As a result, if higher-order entities are used in a model, the functions returning connectivity and vertex adjacencies always return all vertices, rather than providing an option to return just corner vertices.
01258 
01259 - Self-adjacencies: iMesh specifies that entities are not self-adjacent; that is, requesting adjacencies of the same dimension/type results in an error.  MOAB does not consider this an error, returning the entity itself.
01260 
01261 - Adjacency table and AEntities: iMesh uses the concept of an “adjacency table” to determine which AEntities are available and created by default.  MOAB uses input arguments to the get_adjacencies functions to control whether AEntities are created.  These flags provide finer-grained control over AEntities, but make it slightly less convenient to ensure that AEntities of a given dimension are always created.
01262 .
01263 
01264  \page figures List of Figures
01265  
01266   This page is intended to be empty.
01267  
01268   \page tables List of Tables
01269  
01270      \ref tableone
01271 
01272      \ref tabletwo
01273 
01274      \ref tablethree
01275 
01276      \ref tablefour
01277 
01278      \ref tablefive
01279 
01280      \ref tablesix
01281 
01282      \ref tableseven
01283  
01284   \page building Building & Installing
01285  
01286 <h1> Standard Installation Instructions </h1>
01287 
01288   MOAB uses an autoconf and libtool-based build process by default.  The procedure used to build MOAB from scratch depends on whether the source code was obtained from a “tarball” or directly from the Subversion repository.  Assuming the latter, the following steps should be executed for building and installing MOAB:
01289   - Locate and build any required dependencies.  MOAB can be built with no dependencies on other libraries; this may be useful for applications only needing basic mesh representation and not needing to export mesh to formats implemented in other libraries.  MOAB’s native save/restore capability is built on HDF5-based files; applications needing to save and restore files from MOAB reliably should use this library.  MOAB also uses ExodusII, a netCDF-based file format developed at Sandia National Laboratories [10].  Applications needing to execute these tests should also build netCDF.  Note that MOAB uses netCDF’s C++ interface, which is not enabled by default in netCDF but can be enabled using the “–enable-cxx” option to netCDF’s configure script.
01290   - Unpack source code into <moab>, and change current working directory to that location.
01291   - Execute “autoreconf –fi”.
01292   - Run configure script, by executing “./configure <options>”.  Recommended options:
01293        -# –prefix=<install_dir>: directory below which MOAB library and include files will be installed; can either be the directory used for MOAB source (<moab> from step 1), or a different directory.
01294        -# –hdf5-dir=<hdf5_dir>: directory whose “include” and “lib” subdirectories hold HDF5 include and library, respectively.  MOAB uses HDF5 for its native save/restore format (see Section 4.6.1).
01295        -# –netcdf-dir=: directory whose “include” and “lib” subdirectories hold netCDF include and library, respectively.  MOAB uses netCDF-based files for many of its build tests.  If the location of netCDF cannot be found, MOAB’s build tests will not function properly, but MOAB will still be usable.
01296        .
01297   - Run “make check”; this runs a series of build tests, to verify that the MOAB build was successful.  Note this check will fail if netCDF is not used, but MOAB itself will still be usable from applications.
01298   - Run “make install”; this copies include files and libraries to subdirectories of the directory specified in the “prefix” option.
01299   .
01300 
01301 These steps are sufficient for building MOAB against HDF5 and netCDF.  By default, a small number of standard MOAB-based applications are also built, including mbconvert (a utility for reading and writing files), mbsize (for querying basic information about a mesh), and the iMesh interface (see Section 7).  Other utilities can be enabled using various other options to the configure script; for a complete list of build options, execute “./configure –help”.
01302 
01303 <h1> Other MOAB Components </h1>
01304 
01305 <h2> PyMOAB </h2>
01306 
01307   To install the PyMOAB module, add "--enable-pymoab" in the configuration step described above. PyMOAB will build using whatever version of Python is available on the system when the configuration command is executed. It requires the following packages to operate properly:
01308     - Python's <a href=https://pypi.python.org/pypi/setuptools> setuptools </a> and <a href=https://docs.python.org/3/library/distutils.html> distutils </a> for building, linking, and distribution
01309     - <a href=https://docs.python.org/3/library/distutils.html> NumPy </a> for management of data being passed into and out of the MOAB database
01310 
01311   PyMOAB can be installed with a custom prefix path which can be set adding "PYMOAB_PREFIX=/path/to/custom/install/location" to the configuration command.
01312 
01313  */
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines