.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/structured/plot_2_boundary.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_examples_structured_plot_2_boundary.py: Contour and scatter from a boundary/patch ========================================= This example reads and plots the first component of an OpenFoam vector field from a boundary (patch) of a structured mesh .. GENERATED FROM PYTHON SOURCE LINES 10-15 First reads the mesh and print the shape/size of the mesh boundary ------------------------------------------------------------------ .. note:: It reads the mesh coordinates of a boundary for a structured mesh and stores them in variables x, y and z .. GENERATED FROM PYTHON SOURCE LINES 15-27 .. code-block:: Python # import readmesh function from fluidfoam package from fluidfoam import readmesh sol = "../../output_samples/box/" x, y, z = readmesh(path=sol, structured=True, boundary="topWall") nface = x.shape print("Boundary shape = ", nface) .. rst-class:: sphx-glr-script-out .. code-block:: none Reading file ../../output_samples/box//constant/polyMesh/faces Reading file ../../output_samples/box//constant/polyMesh/points Reading file ../../output_samples/box//constant/polyMesh/boundary Boundary shape = (20, 1, 15) .. GENERATED FROM PYTHON SOURCE LINES 28-33 Reads a vector field -------------------- .. note:: It reads a vector field of a boundary from a structured mesh and stores it in vel variable .. GENERATED FROM PYTHON SOURCE LINES 33-40 .. code-block:: Python # import readvector function from fluidfoam package from fluidfoam import readvector timename = "0" vel = readvector(sol, timename, "U", structured=True, boundary="topWall") .. rst-class:: sphx-glr-script-out .. code-block:: none Reading file ../../output_samples/box/0/U Warning : No data on boundary/patch Using the values of the nearest cells Reading file ../../output_samples/box//constant/polyMesh/boundary Reading file ../../output_samples/box//constant/polyMesh/owner Reading file ../../output_samples/box//constant/polyMesh/faces Reading file ../../output_samples/box//constant/polyMesh/points Reading file ../../output_samples/box//constant/polyMesh/boundary .. GENERATED FROM PYTHON SOURCE LINES 41-45 Now plots the contour of the first velocity component on the topWall boundary ----------------------------------------------------------------------------- .. note:: Here the topWall boundary is in (x, z) plane .. GENERATED FROM PYTHON SOURCE LINES 45-60 .. code-block:: Python import matplotlib.pyplot as plt import numpy as np plt.figure() levels = np.arange(0, np.max(vel[0]), 0.001) ax = plt.contourf(x[:, 0, :], z[:, 0, :], vel[0, :, 0, :], levels=levels) cbar = plt.colorbar(ax) cbar.set_label("Ux (m/s)") # Setting axis labels plt.xlabel("x (m)") plt.ylabel("z (m)") .. image-sg:: /auto_examples/structured/images/sphx_glr_plot_2_boundary_001.png :alt: plot 2 boundary :srcset: /auto_examples/structured/images/sphx_glr_plot_2_boundary_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Text(33.972222222222214, 0.5, 'z (m)') .. GENERATED FROM PYTHON SOURCE LINES 61-64 If you don't know the plane --------------------------- .. GENERATED FROM PYTHON SOURCE LINES 64-76 .. code-block:: Python fig = plt.figure() ax = fig.add_subplot(111, projection="3d") ax3d = ax.scatter(x, y, z, c=vel[0, :, :, :]) # Setting axis labels ax.set_xlabel("x (m)") ax.set_ylabel("y (m)") ax.set_zlabel("z (m)") cbar = plt.colorbar(ax3d) cbar.set_label("Ux (m/s)") .. image-sg:: /auto_examples/structured/images/sphx_glr_plot_2_boundary_002.png :alt: plot 2 boundary :srcset: /auto_examples/structured/images/sphx_glr_plot_2_boundary_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.402 seconds) .. _sphx_glr_download_auto_examples_structured_plot_2_boundary.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_2_boundary.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_2_boundary.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: plot_2_boundary.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_