:orphan: .. _tutorials: 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 :class:`cube ` with hyperelastic material behaviour is subjected to a :func:`~felupe.dof.uniaxial` elongation applied at a clamped end-face. First, let’s import FElupe and create a meshed :class:`cube ` out of :class:`hexahedron ` cells with a given number of points per axis. A numeric :class:`region `, pre-defined for hexahedrons, is created on the mesh. A vector-valued displacement :class:`field ` is initiated on the region. Next, a :class:`field container ` is created on top of this field. A :func:`~felupe.dof.uniaxial` load case is applied on the displacement :class:`field ` stored inside the :class:`field container `. This involves setting up :func:`~felupe.dof.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 :class:`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 :class:`Ogden-Roxburgh ` Mullins-softening model formulation in combination with an isotropic hyperelastic :class:`Neo-Hookean ` material formulation is applied on the displacement :class:`field ` of a :class:`nearly-incompressible solid body `. A :class:`step ` generates the consecutive substep-movements of a given :class:`boundary ` condition. The :class:`step ` is further added to a list of steps of a :class:`job ` πŸ‘©β€πŸ’» (here, a :class:`characteristic curve ` πŸ“ˆ job is used). During :meth:`evaluation ` ⏳, each substep of each :class:`step ` is solved by an iterative :func:`Newton-Rhapson ` procedure βš–οΈ. The :func:`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. .. tab:: 3D .. tab:: Hexahedron .. code-block:: python 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() .. tab:: Quadratic Hexahedron .. code-block:: python 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() .. tab:: Lagrange Hexahedron .. code-block:: python 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() .. tab:: Plane Strain .. tab:: Quad .. code-block:: python 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() .. tab:: Axisymmetric .. tab:: Quad .. code-block:: python 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. .. raw:: html
.. thumbnail-parent-div-open .. raw:: html
.. only:: html .. image:: /tutorial/images/thumb/sphx_glr_extut01_getting_started_thumb.png :alt: :ref:`sphx_glr_tutorial_extut01_getting_started.py` .. raw:: html
Getting Started
.. raw:: html
.. only:: html .. image:: /tutorial/images/thumb/sphx_glr_extut02_job_thumb.png :alt: :ref:`sphx_glr_tutorial_extut02_job.py` .. raw:: html
Run a Job
.. raw:: html
.. only:: html .. image:: /tutorial/images/thumb/sphx_glr_extut03_building_blocks_thumb.png :alt: :ref:`sphx_glr_tutorial_extut03_building_blocks.py` .. raw:: html
Building Blocks
.. raw:: html
.. only:: html .. image:: /tutorial/images/thumb/sphx_glr_extut04_elementary_deformations_thumb.png :alt: :ref:`sphx_glr_tutorial_extut04_elementary_deformations.py` .. raw:: html
Non-homogeneous Deformations
.. thumbnail-parent-div-close .. raw:: html
.. toctree:: :hidden: /tutorial/extut01_getting_started /tutorial/extut02_job /tutorial/extut03_building_blocks /tutorial/extut04_elementary_deformations .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_