To represent conventional, bending moment only sections numerically, only geometry information, such as location and area, is required. This allows the usage of predefined basic shapes, such as rectangles, circles. And complex sections can be built up by using those basic shapes. However, such an approach is not practical for sections involving torsion and warping, as some important sectional properties can only be obtained by analysing the section as a whole. Hence, a general approach that only focuses on section discretisation is used. The discretisation refers to the division of the section into a number of fibres/cells. Those fibres are the smallest units that are used in numerical integration. The users shall provide the corresponding sectional properties.
To obtain necessary sectional properties, a section analysis shall be performed. For arbitrary sections, often a general purpose FEA based method shall be used. This book provides a thorough discussion on relevant numerical methods.
The sectionproperties library implements such a method. It can be used to prepare data to be used in suanPan.
Theory
We first discuss what sectional properties are required and provide a theoretical background.
section Cell3DOS (1) (2) (3) (4) (5) (6) [7] [8]
# (1) int, unique section tag
# (2) double, area
# (3) double, sectional coordiante
# (4) double, py
# (5) double, pz
# (6) int, material tag
# [7] double, eccentricity/location along y axis, default: 0.0
# [8] double, eccentricity/location along z axis, default: 0.0
Apart from the location (coordinates/eccentricities) and area, it requires sectional coordinate (warping function) (3) and its derivatives (4) and (5).
The warping function is the conventional one that can be expressed as a function of locations.
ω=ω(y,z).
The integration shall result in a net-zero warping.
∫Aω(y,z)dA=0.
Its derivatives along two axes are also required.
∂y∂ω,∂z∂ω.
Alemdar's thesis, and some other derivatives, uses γ=2nϕ′ to approximate shear strain. This is only suitable for thin-walled open sections. Here we use general formulation instead for complex, universal, arbitrary sections. Interested readers can refer to, for example, Pilkey's book, for more details.
Two shear strains can be expressed as follows.
γxy=ϕ′(∂y∂ω−z),γxz=ϕ′(∂z∂ω+y).
A Simple Example
We start with a simple circular example. Let's assume it is centred and has a diameter of 10.
We use built-in functionality to create such a geometry and assign a mesh size of 0.1.
import matplotlib.pyplot as pltimport sectionproperties.pre.library.primitive_sections as sectionsfrom sectionproperties.analysis.section import Section# approximate a circle by a polygongeometry = sections.circular_section(d=10, n=64)geometry.create_mesh(mesh_sizes=[0.1])section =Section(geometry)section.plot_mesh(title="Circular Section Mesh")
<Axes: title={'center': 'Circular Section Mesh'}>
The mesh is generated internally using second-order triangular elements (CPS6 in ABAQUS notation). One can check the details of the created elements. The following shows the first element.
Alternatively, the area can also be computed as the determinant of the Jacobian matrix of the isoparametric mapping.
Locations
We assume the triangle element is small enough so that the centroid can represent the location of the element. The centroid can be obtained by averaging the coordinates of the first three nodes.
As we are using a circular section as the example, the warping function is zero everywhere. The above method computes the warping in the local coordinate system. If one wants to use the centroid as the origin, one shall translate the coordinates first.
Partial Derivatives
To obtain ∂y∂ω and ∂z∂ω, we use the shape functions.
from sectionproperties.analysis.fea import shape_function# use centroidN, B, _, _, _ =shape_function(ele.coords, [0, 1.0/3.0, 1.0/3.0, 1.0/3.0])centre_omega = omega[ele.node_ids].dot(N)centre_partial = B @ omega[ele.node_ids]print(f"warping function = {centre_omega}")print(f"derivatives = {centre_partial}")
warping function = 1.1776217629584232e-12
derivatives = [-2.11447839e-12 -1.18613467e-12]
Since for circular sections, there is no warping, the corresponding derivatives are trivial.
A Universal Function
It appears to be easy to obtain necessary properties for each element. Here we provide a universal function to generate the required data.
One shall adjust it accordingly for specific purposes.
from numpy import ndarray, arrayfrom sectionproperties.analysis.fea import Tri6classCell3DOS: tag:int area:float y:float z:float omega:float py:float pz:floatdef__init__(self,ele: Tri6,omega: ndarray =None):""" If `omega` is None, no warping is considered. """ N, B, self.area, _, _ =shape_function(ele.coords, array([0, 1.0/3.0, 1.0/3.0, 1.0/3.0])) self.tag = ele.el_id +1# el_id is zero-based self.y, self.z = ele.coords[:,:3].mean(axis=1) self.omega = omega.dot(N)if omega isnotNoneelse0 self.py, self.pz = (B @ omega) if omega isnotNoneelse (0,0)defexport(self):return (f"section Cell3DOS"f" {self.tag}"f" {self.area:+E}"f" {self.omega:+E}"f" {self.py:+E}"f" {self.pz:+E}"# f" $mat_tag"f" 1"# here we just user 1 for brevityf" {self.y:+E}"f" {self.z:+E}" )defto_cell3dos(geometry,*, no_warping:bool=False): section =Section(geometry)# if shift is required, one can, for example, do the following# section.calculate_geometric_properties()# for ele in section.elements:# ele.coords[0, :] -= section.section_props.cx# ele.coords[1, :] -= section.section_props.cyif no_warping:return [Cell3DOS(ele)for ele in section.elements] stiffness, force = section.assemble_torsion() omega = solver.solve_direct_lagrange(stiffness, force)return [Cell3DOS(ele, omega[ele.node_ids])for ele in section.elements]defto_file(cells,file_name):withopen(file_name, "w")as f:for cell in cells: f.write(cell.export() +"\n") f.write(f"section Fibre3DOS {1+len(cells)} ") f.write(" ".join([str(cell.tag) for cell in cells]) +"\n")
Now it is possible to pass in a geometry, and call the export function to generate cell data.
to_file(to_cell3dos(geometry), "circular.sp")
!!! Note The above functionality has been added to the latest sectionproperties library, see this page.
Validation
Now it is possible to create a simple single element model to validate the data. Remember to change the material tag accordingly.
node 1 0 0 0
node 2 1 0 0
material ElasticOS 1 2. 0. ! <--- so that G=1
file circular.sp
orientation B3DOSL 1 0. 0. 1.
element B31OS 1 1 2 1203 1 6 ! <--- 1203 is the tag for Fibre3DOS section
fix2 1 E 1
displacement 1 0 1E-3 4 2
step static 1
set ini_step_size 2E-1
set fixed_step_size true
analyze
peek node 2
exit
We shall see the initial stiffness is 74.93, which is close to the torsional constant J. Due to the existence of bi-moments, the torque is not linearly proportional to the rotation. If one records the S output, the St. Venant torsion would be GJϕ′.
Flat Bar Example
We show another example with analytical solution. This example is taken from this paper.
# flat bar
node 1 0 0 0
node 2 1000 0 0
material ElasticOS 1 200 .25
file flat.sp
orientation B3DOSL 1 0. 0. 1.
element B31OS 1 1 2 1574 1 6
fix2 1 E 1
displacement 1 0 1.4 4 2
hdf5recorder 1 Node RF4 2
step static 1
set ini_step_size 1E-2
set fixed_step_size true
converger RelIncreDisp 1 1E-10 5 1
analyze
save recorder 1
exit
Use the following script to plot the results.
with h5py.File('R1-RF4.h5', 'r')as f: data = f['/R1-RF4/R1-RF42'] x = data[:,0]*1.4 plt.plot(x, data[:, 1] /1000, label='B31OS') ref = [z *.080*66.667+.100*17.778* z **3for z in x] plt.plot(x, ref, label='theoretical') plt.legend() plt.xlabel('rotation (rad)') plt.ylabel('torque (kNm)') plt.tight_layout(pad=.2) plt.show()
The discretisation introduces some errors. With mesh refinement, the numerical result approaches the analytical solution.
Closing Remarks
Cell3DOS is the basic building block. Fibre3DOS collects all cells and defines the section. B31OS is the corresponding element.
Sectional integration is always about the origin of the local coordinate system. Shifting of axis can be accounted for by directly defining the section in the shifted position.