:py:mod:`geomesher.mesher`
==========================

.. py:module:: geomesher.mesher

.. autoapi-nested-parse::

   Generate mesh using Gmsh.

   ..
       !! processed by numpydoc !!


Module Contents
---------------

.. py:class:: Mesher(gdf)


   
   Wrapper for the python bindings to Gmsh.

   This class must be initialized
   with a ``geopandas.GeoDataFrame`` containing at least one polygon, and a column
   named ``cellsize``.

   Optionally, multiple polygons with different cell sizes can be included in
   the geodataframe. These can be used to achieve local mesh remfinement.

   Linestrings and points may also be included. The segments of linestrings
   will be directly forced into the triangulation. Points can also be forced
   into the triangulation. Unlike Triangle, the cell size values associated
   with these geometries **will** be used.

   Gmsh cannot automatically resolve overlapping polygons, or points
   located exactly on segments. During initialization, the geometries of
   the geodataframe are checked:

   * Polygons should not have any overlap with each other.
   * Linestrings should not intersect each other.
   * Every linestring should be fully contained by a single polygon;
     a linestring may not intersect two or more polygons.
   * Linestrings and points should not "touch" or be located on
     polygon borders.
   * Holes in polygons are fully supported, but they must not contain
     any linestrings or points.

   If such cases are detected, the initialization will error.

   For more details on Gmsh, see:
   https://gmsh.info/doc/texinfo/gmsh.html

   A helpful index can be found near the bottom:
   https://gmsh.info/doc/texinfo/gmsh.html#Syntax-index















   ..
       !! processed by numpydoc !!
   .. py:property:: field_combination
      :type: Combinations

      
      Control how cell size fields are combined.

      When they are found at the same location. Accepted values are:

      - ``MIN``
      - ``MAX``
      - ``MEAN``















      ..
          !! processed by numpydoc !!

   .. py:property:: force_geometry
      :type: bool

      
      Force the geometry to be used in the mesh, or not per surface.
















      ..
          !! processed by numpydoc !!

   .. py:property:: mesh_algorithm
      :type: Algorithms

      
      Meshing algorithm to use.

      Available algorithms are:

      * ``MESH_ADAPT``
      * ``AUTOMATIC``
      * ``INITIAL_MESH_ONLY``
      * ``FRONTAL_DELAUNAY``
      * ``BAMG``
      * ``FRONTAL_DELAUNAY_FOR_QUADS``
      * ``PACKING_OF_PARALLELLOGRAMS``

      .. rubric:: Notes

      Each algorithm has its own advantages and disadvantages.

      For all 2D unstructured algorithms a Delaunay mesh that contains all
      the points of the 1D mesh is initially constructed using a
      divide-and-conquer algorithm. Missing edges are recovered using edge
      swaps. After this initial step several algorithms can be applied to
      generate the final mesh:

      * The MeshAdapt algorithm is based on local mesh modifications. This
        technique makes use of edge swaps, splits, and collapses: long edges
        are split, short edges are collapsed, and edges are swapped if a
        better geometrical configuration is obtained.
      * The Delaunay algorithm is inspired by the work of the GAMMA team at
        INRIA. New points are inserted sequentially at the circumcenter of
        the element that has the largest adimensional circumradius. The mesh
        is then reconnected using an anisotropic Delaunay criterion.
      * The Frontal-Delaunay algorithm is inspired by the work of S. Rebay.
      * Other experimental algorithms with specific features are also
        available. In particular, Frontal-Delaunay for Quads is a variant of
        the Frontal-Delaunay algorithm aiming at generating right-angle
        triangles suitable for recombination; and BAMG allows to generate
        anisotropic triangulations.

      For very complex curved surfaces the MeshAdapt algorithm is the most robust.
      When high element quality is important, the Frontal-Delaunay algorithm should
      be tried. For very large meshes of plane surfaces the Delaunay algorithm is
      the fastest; it usually also handles complex mesh size fields better than the
      Frontal-Delaunay. When the Delaunay or Frontal-Delaunay algorithms fail,
      MeshAdapt is automatically triggered. The Automatic algorithm uses
      Delaunay for plane surfaces and MeshAdapt for all other surfaces.















      ..
          !! processed by numpydoc !!

   .. py:property:: mesh_size_extend_from_boundary
      :type: bool

      
      Forces the mesh size to be extended from the boundary, or not per surface.
















      ..
          !! processed by numpydoc !!

   .. py:property:: mesh_size_from_curvature
      :type: bool

      
      Automatically compute mesh element sizes from curvature.

      It uses the value as the target number of elements per
      :math:`2 \pi` radians.















      ..
          !! processed by numpydoc !!

   .. py:property:: mesh_size_from_points
      :type: bool

      
      Compute mesh element sizes from values given at geometry points.
















      ..
          !! processed by numpydoc !!

   .. py:property:: recombine_all
      :type: bool

      
      Apply recombination algorithm to all surfaces, ignoring per-surface spec.
















      ..
          !! processed by numpydoc !!

   .. py:property:: subdivision_algorithm
      :type: Subdivisions

      
      Subdivision algorithm.

      All meshes can be subdivided to generate fully quadrangular cells.
      Available algorithms are:

      - ``NONE``
      - ``ALL_QUADRANGLES``
      - ``BARYCENTRIC``















      ..
          !! processed by numpydoc !!

   .. py:method:: add_distance_field(gdf, minimum_cellsize)

      
      Add a distance field to the mesher.

      The of geometry of these fields are not forced into the mesh, but they
      can be used to specify zones of with cell sizes.

      :Parameters: * **gdf** (:class:`geopandas.GeoDataFrame`) -- Location and cell size of the fields, as vector data.
                   * **minimum_cellsize** (:class:`float`) -- Minimum cell size.















      ..
          !! processed by numpydoc !!

   .. py:method:: add_structured_field(cellsize, xmin, ymin, dx, dy)

      
      Add a structured field specifying cell sizes.

      Gmsh will interpolate between the points to determine the desired cell size.

      :Parameters: * **cellsize** (:class:`FloatArray`) -- Specifies the cell size on a structured grid. The location of this grid
                     is determined by ``xmin``, ``ymin``, ``dx``, ``dy``.
                   * **xmin** (:class:`float`) -- x-origin.
                   * **ymin** (:class:`float`) -- y-origin.
                   * **dx** (:class:`float`) -- Spacing along the x-axis.
                   * **dy** (:class:`float`) -- Spacing along the y-axis.















      ..
          !! processed by numpydoc !!

   .. py:method:: clear_fields()

      
      Clear all cell size fields from the mesher.
















      ..
          !! processed by numpydoc !!

   .. py:method:: finalize_gmsh()
      :staticmethod:

      
      Finalize Gmsh.
















      ..
          !! processed by numpydoc !!

   .. py:method:: generate()

      
      Generate the mesh and return it as a geopandas.GeoDataFrame.
















      ..
          !! processed by numpydoc !!

   .. py:method:: run_gmsh()

      
      Generate a mesh of triangles or quadrangles.

      :returns: * **vertices** (:class:`numpy.ndarray`) -- Coordinates of mesh vertices with shape ``(n_vertex, 2)``
                * **faces** (:class:`numpy.ndarray`) -- Index number of mesh faces  with shape ``(n_face, nmax_per_face)``.
                  Where ``nmax_per_face`` is 3 for exclusively triangles and 4 if
                  quadrangles are included. A fill value of -1 is used as a last
                  entry for triangles in that case.















      ..
          !! processed by numpydoc !!

   .. py:method:: write(path)

      
      Write a gmsh .msh file.

      :Parameters: **path** (:class:`str` or :class:`pathlib.Path`)















      ..
          !! processed by numpydoc !!


