Quad4 - Quadrilateral element with mixed integration (pyfe3d.quad4)#

The Quad4 is the recommended quadrilateral plane stress finite element.

Another option is the pyfe3d.quad4r with full reduced integration, more efficient, but with an hourglass control to compensate the reduced integration that is not robust and creates significant artificial stiffness.

The Quad4 element has 6 degrees-of-freedom (DOF): \(u\), \(v\), \(w\), \(r_x\), \(r_y\), \(r_z\). All DOF are interpolated bi-linearly between the nodes, such that any of the DOF gradients can be constant over the element when the element is rectangular.

The stiffness for the degrees of freedom \(w\), \(r_x\) and \(r_y\) is based on the paper below, where \(r_x = \theta_1\) and \(r_y = \theta_2\):

Hughes T.J.R., Taylor R.L., Kanoknukulchai W. “A simple and efficient finite element for plate bending”. International Journal of Numerical Methods in Engineering, Volume 11, 1977. https://doi.org/10.1002/nme.1620111005

Hughes et al. (1977) proposed the following integration scheme:

  • For thin plates, when \(h/\ell < 1\), where \(\ell\) is the element characteristic length, here calculated as the square root of the element area \(\ell = \sqrt{\text{area}}\)

– two-by-two quadrature for the bending energy terms

– one-point quadrature for the transverse shear energy terms

  • For thick plates, when \(h/\ell >= 1\)

– two-by-two quadrature for the bending energy terms

– two-by-two quadrature for the transverse shear terms with gradients

– one-point quadrature for the transverse shear terms without gradients

The in-plane stiffness terms are integrated with 2 quadrature points, and the drilling stiffness is integrated with 1 quadrature point. These are not specified in the paper of Hughes et al. (1977).

class pyfe3d.quad4.Quad4#

Nodal connectivity for the plate element similar to Nastran’s CQUAD4:

 ^ y
 |
4 ________ 3
 |       |
 |       |   --> x
 |       |
 |_______|
1         2

The element coordinate system is determined identically what is explained in Nastran’s quick reference guide for the CQUAD4 element, as illustrated below.

_images/nastran_cquad4.svg
Attributes:
eid,int

Element identification number.

area,double

Element area.

K6ROT,double

Element drilling stiffness, added only to the diagonal of the local stiffness matrix. The default value is according to AUTODESK NASTRAN’s quick reference guide is K6ROT = 100. for static analysis. For modal solutions, this value should be K6ROT=1.e4. MSC NASTRAN’s quick reference guide states that K6ROT > 100. should not be used, but this is controversial, already being controversial to what AUTODESK NASTRAN’s manual says.

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

Rotation matrix from local to global coordinates.

m11, m12, m21, m22double

Rotation matrix only for the constitutive relations. Used when a material direction is used instead of the element local coordinates.

c1, c2, c3, c4int

Position of each node in the global stiffness matrix.

n1, n2, n3, n4int

Node identification number.

init_k_KC0, init_k_KG, init_k_Mint

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

init_k_KA_beta, init_k_KA_gamma, init_k_CAint

Position in the arrays storing the sparse data for the aerodynamic matrices based on the Piston theory.

probe,Quad4Probe object

Pointer to the probe.

Methods

