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

Background

The UnxBoucWen2DLink class creates a uniaxial Bouc-Wen element with two-nodes in two-dimensions. The number of degree-of-freedom depend on the node. Currently node with 2 (two translation) and 3 (two translation and one rotation) dofs are supported. The main purpose of this element is to force uniaxial behavior along a certain direction. Figure provides a simple representation of this element where (1) and (2) represent the start and end nodes.

UnxBoucWen2DLink.png

The UnxBoucWen2DLink mass matrix in global coordinates is \(\textbf{M}^\textrm{e} = 0\). On the other hand, the UnxBoucWen2DLink solves the differential equation:

\[ \dot{z}(t)=(1 - (\beta \, \textrm{sign}[z(t)\dot{u}(t)] + \gamma)\|z(t)\|^\eta )\dot{u}(t) \,, \]

and the force is computed as

\[F(z,u) = q_y\,z + \alpha\,K_0\,U \,, \]

where \(q_y\) is the yield force of hysteretic component.

BoucWenParameters.png

The local to global transformation matrix is

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

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

\[ \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} \,, \]

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} = 4, \text{ or } 6\) (solid or structural).

REFERENCE:

  • Bouc, R. (1971). "Mathematical model for hysteresis." Report to the Centre de Recherches Physiques, pp16-25, Marseille, France.
  • Wen, Y.-K. (1976). for random vibration of hysteretic systems." Journal of Engineering Mechanics Division, 102(EM2), 249-263.

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 UnxBoucWen2DLink, using json format, use:

  • addElement(tag, name='UnxBoucWen2DLink', 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
      • 'alpha' Post-yield ratio
      • 'eta' Yielding exponent (sharpness of hysteresis loop corners)
      • 'beta' First hysteretic shape parameter
      • 'gamma' Second hysteretic shape parameter
      • 'tol' Newton-Raphson Tolerance.
      • 'Fy' Yield force
      • 'k0' Initial stiffness of hysteretic component
      • 'dir' The direction of action of this link

    Example

    A UNXBOUCWEN2DLINK element can be defined using the python interface as follows:
    SVL.addElement(tag=1, name='UnxBoucWen2DLink', conn=[1,2], attributes={'alpha': 0.1, 'eta': 1.0, 'beta': 0.5, 'gamma': 0.5, 'Fy': 250, 'k0': 250, 'dir': 1})

    Application
    Please refer to the A08-DY_2D_UniAxial_BoucWen_Link_Param1.py file located at 03-Validations/01-Debugging/ to see an example on how to define a UnxBoucWen2DLink 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/02-Link/UnxBoucWen2DLink.cpp file provides the class implementation. A UnxBoucWen2DLink element is created using the built-in json parse-structure provided in the Driver.hpp. A UnxBoucWen2DLink is defined inside the "Elements" json field indicating its "Tag" as follows,

  • {
        "Elements": {
            "Tag": {
                "name" : "UNXBOUCWEN2DLINK",
                "conn" : [ ],
                "attributes": {
                    "nmax": int,
                    "Fy": double,
                    "k0": double,
                    "alpha": double,
                    "eta": double,
                    "beta": double,
                    "gamma": double,
                    "tol": double,
                    "dim": double,
                    "dir": int
                }
            }
        }
    }
    
    Variable Description
    Tag Unique positive element identifier.
    conn The element connectivity node array.
    Fy The yield force.
    k0 The initial stiffness.
    alpha The stiffness ratio \(\alpha \in [0,1]\).
    eta The hysteresis shape parameter.
    beta The hysteresis shape parameter.
    gamma The hysteresis shape parameter.
    nmax Number of maximum iteration Newton-Raphson.
    tol Tolerance to accept convergence of solution.
    dim The number of degree of freedom of the nodes.
    dir The direction of action in local coordinates (1, or 2).

    Attention
    The element acts only in one direction (local coordinates).
    The element provides no-contribution to the mass matrix.
    The element only provides contribution to the stiffness matrix.
    The length of the link is not used.
    Strain and Stress link responses are relative deformation and non-linear force and in local coordinates respectively.
    Internal force responses are provided in global coordinates.
    Example

    A UNXBOUCWEN2DLINK element between nodes 1 and 2 with 2 dofs, acting on 1-direction in a 2D problem is constructed as:
    { "Elements": { "1": { "name" : "UNXBOUCWEN2DLINK", "conn" : [1,2], "attributes": { "nmax": 50, "Fy": 250.0, "k0": 250.0, "alpha": 0.1, "eta": 1.0, "beta": 0.5, "gamma": 0.5, "tol": 1.0E-08, "dim": 2, "dir": 0 } } } }