Layer results

Once the laminate load-deformation relations are solved, the strains in the laminate coordinate system at any position $z$ through the thickness of the laminate can be computed by:

\begin{equation} \begin{bmatrix} \varepsilon_x \\ \varepsilon_y \\ \gamma_{xy} \end{bmatrix}= \begin{bmatrix} \varepsilon_x^0 \\ \varepsilon_y^0 \\ \gamma_{xy}^0 \end{bmatrix}+ z\begin{bmatrix} \kappa_x \\ \kappa_y \\ \kappa_{xy} \end{bmatrix} \end{equation}

Example:

In [1]:
import numpy as np
from laminatelib import laminateStiffnessMatrix

m1={'E1':140000, 'E2':10000, 'v12':0.3, 'G12':5000}

layup  =[ {'mat':m1 , 'ori':  0  , 'thi':1} , 
          {'mat':m1 , 'ori': 90  , 'thi':1} ,
          {'mat':m1 , 'ori':  0  , 'thi':1}  ]

ABD = laminateStiffnessMatrix(layup)

deforms = np.linalg.solve(ABD, [100,0,0,200,0,0])

print(deforms[0:3])  # mid-plane strains
print(deforms[3:6])  # curvatures
[ 3.43209976e-04 -1.93055611e-05  5.93885692e-21]
[ 6.56256349e-04 -1.32891911e-04  6.82900397e-21]

A structure for the layer results

The derived results for one layer includes stresses and strains in two coordinate systems, as well as stress exposure factors and other relevant data such as the z-coordinates. All of these results should generally be computed at two locations; at the bottom and at the top of the layer. A suggestion for a useful datastructure is illustrated:

                                    Layer k
                                       |
           ------------------------------------------------------------------
           |                   |                      |                     |
        stress              strain                   fail                   h
     -----------         -----------         -------------------          -----
     |         |         |         |         |        |        |          |   |    
    xyz       123       xyz       123        MS       ME       TW        bot top
   -----     -----     -----     -----      -----    -----    -----
   |   |     |   |     |   |     |   |      |   |    |   |    |   |
  bot top   bot top   bot top   bot top    bot top  bot top  bot top

The data structure can easily be created using dictionaries and then be organized into a list containing the results from all layers:

In [2]:
from laminatelib import laminateThickness, T2Ds, T2De, Q2D, fE2DMS, fE2DME, fE2DTW

def layerResults(layup,deformations):
    # An empty list that shall be populated with results from all layers:
    results=[]
    # The bottom coordinate of the laminate:
    bot = -laminateThickness(layup)/2
    for layer in layup:
        # the top of current layer is just the bottom + the thickness:
        top=bot+layer['thi']
        # creating a dictionary with the two coordinates:
        h={'bot':bot, 'top':top}
        # computing the strains in laminate coordinate system:
        strnXYZbot=deformations[0:3]+bot*deformations[3:6]
        strnXYZtop=deformations[0:3]+top*deformations[3:6]
        # put the strains into a dictionary:
        strnXYZ={'bot':strnXYZbot, 'top':strnXYZtop}
        # the transformation matrix for strains:
        Te=T2De(layer['ori'])
        # transforming the strains into the material coordinate system:
        strn123bot=np.dot(Te, strnXYZbot)
        strn123top=np.dot(Te, strnXYZtop)
        # the strains at top and bottom as a dictionary:
        strn123={'bot':strn123bot, 'top':strn123top}
        # all strains as a new dictionary:
        strn={'xyz':strnXYZ, '123':strn123}
        # stiffness matrix of the material:
        Q=Q2D(layer['mat'])
        # Hooke's law: finding the stresses in material system:
        strs123bot=np.dot(Q,strn123bot)
        strs123top=np.dot(Q,strn123top)
        # the material stresses as a dictionary having bottom and topp results:
        strs123={'bot':strs123bot, 'top':strs123top}
        # transformation matrix for stress, negativ rotation applies now,
        # since this is from material system to laminat system (reversed..)
        Ts=T2Ds(-layer['ori'])
        # transforming stresses into laminate system:
        strsXYZbot=np.dot(Ts,strs123bot)
        strsXYZtop=np.dot(Ts,strs123top)
        # organizing the stresses for the two location in a dictionary
        strsXYZ={'bot':strsXYZbot, 'top':strsXYZtop}
        # all stresses as a new nested dictionary:
        strs={'xyz':strsXYZ, '123':strs123}
        # Failure criteria...follows the same logic as above
        MSbot=fE2DMS(strs123bot,layer['mat'])
        MStop=fE2DMS(strs123top,layer['mat'])
        MS={'bot':MSbot, 'top':MStop}
        MEbot=fE2DME(strs123bot,layer['mat'])
        MEtop=fE2DME(strs123top,layer['mat'])
        ME={'bot':MEbot, 'top':MEtop}
        TWbot=fE2DTW(strs123bot,layer['mat'])
        TWtop=fE2DTW(strs123top,layer['mat'])
        TW={'bot':TWbot, 'top':TWtop}
        fail={'MS':MS, 'ME':ME, 'TW':TW}
        # and finally put everything into a top level dictionary,
        # and add that to the list
        results.append( {'h':h , 'strain':strn, 'stress':strs, 'fail':fail }  )
        bot=top
    return results

