Python API Documentation

Pre-Processor Package

sections Module

Geometry Class

class sectionproperties.pre.sections.Geometry(control_points, shift)[source]

Parent class for a cross-section geometry input.

Provides an interface for the user to specify the geometry defining a cross-section. A method is provided for generating a triangular mesh, for translating the cross-section by (x, y) and for plotting the geometry.

Variables:
  • points (list[list[float, float]]) – List of points (x, y) defining the vertices of the cross-section
  • facets (list[list[int, int]]) – List of point index pairs (p1, p2) defining the edges of the cross-section
  • holes (list[list[float, float]]) – List of points (x, y) defining the locations of holes within the cross-section. If there are no holes, provide an empty list [].
  • control_points (list[list[float, float]]) – A list of points (x, y) that define different regions of the cross-section. A control point is an arbitrary point within a region enclosed by facets.
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)
  • perimeter (list[int]) – List of facet indices defining the perimeter of the cross-section
add_control_point(control_point)[source]

Adds a control point to the geometry and returns the added control point id.

Parameters:hole (list[float, float]) – Location of the control point
Returns:Control point id
Return type:int
add_facet(facet)[source]

Adds a facet to the geometry and returns the added facet id.

Parameters:facet (list[float, float]) – Point indices of the facet
Returns:Facet id
Return type:int
add_hole(hole)[source]

Adds a hole location to the geometry and returns the added hole id.

Parameters:hole (list[float, float]) – Location of the hole
Returns:Hole id
Return type:int
add_point(point)[source]

Adds a point to the geometry and returns the added point id.

Parameters:point (list[float, float]) – Location of the point
Returns:Point id
Return type:int
calculate_extents()[source]

Calculates the minimum and maximum x and y-values amongst the list of points.

Returns:Minimum and maximum x and y-values (x_min, x_max, y_min, y_max)
Return type:tuple(float, float, float, float)
calculate_facet_length(facet)[source]

Calculates the length of the facet.

Parameters:facet – Point index pair (p1, p2) defining a facet
Returns:Facet length
Return type:float
calculate_perimeter()[source]

Calculates the perimeter of the cross-section by summing the length of all facets in the perimeter class variable.

Returns:Cross-section perimeter, returns 0 if there is no perimeter defined
Return type:float
clean_geometry(verbose=False)[source]

Peforms a full clean on the geometry.

Parameters:verbose (bool) – If set to true, information related to the geometry cleaning process is printed to the terminal.

Note

Cleaning the geometry is always recommended when creating a merged section, which may result in overlapping or intersecting facets, or duplicate nodes.

create_mesh(mesh_sizes)[source]

Creates a quadratic triangular mesh from the Geometry object.

Parameters:mesh_sizes – A list of maximum element areas corresponding to each region within the cross-section geometry.
Returns:Object containing generated mesh data
Return type:meshpy.triangle.MeshInfo
Raises:AssertionError – If the number of mesh sizes does not match the number of regions

The following example creates a circular cross-section with a diameter of 50 with 64 points, and generates a mesh with a maximum triangular area of 2.5:

import sectionproperties.pre.sections as sections

geometry = sections.CircularSection(d=50, n=64)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
../_images/circle_mesh.png

Mesh generated from the above geometry.

draw_radius(pt, r, theta, n, anti=True)[source]

Adds a quarter radius of points to the points list - centered at point pt, with radius r, starting at angle theta, with n points. If r = 0, adds pt only.

Parameters:
  • pt (list[float, float]) – Centre of radius (x,y)
  • r (float) – Radius
  • theta (float) – Initial angle
  • n (int) – Number of points
  • anti (bool) – Anticlockwise rotation?
mirror_section(axis='x', mirror_point=None)[source]

Mirrors the geometry about a point on either the x or y-axis. If no point is provided, mirrors the geometry about the first control point in the list of control points of the Geometry object.

Parameters:
  • axis (string) – Axis about which to mirror the geometry, ‘x’ or ‘y’
  • mirror_point (list[float, float]) – Point about which to mirror the geometry (x, y)

The following example mirrors a 200PFC section about the y-axis and the point (0, 0):

import sectionproperties.pre.sections as sections

geometry = sections.PfcSection(d=200, b=75, t_f=12, t_w=6, r=12, n_r=8)
geometry.mirror_section(axis='y', mirror_point=[0, 0])
plot_geometry(ax=None, pause=True, labels=False, perimeter=False)[source]

Plots the geometry defined by the input section. If no axes object is supplied a new figure and axis is created.

Parameters:
  • ax (matplotlib.axes.Axes) – Axes object on which the mesh is plotted
  • pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
  • labels (bool) – If set to true, node and facet labels are displayed
  • perimeter (bool) – If set to true, boldens the perimeter of the cross-section
Returns:

Matplotlib figure and axes objects (fig, ax)

Return type:

(matplotlib.figure.Figure, matplotlib.axes)

The following example creates a CHS discretised with 64 points, with a diameter of 48 and thickness of 3.2, and plots the geometry:

import sectionproperties.pre.sections as sections

geometry = sections.Chs(d=48, t=3.2, n=64)
geometry.plot_geometry()
../_images/chs_geometry.png

Geometry generated by the above example.

rotate_section(angle, rot_point=None)[source]

Rotates the geometry and specified angle about a point. If the rotation point is not provided, rotates the section about the first control point in the list of control points of the Geometry object.

Parameters:
  • angle (float) – Angle (degrees) by which to rotate the section. A positive angle leads to a counter-clockwise rotation.
  • rot_point (list[float, float]) – Point (x, y) about which to rotate the section

The following example rotates a 200UB25 section clockwise by 30 degrees:

import sectionproperties.pre.sections as sections

geometry = sections.ISection(d=203, b=133, t_f=7.8, t_w=5.8, r=8.9, n_r=8)
geometry.rotate_section(angle=-30)
shift_section()[source]

Shifts the cross-section parameters by the class variable vector shift.

CustomSection Class

class sectionproperties.pre.sections.CustomSection(points, facets, holes, control_points, shift=[0, 0], perimeter=[])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a cross-section from a list of points, facets, holes and a user specified control point.

Parameters:
  • points (list[list[float, float]]) – List of points (x, y) defining the vertices of the cross-section
  • facets (list[list[int, int]]) – List of point index pairs (p1, p2) defining the edges of the cross-section
  • holes (list[list[float, float]]) – List of points (x, y) defining the locations of holes within the cross-section. If there are no holes, provide an empty list [].
  • control_points (list[list[float, float]]) – A list of points (x, y) that define different regions of the cross-section. A control point is an arbitrary point within a region enclosed by facets.
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)
  • perimeter – List of facet indices defining the perimeter of the cross-section

The following example creates a hollow trapezium with a base width of 100, top width of 50, height of 50 and a wall thickness of 10. A mesh is generated with a maximum triangular area of 2.0:

import sectionproperties.pre.sections as sections

points = [[0, 0], [100, 0], [75, 50], [25, 50], [15, 10], [85, 10], [70, 40], [30, 40]]
facets = [[0, 1], [1, 2], [2, 3], [3, 0], [4, 5], [5, 6], [6, 7], [7, 4]]
holes = [[50, 25]]
control_points = [[5, 5]]
perimeter = [0, 1, 2, 3]

geometry = sections.CustomSection(
    points, facets, holes, control_points, perimeter=perimeter
)
mesh = geometry.create_mesh(mesh_sizes=[2.0])
../_images/custom_geometry.png

Custom section geometry.

../_images/custom_mesh.png

Mesh generated from the above geometry.

RectangularSection Class

