Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Bump versions and add streamlines #5

Merged
merged 2 commits into from
Feb 2, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion .github/workflows/run_tests.yml
Original file line number Diff line number Diff line change
Expand Up @@ -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

Expand Down
5 changes: 4 additions & 1 deletion febug/__init__.py
Original file line number Diff line number Diff line change
@@ -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)


Expand Down
50 changes: 50 additions & 0 deletions febug/plot.py
Original file line number Diff line number Diff line change
Expand Up @@ -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()
Expand Down
8 changes: 4 additions & 4 deletions setup.py
Original file line number Diff line number Diff line change
Expand Up @@ -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):
Expand Down
29 changes: 28 additions & 1 deletion test/test_plot.py
Original file line number Diff line number Diff line change
@@ -1,5 +1,6 @@
import numpy as np
import pytest
import pyvista
from mpi4py import MPI
import dolfinx
import dolfinx.mesh
Expand Down Expand Up @@ -44,14 +45,40 @@ 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):
u = dolfinx.fem.Function(dolfinx.fem.VectorFunctionSpace(mesh, ("CG", p)))
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):
Expand Down
Loading