Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations

Background

The lin3DFrame2 class creates a linearized three-dimensional frame element with two-nodes. Each node has six degrees-of-freedom, for which three degree-of-freedom represent three translation in the X-, Y-, and Z-directions, and the other three rotation about the X-, Y-, and Z-directions. Figure provides a simple representation of this element in which the local coordinates are represented by r-, s-, t-axis, and (1) and (2) represent the start and end nodes of the element respectively. The element mass, damping and stiffness matrices, as well as internal force vector are computed in local-coordinate and then transformed into global coordinates. No numerical integration is employed for this element.

lin3DFrame2.png

The coordinates of any point \(x_1,x_2,x_3\) in the element is written in terms of the shape functions \(\textbf{N}_\textrm{i}(r)\) in the reference coordinates and the values of the i-th element nodal coordinates \(x^\textrm{(i)}_\textrm{j}\) as follows

\[ x_\textrm{j}^\textrm{e}(r ) = \sum_\textrm{i = 1}^{\textrm{N}_\textrm{n}} x^\textrm{(i)}_\textrm{j} \, \textbf{N}_\textrm{i}^\textrm{e}(r) \,. \]

Computation of the strain-displacement matrix \(\textbf{B}^\textrm{e}(x_1)\) now requires to be mapped into the reference coordinate systems. Application of the chain rule for differentiation to the shape functions yields

\[ \frac{\partial N_\textrm{i}}{\partial r} = \frac{\partial x}{\partial r} \frac{\partial N_\textrm{i}}{\partial x} = \mathbf{J} \frac{\partial N_\textrm{i}}{\partial x} \,, \nonumber \]

where \(\mathbf{J}: \mathbb{R} \to \mathbb{R}\) is called the Jacobian matrix, here \(\mathbf{J} = \frac{L}{2}\). Hence, derivatives in \(\textbf{B}^\textrm{e}(x_1)\) can be transformed into \(\textbf{B}^\textrm{e}(r)\) in the reference coordinates using the Jacobian matrix.

The shape functions \(\textbf{N}_\textrm{i}^\textrm{e}(r)\) required to evaluate the shape function matrix \(\textbf{N}^\textrm{e}(r)\), the strain-displacement matrix \(\textbf{B}^\textrm{e}(r)\), and Jacobian matrix \(\textbf{J}^\textrm{e}(r)\) are:

\[ \textbf{N}_\textrm{1}^\textrm{e}(r ) = \frac{1}{2} (1 - r) \,, \\ \textbf{N}_\textrm{2}^\textrm{e}(r ) = \frac{1}{2} (1 + r) \,, \\ \textbf{N}_\textrm{3}^\textrm{e}(r ) = \frac{1}{4(1 + \phi_2)} (1 - r)^2(2 + r) + \frac{\phi_2}{2(1 + \phi_2)}(1 - r) \,, \\ \textbf{N}_\textrm{4}^\textrm{e}(r ) = \frac{L}{8(1 + \phi_2)} (1 - r)^2(1 + r) + \frac{L\, \phi_2}{8(1 + \phi_2)}(1 - r^2) \,, \\ \textbf{N}_\textrm{5}^\textrm{e}(r ) = \frac{1}{4(1 + \phi_2)} (1 + r)^2(2 - r) + \frac{\phi_2}{2(1 + \phi_2)}(1 + r) \,, \\ \textbf{N}_\textrm{6}^\textrm{e}(r ) = \frac{L}{8(1 + \phi_2)} (1 + r)^2(r - 1) - \frac{L\, \phi_2}{8(1 + \phi_2)}(1 - r^2) \,, \\ \textbf{N}_\textrm{7}^\textrm{e}(r ) = \frac{1}{4(1 + \phi_3)} (1 - r)^2(2 + r) + \frac{\phi_3}{2(1 + \phi_2)}(1 - r) \,, \\ \textbf{N}_\textrm{8}^\textrm{e}(r ) = \frac{L}{8(1 + \phi_3)} (1 - r)^2(1 + r) + \frac{L\, \phi_3}{8(1 + \phi_3)}(1 - r^2) \,, \\ \textbf{N}_\textrm{9}^\textrm{e}(r ) = \frac{1}{4(1 + \phi_3)} (1 + r)^2(2 - r) + \frac{\phi_3}{2(1 + \phi_3)}(1 + r) \,, \\ \textbf{N}_\textrm{10}^\textrm{e}(r ) = \frac{L}{8(1 + \phi_3)} (1 + r)^2(r - 1) - \frac{L\, \phi_3}{8(1 + \phi_3)}(1 - r^2) \,, \]

where the section shear factor is defined as \(\phi_2 = \frac{12 \textrm{EI}_\textrm{33}}{\textrm{GA}_\textrm{s2}}\), \(\phi_3 = \frac{12 \textrm{EI}_\textrm{22}}{\textrm{GA}_\textrm{s3}}\), with \(G\) the shear stiffness, and \(\textrm{A}_\textrm{s2}\), \(\textrm{A}_\textrm{s3}\) are the shear area along the 2 and 3-axes respectively.

The isoparametric transformation has the advantage that every element is treated in the same manner. In this regard, mass and stiffness matrices as well as force vector are evaluated using the same reference coordinates as follows

