.. DO NOT EDIT.
.. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY.
.. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE:
.. "examples/02_gmsh-basic.py"
.. LINE NUMBERS ARE GIVEN BELOW.

.. only:: html

    .. note::
        :class: sphx-glr-download-link-note

        :ref:`Go to the end <sphx_glr_download_examples_02_gmsh-basic.py>`
        to download the full example code.

.. rst-class:: sphx-glr-example-title

.. _sphx_glr_examples_02_gmsh-basic.py:


Basic Gmsh Example
==================

In this example we'll create some basic geometries and turn them into meshes.
to illustrate some of the mesh generation features that Gmsh provides in
combination with polygon, point, and linestring geometries represented by
geopandas.

The :py:class:`GmshMesher` supports the geometry show in the basic Triangle
example and has a number of additional features.

.. GENERATED FROM PYTHON SOURCE LINES 14-22

.. code-block:: Python

    import geopandas as gpd
    import matplotlib.pyplot as plt
    import numpy as np
    import shapely.geometry as sg

    import pandamesh as pm









.. GENERATED FROM PYTHON SOURCE LINES 26-31

A simple rectangular mesh
-------------------------

The most simple example is perhaps a rectangle. We'll create a vector
geometry, store this in a geodataframe, and associate a cell size.

.. GENERATED FROM PYTHON SOURCE LINES 31-43

.. code-block:: Python


    polygon = sg.Polygon(
        [
            [0.0, 0.0],
            [10.0, 0.0],
            [10.0, 10.0],
            [0.0, 10.0],
        ]
    )
    gdf = gpd.GeoDataFrame(geometry=[polygon])
    gdf["cellsize"] = 2.0








.. GENERATED FROM PYTHON SOURCE LINES 44-48

We'll use this polygon to generate a mesh. We start by initializing a
TriangleMesher, which is a simple wrapper around the Python bindings to the
Gmsh C++-library. This wrapper extracts the coordinates and presents them
in the appropriate manner for Gmsh.

.. GENERATED FROM PYTHON SOURCE LINES 48-53

.. code-block:: Python


    mesher = pm.GmshMesher(gdf)
    vertices, triangles = mesher.generate()
    pm.plot(vertices, triangles)




.. image-sg:: /examples/images/sphx_glr_02_gmsh-basic_001.png
   :alt: 02 gmsh basic
   :srcset: /examples/images/sphx_glr_02_gmsh-basic_001.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 54-56

Before we can instantiate another GmshMesher, we need to ``finalize`` the old
one.

.. GENERATED FROM PYTHON SOURCE LINES 56-59

.. code-block:: Python


    mesher.finalize()








.. GENERATED FROM PYTHON SOURCE LINES 60-63

As the name suggests, Triangle only generates triangular meshes. Gmsh is
capable of generating quadrilateral-dominant meshes, and has a lot more bells
and whistles for defining cellsizes.

.. GENERATED FROM PYTHON SOURCE LINES 63-78

.. code-block:: Python


    line = sg.LineString([(2.0, 8.0), (8.0, 2.0)])
    gdf = gpd.GeoDataFrame(geometry=[polygon, line])
    gdf["cellsize"] = [2.0, 0.5]

    fig, (ax0, ax1) = plt.subplots(ncols=2)

    mesher = pm.TriangleMesher(gdf)
    vertices, triangles = mesher.generate()
    pm.plot(vertices, triangles, ax=ax0)

    mesher = pm.GmshMesher(gdf)
    vertices, triangles = mesher.generate()
    pm.plot(vertices, triangles, ax=ax1)




.. image-sg:: /examples/images/sphx_glr_02_gmsh-basic_002.png
   :alt: 02 gmsh basic
   :srcset: /examples/images/sphx_glr_02_gmsh-basic_002.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 79-88

Gmsh allows for specifying cell sizes in a more flexible way. Triangle (left)
only supports polygons (regions) with fixed cell sizes and explicitly placed
vertices. Gmsh is capable of forcing refinement in a larger zone around
features as is visible around the diagonal (right).

Defaults
--------

The GmshMesher class is initialized with a number of default parameters:

.. GENERATED FROM PYTHON SOURCE LINES 88-93

.. code-block:: Python


    print(mesher)

    mesher.finalize()





.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    GmshMesher
        initialized = True
        xoff = 5.0
        yoff = 5.0
        fields = []
        combination_field = CombinationField(fields=[], combination=<FieldCombination.MIN: 'Min'>, id=1, fields_list=[])
        tmpdir = <TemporaryDirectory '/tmp/tmpt4p7eiyc'>
        recombine_all = False
        mesh_size_extend_from_boundary = True
        mesh_size_from_points = True
        mesh_size_from_curvature = False
        field_combination = FieldCombination.MIN




