Modules

Module Basic

Module pylabfea.basic introduces basic methods and attributes like calculation of equivalent stresses and strain, conversions of cylindrical to principal stresses and vice versa, and the classes Stress and Strain for efficient operations with these quantities.

uses NumPy, pickle

Version: 4.0 (2021-11-27) Author: Alexander Hartmaier, ICAMS/Ruhr University Bochum, Germany Email: alexander.hartmaier@rub.de distributed under GNU General Public License (GPLv3)

class pylabfea.basic.Strain(sv)[source]

Bases: object

Stores and converts Voigt strain tensors into different formats, calculates principle strain and equivalent strain.

Parameters:

sv (list-like object, must be 1D with length 6) – Voigt-strain components

voigt, v

Strain tensor in Voigt notation

Type:

1d-array (size 6)

tens, t

Strain tensor in matrix notation

Type:

3x3 array

princ, p

Principal strains

Type:

1d-array (size 3)

eeq()[source]

Calculate equivalent strain

Returns:

equiv – Equivalent strain

Return type:

float

inv()[source]

Calculate inverse of strain tensor ignoring zeros.

Returns:

inv

Return type:

(6,) array

class pylabfea.basic.Stress(sv)[source]

Bases: object

Stores and converts Voigt stress tensors into different formats, calculates principle stresses, equivalent stresses and transforms into cylindrical coordinates.

Parameters:

sv (list-like object, must be 1D with length 6) – Voigt-stress components

voigt, v

Stress tensor in Voigt notation

Type:

1d-array (size 6)

tens, t

Stress tensor in matrix notation

Type:

3x3 array

princ, p

Principal stresses

Type:

1d-array (size 3)

hydrostatic, h

Hydrostatic stress component

Type:

float

cyl()[source]

Calculate cylindrical stress tensor

Returns:

cyl – stress in cylindrical form: (J2 eqiv. stress, polar angle, hydrostatic)

Return type:

(3,) array

lode_ang(arg)[source]

Calculate Lode angle: Transforms principal stress space into hydrostatic stress, eqiv. stress, and Lode angle; definition of positive cosine for Lode angle is applied

Parameters:

arg (either float or object of class Material) – if type is float: interpreted as equivalent stress if type is Material: method of that class is used to calculate equivalent stress

Returns:

la – Lode angle

Return type:

float

seq(mat=None)[source]

calculate Hill-type equivalent stress, invokes corresponding method of class Material

Parameters:

mat (object of class Material) – contains Hill parameters and method needed for Hill-type equivalent stress (optional, default=None)

Returns:

seq – equivalent stress; material specific (J2, Hill, Tresca, Barlat) if mat is provided, J2 equivalent stress otherwise

Return type:

float

seq_j2()[source]

Calculate J2 principal stress

Returns:

sj2 – equivalent stress

Return type:

float

theta()[source]

Calculate polar angle in deviatoric plane

Returns:

ang – polar angle of stress in deviatoric plane

Return type:

float

pylabfea.basic.a_vec = array([ 0.81649658, -0.40824829, -0.40824829])

First unit vector spanning deviatoric stress plane (real axis)

pylabfea.basic.b_vec = array([ 0.        ,  0.70710678, -0.70710678])

Second unit vector spanning deviatoric stress plane (imaginary axis)

pylabfea.basic.eps_eq(eps: ndarray)[source]

Calculate equivalent strain

Parameters:

eps ((3,), (6,), (nsc,3) or (nsc,6) array) – (3,) or (nsc,3): Principal strains; (6,) or (nsc,6): Voigt strains

Returns:

eeq – equivalent strains

Return type:

float or (nsc,) array

pylabfea.basic.pickle2mat(name, path='./')[source]

Read pickled material file.

Parameters:
  • name (string) – File name of pickled material to be read.

  • path (string) – Path under which pickle-files are stored (optional, default: ‘./’)

Returns:

pcl – unpickled material

Return type:

Material object

pylabfea.basic.polar_ang(sig)[source]
pylabfea.basic.s_cyl(sig, mat=None)[source]
pylabfea.basic.sdev(sig)[source]
pylabfea.basic.seq_J2(sig)[source]
pylabfea.basic.sig_cyl2princ(s_cyl) ndarray[source]

Convert cylindrical stress into 3D Cartesian principle stress

Parameters:

s_cyl ((2,), (3,), (N,2) or (N,3) array) – Cylindrical stress in form (seq, theta, (optional: p))

Returns:

s_princ – principle deviatoric stresses

Return type:

(3,) or (N,3) array

pylabfea.basic.sig_cyl2voigt(sig_cyl: ndarray, eigen_vector: ndarray) ndarray[source]

Convert cylindrical stress and eigenvectors into Voigt stress tensor

Parameters:
  • sig_cyl ((3,) or (N,3) array) – Cylindrical stress in form (seq, theta, p)

  • eigen_vector ((3,3) or (N,3,3) array) – Eigenvectors of stress tensor

Returns:

sig_voigt – Voigt stress tensor

Return type:

(6,) or (N,6) array

pylabfea.basic.sig_dev(sig: ndarray) ndarray[source]

Calculate deviatoric stress component from given stress tensor

Parameters:

sig ((3,), (6,) (N,3) or (N,6) array)

Returns:

sd – deviatoric stresses

Return type:

float or (N,) array

pylabfea.basic.sig_eq_j2(sig: ndarray)[source]

Calculate sj2 equivalent stress from any stress tensor

Parameters:

sig ((3,), (6,) (N,3) or (N,6) array) – (3,), (N,3): Principal stress or list of principal stresses; (6,), (N,6): Voigt stress

Returns:

seq – sj2 equivalent stresses

Return type:

float or (N,) array

pylabfea.basic.sig_polar_ang(sig: ndarray)[source]

Transform stresses into polar angle on deviatoric plane spanned by a_vec and b_vec

Parameters:

sig ((3,), (6,) (N,3) or (N,6) array) – (3,), (N,3): Principal stresses; (6,), (N,6): Voigt stress

Returns:

theta – polar angles in deviatoric plane as positive angle between sig and a_vec in range [-pi,+p]

Return type:

float or (N,) array

pylabfea.basic.sig_princ(sig: ndarray)[source]

Convert Voigt stress tensors into principal stresses and eigenvectors.

Parameters:

sig ((6,), (nsc,6), (3,3), or (nsc,3,3) array) – Voigt stress tensor (dim=6) or Cartesian stress tensor (dim=3x3)

Returns:

  • spa ((3,) or (nsc,3) array) – Principal stresses

  • eva ((3,3) or (nsc,3,3) array) – Eigenvectors/rotation matrices of stress tensor

pylabfea.basic.sig_princ2cyl(sig: ndarray, mat=None) ndarray[source]

Convert principal stress into cylindrical stress vector

Parameters:
  • sig ((3,), (6,), (N,3) or (N,6) array) – stress to be converted, if (3,) or (N,3) principal stress is assumed

  • mat (object of class Material) – Material for Hill-type principal stress (optional)

Returns:

sc – stress in cylindrical coordinates (seq, theta, p)

Return type:

(3,) or (N,3) array

pylabfea.basic.sp_cart(scyl)[source]
pylabfea.basic.sprinc(sig)[source]
pylabfea.basic.svoigt(scyl, evec)[source]
pylabfea.basic.yf_tolerance = 0.005

Plastic yielding if yield function > yf_tolerance

Type:

Tolerance

Module Model

Module pylabefea.model introduces global functions for mechanical quantities and class Model that contains that attributes and methods needed in FEA. Materials are defined in module pylabfea.material

uses NumPy, SciPy, MatPlotLib

Version: 4.1 (2022-02-22) Author: Alexander Hartmaier, ICAMS/Ruhr University Bochum, Germany Email: alexander.hartmaier@rub.de distributed under GNU General Public License (GPLv3)

class pylabfea.model.Model(dim=1, planestress=False)[source]

Bases: object

Class for finite element model. Defines necessary attributes and methods for pre-processing (defining geometry, material assignments, mesh and boundary conditions); solution (invokes numerical solver for non-linear problems); and post-processing (visualization of results on meshed geometry, homogenization of results into global quantities).

Geometry and sections can be defined with the methods geom and assign. There are pre-defined options that make the generation of laminate structures particularly easy. But also more general section can be defined and each section can be associated with a different material.