class sectionproperties.pre.sections.RectangularSection(d, b, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a rectangular section with the bottom left corner at the origin (0, 0), with depth d and width b.

Parameters:
  • d (float) – Depth (y) of the rectangle
  • b (float) – Width (x) of the rectangle
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a rectangular cross-section with a depth of 100 and width of 50, and generates a mesh with a maximum triangular area of 5:

import sectionproperties.pre.sections as sections

geometry = sections.RectangularSection(d=100, b=50)
mesh = geometry.create_mesh(mesh_sizes=[5])
../_images/rectangle_geometry.png

Rectangular section geometry.

../_images/rectangle_mesh.png

Mesh generated from the above geometry.

CircularSection Class

class sectionproperties.pre.sections.CircularSection(d, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a solid circle centered at the origin (0, 0) with diameter d and using n points to construct the circle.

Parameters:
  • d (float) – Diameter of the circle
  • n (int) – Number of points discretising the circle
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a circular cross-section with a diameter of 50 with 64 points, and generates a mesh with a maximum triangular area of 2.5:

import sectionproperties.pre.sections as sections

geometry = sections.CircularSection(d=50, n=64)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
../_images/circle_geometry.png

Circular section geometry.

../_images/circle_mesh.png

Mesh generated from the above geometry.

Chs Class

class sectionproperties.pre.sections.Chs(d, t, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a circular hollow section centered at the origin (0, 0), with diameter d and thickness t, using n points to construct the inner and outer circles.

Parameters:
  • d (float) – Outer diameter of the CHS
  • t (float) – Thickness of the CHS
  • n (int) – Number of points discretising the inner and outer circles
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a CHS discretised with 64 points, with a diameter of 48 and thickness of 3.2, and generates a mesh with a maximum triangular area of 1.0:

import sectionproperties.pre.sections as sections

geometry = sections.Chs(d=48, t=3.2, n=64)
mesh = geometry.create_mesh(mesh_sizes=[1.0])
../_images/chs_geometry.png

CHS geometry.

../_images/chs_mesh.png

Mesh generated from the above geometry.

EllipticalSection Class

class sectionproperties.pre.sections.EllipticalSection(d_y, d_x, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a solid ellipse centered at the origin (0, 0) with vertical diameter d_y and horizontal diameter d_x, using n points to construct the ellipse.

Parameters:
  • d_y (float) – Diameter of the ellipse in the y-dimension
  • d_x (float) – Diameter of the ellipse in the x-dimension
  • n (int) – Number of points discretising the ellipse
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates an elliptical cross-section with a vertical diameter of 25 and horizontal diameter of 50, with 40 points, and generates a mesh with a maximum triangular area of 1.0:

import sectionproperties.pre.sections as sections

geometry = sections.EllipticalSection(d_y=25, d_x=50, n=40)
mesh = geometry.create_mesh(mesh_sizes=[1.0])
../_images/ellipse_geometry.png

Elliptical section geometry.

../_images/ellipse_mesh.png

Mesh generated from the above geometry.

Ehs Class

class sectionproperties.pre.sections.Ehs(d_y, d_x, t, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs an elliptical hollow section centered at the origin (0, 0), with outer vertical diameter d_y, outer horizontal diameter d_x, and thickness t, using n points to construct the inner and outer ellipses.

Parameters:
  • d_y (float) – Diameter of the ellipse in the y-dimension
  • d_x (float) – Diameter of the ellipse in the x-dimension
  • t (float) – Thickness of the EHS
  • n (int) – Number of points discretising the inner and outer ellipses
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a EHS discretised with 30 points, with a outer vertical diameter of 25, outer horizontal diameter of 50, and thickness of 2.0, and generates a mesh with a maximum triangular area of 0.5:

import sectionproperties.pre.sections as sections

geometry = sections.Ehs(d_y=25, d_x=50, t=2.0, n=64)
mesh = geometry.create_mesh(mesh_sizes=[0.5])
../_images/ehs_geometry.png

EHS geometry.

../_images/ehs_mesh.png

Mesh generated from the above geometry.

Rhs Class

class sectionproperties.pre.sections.Rhs(d, b, t, r_out, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a rectangular hollow section centered at (b/2, d/2), with depth d, width b, thickness t and outer radius r_out, using n_r points to construct the inner and outer radii. If the outer radius is less than the thickness of the RHS, the inner radius is set to zero.

Parameters:
  • d (float) – Depth of the RHS
  • b (float) – Width of the RHS
  • t (float) – Thickness of the RHS
  • r_out (float) – Outer radius of the RHS
  • n_r (int) – Number of points discretising the inner and outer radii
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates an RHS with a depth of 100, a width of 50, a thickness of 6 and an outer radius of 9, using 8 points to discretise the inner and outer radii. A mesh is generated with a maximum triangular area of 2.0:

import sectionproperties.pre.sections as sections

geometry = sections.Rhs(d=100, b=50, t=6, r_out=9, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.0])
../_images/rhs_geometry.png

RHS geometry.

../_images/rhs_mesh.png

Mesh generated from the above geometry.

ISection Class

class sectionproperties.pre.sections.ISection(d, b, t_f, t_w, r, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs an I-section centered at (b/2, d/2), with depth d, width b, flange thickness t_f, web thickness t_w, and root radius r, using n_r points to construct the root radius.

Parameters:
  • d (float) – Depth of the I-section
  • b (float) – Width of the I-section
  • t_f (float) – Flange thickness of the I-section
  • t_w (float) – Web thickness of the I-section
  • r (float) – Root radius of the I-section
  • n_r (int) – Number of points discretising the root radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates an I-section with a depth of 203, a width of 133, a flange thickness of 7.8, a web thickness of 5.8 and a root radius of 8.9, using 16 points to discretise the root radius. A mesh is generated with a maximum triangular area of 3.0:

import sectionproperties.pre.sections as sections

geometry = sections.ISection(d=203, b=133, t_f=7.8, t_w=5.8, r=8.9, n_r=16)
mesh = geometry.create_mesh(mesh_sizes=[3.0])
../_images/isection_geometry.png

I-section geometry.

../_images/isection_mesh.png

Mesh generated from the above geometry.

MonoISection Class

class sectionproperties.pre.sections.MonoISection(d, b_t, b_b, t_fb, t_ft, t_w, r, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a monosymmetric I-section centered at (max(b_t, b_b)/2, d/2), with depth d, top flange width b_t, bottom flange width b_b, top flange thickness t_ft, top flange thickness t_fb, web thickness t_w, and root radius r, using n_r points to construct the root radius.

Parameters:
  • d (float) – Depth of the I-section
  • b_t (float) – Top flange width
  • b_b (float) – Bottom flange width
  • t_ft (float) – Top flange thickness of the I-section
  • t_fb (float) – Bottom flange thickness of the I-section
  • t_w (float) – Web thickness of the I-section
  • r (float) – Root radius of the I-section
  • n_r (int) – Number of points discretising the root radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a monosymmetric I-section with a depth of 200, a top flange width of 50, a top flange thickness of 12, a bottom flange width of 130, a bottom flange thickness of 8, a web thickness of 6 and a root radius of 8, using 16 points to discretise the root radius. A mesh is generated with a maximum triangular area of 3.0:

import sectionproperties.pre.sections as sections

geometry = sections.MonoISection(
    d=200, b_t=50, b_b=130, t_ft=12, t_fb=8, t_w=6, r=8, n_r=16
)
mesh = geometry.create_mesh(mesh_sizes=[3.0])
../_images/monoisection_geometry.png

I-section geometry.

../_images/monoisection_mesh.png

Mesh generated from the above geometry.

TaperedFlangeISection Class

class sectionproperties.pre.sections.TaperedFlangeISection(d, b, t_f, t_w, r_r, r_f, alpha, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Tapered Flange I-section centered at (b/2, d/2), with depth d, width b, mid-flange thickness t_f, web thickness t_w, root radius r_r, flange radius r_f and flange angle alpha, using n_r points to construct the radii.

Parameters:
  • d (float) – Depth of the Tapered Flange I-section
  • b (float) – Width of the Tapered Flange I-section
  • t_f (float) – Mid-flange thickness of the Tapered Flange I-section (measured at the point equidistant from the face of the web to the edge of the flange)
  • t_w (float) – Web thickness of the Tapered Flange I-section
  • r_r (float) – Root radius of the Tapered Flange I-section
  • r_f (float) – Flange radius of the Tapered Flange I-section
  • alpha (float) – Flange angle of the Tapered Flange I-section (degrees)
  • n_r (int) – Number of points discretising the radii
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a Tapered Flange I-section with a depth of 588, a width of 191, a mid-flange thickness of 27.2, a web thickness of 15.2, a root radius of 17.8, a flange radius of 8.9 and a flange angle of 8°, using 16 points to discretise the radii. A mesh is generated with a maximum triangular area of 20.0:

import sectionproperties.pre.sections as sections

geometry = sections.TaperedFlangeISection(
    d=588, b=191, t_f=27.2, t_w=15.2, r_r=17.8, r_f=8.9, alpha=8, n_r=16
)
mesh = geometry.create_mesh(mesh_sizes=[20.0])
../_images/taperedisection_geometry.png

I-section geometry.

../_images/taperedisection_mesh.png

Mesh generated from the above geometry.

PfcSection Class

class sectionproperties.pre.sections.PfcSection(d, b, t_f, t_w, r, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a PFC section with the bottom left corner at the origin (0, 0), with depth d, width b, flange thickness t_f, web thickness t_w and root radius r, using n_r points to construct the root radius.

Parameters:
  • d (float) – Depth of the PFC section
  • b (float) – Width of the PFC section
  • t_f (float) – Flange thickness of the PFC section
  • t_w (float) – Web thickness of the PFC section
  • r (float) – Root radius of the PFC section
  • n_r (int) – Number of points discretising the root radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a PFC section with a depth of 250, a width of 90, a flange thickness of 15, a web thickness of 8 and a root radius of 12, using 8 points to discretise the root radius. A mesh is generated with a maximum triangular area of 5.0:

import sectionproperties.pre.sections as sections

geometry = sections.PfcSection(d=250, b=90, t_f=15, t_w=8, r=12, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[5.0])
../_images/pfc_geometry.png

PFC geometry.

../_images/pfc_mesh.png

Mesh generated from the above geometry.

TaperedFlangeChannel Class

class sectionproperties.pre.sections.TaperedFlangeChannel(d, b, t_f, t_w, r_r, r_f, alpha, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Tapered Flange Channel section with the bottom left corner at the origin (0, 0), with depth d, width b, mid-flange thickness t_f, web thickness t_w, root radius r_r, flange radius r_f and flange angle alpha, using n_r points to construct the radii.

Parameters:
  • d (float) – Depth of the Tapered Flange Channel section
  • b (float) – Width of the Tapered Flange Channel section
  • t_f (float) – Mid-flange thickness of the Tapered Flange Channel section (measured at the point equidistant from the face of the web to the edge of the flange)
  • t_w (float) – Web thickness of the Tapered Flange Channel section
  • r_r (float) – Root radius of the Tapered Flange Channel section
  • r_f (float) – Flange radius of the Tapered Flange Channel section
  • alpha (float) – Flange angle of the Tapered Flange Channel section (degrees)
  • n_r (int) – Number of points discretising the radii
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a Tapered Flange Channel section with a depth of 10, a width of 3.5, a mid-flange thickness of 0.575, a web thickness of 0.475, a root radius of 0.575, a flange radius of 0.4 and a flange angle of 8°, using 16 points to discretise the radii. A mesh is generated with a maximum triangular area of 0.02:

import sectionproperties.pre.sections as sections

geometry = sections.TaperedFlangeChannel(
    d=10, b=3.5, t_f=0.575, t_w=0.475, r_r=0.575, r_f=0.4, alpha=8, n_r=16
)
mesh = geometry.create_mesh(mesh_sizes=[0.02])
../_images/taperedchannel_geometry.png

I-section geometry.

../_images/taperedchannel_mesh.png

Mesh generated from the above geometry.

TeeSection Class

class sectionproperties.pre.sections.TeeSection(d, b, t_f, t_w, r, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Tee section with the top left corner at (0, d), with depth d, width b, flange thickness t_f, web thickness t_w and root radius r, using n_r points to construct the root radius.

Parameters:
  • d (float) – Depth of the Tee section
  • b (float) – Width of the Tee section
  • t_f (float) – Flange thickness of the Tee section
  • t_w (float) – Web thickness of the Tee section
  • r (float) – Root radius of the Tee section
  • n_r (int) – Number of points discretising the root radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a Tee section with a depth of 200, a width of 100, a flange thickness of 12, a web thickness of 6 and a root radius of 8, using 8 points to discretise the root radius. A mesh is generated with a maximum triangular area of 3.0:

import sectionproperties.pre.sections as sections

geometry = sections.TeeSection(d=200, b=100, t_f=12, t_w=6, r=8, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[3.0])
../_images/tee_geometry.png

Tee section geometry.

../_images/tee_mesh.png

Mesh generated from the above geometry.

AngleSection Class

class sectionproperties.pre.sections.AngleSection(d, b, t, r_r, r_t, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs an angle section with the bottom left corner at the origin (0, 0), with depth d, width b, thickness t, root radius r_r and toe radius r_t, using n_r points to construct the radii.

Parameters:
  • d (float) – Depth of the angle section
  • b (float) – Width of the angle section
  • t (float) – Thickness of the angle section
  • r_r (float) – Root radius of the angle section
  • r_t (float) – Toe radius of the angle section
  • n_r (int) – Number of points discretising the radii
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates an angle section with a depth of 150, a width of 100, a thickness of 8, a root radius of 12 and a toe radius of 5, using 16 points to discretise the radii. A mesh is generated with a maximum triangular area of 2.0:

import sectionproperties.pre.sections as sections

geometry = sections.AngleSection(d=150, b=100, t=8, r_r=12, r_t=5, n_r=16)
mesh = geometry.create_mesh(mesh_sizes=[2.0])
../_images/angle_geometry.png

Angle section geometry.

../_images/angle_mesh.png

CeeSection Class

class sectionproperties.pre.sections.CeeSection(d, b, l, t, r_out, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Cee section with the bottom left corner at the origin (0, 0), with depth d, width b, lip l, thickness t and outer radius r_out, using n_r points to construct the radius. If the outer radius is less than the thickness of the Cee Section, the inner radius is set to zero.

Parameters:
  • d (float) – Depth of the Cee section
  • b (float) – Width of the Cee section
  • l (float) – Lip of the Cee section
  • t (float) – Thickness of the Cee section
  • r_out (float) – Outer radius of the Cee section
  • n_r (int) – Number of points discretising the outer radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)
Raises:

Exception – Lip length must be greater than the outer radius

The following example creates a Cee section with a depth of 125, a width of 50, a lip of 30, a thickness of 1.5 and an outer radius of 6, using 8 points to discretise the radius. A mesh is generated with a maximum triangular area of 0.25:

import sectionproperties.pre.sections as sections

geometry = sections.CeeSection(d=125, b=50, l=30, t=1.5, r_out=6, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[0.25])
../_images/cee_geometry.png

Cee section geometry.

../_images/cee_mesh.png

ZedSection Class

class sectionproperties.pre.sections.ZedSection(d, b_l, b_r, l, t, r_out, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Zed section with the bottom left corner at the origin (0, 0), with depth d, left flange width b_l, right flange width b_r, lip l, thickness t and outer radius r_out, using n_r points to construct the radius. If the outer radius is less than the thickness of the Zed Section, the inner radius is set to zero.

Parameters:
  • d (float) – Depth of the Zed section
  • b_l (float) – Left flange width of the Zed section
  • b_r (float) – Right flange width of the Zed section
  • l (float) – Lip of the Zed section
  • t (float) – Thickness of the Zed section
  • r_out (float) – Outer radius of the Zed section
  • n_r (int) – Number of points discretising the outer radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)
Raises:

Exception – Lip length must be greater than the outer radius

The following example creates a Zed section with a depth of 100, a left flange width of 40, a right flange width of 50, a lip of 20, a thickness of 1.2 and an outer radius of 5, using 8 points to discretise the radius. A mesh is generated with a maximum triangular area of 0.15:

import sectionproperties.pre.sections as sections

geometry = sections.ZedSection(d=100, b_l=40, b_r=50, l=20, t=1.2, r_out=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[0.15])
../_images/zed_geometry.png

Zed section geometry.

../_images/zed_mesh.png

CruciformSection Class

class sectionproperties.pre.sections.CruciformSection(d, b, t, r, n_r, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a cruciform section centered at the origin (0, 0), with depth d, width b, thickness t and root radius r, using n_r points to construct the root radius.

Parameters:
  • d (float) – Depth of the cruciform section
  • b (float) – Width of the cruciform section
  • t (float) – Thickness of the cruciform section
  • r (float) – Root radius of the cruciform section
  • n_r (int) – Number of points discretising the root radius
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a cruciform section with a depth of 250, a width of 175, a thickness of 12 and a root radius of 16, using 16 points to discretise the radius. A mesh is generated with a maximum triangular area of 5.0:

import sectionproperties.pre.sections as sections

geometry = sections.CruciformSection(d=250, b=175, t=12, r=16, n_r=16)
mesh = geometry.create_mesh(mesh_sizes=[5.0])
../_images/cruciform_geometry.png

Cruciform section geometry.

../_images/cruciform_mesh.png

PolygonSection Class

class sectionproperties.pre.sections.PolygonSection(d, t, n_sides, r_in=0, n_r=1, rot=0, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a regular hollow polygon section centered at (0, 0), with a pitch circle diameter of bounding polygon d, thickness t, number of sides n_sides and an optional inner radius r_in, using n_r points to construct the inner and outer radii (if radii is specified).

Parameters:
  • d (float) – Pitch circle diameter of the outer bounding polygon (i.e. diameter of circle that passes through all vertices of the outer polygon)
  • t (float) – Thickness of the polygon section wall
  • r_in (float) – Inner radius of the polygon corners. By default, if not specified, a polygon with no corner radii is generated.
  • n_r (int) – Number of points discretising the inner and outer radii, ignored if no inner radii is specified
  • rot – Initial counterclockwise rotation in degrees. By default bottom face is aligned with x axis.
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)
Raises:

Exception – Number of sides in polygon must be greater than or equal to 3

The following example creates an Octagonal section (8 sides) with a diameter of 200, a thickness of 6 and an inner radius of 20, using 12 points to discretise the inner and outer radii. A mesh is generated with a maximum triangular area of 5:

import sectionproperties.pre.sections as sections

geometry = sections.PolygonSection(d=200, t=6, n_sides=8, r_in=20, n_r=12)
mesh = geometry.create_mesh(mesh_sizes=[5])
../_images/polygon_geometry.png

Octagonal section geometry.

../_images/polygon_mesh.png

Mesh generated from the above geometry.

BoxGirderSection Class

class sectionproperties.pre.sections.BoxGirderSection(d, b_t, b_b, t_ft, t_fb, t_w, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Box Girder section centered at at (max(b_t, b_b)/2, d/2), with depth d, top width b_t, bottom width b_b, top flange thickness t_ft, bottom flange thickness t_fb and web thickness t_w.

Parameters:
  • d (float) – Depth of the Box Girder section
  • b_t (float) – Top width of the Box Girder section
  • b_b (float) – Bottom width of the Box Girder section
  • t_ft (float) – Top lange thickness of the Box Girder section
  • t_fb (float) – Bottom flange thickness of the Box Girder section
  • t_w (float) – Web thickness of the Box Girder section
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a Box Gider section with a depth of 1200, a top width of 1200, a bottom width of 400, a top flange thickness of 16, a bottom flange thickness of 12 and a web thickness of 8. A mesh is generated with a maximum triangular area of 5.0:

import sectionproperties.pre.sections as sections

geometry = sections.BoxGirderSection(d=1200, b_t=1200, b_b=400, t_ft=100, t_fb=80, t_w=50)
mesh = geometry.create_mesh(mesh_sizes=[200.0])
../_images/box_girder_geometry.png

Box Girder geometry.

../_images/box_girder_mesh.png

Mesh generated from the above geometry.

MergedSection Class

class sectionproperties.pre.sections.MergedSection(sections)[source]

Bases: sectionproperties.pre.sections.Geometry

Merges a number of section geometries into one geometry. Note that for the meshing algorithm to work, there needs to be connectivity between all regions of the provided geometries. Overlapping of geometries is permitted.

Parameters:sections (list[Geometry]) – A list of geometry objects to merge into one Geometry object

The following example creates a combined cross-section with a 150x100x6 RHS placed on its side on top of a 200UB25.4. A mesh is generated with a maximum triangle size of 5.0 for the I-section and 2.5 for the RHS:

import sectionproperties.pre.sections as sections

isection = sections.ISection(d=203, b=133, t_f=7.8, t_w=5.8, r=8.9, n_r=8)
box = sections.Rhs(d=100, b=150, t=6, r_out=15, n_r=8, shift=[-8.5, 203])

geometry = sections.MergedSection([isection, box])
geometry.clean_geometry()
mesh = geometry.create_mesh(mesh_sizes=[5.0, 2.5])
../_images/merged_geometry.png

Merged section geometry.

../_images/merged_mesh.png

pre Module

Material Class

class sectionproperties.pre.pre.Material(name, elastic_modulus, poissons_ratio, yield_strength, color='w')[source]

Bases: object

Class for structural materials.

Provides a way of storing material properties related to a specific material. The color can be a multitude of different formats, refer to https://matplotlib.org/api/colors_api.html and https://matplotlib.org/examples/color/named_colors.html for more information.

Parameters:
  • name (string) – Material name
  • elastic_modulus (float) – Material modulus of elasticity
  • poissons_ratio (float) – Material Poisson’s ratio
  • yield_strength (float) – Material yield strength
  • color (matplotlib.colors) – Material color for rendering
Variables:
  • name (string) – Material name
  • elastic_modulus (float) – Material modulus of elasticity
  • poissons_ratio (float) – Material Poisson’s ratio
  • shear_modulus (float) – Material shear modulus, derived from the elastic modulus and Poisson’s ratio assuming an isotropic material
  • yield_strength (float) – Material yield strength
  • color (matplotlib.colors) – Material color for rendering

The following example creates materials for concrete, steel and timber:

from sectionproperties.pre.pre import Material

concrete = Material(
    name='Concrete', elastic_modulus=30.1e3, poissons_ratio=0.2, yield_strength=32,
        color='lightgrey'
)
steel = Material(
    name='Steel', elastic_modulus=200e3, poissons_ratio=0.3, yield_strength=500,
        color='grey'
)
timber = Material(
    name='Timber', elastic_modulus=8e3, poissons_ratio=0.35, yield_strength=20,
        color='burlywood'
)

GeometryCleaner Class

class sectionproperties.pre.pre.GeometryCleaner(geometry, verbose)[source]

Bases: object

Class for cleaning Geometry objects.

Parameters:
  • geometry (Geometry) – Geometry object to clean
  • verbose (bool) – If set to true, information related to the geometry cleaning process is printed to the terminal.

Provides methods to clean various aspects of the geometry including:

  • Zipping nodes - Find nodes that are close together (relative and absolute tolerance) and deletes one of the nodes and rejoins the facets to the remaining node.
  • Removing zero length facets - Removes facets that start and end at the same point.
  • Remove duplicate facets - Removes facets that have the same starting and ending point as an existing facet.
  • Removing overlapping facets - Searches for facets that overlap each other, given a tolerance angle, and reconstructs a unique set of facets along the overlapping region.
  • Remove unused points - Removes points that are not connected to any facets.
  • Intersect facets - Searches for intersections between two facets and adds the intersection point to the points list and splits the intersected facets.

Note that a geometry cleaning method is provided to all Geometry objects.

Variables:
  • geometry (Geometry) – Geometry object to clean
  • verbose (bool) – If set to true, information related to the geometry cleaning process is printed to the terminal.

The following example creates a back-to-back 200PFC geometry, rotates the geometry by 30 degrees, and cleans the geometry before meshing:

import sectionproperties.pre.sections as sections

pfc_right = sections.PfcSection(d=203, b=133, t_f=7.8, t_w=5.8, r=8.9, n_r=8)
pfc_left = sections.PfcSection(d=203, b=133, t_f=7.8, t_w=5.8, r=8.9, n_r=8)
pfc_left.mirror_section(axis='y', mirror_point=[0, 0])
geometry = sections.MergedSection([pfc_left, pfc_right])
geometry.rotate_section(angle=30)
geometry.clean_geometry(verbose=True)
mesh = geometry.create_mesh(mesh_sizes=[5, 5])

Warning

If the geometry were not cleaned in the previous example, the meshing algorithm would crash (most likely return a segment error). Cleaning the geometry is always recommended when creating a merged section, which may result in overlapping or intersecting facets, or duplicate nodes.

clean_geometry()[source]

Performs a full geometry clean on the geometry object.

intersect_facets()[source]

Searches through all facet combinations and finds facets that intersect each other. The intersection point is added and the facets rebuilt.

is_duplicate_facet(fct1, fct2)[source]

Checks to see if to facets are duplicates.

Parameters:
  • fct1 (list[int, int]) – First facet to compare
  • fct2 (list[int, int]) – Second facet to compare
Returns:

Whether or not the facets are identical

Return type:

bool

is_intersect(p, q, r, s)[source]

Determines if the line segment p->p+r intersects q->q+s. Implements Gareth Rees’s answer: https://stackoverflow.com/questions/563198.

Parameters:
  • p (numpy.ndarray [float, float]) – Starting point of the first line segment
  • q (numpy.ndarray [float, float]) – Starting point of the second line segment
  • r (numpy.ndarray [float, float]) – Vector of the first line segment
  • s (numpy.ndarray [float, float]) – Vector of the second line segment
Returns:

The intersection point of the line segments. If there is no intersection, returns None.

Return type:

numpy.ndarray [float, float]

is_overlap(p, q, r, s, fct1, fct2)[source]

Determines if the line segment p->p+r overlaps q->q+s. Implements Gareth Rees’s answer: https://stackoverflow.com/questions/563198.

Parameters:
  • p (numpy.ndarray [float, float]) – Starting point of the first line segment
  • q (numpy.ndarray [float, float]) – Starting point of the second line segment
  • r (numpy.ndarray [float, float]) – Vector of the first line segment
  • s (numpy.ndarray [float, float]) – Vector of the second line segment
  • fct1 – sadkjas;dkas;dj
Returns:

A list containing the points required for facet rebuilding. If there is no rebuild to be done, returns None.

Return type:

list[list[float, float]]

remove_duplicate_facets()[source]

Searches through all facets and removes facets that are duplicates, independent of the point order.

remove_overlapping_facets()[source]

Searches through all facet combinations and fixes facets that overlap within a tolerance.

remove_point_id(point_id)[source]

Removes point point_id from the points list and renumbers the references to points after point_id in the facet list.

Parameters:point_id (int) – Index of point to be removed
remove_unused_points()[source]

Searches through all facets and removes points that are not connected to any facets.

remove_zero_length_facets()[source]

Searches through all facets and removes those that have the same starting and ending point.

replace_point_id(id_old, id_new)[source]

Searches all facets and replaces references to point id_old with id_new.

Parameters:
  • id_old (int) – Point index to be replaced
  • id_new (int) – Point index to replace point id_old
zip_points(atol=1e-08, rtol=1e-05)[source]

Zips points that are close to each other. Searches through the point list and merges two points if there are deemed to be sufficiently close. The average value of the coordinates is used for the new point. One of the points is deleted from the point list and the facet list is updated to remove references to the old points and renumber the remaining point indices in the facet list.

Parameters:
  • atol (float) – Absolute tolerance for point zipping
  • rtol (float) – Relative tolerance (to geometry extents) for point zipping

pre Functions

sectionproperties.pre.pre.create_mesh(points, facets, holes, control_points, mesh_sizes)[source]

Creates a quadratic triangular mesh using the meshpy module, which utilises the code ‘Triangle’, by Jonathan Shewchuk.

Parameters:
  • points (list[list[int, int]]) – List of points (x, y) defining the vertices of the cross-section
  • facets – List of point index pairs (p1, p2) defining the edges of the cross-section
  • holes (list[list[float, float]]) – List of points (x, y) defining the locations of holes within the cross-section. If there are no holes, provide an empty list [].
  • control_points (list[list[float, float]]) – A list of points (x, y) that define different regions of the cross-section. A control point is an arbitrary point within a region enclosed by facets.
  • mesh_sizes (list[float]) – List of maximum element areas for each region defined by a control point
Returns:

Object containing generated mesh data

Return type:

meshpy.triangle.MeshInfo

offset Module

sectionproperties.pre.offset.offset_perimeter(geometry, offset, side='left', plot_offset=False)[source]

Offsets the perimeter of a geometry of a Geometry object by a certain distance. Note that the perimeter facet list must be entered in a consecutive order.

Parameters:
  • geometry (Geometry) – Cross-section geometry object
  • offset (float) – Offset distance for the perimeter
  • side (string) – Side of the perimeter offset, either ‘left’ or ‘right’. E.g. ‘left’ for a counter-clockwise offsets the perimeter inwards.
  • plot_offset (bool) – If set to True, generates a plot comparing the old and new geometry

The following example ‘corrodes’ a 200UB25 I-section by 1.5 mm and compares a few of the section properties:

import sectionproperties.pre.sections as sections
from sectionproperties.pre.offset import offset_perimeter
from sectionproperties.analysis.cross_section import CrossSection

# calculate original section properties
original_geometry = sections.ISection(d=203, b=133, t_f=7.8, t_w=5.8, r=8.9, n_r=16)
original_mesh = original_geometry.create_mesh(mesh_sizes=[3.0])
original_section = CrossSection(original_geometry, original_mesh)
original_section.calculate_geometric_properties()
original_area = original_section.get_area()
(original_ixx, _, _) = original_section.get_ic()

# calculate corroded section properties
corroded_geometry = offset_perimeter(original_geometry, 1.5, plot_offset=True)
corroded_mesh = corroded_geometry.create_mesh(mesh_sizes=[3.0])
corroded_section = CrossSection(corroded_geometry, corroded_mesh)
corroded_section.calculate_geometric_properties()
corroded_area = corroded_section.get_area()
(corroded_ixx, _, _) = corroded_section.get_ic()

# compare section properties
print("Area reduction = {0:.2f}%".format(
    100 * (original_area - corroded_area) / original_area))
print("Ixx reduction = {0:.2f}%".format(
    100 *(original_ixx - corroded_ixx) / original_ixx))

The following plot is generated by the above example:

../_images/offset_example.png

200UB25 with 1.5 mm corrosion.

The following is printed to the terminal:

Area reduction = 41.97%
Ixx reduction = 39.20%

nastran_sections Module

This module contains cross-sections as defined by Nastran and Nastran-based programs, such as MYSTRAN and ASTROS.

BARSection Class

class sectionproperties.pre.nastran_sections.BARSection(DIM1, DIM2, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a BAR section with the center at the origin (0, 0), with two parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] [5] for definition of parameters. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of bar
  • DIM2 (float) – Depth (y) of bar
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a BAR cross-section with a depth of 1.5 and width of 2.0, and generates a mesh with a maximum triangular area of 0.001:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.BARSection(DIM1=2.0, DIM2=1.5)
mesh = geometry.create_mesh(mesh_sizes=[0.001])
../_images/bar_geometry.png

BAR section geometry.

../_images/bar_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

BOXSection Class

class sectionproperties.pre.nastran_sections.BOXSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a BOX section with the center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] [5] for definition of parameters. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of box
  • DIM2 (float) – Depth (y) of box
  • DIM3 (float) – Thickness of box in y direction
  • DIM4 (float) – Thickness of box in x direction
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a BOX cross-section with a depth of 3.0 and width of 4.0, and generates a mesh with a maximum triangular area of 0.001:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.BOXSection(DIM1=4.0, DIM2=3.0, DIM3=0.375, DIM4=0.5)
mesh = geometry.create_mesh(mesh_sizes=[0.001])
../_images/box_geometry.png

BOX section geometry.

../_images/box_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

BOX1Section Class

class sectionproperties.pre.nastran_sections.BOX1Section(DIM1, DIM2, DIM3, DIM4, DIM5, DIM6, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a BOX1 section with the center at the origin (0, 0), with six parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of box
  • DIM2 (float) – Depth (y) of box
  • DIM3 (float) – Thickness of top wall
  • DIM4 (float) – Thickness of bottom wall
  • DIM5 (float) – Thickness of left wall
  • DIM6 (float) – Thickness of right wall
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a BOX1 cross-section with a depth of 3.0 and width of 4.0, and generates a mesh with a maximum triangular area of 0.007:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.BOX1Section(
    DIM1=4.0, DIM2=3.0, DIM3=0.375, DIM4=0.5, DIM5=0.25, DIM6=0.75
)
mesh = geometry.create_mesh(mesh_sizes=[0.007])
../_images/box1_geometry.png

BOX1 section geometry.

../_images/box1_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

CHANSection Class

class sectionproperties.pre.nastran_sections.CHANSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a CHAN (C-Channel) section with the web’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of the CHAN-section
  • DIM2 (float) – Depth (y) of the CHAN-section
  • DIM3 (float) – Thickness of web (vertical portion)
  • DIM4 (float) – Thickness of flanges (top/bottom portion)
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a CHAN cross-section with a depth of 4.0 and width of 2.0, and generates a mesh with a maximum triangular area of 0.008:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.CHANSection(DIM1=2.0, DIM2=4.0, DIM3=0.25, DIM4=0.5)
mesh = geometry.create_mesh(mesh_sizes=[0.008])
../_images/chan_geometry.png

CHAN section geometry.

../_images/chan_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

CHAN1Section Class

class sectionproperties.pre.nastran_sections.CHAN1Section(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a CHAN1 (C-Channel) section with the web’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of channels
  • DIM2 (float) – Thicknesss (x) of web
  • DIM3 (float) – Spacing between channels (length of web)
  • DIM4 (float) – Depth (y) of CHAN1-section
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a CHAN1 cross-section with a depth of 4.0 and width of 1.75, and generates a mesh with a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.CHAN1Section(DIM1=0.75, DIM2=1.0, DIM3=3.5, DIM4=4.0)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/chan1_geometry.png

CHAN1 section geometry.

../_images/chan1_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

CHAN2Section Class

class sectionproperties.pre.nastran_sections.CHAN2Section(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a CHAN2 (C-Channel) section with the bottom web’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Thickness of channels
  • DIM2 (float) – Thickness of web
  • DIM3 (float) – Depth (y) of CHAN2-section
  • DIM4 (float) – Width (x) of CHAN2-section
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a CHAN2 cross-section with a depth of 2.0 and width of 4.0, and generates a mesh with a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.CHAN2Section(DIM1=0.375, DIM2=0.5, DIM3=2.0, DIM4=4.0)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/chan2_geometry.png

CHAN2 section geometry.

../_images/chan2_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

CROSSSection Class

class sectionproperties.pre.nastran_sections.CROSSSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs Nastran’s cruciform/cross section with the intersection’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Twice the width of horizontal member protruding from the vertical center member
  • DIM2 (float) – Thickness of the vertical member
  • DIM3 (float) – Depth (y) of the CROSS-section
  • DIM4 (float) – Thickness of the horizontal members
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a rectangular cross-section with a depth of 3.0 and width of 1.875, and generates a mesh with a maximum triangular area of 0.008:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.CROSSSection(DIM1=1.5, DIM2=0.375, DIM3=3.0, DIM4=0.25)
mesh = geometry.create_mesh(mesh_sizes=[0.008])
../_images/cross_geometry.png

Cruciform/cross section geometry.

../_images/cross_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

DBOXSection Class

class sectionproperties.pre.nastran_sections.DBOXSection(DIM1, DIM2, DIM3, DIM4, DIM5, DIM6, DIM7, DIM8, DIM9, DIM10, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a DBOX section with the center at the origin (0, 0), with ten parameters defining dimensions. See MSC Nastran documentation [1] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of the DBOX-section
  • DIM2 (float) – Depth (y) of the DBOX-section
  • DIM3 (float) – Width (x) of left-side box
  • DIM4 (float) – Thickness of left wall
  • DIM5 (float) – Thickness of center wall
  • DIM6 (float) – Thickness of right wall
  • DIM7 (float) – Thickness of top left wall
  • DIM8 (float) – Thickness of bottom left wall
  • DIM9 (float) – Thickness of top right wall
  • DIM10 (float) – Thickness of bottom right wall
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a DBOX cross-section with a depth of 3.0 and width of 8.0, and generates a mesh with a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.DBOXSection(
    DIM1=8.0, DIM2=3.0, DIM3=3.0, DIM4=0.5, DIM5=0.625, DIM6=0.75, DIM7=0.375, DIM8=0.25,
    DIM9=0.5, DIM10=0.375
)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/dbox_geometry.png

DBOX section geometry.

../_images/dbox_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

FCROSSSection Class

class sectionproperties.pre.nastran_sections.FCROSSSection(DIM1, DIM2, DIM3, DIM4, DIM5, DIM6, DIM7, DIM8, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a flanged cruciform/cross section with the intersection’s middle center at the origin (0, 0), with eight parameters defining dimensions. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Depth (y) of flanged cruciform
  • DIM2 (float) – Width (x) of flanged cruciform
  • DIM3 (float) – Thickness of vertical web
  • DIM4 (float) – Thickness of horizontal web
  • DIM5 (float) – Length of flange attached to vertical web
  • DIM6 (float) – Thickness of flange attached to vertical web
  • DIM7 (float) – Length of flange attached to horizontal web
  • DIM8 (float) – Thickness of flange attached to horizontal web
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example demonstrates the creation of a flanged cross section:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.FCROSSSection(
    DIM1=9.0, DIM2=6.0, DIM3=0.75, DIM4=0.625, DIM5=2.1, DIM6=0.375, DIM7=4.5, DIM8=0.564
)
mesh = geometry.create_mesh(mesh_sizes=[0.03])
../_images/fcross_geometry.png

Flanged Cruciform/cross section geometry.

../_images/fcross_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (list[float, float]) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

GBOXSection Class

class sectionproperties.pre.nastran_sections.GBOXSection(DIM1, DIM2, DIM3, DIM4, DIM5, DIM6, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a GBOX section with the center at the origin (0, 0), with six parameters defining dimensions. See ASTROS documentation [5] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of the GBOX-section
  • DIM2 (float) – Depth (y) of the GBOX-section
  • DIM3 (float) – Thickness of top flange
  • DIM4 (float) – Thickness of bottom flange
  • DIM5 (float) – Thickness of webs
  • DIM6 (float) – Spacing between webs
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a GBOX cross-section with a depth of 2.5 and width of 6.0, and generates a mesh with a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.GBOXSection(
    DIM1=6.0, DIM2=2.5, DIM3=0.375, DIM4=0.25, DIM5=0.625, DIM6=1.0
)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/gbox_geometry.png

GBOX section geometry.

../_images/gbox_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

HSection Class

class sectionproperties.pre.nastran_sections.HSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a H section with the middle web’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Spacing between vertical flanges (length of web)
  • DIM2 (float) – Twice the thickness of the vertical flanges
  • DIM3 (float) – Depth (y) of the H-section
  • DIM4 (float) – Thickness of the middle web
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a H cross-section with a depth of 3.5 and width of 2.75, and generates a mesh with a maximum triangular area of 0.005:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.HSection(DIM1=2.0, DIM2=0.75, DIM3=3.5, DIM4=0.25)
mesh = geometry.create_mesh(mesh_sizes=[0.005])
../_images/h_geometry.png

H section geometry.

../_images/h_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

HATSection Class

class sectionproperties.pre.nastran_sections.HATSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Hat section with the top most section’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Note that HAT in ASTROS is actually HAT1 in this code. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Depth (y) of HAT-section
  • DIM2 (float) – Thickness of HAT-section
  • DIM3 (float) – Width (x) of top most section
  • DIM4 (float) – Width (x) of bottom sections
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a HAT cross-section with a depth of 1.25 and width of 2.5, and generates a mesh with a maximum triangular area of 0.001:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.HATSection(DIM1=1.25, DIM2=0.25, DIM3=1.5, DIM4=0.5)
mesh = geometry.create_mesh(mesh_sizes=[0.001])
../_images/hat_geometry.png

HAT section geometry.

../_images/hat_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the origin by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

HAT1Section Class

class sectionproperties.pre.nastran_sections.HAT1Section(DIM1, DIM2, DIM3, DIM4, DIM5, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a HAT1 section with the bottom plate’s bottom center at the origin (0, 0), with five parameters defining dimensions. See Nastran documentation [1] [2] [3] [5] for definition of parameters. Note that in ASTROS, HAT1 is called HAT. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width(x) of the HAT1-section
  • DIM2 (float) – Depth (y) of the HAT1-section
  • DIM3 (float) – Width (x) of hat’s top flange
  • DIM4 (float) – Thickness of hat stiffener
  • DIM5 (float) – Thicknesss of bottom plate
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a HAT1 cross-section with a depth of 2.0 and width of 4.0, and generates a mesh with a maximum triangular area of 0.005:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.HAT1Section(DIM1=4.0, DIM2=2.0, DIM3=1.5, DIM4=0.1875, DIM5=0.375)
mesh = geometry.create_mesh(mesh_sizes=[0.005])
../_images/hat1_geometry.png

HAT1 section geometry.

../_images/hat1_mesh.png

Mesh generated from the above geometry.

create_mesh(mesh_sizes)[source]

Creates a quadratic triangular mesh from the Geometry object. This is overloaded here to allow specifying only one mesh_size which is used for both regions in the Hat1 section.

Parameters:mesh_sizes – A list of maximum element areas corresponding to each region within the cross-section geometry.
Returns:Object containing generated mesh data
Return type:meshpy.triangle.MeshInfo
Raises:AssertionError – If the number of mesh sizes does not match the number of regions
getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the origin by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

HEXASection Class

class sectionproperties.pre.nastran_sections.HEXASection(DIM1, DIM2, DIM3, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a HEXA (hexagon) section with the center at the origin (0, 0), with three parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Spacing between bottom right point and right most point
  • DIM2 (float) – Width (x) of hexagon
  • DIM3 (float) – Depth (y) of hexagon
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a rectangular cross-section with a depth of 1.5 and width of 2.0, and generates a mesh with a maximum triangular area of 0.005:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.HEXASection(DIM1=0.5, DIM2=2.0, DIM3=1.5)
mesh = geometry.create_mesh(mesh_sizes=[0.005])
../_images/hexa_geometry.png

HEXA section geometry.

../_images/hexa_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

NISection Class

class sectionproperties.pre.nastran_sections.NISection(DIM1, DIM2, DIM3, DIM4, DIM5, DIM6, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs Nastran’s I section with the bottom flange’s middle center at the origin (0, 0), with six parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for definition of parameters. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Depth(y) of the I-section
  • DIM2 (float) – Width (x) of bottom flange
  • DIM3 (float) – Width (x) of top flange
  • DIM4 (float) – Thickness of web
  • DIM5 (float) – Thickness of bottom web
  • DIM6 (float) – Thickness of top web
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a Nastran I cross-section with a depth of 5.0, and generates a mesh with a maximum triangular area of 0.008:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.NISection(
    DIM1=5.0, DIM2=2.0, DIM3=3.0, DIM4=0.25, DIM5=0.375, DIM6=0.5
)
mesh = geometry.create_mesh(mesh_sizes=[0.008])
../_images/ni_geometry.png

Nastran’s I section geometry.

../_images/ni_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

I1Section Class

class sectionproperties.pre.nastran_sections.I1Section(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a I1 section with the web’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Twice distance from web end to flange end
  • DIM2 (float) – Thickness of web
  • DIM3 (float) – Length of web (spacing between flanges)
  • DIM4 (float) – Depth (y) of the I1-section
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a I1 cross-section with a depth of 5.0 and width of 1.75, and generates a mesh with a maximum triangular area of 0.02:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.I1Section(DIM1=1.0, DIM2=0.75, DIM3=4.0, DIM4=5.0)
mesh = geometry.create_mesh(mesh_sizes=[0.02])
../_images/i1_geometry.png

I1 section geometry.

../_images/i1_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

LSection Class

class sectionproperties.pre.nastran_sections.LSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a L section with the intersection’s center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of the L-section
  • DIM2 (float) – Depth (y) of the L-section
  • DIM3 (float) – Thickness of flange (horizontal portion)
  • DIM4 (float) – Thickness of web (vertical portion)
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a L cross-section with a depth of 6.0 and width of 3.0, and generates a mesh with a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.LSection(DIM1=3.0, DIM2=6.0, DIM3=0.375, DIM4=0.625)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/l_geometry.png

L section geometry.

../_images/l_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

RODSection Class

class sectionproperties.pre.nastran_sections.RODSection(DIM1, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a circular rod section with the center at the origin (0, 0), with one parameter defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Radius of the circular rod section
  • n (int) – Number of points discretising the circle
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a circular rod with a radius of 3.0 and 50 points discretizing the boundary, and generates a mesh with a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.RODSection(DIM1=3.0, n=50)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/rod_geometry.png

Rod section geometry.

../_images/rod_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:
  • DIM1 (float) – Radius of the circular rod section
  • shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:

Stress evaluation points relative to shifted origin - C, D, E, F

TSection Class

class sectionproperties.pre.nastran_sections.TSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a T section with the top flange’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] [5] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of top flange
  • DIM2 (float) – Depth (y) of the T-section
  • DIM3 (float) – Thickness of top flange
  • DIM4 (float) – Thickness of web
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a T cross-section with a depth of 4.0 and width of 3.0, and generates a mesh with a maximum triangular area of 0.001:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.TSection(DIM1=3.0, DIM2=4.0, DIM3=0.375, DIM4=0.25)
mesh = geometry.create_mesh(mesh_sizes=[0.001])
../_images/t_geometry.png

T section geometry.

../_images/t_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

T1Section Class

class sectionproperties.pre.nastran_sections.T1Section(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a T1 section with the right flange’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Depth (y) of T1-section
  • DIM2 (float) – Length (x) of web
  • DIM3 (float) – Thickness of right flange
  • DIM4 (float) – Thickness of web
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a T1 cross-section with a depth of 3.0 and width of 3.875, and generates a mesh with a maximum triangular area of 0.001:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.T1Section(DIM1=3.0, DIM2=3.5, DIM3=0.375, DIM4=0.25)
mesh = geometry.create_mesh(mesh_sizes=[0.001])
../_images/t1_geometry.png

T1 section geometry.

../_images/t1_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

T2Section Class

class sectionproperties.pre.nastran_sections.T2Section(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a T2 section with the bottom flange’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of T2-section
  • DIM2 (float) – Depth (y) of T2-section
  • DIM3 (float) – Thickness of bottom flange
  • DIM4 (float) – Thickness of web
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a T2 cross-section with a depth of 4.0 and width of 3.0, and generates a mesh with a maximum triangular area of 0.005:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.T2Section(DIM1=3.0, DIM2=4.0, DIM3=0.375, DIM4=0.5)
mesh = geometry.create_mesh(mesh_sizes=[0.005])
../_images/t2_geometry.png

T2 section geometry.

../_images/t2_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

TUBESection Class

class sectionproperties.pre.nastran_sections.TUBESection(DIM1, DIM2, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a circular tube section with the center at the origin (0, 0), with two parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Outer radius of the circular tube section
  • DIM2 (float) – Inner radius of the circular tube section
  • n (int) – Number of points discretising the circle
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a circular tube cross-section with an outer radius of 3.0 and an inner radius of 2.5, and generates a mesh with 37 points discretizing the boundaries and a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.TUBESection(DIM1=3.0, DIM2=2.5, n=37)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/tube_geometry.png

TUBE section geometry.

../_images/tube_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

TUBE2Section Class

class sectionproperties.pre.nastran_sections.TUBE2Section(DIM1, DIM2, n, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a circular TUBE2 section with the center at the origin (0, 0), with two parameters defining dimensions. See MSC Nastran documentation [1] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Outer radius of the circular tube section
  • DIM2 (float) – Thickness of wall
  • n (int) – Number of points discretising the circle
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a ciruclar TUBE2 cross-section with an outer radius of 3.0 and a wall thickness of 0.5, and generates a mesh with 37 point discritizing the boundary and a maximum triangular area of 0.01:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.TUBE2Section(DIM1=3.0, DIM2=0.5, n=37)
mesh = geometry.create_mesh(mesh_sizes=[0.01])
../_images/tube2_geometry.png

TUBE2 section geometry.

../_images/tube2_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

ZSection Class

class sectionproperties.pre.nastran_sections.ZSection(DIM1, DIM2, DIM3, DIM4, shift=[0, 0])[source]

Bases: sectionproperties.pre.sections.Geometry

Constructs a Z section with the web’s middle center at the origin (0, 0), with four parameters defining dimensions. See Nastran documentation [1] [2] [3] [4] for more details. Added by JohnDN90.

Parameters:
  • DIM1 (float) – Width (x) of horizontal members
  • DIM2 (float) – Thickness of web
  • DIM3 (float) – Spacing between horizontal members (length of web)
  • DIM4 (float) – Depth (y) of Z-section
  • shift (list[float, float]) – Vector that shifts the cross-section by (x, y)

The following example creates a rectangular cross-section with a depth of 4.0 and width of 2.75, and generates a mesh with a maximum triangular area of 0.005:

import sectionproperties.pre.nastran_sections as nsections

geometry = nsections.ZSection(DIM1=1.125, DIM2=0.5, DIM3=3.5, DIM4=4.0)
mesh = geometry.create_mesh(mesh_sizes=[0.005])
../_images/z_geometry.png

Z section geometry.

../_images/z_mesh.png

Mesh generated from the above geometry.

getStressPoints(shift=(0.0, 0.0))[source]

Returns the coordinates of the stress evaluation points relative to the origin of the cross-section. The shift parameter can be used to make the coordinates relative to the centroid or the shear center.

Parameters:shift (tuple(float, float)) – Vector that shifts the cross-section by (x, y)
Returns:Stress evaluation points relative to shifted origin - C, D, E, F

References

[1](1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22) MSC Nastran Quick Reference Guide 2012, PBEAML - Simple Beam Cross-Section Property, pp. 2890-2894 https://simcompanion.mscsoftware.com/infocenter/index?page=content&id=DOC10351
[2](1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20) Siemens NX Nastran 12 Quick Reference Guide, PBEAML, pp. 16-59 - 16-62 https://docs.plm.automation.siemens.com/data_services/resources/nxnastran/12/help/tdoc/en_US/pdf/QRG.pdf
[3](1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20) AutoDesk Nastran Online Documentation, Nastran Reference Guide, Section 4 - Bulk Data, PBEAML http://help.autodesk.com/view/NSTRN/2018/ENU/?guid=GUID-B7044BA7-3C26-49DA-9EE7-DA7505FD4B2C
[4](1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18) Users Reference Manual for the MYSTRAN General Purpose Finite Element Structural Analysis Computer Program, Jan. 2019, Section 6.4.1.53 - PBARL, pp. 154-156 https://www.mystran.com/Executable/MYSTRAN-Users-Manual.pdf
[5](1, 2, 3, 4, 5) Astros Enhancements - Volume III - Astros Theoretical Manual, Section 5.1.3.2, pp. 56 https://apps.dtic.mil/dtic/tr/fulltext/u2/a308134.pdf

Analysis Package

cross_section Module

CrossSection Class

class sectionproperties.analysis.cross_section.CrossSection(geometry, mesh, materials=None, time_info=False)[source]

Bases: object

Class for structural cross-sections.

Stores the finite element geometry, mesh and material information and provides methods to compute the cross-section properties. The element type used in this program is the six-noded quadratic triangular element.

The constructor extracts information from the provided mesh object and creates and stores the corresponding Tri6 finite element objects.

Parameters:
  • geometry (Geometry) – Cross-section geometry object used to generate the mesh
  • mesh (meshpy.triangle.MeshInfo) – Mesh object returned by meshpy
  • materials (list[Material]) – A list of material properties corresponding to various regions in the geometry and mesh. Note that if materials are specified, the number of material objects ust equal the number of regions in the geometry. If no materials are specified, only a purely geometric analysis can take place, and all regions will be assigned a default material with an elastic modulus and yield strength equal to 1, and a Poisson’s ratio equal to 0.
  • time_info (bool) – If set to True, a detailed description of the computation and the time cost is printed to the terminal.

The following example creates a CrossSection object of a 100D x 50W rectangle using a mesh size of 5:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.RectangularSection(d=100, b=50)
mesh = geometry.create_mesh(mesh_sizes=[5])
section = CrossSection(geometry, mesh)

The following example creates a 100D x 50W rectangle, with the top half of the section comprised of timber and the bottom half steel. The timber section is meshed with a maximum area of 10 and the steel section mesh with a maximum area of 5:

import sectionproperties.pre.sections as sections
from sectionproperties.pre.pre import Material
from sectionproperties.analysis.cross_section import CrossSection

geom_steel = sections.RectangularSection(d=50, b=50)
geom_timber = sections.RectangularSection(d=50, b=50, shift=[0, 50])
geometry = sections.MergedSection([geom_steel, geom_timber])
geometry.clean_geometry()

mesh = geometry.create_mesh(mesh_sizes=[5, 10])

steel = Material(
    name='Steel', elastic_modulus=200e3, poissons_ratio=0.3, yield_strength=250,
    color='grey'
)
timber = Material(
    name='Timber', elastic_modulus=8e3, poissons_ratio=0.35, yield_strength=20,
    color='burlywood'
)

section = CrossSection(geometry, mesh, [steel, timber])
section.plot_mesh(materials=True, alpha=0.5)
Variables:
  • elements (list[Tri6]) – List of finite element objects describing the cross-section mesh
  • num_nodes (int) – Number of nodes in the finite element mesh
  • geometry (Geometry) – Cross-section geometry object used to generate the mesh
  • mesh (meshpy.triangle.MeshInfo) – Mesh object returned by meshpy
  • mesh_nodes (numpy.ndarray) – Array of node coordinates from the mesh
  • mesh_elements (numpy.ndarray) – Array of connectivities from the mesh
  • mesh_attributes (numpy.ndarray) – Array of attributes from the mesh
  • materials – List of materials
  • material_groups – List of objects containing the elements in each defined material
  • section_props (SectionProperties) – Class to store calculated section properties
Raises:

AssertionError – If the number of materials does not equal the number of regions

assemble_torsion(lg=True)[source]

Assembles stiffness matrices to be used for the computation of warping properties and the torsion load vector (f_torsion). Both a regular (k) and Lagrangian multiplier (k_lg) stiffness matrix are returned. The stiffness matrices are assembled using the sparse COO format and returned in the sparse CSC format.

Parameters:lg (bool) – Whether or not to calculate the Lagrangian multiplier stiffness matrix
Returns:Regular stiffness matrix, Lagrangian multiplier stiffness matrix and torsion load vector (k, k_lg, f_torsion)
Return type:tuple(scipy.sparse.csc_matrix, scipy.sparse.csc_matrix, numpy.ndarray)
calculate_frame_properties(time_info=False, solver_type='direct')[source]

Calculates and returns the properties required for a frame analysis. The properties are also stored in the SectionProperties object contained in the section_props class variable.

Parameters:
  • time_info (bool) – If set to True, a detailed description of the computation and the time cost is printed to the terminal.
  • solver_type (string) – Solver used for solving systems of linear equations, either using the ‘direct’ method or ‘cgs’ iterative method
Returns:

Cross-section properties to be used for a frame analysis (area, ixx, iyy, ixy, j, phi)

Return type:

tuple(float, float, float, float, float, float)

The following section properties are calculated:

  • Cross-sectional area (area)
  • Second moments of area about the centroidal axis (ixx, iyy, ixy)
  • Torsion constant (j)
  • Principal axis angle (phi)

If materials are specified for the cross-section, the area, second moments of area and torsion constant are elastic moulus weighted.

The following example demonstrates the use of this method:

section = CrossSection(geometry, mesh)
(area, ixx, iyy, ixy, j, phi) = section.calculate_frame_properties()
calculate_geometric_properties(time_info=False)[source]

Calculates the geometric properties of the cross-section and stores them in the SectionProperties object contained in the section_props class variable.

Parameters:time_info (bool) – If set to True, a detailed description of the computation and the time cost is printed to the terminal.

The following geometric section properties are calculated:

  • Cross-sectional area
  • Cross-sectional perimeter
  • Modulus weighted area (axial rigidity)
  • First moments of area
  • Second moments of area about the global axis
  • Second moments of area about the centroidal axis
  • Elastic centroid
  • Centroidal section moduli
  • Radii of gyration
  • Principal axis properties

If materials are specified for the cross-section, the moments of area and section moduli are elastic modulus weighted.

The following example demonstrates the use of this method:

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
calculate_plastic_properties(time_info=False, verbose=False, debug=False)[source]

Calculates the plastic properties of the cross-section and stores the, in the SectionProperties object contained in the section_props class variable.

Parameters:
  • time_info (bool) – If set to True, a detailed description of the computation and the time cost is printed to the terminal.
  • verbose (bool) – If set to True, the number of iterations required for each plastic axis is printed to the terminal.
  • debug (bool) – If set to True, the geometry is plotted each time a new mesh is generated by the plastic centroid algorithm.

The following warping section properties are calculated:

  • Plastic centroid for bending about the centroidal and principal axes
  • Plastic section moduli for bending about the centroidal and principal axes
  • Shape factors for bending about the centroidal and principal axes

If materials are specified for the cross-section, the plastic section moduli are displayed as plastic moments (i.e \(M_p = f_y S\)) and the shape factors are not calculated.

Note that the geometric properties must be calculated before the plastic properties are calculated:

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
Raises:RuntimeError – If the geometric properties have not been calculated prior to calling this method
calculate_stress(N=0, Vx=0, Vy=0, Mxx=0, Myy=0, M11=0, M22=0, Mzz=0, time_info=False)[source]

Calculates the cross-section stress resulting from design actions and returns a StressPost object allowing post-processing of the stress results.

Parameters:
  • N (float) – Axial force
  • Vx (float) – Shear force acting in the x-direction
  • Vy (float) – Shear force acting in the y-direction
  • Mxx (float) – Bending moment about the centroidal xx-axis
  • Myy (float) – Bending moment about the centroidal yy-axis
  • M11 (float) – Bending moment about the centroidal 11-axis
  • M22 (float) – Bending moment about the centroidal 22-axis
  • Mzz (float) – Torsion moment about the centroidal zz-axis
  • time_info (bool) – If set to True, a detailed description of the computation and the time cost is printed to the terminal.
Returns:

Object for post-processing cross-section stresses

Return type:

StressPost

Note that a geometric and warping analysis must be performed before a stress analysis is carried out:

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(N=1e3, Vy=3e3, Mxx=1e6)
Raises:RuntimeError – If a geometric and warping analysis have not been performed prior to calling this method
calculate_warping_properties(time_info=False, solver_type='direct')[source]

Calculates all the warping properties of the cross-section and stores them in the SectionProperties object contained in the section_props class variable.

Parameters:
  • time_info (bool) – If set to True, a detailed description of the computation and the time cost is printed to the terminal.
  • solver_type (string) – Solver used for solving systems of linear equations, either using the ‘direct’ method or ‘cgs’ iterative method

The following warping section properties are calculated:

  • Torsion constant
  • Shear centre
  • Shear area
  • Warping constant
  • Monosymmetry constant

If materials are specified, the values calculated for the torsion constant, warping constant and shear area are elastic modulus weighted.

Note that the geometric properties must be calculated first for the calculation of the warping properties to be correct:

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
Raises:RuntimeError – If the geometric properties have not been calculated prior to calling this method
display_mesh_info()[source]

Prints mesh statistics (number of nodes, elements and regions) to the command window.

The following example displays the mesh statistics for a Tee section merged from two rectangles:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

rec1 = sections.RectangularSection(d=100, b=25, shift=[-12.5, 0])
rec2 = sections.RectangularSection(d=25, b=100, shift=[-50, 100])
geometry = sections.MergedSection([rec1, rec2])
mesh = geometry.create_mesh(mesh_sizes=[5, 2.5])
section = CrossSection(geometry, mesh)
section.display_mesh_info()

>>>Mesh Statistics:
>>>--4920 nodes
>>>--2365 elements
>>>--2 regions
display_results(fmt='8.6e')[source]

Prints the results that have been calculated to the terminal.

Parameters:fmt (string) – Number formatting string

The following example displays the geometric section properties for a 100D x 50W rectangle with three digits after the decimal point:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.RectangularSection(d=100, b=50)
mesh = geometry.create_mesh(mesh_sizes=[5])

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()

section.display_results(fmt='.3f')
get_As()[source]
Returns:Shear area for loading about the centroidal axis (A_sx, A_sy)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(A_sx, A_sy) = section.get_As()
get_As_p()[source]
Returns:Shear area for loading about the principal bending axis (A_s11, A_s22)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(A_s11, A_s22) = section.get_As_p()
get_area()[source]
Returns:Cross-section area
Return type:float
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
area = section.get_area()
get_beta()[source]
Returns:Monosymmetry constant for bending about both global axes (beta_x_plus, beta_x_minus, beta_y_plus, beta_y_minus). The plus value relates to the top flange in compression and the minus value relates to the bottom flange in compression.
Return type:tuple(float, float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(beta_x_plus, beta_x_minus, beta_y_plus, beta_y_minus) = section.get_beta()
get_beta_p()[source]
Returns:Monosymmetry constant for bending about both principal axes (beta_11_plus, beta_11_minus, beta_22_plus, beta_22_minus). The plus value relates to the top flange in compression and the minus value relates to the bottom flange in compression.
Return type:tuple(float, float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(beta_11_plus, beta_11_minus, beta_22_plus, beta_22_minus) = section.get_beta_p()
get_c()[source]
Returns:Elastic centroid (cx, cy)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(cx, cy) = section.get_c()
get_ea()[source]
Returns:Modulus weighted area (axial rigidity)
Return type:float
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
ea = section.get_ea()
get_gamma()[source]
Returns:Warping constant
Return type:float
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
gamma = section.get_gamma()
get_ic()[source]
Returns:Second moments of area centroidal axis (ixx_c, iyy_c, ixy_c)
Return type:tuple(float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(ixx_c, iyy_c, ixy_c) = section.get_ic()
get_ig()[source]
Returns:Second moments of area about the global axis (ixx_g, iyy_g, ixy_g)
Return type:tuple(float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(ixx_g, iyy_g, ixy_g) = section.get_ig()
get_ip()[source]
Returns:Second moments of area about the principal axis (i11_c, i22_c)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(i11_c, i22_c) = section.get_ip()
get_j()[source]
Returns:St. Venant torsion constant
Return type:float
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
j = section.get_j()
get_pc()[source]
Returns:Centroidal axis plastic centroid (x_pc, y_pc)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
(x_pc, y_pc) = section.get_pc()
get_pc_p()[source]
Returns:Principal bending axis plastic centroid (x11_pc, y22_pc)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
(x11_pc, y22_pc) = section.get_pc_p()
get_perimeter()[source]
Returns:Cross-section perimeter
Return type:float
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
perimeter = section.get_perimeter()
get_phi()[source]
Returns:Principal bending axis angle
Return type:float
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
phi = section.get_phi()
get_q()[source]
Returns:First moments of area about the global axis (qx, qy)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(qx, qy) = section.get_q()
get_rc()[source]
Returns:Radii of gyration about the centroidal axis (rx, ry)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(rx, ry) = section.get_rc()
get_rp()[source]
Returns:Radii of gyration about the principal axis (r11, r22)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(r11, r22) = section.get_rp()
get_s()[source]
Returns:Plastic section moduli about the centroidal axis (sxx, syy)
Return type:tuple(float, float)

If material properties have been specified, returns the plastic moment \(M_p = f_y S\).

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
(sxx, syy) = section.get_s()
get_sc()[source]
Returns:Centroidal axis shear centre (elasticity approach) (x_se, y_se)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(x_se, y_se) = section.get_sc()
get_sc_p()[source]
Returns:Principal axis shear centre (elasticity approach) (x11_se, y22_se)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(x11_se, y22_se) = section.get_sc_p()
get_sc_t()[source]
Returns:Centroidal axis shear centre (Trefftz’s approach) (x_st, y_st)
Return type:tuple(float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
(x_st, y_st) = section.get_sc_t()
get_sf()[source]
Returns:Centroidal axis shape factors with respect to the top and bottom fibres (sf_xx_plus, sf_xx_minus, sf_yy_plus, sf_yy_minus)
Return type:tuple(float, float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
(sf_xx_plus, sf_xx_minus, sf_yy_plus, sf_yy_minus) = section.get_sf()
get_sf_p()[source]
Returns:Principal bending axis shape factors with respect to the top and bottom fibres (sf_11_plus, sf_11_minus, sf_22_plus, sf_22_minus)
Return type:tuple(float, float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
(sf_11_plus, sf_11_minus, sf_22_plus, sf_22_minus) = section.get_sf_p()
get_sp()[source]
Returns:Plastic section moduli about the principal bending axis (s11, s22)
Return type:tuple(float, float)

If material properties have been specified, returns the plastic moment \(M_p = f_y S\).

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_plastic_properties()
(s11, s22) = section.get_sp()
get_z()[source]
Returns:Elastic section moduli about the centroidal axis with respect to the top and bottom fibres (zxx_plus, zxx_minus, zyy_plus, zyy_minus)
Return type:tuple(float, float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(zxx_plus, zxx_minus, zyy_plus, zyy_minus) = section.get_z()
get_zp()[source]
Returns:Elastic section moduli about the principal axis with respect to the top and bottom fibres (z11_plus, z11_minus, z22_plus, z22_minus)
Return type:tuple(float, float, float, float)
section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
(z11_plus, z11_minus, z22_plus, z22_minus) = section.get_zp()
plot_centroids(pause=True)[source]

Plots the elastic centroid, the shear centre, the plastic centroids and the principal axis, if they have been calculated, on top of the finite element mesh.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example analyses a 200 PFC section and displays a plot of the centroids:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.PfcSection(d=200, b=75, t_f=12, t_w=6, r=12, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
section.calculate_plastic_properties()

section.plot_centroids()
../_images/pfc_centroids.png

Plot of the centroids generated by the above example.

The following example analyses a 150x90x12 UA section and displays a plot of the centroids:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])

section = CrossSection(geometry, mesh)
section.calculate_geometric_properties()
section.calculate_warping_properties()
section.calculate_plastic_properties()

section.plot_centroids()
../_images/angle_centroids.png

Plot of the centroids generated by the above example.

plot_mesh(ax=None, pause=True, alpha=1, materials=False, mask=None)[source]

Plots the finite element mesh. If no axes object is supplied a new figure and axis is created.

Parameters:
  • ax (matplotlib.axes.Axes) – Axes object on which the mesh is plotted
  • pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
  • alpha (float) – Transparency of the mesh outlines: \(0 \leq \alpha \leq 1\)
  • materials (bool) – If set to true and material properties have been provided to the CrossSection object, shades the elements with the specified material colours
  • mask (list[bool]) – Mask array, of length num_nodes, to mask out triangles
Returns:

Matplotlib figure and axes objects (fig, ax)

Return type:

(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the mesh generated for the second example listed under the CrossSection object definition:

import sectionproperties.pre.sections as sections
from sectionproperties.pre.pre import Material
from sectionproperties.analysis.cross_section import CrossSection

geom_steel = sections.RectangularSection(d=50, b=50)
geom_timber = sections.RectangularSection(d=50, b=50, shift=[50, 0])
geometry = sections.MergedSection([geom_steel, geom_timber])
geometry.clean_geometry()

mesh = geometry.create_mesh(mesh_sizes=[5, 10])

steel = Material(
    name='Steel', elastic_modulus=200e3, poissons_ratio=0.3, yield_strength=250,
    color='grey'
)
timber = Material(
    name='Timber', elastic_modulus=8e3, poissons_ratio=0.35, yield_strength=20,
    color='burlywood'
)

section = CrossSection(geometry, mesh, [steel, timber])
section.plot_mesh(materials=True, alpha=0.5)
../_images/composite_mesh.png

Finite element mesh generated by the above example.

PlasticSection Class

class sectionproperties.analysis.cross_section.PlasticSection(geometry, materials, debug)[source]

Bases: object

Class for the plastic analysis of cross-sections.

Stores the finite element geometry and material information and provides methods to compute the plastic section properties.

Parameters:
  • geometry (Geometry) – Cross-section geometry object
  • materials (list[Material]) – A list of material properties corresponding to various regions in the geometry and mesh.
  • debug (bool) – If set to True, the geometry is plotted each time a new mesh is generated by the plastic centroid algorithm.
Variables:
  • geometry (Geometry) – Deep copy of the cross-section geometry object provided to the constructor
  • materials (list[Material]) – A list of material properties corresponding to various regions in the geometry and mesh.
  • debug (bool) – If set to True, the geometry is plotted each time a new mesh is generated by the plastic centroid algorithm.
  • mesh (meshpy.triangle.MeshInfo) – Mesh object returned by meshpy
  • mesh_nodes (numpy.ndarray) – Array of node coordinates from the mesh
  • mesh_elements (numpy.ndarray) – Array of connectivities from the mesh
  • elements (list[Tri6]) – List of finite element objects describing the cross-section mesh
  • f_top (float) – Current force in the top region
  • c_top – Centroid of the force in the top region (c_top_x, c_top_y)
  • c_bot – Centroid of the force in the bottom region (c_bot_x, c_bot_y)
add_line(geometry, line)[source]

Adds a line a geometry object. Finds the intersection points of the line with the current facets and splits the existing facets to accomodate the new line.

Parameters:
  • geometry (Geometry) – Cross-section geometry object used to generate the mesh
  • line (list[numpy.ndarray, numpy.ndarray]) – A point p and a unit vector u defining a line to add to the mesh (line: p -> p + u)
calculate_centroid(elements)[source]

Calculates the elastic centroid from a list of finite elements.

Parameters:elements (list[Tri6]) – A list of Tri6 finite elements.
Returns:A tuple containing the x and y location of the elastic centroid.
Return type:tuple(float, float)
calculate_extreme_fibres(angle)[source]

Calculates the locations of the extreme fibres along and perpendicular to the axis specified by ‘angle’ using the elements stored in self.elements.

Parameters:angle (float) – Angle (in radians) along which to calculate the extreme fibre locations
Returns:The location of the extreme fibres parallel (u) and perpendicular (v) to the axis (u_min, u_max, v_min, v_max)
Return type:tuple(float, float, float, float)
calculate_plastic_force(elements, u, p)[source]

Sums the forces above and below the axis defined by unit vector u and point p. Also returns the force centroid of the forces above and below the axis.

Parameters:
  • elements (list[Tri6]) – A list of Tri6 finite elements.
  • u (numpy.ndarray) – Unit vector defining the direction of the axis
  • p (numpy.ndarray) – Point on the axis
Returns:

Force in the top and bottom areas (f_top, f_bot)

Return type:

tuple(float, float)

calculate_plastic_properties(cross_section, verbose)[source]

Calculates the location of the plastic centroid with respect to the centroidal and principal bending axes, the plastic section moduli and shape factors and stores the results to the supplied CrossSection object.

Parameters:
  • cross_section (CrossSection) – Cross section object that uses the same geometry and materials specified in the class constructor
  • verbose (bool) – If set to True, the number of iterations required for each plastic axis is printed to the terminal.
check_convergence(root_result, axis)[source]

Checks that the function solver converged and if not, raises a helpful error.

Parameters:
  • root_result (scipy.optimize.RootResults) – Result object from the root finder
  • axis (string) – Axis being considered by the function sovler
Raises:

RuntimeError – If the function solver did not converge

create_plastic_mesh(new_line=None)[source]

Generates a triangular mesh of a deep copy of the geometry stored in self.geometry. Optionally, a line can be added to the copied geometry, which is defined by a point p and a unit vector u.

Parameters:
  • new_line (list[numpy.ndarray, numpy.ndarray]) – A point p and a unit vector u defining a line to add to the mesh (new_line: p -> p + u) [p, u]
  • mesh (meshpy.triangle.MeshInfo) – Mesh object returned by meshpy
evaluate_force_eq(d, u, u_p, verbose)[source]

Given a distance d from the centroid to an axis (defined by unit vector u), creates a mesh including the new and axis and calculates the force equilibrium. The resultant force, as a ratio of the total force, is returned.

Parameters:
  • d (float) – Distance from the centroid to current axis
  • u (numpy.ndarray) – Unit vector defining the direction of the axis
  • u_p (numpy.ndarray) – Unit vector perpendicular to the direction of the axis
  • verbose (bool) – If set to True, the number of iterations required for each plastic axis is printed to the terminal.
Returns:

The force equilibrium norm

Return type:

float

get_elements(mesh)[source]

Extracts finite elements from the provided mesh and returns Tri6 finite elements with their associated material properties.

Parameters:mesh (meshpy.triangle.MeshInfo) – Mesh object returned by meshpy
Returns:A tuple containing an array of the nodes locations, element indicies and a list of the finite elements.
Return type:tuple(numpy.ndarray, numpy.ndarray, list[Tri6])
pc_algorithm(u, dlim, axis, verbose)[source]

An algorithm used for solving for the location of the plastic centroid. The algorithm searches for the location of the axis, defined by unit vector u and within the section depth, that satisfies force equilibrium.

Parameters:
  • u (numpy.ndarray) – Unit vector defining the direction of the axis
  • dlim (list[float, float]) – List [dmax, dmin] containing the distances from the centroid to the extreme fibres perpendicular to the axis
  • axis (int) – The current axis direction: 1 (e.g. x or 11) or 2 (e.g. y or 22)
  • verbose (bool) – If set to True, the number of iterations required for each plastic axis is printed to the terminal.
Returns:

The distance to the plastic centroid axis d, the result object r, the force in the top of the section f_top and the location of the centroids of the top and bottom areas c_top and c_bottom

Return type:

tuple(float, scipy.optimize.RootResults, float, list[float, float], list[float, float])

plot_mesh(nodes, elements, element_list, materials)[source]

Watered down implementation of the CrossSection method to plot the finite element mesh, showing material properties.

point_within_element(pt)[source]

Determines whether a point lies within an element in the mesh stored in self.mesh_elements.

Parameters:pt (numpy.ndarray) – Point to check
Returns:Whether the point lies within an element
Return type:bool
print_verbose(d, root_result, axis)[source]

Prints information related to the function solver convergence to the terminal.

Parameters:
  • d (float) – Location of the plastic centroid axis
  • root_result (scipy.optimize.RootResults) – Result object from the root finder
  • axis (string) – Axis being considered by the function sovler
rebuild_parent_facet(geometry, fct_idx, pt_idx)[source]

Splits and rebuilds a facet at a given point.

Parameters:
  • geometry (Geometry) – Cross-section geometry object used to generate the mesh
  • fct_idx (int) – Index of the facet to be split
  • pt_idx (int) – Index of the point to insert into the facet

StressPost Class

class sectionproperties.analysis.cross_section.StressPost(cross_section)[source]

Bases: object

Class for post-processing finite element stress results.

A StressPost object is created when a stress analysis is carried out and is returned as an object to allow post-processing of the results. The StressPost object creates a deep copy of the MaterialGroups within the cross-section to allow the calculation of stresses for each material. Methods for post-processing the calculated stresses are provided.

Parameters:

cross_section (CrossSection) – Cross section object for stress calculation

Variables:
  • cross_section (CrossSection) – Cross section object for stress calculation
  • material_groups (list[MaterialGroup]) – A deep copy of the cross_section material groups to allow a new stress analysis
get_stress()[source]

Returns the stresses within each material belonging to the current StressPost object.

Returns:A list of dictionaries containing the cross-section stresses for each material.
Return type:list[dict]

A dictionary is returned for each material in the cross-section, containing the following keys and values:

  • ‘Material’: Material name
  • ‘sig_zz_n’: Normal stress \(\sigma_{zz,N}\) resulting from the axial load \(N\)
  • ‘sig_zz_mxx’: Normal stress \(\sigma_{zz,Mxx}\) resulting from the bending moment \(M_{xx}\)
  • ‘sig_zz_myy’: Normal stress \(\sigma_{zz,Myy}\) resulting from the bending moment \(M_{yy}\)
  • ‘sig_zz_m11’: Normal stress \(\sigma_{zz,M11}\) resulting from the bending moment \(M_{11}\)
  • ‘sig_zz_m22’: Normal stress \(\sigma_{zz,M22}\) resulting from the bending moment \(M_{22}\)
  • ‘sig_zz_m’: Normal stress \(\sigma_{zz,\Sigma M}\) resulting from all bending moments
  • ‘sig_zx_mzz’: x-component of the shear stress \(\sigma_{zx,Mzz}\) resulting from the torsion moment
  • ‘sig_zy_mzz’: y-component of the shear stress \(\sigma_{zy,Mzz}\) resulting from the torsion moment
  • ‘sig_zxy_mzz’: Resultant shear stress \(\sigma_{zxy,Mzz}\) resulting from the torsion moment
  • ‘sig_zx_vx’: x-component of the shear stress \(\sigma_{zx,Vx}\) resulting from the shear force \(V_{x}\)
  • ‘sig_zy_vx’: y-component of the shear stress \(\sigma_{zy,Vx}\) resulting from the shear force \(V_{x}\)
  • ‘sig_zxy_vx’: Resultant shear stress \(\sigma_{zxy,Vx}\) resulting from the shear force \(V_{x}\)
  • ‘sig_zx_vy’: x-component of the shear stress \(\sigma_{zx,Vy}\) resulting from the shear force \(V_{y}\)
  • ‘sig_zy_vy’: y-component of the shear stress \(\sigma_{zy,Vy}\) resulting from the shear force \(V_{y}\)
  • ‘sig_zxy_vy’: Resultant shear stress \(\sigma_{zxy,Vy}\) resulting from the shear force \(V_{y}\)
  • ‘sig_zx_v’: x-component of the shear stress \(\sigma_{zx,\Sigma V}\) resulting from all shear forces
  • ‘sig_zy_v’: y-component of the shear stress \(\sigma_{zy,\Sigma V}\) resulting from all shear forces
  • ‘sig_zxy_v’: Resultant shear stress \(\sigma_{zxy,\Sigma V}\) resulting from all shear forces
  • ‘sig_zz’: Combined normal stress \(\sigma_{zz}\) resulting from all actions
  • ‘sig_zx’: x-component of the shear stress \(\sigma_{zx}\) resulting from all actions
  • ‘sig_zy’: y-component of the shear stress \(\sigma_{zy}\) resulting from all actions
  • ‘sig_zxy’: Resultant shear stress \(\sigma_{zxy}\) resulting from all actions
  • ‘sig_vm’: von Mises stress \(\sigma_{vM}\) resulting from all actions

The following example returns the normal stress within a 150x90x12 UA section resulting from an axial force of 10 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(N=10e3)

stresses = stress_post.get_stress()
print('Material: {0}'.format(stresses[0]['Material']))
print('Axial Stresses: {0}'.format(stresses[0]['sig_zz_n']))

$ Material: default
$ Axial Stresses: [3.6402569 3.6402569 3.6402569 ... 3.6402569 3.6402569 3.6402569]
plot_stress_contour(sigs, title, pause)[source]

Plots filled stress contours over the finite element mesh.

Parameters:
  • sigs (list[numpy.ndarray]) – List of nodal stress values for each material
  • title (string) – Plot title
  • pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:

Matplotlib figure and axes objects (fig, ax)

Return type:

(matplotlib.figure.Figure, matplotlib.axes)

plot_stress_m11_zz(pause=True)[source]

Produces a contour plot of the normal stress \(\sigma_{zz,M11}\) resulting from the bending moment \(M_{11}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from a bending moment about the 11-axis of 5 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(M11=5e6)

stress_post.plot_stress_m11_zz()
../_images/stress_m11_zz.png

Contour plot of the bending stress.

plot_stress_m22_zz(pause=True)[source]

Produces a contour plot of the normal stress \(\sigma_{zz,M22}\) resulting from the bending moment \(M_{22}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from a bending moment about the 22-axis of 2 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(M22=5e6)

stress_post.plot_stress_m22_zz()
../_images/stress_m22_zz.png

Contour plot of the bending stress.

plot_stress_m_zz(pause=True)[source]

Produces a contour plot of the normal stress \(\sigma_{zz,\Sigma M}\) resulting from all bending moments \(M_{xx} + M_{yy} + M_{11} + M_{22}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from a bending moment about the x-axis of 5 kN.m, a bending moment about the y-axis of 2 kN.m and a bending moment of 3 kN.m about the 11-axis:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mxx=5e6, Myy=2e6, M11=3e6)

stress_post.plot_stress_m_zz()
../_images/stress_m_zz.png

Contour plot of the bending stress.

plot_stress_mxx_zz(pause=True)[source]

Produces a contour plot of the normal stress \(\sigma_{zz,Mxx}\) resulting from the bending moment \(M_{xx}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from a bending moment about the x-axis of 5 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mxx=5e6)

stress_post.plot_stress_mxx_zz()
../_images/stress_mxx_zz.png

Contour plot of the bending stress.

plot_stress_myy_zz(pause=True)[source]

Produces a contour plot of the normal stress \(\sigma_{zz,Myy}\) resulting from the bending moment \(M_{yy}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from a bending moment about the y-axis of 2 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Myy=2e6)

stress_post.plot_stress_myy_zz()
../_images/stress_myy_zz.png

Contour plot of the bending stress.

plot_stress_mzz_zx(pause=True)[source]

Produces a contour plot of the x-component of the shear stress \(\sigma_{zx,Mzz}\) resulting from the torsion moment \(M_{zz}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the x-component of the shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6)

stress_post.plot_stress_mzz_zx()
../_images/stress_mzz_zx.png

Contour plot of the shear stress.

plot_stress_mzz_zxy(pause=True)[source]

Produces a contour plot of the resultant shear stress \(\sigma_{zxy,Mzz}\) resulting from the torsion moment \(M_{zz}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots a contour of the resultant shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6)

stress_post.plot_stress_mzz_zxy()
../_images/stress_mzz_zxy.png

Contour plot of the shear stress.

plot_stress_mzz_zy(pause=True)[source]

Produces a contour plot of the y-component of the shear stress \(\sigma_{zy,Mzz}\) resulting from the torsion moment \(M_{zz}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the y-component of the shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6)

stress_post.plot_stress_mzz_zy()
../_images/stress_mzz_zy.png

Contour plot of the shear stress.

plot_stress_n_zz(pause=True)[source]

Produces a contour plot of the normal stress \(\sigma_{zz,N}\) resulting from the axial load \(N\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from an axial force of 10 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(N=10e3)

stress_post.plot_stress_n_zz()
../_images/stress_n_zz.png

Contour plot of the axial stress.

plot_stress_v_zx(pause=True)[source]

Produces a contour plot of the x-component of the shear stress \(\sigma_{zx,\Sigma V}\) resulting from the sum of the applied shear forces \(V_{x} + V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the x-component of the shear stress within a 150x90x12 UA section resulting from a shear force of 15 kN in the x-direction and 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3, Vy=30e3)

stress_post.plot_stress_v_zx()
../_images/stress_v_zx.png

Contour plot of the shear stress.

plot_stress_v_zxy(pause=True)[source]

Produces a contour plot of the resultant shear stress \(\sigma_{zxy,\Sigma V}\) resulting from the sum of the applied shear forces \(V_{x} + V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots a contour of the resultant shear stress within a 150x90x12 UA section resulting from a shear force of 15 kN in the x-direction and 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3, Vy=30e3)

stress_post.plot_stress_v_zxy()
../_images/stress_v_zxy.png

Contour plot of the shear stress.

plot_stress_v_zy(pause=True)[source]

Produces a contour plot of the y-component of the shear stress \(\sigma_{zy,\Sigma V}\) resulting from the sum of the applied shear forces \(V_{x} + V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the y-component of the shear stress within a 150x90x12 UA section resulting from a shear force of 15 kN in the x-direction and 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3, Vy=30e3)

stress_post.plot_stress_v_zy()
../_images/stress_v_zy.png

Contour plot of the shear stress.

plot_stress_vector(sigxs, sigys, title, pause)[source]

Plots stress vectors over the finite element mesh.

Parameters:
  • sigxs (list[numpy.ndarray]) – List of x-components of the nodal stress values for each material
  • sigys (list[numpy.ndarray]) – List of y-components of the nodal stress values for each material
  • title (string) – Plot title
  • pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:

Matplotlib figure and axes objects (fig, ax)

Return type:

(matplotlib.figure.Figure, matplotlib.axes)

plot_stress_vm(pause=True)[source]

Produces a contour plot of the von Mises stress \(\sigma_{vM}\) resulting from all actions.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots a contour of the von Mises stress within a 150x90x12 UA section resulting from the following actions:

  • \(N = 50\) kN
  • \(M_{xx} = -5\) kN.m
  • \(M_{22} = 2.5\) kN.m
  • \(M_{zz} = 1.5\) kN.m
  • \(V_{x} = 10\) kN
  • \(V_{y} = 5\) kN
import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(
    N=50e3, Mxx=-5e6, M22=2.5e6, Mzz=0.5e6, Vx=10e3, Vy=5e3
)

stress_post.plot_stress_vm()
../_images/stress_vm.png

Contour plot of the von Mises stress.

plot_stress_vx_zx(pause=True)[source]

Produces a contour plot of the x-component of the shear stress \(\sigma_{zx,Vx}\) resulting from the shear force \(V_{x}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the x-component of the shear stress within a 150x90x12 UA section resulting from a shear force in the x-direction of 15 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3)

stress_post.plot_stress_vx_zx()
../_images/stress_vx_zx.png

Contour plot of the shear stress.

plot_stress_vx_zxy(pause=True)[source]

Produces a contour plot of the resultant shear stress \(\sigma_{zxy,Vx}\) resulting from the shear force \(V_{x}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots a contour of the resultant shear stress within a 150x90x12 UA section resulting from a shear force in the x-direction of 15 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3)

stress_post.plot_stress_vx_zxy()
../_images/stress_vx_zxy.png

Contour plot of the shear stress.

plot_stress_vx_zy(pause=True)[source]

Produces a contour plot of the y-component of the shear stress \(\sigma_{zy,Vx}\) resulting from the shear force \(V_{x}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the y-component of the shear stress within a 150x90x12 UA section resulting from a shear force in the x-direction of 15 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3)

stress_post.plot_stress_vx_zy()
../_images/stress_vx_zy.png

Contour plot of the shear stress.

plot_stress_vy_zx(pause=True)[source]

Produces a contour plot of the x-component of the shear stress \(\sigma_{zx,Vy}\) resulting from the shear force \(V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the x-component of the shear stress within a 150x90x12 UA section resulting from a shear force in the y-direction of 30 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vy=30e3)

stress_post.plot_stress_vy_zx()
../_images/stress_vy_zx.png

Contour plot of the shear stress.

plot_stress_vy_zxy(pause=True)[source]

Produces a contour plot of the resultant shear stress \(\sigma_{zxy,Vy}\) resulting from the shear force \(V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots a contour of the resultant shear stress within a 150x90x12 UA section resulting from a shear force in the y-direction of 30 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vy=30e3)

stress_post.plot_stress_vy_zxy()
../_images/stress_vy_zxy.png

Contour plot of the shear stress.

plot_stress_vy_zy(pause=True)[source]

Produces a contour plot of the y-component of the shear stress \(\sigma_{zy,Vy}\) resulting from the shear force \(V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the y-component of the shear stress within a 150x90x12 UA section resulting from a shear force in the y-direction of 30 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vy=30e3)

stress_post.plot_stress_vy_zy()
../_images/stress_vy_zy.png

Contour plot of the shear stress.

plot_stress_zx(pause=True)[source]

Produces a contour plot of the x-component of the shear stress \(\sigma_{zx}\) resulting from all actions.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the x-component of the shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m and a shear force of 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6, Vy=30e3)

stress_post.plot_stress_zx()
../_images/stress_zx.png

Contour plot of the shear stress.

plot_stress_zxy(pause=True)[source]

Produces a contour plot of the resultant shear stress \(\sigma_{zxy}\) resulting from all actions.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots a contour of the resultant shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m and a shear force of 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6, Vy=30e3)

stress_post.plot_stress_zxy()
../_images/stress_zxy.png

Contour plot of the shear stress.

plot_stress_zy(pause=True)[source]

Produces a contour plot of the y-component of the shear stress \(\sigma_{zy}\) resulting from all actions.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the y-component of the shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m and a shear force of 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6, Vy=30e3)

stress_post.plot_stress_zy()
../_images/stress_zy.png

Contour plot of the shear stress.

plot_stress_zz(pause=True)[source]

Produces a contour plot of the combined normal stress \(\sigma_{zz}\) resulting from all actions.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example plots the normal stress within a 150x90x12 UA section resulting from an axial force of 100 kN, a bending moment about the x-axis of 5 kN.m and a bending moment about the y-axis of 2 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(N=100e3, Mxx=5e6, Myy=2e6)

stress_post.plot_stress_zz()
../_images/stress_zz.png

Contour plot of the normal stress.

plot_vector_mzz_zxy(pause=True)[source]

Produces a vector plot of the resultant shear stress \(\sigma_{zxy,Mzz}\) resulting from the torsion moment \(M_{zz}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example generates a vector plot of the shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6)

stress_post.plot_vector_mzz_zxy()
../_images/vector_mzz_zxy.png

Vector plot of the shear stress.

plot_vector_v_zxy(pause=True)[source]

Produces a vector plot of the resultant shear stress \(\sigma_{zxy,\Sigma V}\) resulting from the sum of the applied shear forces \(V_{x} + V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example generates a vector plot of the shear stress within a 150x90x12 UA section resulting from a shear force of 15 kN in the x-direction and 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3, Vy=30e3)

stress_post.plot_vector_v_zxy()
../_images/vector_v_zxy.png

Vector plot of the shear stress.

plot_vector_vx_zxy(pause=True)[source]

Produces a vector plot of the resultant shear stress \(\sigma_{zxy,Vx}\) resulting from the shear force \(V_{x}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example generates a vector plot of the shear stress within a 150x90x12 UA section resulting from a shear force in the x-direction of 15 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vx=15e3)

stress_post.plot_vector_vx_zxy()
../_images/vector_vx_zxy.png

Vector plot of the shear stress.

plot_vector_vy_zxy(pause=True)[source]

Produces a vector plot of the resultant shear stress \(\sigma_{zxy,Vy}\) resulting from the shear force \(V_{y}\).

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example generates a vector plot of the shear stress within a 150x90x12 UA section resulting from a shear force in the y-direction of 30 kN:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Vy=30e3)

stress_post.plot_vector_vy_zxy()
../_images/vector_vy_zxy.png

Vector plot of the shear stress.

plot_vector_zxy(pause=True)[source]

Produces a vector plot of the resultant shear stress \(\sigma_{zxy}\) resulting from all actions.

Parameters:pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
Returns:Matplotlib figure and axes objects (fig, ax)
Return type:(matplotlib.figure.Figure, matplotlib.axes)

The following example generates a vector plot of the shear stress within a 150x90x12 UA section resulting from a torsion moment of 1 kN.m and a shear force of 30 kN in the y-direction:

import sectionproperties.pre.sections as sections
from sectionproperties.analysis.cross_section import CrossSection

geometry = sections.AngleSection(d=150, b=90, t=12, r_r=10, r_t=5, n_r=8)
mesh = geometry.create_mesh(mesh_sizes=[2.5])
section = CrossSection(geometry, mesh)

section.calculate_geometric_properties()
section.calculate_warping_properties()
stress_post = section.calculate_stress(Mzz=1e6, Vy=30e3)

stress_post.plot_vector_zxy()
../_images/vector_zxy.png

Vector plot of the shear stress.

MaterialGroup Class

class sectionproperties.analysis.cross_section.MaterialGroup(material, num_nodes)[source]

Bases: object

Class for storing elements of different materials.

A MaterialGroup object contains the finite element objects for a specified material. The stress_result variable provides storage for stresses related each material.

Parameters:
  • material (Material) – Material object for the current MaterialGroup
  • num_nods (int) – Number of nodes for the entire cross-section
Variables:
  • material (Material) – Material object for the current MaterialGroup
  • stress_result (StressResult) – A StressResult object for saving the stresses of the current material
  • elements (list[Tri6]) – A list of finite element objects that are of the current material type
  • el_ids (list[int]) – A list of the element IDs of the elements that are of the current material type
add_element(element)[source]

Adds an element and its element ID to the MaterialGroup.

Parameters:element (Tri6) – Element to add to the MaterialGroup

StressResult Class

class sectionproperties.analysis.cross_section.StressResult(num_nodes)[source]

Bases: object

Class for storing a stress result.

Provides variables to store the results from a cross-section stress analysis. Also provides a method to calculate combined stresses.

Parameters:

num_nodes (int) – Number of nodes in the finite element mesh

Variables:
  • sig_zz_n (numpy.ndarray) – Normal stress (\(\sigma_{zz,N}\)) resulting from an axial force
  • sig_zz_mxx (numpy.ndarray) – Normal stress (\(\sigma_{zz,Mxx}\)) resulting from a bending moment about the xx-axis
  • sig_zz_myy (numpy.ndarray) – Normal stress (\(\sigma_{zz,Myy}\)) resulting from a bending moment about the yy-axis
  • sig_zz_m11 (numpy.ndarray) – Normal stress (\(\sigma_{zz,M11}\)) resulting from a bending moment about the 11-axis
  • sig_zz_m22 (numpy.ndarray) – Normal stress (\(\sigma_{zz,M22}\)) resulting from a bending moment about the 22-axis
  • sig_zx_mzz (numpy.ndarray) – Shear stress (\(\sigma_{zx,Mzz}\)) resulting from a torsion moment about the zz-axis
  • sig_zy_mzz (numpy.ndarray) – Shear stress (\(\sigma_{zy,Mzz}\)) resulting from a torsion moment about the zz-axis
  • sig_zx_vx (numpy.ndarray) – Shear stress (\(\sigma_{zx,Vx}\)) resulting from a shear force in the x-direction
  • sig_zy_vx (numpy.ndarray) – Shear stress (\(\sigma_{zy,Vx}\)) resulting from a shear force in the x-direction
  • sig_zx_vy (numpy.ndarray) – Shear stress (\(\sigma_{zx,Vy}\)) resulting from a shear force in the y-direction
  • sig_zy_vy (numpy.ndarray) – Shear stress (\(\sigma_{zy,Vy}\)) resulting from a shear force in the y-direction
  • sig_zz_m (numpy.ndarray) – Normal stress (\(\sigma_{zz,\Sigma M}\)) resulting from all bending moments
  • sig_zxy_mzz (numpy.ndarray) – Resultant shear stress (\(\sigma_{zxy,Mzz}\)) resulting from a torsion moment in the zz-direction
  • sig_zxy_vx (numpy.ndarray) – Resultant shear stress (\(\sigma_{zxy,Vx}\)) resulting from a a shear force in the x-direction
  • sig_zxy_vy (numpy.ndarray) – Resultant shear stress (\(\sigma_{zxy,Vy}\)) resulting from a a shear force in the y-direction
  • sig_zx_v (numpy.ndarray) – Shear stress (\(\sigma_{zx,\Sigma V}\)) resulting from all shear forces
  • sig_zy_v (numpy.ndarray) – Shear stress (\(\sigma_{zy,\Sigma V}\)) resulting from all shear forces
  • sig_zxy_v (numpy.ndarray) – Resultant shear stress (\(\sigma_{zxy,\Sigma V}\)) resulting from all shear forces
  • sig_zz (numpy.ndarray) – Combined normal force (\(\sigma_{zz}\)) resulting from all actions
  • sig_zx (numpy.ndarray) – Combined shear stress (\(\sigma_{zx}\)) resulting from all actions
  • sig_zy (numpy.ndarray) – Combined shear stress (\(\sigma_{zy}\)) resulting from all actions
  • sig_zxy (numpy.ndarray) – Combined resultant shear stress (\(\sigma_{zxy}\)) resulting from all actions
  • sig_vm (numpy.ndarray) – von Mises stress (\(\sigma_{VM}\)) resulting from all actions
calculate_combined_stresses()[source]

Calculates the combined cross-section stresses.

SectionProperties Class

class sectionproperties.analysis.cross_section.SectionProperties[source]

Bases: object

Class for storing section properties.

Stores calculated section properties. Also provides methods to calculate section properties entirely derived from other section properties.

Variables:
  • area (float) – Cross-sectional area
  • perimeter (float) – Cross-sectional perimeter
  • ea (float) – Modulus weighted area (axial rigidity)
  • ga (float) – Modulus weighted product of shear modulus and area
  • nu_eff (float) – Effective Poisson’s ratio
  • qx (float) – First moment of area about the x-axis
  • qy (float) – First moment of area about the y-axis
  • ixx_g (float) – Second moment of area about the global x-axis
  • iyy_g (float) – Second moment of area about the global y-axis
  • ixy_g (float) – Second moment of area about the global xy-axis
  • cx (float) – X coordinate of the elastic centroid
  • cy (float) – Y coordinate of the elastic centroid
  • ixx_c (float) – Second moment of area about the centroidal x-axis
  • iyy_c (float) – Second moment of area about the centroidal y-axis
  • ixy_c (float) – Second moment of area about the centroidal xy-axis
  • zxx_plus (float) – Section modulus about the centroidal x-axis for stresses at the positive extreme value of y
  • zxx_minus (float) – Section modulus about the centroidal x-axis for stresses at the negative extreme value of y
  • zyy_plus (float) – Section modulus about the centroidal y-axis for stresses at the positive extreme value of x
  • zyy_minus (float) – Section modulus about the centroidal y-axis for stresses at the negative extreme value of x
  • rx_c (float) – Radius of gyration about the centroidal x-axis.
  • ry_c (float) – Radius of gyration about the centroidal y-axis.
  • i11_c (float) – Second moment of area about the centroidal 11-axis
  • i22_c (float) – Second moment of area about the centroidal 22-axis
  • phi (float) – Principal axis angle
  • z11_plus (float) – Section modulus about the principal 11-axis for stresses at the positive extreme value of the 22-axis
  • z11_minus (float) – Section modulus about the principal 11-axis for stresses at the negative extreme value of the 22-axis
  • z22_plus (float) – Section modulus about the principal 22-axis for stresses at the positive extreme value of the 11-axis
  • z22_minus (float) – Section modulus about the principal 22-axis for stresses at the negative extreme value of the 11-axis
  • r11_c (float) – Radius of gyration about the principal 11-axis.
  • r22_c (float) – Radius of gyration about the principal 22-axis.
  • j (float) – Torsion constant
  • omega (numpy.ndarray) – Warping function
  • psi_shear (numpy.ndarray) – Psi shear function
  • phi_shear (numpy.ndarray) – Phi shear function
  • Delta_s (float) – Shear factor
  • x_se (float) – X coordinate of the shear centre (elasticity approach)
  • y_se (float) – Y coordinate of the shear centre (elasticity approach)
  • x11_se (float) – 11 coordinate of the shear centre (elasticity approach)
  • y22_se (float) – 22 coordinate of the shear centre (elasticity approach)
  • x_st (float) – X coordinate of the shear centre (Trefftz’s approach)
  • y_st (float) – Y coordinate of the shear centre (Trefftz’s approach)
  • gamma (float) – Warping constant
  • A_sx (float) – Shear area about the x-axis
  • A_sy (float) – Shear area about the y-axis
  • A_sxy (float) – Shear area about the xy-axis
  • A_s11 (float) – Shear area about the 11 bending axis
  • A_s22 (float) – Shear area about the 22 bending axis
  • beta_x_plus (float) – Monosymmetry constant for bending about the x-axis with the top flange in compression
  • beta_x_minus (float) – Monosymmetry constant for bending about the x-axis with the bottom flange in compression
  • beta_y_plus (float) – Monosymmetry constant for bending about the y-axis with the top flange in compression
  • beta_y_minus (float) – Monosymmetry constant for bending about the y-axis with the bottom flange in compression
  • beta_11_plus (float) – Monosymmetry constant for bending about the 11-axis with the top flange in compression
  • beta_11_minus (float) – Monosymmetry constant for bending about the 11-axis with the bottom flange in compression
  • beta_22_plus (float) – Monosymmetry constant for bending about the 22-axis with the top flange in compression
  • beta_22_minus (float) – Monosymmetry constant for bending about the 22-axis with the bottom flange in compression
  • x_pc (float) – X coordinate of the global plastic centroid
  • y_pc (float) – Y coordinate of the global plastic centroid
  • x11_pc (float) – 11 coordinate of the principal plastic centroid
  • y22_pc (float) – 22 coordinate of the principal plastic centroid
  • sxx (float) – Plastic section modulus about the centroidal x-axis
  • syy (float) – Plastic section modulus about the centroidal y-axis
  • sf_xx_plus (float) – Shape factor for bending about the x-axis with respect to the top fibre
  • sf_xx_minus (float) – Shape factor for bending about the x-axis with respect to the bottom fibre
  • sf_yy_plus (float) – Shape factor for bending about the y-axis with respect to the top fibre
  • sf_yy_minus (float) – Shape factor for bending about the y-axis with respect to the bottom fibre
  • s11 (float) – Plastic section modulus about the 11-axis
  • s22 (float) – Plastic section modulus about the 22-axis
  • sf_11_plus (float) – Shape factor for bending about the 11-axis with respect to the top fibre
  • sf_11_minus (float) – Shape factor for bending about the 11-axis with respect to the bottom fibre
  • sf_22_plus (float) – Shape factor for bending about the 22-axis with respect to the top fibre
  • sf_22_minus (float) – Shape factor for bending about the 22-axis with respect to the bottom fibre
calculate_centroidal_properties(mesh)[source]

Calculates the geometric section properties about the centroidal and principal axes based on the results about the global axis.

calculate_elastic_centroid()[source]

Calculates the elastic centroid based on the cross-section area and first moments of area.

fea Module

Tri6 Class

class sectionproperties.analysis.fea.Tri6(el_id, coords, node_ids, material)[source]

Bases: object

Class for a six noded quadratic triangular element.

Provides methods for the calculation of section properties based on the finite element method.

Parameters:
  • el_id (int) – Unique element id
  • coords (numpy.ndarray) – A 2 x 6 array of the coordinates of the tri-6 nodes. The first three columns relate to the vertices of the triangle and the last three columns correspond to the mid-nodes.
  • node_ids (list[int]) – A list of the global node ids for the current element
  • material (Material) – Material object for the current finite element.
Variables:
  • el_id (int) – Unique element id
  • coords (numpy.ndarray) – A 2 x 6 array of the coordinates of the tri-6 nodes. The first three columns relate to the vertices of the triangle and the last three columns correspond to the mid-nodes.
  • node_ids (list[int]) – A list of the global node ids for the current element
  • material (Material) – Material of the current finite element.
element_stress(N, Mxx, Myy, M11, M22, Mzz, Vx, Vy, ea, cx, cy, ixx, iyy, ixy, i11, i22, phi, j, nu, omega, psi_shear, phi_shear, Delta_s)[source]

Calculates the stress within an element resulting from a specified loading. Also returns the shape function weights.

Parameters:
  • N (float) – Axial force
  • Mxx (float) – Bending moment about the centroidal xx-axis
  • Myy (float) – Bending moment about the centroidal yy-axis
  • M11 (float) – Bending moment about the centroidal 11-axis
  • M22 (float) – Bending moment about the centroidal 22-axis
  • Mzz (float) – Torsion moment about the centroidal zz-axis
  • Vx (float) – Shear force acting in the x-direction
  • Vy (float) – Shear force acting in the y-direction
  • ea (float) – Modulus weighted area
  • cx (float) – x position of the elastic centroid
  • cy (float) – y position of the elastic centroid
  • ixx (float) – Second moment of area about the centroidal x-axis
  • iyy (float) – Second moment of area about the centroidal y-axis
  • ixy (float) – Second moment of area about the centroidal xy-axis
  • i11 (float) – Second moment of area about the principal 11-axis
  • i22 (float) – Second moment of area about the principal 22-axis
  • phi (float) – Principal bending axis angle
  • j (float) – St. Venant torsion constant
  • nu (float) – Effective Poisson’s ratio for the cross-section
  • omega (numpy.ndarray) – Values of the warping function at the element nodes
  • psi_shear (numpy.ndarray) – Values of the psi shear function at the element nodes
  • phi_shear (numpy.ndarray) – Values of the phi shear function at the element nodes
  • Delta_s (float) – Cross-section shear factor
Returns:

Tuple containing element stresses and integration weights (\(\sigma_{zz,n}\), \(\sigma_{zz,mxx}\), \(\sigma_{zz,myy}\), \(\sigma_{zz,m11}\), \(\sigma_{zz,m22}\), \(\sigma_{zx,mzz}\), \(\sigma_{zy,mzz}\), \(\sigma_{zx,vx}\), \(\sigma_{zy,vx}\), \(\sigma_{zx,vy}\), \(\sigma_{zy,vy}\), \(w_i\))

Return type:

tuple(numpy.ndarray, numpy.ndarray, …)

geometric_properties()[source]

Calculates the geometric properties for the current finite element.

Returns:Tuple containing the geometric properties and the elastic and shear moduli of the element: (area, qx, qy, ixx, iyy, ixy, e, g)
Return type:tuple(float)
monosymmetry_integrals(phi)[source]

Calculates the integrals used to evaluate the monosymmetry constant about both global axes and both prinicipal axes.

Parameters:phi (float) – Principal bending axis angle
Returns:Integrals used to evaluate the monosymmetry constants (int_x, int_y, int_11, int_22)
Return type:tuple(float, float, float, float)
plastic_properties(u, p)[source]

Calculates total force resisted by the element when subjected to a stress equal to the yield strength. Also returns the modulus weighted area and first moments of area, and determines whether or not the element is above or below the line defined by the unit vector u and point p.

Parameters:
  • u (numpy.ndarray) – Unit vector in the direction of the line
  • p (numpy.ndarray) – Point on the line
Returns:

Element force (force), modulus weighted area properties (ea, e.qx, e.qy) and whether or not the element is above the line

Return type:

tuple(float, float, float, float, bool)

point_within_element(pt)[source]

Determines whether a point lies within the current element.

Parameters:pt (list[float, float]) – Point to check (x, y)
Returns:Whether the point lies within an element
Return type:bool
shear_coefficients(ixx, iyy, ixy, psi_shear, phi_shear, nu)[source]

Calculates the variables used to determine the shear deformation coefficients.

Parameters:
  • ixx (float) – Second moment of area about the centroidal x-axis
  • iyy (float) – Second moment of area about the centroidal y-axis
  • ixy (float) – Second moment of area about the centroidal xy-axis
  • psi_shear (numpy.ndarray) – Values of the psi shear function at the element nodes
  • phi_shear (numpy.ndarray) – Values of the phi shear function at the element nodes
  • nu (float) – Effective Poisson’s ratio for the cross-section
Returns:

Shear deformation variables (kappa_x, kappa_y, kappa_xy)

Return type:

tuple(float, float, float)

shear_load_vectors(ixx, iyy, ixy, nu)[source]

Calculates the element shear load vectors used to evaluate the shear functions.

Parameters:
  • ixx (float) – Second moment of area about the centroidal x-axis
  • iyy (float) – Second moment of area about the centroidal y-axis
  • ixy (float) – Second moment of area about the centroidal xy-axis
  • nu (float) – Effective Poisson’s ratio for the cross-section
Returns:

Element shear load vector psi (f_psi) and phi (f_phi)

Return type:

tuple(numpy.ndarray, numpy.ndarray)

shear_warping_integrals(ixx, iyy, ixy, omega)[source]

Calculates the element shear centre and warping integrals required for shear analysis of the cross-section.

Parameters:
  • ixx (float) – Second moment of area about the centroidal x-axis
  • iyy (float) – Second moment of area about the centroidal y-axis
  • ixy (float) – Second moment of area about the centroidal xy-axis
  • omega (numpy.ndarray) – Values of the warping function at the element nodes
Returns:

Shear centre integrals about the x and y-axes (sc_xint, sc_yint), warping integrals (q_omega, i_omega, i_xomega, i_yomega)

Return type:

tuple(float, float, float, float, float, float)

torsion_properties()[source]

Calculates the element stiffness matrix used for warping analysis and the torsion load vector.

Returns:Element stiffness matrix (k_el) and element torsion load vector (f_el)
Return type:tuple(numpy.ndarray, numpy.ndarray)

fea Functions

sectionproperties.analysis.fea.gauss_points(n)[source]

Returns the Gaussian weights and locations for n point Gaussian integration of a quadratic triangular element.

Parameters:n (int) – Number of Gauss points (1, 3 or 6)
Returns:An n x 4 matrix consisting of the integration weight and the eta, xi and zeta locations for n Gauss points
Return type:numpy.ndarray
sectionproperties.analysis.fea.shape_function(coords, gauss_point)[source]

Computes shape functions, shape function derivatives and the determinant of the Jacobian matrix for a tri 6 element at a given Gauss point.

Parameters:
  • coords (numpy.ndarray) – Global coordinates of the quadratic triangle vertices [2 x 6]
  • gauss_point (numpy.ndarray) – Gaussian weight and isoparametric location of the Gauss point
Returns:

The value of the shape functions N(i) at the given Gauss point [1 x 6], the derivative of the shape functions in the j-th global direction B(i,j) [2 x 6] and the determinant of the Jacobian matrix j

Return type:

tuple(numpy.ndarray, numpy.ndarray, float)

sectionproperties.analysis.fea.extrapolate_to_nodes(w)[source]

Extrapolates results at six Gauss points to the six noes of a quadratic triangular element.

Parameters:w (numpy.ndarray) – Result at the six Gauss points [1 x 6]
Returns:Extrapolated nodal values at the six nodes [1 x 6]
Return type:numpy.ndarray
sectionproperties.analysis.fea.principal_coordinate(phi, x, y)[source]

Determines the coordinates of the cartesian point (x, y) in the principal axis system given an axis rotation angle phi.

Parameters:
  • phi (float) – Prinicpal bending axis angle (degrees)
  • x (float) – x coordinate in the global axis
  • y (float) – y coordinate in the global axis
Returns:

Principal axis coordinates (x1, y2)

Return type:

tuple(float, float)

sectionproperties.analysis.fea.global_coordinate(phi, x11, y22)[source]

Determines the global coordinates of the principal axis point (x1, y2) given principal axis rotation angle phi.

Parameters:
  • phi (float) – Prinicpal bending axis angle (degrees)
  • x11 (float) – 11 coordinate in the principal axis
  • y22 (float) – 22 coordinate in the principal axis
Returns:

Global axis coordinates (x, y)

Return type:

tuple(float, float)

sectionproperties.analysis.fea.point_above_line(u, px, py, x, y)[source]

Determines whether a point (x, y) is a above or below the line defined by the parallel unit vector u and the point (px, py).

Parameters:
  • u (numpy.ndarray) – Unit vector parallel to the line [1 x 2]
  • px (float) – x coordinate of a point on the line
  • py (float) – y coordinate of a point on the line
  • x (float) – x coordinate of the point to be tested
  • y (float) – y coordinate of the point to be tested
Returns:

This method returns True if the point is above the line or False if the point is below the line

Return type:

bool

solver Module

solver Functions

sectionproperties.analysis.solver.solve_cgs(k, f, m=None, tol=1e-05)[source]

Solves a linear system of equations (Ku = f) using the CGS iterative method.

Parameters:
  • k (scipy.sparse.csc_matrix) – N x N matrix of the linear system
  • f (numpy.ndarray) – N x 1 right hand side of the linear system
  • tol (float) – Tolerance for the solver to acheieve. The algorithm terminates when either the relative or the absolute residual is below tol.
  • m (scipy.linalg.LinearOperator) – Preconditioner for the linear matrix approximating the inverse of k
Returns:

The solution vector to the linear system of equations

Return type:

numpy.ndarray

Raises:

RuntimeError – If the CGS iterative method does not converge

sectionproperties.analysis.solver.solve_cgs_lagrange(k_lg, f, tol=1e-05, m=None)[source]

Solves a linear system of equations (Ku = f) using the CGS iterative method and the Lagrangian multiplier method.

Parameters:
  • k (scipy.sparse.csc_matrix) – (N+1) x (N+1) Lagrangian multiplier matrix of the linear system
  • f (numpy.ndarray) – N x 1 right hand side of the linear system
  • tol (float) – Tolerance for the solver to acheieve. The algorithm terminates when either the relative or the absolute residual is below tol.
  • m (scipy.linalg.LinearOperator) – Preconditioner for the linear matrix approximating the inverse of k
Returns:

The solution vector to the linear system of equations

Return type:

numpy.ndarray

Raises:

RuntimeError – If the CGS iterative method does not converge or the error from the Lagrangian multiplier method exceeds the tolerance

sectionproperties.analysis.solver.solve_direct(k, f)[source]

Solves a linear system of equations (Ku = f) using the direct solver method.

Parameters:
  • k (scipy.sparse.csc_matrix) – N x N matrix of the linear system
  • f (numpy.ndarray) – N x 1 right hand side of the linear system
Returns:

The solution vector to the linear system of equations

Return type:

numpy.ndarray

sectionproperties.analysis.solver.solve_direct_lagrange(k_lg, f)[source]

Solves a linear system of equations (Ku = f) using the direct solver method and the Lagrangian multiplier method.

Parameters:
  • k (scipy.sparse.csc_matrix) – (N+1) x (N+1) Lagrangian multiplier matrix of the linear system
  • f (numpy.ndarray) – N x 1 right hand side of the linear system
Returns:

The solution vector to the linear system of equations

Return type:

numpy.ndarray

Raises:

RuntimeError – If the Lagrangian multiplier method exceeds a tolerance of 1e-5

sectionproperties.analysis.solver.function_timer(text, function, *args)[source]

Displays the message text and returns the time taken for a function, with arguments args, to execute. The value returned by the timed function is also returned.

Parameters:
  • text (string) – Message to display
  • function (function) – Function to time and execute
  • args – Function arguments
Returns:

Value returned from the function

Post-Processor Package

post Module

post Functions

sectionproperties.post.post.setup_plot(ax, pause)[source]

Executes code required to set up a matplotlib figure.

Parameters:
  • ax (matplotlib.axes.Axes) – Axes object on which to plot
  • pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
sectionproperties.post.post.finish_plot(ax, pause, title='')[source]

Executes code required to finish a matplotlib figure.

Parameters:
  • ax (matplotlib.axes.Axes) – Axes object on which to plot
  • pause (bool) – If set to true, the figure pauses the script until the window is closed. If set to false, the script continues immediately after the window is rendered.
  • title (string) – Plot title
sectionproperties.post.post.draw_principal_axis(ax, phi, cx, cy)[source]

Draws the principal axis on a plot.

Parameters:
  • ax (matplotlib.axes.Axes) – Axes object on which to plot
  • phi (float) – Principal axis angle in radians
  • cx (float) – x-location of the centroid
  • cy (float) – y-location of the centroid
sectionproperties.post.post.print_results(cross_section, fmt)[source]

Prints the results that have been calculated to the terminal.

Parameters:
  • cross_section (CrossSection) – Structural cross-section object
  • fmt (string) – Number format