.. GENERATED FROM PYTHON SOURCE LINES 94-101

The parameters of Gmsh differ from Triangle, but they work the same: they can
be altered after initialization to control the triangulation.

Forcing points, lines, local refinement
---------------------------------------

We can force points and lines into the triangulation:

.. GENERATED FROM PYTHON SOURCE LINES 101-145

.. code-block:: Python


    outer = [(0.0, 0.0), (10.0, 0.0), (10.0, 10.0), (0.0, 10.0)]
    inner = [(3.0, 3.0), (7.0, 3.0), (7.0, 7.0), (3.0, 7.0)]
    donut = sg.Polygon(shell=outer, holes=[inner])
    refined = sg.Polygon(inner)

    y = np.arange(0.5, 10.0, 0.5)
    x = np.full(y.size, 1.0)
    points = gpd.points_from_xy(x, y)

    line = sg.LineString(
        [
            [9.0, 2.0],
            [9.0, 8.0],
        ]
    )

    gdf = gpd.GeoDataFrame(geometry=[donut, refined, line, *points])
    gdf["cellsize"] = [2.0, 0.5, 2.0] + (len(points) * [2.0])

    mesher = pm.GmshMesher(gdf)
    vertices, triangles = mesher.generate()
    mesher.finalize()

    fig, ax = plt.subplots()
    pm.plot(vertices, triangles, ax=ax)
    gdf.plot(facecolor="none", edgecolor="red", ax=ax)


    # Quadrilateral meshes
    # --------------------
    #
    # One of the features of Gmsh is that it is also capable of generating
    # quadrilateral (dominant) meshes, by recombining triangles. We can achieve
    # this by changing a parameter on the mesher:

    gdf = gpd.GeoDataFrame(geometry=[polygon])
    gdf["cellsize"] = 2.0
    mesher = pm.GmshMesher(gdf)
    mesher.recombine_all = True
    vertices, faces = mesher.generate()

    pm.plot(vertices, faces)




.. rst-class:: sphx-glr-horizontal


    *

      .. image-sg:: /examples/images/sphx_glr_02_gmsh-basic_003.png
         :alt: 02 gmsh basic
         :srcset: /examples/images/sphx_glr_02_gmsh-basic_003.png
         :class: sphx-glr-multi-img

    *

      .. image-sg:: /examples/images/sphx_glr_02_gmsh-basic_004.png
         :alt: 02 gmsh basic
         :srcset: /examples/images/sphx_glr_02_gmsh-basic_004.png
         :class: sphx-glr-multi-img





.. GENERATED FROM PYTHON SOURCE LINES 146-151

Writing to file
---------------
It's also possible to use the Python bindings to write a Gmsh ``.msh`` file.
This file can be opened using the Gmsh GUI to e.g. inspect the generated
mesh.

.. GENERATED FROM PYTHON SOURCE LINES 151-154

.. code-block:: Python


    mesher.write("my-mesh.msh")








.. GENERATED FROM PYTHON SOURCE LINES 155-164

Conclusion
----------

In real use, the vector geometries will be more complex, and not based on
just a few coordinate pairs. Such cases are presented in the other examples,
but the same principles apply: we may use polygons, linestrings and points
with associated cell sizes to steer the triangulation; unlike Triangle,
for Gmsh cell sizes can associated to linestrings and points, not just
polygons.


.. rst-class:: sphx-glr-timing

   **Total running time of the script:** (0 minutes 0.362 seconds)


.. _sphx_glr_download_examples_02_gmsh-basic.py:

.. only:: html

  .. container:: sphx-glr-footer sphx-glr-footer-example

    .. container:: sphx-glr-download sphx-glr-download-jupyter

      :download:`Download Jupyter notebook: 02_gmsh-basic.ipynb <02_gmsh-basic.ipynb>`

    .. container:: sphx-glr-download sphx-glr-download-python

      :download:`Download Python source code: 02_gmsh-basic.py <02_gmsh-basic.py>`

    .. container:: sphx-glr-download sphx-glr-download-zip

      :download:`Download zipped: 02_gmsh-basic.zip <02_gmsh-basic.zip>`


.. only:: html

 .. rst-class:: sphx-glr-signature

    `Gallery generated by Sphinx-Gallery <https://sphinx-gallery.github.io>`_