Boundary conditions on left-hand-side and bottom nodes are always assumed to be static (typically with values of zero); Boundary conditions on right-hand-side and top nodes are considered as load steps starting from zero. Default boundary conditions are:

  • lhs: fixed in x-direction (ux=0), free in y-direction (fy=0)

  • bot: fixed in y-direction (uy=0), free in x-direction (fx=0)

  • rhs: free (fx=fy=0)

  • top: free (fx=fy=0)

There are pre-defined sets of nodes to top, bottom, left and right boundaries to which either force or displacement controlled loads can be applied in x or y-direction with the methods bctop, bcbot, bcleft, bcright. Furthermore, it is possible to define boundary conditions for a freely defined set of nodes with the method bcnode. The latter option if useful, to fix only one corner node when laterally free boundaries shall be implemented. It can also be used to simulate loads on parts of boundaries, as it occurs for indentations, etc.

Visualization is performed with the method plot, which can display various mechanical quantities on the deformed mesh. Homogenization of boundary loads is simply performed with the method calc_glob by which all global stresses and strains are obtained by averaging over the element values and by summing up the boundary loads for comparison. This is particularly useful when calculating the stress-strain behavior of a structure.

Parameters:
  • dim (int) – Dimensionality of model (optional, default: 1)

  • planestress (Boolean) – Sets plane-stress condition (optional, default: False)

dim

Dimensionality of model (1 and 2 are supported in this version)

Type:

integer

planestress

Sets plane-stress condition

Type:

Boolean

Nsec

Number of sections (defined in geom)

Type:

int

LS

Absolute lengths of sections (defined in geom)

Type:

1-d array

lenx

Size of model in x-direction (defined in geom)

Type:

float

leny

Size of model in y-direction (defined in geom, default leny=1)

Type:

float

thick

Thickness of 2d-model (defined in geom, default thick=1)

Type:

float

nonlin

Indicates non-linearity of model (defined in assign)

Type:

Boolean

mat

List of materials assigned to sections of model, same dimensions as LS (defined in assign)

Type:

list of objects to class Material

ubcbot

True: displacement BC on rhs nodes; False: force BC on rhs nodes (defined in bcright)

Type:

dim-array of Boolean

bcb

Nodal displacements in x or y-direction for lhs nodes (defined in bcbot)

Type:

dim-array

ubcleft

True: displacement BC on rhs nodes; False: force BC on rhs nodes (defined in bcright)

Type:

dim-array of Boolean

bcl

Nodal displacement in x or y-direction for lhs nodes (defined in bcleft)

Type:

dim-array

ubcright

True: displacement BC on rhs nodes; False: force BC on rhs nodes (defined in bcright)

Type:

dim-array of Boolean

bcr

Nodal displacements/forces in x or y-direction for rhs nodes (defined in bcright)

Type:

dim-array

ubctop

True: displacement BC on top nodes; False: force BC on top nodes (defined in bctop)

Type:

Boolean

bct

Nodal displacements/forces in x or y-direction for top nodes (defined in bctop)

Type:

dim-array

Nnode

Total number of nodes in Model (defined in mesh)

Type:

int

NnodeX, NnodeY

Numbers of nodes in x and y-direction (defined in mesh)

Type:

int

Nel

Number of elements (defined in mesh)

Type:

int

Ndof

Number of degrees of freedom (defined in mesh)

Type:

int

npos
Nodal positions, dimension Ndof, same structure as u,f-arrays:

(x1, y1, x2, y1, …) (defined in mesh)

Type:

1d-array

noleft, noright, nobot, notop

Lists of nodes on boundaries (defined in mesh)

Type:

list

noinner

List of inner nodes (defined in mesh)

Type:

list

element

List of objects of class Element, dimension Nel (defined in mesh)

Type:

list

u

List of nodal displacements (defined in solve)

Type:

(Ndof,) array

f

List of nodal forces (defined in solve)

Type:

(Ndof,) array

sgl

Time evolution of global stress tensor with incremental load steps (defined in solve)

Type:

(N,6) array

egl

Time evolution of global total strain tensor with incremental load steps (defined in solve)

Type:

(N,6) array

epgl

Time evolution of global plastic strain tensor with incremental load steps (defined in solve)

Type:

(N,6) array

glob

Global values homogenized from BC or element solutions, contains the elements: ‘ebc1’, ‘ebc2’, ‘sbc1’, ‘sbc2’ : global strain and stress from BC (type: float) ‘eps’, ‘epl’, ‘sig’, : global strain, plastic strain, and stress tensors homogenized from element solutions (type: Voigt tensor) (defined in calc_global)

Type:

python dictionary

class Element(model, nodes, lx, ly, mat)[source]

Bases: object

Class for isoparametric elements; supports 1-d elements with linear or quadratic shape function and full integration; 2-d quadrilateral elements with linear shape function and full integration

Parameters:
  • model (object of class Model) – Refers to parent class Model to inherit all attributes

  • nodes (list) – List of nodes belonging to element

  • lx (float) – Element size in x-direction

  • ly (float) – Element size in y-direction

  • sect (int) – Section in which element is located

  • mat (object of class Material) – Material card with attributes and methods to be applied in element

Model

Refers to parent object (FE model)

Type:

object of class Model

nodes

List of nodes of this element

Type:

list

Lelx

Size of element in x-direction

Type:

float

Lely

Size of element in y-direction

Type:

float

ngp

number of Gauss points

Type:

int

gpx

x-locations of Gauss points in element

Type:

1-d array

gpy

y-locations of Gauss points in element

Type:

1-d array

Bmat

List of B-matrices at Gauss points

Type:

list

wght

List of weight of each Gauss point in integration

Type:

1-d array

Vel

Volume of element

Type:

float

Jac

Determinant of Jacobian of iso-parametric formulation

Type:

float

CV

Voigt stiffness matrix of element material

Type:

(6,6) array

elstiff

Tangent stiffness tensor of element

Type:

2-d array

Kel

Element stiffness matrix

Type:

2-d array

Mat

Material model to be applied

Type:

object of class Material

eps

average (total) element strain (defined in Model.solve)

Type:

Voigt tensor

sig

average element stress (defined in Model.solve)

Type:

Voigt tensor

epl

average plastic element strain (defined in Model.solve)

Type:

Voigt tensor

calc_Bmat(x=0.0, y=0.0)[source]

Calculate B matrix at position x in element

Parameters:
  • x (float) – absolute x-position in element (optional, default: 0)

  • y (float) – absolute y-position in element (optional, default: 0)

Returns:

B – matrix (N=dim^2*(SF+1))

Return type:

6xN array

calc_Kel()[source]

Calculate element stiffness matrix by Gaussian integration

depl()[source]

Calculate plastic strain increment

Returns:

depl – Plastic strain increment

Return type:

Voigt tensor

deps()[source]

Calculate strain increment in element

Returns:

deps – Strain increment

Return type:

Voigt tensor

dsig()[source]

Calculate stress increment

Returns:

dsig – Stress increment

Return type:

Voigt tensor

eps_t()[source]

Calculate total strain in element

Returns:

eps_t – Total strain

Return type:

Voigt tensor

node_num()[source]

Calculate indices of DOF associated with element

Returns:

ind – List of indices

Return type:

list of int

assign(mats)[source]

Assigns an object of class Material to each section.

Parameters:

mats (list) – List of materials, dimension must be equal to number of sections

material.mat

Internal variable for materials assigned to each section of the geometry

Type:

List of material objects

material.nonlin

Indicate if material non-linearity must be considered

Type:

bool

bcbot(val=0.0, bctype='disp', bcdir='y')[source]

Define boundary conditions on bottom nodes, either force or displacement type; static boundary conditions are assumed for bottom boundary. The default is freezing all y-displacements to zero.

Parameters:
  • val (float) – Displacement in bcdir direction (optional, default: 0)

  • bctype (str) – Type of boundary condition (‘disp’ or ‘force’) (optional, default: ‘disp’)

  • bcdir (str or int) – Direction of boundary load (optional, default: ‘y’)

bcleft(val=0.0, bctype='disp', bcdir='x')[source]

Define boundary conditions on lhs nodes, either force or displacement type; static boundary conditions are assumed for lhs boundary. The default is freezing all x-displacements to zero.

Parameters:
  • val (float) – Displacement or force of lhs nodes in bc_dir direction (optional, default: 0)

  • bctype (str) – Type of boundary condition (‘disp’ or ‘force’) (optional, default: ‘disp’)

  • bcdir (str or int) – Direction of boundary load (optional, default: ‘x’)

