Beginner’s Guide#

This minimal code-block covers the essential high-level parts of creating and solving problems with FElupe. As an introductory example 👨‍🏫, a quarter model of a solid cube with hyperelastic material behaviour is subjected to a uniaxial() elongation applied at a clamped end-face.

First, let’s import FElupe and create a meshed cube out of hexahedron cells with a given number of points per axis. A numeric region, pre-defined for hexahedrons, is created on the mesh. A vector-valued displacement field is initiated on the region. Next, a field container is created on top of this field.

A uniaxial() load case is applied on the displacement field stored inside the field container. This involves setting up symmetry() planes as well as the absolute value of the prescribed displacement at the mesh-points on the right-end face of the cube. The right-end face is clamped 🛠️: only displacements in direction x are allowed. The dict of boundary conditions for this pre-defined load case are returned as boundaries and the partitioned degrees of freedom as well as the external displacements are stored within the returned dict loadcase.

An isotropic pseudo-elastic Ogden-Roxburgh Mullins-softening model formulation in combination with an isotropic hyperelastic Neo-Hookean material formulation is applied on the displacement field of a nearly-incompressible solid body.

A step generates the consecutive substep-movements of a given boundary condition. The step is further added to a list of steps of a job 👩‍💻 (here, a characteristic curve 📈 job is used). During evaluation ⏳, each substep of each step is solved by an iterative Newton-Rhapson procedure ⚖️. The solution is exported after each completed substep as a time-series ⌚ XDMF file. Finally, the result of the last completed substep is plotted.

Slightly modified code-blocks are provided for different kind of analyses and element formulations.

import felupe as fem

mesh = fem.Cube(n=6)
region = fem.RegionHexahedron(mesh)
field = fem.FieldContainer([fem.Field(region, dim=3)])

boundaries, loadcase = fem.dof.uniaxial(field, clamped=True)

umat = fem.NeoHooke(mu=1)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)

move = fem.math.linsteps([0, 1], num=5)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)

job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(filename="result.xdmf")
fig, ax = job.plot(
    xlabel="Displacement $d_1$ in mm $\longrightarrow$",
    ylabel="Normal Force $F_1$ in N $\longrightarrow$",
)
solid.plot(
    "Principal Values of Cauchy Stress"
).show()
import felupe as fem

mesh = fem.Cube(n=(9, 5, 5)).add_midpoints_edges()
region = fem.RegionQuadraticHexahedron(mesh)
field = fem.FieldContainer([fem.Field(region, dim=3)])

boundaries, loadcase = fem.dof.uniaxial(field, clamped=True)

umat = fem.NeoHooke(mu=1, bulk=50)
solid = fem.SolidBody(umat, field)

move = fem.math.linsteps([0, 1], num=5)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)

job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(parallel=True)
fig, ax = job.plot(
    xlabel="Displacement $u$ in mm $\longrightarrow$",
    ylabel="Normal Force $F$ in N $\longrightarrow$",
)
solid.plot(
    "Principal Values of Cauchy Stress", project=fem.topoints, nonlinear_subdivision=4
).show()
import felupe as fem

mesh = fem.mesh.CubeArbitraryOrderHexahedron(order=6)
region = fem.RegionLagrange(mesh, order=6, dim=3)
field = fem.FieldContainer([fem.Field(region, dim=3)])

boundaries, loadcase = fem.dof.uniaxial(field, clamped=True)

umat = fem.NeoHooke(mu=1, bulk=50)
solid = fem.SolidBody(umat, field)

move = fem.math.linsteps([0, 1], num=5)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)

job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(parallel=True)
fig, ax = job.plot(
    xlabel="Displacement $u$ in mm $\longrightarrow$",
    ylabel="Normal Force $F$ in N $\longrightarrow$",
)
solid.plot(
    "Principal Values of Cauchy Stress", project=fem.topoints, nonlinear_subdivision=4
).show()
import felupe as fem

mesh = fem.Rectangle(n=6)
region = fem.RegionQuad(mesh)
field = fem.FieldContainer([fem.FieldPlaneStrain(region, dim=2)])

boundaries, loadcase = fem.dof.uniaxial(field, clamped=True)

umat = fem.NeoHooke(mu=1)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)

move = fem.math.linsteps([0, 1], num=5)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)

job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(filename="result.xdmf")
fig, ax = job.plot(
    xlabel="Displacement $d_1$ in mm $\longrightarrow$",
    ylabel="Normal Force $F_1$ in N $\longrightarrow$",
)
solid.plot(
    "Principal Values of Cauchy Stress"
).show()
import felupe as fem

mesh = fem.Rectangle(n=6)
region = fem.RegionQuad(mesh)
field = fem.FieldContainer([fem.FieldAxisymmetric(region, dim=2)])

boundaries, loadcase = fem.dof.uniaxial(field, clamped=True)

umat = fem.NeoHooke(mu=1)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)

move = fem.math.linsteps([0, 1], num=5)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)

job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(filename="result.xdmf")
fig, ax = job.plot(
    xlabel="Displacement $d_1$ in mm $\longrightarrow$",
    ylabel="Normal Force $F_1$ in N $\longrightarrow$",
)
solid.plot(
    "Principal Values of Cauchy Stress"
).show()

Tutorials#

This section is all about learning. Each tutorial focusses on some lessons to learn.

Getting Started

Getting Started

Run a Job

Run a Job

Building Blocks

Building Blocks

Non-homogeneous Deformations

Non-homogeneous Deformations

Gallery generated by Sphinx-Gallery