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

Add BoundaryDict #915

Merged
merged 6 commits into from
Dec 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
1 change: 1 addition & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@ All notable changes to this project will be documented in this file. The format
- Add `SolidBody.evaluate.stress(field)` to evaluate the (first Piola-Kirchhoff) stress tensor (engineering stress in linear elasticity).
- Add a free-vibration modal analysis Step/Job `FreeVibration(items, boundaries)` with methods to evaluate `FreeVibration.evaluate()` and to extract `field, frequency = FreeVibration.extract(n)` its n-th result.
- Add `Boundary.plot()` to plot the points and prescribed directions of a boundary.
- Add `BoundaryDict` as a subclassed dict with methods to `plot()`, `screenshot()` and `imshow()`.

### Changed
- The first Piola-Kirchhoff stress tensor is evaluated if `ViewSolid(stress_type=None)`.
Expand Down
6 changes: 6 additions & 0 deletions docs/felupe/dof.rst
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,7 @@ This module contains the definition of a boundary condition, tools related to th
.. autosummary::

Boundary
BoundaryDict


**Tools**
Expand Down Expand Up @@ -37,6 +38,11 @@ This module contains the definition of a boundary condition, tools related to th
:undoc-members:
:inherited-members:

.. autoclass:: felupe.BoundaryDict
:members:
:undoc-members:
:inherited-members:

.. autofunction:: felupe.dof.partition

.. autofunction:: felupe.dof.apply
Expand Down
7 changes: 3 additions & 4 deletions docs/tutorial/examples/extut03_building_blocks.py
Original file line number Diff line number Diff line change
Expand Up @@ -171,14 +171,13 @@ def W(C, mu, bulk):
f0 = lambda x: np.isclose(x, 0)
f1 = lambda x: np.isclose(x, 1)

boundaries = {}
boundaries = fem.BoundaryDict()

boundaries["left"] = fem.Boundary(displacement, fx=f0)
boundaries["right"] = fem.Boundary(displacement, fx=f1, skip=(1, 0, 0))
boundaries["move"] = fem.Boundary(displacement, fx=f1, skip=(0, 1, 1), value=0.5)

plotter = boundaries["left"].plot(color="green")
plotter = boundaries["right"].plot(plotter=plotter, color="red")
boundaries["move"].plot(plotter=plotter, color="blue", point_size=5).show()
boundaries.plot().show()