bcnode(node, val, bctype, bcdir)[source]

Define boundary conditions on a set of nodes defined in node, either force or displacement type in x or y-direction are accepted. If non.zero, the boundary loads will be incremented step-wise until the given boundary conditions are fulfilled.

Since nodes must be given, this subroutine can only be called after meshing.

Parameters:
  • node (int or list of int) – Node or set of nodes to which BC shall be applied

  • val (float) – Displacement or force in bcdir direction

  • bctype (str) – Type of boundary condition (‘disp’ or ‘force’)

  • bcdir (str or int) – Direction of boundary load (‘x’ or ‘y’; 0 or 1)

bcright(val, bctype, bcdir='x')[source]

Define boundary conditions on rhs nodes, either force or displacement type. If non-zero, the boundary loads will be incremented step-wise until the given boundary conditions are fulfilled.

Parameters:
  • val (float) – Displacement or force on rhs nodes in bc_dir direction

  • bctype (str) – Type of boundary condition (‘disp’ or ‘force’)

  • bcdir (str or int) – Direction of boundary load (optional, default: ‘x’)

bctop(val, bctype, bcdir='y')[source]

Define boundary conditions on top nodes, either force or displacement type. If non-zero, the boundary loads will be incremented step-wise until the given boundary conditions are fulfilled.

Parameters:
  • val (float) – Displacement or force in bcdir direction

  • bctype (str) – Type of boundary condition (‘disp’ or ‘force’)

  • bcdir (str or int) – Direction of boundary load (optional, default: ‘y’)

bcval(nodes)[source]

Calculate average displacement and total force at (boundary) nodes

Parameters:

nodes (list) – List of nodes

calc_global()[source]

Calculate global quantities and store in Model.glob; homogenization done by averaging residual forces (sbc1/2) and displacements (ebc1/2) at boundary nodes or by averaging element quantities (sig, eps, epl)

Yields:

Model.glob (dictionary) – Values for stress (‘sig’), total strain (‘eps’) and plastic strain (‘epl’) as homogenized element solutions (Voigt tensors); values for (11)- and (22)-components of stress (‘sbc1’,’scb2’) and total strain (‘ebc1’, ‘ebc2’) as homogenized values of boundary nodes.

geom(sect=1, LX=None, LY=1.0, LZ=1.0)[source]

Specify geometry of FE model with dimensions LX, LY and LZ and its subdivision into a number of sections; for 2-d model a laminate structure normal to x-direction can be created easily with the in-built option to pass a list with the absolute lengths of each section in the parameter sect. When this parameter is an integer it merely reserves space for sections with arbitrary geometries; adds attributes to class Model

Parameters:
  • sect (list or int) – Either number of sections (int) or list with with absolute length of each section

  • LX (float) – Length of model in x-direction (optional, default None in which case LX is calculated as the sum of the lengths of all sections)

  • LY (float) – Length in y direction (optional, default: 1)

  • LZ (float) – Thickness of model in z-direction (optional, default: 1)

mesh(elmts=None, nodes=None, NX=10, NY=1, SF=1)[source]

Import mesh or generate structured mesh with quadrilateral elements (2d models). First, nodal positions Model.npos are defined such that nodes lie at corners (linear shape function) and edges (quadratic shape function) of elements. Then, elements are initialized as object of class Model.Element, which requires the list of nodes associated with the element, the dimensions of the element, and the material of the section in which the element is situated to be passed.

Parameters:
  • elmts ((NX, NY) array) – Represents number of material as defined by list Model.mat

  • nodes ((2,) array) – Defines positions of nodes on regular grid.

  • NX (int) – Number of elements in x-direction (optional, default: 10)

  • NY (int) – Number of elements in y-direction (optional, default: 1)

  • SF (int) – Degree of shape functions: 1=linear, 2=quadratic (optional, default: 1)

plot(fsel, mag=10, colormap='viridis', cdepth=20, showmesh=True, shownodes=True, vmin=None, vmax=None, annot=True, file=None, show_fig=True, pos_bar=0.83, fig=None, ax=None, show_bar=True)[source]

Produce graphical output: draw elements in deformed shape with color according to field variable ‘fsel’; uses matplotlib

Parameters:
  • fsel (str) – Field selector for library field, see Keyword Arguments for possible values

  • mag (float) – Magnification factor for displacements (optional, default: 10)

  • vmin (float) – Start value for range of plotted values

  • vmax (float) – End value for range of plotted values

  • cdepth (int) – Number of colors in colormap (optional, default: 20)

  • showmesh (Boolean) – Set/unset plotting of lines for element edges (optional, default: True)

  • shownodes (Boolean) – Set/unset plotting of nodes (optional, default: True)

  • colormap (str) – Name of colormap to be used (optional, default: viridis)

  • annot (Boolean) – Show annotations for x and y-axis (optional, default: True)

  • file (str) – If a filename is provided, plot is exported as PDF (optional, default: None)

  • show_fig (bool) – True: show figure, False: return figure handle (optional, default: True)

  • pos_bar (float) – Position of color bar, (optional; for use in Jupyter notebook: 1.01, for python: 0.83)

Keyword Arguments:
  • strain1 – total strain in horizontal direction

  • strain2 – total strain in vertical direction

  • strain12 – total shear strain, xy-component

  • stress1 – horizontal stress component

  • stress2 – vertical stress component

  • stress12 – xy-shear stress component

  • plastic1 – plastic strain in horizontal direction

  • plastic2 – plastic strain in vertical direction

  • plastic12 – plastic strain, xy-component

  • seq – equivalent stress (Hill-formulation for anisotropic plasticity)

  • seqJ2 – equivalent J2 stress

  • peeq – equivalent plastic strain

  • etot – equivalent total strain

  • ux – horizontal displacement

  • uy – vertical displacement

  • mat – materials and sections of model

setupK()[source]

Calculate and assemble system stiffness matrix based on element stiffness matrices.

Returns:

K – System stiffness matrix

Return type:

2d-array

solve(min_step=None, verb=False)[source]

Solve linear system of equations K.u = f with respect to u, to obtain distortions of the system under the applied boundary conditions for mechanical equilibrium, i.e., when the total force on internal nodes is zero. In the first step, the stiffness matrix K, the nodal displacements u and the nodal forces f are modified to conform with the boundary conditions (calc_BC), then an elastic predictor step is calculated, which is already the final solution for linear problems. For non-linear problems, i.e. for plastic materials, the load step must be controlled and subdivided to fulfill the side conditions of plastic materials, i.e. the equivalent stress must remain smaller than the yield strength, or the flow stress during plastic yielding. This is ensured in a self.consistency loop for non-linear models. The system of equations is solved by invoking the subroutine numpy.linalg.solve. The method yields the final solution for nodal displacements u and nodal forces f as attributes of class Model; global stress, global total strain and global plastic strain are evaluated and stored as attributes, too. Element solutions for stresses and strains are stored as attributes of class Element, see documentation of this class.

Parameters:
  • min_step (int) – Minimum number of load steps (optional)

  • verb (Boolean) – Be verbose in text output (optional, default: False)

Yields:
  • Model.u ((Model.Ndof,) array) – Nodal displacements

  • Model.f ((Model.Ndof,) array) – Nodal forces

  • Model.sgl ((N,6) array) – Global stress as as Voigt tensor for each incremental load step (homogenized element solution)

  • Model.egl ((N,6) array) – Global total strain as as Voigt tensor for each incremental load step (homogenized element solution)

  • Model.epgl ((N,6) array) – Global plastic strain as as Voigt tensor for each incremental load step (homogenized element solution)

  • Element.sig ((6,) array) – Element solution for stress tensor

  • Element.eps ((6,) array) – Element solution for total strain tensor

  • Element.epl ((6,) array) – Element solution for plastic strain tensor

Module Material

Module pylabfea.material introduces class Material that contains attributes and methods needed for elastic-plastic material definitions in FEA. It also enables the training of machine learning algorithms as yield functions for plasticity. The module pylabfea.model is used to calculate mechanical properties of a defined material under various loading conditions.

uses NumPy, ScipPy, MatPlotLib, sklearn, pickle, and pyLabFEA.model

Version: 4.1 (2022-01-23) Authors: Alexander Hartmaier, Ronak Shoghi, ICAMS/Ruhr University Bochum, Germany Email: alexander.hartmaier@rub.de