\[ \mathcal{M}^\textrm{e} = \int_{-1}^1 {\textbf{N}^\textrm{e}} (r)^\top \, \rho^\textrm{e} \, \textbf{N}^\textrm{e}(r) \det \mathbf{J}^\textrm{e} \, dr \,, \\ \mathcal{K}^\textrm{e} = \int_{-1}^1 {\textbf{B}^\textrm{e}} (r)^\top \, \mathbb{C}^\textrm{e} \, \textbf{B}^\textrm{e}(r) \det \mathbf{J}^\textrm{e} \, dr \,, \\ \mathcal{F}^\textrm{e} = \int_{-1}^1 {\textbf{N}^\textrm{e}} (r)^\top \, \textbf{F}_b^\textrm{e} \det \mathbf{J}^\textrm{e} \, dr \,. \]

and the mass and stiffness matrices as well as the force vector in global coordinates are defined as:

\[ \textbf{M}^\textrm{e} = \textbf{T}_\textrm{e}^\top \, \mathcal{M}^\textrm{e} \, \textbf{T}_\textrm{e} \,, \; \textbf{K}^\textrm{e} = \textbf{T}_\textrm{e}^\top \, \mathcal{K}^\textrm{e} \, \textbf{T}_\textrm{e} \,, \; \textbf{F}^\textrm{e} = \textbf{T}_\textrm{e}^\top \, \mathcal{F}^\textrm{e} \,, \]

In addition, the transformation matrix is

\[ \textbf{T}_\textrm{e} = \begin{bmatrix} \textbf{N} & 0 \\ 0 & \textbf{N} \end{bmatrix} \,, \; \text{ with } \textbf{N} = \begin{bmatrix} \hat{\textbf{N}} & 0 \\ 0 & \hat{\textbf{N}} \end{bmatrix} \,, \]

and \(\hat{\textbf{N}} = [\hat{\textrm{n}}_1, \hat{\textrm{n}}_2, \hat{\textrm{n}}_3]^\top\) is the unit matrix of local coordinate axes, and \(\hat{\textrm{n}}_\textrm{i}\) is the unit vector of the i-th local coordinate axis. Note that \(\hat{\textbf{N}} \in \mathbb{R}^{3 \times 3}\), and \(\hat{\textrm{n}}_\textrm{i} \in \mathbb{R}^3\). Note that \(\textbf{M}^\textrm{e}, \textbf{K}^\textrm{e} \in \mathbb{R}^{\textrm{N}_\textrm{dof}^\textrm{e} \times \textrm{N}_\textrm{dof}^\textrm{e}}\) with \(\textrm{N}_\textrm{dof}^\textrm{e} = 12\) in this case.

REFERENCE:

  • Bathe, K. Jurgen, "Finite Element Procedures", Chapter 5: pages 400-408. Prentice-Hall, 1996.
  • Lars Andersen and Søren R.K. Nielsen, "Elastic Beams in Three Dimensions", Aalborg University: Department of Civil Engineering, ISSN 1901-7286, DCE Lecture Notes No. 23.

Pre-Analysis

The python Pre-Analysis in the 01-Pre_Process/Method/Attach.py file provides with an interface to populate the Entities dictionary. This file contains several functions to populate specific fields. For example, to create a lin3DFrame2, using json format, use:

  • addElement(tag, name='lin3DFrame2', conn, attributes):

    • tag : The identifier of this element, i.e., tag > -1
    • name : Seismo-VLAB element class name
    • conn : Connectivity array of this element
    • attributes : Specific properties for the created element, for example
      • 'section' is section identifier.
      • 'form' is the beam formulation (Bernoulli, Timoshenko).
      • 'rule' The integration rule to apply: GAUSS, LOBATTO.
      • 'np' Number of quadrature-points used for the integration.

    Example

    A LIN3DFRAME2 element can be defined using the python interface as follows:
    SVL.addElement(tag=1, name='lin3DFrame2', conn=[1,2], attributes={'form': 'Timoshenko', 'section': 1, 'rule': 'LOBATTO', 'np': 5})

    Application
    Please refer to the D17-DY_Damped_WideFlange_3DPointLoad_Elastic_Frame2.py file located at 03-Validations/01-Debugging/ to see an example on how to define a lin3DFrame2 element.

On the contrary, the 01-Pre_Process/Method/Remove.py file provides with an interface to depopulate the Entities dictionary. For example, to remove an already define Element, use:

  • delElement(tag):
    • tag : The identifier of the element to be removed, i.e., tag > -1

Run-Analysis

The C++ Run-Analysis in the 02-Run_Process/04-Elements/04-Frame/lin3DFrame2.cpp file provides the class implementation. A lin3DFrame2 element is created using the built-in json parse-structure provided in the Driver.hpp. A lin3DFrame2 is defined inside the "Elements" json field indicating its "Tag" as follows,

  • {
        "Elements": {
            "Tag": {
                "name" : "LIN3DFRAME2",
                "conn" : [ ],
                "attributes": {
                    "np": int,
                    "rule": str,
                    "section": int,
                    "formulation": str
                }
            }
        }
    }
    
    Variable Description
    Tag Unique element object identifier.
    conn The element connectivity node array.
    np Number of quadrature-points used for the integration.
    rule The quadrature integration name rule=GAUSS, LOBATTO.
    section The section identifier.
    formulationThe frame model, it can be Bernoulli, or Timoshenko.

    Attention
    This formulation only allows small deformations.
    Only one section can be used to define this element.
    The section is assumed to be constant along the element.
    No integration is performed to evaluate mass and stiffness matrices.
    The domain reduction force method is incompatible with this element.
    In order to includes the effect of the shear stresses on the deformation use Timoshenko model.
    Example

    A LIN3DFRAME2 Bernoulli element between nodes 1 and 2, made of section 1 is constructed as:
    { "Elements": { "1": { "name" : "LIN3DFRAME2", "conn" : [1,2], "attributes": { "np": 5, "rule": "Lobatto", "section": 1, "formulation": "Bernoulli" } } } }