BeamLR - Linear Timoshenko 3D beam element with reduced integration (pyfe3d.beamlr)#

class pyfe3d.beamlr.BeamLR#

Timoshenko 3D beam element with linear shape functions

Formulation based on reference, replacing the consistent shape functions by linear functions and performing numerical integration with just 1 point at the beam center:

Luo, Y., 2008, “An Efficient 3D Timoshenko Beam Element with Consistent Shape Functions,” Adv. Theor. Appl. Mech., 1(3), pp. 95–106.

Nodal connectivity for the beam element:

^ y axis
|
|
______   --> x axis
1    2
Attributes:
eid,int

Element identification number.

pid,int

Property identification number.

length,double

Element length.

r11, r12, r13, r21, r22, r23, r31, r32, r33double

Rotation matrix from local to global coordinates.

vxyi, vxyj, vxykdouble

Components of a vector on the \(XY\) plane of the element coordinate system, defined using global coordinates.

c1, c2int

Position of each node in the global stiffness matrix.

n1, n2int

Node identification number.

init_k_KC0, init_k_KG, init_k_Mint

Position in the arrays storing the sparse data for the structural matrices.

probeBeamLRProbe object

probe: pyfe3d.beamlr.BeamLRProbe

Methods

update_KC0(self, long[, long[, double[, ...)

Update sparse vectors for linear constitutive stiffness matrix KC0

update_KG(self, long[, long[, double[, ...)

Update sparse vectors for geometric stiffness matrix KG

update_M(self, long[, long[, double[, ...)

Update sparse vectors for mass matrix M

update_length(self)

Update element length

update_probe_ue(self, double[)

Update the local displacement vector of the probe of the element

update_probe_xe(self, double[)

Update the 3D coordinates of the probe of the element

update_rotation_matrix(self, double vxyi, ...)

Update the rotation matrix of the element

c1#

c1: ‘int’

c2#

c2: ‘int’

eid#

eid: ‘int’

init_k_KC0#

init_k_KC0: ‘int’

init_k_KG#

init_k_KG: ‘int’

init_k_M#

init_k_M: ‘int’

length#

length: ‘double’

n1#

n1: ‘int’

n2#

n2: ‘int’

pid#

pid: ‘int’

probe#

probe: pyfe3d.beamlr.BeamLRProbe

r11#

r11: ‘double’

r12#

r12: ‘double’

r13#

r13: ‘double’

r21#

r21: ‘double’

r22#

r22: ‘double’

r23#

r23: ‘double’

r31#

r31: ‘double’

r32#

r32: ‘double’

r33#

r33: ‘double’

update_KC0(self, long[::1] KC0r, long[::1] KC0c, double[::1] KC0v, BeamProp prop, int update_KC0v_only=0) void#

Update sparse vectors for linear constitutive stiffness matrix KC0

Parameters:
KC0rnp.array

Array to store row positions of sparse values

KC0cnp.array

Array to store column positions of sparse values

KC0vnp.array

Array to store sparse values

propBeamProp object

Beam property object from where the stiffness and mass attributes are read from.

update_KC0v_onlyint

The default 0 means that the row and column indices KC0r and KC0c should also be updated. Any other value will only update the stiffness matrix values KC0v.

update_KG(self, long[::1] KGr, long[::1] KGc, double[::1] KGv, BeamProp prop, int update_KGv_only=0) void#

Update sparse vectors for geometric stiffness matrix KG

Parameters:
KGrnp.array

Array to store row positions of sparse values

KGcnp.array

Array to store column positions of sparse values

KGvnp.array

Array to store sparse values

propBeamProp object

Beam property object from where the stiffness and mass attributes are read from.

update_KGv_onlyint

The default \(0\) means that only \(KGv\) is updated. Any other value will lead to \(KGr\) and \(KGc\) also being updated.

update_M(self, long[::1] Mr, long[::1] Mc, double[::1] Mv, BeamProp prop, int mtype=0) void#

Update sparse vectors for mass matrix M

Parameters:
Mrnp.array

Array to store row positions of sparse values

Mcnp.array

Array to store column positions of sparse values

Mvnp.array

Array to store sparse values

mtypeint, optional

0 for consistent mass matrix using method from Brockman 1987 1 for reduced integration mass matrix using method from Brockman 1987 2 for lumped mass matrix using method from Brockman 1987

update_length(self) void#

Update element length

update_probe_ue(self, double[::1] u) void#

Update the local displacement vector of the probe of the element

Note

The probe attribute object BeamLRProbe is updated, not the element object.

Parameters:
uarray-like

Array with global displacements, for a total of \(M\) nodes in the model, this array will be arranged as: \(u_1, v_1, w_1, {r_x}_1, {r_y}_1, {r_z}_1, u_2, v_2, w_2, {r_x}_2, {r_y}_2, {r_z}_2, ..., u_M, v_M, w_M, {r_x}_M, {r_y}_M, {r_z}_M\).

update_probe_xe(self, double[::1] x) void#

Update the 3D coordinates of the probe of the element

Note

The probe attribute object BeamLRProbe is updated, not the element object.

Parameters:
xarray-like

Array with global nodal coordinates, for a total of \(M\) nodes in the model, this array will be arranged as: \(x_1, y_1, z_1, x_2, y_2, z_2, ..., x_M, y_M, z_M\).

update_rotation_matrix(self, double vxyi, double vxyj, double vxyk, double[::1] x) void#

Update the rotation matrix of the element

Attributes r11,r12,r13,r21,r22,r23,r31,r32,r33 are updated, corresponding to the rotation matrix from local to global coordinates.

The element attributes vxyi, vxyj and vxyk are also updated when this function is called.

The element coordinate system is determined, identifying the \(ijk\) components of each axis: \({x_e}_i, {x_e}_j, {x_e}_k\); \({y_e}_i, {y_e}_j, {y_e}_k\); \({z_e}_i, {z_e}_j, {z_e}_k\).

Parameters:
vxyi, vxyj, vxykdouble

Components of a vector on the \(XY\) plane of the element coordinate system, defined using global coordinates.

xarray-like

Array with global nodal coordinates, for a total of \(M\) nodes in the model, this array will be arranged as: \(x_1, y_1, z_1, x_2, y_2, z_2, ..., x_M, y_M, z_M\).

vxyi#

vxyi: ‘double’

vxyj#

vxyj: ‘double’

vxyk#

vxyk: ‘double’

class pyfe3d.beamlr.BeamLRData#

Used to allocate memory for the sparse matrices.

Attributes:
KC0_SPARSE_SIZE,int

KC0_SPARSE_SIZE = 144

KG_SPARSE_SIZE,int

KG_SPARSE_SIZE = 36

M_SPARSE_SIZE,int

M_SPARSE_SIZE = 144

KC0_SPARSE_SIZE#

KC0_SPARSE_SIZE: ‘int’

KG_SPARSE_SIZE#

KG_SPARSE_SIZE: ‘int’

M_SPARSE_SIZE#

M_SPARSE_SIZE: ‘int’

class pyfe3d.beamlr.BeamLRProbe#

Probe used for local coordinates, local displacements, local stresses etc

Attributes:
xe,array-like

Array of size NUM_NODES*DOF//2=6 containing the nodal coordinates in the element coordinate system, in the following order \({x_e}_1, {y_e}_1, {z_e}_1, {x_e}_2, {y_e}_2, {z_e}_2\).

ue,array-like

Array of size NUM_NODES*DOF=12 containing the element displacements in the following order \({u_e}_1, {v_e}_1, {w_e}_1, {{r_x}_e}_1, {{r_y}_e}_1, {{r_z}_e}_1, {u_e}_2, {v_e}_2, {w_e}_2, {{r_x}_e}_2, {{r_y}_e}_2, {{r_z}_e}_2\).

ue#

ue: ‘double[::1]’

xe#

xe: ‘double[::1]’