Example:

In [3]:
from laminatelib import solveLaminateLoadCase

m1={'E1':140000, 'E2':10000, 'v12':0.3, 'G12':5000, 'XT':1200, 'XC':800, 'YT':50, 'YC':120, 'S12':75, 'f12':-0.5}

layup1  =[ {'mat':m1 , 'ori':   0  , 'thi':1} , 
           {'mat':m1 , 'ori': -45  , 'thi':1} ,
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori': +45  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1}  ]

ABD1=laminateStiffnessMatrix(layup1)

load,defs = solveLaminateLoadCase(ABD=ABD1, Nx=500)

res = layerResults(layup1,defs)

# Accessing, examples:

print(  res[0]['h']['bot']        )           # z-coordinate at the bottom of the first layer
print(  res[0]['stress']['xyz']['bot']   )    # stress in the global c-sys at the bottom of the first layer
print(  res[1]['stress']['123']['top']   )    # stress in the material cooridnate system at the top of the second layer
-3.0
[203.66497762   2.71737626   8.01833861]
[55.56324148 10.28461825  8.06807109]

Visualization

The function plotLayerStresses(results) found in Plot Gallery makes graphs of layer results through the thickness of a laminate where results are obtained by the function layerResults(layup,deformations).

In [4]:
from plotlib import plotLayerStresses

Examples:

In [5]:
m1={'E1':140000, 'E2':10000, 'v12':0.3, 'G12':5000, 'XT':1200, 'XC':800, 'YT':50, 'YC':120, 'S12':75, 'f12':-0.5}

layup1  =[ {'mat':m1 , 'ori':   0  , 'thi':1} , 
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1}   ]

ABD1=laminateStiffnessMatrix(layup1)
load,defs = solveLaminateLoadCase(ABD=ABD1, Nx=500, Nxy=200)
res = layerResults(layup1,defs)    
plotLayerStresses(res)
In [6]:
m1={'E1':140000, 'E2':10000, 'v12':0.3, 'G12':5000, 'XT':1200, 'XC':800, 'YT':50, 'YC':120, 'S12':75, 'f12':-0.5}

layup1  =[ {'mat':m1 , 'ori':   0  , 'thi':1} , 
           {'mat':m1 , 'ori': -45  , 'thi':1} ,
           {'mat':m1 , 'ori': +45  , 'thi':1} ,
           {'mat':m1 , 'ori': +45  , 'thi':1} ,
           {'mat':m1 , 'ori': -45  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1}  ]

ABD1=laminateStiffnessMatrix(layup1)
load,defs = solveLaminateLoadCase(ABD=ABD1, Nx=500)
res = layerResults(layup1,defs)    
plotLayerStresses(res)
In [7]:
m1={'E1':140000, 'E2':10000, 'v12':0.3, 'G12':5000, 'XT':1200, 'XC':800, 'YT':50, 'YC':120, 'S12':75, 'f12':-0.5}

layup1  =[ {'mat':m1 , 'ori':   0  , 'thi':1} , 
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1} ,
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1}  ]

ABD1=laminateStiffnessMatrix(layup1)
load,defs = solveLaminateLoadCase(ABD=ABD1, Ny=500, Mx=200)
res = layerResults(layup1,defs)    
plotLayerStresses(res)

The function plotLayerFailure(results) found in Plot Gallery makes graphs of layer results through the thickness of a laminate where results are obtained by the function layerResults(layup,deformations).

In [8]:
from plotlib import plotLayerFailure

Examples:

In [9]:
m1={'E1':140000, 'E2':10000, 'v12':0.3, 'G12':5000, 'XT':1200, 'XC':800, 'YT':50, 'YC':120, 'S12':75, 'f12':-0.5}

layup1  =[ {'mat':m1 , 'ori':   0  , 'thi':1} , 
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1} ,
           {'mat':m1 , 'ori':  90  , 'thi':1} ,
           {'mat':m1 , 'ori':   0  , 'thi':1}  ]

ABD1=laminateStiffnessMatrix(layup1)
load,defs = solveLaminateLoadCase(ABD=ABD1, Ny=500, Mx=200)
res = layerResults(layup1,defs)    
plotLayerFailure(res)
TOC Next Prev

Disclaimer:This site is designed for educational purposes only. There are most likely errors, mistakes, typos, and poorly crafted statements that are not detected yet... www.ntnu.edu/employees/nils.p.vedvik

Copyright 2024, All rights reserved