+
Skip to content

Remove SolidBody._area_change #970

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

Merged
merged 1 commit into from
Jul 14, 2025
Merged
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
54 changes: 44 additions & 10 deletions src/felupe/mechanics/_solidbody.py
Original file line number Diff line number Diff line change
Expand Up @@ -22,7 +22,6 @@
import numpy as np

from ..assembly import IntegralForm
from ..constitution import AreaChange
from ..math import det, dot, transpose
from ..view import ViewSolid
from ._helpers import Assemble, Evaluate, Results
Expand Down Expand Up @@ -174,8 +173,7 @@ class SolidBody(Solid):

Notes
-----
The total potential energy of internal forces is given in Eq.
:eq:`solidbody`
The total potential energy of internal forces is given in Eq. :eq:`solidbody`

.. math::
:label: solidbody
Expand Down Expand Up @@ -269,11 +267,20 @@ def __init__(
self.results.kinematics = self._extract(self.field)

if statevars is not None:
# state variables are located in the results
self.results.statevars = statevars

else:
# init the shape of the state variables
statevars_shape = (0,)

# take the shape of the last item of the list if the provided umat has an
# x-attribute
if hasattr(umat, "x"):
statevars_shape = umat.x[-1].shape

# create and fill the array of state variables with zeros
# state variables are located in the results
self.results.statevars = np.zeros(
(
*statevars_shape,
Expand All @@ -296,10 +303,6 @@ def __init__(
kirchhoff_stress=self._kirchhoff_stress,
)

self._area_change = AreaChange()

self._form = IntegralForm

def _vector(
self,
field=None,
Expand All @@ -322,13 +325,18 @@ def _vector(
if apply is None:
apply = self.apply

# evaluate the (first Piola-Kirchhoff) stress tensor and store it in the results
self.results.stress = self._gradient(field, args=args, kwargs=kwargs)
self.results.force = self._form(

# assemble the internal force vector
# optionally, use only the first n items for mixed-field formulations
self.results.force = IntegralForm(
fun=self.results.stress[slice(items)],
v=self.field,
dV=self.field.region.dV,
).assemble(parallel=parallel, block=block)

# apply a callback on the assembled internal force vector
if apply is not None:
self.results.force = apply(self.results.force)

Expand Down Expand Up @@ -356,29 +364,43 @@ def _matrix(
if apply is None:
apply = self.apply

# evaluate the fourth-order elasticity tensor and store it in the results
# (associated to the first Piola-Kirchhoff stress tensor, i.e. the partial
# derivative of the first Piola-Kirchhoff stress tensor w.r.t. the deformation
# gradient tensor)
self.results.elasticity = self._hessian(field, args=args, kwargs=kwargs)

form = self._form(
# assemble the (sparse) tangent stiffness matrix
# optionally, use only the first n items for mixed-field formulations
form = IntegralForm(
fun=self.results.elasticity[slice(items)],
v=self.field,
u=self.field,
dV=self.field.region.dV,
)

# in a first step, integrate the weak-form and store the stiffness values
# (this dense array is only allocated once and will be re-used in all following
# evaluations)
self.results.stiffness_values = form.integrate(
parallel=parallel, out=self.results.stiffness_values
)

# finally, the dense array of stiffness-values is assembled into the sparse
# matrix
self.results.stiffness = form.assemble(
values=self.results.stiffness_values, block=block
)

# apply a callback on the assembled tangent stiffness matrix
if apply is not None:
self.results.stiffness = apply(self.results.stiffness)

return self.results.stiffness

def _extract(self, field):
"Evaluate and return the kinematics (the deformation gradient tensor)."

self.field = field
self.results.kinematics = self.field.extract(out=self.results.kinematics)

Expand All @@ -388,18 +410,25 @@ def _gradient(self, field=None, args=(), kwargs=None):
if kwargs is None:
kwargs = {}

# update the deformation gradient
if field is not None:
self.field = field
self.results.kinematics = self._extract(self.field)

if "out" in inspect.signature(self.umat.gradient).parameters:
kwargs["out"] = self.results.gradient

# evaluate the partial derivative of the strain energy density w.r.t. the
# deformation gradient
gradient = self.umat.gradient(
[*self.results.kinematics, self.results.statevars], *args, **kwargs
)
# store the first Piola-Kirchhoff stress tensor in the results
self.results.gradient = gradient[0]

# store all gradients and the temporary state variables in the results
# the state variables are only updated to the public results-attribute if the
# solution is valid
self.results.stress, self.results._statevars = gradient[:-1], gradient[-1]

return self.results.stress
Expand All @@ -415,9 +444,14 @@ def _hessian(self, field=None, args=(), kwargs=None):
if "out" in inspect.signature(self.umat.hessian).parameters:
kwargs["out"] = self.results.hessian

# evaluate the second partial derivative of the strain energy density w.r.t. the
# deformation gradient
self.results.elasticity = self.umat.hessian(
[*self.results.kinematics, self.results.statevars], *args, **kwargs
)

# store the partial derivative of the first Piola-Kirchhoff stress tensor w.r.t.
# the deformation gradient in the results
self.results.hessian = self.results.elasticity[0]

return self.results.elasticity
Expand Down Expand Up @@ -458,7 +492,7 @@ def _mass(self, density=None):
field = self.field[0].as_container()
dim = field[0].dim

form = self._form(
form = IntegralForm(
fun=[density * np.eye(dim).reshape(dim, dim, 1, 1)],
v=field,
u=field,
Expand Down
点击 这是indexloc提供的php浏览器服务,不要输入任何密码和下载