distributed under GNU General Public License (GPLv3)

class pylabfea.material.Material(name='Material', num=1)[source]

Bases: object

Define class for Materials including material parameters (attributes), constitutive relations (methods) and derived properties und various loading conditions (dictionary).

Parameters:
  • name (str) – Name of material (optional, default: ‘Material’)

  • num (int) – Material number (optional, default: 1)

name

Name of material

Type:

str

num

Material number

Type:

int

sy

Yield strength

Type:

float

ML_yf

Existence of trained machine learning (ML) yield function (default: False)

Type:

Boolean

ML_grad

Existence of trained ML gradient (default: False)

Type:

Boolean

dev_only

Consider only deviatoric stress tensor (default: False)

Type:

Boolean

tresca

Indicate if Tresca equivalent stress should be used (default: False)

Type:

Boolean

hill_3p

Indicates whether 3-paramater Hill model should be used (default: False)

Type:

Boolean

hill_6p

Indicates whether 6-paramater Hill model should be used (default: False)

Type:

Boolean

barlat

Indicate if Barlat equivalent stress should be used (default: False)

Type:

Boolean

msg

Messages returned

Type:

dictionary

prop

Derived properties under defined load paths

Type:

dictionary

propJ2

Derived properties in form of J2 equivalent stress

Type:

dictionary

sigeps

Data of stress-strain curves under defined load paths

Type:

dictionary

C11, C12, C44

Anisotropic elastic constants

Type:

float

E, nu

Isotropic elastic constants, Young modulus and Poisson number

Type:

float

msparam

Dicitionary with microstructural parameters assigned to this material

Type:

ditionary

whdat

Indicates existence of work hardening data

Type:

Boolean

txdat

Indicates existance of data for different textures

Type:

Boolean

Ndof

degress of freedom for yield function, mandatory: 1:seq, 2:theta; optional: 3:work_hard, 4:texture)

Type:

int

Keyword Arguments:
  • prop-propJ2 – Store properties of material (Hill-formulation or J2 formulation) in sub-dictonaries: ‘stx’ (tensile horiz. stress), ‘sty’ (tensile vert. stress), ‘et2’ (equibiaxial tensile strain), ‘ect’ (pure shear)

  • stx-sty-et2-ect (sub-dictionaries) – Store data for ‘ys’ (float - yield strength), seq (array - eqiv. stress), ‘eeq’ (array - equiv. total strain), ‘peeq’ (array - equiv. plastic strain), ‘sytel’ (str - line style for plot), ‘name’ (str - name in legend)

  • sigeps – Store tensorial stress strain data in sub-directories; Contains data for ‘sig’ (2d-array - stress), ‘eps’ (2d-array - strain), ‘epl’ (2d-array - plastic strain)

  • msparam – Store data on microstructure parameters: ‘Npl’, ‘Nlc, ‘Ntext’, ‘texture’, ‘peeq_max’, ‘work_hard’, ‘flow_stress’ are obtained from data analysis module. Other parameters can be added.

  • msg – Messages that can be retrieved: ‘yield_fct’, ‘gradient’, ‘nsteps’, ‘equiv’

C_tan(sig, Cel, epl=None)[source]

Calculate tangent stiffness relaxing stress back to yield locus; Reference: M.A. Crisfield, Non-linear finite element analysis of solids and structures, Chapter 6, Eqs. (6.9) and (6.18)

Parameters:
  • sig (Voigt tensor) – Stress

  • Cel ((6,6) array) – Elastic stiffness tensor used for predictor step

  • epl ((sdim,) array) – Equivalent plastic strain tensor (optional, default: 0.)

Returns:

Ct – Tangent stiffness tensor

Return type:

(6,6) array

ML_full_yf(sig, epl=None, ld=None, verb=True)[source]

Calculate full ML yield function as distance of a single given stress tensor to the yield locus in loading direction.

Parameters:
  • sig ((sdim,) array) – Voigt stress tensor

  • epl ((sdim,) array) – Equivalent plastic strain (optional, default: 0)

  • ld ((6,) array) – Vector of loading direction in princ. stress space (optional)

  • verb (Boolean) – Indicate whether to be verbose in text output (optional, default: False)

Returns:

yf – Full ML yield function, i.e. distance of sig to yield locus in ld-direction

Return type:

float

calc_fgrad(sig, epl=None, seq=None, accumulated_strain=0.0, max_stress=0.0, ana=False)[source]

Calculate gradient to yield surface. Three different methods can be used: (i) analytical gradient to Hill-like yield function (default if no ML yield function exists - ML_yf=False), (ii) gradient to ML yield function (default if ML yield function exists - ML_yf=True; can be overwritten if ana=True), (iii) ML gradient fitted seperately from ML yield function (activated if ML_grad=True and ana=False)

Parameters:
  • sig ((sdim,) or (N,sdim) array) – Stress value (Pricipal stress or full stress tensor)

  • epl ((sdim,) array) – Plastic strain tensor (optional, default = None)

  • accumulated_strain

  • max_stress

  • seq (float or (N,) array) – Equivalent stresses (optional)

  • ana (Boolean) – Indicator if analytical solution should be used, rather than ML yield fct (optional, default: False)

Returns:

fgrad – Gradient to yield surface at given position in stress space, same dimension as sdim

Return type:

(sdim,), (N,sdim) array

calc_properties(size=2, Nel=2, verb=False, eps=0.005, min_step=None, sigeps=False, load_cases=['stx', 'sty', 'et2', 'ect'])[source]

Use pylabfea.model to calculate material strength and stress-strain data along a given load path.

Parameters:
  • size (int) – Size of FE model (optional, defaul: 2)

  • Nel (int) – Number of elements per axis (optional, defaul: 2)

  • verb (Boolean) – Be verbose with output (optional, default: False)

  • eps (float) – Maximum total strain (optional, default:0.005)

  • min_step (int) – Minumum number of load steps (optional)

  • sigeps (Boolean) – Decide if data for stress and strain tensors is stored in dictionary Material.sigeps (optional, default: False)

  • load_cases (list) –

    List of load cases to be performed (optional, default:

    [‘stx’,’sty’,’et2’,’ect’]);

    ’stx’: uniaxial tensile yield stress in horizontal (x-)direction; ‘sty’: uniaxial tensile yield stress in vertical (y-)direction; ‘et2’: plane stress, equibiaxial strain in x and y direction; ‘ect’: pure shear strain (x-compression, y-tension), plane stress

calc_seq(sig)[source]

Calculate generalized equivalent stress from stress tensor; equivalent J2 stress for isotropic flow behavior and tension compression invariance; Hill-type approach for anisotropic plastic yielding; Drucker-like approach for tension-compression asymmetry; Barlat 2004-18p model for plastic anisotropy; Tresca equivalent stress

Step 1: transform input into
  1. sig: (N,6)-array of Voigt stresses (zeros added if input is princ. stresses)

  2. sp: (N,3)-array of principal stresses

N=1 if input is single stress in which case return value is of type float

Step 2: Call appropriate subroutines for evaluation of equiv. stress. Currently supported are:
  1. Tresca

  2. Barlat Yld2004-18p

  3. Hill 3-parameter (hill_3p) or 6-parameter (hill_6p)

  4. von Mises/J2 (special case of Hill with all coefficients equal 1,

    material independent, works also for elastic and ML materials)

Parameters:

sig ((sdim,) or (N,sdim) array) – Stress values (for dim=3 principal stresses are assumed, otherwise Voigt stress)

Returns:

seq – Hill-Drucker-type equivalent stress

Return type:

float or (N,) array

calc_seqB(sv)[source]

Calculate equivalent stress based on Yld2004-18p yield function proposed by Barlat et al, Int. J. Plast. 21 (2005) 1009

Parameters:

sv ((6,) array) – Voigt stress tensor

Returns:

seq – Equivalent stress

Return type:

float

calc_yf(sig, epl=None, accumulated_strain=0.0, max_stress=0.0, ana=False, pred=False)[source]

Calculate yield function

Parameters:
  • sig ((sdim,) or (N, sdim) array) – Stresses (arrays of Voigt or principal stresses)

  • epl ((sdim, ) array) – Equivalent plastic strain tensor (optional, default: 0)

  • max_stress

  • accumulated_strain

  • ana (Boolean) – Indicator if analytical solution should be used, rather than ML yield fct (optional, default: False)

  • pred (Boolean) – Indicator if prediction value should be returned, rather than decision function (optional, default: False)