# %%
# Partition of deegrees of freedom
Expand Down
3 changes: 2 additions & 1 deletion src/felupe/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -40,7 +40,7 @@
linear_elastic,
linear_elastic_plastic_isotropic_hardening,
)
from .dof import Boundary
from .dof import Boundary, BoundaryDict
from .element import ArbitraryOrderLagrange as ArbitraryOrderLagrangeElement
from .element import (
BiQuadraticQuad,
Expand Down Expand Up @@ -168,6 +168,7 @@
"linear_elastic",
"linear_elastic_plastic_isotropic_hardening",
"Boundary",
"BoundaryDict",
"ArbitraryOrderLagrangeElement",
"BiQuadraticQuad",
"ConstantHexahedron",
Expand Down
2 changes: 2 additions & 0 deletions src/felupe/dof/__init__.py
Original file line number Diff line number Diff line change
@@ -1,9 +1,11 @@
from ._boundary import Boundary
from ._dict import BoundaryDict
from ._loadcase import biaxial, shear, symmetry, uniaxial
from ._tools import apply, get_dof0, get_dof1, partition

__all__ = [
"Boundary",
"BoundaryDict",
"biaxial",
"shear",
"symmetry",
Expand Down
53 changes: 35 additions & 18 deletions src/felupe/dof/_boundary.py
Original file line number Diff line number Diff line change
Expand Up @@ -263,7 +263,16 @@ def update(self, value):
self.value = value #

def plot(
self, plotter=None, color="black", scale=0.125, point_size=10, width=3, **kwargs
self,
plotter=None,
color="black",
scale=0.125,
point_size=10,
width=3,
label=None,
show_points=True,
show_lines=True,
**kwargs,
):
"Plot the points and their prescribed directions of a boundary condition."

Expand All @@ -283,24 +292,32 @@ def plot(
)

if len(self.points) > 0:
magnitude = min(mesh.points.max(axis=0) - mesh.points.min(axis=0)) * scale
if label is None:
label = self.name

if show_points or show_lines:
points = np.pad(mesh.points, ((0, 0), (0, 3 - mesh.dim)))

if show_points:
_ = plotter.add_points(
points[self.points],
color=color,
point_size=point_size,
label=label,
)

_ = plotter.add_points(
mesh.points[self.points],
color=color,
point_size=point_size,
label=self.name,
)
if show_lines:
magnitude = (
min(mesh.points.max(axis=0) - mesh.points.min(axis=0)) * scale
)

for skip, direction in zip(self.skip, np.eye(mesh.dim)):
if not skip:
end = mesh.points[self.points] + direction * magnitude
_ = plotter.add_lines(
np.hstack([mesh.points[self.points], end]).reshape(
-1, mesh.dim
),
color=color,
width=width,
)
for skip, direction in zip(self.skip, np.eye(3)):
if not skip:
end = points[self.points] + direction * magnitude
_ = plotter.add_lines(
np.hstack([points[self.points], end]).reshape(-1, 3),
color=color,
width=width,
)

return plotter
79 changes: 79 additions & 0 deletions src/felupe/dof/_dict.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,79 @@
class BoundaryDict(dict):
"A dict of boundary conditions."

def plot(
self,
plotter=None,
colors=None,
size=(0.1, 0.1),
show_points=True,
show_lines=True,
**kwargs,
):
"Plot the boundary conditions."

if colors is None:
import matplotlib.colors as mcolors

colors = list(mcolors.TABLEAU_COLORS.values())

colors_list = []
while len(colors_list) < len(self.keys()):
colors_list = [*colors_list, *colors]

for (key, boundary), color in zip(self.items(), colors_list):
label = key
if boundary.name != "default":
label = boundary.name

plotter = boundary.plot(
label=label,
color=color,
plotter=plotter,
show_points=show_points,
show_lines=show_lines,
**kwargs,
)

plotter.add_legend(size=size, bcolor="white")

return plotter

def screenshot(
self,
*args,
filename="boundaries.png",
transparent_background=None,
scale=None,
colors=None,
plotter=None,
**kwargs,
):
"Take a screenshot of the boundary conditions."

if plotter is None:
mesh = self[list(self.keys())[0]].field.region.mesh
plotter = mesh.plot(off_screen=True)

Check warning on line 56 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L54-L56

Added lines #L54 - L56 were not covered by tests

plotter = self.plot(plotter=plotter, colors=colors, **kwargs)

Check warning on line 58 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L58

Added line #L58 was not covered by tests

return plotter.screenshot(

Check warning on line 60 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L60

Added line #L60 was not covered by tests
filename=filename,
transparent_background=transparent_background,
scale=scale,
)

def imshow(self, *args, ax=None, dpi=None, **kwargs):
"""Take a screenshot of the boundary conditions, show the image data in
a figure and return the ax.
"""

if ax is None:
import matplotlib.pyplot as plt

Check warning on line 72 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L71-L72

Added lines #L71 - L72 were not covered by tests

fig, ax = plt.subplots(dpi=dpi)

Check warning on line 74 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L74

Added line #L74 was not covered by tests

ax.imshow(self.screenshot(*args, filename=None, **kwargs))
ax.set_axis_off()

Check warning on line 77 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L76-L77

Added lines #L76 - L77 were not covered by tests

return ax

Check warning on line 79 in src/felupe/dof/_dict.py

View check run for this annotation

Codecov / codecov/patch

src/felupe/dof/_dict.py#L79

Added line #L79 was not covered by tests
14 changes: 14 additions & 0 deletions tests/test_dof.py
Original file line number Diff line number Diff line change
Expand Up @@ -180,8 +180,22 @@ def test_boundary_plot():
_ = boundaries["left"].plot()


def test_boundary_dict():
field = fem.FieldPlaneStrain(
fem.RegionQuad(fem.Rectangle(b=(3, 1), n=5)), dim=2
).as_container()
boundaries = fem.BoundaryDict(
left=fem.Boundary(field[0], fx=0, skip=(0, 0)),
right=fem.Boundary(field[0], name="my_label", fx=3, skip=(0, 1)),
)
plotter = boundaries.plot()
# img = boundaries.screenshot()
# ax = boundaries.imshow()


if __name__ == "__main__":
test_boundary()
test_boundary_dict()
test_boundary_multiaxial()
test_boundary_plot()
test_loadcase()