From f85f581e95fb71125c6f3dd590e37c0be8f5bcbd Mon Sep 17 00:00:00 2001 From: Andreas Dutzler Date: Fri, 3 Jan 2025 23:01:02 +0100 Subject: [PATCH] Add a docstring to `Element` --- docs/felupe/element.rst | 13 ++++++++++-- src/felupe/element/_base.py | 41 +++++++++++++++++++++++++++++++++++++ 2 files changed, 52 insertions(+), 2 deletions(-) diff --git a/docs/felupe/element.rst b/docs/felupe/element.rst index d65c2528..ec91e7c2 100644 --- a/docs/felupe/element.rst +++ b/docs/felupe/element.rst @@ -5,10 +5,14 @@ Element This module provides classes for the finite element formulations. -**Linear Elements** - .. currentmodule:: felupe +.. autosummary:: + + Element + +**Linear Elements** + .. autosummary:: Line @@ -50,6 +54,11 @@ This module provides classes for the finite element formulations. **Detailed API Reference** +.. autoclass:: felupe.Element + :members: + :undoc-members: + :show-inheritance: + .. autoclass:: felupe.Line :members: :undoc-members: diff --git a/src/felupe/element/_base.py b/src/felupe/element/_base.py index 0fd3c1af..877e5d8c 100644 --- a/src/felupe/element/_base.py +++ b/src/felupe/element/_base.py @@ -20,6 +20,47 @@ class Element: + """ + Base-class for a finite element which provides methods for plotting. + + Examples + -------- + This example shows how to implement a hexahedron element. + + .. pyvista-plot:: + :force_static: + + >>> import felupe as fem + >>> import numpy as np + >>> + >>> class Hexahedron(fem.Element): + ... def __init__(self): + ... a = [-1, 1, 1, -1, -1, 1, 1, -1] + ... b = [-1, -1, 1, 1, -1, -1, 1, 1] + ... c = [-1, -1, -1, -1, 1, 1, 1, 1] + ... self.points = np.vstack([a, b, c]).T + ... + ... # additional attributes for plotting, optional + ... self.cells = np.array([[0, 1, 2, 3, 4, 5, 6, 7]]) + ... self.cell_type = "hexahedron" + ... + ... def function(self, rst): + ... r, s, t = rst + ... a, b, c = self.points.T + ... ar, bs, ct = 1 + a * r, 1 + b * s, 1 + c * t + ... return (ar * bs * ct) / 8 + ... + ... def gradient(self, rst): + ... r, s, t = rst + ... a, b, c = self.points.T + ... ar, bs, ct = 1 + a * r, 1 + b * s, 1 + c * t + ... return np.stack([a * bs * ct, ar * b * ct, ar * bs * c], axis=1) + >>> + >>> mesh = fem.Cube(n=6) + >>> element = Hexahedron() + >>> quadrature = fem.GaussLegendre(order=1, dim=3) + >>> region = fem.Region(mesh, element, quadrature) + """ def view(self, point_data=None, cell_data=None, cell_type=None): """View the element with optional given dicts of point- and cell-data items.