Returns:

f – Yield function for given stress (same length as sig)

Return type:

flot or 1d-array

create_sig_data(N=None, mat_ref=None, sdata=None, Nseq=2, sflow=None, offs=0.01, extend=False, rand=False, Fe=0.1, Ce=0.99)[source]

Function to create consistent data sets on the deviatoric stress plane for training or testing of ML yield function. Either the number “N” of raw data points, i.e. load angles, to be generated and a reference material “mat_ref” has to be provided, or a list of raw data points “sdata” with Cartesian stress tensors lying on the yield surface serves as input. Based on the raw data, stress tensors from the yield locus are distributed into the entire deviatoric space, by linear downscaling into the elastic region and upscaling into the plastic region. Data is created in form of cylindrical stresses that lie densly around the expected yield locus and more sparsely in the outer plastic region.

Parameters:
  • N (int) – Number of load cases (polar angles) to be created (optional, either N and mat_ref or sdata must be provided)

  • mat_ref (object of class Material) – reference material needed to calculate yield function if only N is provided (optional, ignored if sdata is given)

  • sdata ((N, sdim) array) – List of Cartesian stress tensors lying on yield locus. Based on these yield stresses, training data in entire deviatoric stress space is created (optional, either sdata or N and mat_ref must be provided)

  • Nseq (int) – Number of training stresses to be generated in the range ‘offs’ to the yield strength (optional, default: 12)

  • sflow (float) – Expected flow stress of data set (optional, default: self.sy)

  • offs (float) – Start of range for equiv. stress (optional, default: 0.01)

  • extend (Boolean) – Create additional data in plastic regime (optional, default: False)

  • rand (Boolean) – Chose random load cases (polar angles) (optional, default: False)

  • Fe (float) – Relative value of lowest stress in elastic regime

  • Ce (float) – Relative value of largest stress in elastic regime

Returns:

  • st ((M, sdim) array) – Cartesian training stresses, M = N (2 Nseq + Nextend)

  • yt ((M,) array) – Result vector of categorial yield function (-1 or +1) for supervised training

elasticity(C11=None, C12=None, C44=None, CV=None, E=None, nu=None)[source]

Define elastic material properties

Parameters:
  • C11 (float)

  • C12 (float)

  • C44 (float) – Anisoptropic elastic constants of material (optional, if (C11,C12,C44) not given, either (E,nu) or CV must be specified)

  • E (float)

  • nu (float) – Isotropic parameters Young’s modulus and Poisson’s number (optional, if (E,nu) not given, either (C11,C12,C44) or CV must be specified)

  • CV ((6,6) array) – Voigt matrix of elastic constants (optional, if CV not given, either (C11,C12,C44) or (E,nu) must be specified)

Return type:

None.

ellipsis(a=1.0, b=np.float64(0.5773502691896258), n=72)[source]

Create ellipsis with main axis along 45° axis, used for graphical representation of isotropic yield locus.

Parameters:
  • a (float) – Long half-axis of ellipsis (optional, default: 1)

  • a – Short half-axis of ellipsis (optional, default: 1/sqrt(3))

  • n (int) – Number of points on ellipsis to be calculated

Returns:

x, y – x and y coordinates of points on ellipsis

Return type:

(n,) arrayy

epl_dot(sig, epl, Cel, deps)[source]

Calculate plastic strain increment relaxing stress back to yield locus; Reference: M.A. Crisfield, Non-linear finite element analysis of solids and structures, Chapter 6, Eqs. (6.4), (6.8) and (6.17)

Parameters:
  • sig ((6,)-array) – Voigt stress tensor

  • epl ((6,)-array) – Voigt plastic strain tensor

  • Cel ((6,6) array) – Elastic stiffnes tensor

  • deps (Voigt tensor) – Strain increment from predictor step

Returns:

pdot – Plastic strain increment

Return type:

Voigt tensor

export_MLparam(sname, source=None, file=None, path='../../models/', descr=None, param=None)[source]

The parameters of the trained Ml flow rule (support vectors, dual coefficients, offset and scaling parameters) are written to a csv file that is readable to Abaqus (8 numbers per line).

Parameters:
  • sname (str) – Name of script that created this material

  • source (str) – Source of parameters (optional, default: None)

  • file (str) – Trunk of filename to which CSV flies are written (optional, default: None)

  • path (str) – Path to which files are written (optional: default: ‘../../models/’)

  • descr (list) – List of names of model parameters used for generating this ML material (optional, default: [])

  • param (list) – List of values of parameters used for generating this ML material (optional, default: []); descr and param must be of the same size

Yields:
  • CSV file with name path+file+’-svm.csv’ containing

  • support vectors, dual coefficients, and (Ndof, elastic parameters, offset,

  • gamma value, scaling factors) in Abaqus format; and

  • JSON file with name path+file+’-svm_meta.json’ containing meta data in given format.

find_yloc(x, su, epl=None)[source]

Function to expand unit stresses by factor and calculate yield function; used by search algorithm to find zeros of yield function.

Parameters:
  • x ((N,)-array) – Multiplyer for stress

  • su ((N,sdim) array) – Unit stress

  • epl ((sdim,) array)) – Plastic strain tensor (optional, default: None)

Returns:

f – Yield function evaluated at sig=x.sp

Return type:

(N,)-array

find_yloc_scalar(x, su, epl=None)[source]

Function to expand unit stresses by factor and calculate yield function; used by search algorithm to find zeros of yield function.

Parameters:
  • x (float) – Multiplier for stress

  • su ((sdim,) array) – Unit stress

  • epl ((sdim,) array)) – Plastic strain tensor (optional, default: None)

Returns:

f – Yield function evaluated at sig=x.sp

Return type:

float

from_MLparam(name, path='../../models/')[source]

Define material properties from parameters of trained machine learning models that have been written with Material.export_MLparam. Will invoke definition of elastic parameters by calls to the methods Material.elasticity with the parameters provided in the data set. Also initializes current texture to first one in list and resets work hardening parameters.

Parameters:
  • name (string) – Name of parameter files (name.csv file and metadata file name_meta.json)

  • path (string) – Path in which files are stored (optional, default: ‘../../models/’)

from_data(param)[source]

Define material properties from data sets generated in module Data: contains data on elastic and plastic behavior, including work hardening, for different crystallographic textures. Possible to extend to grain sizes, grain shapes and porosities. Will invoke definition of elastic and plastic parameters by calls to the methods Material.elasticity and Material.plasticity with the parameters provided in the data set. Also initializes current texture to first one in list and resets work hardening parameters.

Parameters:

param (list of directories) – Data.mat_param directories containing material data sets

get_sflow(epl)[source]

Calculate an estimate of the scalar flow stress (strength) of the material for a given plastic strain.

NOTE: Currently assumes only linear isotropic hardening with the current hardening rate, does not include texture information and needs to be adapted to data contained in ms.param

Parameters:

epl (float or (sdim,) array) – Current value of equiv. plastic strain (float) or plastic strain tensor

Yields:

sflow (float) – Average flow stress

pckl(name=None, path='../../materials/')[source]

Write material into pickle file. Usefull for materials with trained machine learning flow rules to avoid time-consuming re-training.

Parameters:
  • name (string (optional, default: None)) – File name for pickled material. The default is None, in which case the filename will be the material name + ‘.pckl’.

  • path (string) – Path to location for pickles

Return type:

None.

plasticity(sy=None, sdim=6, drucker=0.0, khard=0.0, tresca=False, barlat=None, barlat_exp=None, hill=None, hill_3p=None, hill_6p=None, rv=None)[source]

Define plastic material parameters; anisotropic Hill-like and Drucker-like behavior is supported

Parameters:
  • sy (float) – Yield strength

  • hill ((3,) or (6,) array) – Parameters for Hill-like orthotropic anisotropy (optional, default: isotropy)

  • drucker (float) – Parameter for Drucker-like tension-compression asymmetry (optional, default: 0)

  • khard (float) – Linear strain hardening slope (ds/de_p) (optional, default: 0)

  • tresca (Boolean) – Indicate if Tresca equivalent stress should be used (optional, default: False)

  • barlat ((18,) array) – Array with parameters for Barlat Yld2004-18p yield function (optional)

  • barlat_exp (int) – Exponent for Barlat Yld2004-18p yield function (optional)

  • hill_3p (Boolean) – Indicate if 3-parameter Hill model shall be applied (optional, default: None) Will be set True automatically if 3 Hill paramaters != 1 are provided

  • hill_6p (Boolean) – Indicate if 6-parameter Hill model shall be applied (optional, default: None) Will be set True automatically if 6 Hill parameters are provided

  • rv ((6,) array) – Parameters for anisotropic flow aspect ratios that can be given alternatively to Hill parameters (optional)

  • sdim (int) – Dimensionality of stress tensor to be used for plasticity, must be either 3 (only principal stresses are considered) or 6 (full stress tensor is considered), (optional, default: 6)