update_CA(self, long[, long[, double[)

Update sparse vectors for piston-theory aerodynamic damping matrix \(CA\)

update_KA_beta(self, long[, long[, double[)

Update sparse vectors for piston-theory aerodynamic matrix \(KA_{\beta}\)

update_KA_gamma(self, long[, long[, double[)

Update sparse vectors for piston-theory aerodynamic matrix \(KA_{\gamma}\)

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_KG_given_stress(self, double Nxx, ...)

Update sparse vectors for geometric stiffness matrix KG

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

Update sparse vectors for mass matrix M

update_area(self)

Update element area

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[, ...)

Update the rotation matrix of the element

K6ROT#

K6ROT: ‘double’

area#

area: ‘double’

c1#

c1: ‘int’

c2#

c2: ‘int’

c3#

c3: ‘int’

c4#

c4: ‘int’

eid#

eid: ‘int’

init_k_CA#

init_k_CA: ‘int’

init_k_KA_beta#

init_k_KA_beta: ‘int’

init_k_KA_gamma#

init_k_KA_gamma: ‘int’

init_k_KC0#

init_k_KC0: ‘int’

init_k_KG#

init_k_KG: ‘int’

init_k_M#

init_k_M: ‘int’

m11#

m11: ‘double’

m12#

m12: ‘double’

m21#

m21: ‘double’

m22#

m22: ‘double’

n1#

n1: ‘int’

n2#

n2: ‘int’

n3#

n3: ‘int’

n4#

n4: ‘int’

probe#

probe: pyfe3d.quad4.Quad4Probe

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_CA(self, long[::1] CAr, long[::1] CAc, double[::1] CAv) void#

Update sparse vectors for piston-theory aerodynamic damping matrix \(CA\)

Parameters:
CArnp.array

Array to store row positions of sparse values

CAcnp.array

Array to store column positions of sparse values

CAvnp.array

Array to store sparse values

update_KA_beta(self, long[::1] KA_betar, long[::1] KA_betac, double[::1] KA_betav) void#

Update sparse vectors for piston-theory aerodynamic matrix \(KA_{\beta}\)

Parameters:
KA_betarnp.array

Array to store row positions of sparse values

KA_betacnp.array

Array to store column positions of sparse values

KA_betavnp.array

Array to store sparse values

update_KA_gamma(self, long[::1] KA_gammar, long[::1] KA_gammac, double[::1] KA_gammav) void#

Update sparse vectors for piston-theory aerodynamic matrix \(KA_{\gamma}\)

Parameters:
KA_gammarnp.array

Array to store row positions of sparse values

KA_gammacnp.array

Array to store column positions of sparse values

KA_gammavnp.array

Array to store sparse values

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

Update sparse vectors for linear constitutive stiffness matrix KC0

Drilling stiffness is used according to Adam et al. 2013:

Adam, F. M., Mohamed, A. E., and Hassaballa, A. E., 2013, “Degenerated Four Nodes Shell Element with Drilling Degree of Freedom,” IOSR J. Eng., 3(8), pp. 10–20.

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

propShellProp object

Shell 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, ShellProp prop) void#

Update sparse vectors for geometric stiffness matrix KG

Two-point Gauss-Legendre quadrature is used, which showed more accuracy for linear buckling load predictions.

Before this function is called, the probe Quad4Probe attribute of the Quad4 object must be updated using update_probe_ue() with the correct pre-buckling (fundamental state) displacements; and update_probe_xe() with the node coordinates.

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

propShellProp object

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

update_KG_given_stress(self, double Nxx, double Nyy, double Nxy, long[::1] KGr, long[::1] KGc, double[::1] KGv) void#

Update sparse vectors for geometric stiffness matrix KG

Note

A constant stress state is assumed within the element, according to the given values of \(N_{xx}, N_{yy}, N_{xy}\).

Two-point Gauss-Legendre quadrature is used, which showed more accuracy for linear buckling load predictions.

Before this function is called, the probe Quad4Probe attribute of the Quad4 object must be updated using update_probe_xe() with the node coordinates.

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

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

Update sparse vectors for mass matrix M

Different integration schemes are available by means of the mtype parameter.

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_area(self) void#

Update element area

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

Update the local displacement vector of the probe of the element

Note

The probe attribute object Quad4Probe 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 Quad4Probe 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[::1] x, double xmati=0., double xmatj=0., double xmatk=0.) 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 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\).

The rotation matrix terms are calculated after solving 9 equations.

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\).

xmati, xmatj, xmatk: array-like

Vector in global coordinates representing the material direction. This vector is projected onto the plate element, thus becoming the material direction. The \(ABD\) matrix defining the constitutive behavior of the element is rotated from the material direction to the element \(x\) axis while calculating the stiffness matrices.

class pyfe3d.quad4.Quad4Data#

Used to allocate memory for the sparse matrices.

Attributes:
KC0_SPARSE_SIZE,int

KC0_SPARSE_SIZE = 576

KG_SPARSE_SIZE,int

KG_SPARSE_SIZE = 144

M_SPARSE_SIZE,int

M_SPARSE_SIZE = 480

KA_BETA_SPARSE_SIZE,int

KA_BETA_SPARSE_SIZE = 144

KA_GAMMA_SPARSE_SIZE,int

KA_GAMMA_SPARSE_SIZE = 144

CA_SPARSE_SIZE,int

CA_SPARSE_SIZE = 144

CA_SPARSE_SIZE#

CA_SPARSE_SIZE: ‘int’

KA_BETA_SPARSE_SIZE#

KA_BETA_SPARSE_SIZE: ‘int’

KA_GAMMA_SPARSE_SIZE#

KA_GAMMA_SPARSE_SIZE: ‘int’

KC0_SPARSE_SIZE#

KC0_SPARSE_SIZE: ‘int’

KG_SPARSE_SIZE#

KG_SPARSE_SIZE: ‘int’

M_SPARSE_SIZE#

M_SPARSE_SIZE: ‘int’

class pyfe3d.quad4.Quad4Probe#

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

The idea behind using a probe is to avoid allocating larger memory buffers per finite element. The memory buffers are allocated per probe, and one probe can be shared amongst many finite elements, with the information being updated and retrieved on demand.

Attributes:
xe,array-like

Array of size NUM_NODES*DOF//2=12 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\), \({x_e}_3, {y_e}_3, {z_e}_3\), \({x_e}_4, {y_e}_4, {z_e}_4\).

ue,array-like

Array of size NUM_NODES*DOF=24 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\), \({u_e}_3, {v_e}_3, {w_e}_3, {{r_x}_e}_3, {{r_y}_e}_3, {{r_z}_e}_3\), \({u_e}_4, {v_e}_4, {w_e}_4, {{r_x}_e}_4, {{r_y}_e}_4, {{r_z}_e}_4\).

KC0ve,array-like

Local stiffness matrix stored as a 1D array of size (NUM_NODES*DOF)**2.

BLexx, BLeyy, BLgxyarray-like

Arrays of size NUM_NODES*DOF=24 containing the in-plane strain interpolation functions evaluated at a given natural coordinate point \(\xi\), \(\eta\).

BLkxx, BLkyy, BLkxyarray-like

Arrays of size NUM_NODES*DOF=24 containing the bending strain interpolation functions evaluated at a given natural coordinate point \(\xi\), \(\eta\).

BLgyz_grad, BLgyz_rot, BLgxz_grad, BLgxz_rotarray-like

Arrays of size NUM_NODES*DOF=24 containing the transverse shear strain interpolation functions evaluated at a given natural coordinate point \(\xi\), \(\eta\).

Methods

update_BL(self, double xi, double eta)

Update all components of the interpolation matrix \(\pmb{B_L}\) at a given natural coordinate point \(\xi\), \(\eta\).

BLexx#

BLexx: ‘double[::1]’

BLeyy#

BLeyy: ‘double[::1]’

BLgxy#

BLgxy: ‘double[::1]’

BLgxz_grad#

BLgxz_grad: ‘double[::1]’

BLgxz_rot#

BLgxz_rot: ‘double[::1]’

BLgyz_grad#

BLgyz_grad: ‘double[::1]’

BLgyz_rot#

BLgyz_rot: ‘double[::1]’

BLkxx#

BLkxx: ‘double[::1]’

BLkxy#

BLkxy: ‘double[::1]’

BLkyy#

BLkyy: ‘double[::1]’

KC0ve#

KC0ve: ‘double[::1]’

ue#

ue: ‘double[::1]’

update_BL(self, double xi, double eta) void#

Update all components of the interpolation matrix \(\pmb{B_L}\) at a given natural coordinate point \(\xi\), \(\eta\).

xe#

xe: ‘double[::1]’