.. py:function:: gdf_mesher(gdf, meshing_algorithm = 'AUTOMATIC', extensive_variables = None, intensive_variables = None, categorical_variables = None)

   
   Generate a mesh from a geodataframe.

   This function uses default Gmsh parameters. For more control, use
   :class:`Mesher`.

   :Parameters: * **gdf** (:class:`geopandas.GeoDataFrame`) -- The input must have at least one polygon and a column named ``cellsize``.
                  Optionally, multiple polygons with different cell sizes can be included in
                  the geodataframe. These can be used to achieve local mesh remfinement.
                * **meshing_algorithm** (:class:`str`, *optional*) -- Meshing algorithm to use. Available algorithms are:

                  * ``MESH_ADAPT``
                  * ``AUTOMATIC``
                  * ``INITIAL_MESH_ONLY``
                  * ``FRONTAL_DELAUNAY``
                  * ``BAMG``
                  * ``FRONTAL_DELAUNAY_FOR_QUADS``
                  * ``PACKING_OF_PARALLELLOGRAMS``

                  For more details, see :attr:`Mesher.mesh_algorithm`.
                * **extensive_variables** (:class:`list`, *optional*) -- Columns in dataframes for extensive variables for remapping,
                  defaults to ``None``.
                * **intensive_variables** (:class:`list`, *optional*) -- Columns in dataframes for intensive variables for remapping,
                  defaults to ``None``.
                * **categorical_variables** (:class:`list`, *optional*) -- Columns in dataframes for categorical variables for remapping,
                  defaults to ``None``.

   :returns: **mesh** (:class:`geopandas.GeoDataFrame`) -- The mesh as a geopandas.GeoDataFrame.

   .. rubric:: Notes

   Linestrings and points may also be included. The segments of linestrings
   will be directly forced into the triangulation. Points can also be forced
   into the triangulation. Unlike Triangle, the cell size values associated
   with these geometries **will** be used.

   Gmsh cannot automatically resolve overlapping polygons, or points
   located exactly on segments. During initialization, the geometries of
   the geodataframe are checked:

   * Polygons should not have any overlap with each other.
   * Linestrings should not intersect each other.
   * Every linestring should be fully contained by a single polygon;
     a linestring may not intersect two or more polygons.
   * Linestrings and points should not "touch" or be located on polygon borders.
   * Holes in polygons are fully supported, but they must not contain
     any linestrings or points.

   If such cases are detected, the initialization will throw an error.

   For more details on Gmsh, see:
   https://gmsh.info/doc/texinfo/gmsh.html

   A helpful index can be found near the bottom:
   https://gmsh.info/doc/texinfo/gmsh.html#Syntax-index















   ..
       !! processed by numpydoc !!

.. py:function:: gmsh_env()

   
   Context manager for Gmsh initialization and finalization.
















   ..
       !! processed by numpydoc !!

