Skip to content
Snippets Groups Projects
Select Git revision
  • 3b1b0ca2c775554184952dc0b23c89905c53232c
  • master default protected
2 results

main.c

Blame
  • gmsh.h 109.16 KiB
    // Gmsh - Copyright (C) 1997-2019 C. Geuzaine, J.-F. Remacle
    //
    // See the LICENSE.txt file for license information. Please report all
    // issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
    
    #ifndef GMSH_H
    #define GMSH_H
    
    // This file defines the Gmsh C++ API (v4.4).
    //
    // Do not edit it directly: it is automatically generated by `api/gen.py'.
    //
    // By design, the Gmsh C++ API is purely functional, and only uses elementary
    // types from the standard library. See `demos/api' for examples.
    
    #if defined(_MSC_VER) && !defined(_USE_MATH_DEFINES)
    #define _USE_MATH_DEFINES
    #endif
    
    #include <cmath>
    #include <vector>
    #include <string>
    #include <utility>
    
    #define GMSH_API_VERSION "4.4"
    #define GMSH_API_VERSION_MAJOR 4
    #define GMSH_API_VERSION_MINOR 4
    
    #if defined(GMSH_DLL)
    #if defined(GMSH_DLL_EXPORT)
    #define GMSH_API __declspec(dllexport)
    #else
    #define GMSH_API __declspec(dllimport)
    #endif
    #else
    #define GMSH_API
    #endif
    
    namespace gmsh {
    
      // A geometrical entity in the Gmsh API is represented by two integers: its
      // dimension (dim = 0, 1, 2 or 3) and its tag (its unique, strictly positive
      // identifier). When dealing with multiple geometrical entities of possibly
      // different dimensions, the entities are packed as a vector of (dim, tag)
      // integer pairs.
      typedef std::vector<std::pair<int, int> > vectorpair;
    
    }
    
    namespace gmsh { // Top-level functions
    
      // Initialize Gmsh. This must be called before any call to the other functions in
      // the API. If `argc' and `argv' (or just `argv' in Python or Julia) are
      // provided, they will be handled in the same way as the command line arguments
      // in the Gmsh app. If `readConfigFiles' is set, read system Gmsh configuration
      // files (gmshrc and gmsh-options).
      GMSH_API void initialize(int argc = 0, char ** argv = 0,
                               const bool readConfigFiles = true);
    
      // Finalize Gmsh. This must be called when you are done using the Gmsh API.
      GMSH_API void finalize();
    
      // Open a file. Equivalent to the `File->Open' menu in the Gmsh app. Handling of
      // the file depends on its extension and/or its contents: opening a file with
      // model data will create a new model.
      GMSH_API void open(const std::string & fileName);
    
      // Merge a file. Equivalent to the `File->Merge' menu in the Gmsh app. Handling
      // of the file depends on its extension and/or its contents. Merging a file with
      // model data will add the data to the current model.
      GMSH_API void merge(const std::string & fileName);
    
      // Write a file. The export format is determined by the file extension.
      GMSH_API void write(const std::string & fileName);
    
      // Clear all loaded models and post-processing data, and add a new empty model.
      GMSH_API void clear();
    
      namespace option { // Option handling functions
    
        // Set a numerical option to `value'. `name' is of the form "category.option"
        // or "category[num].option". Available categories and options are listed in
        // the Gmsh reference manual.
        GMSH_API void setNumber(const std::string & name,
                                const double value);
    
        // Get the `value' of a numerical option. `name' is of the form
        // "category.option" or "category[num].option". Available categories and
        // options are listed in the Gmsh reference manual.
        GMSH_API void getNumber(const std::string & name,
                                double & value);
    
        // Set a string option to `value'. `name' is of the form "category.option" or
        // "category[num].option". Available categories and options are listed in the
        // Gmsh reference manual.
        GMSH_API void setString(const std::string & name,
                                const std::string & value);
    
        // Get the `value' of a string option. `name' is of the form "category.option"
        // or "category[num].option". Available categories and options are listed in
        // the Gmsh reference manual.
        GMSH_API void getString(const std::string & name,
                                std::string & value);
    
        // Set a color option to the RGBA value (`r', `g', `b', `a'), where where `r',
        // `g', `b' and `a' should be integers between 0 and 255. `name' is of the form
        // "category.option" or "category[num].option". Available categories and
        // options are listed in the Gmsh reference manual, with the "Color." middle
        // string removed.
        GMSH_API void setColor(const std::string & name,
                               const int r,
                               const int g,
                               const int b,
                               const int a = 0);
    
        // Get the `r', `g', `b', `a' value of a color option. `name' is of the form
        // "category.option" or "category[num].option". Available categories and
        // options are listed in the Gmsh reference manual, with the "Color." middle
        // string removed.
        GMSH_API void getColor(const std::string & name,
                               int & r,
                               int & g,
                               int & b,
                               int & a);
    
      } // namespace option
    
      namespace model { // Model functions
    
        // Add a new model, with name `name', and set it as the current model.
        GMSH_API void add(const std::string & name);
    
        // Remove the current model.
        GMSH_API void remove();
    
        // List the names of all models.
        GMSH_API void list(std::vector<std::string> & names);
    
        // Set the current model to the model with name `name'. If several models have
        // the same name, select the one that was added first.
        GMSH_API void setCurrent(const std::string & name);
    
        // Get all the entities in the current model. If `dim' is >= 0, return only the
        // entities of the specified dimension (e.g. points if `dim' == 0). The
        // entities are returned as a vector of (dim, tag) integer pairs.
        GMSH_API void getEntities(gmsh::vectorpair & dimTags,
                                  const int dim = -1);
    
        // Set the name of the entity of dimension `dim' and tag `tag'.
        GMSH_API void setEntityName(const int dim,
                                    const int tag,
                                    const std::string & name);
    
        // Get the name of the entity of dimension `dim' and tag `tag'.
        GMSH_API void getEntityName(const int dim,
                                    const int tag,
                                    std::string & name);
    
        // Get all the physical groups in the current model. If `dim' is >= 0, return
        // only the entities of the specified dimension (e.g. physical points if `dim'
        // == 0). The entities are returned as a vector of (dim, tag) integer pairs.
        GMSH_API void getPhysicalGroups(gmsh::vectorpair & dimTags,
                                        const int dim = -1);
    
        // Get the tags of the model entities making up the physical group of dimension
        // `dim' and tag `tag'.
        GMSH_API void getEntitiesForPhysicalGroup(const int dim,
                                                  const int tag,
                                                  std::vector<int> & tags);
    
        // Get the tags of the physical groups (if any) to which the model entity of
        // dimension `dim' and tag `tag' belongs.
        GMSH_API void getPhysicalGroupsForEntity(const int dim,
                                                 const int tag,
                                                 std::vector<int> & physicalTags);
    
        // Add a physical group of dimension `dim', grouping the model entities with
        // tags `tags'. Return the tag of the physical group, equal to `tag' if `tag'
        // is positive, or a new tag if `tag' < 0.
        GMSH_API int addPhysicalGroup(const int dim,
                                      const std::vector<int> & tags,
                                      const int tag = -1);
    
        // Set the name of the physical group of dimension `dim' and tag `tag'.
        GMSH_API void setPhysicalName(const int dim,
                                      const int tag,
                                      const std::string & name);
    
        // Get the name of the physical group of dimension `dim' and tag `tag'.
        GMSH_API void getPhysicalName(const int dim,
                                      const int tag,
                                      std::string & name);
    
        // Get the boundary of the model entities `dimTags'. Return in `outDimTags' the
        // boundary of the individual entities (if `combined' is false) or the boundary
        // of the combined geometrical shape formed by all input entities (if
        // `combined' is true). Return tags multiplied by the sign of the boundary
        // entity if `oriented' is true. Apply the boundary operator recursively down
        // to dimension 0 (i.e. to points) if `recursive' is true.
        GMSH_API void getBoundary(const gmsh::vectorpair & dimTags,
                                  gmsh::vectorpair & outDimTags,
                                  const bool combined = true,
                                  const bool oriented = true,
                                  const bool recursive = false);
    
        // Get the model entities in the bounding box defined by the two points
        // (`xmin', `ymin', `zmin') and (`xmax', `ymax', `zmax'). If `dim' is >= 0,
        // return only the entities of the specified dimension (e.g. points if `dim' ==
        // 0).
        GMSH_API void getEntitiesInBoundingBox(const double xmin,
                                               const double ymin,
                                               const double zmin,
                                               const double xmax,
                                               const double ymax,
                                               const double zmax,
                                               gmsh::vectorpair & tags,
                                               const int dim = -1);
    
        // Get the bounding box (`xmin', `ymin', `zmin'), (`xmax', `ymax', `zmax') of
        // the model entity of dimension `dim' and tag `tag'. If `dim' and `tag' are
        // negative, get the bounding box of the whole model.
        GMSH_API void getBoundingBox(const int dim,
                                     const int tag,
                                     double & xmin,
                                     double & ymin,
                                     double & zmin,
                                     double & xmax,
                                     double & ymax,
                                     double & zmax);
    
        // Get the geometrical dimension of the current model.
        GMSH_API int getDimension();
    
        // Add a discrete model entity (defined by a mesh) of dimension `dim' in the
        // current model. Return the tag of the new discrete entity, equal to `tag' if
        // `tag' is positive, or a new tag if `tag' < 0. `boundary' specifies the tags
        // of the entities on the boundary of the discrete entity, if any. Specifying
        // `boundary' allows Gmsh to construct the topology of the overall model.
        GMSH_API int addDiscreteEntity(const int dim,
                                       const int tag = -1,
                                       const std::vector<int> & boundary = std::vector<int>());
    
        // Remove the entities `dimTags' of the current model. If `recursive' is true,
        // remove all the entities on their boundaries, down to dimension 0.
        GMSH_API void removeEntities(const gmsh::vectorpair & dimTags,
                                     const bool recursive = false);
    
        // Remove the entity name `name' from the current model.
        GMSH_API void removeEntityName(const std::string & name);
    
        // Remove the physical groups `dimTags' of the current model. If `dimTags' is
        // empty, remove all groups.
        GMSH_API void removePhysicalGroups(const gmsh::vectorpair & dimTags = gmsh::vectorpair());
    
        // Remove the physical name `name' from the current model.
        GMSH_API void removePhysicalName(const std::string & name);
    
        // Get the type of the entity of dimension `dim' and tag `tag'.
        GMSH_API void getType(const int dim,
                              const int tag,
                              std::string & entityType);
    
        // In a partitioned model, get the parent of the entity of dimension `dim' and
        // tag `tag', i.e. from which the entity is a part of, if any. `parentDim' and
        // `parentTag' are set to -1 if the entity has no parent.
        GMSH_API void getParent(const int dim,
                                const int tag,
                                int & parentDim,
                                int & parentTag);
    
        // In a partitioned model, return the tags of the partition(s) to which the
        // entity belongs.
        GMSH_API void getPartitions(const int dim,
                                    const int tag,
                                    std::vector<int> & partitions);
    
        // Evaluate the parametrization of the entity of dimension `dim' and tag `tag'
        // at the parametric coordinates `parametricCoord'. Only valid for `dim' equal
        // to 0 (with empty `parametricCoord'), 1 (with `parametricCoord' containing
        // parametric coordinates on the curve) or 2 (with `parametricCoord' containing
        // pairs of u, v parametric coordinates on the surface, concatenated: [p1u,
        // p1v, p2u, ...]). Return triplets of x, y, z coordinates in `points',
        // concatenated: [p1x, p1y, p1z, p2x, ...].
        GMSH_API void getValue(const int dim,
                               const int tag,
                               const std::vector<double> & parametricCoord,
                               std::vector<double> & points);
    
        // Evaluate the derivative of the parametrization of the entity of dimension
        // `dim' and tag `tag' at the parametric coordinates `parametricCoord'. Only
        // valid for `dim' equal to 1 (with `parametricCoord' containing parametric
        // coordinates on the curve) or 2 (with `parametricCoord' containing pairs of
        // u, v parametric coordinates on the surface, concatenated: [p1u, p1v, p2u,
        // ...]). For `dim' equal to 1 return the x, y, z components of the derivative
        // with respect to u [d1ux, d1uy, d1uz, d2ux, ...]; for `dim' equal to 2 return
        // the x, y, z components of the derivate with respect to u and v: [d1ux, d1uy,
        // d1uz, d1vx, d1vy, d1vz, d2ux, ...].
        GMSH_API void getDerivative(const int dim,
                                    const int tag,
                                    const std::vector<double> & parametricCoord,
                                    std::vector<double> & derivatives);
    
        // Evaluate the (maximum) curvature of the entity of dimension `dim' and tag
        // `tag' at the parametric coordinates `parametricCoord'. Only valid for `dim'
        // equal to 1 (with `parametricCoord' containing parametric coordinates on the
        // curve) or 2 (with `parametricCoord' containing pairs of u, v parametric
        // coordinates on the surface, concatenated: [p1u, p1v, p2u, ...]).
        GMSH_API void getCurvature(const int dim,
                                   const int tag,
                                   const std::vector<double> & parametricCoord,
                                   std::vector<double> & curvatures);
    
        // Evaluate the principal curvatures of the surface with tag `tag' at the
        // parametric coordinates `parametricCoord', as well as their respective
        // directions. `parametricCoord' are given by pair of u and v coordinates,
        // concatenated: [p1u, p1v, p2u, ...].
        GMSH_API void getPrincipalCurvatures(const int tag,
                                             const std::vector<double> & parametricCoord,
                                             std::vector<double> & curvatureMax,
                                             std::vector<double> & curvatureMin,
                                             std::vector<double> & directionMax,
                                             std::vector<double> & directionMin);
    
        // Get the normal to the surface with tag `tag' at the parametric coordinates
        // `parametricCoord'. `parametricCoord' are given by pairs of u and v
        // coordinates, concatenated: [p1u, p1v, p2u, ...]. `normals' are returned as
        // triplets of x, y, z components, concatenated: [n1x, n1y, n1z, n2x, ...].
        GMSH_API void getNormal(const int tag,
                                const std::vector<double> & parametricCoord,
                                std::vector<double> & normals);
    
        // Set the visibility of the model entities `dimTags' to `value'. Apply the
        // visibility setting recursively if `recursive' is true.
        GMSH_API void setVisibility(const gmsh::vectorpair & dimTags,
                                    const int value,
                                    const bool recursive = false);
    
        // Get the visibility of the model entity of dimension `dim' and tag `tag'.
        GMSH_API void getVisibility(const int dim,
                                    const int tag,
                                    int & value);
    
        // Set the color of the model entities `dimTags' to the RGBA value (`r', `g',
        // `b', `a'), where `r', `g', `b' and `a' should be integers between 0 and 255.
        // Apply the color setting recursively if `recursive' is true.
        GMSH_API void setColor(const gmsh::vectorpair & dimTags,
                               const int r,
                               const int g,
                               const int b,
                               const int a = 0,
                               const bool recursive = false);
    
        // Get the color of the model entity of dimension `dim' and tag `tag'.
        GMSH_API void getColor(const int dim,
                               const int tag,
                               int & r,
                               int & g,
                               int & b,
                               int & a);
    
        // Set the `x', `y', `z' coordinates of a geometrical point.
        GMSH_API void setCoordinates(const int tag,
                                     const double x,
                                     const double y,
                                     const double z);
    
        namespace mesh { // Mesh functions
    
          // Generate a mesh of the current model, up to dimension `dim' (0, 1, 2 or
          // 3).
          GMSH_API void generate(const int dim = 3);
    
          // Partition the mesh of the current model into `numPart' partitions.
          GMSH_API void partition(const int numPart);
    
          // Unpartition the mesh of the current model.
          GMSH_API void unpartition();
    
          // Optimize the mesh of the current model using `method' (empty for default
          // tetrahedral mesh optimizer, "Netgen" for Netgen optimizer, "HighOrder" for
          // direct high-order mesh optimizer, "HighOrderElastic" for high-order
          // elastic smoother).
          GMSH_API void optimize(const std::string & method);
    
          // Recombine the mesh of the current model.
          GMSH_API void recombine();
    
          // Refine the mesh of the current model by uniformly splitting the elements.
          GMSH_API void refine();
    
          // Smooth the mesh of the current model.
          GMSH_API void smooth();
    
          // Set the order of the elements in the mesh of the current model to `order'.
          GMSH_API void setOrder(const int order);
    
          // Get the last entities (if any) where a meshing error occurred. Currently
          // only populated by the new 3D meshing algorithms.
          GMSH_API void getLastEntityError(gmsh::vectorpair & dimTags);
    
          // Get the last nodes (if any) where a meshing error occurred. Currently only
          // populated by the new 3D meshing algorithms.
          GMSH_API void getLastNodeError(std::vector<std::size_t> & nodeTags);
    
          // Clear the mesh, i.e. delete all the nodes and elements.
          GMSH_API void clear();
    
          // Get the nodes classified on the entity of dimension `dim' and tag `tag'.
          // If `tag' < 0, get the nodes for all entities of dimension `dim'. If `dim'
          // and `tag' are negative, get all the nodes in the mesh. `nodeTags' contains
          // the node tags (their unique, strictly positive identification numbers).
          // `coord' is a vector of length 3 times the length of `nodeTags' that
          // contains the x, y, z coordinates of the nodes, concatenated: [n1x, n1y,
          // n1z, n2x, ...]. If `dim' >= 0 and `returnParamtricCoord' is set,
          // `parametricCoord' contains the parametric coordinates ([u1, u2, ...] or
          // [u1, v1, u2, ...]) of the nodes, if available. The length of
          // `parametricCoord' can be 0 or `dim' times the length of `nodeTags'. If
          // `includeBoundary' is set, also return the nodes classified on the boundary
          // of the entity (which will be reparametrized on the entity if `dim' >= 0 in
          // order to compute their parametric coordinates).
          GMSH_API void getNodes(std::vector<std::size_t> & nodeTags,
                                 std::vector<double> & coord,
                                 std::vector<double> & parametricCoord,
                                 const int dim = -1,
                                 const int tag = -1,
                                 const bool includeBoundary = false,
                                 const bool returnParametricCoord = true);
    
          // Get the nodes classified on the entity of tag `tag', for all the elements
          // of type `elementType'. The other arguments are treated as in `getNodes'.
          GMSH_API void getNodesByElementType(const int elementType,
                                              std::vector<std::size_t> & nodeTags,
                                              std::vector<double> & coord,
                                              std::vector<double> & parametricCoord,
                                              const int tag = -1,
                                              const bool returnParametricCoord = true);
    
          // Get the coordinates and the parametric coordinates (if any) of the node
          // with tag `tag'. This is a sometimes useful but inefficient way of
          // accessing nodes, as it relies on a cache stored in the model. For large
          // meshes all the nodes in the model should be numbered in a continuous
          // sequence of tags from 1 to N to maintain reasonable performance (in this
          // case the internal cache is based on a vector; otherwise it uses a map).
          GMSH_API void getNode(const std::size_t nodeTag,
                                std::vector<double> & coord,
                                std::vector<double> & parametricCoord);
    
          // Rebuild the node cache.
          GMSH_API void rebuildNodeCache(const bool onlyIfNecessary = true);
    
          // Get the nodes from all the elements belonging to the physical group of
          // dimension `dim' and tag `tag'. `nodeTags' contains the node tags; `coord'
          // is a vector of length 3 times the length of `nodeTags' that contains the
          // x, y, z coordinates of the nodes, concatenated: [n1x, n1y, n1z, n2x, ...].
          GMSH_API void getNodesForPhysicalGroup(const int dim,
                                                 const int tag,
                                                 std::vector<std::size_t> & nodeTags,
                                                 std::vector<double> & coord);
    
          // Add nodes classified on the model entity of dimension `dim' and tag `tag'.
          // `nodeTags' contains the node tags (their unique, strictly positive
          // identification numbers). `coord' is a vector of length 3 times the length
          // of `nodeTags' that contains the x, y, z coordinates of the nodes,
          // concatenated: [n1x, n1y, n1z, n2x, ...]. The optional `parametricCoord'
          // vector contains the parametric coordinates of the nodes, if any. The
          // length of `parametricCoord' can be 0 or `dim' times the length of
          // `nodeTags'. If the `nodeTags' vector is empty, new tags are automatically
          // assigned to the nodes.
          GMSH_API void addNodes(const int dim,
                                 const int tag,
                                 const std::vector<std::size_t> & nodeTags,
                                 const std::vector<double> & coord,
                                 const std::vector<double> & parametricCoord = std::vector<double>());
    
          // Reclassify all nodes on their associated model entity, based on the
          // elements. Can be used when importing nodes in bulk (e.g. by associating
          // them all to a single volume), to reclassify them correctly on model
          // surfaces, curves, etc. after the elements have been set.
          GMSH_API void reclassifyNodes();
    
          // Relocate the nodes classified on the entity of dimension `dim' and tag
          // `tag' using their parametric coordinates. If `tag' < 0, relocate the nodes
          // for all entities of dimension `dim'. If `dim' and `tag' are negative,
          // relocate all the nodes in the mesh.
          GMSH_API void relocateNodes(const int dim = -1,
                                      const int tag = -1);
    
          // Get the elements classified on the entity of dimension `dim' and tag
          // `tag'. If `tag' < 0, get the elements for all entities of dimension `dim'.
          // If `dim' and `tag' are negative, get all the elements in the mesh.
          // `elementTypes' contains the MSH types of the elements (e.g. `2' for 3-node
          // triangles: see `getElementProperties' to obtain the properties for a given
          // element type). `elementTags' is a vector of the same length as
          // `elementTypes'; each entry is a vector containing the tags (unique,
          // strictly positive identifiers) of the elements of the corresponding type.
          // `nodeTags' is also a vector of the same length as `elementTypes'; each
          // entry is a vector of length equal to the number of elements of the given
          // type times the number N of nodes for this type of element, that contains
          // the node tags of all the elements of the given type, concatenated: [e1n1,
          // e1n2, ..., e1nN, e2n1, ...].
          GMSH_API void getElements(std::vector<int> & elementTypes,
                                    std::vector<std::vector<std::size_t> > & elementTags,
                                    std::vector<std::vector<std::size_t> > & nodeTags,
                                    const int dim = -1,
                                    const int tag = -1);
    
          // Get the type and node tags of the element with tag `tag'. This is a
          // sometimes useful but inefficient way of accessing elements, as it relies
          // on a cache stored in the model. For large meshes all the elements in the
          // model should be numbered in a continuous sequence of tags from 1 to N to
          // maintain reasonable performance (in this case the internal cache is based
          // on a vector; otherwise it uses a map).
          GMSH_API void getElement(const std::size_t elementTag,
                                   int & elementType,
                                   std::vector<std::size_t> & nodeTags);
    
          // Search the mesh for an element located at coordinates (`x', `y', `z').
          // This is a sometimes useful but inefficient way of accessing elements, as
          // it relies on a search in a spatial octree. If an element is found, return
          // its tag, type and node tags, as well as the local coordinates (`u', `v',
          // `w') within the element corresponding to search location. If `dim' is >=
          // 0, only search for elements of the given dimension. If `strict' is not
          // set, use a tolerance to find elements near the search location.
          GMSH_API void getElementByCoordinates(const double x,
                                                const double y,
                                                const double z,
                                                std::size_t & elementTag,
                                                int & elementType,
                                                std::vector<std::size_t> & nodeTags,
                                                double & u,
                                                double & v,
                                                double & w,
                                                const int dim = -1,
                                                const bool strict = false);
    
          // Get the types of elements in the entity of dimension `dim' and tag `tag'.
          // If `tag' < 0, get the types for all entities of dimension `dim'. If `dim'
          // and `tag' are negative, get all the types in the mesh.
          GMSH_API void getElementTypes(std::vector<int> & elementTypes,
                                        const int dim = -1,
                                        const int tag = -1);
    
          // Return an element type given its family name `familyName' ("point",
          // "line", "triangle", "quadrangle", "tetrahedron", "pyramid", "prism",
          // "hexahedron") and polynomial order `order'. If `serendip' is true, return
          // the corresponding serendip element type (element without interior nodes).
          GMSH_API int getElementType(const std::string & familyName,
                                      const int order,
                                      const bool serendip = false);
    
          // Get the properties of an element of type `elementType': its name
          // (`elementName'), dimension (`dim'), order (`order'), number of nodes
          // (`numNodes') and coordinates of the nodes in the reference element
          // (`nodeCoord' vector, of length `dim' times `numNodes').
          GMSH_API void getElementProperties(const int elementType,
                                             std::string & elementName,
                                             int & dim,
                                             int & order,
                                             int & numNodes,
                                             std::vector<double> & nodeCoord);
    
          // Get the elements of type `elementType' classified on the entity of tag
          // `tag'. If `tag' < 0, get the elements for all entities. `elementTags' is a
          // vector containing the tags (unique, strictly positive identifiers) of the
          // elements of the corresponding type. `nodeTags' is a vector of length equal
          // to the number of elements of the given type times the number N of nodes
          // for this type of element, that contains the node tags of all the elements
          // of the given type, concatenated: [e1n1, e1n2, ..., e1nN, e2n1, ...]. If
          // `numTasks' > 1, only compute and return the part of the data indexed by
          // `task'.
          GMSH_API void getElementsByType(const int elementType,
                                          std::vector<std::size_t> & elementTags,
                                          std::vector<std::size_t> & nodeTags,
                                          const int tag = -1,
                                          const std::size_t task = 0,
                                          const std::size_t numTasks = 1);
    
          // Preallocate data before calling `getElementsByType' with `numTasks' > 1.
          // For C and C++ only.
          GMSH_API void preallocateElementsByType(const int elementType,
                                                  const bool elementTag,
                                                  const bool nodeTag,
                                                  std::vector<std::size_t> & elementTags,
                                                  std::vector<std::size_t> & nodeTags,
                                                  const int tag = -1);
    
          // Add elements classified on the entity of dimension `dim' and tag `tag'.
          // `types' contains the MSH types of the elements (e.g. `2' for 3-node
          // triangles: see the Gmsh reference manual). `elementTags' is a vector of
          // the same length as `types'; each entry is a vector containing the tags
          // (unique, strictly positive identifiers) of the elements of the
          // corresponding type. `nodeTags' is also a vector of the same length as
          // `types'; each entry is a vector of length equal to the number of elements
          // of the given type times the number N of nodes per element, that contains
          // the node tags of all the elements of the given type, concatenated: [e1n1,
          // e1n2, ..., e1nN, e2n1, ...].
          GMSH_API void addElements(const int dim,
                                    const int tag,
                                    const std::vector<int> & elementTypes,
                                    const std::vector<std::vector<std::size_t> > & elementTags,
                                    const std::vector<std::vector<std::size_t> > & nodeTags);
    
          // Add elements of type `elementType' classified on the entity of tag `tag'.
          // `elementTags' contains the tags (unique, strictly positive identifiers) of
          // the elements of the corresponding type. `nodeTags' is a vector of length
          // equal to the number of elements times the number N of nodes per element,
          // that contains the node tags of all the elements, concatenated: [e1n1,
          // e1n2, ..., e1nN, e2n1, ...]. If the `elementTag' vector is empty, new tags
          // are automatically assigned to the elements.
          GMSH_API void addElementsByType(const int tag,
                                          const int elementType,
                                          const std::vector<std::size_t> & elementTags,
                                          const std::vector<std::size_t> & nodeTags);
    
          // Get the numerical quadrature information for the given element type
          // `elementType' and integration rule `integrationType' (e.g. "Gauss4" for a
          // Gauss quadrature suited for integrating 4th order polynomials).
          // `integrationPoints' contains the u, v, w coordinates of the G integration
          // points in the reference element: [g1u, g1v, g1w, ..., gGu, gGv, gGw].
          // `integrationWeigths' contains the associated weights: [g1q, ..., gGq].
          GMSH_API void getIntegrationPoints(const int elementType,
                                             const std::string & integrationType,
                                             std::vector<double> & integrationPoints,
                                             std::vector<double> & integrationWeights);
    
          // Get the Jacobians of all the elements of type `elementType' classified on
          // the entity of tag `tag', at the G integration points `integrationPoints'
          // given as concatenated triplets of coordinates in the reference element
          // [g1u, g1v, g1w, ..., gGu, gGv, gGw]. Data is returned by element, with
          // elements in the same order as in `getElements' and `getElementsByType'.
          // `jacobians' contains for each element the 9 entries of the 3x3 Jacobian
          // matrix at each integration point. The matrix is returned by column:
          // [e1g1Jxu, e1g1Jyu, e1g1Jzu, e1g1Jxv, ..., e1g1Jzw, e1g2Jxu, ..., e1gGJzw,
          // e2g1Jxu, ...], with Jxu=dx/du, Jyu=dy/du, etc. `determinants' contains for
          // each element the determinant of the Jacobian matrix at each integration
          // point: [e1g1, e1g2, ... e1gG, e2g1, ...]. `points' contains for each
          // element the x, y, z coordinates of the integration points. If `tag' < 0,
          // get the Jacobian data for all entities. If `numTasks' > 1, only compute
          // and return the part of the data indexed by `task'.
          GMSH_API void getJacobians(const int elementType,
                                     const std::vector<double> & integrationPoints,
                                     std::vector<double> & jacobians,
                                     std::vector<double> & determinants,
                                     std::vector<double> & points,
                                     const int tag = -1,
                                     const std::size_t task = 0,
                                     const std::size_t numTasks = 1);
    
          // Preallocate data before calling `getJacobians' with `numTasks' > 1. For C
          // and C++ only.
          GMSH_API void preallocateJacobians(const int elementType,
                                             const int numIntegrationPoints,
                                             const bool jacobian,
                                             const bool determinant,
                                             const bool point,
                                             std::vector<double> & jacobians,
                                             std::vector<double> & determinants,
                                             std::vector<double> & points,
                                             const int tag = -1);
    
          // Get the basis functions of the element of type `elementType' at the
          // integration points `integrationPoints' (given as concatenated triplets of
          // coordinates in the reference element [g1u, g1v, g1w, ..., gGu, gGv, gGw]),
          // for the function space `functionSpaceType' (e.g. "Lagrange" or
          // "GradLagrange" for Lagrange basis functions or their gradient, in the u,
          // v, w coordinates of the reference element). `numComponents' returns the
          // number C of components of a basis function. `basisFunctions' returns the
          // value of the N basis functions at the integration points, i.e. [g1f1,
          // g1f2, ..., g1fN, g2f1, ...] when C == 1 or [g1f1u, g1f1v, g1f1w, g1f2u,
          // ..., g1fNw, g2f1u, ...] when C == 3.
          GMSH_API void getBasisFunctions(const int elementType,
                                          const std::vector<double> & integrationPoints,
                                          const std::string & functionSpaceType,
                                          int & numComponents,
                                          std::vector<double> & basisFunctions);
    
          // Get the element-dependent basis functions of the elements of type
          // `elementType' in the entity of tag `tag'at the integration points
          // `integrationPoints' (given as concatenated triplets of coordinates in the
          // reference element [g1u, g1v, g1w, ..., gGu, gGv, gGw]), for the function
          // space `functionSpaceType' (e.g. "H1Legendre3" or "GradH1Legendre3" for 3rd
          // order hierarchical H1 Legendre functions or their gradient, in the u, v, w
          // coordinates of the reference elements). `numComponents' returns the number
          // C of components of a basis function. `numBasisFunctions' returns the
          // number N of basis functions per element. `basisFunctions' returns the
          // value of the basis functions at the integration points for each element:
          // [e1g1f1,..., e1g1fN, e1g2f1,..., e2g1f1, ...] when C == 1 or [e1g1f1u,
          // e1g1f1v,..., e1g1fNw, e1g2f1u,..., e2g1f1u, ...]. Warning: this is an
          // experimental feature and will probably change in a future release.
          GMSH_API void getBasisFunctionsForElements(const int elementType,
                                                     const std::vector<double> & integrationPoints,
                                                     const std::string & functionSpaceType,
                                                     int & numComponents,
                                                     int & numFunctionsPerElements,
                                                     std::vector<double> & basisFunctions,
                                                     const int tag = -1);
    
          // Generate the `keys' for the elements of type `elementType' in the entity
          // of tag `tag', for the `functionSpaceType' function space. Each key
          // uniquely identifies a basis function in the function space. If
          // `returnCoord' is set, the `coord' vector contains the x, y, z coordinates
          // locating basis functions for sorting purposes. Warning: this is an
          // experimental feature and will probably change in a future release.
          GMSH_API void getKeysForElements(const int elementType,
                                           const std::string & functionSpaceType,
                                           gmsh::vectorpair & keys,
                                           std::vector<double> & coord,
                                           const int tag = -1,
                                           const bool returnCoord = true);
    
          // Get information about the `keys'. Warning: this is an experimental feature
          // and will probably change in a future release.
          GMSH_API void getInformationForElements(const gmsh::vectorpair & keys,
                                                  gmsh::vectorpair & info,
                                                  const int order,
                                                  const int elementType);
    
          // Precomputes the basis functions corresponding to `elementType'.
          GMSH_API void precomputeBasisFunctions(const int elementType);
    
          // Get the barycenters of all elements of type `elementType' classified on
          // the entity of tag `tag'. If `primary' is set, only the primary nodes of
          // the elements are taken into account for the barycenter calculation. If
          // `fast' is set, the function returns the sum of the primary node
          // coordinates (without normalizing by the number of nodes). If `tag' < 0,
          // get the barycenters for all entities. If `numTasks' > 1, only compute and
          // return the part of the data indexed by `task'.
          GMSH_API void getBarycenters(const int elementType,
                                       const int tag,
                                       const bool fast,
                                       const bool primary,
                                       std::vector<double> & barycenters,
                                       const std::size_t task = 0,
                                       const std::size_t numTasks = 1);
    
          // Preallocate data before calling `getBarycenters' with `numTasks' > 1. For
          // C and C++ only.
          GMSH_API void preallocateBarycenters(const int elementType,
                                               std::vector<double> & barycenters,
                                               const int tag = -1);
    
          // Get the nodes on the edges of all elements of type `elementType'
          // classified on the entity of tag `tag'. `nodeTags' contains the node tags
          // of the edges for all the elements: [e1a1n1, e1a1n2, e1a2n1, ...]. Data is
          // returned by element, with elements in the same order as in `getElements'
          // and `getElementsByType'. If `primary' is set, only the primary (begin/end)
          // nodes of the edges are returned. If `tag' < 0, get the edge nodes for all
          // entities. If `numTasks' > 1, only compute and return the part of the data
          // indexed by `task'.
          GMSH_API void getElementEdgeNodes(const int elementType,
                                            std::vector<std::size_t> & nodeTags,
                                            const int tag = -1,
                                            const bool primary = false,
                                            const std::size_t task = 0,
                                            const std::size_t numTasks = 1);
    
          // Get the nodes on the faces of type `faceType' (3 for triangular faces, 4
          // for quadrangular faces) of all elements of type `elementType' classified
          // on the entity of tag `tag'. `nodeTags' contains the node tags of the faces
          // for all elements: [e1f1n1, ..., e1f1nFaceType, e1f2n1, ...]. Data is
          // returned by element, with elements in the same order as in `getElements'
          // and `getElementsByType'. If `primary' is set, only the primary (corner)
          // nodes of the faces are returned. If `tag' < 0, get the face nodes for all
          // entities. If `numTasks' > 1, only compute and return the part of the data
          // indexed by `task'.
          GMSH_API void getElementFaceNodes(const int elementType,
                                            const int faceType,
                                            std::vector<std::size_t> & nodeTags,
                                            const int tag = -1,
                                            const bool primary = false,
                                            const std::size_t task = 0,
                                            const std::size_t numTasks = 1);
    
          // Get the ghost elements `elementTags' and their associated `partitions'
          // stored in the ghost entity of dimension `dim' and tag `tag'.
          GMSH_API void getGhostElements(const int dim,
                                         const int tag,
                                         std::vector<std::size_t> & elementTags,
                                         std::vector<int> & partitions);
    
          // Set a mesh size constraint on the model entities `dimTags'. Currently only
          // entities of dimension 0 (points) are handled.
          GMSH_API void setSize(const gmsh::vectorpair & dimTags,
                                const double size);
    
          // Set a transfinite meshing constraint on the curve `tag', with `numNodes'
          // nodes distributed according to `meshType' and `coef'. Currently supported
          // types are "Progression" (geometrical progression with power `coef') and
          // "Bump" (refinement toward both extremities of the curve).
          GMSH_API void setTransfiniteCurve(const int tag,
                                            const int numNodes,
                                            const std::string & meshType = "Progression",
                                            const double coef = 1.);
    
          // Set a transfinite meshing constraint on the surface `tag'. `arrangement'
          // describes the arrangement of the triangles when the surface is not flagged
          // as recombined: currently supported values are "Left", "Right",
          // "AlternateLeft" and "AlternateRight". `cornerTags' can be used to specify
          // the (3 or 4) corners of the transfinite interpolation explicitly;
          // specifying the corners explicitly is mandatory if the surface has more
          // that 3 or 4 points on its boundary.
          GMSH_API void setTransfiniteSurface(const int tag,
                                              const std::string & arrangement = "Left",
                                              const std::vector<int> & cornerTags = std::vector<int>());
    
          // Set a transfinite meshing constraint on the surface `tag'. `cornerTags'
          // can be used to specify the (6 or 8) corners of the transfinite
          // interpolation explicitly.
          GMSH_API void setTransfiniteVolume(const int tag,
                                             const std::vector<int> & cornerTags = std::vector<int>());
    
          // Set a recombination meshing constraint on the model entity of dimension
          // `dim' and tag `tag'. Currently only entities of dimension 2 (to recombine
          // triangles into quadrangles) are supported.
          GMSH_API void setRecombine(const int dim,
                                     const int tag);
    
          // Set a smoothing meshing constraint on the model entity of dimension `dim'
          // and tag `tag'. `val' iterations of a Laplace smoother are applied.
          GMSH_API void setSmoothing(const int dim,
                                     const int tag,
                                     const int val);
    
          // Set a reverse meshing constraint on the model entity of dimension `dim'
          // and tag `tag'. If `val' is true, the mesh orientation will be reversed
          // with respect to the natural mesh orientation (i.e. the orientation
          // consistent with the orientation of the geometry). If `val' is false, the
          // mesh is left as-is.
          GMSH_API void setReverse(const int dim,
                                   const int tag,
                                   const bool val = true);
    
          // Set meshing constraints on the bounding surfaces of the volume of tag
          // `tag' so that all surfaces are oriented with outward pointing normals.
          // Currently only available with the OpenCASCADE kernel, as it relies on the
          // STL triangulation.
          GMSH_API void setOutwardOrientation(const int tag);
    
          // Embed the model entities of dimension `dim' and tags `tags' in the (inDim,
          // inTag) model entity. `inDim' must be strictly greater than `dim'.
          GMSH_API void embed(const int dim,
                              const std::vector<int> & tags,
                              const int inDim,
                              const int inTag);
    
          // Remove embedded entities in the model entities `dimTags'. if `dim' is >=
          // 0, only remove embedded entities of the given dimension (e.g. embedded
          // points if `dim' == 0).
          GMSH_API void removeEmbedded(const gmsh::vectorpair & dimTags,
                                       const int dim = -1);
    
          // Reorder the elements of type `elementType' classified on the entity of tag
          // `tag' according to `ordering'.
          GMSH_API void reorderElements(const int elementType,
                                        const int tag,
                                        const std::vector<std::size_t> & ordering);
    
          // Renumber the node tags in a continuous sequence.
          GMSH_API void renumberNodes();
    
          // Renumber the element tags in a continuous sequence.
          GMSH_API void renumberElements();
    
          // Set the meshes of the entities of dimension `dim' and tag `tags' as
          // periodic copies of the meshes of entities `tagsMaster', using the affine
          // transformation specified in `affineTransformation' (16 entries of a 4x4
          // matrix, by row). Currently only available for `dim' == 1 and `dim' == 2.
          GMSH_API void setPeriodic(const int dim,
                                    const std::vector<int> & tags,
                                    const std::vector<int> & tagsMaster,
                                    const std::vector<double> & affineTransform);
    
          // Get the master entity `tagMaster', the node tags `nodeTags' and their
          // corresponding master node tags `nodeTagsMaster', and the affine transform
          // `affineTransform' for the entity of dimension `dim' and tag `tag'.
          GMSH_API void getPeriodicNodes(const int dim,
                                         const int tag,
                                         int & tagMaster,
                                         std::vector<std::size_t> & nodeTags,
                                         std::vector<std::size_t> & nodeTagsMaster,
                                         std::vector<double> & affineTransform);
    
          // Remove duplicate nodes in the mesh of the current model.
          GMSH_API void removeDuplicateNodes();
    
          // Split (into two triangles) all quadrangles in surface `tag' whose quality
          // is lower than `quality'. If `tag' < 0, split quadrangles in all surfaces.
          GMSH_API void splitQuadrangles(const double quality = 1.,
                                         const int tag = -1);
    
          // Classify ("color") the surface mesh based on the angle threshold `angle'
          // (in radians), and create new discrete surfaces, curves and points
          // accordingly. If `boundary' is set, also create discrete curves on the
          // boundary if the surface is open. If `forReparametrization' is set, create
          // edges and surfaces that can be reparametrized using a single map.
          GMSH_API void classifySurfaces(const double angle,
                                         const bool boundary = true,
                                         const bool forReparametrization = false);
    
          // Create a parametrization for discrete curves and surfaces (i.e. curves and
          // surfaces represented solely by a mesh, without an underlying CAD
          // description), assuming that each can be parametrized with a single map.
          GMSH_API void createGeometry();
    
          // Create a boundary representation from the mesh if the model does not have
          // one (e.g. when imported from mesh file formats with no BRep representation
          // of the underlying model).
          GMSH_API void createTopology();
    
          // Compute a basis representation for homology spaces after a mesh has been
          // generated. The computation domain is given in a list of physical group
          // tags `domainTags'; if empty, the whole mesh is the domain. The computation
          // subdomain for relative homology computation is given in a list of physical
          // group tags `subdomainTags'; if empty, absolute homology is computed. The
          // dimensions homology bases to be computed are given in the list `dim'; if
          // empty, all bases are computed. Resulting basis representation chains are
          // stored as physical groups in the mesh.
          GMSH_API void computeHomology(const std::vector<int> & domainTags = std::vector<int>(),
                                        const std::vector<int> & subdomainTags = std::vector<int>(),
                                        const std::vector<int> & dims = std::vector<int>());
    
          // Compute a basis representation for cohomology spaces after a mesh has been
          // generated. The computation domain is given in a list of physical group
          // tags `domainTags'; if empty, the whole mesh is the domain. The computation
          // subdomain for relative cohomology computation is given in a list of
          // physical group tags `subdomainTags'; if empty, absolute cohomology is
          // computed. The dimensions homology bases to be computed are given in the
          // list `dim'; if empty, all bases are computed. Resulting basis
          // representation cochains are stored as physical groups in the mesh.
          GMSH_API void computeCohomology(const std::vector<int> & domainTags = std::vector<int>(),
                                          const std::vector<int> & subdomainTags = std::vector<int>(),
                                          const std::vector<int> & dims = std::vector<int>());
    
          namespace field { // Mesh size field functions
    
            // Add a new mesh size field of type `fieldType'. If `tag' is positive,
            // assign the tag explicitly; otherwise a new tag is assigned
            // automatically. Return the field tag.
            GMSH_API int add(const std::string & fieldType,
                             const int tag = -1);
    
            // Remove the field with tag `tag'.
            GMSH_API void remove(const int tag);
    
            // Set the numerical option `option' to value `value' for field `tag'.
            GMSH_API void setNumber(const int tag,
                                    const std::string & option,
                                    const double value);
    
            // Set the string option `option' to value `value' for field `tag'.
            GMSH_API void setString(const int tag,
                                    const std::string & option,
                                    const std::string & value);
    
            // Set the numerical list option `option' to value `value' for field `tag'.
            GMSH_API void setNumbers(const int tag,
                                     const std::string & option,
                                     const std::vector<double> & value);
    
            // Set the field `tag' as the background mesh size field.
            GMSH_API void setAsBackgroundMesh(const int tag);
    
            // Set the field `tag' as a boundary layer size field.
            GMSH_API void setAsBoundaryLayer(const int tag);
    
          } // namespace field
    
        } // namespace mesh
    
        namespace geo { // Built-in CAD kernel functions
    
          // Add a geometrical point in the built-in CAD representation, at coordinates
          // (`x', `y', `z'). If `meshSize' is > 0, add a meshing constraint at that
          // point. If `tag' is positive, set the tag explicitly; otherwise a new tag
          // is selected automatically. Return the tag of the point. (Note that the
          // point will be added in the current model only after `synchronize' is
          // called. This behavior holds for all the entities added in the geo module.)
          GMSH_API int addPoint(const double x,
                                const double y,
                                const double z,
                                const double meshSize = 0.,
                                const int tag = -1);
    
          // Add a straight line segment between the two points with tags `startTag'
          // and `endTag'. If `tag' is positive, set the tag explicitly; otherwise a
          // new tag is selected automatically. Return the tag of the line.
          GMSH_API int addLine(const int startTag,
                               const int endTag,
                               const int tag = -1);
    
          // Add a circle arc (strictly smaller than Pi) between the two points with
          // tags `startTag' and `endTag', with center `centertag'. If `tag' is
          // positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. If (`nx', `ny', `nz') != (0,0,0), explicitly set the plane
          // of the circle arc. Return the tag of the circle arc.
          GMSH_API int addCircleArc(const int startTag,
                                    const int centerTag,
                                    const int endTag,
                                    const int tag = -1,
                                    const double nx = 0.,
                                    const double ny = 0.,
                                    const double nz = 0.);
    
          // Add an ellipse arc (strictly smaller than Pi) between the two points
          // `startTag' and `endTag', with center `centertag' and major axis point
          // `majorTag'. If `tag' is positive, set the tag explicitly; otherwise a new
          // tag is selected automatically. If (`nx', `ny', `nz') != (0,0,0),
          // explicitly set the plane of the circle arc. Return the tag of the ellipse
          // arc.
          GMSH_API int addEllipseArc(const int startTag,
                                     const int centerTag,
                                     const int majorTag,
                                     const int endTag,
                                     const int tag = -1,
                                     const double nx = 0.,
                                     const double ny = 0.,
                                     const double nz = 0.);
    
          // Add a spline (Catmull-Rom) curve going through the points `pointTags'. If
          // `tag' is positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. Create a periodic curve if the first and last points are
          // the same. Return the tag of the spline curve.
          GMSH_API int addSpline(const std::vector<int> & pointTags,
                                 const int tag = -1);
    
          // Add a cubic b-spline curve with `pointTags' control points. If `tag' is
          // positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. Creates a periodic curve if the first and last points are
          // the same. Return the tag of the b-spline curve.
          GMSH_API int addBSpline(const std::vector<int> & pointTags,
                                  const int tag = -1);
    
          // Add a Bezier curve with `pointTags' control points. If `tag' is positive,
          // set the tag explicitly; otherwise a new tag is selected automatically.
          // Return the tag of the Bezier curve.
          GMSH_API int addBezier(const std::vector<int> & pointTags,
                                 const int tag = -1);
    
          // Add a curve loop (a closed wire) formed by the curves `curveTags'.
          // `curveTags' should contain (signed) tags of model enties of dimension 1
          // forming a closed loop: a negative tag signifies that the underlying curve
          // is considered with reversed orientation. If `tag' is positive, set the tag
          // explicitly; otherwise a new tag is selected automatically. Return the tag
          // of the curve loop.
          GMSH_API int addCurveLoop(const std::vector<int> & curveTags,
                                    const int tag = -1);
    
          // Add a plane surface defined by one or more curve loops `wireTags'. The
          // first curve loop defines the exterior contour; additional curve loop
          // define holes. If `tag' is positive, set the tag explicitly; otherwise a
          // new tag is selected automatically. Return the tag of the surface.
          GMSH_API int addPlaneSurface(const std::vector<int> & wireTags,
                                       const int tag = -1);
    
          // Add a surface filling the curve loops in `wireTags'. Currently only a
          // single curve loop is supported; this curve loop should be composed by 3 or
          // 4 curves only. If `tag' is positive, set the tag explicitly; otherwise a
          // new tag is selected automatically. Return the tag of the surface.
          GMSH_API int addSurfaceFilling(const std::vector<int> & wireTags,
                                         const int tag = -1,
                                         const int sphereCenterTag = -1);
    
          // Add a surface loop (a closed shell) formed by `surfaceTags'.  If `tag' is
          // positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. Return the tag of the shell.
          GMSH_API int addSurfaceLoop(const std::vector<int> & surfaceTags,
                                      const int tag = -1);
    
          // Add a volume (a region) defined by one or more shells `shellTags'. The
          // first surface loop defines the exterior boundary; additional surface loop
          // define holes. If `tag' is positive, set the tag explicitly; otherwise a
          // new tag is selected automatically. Return the tag of the volume.
          GMSH_API int addVolume(const std::vector<int> & shellTags,
                                 const int tag = -1);
    
          // Extrude the model entities `dimTags' by translation along (`dx', `dy',
          // `dz'). Return extruded entities in `outDimTags'. If `numElements' is not
          // empty, also extrude the mesh: the entries in `numElements' give the number
          // of elements in each layer. If `height' is not empty, it provides the
          // (cumulative) height of the different layers, normalized to 1. If `dx' ==
          // `dy' == `dz' == 0, the entities are extruded along their normal.
          GMSH_API void extrude(const gmsh::vectorpair & dimTags,
                                const double dx,
                                const double dy,
                                const double dz,
                                gmsh::vectorpair & outDimTags,
                                const std::vector<int> & numElements = std::vector<int>(),
                                const std::vector<double> & heights = std::vector<double>(),
                                const bool recombine = false);
    
          // Extrude the model entities `dimTags' by rotation of `angle' radians around
          // the axis of revolution defined by the point (`x', `y', `z') and the
          // direction (`ax', `ay', `az'). Return extruded entities in `outDimTags'. If
          // `numElements' is not empty, also extrude the mesh: the entries in
          // `numElements' give the number of elements in each layer. If `height' is
          // not empty, it provides the (cumulative) height of the different layers,
          // normalized to 1.
          GMSH_API void revolve(const gmsh::vectorpair & dimTags,
                                const double x,
                                const double y,
                                const double z,
                                const double ax,
                                const double ay,
                                const double az,
                                const double angle,
                                gmsh::vectorpair & outDimTags,
                                const std::vector<int> & numElements = std::vector<int>(),
                                const std::vector<double> & heights = std::vector<double>(),
                                const bool recombine = false);
    
          // Extrude the model entities `dimTags' by a combined translation and
          // rotation of `angle' radians, along (`dx', `dy', `dz') and around the axis
          // of revolution defined by the point (`x', `y', `z') and the direction
          // (`ax', `ay', `az'). Return extruded entities in `outDimTags'. If
          // `numElements' is not empty, also extrude the mesh: the entries in
          // `numElements' give the number of elements in each layer. If `height' is
          // not empty, it provides the (cumulative) height of the different layers,
          // normalized to 1.
          GMSH_API void twist(const gmsh::vectorpair & dimTags,
                              const double x,
                              const double y,
                              const double z,
                              const double dx,
                              const double dy,
                              const double dz,
                              const double ax,
                              const double ay,
                              const double az,
                              const double angle,
                              gmsh::vectorpair & outDimTags,
                              const std::vector<int> & numElements = std::vector<int>(),
                              const std::vector<double> & heights = std::vector<double>(),
                              const bool recombine = false);
    
          // Translate the model entities `dimTags' along (`dx', `dy', `dz').
          GMSH_API void translate(const gmsh::vectorpair & dimTags,
                                  const double dx,
                                  const double dy,
                                  const double dz);
    
          // Rotate the model entities `dimTags' of `angle' radians around the axis of
          // revolution defined by the point (`x', `y', `z') and the direction (`ax',
          // `ay', `az').
          GMSH_API void rotate(const gmsh::vectorpair & dimTags,
                               const double x,
                               const double y,
                               const double z,
                               const double ax,
                               const double ay,
                               const double az,
                               const double angle);
    
          // Scale the model entities `dimTag' by factors `a', `b' and `c' along the
          // three coordinate axes; use (`x', `y', `z') as the center of the homothetic
          // transformation.
          GMSH_API void dilate(const gmsh::vectorpair & dimTags,
                               const double x,
                               const double y,
                               const double z,
                               const double a,
                               const double b,
                               const double c);
    
          // Apply a symmetry transformation to the model entities `dimTag', with
          // respect to the plane of equation `a' * x + `b' * y + `c' * z + `d' = 0.
          GMSH_API void symmetrize(const gmsh::vectorpair & dimTags,
                                   const double a,
                                   const double b,
                                   const double c,
                                   const double d);
    
          // Copy the entities `dimTags'; the new entities are returned in
          // `outDimTags'.
          GMSH_API void copy(const gmsh::vectorpair & dimTags,
                             gmsh::vectorpair & outDimTags);
    
          // Remove the entities `dimTags'. If `recursive' is true, remove all the
          // entities on their boundaries, down to dimension 0.
          GMSH_API void remove(const gmsh::vectorpair & dimTags,
                               const bool recursive = false);
    
          // Remove all duplicate entities (different entities at the same geometrical
          // location).
          GMSH_API void removeAllDuplicates();
    
          // Synchronize the built-in CAD representation with the current Gmsh model.
          // This can be called at any time, but since it involves a non trivial amount
          // of processing, the number of synchronization points should normally be
          // minimized.
          GMSH_API void synchronize();
    
          namespace mesh { // Built-in CAD kernel meshing constraints
    
            // Set a mesh size constraint on the model entities `dimTags'. Currently
            // only entities of dimension 0 (points) are handled.
            GMSH_API void setSize(const gmsh::vectorpair & dimTags,
                                  const double size);
    
            // Set a transfinite meshing constraint on the curve `tag', with `numNodes'
            // nodes distributed according to `meshType' and `coef'. Currently
            // supported types are "Progression" (geometrical progression with power
            // `coef') and "Bump" (refinement toward both extremities of the curve).
            GMSH_API void setTransfiniteCurve(const int tag,
                                              const int nPoints,
                                              const std::string & meshType = "Progression",
                                              const double coef = 1.);
    
            // Set a transfinite meshing constraint on the surface `tag'. `arrangement'
            // describes the arrangement of the triangles when the surface is not
            // flagged as recombined: currently supported values are "Left", "Right",
            // "AlternateLeft" and "AlternateRight". `cornerTags' can be used to
            // specify the (3 or 4) corners of the transfinite interpolation
            // explicitly; specifying the corners explicitly is mandatory if the
            // surface has more that 3 or 4 points on its boundary.
            GMSH_API void setTransfiniteSurface(const int tag,
                                                const std::string & arrangement = "Left",
                                                const std::vector<int> & cornerTags = std::vector<int>());
    
            // Set a transfinite meshing constraint on the surface `tag'. `cornerTags'
            // can be used to specify the (6 or 8) corners of the transfinite
            // interpolation explicitly.
            GMSH_API void setTransfiniteVolume(const int tag,
                                               const std::vector<int> & cornerTags = std::vector<int>());
    
            // Set a recombination meshing constraint on the model entity of dimension
            // `dim' and tag `tag'. Currently only entities of dimension 2 (to
            // recombine triangles into quadrangles) are supported.
            GMSH_API void setRecombine(const int dim,
                                       const int tag,
                                       const double angle = 45.);
    
            // Set a smoothing meshing constraint on the model entity of dimension
            // `dim' and tag `tag'. `val' iterations of a Laplace smoother are applied.
            GMSH_API void setSmoothing(const int dim,
                                       const int tag,
                                       const int val);
    
            // Set a reverse meshing constraint on the model entity of dimension `dim'
            // and tag `tag'. If `val' is true, the mesh orientation will be reversed
            // with respect to the natural mesh orientation (i.e. the orientation
            // consistent with the orientation of the geometry). If `val' is false, the
            // mesh is left as-is.
            GMSH_API void setReverse(const int dim,
                                     const int tag,
                                     const bool val = true);
    
          } // namespace mesh
    
        } // namespace geo
    
        namespace occ { // OpenCASCADE CAD kernel functions
    
          // Add a geometrical point in the OpenCASCADE CAD representation, at
          // coordinates (`x', `y', `z'). If `meshSize' is > 0, add a meshing
          // constraint at that point. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically. Return the tag of the
          // point. (Note that the point will be added in the current model only after
          // `synchronize' is called. This behavior holds for all the entities added in
          // the occ module.)
          GMSH_API int addPoint(const double x,
                                const double y,
                                const double z,
                                const double meshSize = 0.,
                                const int tag = -1);
    
          // Add a straight line segment between the two points with tags `startTag'
          // and `endTag'. If `tag' is positive, set the tag explicitly; otherwise a
          // new tag is selected automatically. Return the tag of the line.
          GMSH_API int addLine(const int startTag,
                               const int endTag,
                               const int tag = -1);
    
          // Add a circle arc between the two points with tags `startTag' and `endTag',
          // with center `centerTag'. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically. Return the tag of the
          // circle arc.
          GMSH_API int addCircleArc(const int startTag,
                                    const int centerTag,
                                    const int endTag,
                                    const int tag = -1);
    
          // Add a circle of center (`x', `y', `z') and radius `r'. If `tag' is
          // positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. If `angle1' and `angle2' are specified, create a circle arc
          // between the two angles. Return the tag of the circle.
          GMSH_API int addCircle(const double x,
                                 const double y,
                                 const double z,
                                 const double r,
                                 const int tag = -1,
                                 const double angle1 = 0.,
                                 const double angle2 = 2*M_PI);
    
          // Add an ellipse arc between the major axis point `startTag' and `endTag',
          // with center `centerTag'. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically. Return the tag of the
          // ellipse arc. Note that OpenCASCADE does not allow creating ellipse arcs
          // with the major radius smaller than the minor radius.
          GMSH_API int addEllipseArc(const int startTag,
                                     const int centerTag,
                                     const int endTag,
                                     const int tag = -1);
    
          // Add an ellipse of center (`x', `y', `z') and radii `r1' and `r2' along the
          // x- and y-axes respectively. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically. If `angle1' and `angle2'
          // are specified, create an ellipse arc between the two angles. Return the
          // tag of the ellipse. Note that OpenCASCADE does not allow creating ellipses
          // with the major radius (along the x-axis) smaller than or equal to the
          // minor radius (along the y-axis): rotate the shape or use `addCircle' in
          // such cases.
          GMSH_API int addEllipse(const double x,
                                  const double y,
                                  const double z,
                                  const double r1,
                                  const double r2,
                                  const int tag = -1,
                                  const double angle1 = 0.,
                                  const double angle2 = 2*M_PI);
    
          // Add a spline (C2 b-spline) curve going through the points `pointTags'. If
          // `tag' is positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. Create a periodic curve if the first and last points are
          // the same. Return the tag of the spline curve.
          GMSH_API int addSpline(const std::vector<int> & pointTags,
                                 const int tag = -1);
    
          // Add a b-spline curve of degree `degree' with `pointTags' control points.
          // If `weights', `knots' or `multiplicities' are not provided, default
          // parameters are computed automatically. If `tag' is positive, set the tag
          // explicitly; otherwise a new tag is selected automatically. Create a
          // periodic curve if the first and last points are the same. Return the tag
          // of the b-spline curve.
          GMSH_API int addBSpline(const std::vector<int> & pointTags,
                                  const int tag = -1,
                                  const int degree = 3,
                                  const std::vector<double> & weights = std::vector<double>(),
                                  const std::vector<double> & knots = std::vector<double>(),
                                  const std::vector<int> & multiplicities = std::vector<int>());
    
          // Add a Bezier curve with `pointTags' control points. If `tag' is positive,
          // set the tag explicitly; otherwise a new tag is selected automatically.
          // Return the tag of the Bezier curve.
          GMSH_API int addBezier(const std::vector<int> & pointTags,
                                 const int tag = -1);
    
          // Add a wire (open or closed) formed by the curves `curveTags'. Note that an
          // OpenCASCADE wire can be made of curves that share geometrically identical
          // (but topologically different) points. If `tag' is positive, set the tag
          // explicitly; otherwise a new tag is selected automatically. Return the tag
          // of the wire.
          GMSH_API int addWire(const std::vector<int> & curveTags,
                               const int tag = -1,
                               const bool checkClosed = false);
    
          // Add a curve loop (a closed wire) formed by the curves `curveTags'.
          // `curveTags' should contain tags of curves forming a closed loop. Note that
          // an OpenCASCADE curve loop can be made of curves that share geometrically
          // identical (but topologically different) points. If `tag' is positive, set
          // the tag explicitly; otherwise a new tag is selected automatically. Return
          // the tag of the curve loop.
          GMSH_API int addCurveLoop(const std::vector<int> & curveTags,
                                    const int tag = -1);
    
          // Add a rectangle with lower left corner at (`x', `y', `z') and upper right
          // corner at (`x' + `dx', `y' + `dy', `z'). If `tag' is positive, set the tag
          // explicitly; otherwise a new tag is selected automatically. Round the
          // corners if `roundedRadius' is nonzero. Return the tag of the rectangle.
          GMSH_API int addRectangle(const double x,
                                    const double y,
                                    const double z,
                                    const double dx,
                                    const double dy,
                                    const int tag = -1,
                                    const double roundedRadius = 0.);
    
          // Add a disk with center (`xc', `yc', `zc') and radius `rx' along the x-axis
          // and `ry' along the y-axis. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically. Return the tag of the disk.
          GMSH_API int addDisk(const double xc,
                               const double yc,
                               const double zc,
                               const double rx,
                               const double ry,
                               const int tag = -1);
    
          // Add a plane surface defined by one or more curve loops (or closed wires)
          // `wireTags'. The first curve loop defines the exterior contour; additional
          // curve loop define holes. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically. Return the tag of the
          // surface.
          GMSH_API int addPlaneSurface(const std::vector<int> & wireTags,
                                       const int tag = -1);
    
          // Add a surface filling the curve loops in `wireTags'. If `tag' is positive,
          // set the tag explicitly; otherwise a new tag is selected automatically.
          // Return the tag of the surface. If `pointTags' are provided, force the
          // surface to pass through the given points.
          GMSH_API int addSurfaceFilling(const int wireTag,
                                         const int tag = -1,
                                         const std::vector<int> & pointTags = std::vector<int>());
    
          // Add a surface loop (a closed shell) formed by `surfaceTags'.  If `tag' is
          // positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. Return the tag of the surface loop. Setting `sewing' allows
          // to build a shell made of surfaces that share geometrically identical (but
          // topologically different) curves.
          GMSH_API int addSurfaceLoop(const std::vector<int> & surfaceTags,
                                      const int tag = -1,
                                      const bool sewing = false);
    
          // Add a volume (a region) defined by one or more surface loops `shellTags'.
          // The first surface loop defines the exterior boundary; additional surface
          // loop define holes. If `tag' is positive, set the tag explicitly; otherwise
          // a new tag is selected automatically. Return the tag of the volume.
          GMSH_API int addVolume(const std::vector<int> & shellTags,
                                 const int tag = -1);
    
          // Add a sphere of center (`xc', `yc', `zc') and radius `r'. The optional
          // `angle1' and `angle2' arguments define the polar angle opening (from -Pi/2
          // to Pi/2). The optional `angle3' argument defines the azimuthal opening
          // (from 0 to 2*Pi). If `tag' is positive, set the tag explicitly; otherwise
          // a new tag is selected automatically. Return the tag of the sphere.
          GMSH_API int addSphere(const double xc,
                                 const double yc,
                                 const double zc,
                                 const double radius,
                                 const int tag = -1,
                                 const double angle1 = -M_PI/2,
                                 const double angle2 = M_PI/2,
                                 const double angle3 = 2*M_PI);
    
          // Add a parallelepipedic box defined by a point (`x', `y', `z') and the
          // extents along the x-, y- and z-axes. If `tag' is positive, set the tag
          // explicitly; otherwise a new tag is selected automatically. Return the tag
          // of the box.
          GMSH_API int addBox(const double x,
                              const double y,
                              const double z,
                              const double dx,
                              const double dy,
                              const double dz,
                              const int tag = -1);
    
          // Add a cylinder, defined by the center (`x', `y', `z') of its first
          // circular face, the 3 components (`dx', `dy', `dz') of the vector defining
          // its axis and its radius `r'. The optional `angle' argument defines the
          // angular opening (from 0 to 2*Pi). If `tag' is positive, set the tag
          // explicitly; otherwise a new tag is selected automatically. Return the tag
          // of the cylinder.
          GMSH_API int addCylinder(const double x,
                                   const double y,
                                   const double z,
                                   const double dx,
                                   const double dy,
                                   const double dz,
                                   const double r,
                                   const int tag = -1,
                                   const double angle = 2*M_PI);
    
          // Add a cone, defined by the center (`x', `y', `z') of its first circular
          // face, the 3 components of the vector (`dx', `dy', `dz') defining its axis
          // and the two radii `r1' and `r2' of the faces (these radii can be zero). If
          // `tag' is positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. `angle' defines the optional angular opening (from 0 to
          // 2*Pi). Return the tag of the cone.
          GMSH_API int addCone(const double x,
                               const double y,
                               const double z,
                               const double dx,
                               const double dy,
                               const double dz,
                               const double r1,
                               const double r2,
                               const int tag = -1,
                               const double angle = 2*M_PI);
    
          // Add a right angular wedge, defined by the right-angle point (`x', `y',
          // `z') and the 3 extends along the x-, y- and z-axes (`dx', `dy', `dz'). If
          // `tag' is positive, set the tag explicitly; otherwise a new tag is selected
          // automatically. The optional argument `ltx' defines the top extent along
          // the x-axis. Return the tag of the wedge.
          GMSH_API int addWedge(const double x,
                                const double y,
                                const double z,
                                const double dx,
                                const double dy,
                                const double dz,
                                const int tag = -1,
                                const double ltx = 0.);
    
          // Add a torus, defined by its center (`x', `y', `z') and its 2 radii `r' and
          // `r2'. If `tag' is positive, set the tag explicitly; otherwise a new tag is
          // selected automatically. The optional argument `angle' defines the angular
          // opening (from 0 to 2*Pi). Return the tag of the wedge.
          GMSH_API int addTorus(const double x,
                                const double y,
                                const double z,
                                const double r1,
                                const double r2,
                                const int tag = -1,
                                const double angle = 2*M_PI);
    
          // Add a volume (if the optional argument `makeSolid' is set) or surfaces
          // defined through the open or closed wires `wireTags'. If `tag' is positive,
          // set the tag explicitly; otherwise a new tag is selected automatically. The
          // new entities are returned in `outDimTags'. If the optional argument
          // `makeRuled' is set, the surfaces created on the boundary are forced to be
          // ruled surfaces.
          GMSH_API void addThruSections(const std::vector<int> & wireTags,
                                        gmsh::vectorpair & outDimTags,
                                        const int tag = -1,
                                        const bool makeSolid = true,
                                        const bool makeRuled = false);
    
          // Add a hollowed volume built from an initial volume `volumeTag' and a set
          // of faces from this volume `excludeSurfaceTags', which are to be removed.
          // The remaining faces of the volume become the walls of the hollowed solid,
          // with thickness `offset'. If `tag' is positive, set the tag explicitly;
          // otherwise a new tag is selected automatically.
          GMSH_API void addThickSolid(const int volumeTag,
                                      const std::vector<int> & excludeSurfaceTags,
                                      const double offset,
                                      gmsh::vectorpair & outDimTags,
                                      const int tag = -1);
    
          // Extrude the model entities `dimTags' by translation along (`dx', `dy',
          // `dz'). Return extruded entities in `outDimTags'. If `numElements' is not
          // empty, also extrude the mesh: the entries in `numElements' give the number
          // of elements in each layer. If `height' is not empty, it provides the
          // (cumulative) height of the different layers, normalized to 1.
          GMSH_API void extrude(const gmsh::vectorpair & dimTags,
                                const double dx,
                                const double dy,
                                const double dz,
                                gmsh::vectorpair & outDimTags,
                                const std::vector<int> & numElements = std::vector<int>(),
                                const std::vector<double> & heights = std::vector<double>(),
                                const bool recombine = false);
    
          // Extrude the model entities `dimTags' by rotation of `angle' radians around
          // the axis of revolution defined by the point (`x', `y', `z') and the
          // direction (`ax', `ay', `az'). Return extruded entities in `outDimTags'. If
          // `numElements' is not empty, also extrude the mesh: the entries in
          // `numElements' give the number of elements in each layer. If `height' is
          // not empty, it provides the (cumulative) height of the different layers,
          // normalized to 1.
          GMSH_API void revolve(const gmsh::vectorpair & dimTags,
                                const double x,
                                const double y,
                                const double z,
                                const double ax,
                                const double ay,
                                const double az,
                                const double angle,
                                gmsh::vectorpair & outDimTags,
                                const std::vector<int> & numElements = std::vector<int>(),
                                const std::vector<double> & heights = std::vector<double>(),
                                const bool recombine = false);
    
          // Add a pipe by extruding the entities `dimTags' along the wire `wireTag'.
          // Return the pipe in `outDimTags'.
          GMSH_API void addPipe(const gmsh::vectorpair & dimTags,
                                const int wireTag,
                                gmsh::vectorpair & outDimTags);
    
          // Fillet the volumes `volumeTags' on the curves `curveTags' with radii
          // `radii'. The `radii' vector can either contain a single radius, as many
          // radii as `curveTags', or twice as many as `curveTags' (in which case
          // different radii are provided for the begin and end points of the curves).
          // Return the filleted entities in `outDimTags'. Remove the original volume
          // if `removeVolume' is set.
          GMSH_API void fillet(const std::vector<int> & volumeTags,
                               const std::vector<int> & curveTags,
                               const std::vector<double> & radii,
                               gmsh::vectorpair & outDimTags,
                               const bool removeVolume = true);
    
          // Chamfer the volumes `volumeTags' on the curves `curveTags' with distances
          // `distances' measured on surfaces `surfaceTags'. The `distances' vector can
          // either contain a single distance, as many distances as `curveTags' and
          // `surfaceTags', or twice as many as `curveTags' and `surfaceTags' (in which
          // case the first in each pair is measured on the corresponding surface in
          // `surfaceTags', the other on the other adjacent surface). Return the
          // chamfered entities in `outDimTags'. Remove the original volume if
          // `removeVolume' is set.
          GMSH_API void chamfer(const std::vector<int> & volumeTags,
                                const std::vector<int> & curveTags,
                                const std::vector<int> & surfaceTags,
                                const std::vector<double> & distances,
                                gmsh::vectorpair & outDimTags,
                                const bool removeVolume = true);
    
          // Compute the boolean union (the fusion) of the entities `objectDimTags' and
          // `toolDimTags'. Return the resulting entities in `outDimTags'. If `tag' is
          // positive, try to set the tag explicitly (only valid if the boolean
          // operation results in a single entity). Remove the object if `removeObject'
          // is set. Remove the tool if `removeTool' is set.
          GMSH_API void fuse(const gmsh::vectorpair & objectDimTags,
                             const gmsh::vectorpair & toolDimTags,
                             gmsh::vectorpair & outDimTags,
                             std::vector<gmsh::vectorpair> & outDimTagsMap,
                             const int tag = -1,
                             const bool removeObject = true,
                             const bool removeTool = true);
    
          // Compute the boolean intersection (the common parts) of the entities
          // `objectDimTags' and `toolDimTags'. Return the resulting entities in
          // `outDimTags'. If `tag' is positive, try to set the tag explicitly (only
          // valid if the boolean operation results in a single entity). Remove the
          // object if `removeObject' is set. Remove the tool if `removeTool' is set.
          GMSH_API void intersect(const gmsh::vectorpair & objectDimTags,
                                  const gmsh::vectorpair & toolDimTags,
                                  gmsh::vectorpair & outDimTags,
                                  std::vector<gmsh::vectorpair> & outDimTagsMap,
                                  const int tag = -1,
                                  const bool removeObject = true,
                                  const bool removeTool = true);
    
          // Compute the boolean difference between the entities `objectDimTags' and
          // `toolDimTags'. Return the resulting entities in `outDimTags'. If `tag' is
          // positive, try to set the tag explicitly (only valid if the boolean
          // operation results in a single entity). Remove the object if `removeObject'
          // is set. Remove the tool if `removeTool' is set.
          GMSH_API void cut(const gmsh::vectorpair & objectDimTags,
                            const gmsh::vectorpair & toolDimTags,
                            gmsh::vectorpair & outDimTags,
                            std::vector<gmsh::vectorpair> & outDimTagsMap,
                            const int tag = -1,
                            const bool removeObject = true,
                            const bool removeTool = true);
    
          // Compute the boolean fragments (general fuse) of the entities
          // `objectDimTags' and `toolDimTags'. Return the resulting entities in
          // `outDimTags'. If `tag' is positive, try to set the tag explicitly (only
          // valid if the boolean operation results in a single entity). Remove the
          // object if `removeObject' is set. Remove the tool if `removeTool' is set.
          GMSH_API void fragment(const gmsh::vectorpair & objectDimTags,
                                 const gmsh::vectorpair & toolDimTags,
                                 gmsh::vectorpair & outDimTags,
                                 std::vector<gmsh::vectorpair> & outDimTagsMap,
                                 const int tag = -1,
                                 const bool removeObject = true,
                                 const bool removeTool = true);
    
          // Translate the model entities `dimTags' along (`dx', `dy', `dz').
          GMSH_API void translate(const gmsh::vectorpair & dimTags,
                                  const double dx,
                                  const double dy,
                                  const double dz);
    
          // Rotate the model entities `dimTags' of `angle' radians around the axis of
          // revolution defined by the point (`x', `y', `z') and the direction (`ax',
          // `ay', `az').
          GMSH_API void rotate(const gmsh::vectorpair & dimTags,
                               const double x,
                               const double y,
                               const double z,
                               const double ax,
                               const double ay,
                               const double az,
                               const double angle);
    
          // Scale the model entities `dimTag' by factors `a', `b' and `c' along the
          // three coordinate axes; use (`x', `y', `z') as the center of the homothetic
          // transformation.
          GMSH_API void dilate(const gmsh::vectorpair & dimTags,
                               const double x,
                               const double y,
                               const double z,
                               const double a,
                               const double b,
                               const double c);
    
          // Apply a symmetry transformation to the model entities `dimTag', with
          // respect to the plane of equation `a' * x + `b' * y + `c' * z + `d' = 0.
          GMSH_API void symmetrize(const gmsh::vectorpair & dimTags,
                                   const double a,
                                   const double b,
                                   const double c,
                                   const double d);
    
          // Apply a general affine transformation matrix `a' (16 entries of a 4x4
          // matrix, by row; only the 12 first can be provided for convenience) to the
          // model entities `dimTag'.
          GMSH_API void affineTransform(const gmsh::vectorpair & dimTags,
                                        const std::vector<double> & a);
    
          // Copy the entities `dimTags'; the new entities are returned in
          // `outDimTags'.
          GMSH_API void copy(const gmsh::vectorpair & dimTags,
                             gmsh::vectorpair & outDimTags);
    
          // Remove the entities `dimTags'. If `recursive' is true, remove all the
          // entities on their boundaries, down to dimension 0.
          GMSH_API void remove(const gmsh::vectorpair & dimTags,
                               const bool recursive = false);
    
          // Remove all duplicate entities (different entities at the same geometrical
          // location) after intersecting (using boolean fragments) all highest
          // dimensional entities.
          GMSH_API void removeAllDuplicates();
    
          // Apply various healing procedures to the entities `dimTags' (or to all the
          // entities in the model if `dimTags' is empty). Return the healed entities
          // in `outDimTags'. Available healing options are listed in the Gmsh
          // reference manual.
          GMSH_API void healShapes(gmsh::vectorpair & outDimTags,
                                   const gmsh::vectorpair & dimTags = gmsh::vectorpair(),
                                   const double tolerance = 1e-8,
                                   const bool fixDegenerated = true,
                                   const bool fixSmallEdges = true,
                                   const bool fixSmallFaces = true,
                                   const bool sewFaces = true);
    
          // Import BREP, STEP or IGES shapes from the file `fileName'. The imported
          // entities are returned in `outDimTags'. If the optional argument
          // `highestDimOnly' is set, only import the highest dimensional entities in
          // the file. The optional argument `format' can be used to force the format
          // of the file (currently "brep", "step" or "iges").
          GMSH_API void importShapes(const std::string & fileName,
                                     gmsh::vectorpair & outDimTags,
                                     const bool highestDimOnly = true,
                                     const std::string & format = "");
    
          // Imports an OpenCASCADE `shape' by providing a pointer to a native
          // OpenCASCADE `TopoDS_Shape' object (passed as a pointer to void). The
          // imported entities are returned in `outDimTags'. If the optional argument
          // `highestDimOnly' is set, only import the highest dimensional entities in
          // `shape'. For C and C++ only. Warning: this function is unsafe, as
          // providing an invalid pointer will lead to undefined behavior.
          GMSH_API void importShapesNativePointer(const void * shape,
                                                  gmsh::vectorpair & outDimTags,
                                                  const bool highestDimOnly = true);
    
          // Set a mesh size constraint on the model entities `dimTags'. Currently only
          // entities of dimension 0 (points) are handled.
          GMSH_API void setMeshSize(const gmsh::vectorpair & dimTags,
                                    const double size);
    
          // Get the mass of the model entity of dimension `dim' and tag `tag'.
          GMSH_API void getMass(const int dim,
                                const int tag,
                                double & mass);
    
          // Get the center of mass of the model entity of dimension `dim' and tag
          // `tag'.
          GMSH_API void getCenterOfMass(const int dim,
                                        const int tag,
                                        double & x,
                                        double & y,
                                        double & z);
    
          // Get the matrix of inertia (by row) of the model entity of dimension `dim'
          // and tag `tag'.
          GMSH_API void getMatrixOfInertia(const int dim,
                                           const int tag,
                                           std::vector<double> & mat);
    
          // Synchronize the OpenCASCADE CAD representation with the current Gmsh
          // model. This can be called at any time, but since it involves a non trivial
          // amount of processing, the number of synchronization points should normally
          // be minimized.
          GMSH_API void synchronize();
    
        } // namespace occ
    
      } // namespace model
    
      namespace view { // Post-processing view functions
    
        // Add a new post-processing view, with name `name'. If `tag' is positive use
        // it (and remove the view with that tag if it already exists), otherwise
        // associate a new tag. Return the view tag.
        GMSH_API int add(const std::string & name,
                         const int tag = -1);
    
        // Remove the view with tag `tag'.
        GMSH_API void remove(const int tag);
    
        // Get the index of the view with tag `tag' in the list of currently loaded
        // views. This dynamic index (it can change when views are removed) is used to
        // access view options.
        GMSH_API int getIndex(const int tag);
    
        // Get the tags of all views.
        GMSH_API void getTags(std::vector<int> & tags);
    
        // Add model-based post-processing data to the view with tag `tag'. `modelName'
        // identifies the model the data is attached to. `dataType' specifies the type
        // of data, currently either "NodeData", "ElementData" or "ElementNodeData".
        // `step' specifies the identifier (>= 0) of the data in a sequence. `tags'
        // gives the tags of the nodes or elements in the mesh to which the data is
        // associated. `data' is a vector of the same length as `tags': each entry is
        // the vector of double precision numbers representing the data associated with
        // the corresponding tag. The optional `time' argument associate a time value
        // with the data. `numComponents' gives the number of data components (1 for
        // scalar data, 3 for vector data, etc.) per entity; if negative, it is
        // automatically inferred (when possible) from the input data. `partition'
        // allows to specify data in several sub-sets.
        GMSH_API void addModelData(const int tag,
                                   const int step,
                                   const std::string & modelName,
                                   const std::string & dataType,
                                   const std::vector<std::size_t> & tags,
                                   const std::vector<std::vector<double> > & data,
                                   const double time = 0.,
                                   const int numComponents = -1,
                                   const int partition = 0);
    
        // Get model-based post-processing data from the view with tag `tag' at step
        // `step'. Return the `data' associated to the nodes or the elements with tags
        // `tags', as well as the `dataType' and the number of components
        // `numComponents'.
        GMSH_API void getModelData(const int tag,
                                   const int step,
                                   std::string & dataType,
                                   std::vector<std::size_t> & tags,
                                   std::vector<std::vector<double> > & data,
                                   double & time,
                                   int & numComponents);
    
        // Add list-based post-processing data to the view with tag `tag'. `dataType'
        // identifies the data: "SP" for scalar points, "VP", for vector points, etc.
        // `numEle' gives the number of elements in the data. `data' contains the data
        // for the `numEle' elements.
        GMSH_API void addListData(const int tag,
                                  const std::string & dataType,
                                  const int numEle,
                                  const std::vector<double> & data);
    
        // Get list-based post-processing data from the view with tag `tag'. Return the
        // types `dataTypes', the number of elements `numElements' for each data type
        // and the `data' for each data type.
        GMSH_API void getListData(const int tag,
                                  std::vector<std::string> & dataType,
                                  std::vector<int> & numElements,
                                  std::vector<std::vector<double> > & data);
    
        // Add a post-processing view as an `alias' of the reference view with tag
        // `refTag'. If `copyOptions' is set, copy the options of the reference view.
        // If `tag' is positive use it (and remove the view with that tag if it already
        // exists), otherwise associate a new tag. Return the view tag.
        GMSH_API int addAlias(const int refTag,
                              const bool copyOptions = false,
                              const int tag = -1);
    
        // Copy the options from the view with tag `refTag' to the view with tag `tag'.
        GMSH_API void copyOptions(const int refTag,
                                  const int tag);
    
        // Combine elements (if `what' == "elements") or steps (if `what' == "steps")
        // of all views (`how' == "all"), all visible views (`how' == "visible") or all
        // views having the same name (`how' == "name"). Remove original views if
        // `remove' is set.
        GMSH_API void combine(const std::string & what,
                              const std::string & how,
                              const bool remove = false);
    
        // Probe the view `tag' for its `value' at point (`x', `y', `z'). Return only
        // the value at step `step' is `step' is positive. Return only values with
        // `numComp' if `numComp' is positive. Return the gradient of the `value' if
        // `gradient' is set. Probes with a geometrical tolerance (in the reference
        // unit cube) of `tolerance' if `tolerance' is not zero. Return the result from
        // the element described by its coordinates if `xElementCoord', `yElementCoord'
        // and `zElementCoord' are provided.
        GMSH_API void probe(const int tag,
                            const double x,
                            const double y,
                            const double z,
                            std::vector<double> & value,
                            const int step = -1,
                            const int numComp = -1,
                            const bool gradient = false,
                            const double tolerance = 0.,
                            const std::vector<double> & xElemCoord = std::vector<double>(),
                            const std::vector<double> & yElemCoord = std::vector<double>(),
                            const std::vector<double> & zElemCoord = std::vector<double>());
    
        // Write the view to a file `fileName'. The export format is determined by the
        // file extension. Append to the file if `append' is set.
        GMSH_API void write(const int tag,
                            const std::string & fileName,
                            const bool append = false);
    
      } // namespace view
    
      namespace plugin { // Plugin functions
    
        // Set the numerical option `option' to the value `value' for plugin `name'.
        GMSH_API void setNumber(const std::string & name,
                                const std::string & option,
                                const double value);
    
        // Set the string option `option' to the value `value' for plugin `name'.
        GMSH_API void setString(const std::string & name,
                                const std::string & option,
                                const std::string & value);
    
        // Run the plugin `name'.
        GMSH_API void run(const std::string & name);
    
      } // namespace plugin
    
      namespace graphics { // Graphics functions
    
        // Draw all the OpenGL scenes.
        GMSH_API void draw();
    
      } // namespace graphics
    
      namespace fltk { // FLTK graphical user interface functions
    
        // Create the FLTK graphical user interface. Can only be called in the main
        // thread.
        GMSH_API void initialize();
    
        // Wait at most `time' seconds for user interface events and return. If `time'
        // < 0, wait indefinitely. First automatically create the user interface if it
        // has not yet been initialized. Can only be called in the main thread.
        GMSH_API void wait(const double time = -1.);
    
        // Update the user interface (potentially creating new widgets and windows).
        // First automatically create the user interface if it has not yet been
        // initialized. Can only be called in the main thread: use `awake("update")' to
        // trigger an update of the user interface from another thread.
        GMSH_API void update();
    
        // Awake the main user interface thread and process pending events, and
        // optionally perform an action (currently the only `action' allowed is
        // "update").
        GMSH_API void awake(const std::string & action = "");
    
        // Block the current thread until it can safely modify the user interface.
        GMSH_API void lock();
    
        // Release the lock that was set using lock.
        GMSH_API void unlock();
    
        // Run the event loop of the graphical user interface, i.e. repeatedly calls
        // `wait()'. First automatically create the user interface if it has not yet
        // been initialized. Can only be called in the main thread.
        GMSH_API void run();
    
        // Select entities in the user interface. If `dim' is >= 0, return only the
        // entities of the specified dimension (e.g. points if `dim' == 0).
        GMSH_API int selectEntities(gmsh::vectorpair & dimTags,
                                    const int dim = -1);
    
        // Select elements in the user interface.
        GMSH_API int selectElements(std::vector<std::size_t> & elementTags);
    
        // Select views in the user interface.
        GMSH_API int selectViews(std::vector<int> & viewTags);
    
      } // namespace fltk
    
      namespace onelab { // ONELAB server functions
    
        // Set one or more parameters in the ONELAB database, encoded in `format'.
        GMSH_API void set(const std::string & data,
                          const std::string & format = "json");
    
        // Get all the parameters (or a single one if `name' is specified) from the
        // ONELAB database, encoded in `format'.
        GMSH_API void get(std::string & data,
                          const std::string & name = "",
                          const std::string & format = "json");
    
        // Set the value of the number parameter `name' in the ONELAB database. Create
        // the parameter if it does not exist; update the value if the parameter
        // exists.
        GMSH_API void setNumber(const std::string & name,
                                const std::vector<double> & value);
    
        // Set the value of the string parameter `name' in the ONELAB database. Create
        // the parameter if it does not exist; update the value if the parameter
        // exists.
        GMSH_API void setString(const std::string & name,
                                const std::vector<std::string> & value);
    
        // Get the value of the number parameter `name' from the ONELAB database.
        // Return an empty vector if the parameter does not exist.
        GMSH_API void getNumber(const std::string & name,
                                std::vector<double> & value);
    
        // Get the value of the string parameter `name' from the ONELAB database.
        // Return an empty vector if the parameter does not exist.
        GMSH_API void getString(const std::string & name,
                                std::vector<std::string> & value);
    
        // Clear the ONELAB database, or remove a single parameter if `name' is given.
        GMSH_API void clear(const std::string & name = "");
    
        // Run a ONELAB client. If `name' is provided, create a new ONELAB client with
        // name `name' and executes `command'. If not, try to run a client that might
        // be linked to the processed input files.
        GMSH_API void run(const std::string & name = "",
                          const std::string & command = "");
    
      } // namespace onelab
    
      namespace logger { // Information logging functions
    
        // Write a `message'. `level' can be "info", "warning" or "error".
        GMSH_API void write(const std::string & message,
                            const std::string & level = "info");
    
        // Start logging messages.
        GMSH_API void start();
    
        // Get logged messages.
        GMSH_API void get(std::vector<std::string> & log);
    
        // Stop logging messages.
        GMSH_API void stop();
    
        // Return wall clock time.
        GMSH_API double time();
    
        // Return CPU time.
        GMSH_API double cputime();
    
      } // namespace logger
    
    } // namespace gmsh
    
    #endif