plot_data(Z, axs, xx, yy, field=True, c='red')[source]

Plotting data in stress space to visualize yield loci.

Parameters:
  • Z (array) – Data for field plot

  • axs (handle) – Axes where plot is added

  • xx (meshgrid) – x-coordinates

  • yy (meshgrid) – y-coordinates

  • field (Boolean) – Decide if field is plotted (optional, default: True)

  • c (str) – Color for contour line (optional, default: ‘red’)

Returns:

line – Reference to plotted line

Return type:

handle

plot_stress_strain(Hill=False, file=None, fontsize=14)[source]

Plot stress-strain data and print values for strength. Requires ‘calc_properties’ to be executed beforehand and plots data stored in attribute ‘props’.

Parameters:
  • Hill (Boolean) – Decide if data for Hill-type equivalent stress is presented (optional, default: False)

  • file (str) – Filename to save plot (optional)

  • fontsize (int) – Fontsize for axis annotations (optional, default: 14)

plot_yield_locus(fun=None, label=None, data=None, trange=0.01, peeq=0.0, xstart=None, xend=None, axis1=[0], axis2=[1], iso=False, ref_mat=None, field=False, Nmesh=100, file=None, fontsize=20, scaling=True)[source]

Plot different cuts through yield locus in 3D principal stress space.

Parameters:
  • fun (function handle) – Yield function to be plotted (optional, default: own yield function)

  • label (str) – Label for yield function (optional, default: own name)

  • data ((N,3) array) – principal stress data to be used for scatter plot (optional)

  • trange (float) – Cut-off for data to be plotted on slice (optional, default: 1.e-2)

  • peeq (float) – Level of plastic strain for which yield locus is plotted (isotropic hardening)

  • xstart (float) – Smallest value on x-axis (optional, default: -2)

  • xend (float) – Largest value on x-axis (optional, default: 2)

  • axis1 (list) – Cartesian stress coordinates to be plotted on x-axis of slices (optional, default: [0])

  • axis2 (list) – Cartesian stress coordinates to be plotted on y-axis of slices (optional, default: [1])

  • iso (Boolean) – Decide if reference ellipsis for isotropic material is plotted (optional, default: False)

  • ref_mat=None – Reference material to plot yield locus (optional)

  • field (Boolean) – Decide if field of yield function is plotted (optional, default: False)

  • Nmesh (int) – Number of mesh points per axis on which yield function is evaluated (optional, default:100)

  • file (str) – File name for output of olot (optional)

  • fontsize (int) – Fontsize for axis annotations (optional, default: 20)

  • scaling (Boolean) – Scale stress with yield strength (optional, default: True)

Returns:

axs – Axis of the plot

Return type:

pyplot axis handle

polar_plot_yl(Na=72, cmat=None, data=None, dname='reference', scaling=None, field=False, predict=False, cbar=False, Np=100, file=None, arrow=False, sJ2=False, show=True)[source]

Plot yield locus as polar plot in deviatoric stress plane

Parameters:
  • Na (int) – Number of angles on which yield locus is evaluated (optional, default: 72)

  • cmat (list of materials) – Materials of which YL is plotted in same plot with same scaling (optional)

  • data ((N,3) array) – Array of cylindrical stress added to plot (optional)

  • dname (str) – Label for data (optional, default: reference)

  • scaling (float) – Scaling factor for stresses (optional)

  • field (Boolean) – Field of decision function is plotted, works only together with ML yield function (optional, default: False)

  • predict (Boolean) – Plot ML prediction (-1,1), otherwise decision fucntion is plotted (optional, default: False)

  • Np (int) – Number of points per axis for field plot (optional, default: 100)

  • cbar (Boolean) – Plot colorbar for field (optional, default: False)

  • file (str) – Name of PDF file to which plot is saved

  • arrow (Boolean) – Indicate if arrows for the pricipal stress directions are shown (optional, default: False)

  • sJ2 (Boolean) – Indicate that J2 equivalent stress shall be used instead of material definition of equivalent stress (optional, default: False)

response(sig, epl, deps, CV, maxit=50)[source]

Calculate non-linear material response to deformation defined by load step, corresponds to user material function.

Parameters:
  • sig ((6,) array) – Voigt stress tensor at start of load step (=end of previous load step)

  • epl ((6,) array) – Voigt plastic strain tensor at start of load step

  • deps ((6,) array) – Voigt strain tensor defining deformation (=load step)

  • CV ((6,6) array) – Voigt elastic tensor

  • maxit (int) – Maximum number of iteration steps (optional, default= 5)

Returns:

  • fy1 (real) – Yield function at end of load step (indicates whether convergence is reached)

  • sig ((6,) array) – Voigt stress tensor at end of load step

  • depl ((6,) array) – Voigt tensor of plastic strain increment at end of load step

  • grad_stiff ((6,6) array) – Tangent material stiffness matrix (d_sig/d_eps) at end of load step

set_texture(current, verb=False)[source]

Set parameters for current crystallographic texture of material as defined in microstructure.

Parameters:
  • current (float or list) – Mixture parameter for microstructures in range [0,1] for each defined microstructure, indicates the intensity of the given microstructure. The sum of all mixture parameters must be <=1, the rest will be set to random texture. Must have same dimension as material.msparam.

  • verb (Boolean) – Be verbose

Yields:
  • Material.tx_cur (list) – Current value of microstructural mixture parameter for active microstructure. Has same dimension as material.msparam

  • Material.sy (float) – Yield stength is redefined accoring to texture parameter

  • Material.khard (float) – Work hardening parameter is redefined according to texture parameter

  • Material.epc (float) – Critical PEEQ for which onset of plastic deformation is definied in data

setup_fgrad_SVM(X_grad_train, y_grad_train, C=10.0, gamma=0.1)[source]

Inititalize and train SVM regression for gradient evaluation

Parameters:
  • X_grad_train ((N,3) array)

  • y_grad_train ((N,) array)

  • C (float) – Paramater for training of Support Vector Regression (SVR) (optional, default:10)

  • gamma (float) – Parameter for kernel of SVR (optional, default: 0.1)

setup_yf_SVM(x, y_train, x_test=None, y_test=None, C=15.0, gamma=2.5, fs=0.1, plot=False, cyl=False, gridsearch=False, cvals=None, gvals=None, vlevel=3)[source]

Generic function call to setup and train the SVM yield function, for details see the specific functions setup_yf_SVM_6D and setup_yf_SVM_3D.

setup_yf_SVM_3D(x, y_train, x_test=None, y_test=None, C=10.0, gamma=1.0, fs=0.1, plot=False, cyl=False, gridsearch=False, cvals=None, gvals=None)[source]

Initialize and train Support Vector Classifier (SVC) as machine learning (ML) yield function. Training and test data (features) are accepted as either 3D principal stresses or cylindrical stresses, but principal stresses will be converted to cylindrical stresses, such that training is always performed in cylindrical stress space, with equiv. stress at yield onset and polar angle as degrees of freedom. Graphical output on the trained SVC yield function is possible.

WARNING: 3D data types are no longer actively supported

Parameters:
  • x ((N,2) or (N,3) array) – Training data either as Cartesian princ. stresses (N,3) or cylindrical stresses (N,2)

  • cyl (Boolean) – Indicator for cylindrical stresses if x is has shape (N,3)

  • y_train (1d-array) – Result vector for training data (same size as x)

  • x_test ((N,2) or (N,3) array) – Test data either as Cartesian princ. stresses (N,3) or cylindrical stresses (N,2) (optional)

  • y_test – Result vector for test data (optional)

  • C (float) – Parameter for training of SVC (optional, default: 10)

  • gamma (float) – Parameter for kernel function of SVC (optional, default: 1)

  • fs (float) – Parameters for size of periodic continuation of training data (optional, default:0.1)

  • plot (Boolean) – Indicates if plot of decision function should be generated (optional, default: False)

  • gridsearch (Boolean) – Perform grid search to optimize hyperparameters of ML flow rule (optional, default: False)

  • cvals (array) – Values for SVC training parameter C in gridsearch (optional, default: None)

  • gvals (array) – Values for SVC parameter gamma in gridsearch (optional, default: None)

