Skip to content
Snippets Groups Projects
Select Git revision
  • f9da1a3db25181729c613040fead7e4bf533b4ce
  • master default protected
  • patches-4.14
  • steplayer
  • bl
  • pluginMeshQuality
  • fixBugsAmaury
  • hierarchical-basis
  • alphashapes
  • relaying
  • new_export_boris
  • oras_vs_osm
  • reassign_partitions
  • distributed_fwi
  • rename-classes
  • fix/fortran-api-example-t4
  • robust_partitions
  • reducing_files
  • fix_overlaps
  • 3115-issue-fix
  • 3023-Fillet2D-Update
  • gmsh_4_14_0
  • gmsh_4_13_1
  • gmsh_4_13_0
  • gmsh_4_12_2
  • gmsh_4_12_1
  • gmsh_4_12_0
  • gmsh_4_11_1
  • gmsh_4_11_0
  • gmsh_4_10_5
  • gmsh_4_10_4
  • gmsh_4_10_3
  • gmsh_4_10_2
  • gmsh_4_10_1
  • gmsh_4_10_0
  • gmsh_4_9_5
  • gmsh_4_9_4
  • gmsh_4_9_3
  • gmsh_4_9_2
  • gmsh_4_9_1
  • gmsh_4_9_0
41 results

api.texi

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