Code: Set up some function definitions
# %% Load packages
import numpy as np
import matplotlib.pyplot as plt
from scipy.linalg import block_diag
def beamFE(dim, matProp, freq, numE=100, Method='Direct', pos='', BC=['f', 'f']):
"""
beamFE computes the mobility matrix of a free-free beam.
Parameters
----------
dim : 3x1 numpy array
Dimensions of plate: length x, length y, thickness - [lx, ly, h].
matProp : 3x1 numpy array
Material properties: density, Young's modulus, loss factor -
[rho, E, nu].
freq : Mx1 numpy array
Frequency vector - [f1, f2, ... fM].
numE : Int,
Number of finite elements used in model. The default is 100.
Method : Str,
Method used to calculate mobility matrix: 'Direct' matrix inversion or
'Modal' summation
pos : Nx1 numpy array
Positions of interest when using 'Modal' summation method (by default
all nodal positions are retained)
BC : 2x1 list
Define boundary conditions for two ends (free, simply supported,
or clamped), e.g. ['f', 'f'] ['s', 's'] ['c', 'c'] ['s','c'] etc.
Returns
-------
Y : NxNxM numpy array
Complex mobility matrix for finite element beam model.
"""
if Method == 'Direct' and type(pos) is np.ndarray:
print(" Warning: to use optional argument 'pos' the 'Method' option \
has been changed to 'Modal'. ")
Method = 'Modal'
# Complex Young's modulus (structural damping)
E = matProp[1] * (1+1j*matProp[2])
# a = dim[0]/(2*numE)
L = dim[0]/(numE)
A = dim[2]*dim[1] # Cross-sectional area
I_ = dim[1] * dim[2]**3/12 # Moment of inertia
# Allocate space
M = np.zeros((2*numE+2, 2*numE+2))
K = np.zeros((2*numE+2, 2*numE+2), dtype=complex)
D = np.zeros((2*numE+2, 2*numE+2, freq.shape[0]), dtype=complex)
Y = np.zeros((2*numE+2, 2*numE+2, freq.shape[0]), dtype=complex)
# Element mass matrix
m_e = matProp[0]*A*L/420 * np.array([[156, 22*L, 54, -13*L],
[22*L, 4*L**2, 13*L, -3*L**2],
[54, 13*L, 156, -22*L],
[-13*L, -3*L**2, -22*L, 4*L**2]])
# Element stiffness matrix
k_e = E*I_/(L**3) * np.array([[12, 6*L, -12, 6*L],
[6*L, 4*L**2, -6*L, 2*L**2],
[-12, -6*L, 12, -6*L],
[6*L, 2*L**2, -6*L, 4*L**2]])
# print(k_e)
for elm in range(0, numE*2-1, 2):
M[elm:elm+4, elm:elm+4] = M[elm:elm+4, elm:elm+4] + m_e
K[elm:elm+4, elm:elm+4] = K[elm:elm+4, elm:elm+4] + k_e
if BC[1] == 's': # Simply supported BC on right - delete translational DoF
M = np.delete(M, -2, 0)
M = np.delete(M, -2, 1)
K = np.delete(K, -2, 0)
K = np.delete(K, -2, 1)
Y = np.delete(Y, -2, 0)
Y = np.delete(Y, -2, 1)
D = np.delete(D, -2, 0)
D = np.delete(D, -2, 1)
# Clamped BC on right - delete translational and rotational DoF
if BC[1] == 'c':
M = np.delete(M, (-2, -1), 0)
M = np.delete(M, (-2, -1), 1)
K = np.delete(K, (-2, -1), 0)
K = np.delete(K, (-2, -1), 1)
Y = np.delete(Y, (-2, -1), 0)
Y = np.delete(Y, (-2, -1), 1)
D = np.delete(D, (-2, -1), 0)
D = np.delete(D, (-2, -1), 1)
# Simply supported BC on left - delete translational DoF
if BC[0] == 's':
M = np.delete(M, 0, 0)
M = np.delete(M, 0, 1)
K = np.delete(K, 0, 0)
K = np.delete(K, 0, 1)
Y = np.delete(Y, 0, 0)
Y = np.delete(Y, 0, 1)
D = np.delete(D, 0, 0)
D = np.delete(D, 0, 1)
# Clamped supported BC on left - delete translational and rotational DoF
if BC[0] == 'c':
M = np.delete(M, (1, 0), 0)
M = np.delete(M, (1, 0), 1)
K = np.delete(K, (1, 0), 0)
K = np.delete(K, (1, 0), 1)
Y = np.delete(Y, (1, 0), 0)
Y = np.delete(Y, (1, 0), 1)
D = np.delete(D, (1, 0), 0)
D = np.delete(D, (1, 0), 1)
if Method == 'Direct':
for fi in range(freq.shape[0]):
# Dynamic stiffness matrix
D[:, :, fi] = K - (2*np.pi*freq[fi])**2*M
# Mobility matrix by direct inversion
Y[:, :, fi] = 1j*(2*np.pi*freq[fi])*np.linalg.pinv(D[:, :, fi])
return Y, M, K
elif Method == 'Modal':
# Get natural freqs (squared) and mode shapes
W2, V = sci.linalg.eigh(np.real(K), M)
# Mass normalise mode shapes
V = V @ np.diag(1/np.sqrt(np.diag(V.T @ (M) @ V)))
# Use linear interpolation to get mode shape at arbitary positions
if type(pos) is np.ndarray:
# Interpolation for translation
ft = sci.interpolate.interp1d(np.linspace(0, dim[0], numE+1),
V[0:-1:2], kind='linear', axis=0)
# Interpolation for rotation
fr = sci.interpolate.interp1d(np.linspace(0, dim[0], numE+1),
V[1::2], kind='linear', axis=0)
Vt = ft(pos) # Translational mode shape and defined position
Vr = fr(pos) # Rotational mode shape and defined position
V = np.empty((Vt.shape[0] + Vr.shape[0], Vt.shape[1]),
dtype=Vt.dtype) # Combined mode shape matrix
V[0::2, :] = Vt
V[1::2, :] = Vr
Y = np.zeros((2*pos.shape[0], 2*pos.shape[0], freq.shape[0]),
dtype=complex) # Reallocate space
for fi in range(freq.shape[0]):
wi = 2*np.pi*freq[fi] # Radian frequency
# Modal summation with structural damping
Y[:, :, fi] = 1j*wi * V[:, :] @ np.diag(1/(W2[:]*(1+1j*matProp[2])
- wi**2)) @ V[:, :].T
return Y, M, K
def subStruct(intFlex=np.zeros(1,), **inputs):
"""
Compute coupled admittance using primal or dual substructuring formulation
Primal : Yc = (L^T Y^-1 L)^-1
Dual : Yc = Y - Y B^T (B Y B^T + G)^-1 B Y
Parameters
----------
freq : frequency vector
Y : block diagonal admittance matrix
BL : primal or dual Boolean coupling matrix.
Type : primal or dual specificiation
intFlex : interface flexibility matrix (undefined = rigid coupling)
Returns
-------
Yc : coupled admittance.
"""
f = inputs["freq"]
Y = inputs["Y"]
BL = inputs["BL"]
Type = inputs["Type"]
G = intFlex
if type(f) is np.ndarray: # If frequency array
numF = f.shape[0]
if Type == 'Primal': # Allocate space
Yc = np.zeros((BL.shape[1], BL.shape[1], numF), dtype=complex)
elif Type == 'Dual':
Yc = np.zeros((Y.shape[0], Y.shape[1], numF), dtype=complex)
for fi in range(numF): # Cycle through frequency array
if Type == 'Primal':
# Primal SS:
Yc[:, :, fi] = np.linalg.inv(BL.T @
np.linalg.inv(Y[:, :, fi])@BL)
elif Type == 'Dual':
if G.shape[0] == 1: # If G undefined (rigid)
Yc[:, :, fi] = Y[:, :, fi]-Y[:, :, fi]@BL.T @\
np.linalg.pinv(BL@Y[:, :, fi]@BL.T)@BL@Y[:, :, fi]
else: # If G defined (resilient)
Yc[:,:,fi] = Y[:,:,fi]-Y[:,:,fi]@BL.T@np.linalg.pinv(BL@Y[:,:,fi]@BL.T+G[:,:,fi])@BL@Y[:,:,fi]
else: # If scalar frequency value
if Type =='Primal': # Allocate space
Yc = np.zeros((BL.shape[1],BL.shape[1]),dtype=complex)
Yc = np.linalg.inv(BL.T@np.linalg.inv(Y)@BL) # Primal SS
elif Type == 'Dual':
Yc = np.zeros((Y.shape[0],Y.shape[1]),dtype=complex)
if G.shape[0] == 1: # If G undefined (rigid)
Yc = Y-Y@BL.T@np.linalg.pinv(BL@Y@BL.T)@BL@Y
else: # If G defined (resilient)
Yc = Y-Y@BL.T@np.linalg.pinv(BL@Y@BL.T+G)@BL@Y
return Yc