Returns:

  • train_sc (float) – Training score

  • test_sc (float) – test score

setup_yf_SVM_6D(x, y_train, x_test=None, y_test=None, C=10.0, gamma=1.0, plot=False, gridsearch=False, cvals=None, gvals=None, vlevel=3)[source]

Initialize and train Support Vector Classifier (SVC) as machine learning (ML) yield function. Training and test data (features) are accepted as either 3D principal stresses or cylindrical stresses, but principal stresses will be converted to cylindrical stresses, such that training is always performed in cylindrical stress space, with equiv. stress at yield onset and polar angle as degrees of freedom. Graphical output on the trained SVC yield function is possible.

Parameters:
  • x ((N,self.Ndof) array) – Training data in form of deviatoric Voigt stresses, components s1-s6), s0=-s1-s2. Additional DOF for work hardening and texture if considered.

  • y_train (1d-array) – Result vector for training data (same size as x)

  • x_test ((N,self.Ndof) array) – Test data either as Cartesian princ. stresses (N,3) or cylindrical stresses (N,2) (optional)

  • y_test – Result vector for test data (optional)

  • C (float) – Parameter for training of SVC (optional, default: 10)

  • gamma (float) – Parameter for kernel function of SVC (optional, default: 1)

  • plot (Boolean) – Indicates if plot of decision function should be generated (optional, default: False)

  • gridsearch (Boolean) – Perform grid search to optimize hyper parameters of ML flow rule (optional, default: False)

  • cvals (array) – Values for SVC training parameter C in gridsearch (optional, default: None)

  • gvals (array) – Values for SVC parameter gamma in gridsearch (optional, default: None)

  • vlevel (int) – Value for verbosity of grid search algorithm (optional, defult: 3)

Returns:

  • train_sc (float) – Training score

  • test_sc (float) – test score

test_data_generation(C=10, gamma=4, Nlc=36, Nseq=25, fs=0.3, extend=False, mat_ref=None, sdata=None, fontsize=16, gridsearch=False, cvals=None, gvals=None, Fe=0.1, Ce=0.99)[source]

A function to generate test data to get the scores, which is exactly as we are generating the training data but use those to test and get the score.

Parameters:
  • C (float) – Parameter needed for training process, larger values lead to more flexibility (optional, default: 10)

  • gamma (float) – Parameter of Radial Basis Function of SVC kernel, larger values, lead to faster decay of influence of individual support vectors, i.e., to more short ranged kernels (optional, default: 4)

  • Nlc (int) – Number of load cases to be considered, will be overwritten if material has microstructure information (optional, default: 36)

  • Nseq (int) – Number of training and test stresses to be generated in elastic regime, same number will be produced in plastic regime (optional, default: 25)

  • fs (float) – Parameter to ensure peridicity of yield function wrt. theta

  • extend (Boolean) – Indicate whether training data should be extended further into plastic regime (optional, default: True)

  • mat_ref (object of class Material) – reference material needed to calculate yield function if only N is provided (optional, ignored if sdata is given)

  • sdata ((N, sdim) array) –

    List of Cartsian stresses lying on yield locus. Based on these yield stresses, training data in entire deviatoric stress space is created (optional, f no data in self.msparam is given, either

    sdata or N and mat_ref must be provided)

  • plot (Boolean) – Indicate if graphical output should be performed (optional, default: False)

  • fontsize (int) – Fontsize for graph annotations (optional, default: 16)

  • gridsearch (Boolean) – Perform grid search to optimize hyperparameters of ML flow rule (optional, default: False)

  • Fe (float) – Relative value of lowest stress in elastic regime (optional, default: 0.1)

  • Ce (float) – Relative value of largest stress in elastic regime (optional, default: 0.99)

train_SVC(C=10, gamma=4, Nlc=36, Nseq=25, fs=0.3, extend=False, mat_ref=None, sdata=None, plot=False, fontsize=16, gridsearch=False, cvals=None, gvals=None, vlevel=3, Fe=0.1, Ce=0.99)[source]

Train SVC for all yield functions of the microstructures provided in msparam and for flow stresses to capture work hardening. In first step, the training data for each set is generated by creating stresses on the deviatoric plane and calculating their catgegorial yield function (“-1”: elastic, “+1”: plastic). Furthermore, axes in different dimensions for microstructural features are introduced that describe the relation between the different sets.

Parameters:
  • C (float) – Parameter needed for training process, larger values lead to more flexibility (optional, default: 10)

  • gamma (float) – Parameter of Radial Basis Function of SVC kernel, larger values, lead to faster decay of influence of individual support vectors, i.e., to more short ranged kernels (optional, default: 4)

  • Nlc (int) – Number of load cases to be considered, will be overwritten if material has microstructure information (optional, default: 36)

  • Nseq (int) – Number of training and test stresses to be generated in elastic regime, same number will be produced in plastic regime (optional, default: 25)

  • fs (float) – Parameter to ensure periodicity of yield function wrt. theta

  • extend (Boolean) – Indicate whether training data should be extended further into plastic regime (optional, default: True)

  • mat_ref (object of class Material) – reference material needed to calculate yield function if only N is provided (optional, ignored if sdata is given)

  • sdata ((N, sdim) array) –

    List of Cartsian stresses lying on yield locus. Based on these yield stresses, training data in entire deviatoric stress space is created (optional, f no data in self.msparam is given, either

    sdata or N and mat_ref must be provided)

  • plot (Boolean) – Indicate if graphical output should be performed (optional, default: False)

  • fontsize (int) – Fontsize for graph annotations (optional, default: 16)

  • gridsearch (Boolean) – Perform grid search to optimize hyperparameters of ML flow rule (optional, default: False)

  • Fe (float) – Relative value of lowest stress in elastic regime (optional, default: 0.1)

  • Ce (float) – Relative value of largest stress in elastic regime (optional, default: 0.99)

  • cvals (array) – Values for SVC training parameter C in gridsearch (optional, default: None)

  • gvals (array) – Values for SVC parameter gamma in gridsearch (optional, default: None)

  • vlevel (int) – Value for verbosity of grid search algorithm (optional, defult: 3)

Module Training

Module pylabfea.training introduces methods to create training data for ML flow rule in shape of unit stresses that are evenly distributed in the stress space to define the load cases for which the critical stress tensor at which plastic yielding starts needs to be determined.

uses NumPy, ScipPy, MatPlotLib, sklearn, and pyLabFEA.basic

Version: 4.0 (2021-11-27) Authors: Ronak Shoghi, Alexander Hartmaier, ICAMS/Ruhr University Bochum, Germany Email: alexander.hartmaier@rub.de distributed under GNU General Public License (GPLv3)

Subroutines int_sin_m, primes and uniform_hypersphere have been adapted from code published by Stack Overflow under the CC-BY-SA 4.0 license, see https://stackoverflow.com/questions/57123194/how-to-distribute-points-evenly-on-the-surface-of-hyperspheres-in-higher-dimensi/59279721#59279721 These subroutines are distributed here under the CC-BY-SA 4.0 license, see https://creativecommons.org/licenses/by-sa/4.0/

pylabfea.training.create_test_sig(Json, Number_sig_per_strain=4)[source]

A function to generate test data for micromechanical simulations based on a given material’s stress-strain data. :param Json: Path to the JSON file containing the stress-strain data of the material. :type Json: str :param Number_sig_per_strain: The number of test cases to generate for each strain level in the elastic or in the plastic range.

The total number of load cases per strain level will be 2 * Number_sig_per_strain. Default is 4.

Returns:

  • ts_sig (np.ndarray) – A numpy array of stress values for the generated test cases.

  • epl_tot (np.ndarray) – A numpy array of strain values for the generated test cases.

  • yf_ref (np.ndarray) – A numpy array with the same length as ts_sig and epl_tot, filled with +1 for the first half of the elements and -1 for the second half. This represents the known sign of the yield function: +1 for upscaling and -1 for downscaling. This array will be used as reference in the training_score function.

pylabfea.training.int_sin_m(x, m)[source]

Computes the integral of sin^m(t) dt from 0 to x recursively

Parameters:
  • x (float) – Upper limit of integration

  • m (int) – Power of trigonometric function to be considered

