diff --git a/.github/workflows/run_tests.yml b/.github/workflows/run_tests.yml index e5705c8..ee551b3 100644 --- a/.github/workflows/run_tests.yml +++ b/.github/workflows/run_tests.yml @@ -8,7 +8,7 @@ on: jobs: Run-Tests: runs-on: ubuntu-latest - container: dolfinx/dolfinx:nightly + container: dolfinx/dolfinx:v0.7.3 steps: - uses: actions/checkout@v3 diff --git a/febug/__init__.py b/febug/__init__.py index 00165ef..82e504e 100644 --- a/febug/__init__.py +++ b/febug/__init__.py @@ -1,6 +1,9 @@ import sys from .plot import (plot_mesh, plot_function, plot_meshtags, plot_dofmap, - plot_warp, plot_quiver, plot_entity_indices, + plot_warp, plot_quiver, + plot_streamlines_evenly_spaced_2D, + plot_streamlines_from_source, + plot_entity_indices, plot_mesh_quality, plot_function_dofs) diff --git a/febug/plot.py b/febug/plot.py index c64b1f4..76d2d85 100644 --- a/febug/plot.py +++ b/febug/plot.py @@ -192,6 +192,56 @@ def plot_quiver(u: dolfinx.fem.Function, plotter: pyvista.Plotter=None, return plotter +def plot_streamlines_evenly_spaced_2D( + u: dolfinx.fem.Function, + start_position: typing.Sequence[float], + plotter: pyvista.Plotter=None, + **pv_args): + assert len(u.ufl_shape) == 1 + assert u.ufl_shape[0] in (1, 2, 3) + + if plotter is None: + plotter = pyvista.Plotter() + mesh = u.function_space.mesh + + grid = _to_pyvista_grid(u) + + streamlines = grid.streamlines_evenly_spaced_2D( + vectors=u.name, start_position=start_position, **pv_args) + plotter.add_mesh(streamlines) + + if mesh.geometry.dim == 2: + plotter.enable_parallel_projection() + plotter.view_xy() + + return plotter + + +def plot_streamlines_from_source( + u: dolfinx.fem.Function, + source: pyvista.DataSet, + plotter: pyvista.Plotter=None, + **pv_args): + assert len(u.ufl_shape) == 1 + assert u.ufl_shape[0] in (1, 2, 3) + + if plotter is None: + plotter = pyvista.Plotter() + mesh = u.function_space.mesh + + grid = _to_pyvista_grid(u) + + streamlines = grid.streamlines_from_source( + source, vectors=u.name, **pv_args) + plotter.add_mesh(streamlines) + + if mesh.geometry.dim == 2: + plotter.enable_parallel_projection() + plotter.view_xy() + + return plotter + + def plot_dofmap(V: dolfinx.fem.FunctionSpaceBase, plotter: pyvista.Plotter=None): if plotter is None: plotter = pyvista.Plotter() diff --git a/setup.py b/setup.py index 61203bb..34ea541 100644 --- a/setup.py +++ b/setup.py @@ -9,13 +9,13 @@ from setuptools import Extension, setup from setuptools.command.build_ext import build_ext -if sys.version_info < (3, 7): - print("Python 3.7 or higher required, please upgrade.") +if sys.version_info < (3, 10): + print("Python 3.10 or higher required, please upgrade.") sys.exit(1) -VERSION = "0.0.0" +VERSION = "2024.0.1" -REQUIREMENTS = ["pyvista", "fenics-dolfinx>0.6.0.dev0"] +REQUIREMENTS = ["pyvista", "fenics-dolfinx>0.7.0.dev0"] class CMakeExtension(Extension): diff --git a/test/test_plot.py b/test/test_plot.py index 318686a..6e2f4d9 100644 --- a/test/test_plot.py +++ b/test/test_plot.py @@ -1,5 +1,6 @@ import numpy as np import pytest +import pyvista from mpi4py import MPI import dolfinx import dolfinx.mesh @@ -44,7 +45,7 @@ def test_plot_warp(p, mesh): @pytest.mark.parametrize( "mesh", [pytest.param(msh, marks=pytest.mark.xfail( - msh.topology.dim == 1, reason="No warping 1D meshes")) + msh.topology.dim == 1, reason="No quivering 1D meshes")) for msh in meshes1D + meshes2D + meshes3D]) @pytest.mark.parametrize("p", [1, 2]) def test_plot_quiver(p, mesh): @@ -52,6 +53,32 @@ def test_plot_quiver(p, mesh): febug.plot_quiver(u) +@pytest.mark.parametrize( + "mesh", + [pytest.param(msh, marks=pytest.mark.xfail( + msh.topology.dim != 2, reason="Streamline function for 2D only")) + for msh in meshes1D + meshes2D + meshes3D]) +@pytest.mark.parametrize("p", [1, 2]) +def test_plot_streamlines_evenly_spaced_2D(p, mesh): + u = dolfinx.fem.Function(dolfinx.fem.VectorFunctionSpace(mesh, ("CG", p))) + febug.plot_streamlines_evenly_spaced_2D(u, start_position=(0.5, 0.5, 0.0)) + + +@pytest.mark.parametrize( + "mesh", + [pytest.param(msh, marks=pytest.mark.xfail( + msh.topology.dim == 1, reason="Streamlines not valid for 1D")) + for msh in meshes1D + meshes2D + meshes3D]) +@pytest.mark.parametrize("p", [1, 2]) +def test_plot_streamlines_from_source(p, mesh): + u = dolfinx.fem.Function(dolfinx.fem.VectorFunctionSpace(mesh, ("CG", p))) + pt_cloud_source = pyvista.PolyData(mesh.geometry.x) + febug.plot_streamlines_from_source(u, pt_cloud_source, max_time=1.0) + + mesh_source = febug.plot._to_pyvista_grid(mesh, mesh.topology.dim) + febug.plot_streamlines_from_source(u, mesh_source, max_time=1.0) + + @pytest.mark.parametrize("mesh", meshes1D + meshes2D + meshes3D) def test_plot_mesh(mesh): for tdim in range(mesh.topology.dim + 1):