Skip to content

Commit d482e97

Browse files
Merge pull request #343 from TLCFEM/correct-warping
Add `plot_warping_function()`, enforce Lagrangian Multiplier constraint
2 parents 900b3f4 + ce1d000 commit d482e97

File tree

2 files changed

+77
-13
lines changed

2 files changed

+77
-13
lines changed

src/sectionproperties/analysis/fea.py

Lines changed: 15 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -25,28 +25,33 @@
2525
def _assemble_torsion(
2626
k_el: np.ndarray,
2727
f_el: np.ndarray,
28+
c_el: np.ndarray,
29+
n: np.ndarray,
2830
b: np.ndarray,
2931
weight: float,
3032
nx: float,
3133
ny: float,
32-
) -> tuple[np.ndarray, np.ndarray]:
34+
) -> tuple[np.ndarray, np.ndarray, np.ndarray]:
3335
"""Utility function for calculating the torsion stiffness matrix and load vector.
3436
3537
Args:
3638
k_el: Element stiffness matrix
3739
f_el: Element load vector
40+
c_el: Element constraint vector
41+
n: Shape function
3842
b: Strain matrix
3943
weight: Effective weight
4044
nx: Global x-coordinate
4145
ny: Global y-coordinate
4246
4347
Returns:
44-
Torsion stiffness matrix and load vector (``k_el``, ``f_el``)
48+
Torsion stiffness matrix and load vector (``k_el``, ``f_el``, ``c_el``)
4549
"""
4650
# calculated modulus weighted stiffness matrix and load vector
4751
k_el += weight * b.transpose() @ b
4852
f_el += weight * b.transpose() @ np.array([ny, -nx])
49-
return k_el, f_el
53+
c_el += weight * n
54+
return k_el, f_el, c_el
5055

5156

5257
@njit(cache=True, nogil=True)
@@ -272,30 +277,32 @@ def geometric_properties(
272277
self.material.density,
273278
)
274279

275-
def torsion_properties(self) -> tuple[np.ndarray, np.ndarray]:
280+
def torsion_properties(self) -> tuple[np.ndarray, np.ndarray, np.ndarray]:
276281
"""Calculates the element warping stiffness matrix and the torsion load vector.
277282
278283
Returns:
279-
Element stiffness matrix ``k_el`` and element torsion load vector ``f_el``
284+
Element stiffness matrix ``k_el``, element torsion load vector ``f_el``
285+
and element constraint vector ``c_el``
280286
"""
281287
# initialise stiffness matrix and load vector
282288
k_el = np.zeros(shape=(6, 6), dtype=float)
283289
f_el = np.zeros(shape=6, dtype=float)
290+
c_el = np.zeros(shape=6, dtype=float)
284291

285292
# Gauss points for 4 point Gaussian integration
286293
gps = gauss_points(n=4)
287294

288295
for gp in gps:
289296
# determine shape function, shape function derivative,
290297
# jacobian and global coordinates
291-
_, b, j, nx, ny = shape_function(coords=self.coords, gauss_point=gp)
298+
n, b, j, nx, ny = shape_function(coords=self.coords, gauss_point=gp)
292299

293300
weight = gp[0] * j * self.material.elastic_modulus
294301

295302
# calculated modulus weighted stiffness matrix and load vector
296-
k_el, f_el = _assemble_torsion(k_el, f_el, b, weight, nx, ny)
303+
k_el, f_el, c_el = _assemble_torsion(k_el, f_el, c_el, n, b, weight, nx, ny)
297304

298-
return k_el, f_el
305+
return k_el, f_el, c_el
299306

300307
def shear_load_vectors(
301308
self,

src/sectionproperties/analysis/section.py

Lines changed: 62 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -1411,17 +1411,19 @@ def assemble_torsion(
14111411
col = [] # list holding column indices
14121412
data = [] # list holding stiffness matrix entries
14131413
f_torsion = np.zeros(n_size) # force vector array
1414+
c_torsion = np.zeros(n_size) # constraint vector array
14141415

14151416
# loop through all elements in the mesh
14161417
for el in self.elements:
14171418
# determine number of nodes in the current element
14181419
n = len(el.node_ids)
14191420

14201421
# calculate the element stiffness matrix and torsion load vector
1421-
k_el, f_el = el.torsion_properties()
1422+
k_el, f_el, c_el = el.torsion_properties()
14221423

14231424
# assemble the torsion load vector
14241425
f_torsion[el.node_ids] += f_el
1426+
c_torsion[el.node_ids] += c_el
14251427

14261428
for node_id in el.node_ids:
14271429
row.extend([node_id] * n)
@@ -1432,15 +1434,15 @@ def assemble_torsion(
14321434
progress.update(task_id=task, advance=1)
14331435

14341436
# construct Lagrangian multiplier matrix:
1435-
# column vector of ones
1437+
# column vector
14361438
row.extend(range(n_size))
14371439
col.extend([n_size] * n_size)
1438-
data.extend([1] * n_size)
1440+
data.extend(c_torsion)
14391441

1440-
# row vector of ones
1442+
# row vector
14411443
row.extend([n_size] * n_size)
14421444
col.extend(range(n_size))
1443-
data.extend([1] * n_size)
1445+
data.extend(c_torsion)
14441446

14451447
k_lg = coo_matrix(
14461448
(data, (row, col)), shape=(n_size + 1, n_size + 1), dtype=float
@@ -1673,6 +1675,61 @@ def plot_centroids(
16731675

16741676
return ax
16751677

1678+
def plot_warping_function(
1679+
self,
1680+
title: str = "Warping Function",
1681+
level: int = 20,
1682+
cmap: str = "viridis",
1683+
with_lines: bool = True,
1684+
**kwargs,
1685+
):
1686+
r"""Plots the warping function over the mesh.
1687+
1688+
Args:
1689+
title: Plot title
1690+
level: Number of contour levels
1691+
cmap: Colormap
1692+
with_lines: If set to True, contour lines are displayed
1693+
kwargs: Passed to :func:`~sectionproperties.post.post.plotting_context`
1694+
1695+
Raises:
1696+
RuntimeError: If a warping analysis has not been performed
1697+
1698+
Returns:
1699+
Matplotlib axes object
1700+
"""
1701+
if self.section_props.omega is None:
1702+
raise RuntimeError(
1703+
"Perform a warping analysis before plotting the warping function."
1704+
)
1705+
1706+
# create plot and setup the plot
1707+
with post.plotting_context(title=title, **kwargs) as (fig, ax):
1708+
assert ax
1709+
1710+
loc = self.mesh["vertices"]
1711+
1712+
if with_lines:
1713+
ax.tricontour(
1714+
loc[:, 0],
1715+
loc[:, 1],
1716+
self.section_props.omega,
1717+
colors="k",
1718+
levels=level,
1719+
)
1720+
1721+
ax.tricontourf(
1722+
loc[:, 0],
1723+
loc[:, 1],
1724+
self.section_props.omega,
1725+
cmap=cmap,
1726+
levels=level,
1727+
)
1728+
ax.set_xlabel("X")
1729+
ax.set_ylabel("Y")
1730+
1731+
return ax
1732+
16761733
def display_mesh_info(self) -> None:
16771734
"""Prints mesh statistics to the command line."""
16781735
console = Console()

0 commit comments

Comments
 (0)