Returns:

f – Value of integral

Return type:

float

pylabfea.training.load_cases(number_3d, number_6d, method='brentq')[source]

Generate unit stresses in principal stress space (3d) and full stress space (6d)

Parameters:
  • number_3d (int) – Number of principal unit stresses to be created

  • number_6d (int) – Number of full unit stresses to be created

Returns:

allsig – Unit stresses

Return type:

(number_3d+number6d, 6)-array

pylabfea.training.primes()[source]

Infinite generator of prime numbers

pylabfea.training.training_score(yf_ref, yf_ml, plot=False)[source]

Calculate the accuracy of the training result in form of different measures as compared to given reference values.

Parameters:
  • yf_ref ((N,)-array) – Yield function values of reference material

  • yf_ml ((N,)-array) – Yield function values of ML material at identical sequence of stresses at which reference material is evaluated.

Returns:

  • mae (float) – Mean Average Error

  • precision (float) – Ratio of true positives w.r.t. all positives

  • Accuracy (float) – Ratio of true positives and true negative w.r.t. all results

  • Recall (float) – Ratio of true positives w.r.t. true positives and false negatives

  • F1Score (float) – F1 score

  • MCC (float) – Matthews Correlation Coefficient

pylabfea.training.uniform_hypersphere(d, n, method='brentq')[source]

Generate n usnits stresse on the d dimensional hypersphere representing create load cases in 3D or 6D stress space

Parameters:
  • d (int) – Dimension of stress space in which to create unit stresses

  • n (int) – Number of stresses to be created

Returns:

points – Unit stresses

Return type:

(n,6)-array

Module Data

Module pylabfea.data introduces the class Data for handling of data resulting from virtual or physical mechanical tests in the pyLabFEA package. This class provides the methods necessary for analyzing data. During this processing, all information is gathered that is required to define a material, i.e., all parameters for elasticity, plasticity, and microstructures are provided from the data. Materials are defined in module pylabfea.material based on the analyzed data of this module.

uses NumPy, SciPy, MatPlotLib

Version: 4.0 (2021-11-27) Last Update: (24-04-2023) Authors: Ronak Shoghi, Alexander Hartmaier, ICAMS/Ruhr University Bochum, Germany Email: alexander.hartmaier@rub.de distributed under GNU General Public License (GPLv3)

class pylabfea.data.Data(source, path_data='./', name='Dataset', mat_name='Simulanium', sdim=6, epl_crit=None, epl_start=None, epl_max=None, depl=0.0, plot=False, wh_data=True, texture_name='Random')[source]

Bases: object

Define class for handling data from virtual mechanical tests in micromechanical simulations and data from physical mechanical tests on materials with various microstructures. Data is used to train machine learning flow rules in pyLabFEA.

Parameters:
  • source (str or data array) – Either filename of datafile or array of initial yield stresses

  • path_data (str) – Trunc of pathname for data files (optional, default: ‘./’)

  • name (str) – Name of Dataset (optional, default: ‘Dataset’)

  • mat_name (str) – Name of material (optional, default: ‘Simulanium)

  • sdim (int) – Dimensionality of stresses; if sdim = 3 only principal stresses are considered (optional, default: 6)

  • epl_crit (float) – Critical plastic strain at which yield strength is defined(optional, default: 2.e-3)

  • epl_start (float) – Start value of equiv. plastic strain at which data will be sampled(optional, default: 1.e-3)

  • epl_max (float) – Maximum equiv. plastic strain up to which data is considered (optional, default: 0.05)

  • plot (Boolean) – Plot data (optional, default: False)

  • wh_data (Boolean) – Data for flow stresses and plastic strain tensors in work hardening regime exists (optional, default: True)

name

Name of dataset

Type:

str

sy_av

Average initial yield stress

Type:

float

E_av

Average Young’s modulus

Type:

float

nu_av

Average Poisson ratio

Type:

float

mat_param
Contains available data for ML yield function and for microstructural parameters
Items
epc : critical value of equiv. plastic strain that defines the onset of plastic yielding
ep_start : start value of equiv. plastic strain at which data is acquired
ep_max : maximum value of equiv. plastic strain up to which data is considered
delta_ep : minimum separation of plastic strains used for training, if 0 no separation is enforced
lc_indices : array with start index for each load case
peeq_max : maximum of equiv. plastic strain that occurred in data set, corrected for onset of plasticity at epc
sdim : dimension of stress vector (must be 3 or 6)
Name : material name
Dataset : Name of dataset
wh_data : indicate if strain hardening data exists
Ntext : number of textures
tx_name : name of texture
texture : texture parameters
flow_stress : array of flow stresses correlated to plastic strains
plastic_strain : array of plastic strains corrected for onset of plasticity at epc
E_av : avereage Young’s modulus derived from data
nu_av : average Poisson’s ratio derived from data
sy_av : average yield strength, i.e. flow stress at epc, obtained from data
Nlc : number of load cases in data
sy_list : ???
sig_ideal : interpolated stress tensors at onset of plastic yielding (epc) for each load case
Type:

dictionary

add2mat_data(data_dict, key)[source]
convert_data(sig)[source]

Convert data provided only for stress tensor at yield point into mat_param dictionary

Parameters:

sig (ndarray) – Stress tensors at yield onset

key_parser(key)[source]
parse_data(epl_crit, epl_start, epl_max, depl)[source]

Read data and store in attribute ‘mat_data’ Estimate elastic properties and initial yield strength from data for each load case and form averages.

Parameters:
  • epl_crit (float) – Critical value for onset of yielding

  • epl_start (float) – Start value of equiv. plastic strain at which data acquisition of flow stresses will start

  • epl_max (float) – Maximum equiv. strain up to which data is considered

  • depl (float) – Minimum separation of plastic strains used for training, if 0, no minimum is enforced

plot_data(data, xlabel, ylabel, emax=None)[source]
plot_set(db, mat_param)[source]
plot_training_data(emax=1)[source]
plot_yield_locus(db, mat_data, active, scatter=False, data=None, data_label=None, arrow=False, file=None, title=None, fontsize=18)[source]
read_data(Data_File)[source]
pylabfea.data.find_transition_index(stress)[source]

Calculates the index at which a significant transition in the total stress-strain relationship occurs. The function applies a Savitzky-Golay filter to smooth the stress data and to calculate the first and second derivatives. It identifies the transition index by finding the point where the second derivative starts to change from zero, i.e. the value in the elastic regime. This approach relies on the overall behavior of the stress-strain relationship rather than focusing on a specific plastic strain threshold (e.g., 0.002%).

Parameters:

stress: 1-d numpy array

Array of equivalent stress values along one load path.

Returns:

idx: int

The index within the stress array where a significant transition from linear behavior occurs.

pylabfea.data.get_elastic_coefficients(eps, sig, method='least_square', initial_guess=None)[source]

A function to compute the elastic coefficients (stiffness matrix) for a material based on stress-strain data. This function supports two methods for determining the stiffness matrix: a direct least squares approach and an optimization approach. The least squares method is used when ‘least_square’ is specified, which processes any desired stress-strain pairs( minimum must be 4). The optimization approach, used when ‘decomposition’ method is specified, iteratively adjusts the stiffness matrix to minimize an objective function that measures the fit to the observed data, subject to physical plausibility constraints.

Parameters:
  • eps ((N, 6)-array) – Array with Voigt strain tensors at the end of the linear-elastic regime of different load cases.

  • sig ((N, 6)-array)

  • cases. (Array of corresponding Voigt stress tensors at the end of the linear elastic regime of different load)

  • method (str) – Method to be used for calculating the stiffness matrix. Options include ‘least_square’ and ‘decomposition’. The ‘least_square’ method calculates the stiffness matrix using a least squares fit to the stress-strain data. The ‘decomposition’ method uses an optimization approach with the specified method for interpreting the stiffness matrix from the optimization variables. (Optional: default is ‘least_square’).

  • initial_guess (np.ndarray or None) – Initial guess for the stiffness matrix coefficients used in the optimization approach. If None, a random initial guess is generated. This parameter is ignored when using the ‘least_square’ method. Default is None.

  • Returns

  • C ((6, 6)-array) – The optimized symmetric 6x6 stiffness matrix (elastic coefficients) as a numpy array. If the optimization or calculation fails to converge to a solution within the maximum number of attempts, the function may return the last attempted solution or raise an error, depending on implementation details.